• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

ElasticModulus_Def.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #include <fstream>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 #include "Sacado_ParameterRegistration.hpp"
00011 #include "Albany_Utils.hpp"
00012 
00013 namespace LCM {
00014 
00015 template<typename EvalT, typename Traits>
00016 ElasticModulus<EvalT, Traits>::
00017 ElasticModulus(Teuchos::ParameterList& p) :
00018   elasticModulus(p.get<std::string>("QP Variable Name"),
00019      p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout"))
00020 {
00021   Teuchos::ParameterList* elmd_list = 
00022     p.get<Teuchos::ParameterList*>("Parameter List");
00023 
00024   Teuchos::RCP<PHX::DataLayout> vector_dl =
00025     p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00026   std::vector<PHX::DataLayout::size_type> dims;
00027   vector_dl->dimensions(dims);
00028   numQPs  = dims[1];
00029   numDims = dims[2];
00030 
00031   Teuchos::RCP<ParamLib> paramLib = 
00032     p.get< Teuchos::RCP<ParamLib> >("Parameter Library", Teuchos::null);
00033 
00034   std::string type = elmd_list->get("Elastic Modulus Type", "Constant");
00035   if (type == "Constant") {
00036     is_constant = true;
00037     constant_value = elmd_list->get("Value", 1.0);
00038 
00039     // Add Elastic Modulus as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Elastic Modulus", this, paramLib);
00042   }
00043 //  else if (type == 'Variable') {
00044 //    is_constant = true; // this means no stochastic nature
00045 //    is_field = true;
00046 //    constant_value = elmd_list->get("Value", 1.0);
00047 //
00048 //    // Add Elastic Modulus as a Sacado-ized parameter
00049 //    new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00050 //      "Elastic Modulus", this, paramLib);
00051 //
00052 //  }
00053   else if (type == "Truncated KL Expansion") {
00054     is_constant = false;
00055     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim>
00056       fx(p.get<std::string>("QP Coordinate Vector Name"), vector_dl);
00057     coordVec = fx;
00058     this->addDependentField(coordVec);
00059 
00060     exp_rf_kl = 
00061       Teuchos::rcp(new Stokhos::KL::ExponentialRandomField<MeshScalarT>(*elmd_list));
00062     int num_KL = exp_rf_kl->stochasticDimension();
00063 
00064     // Add KL random variables as Sacado-ized parameters
00065     rv.resize(num_KL);
00066     for (int i=0; i<num_KL; i++) {
00067       std::string ss = Albany::strint("Elastic Modulus KL Random Variable",i);
00068       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00069       rv[i] = elmd_list->get(ss, 0.0);
00070     }
00071   }
00072   else {
00073     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00074              "Invalid elastic modulus type " << type);
00075   } 
00076 
00077   // Optional dependence on Temperature (E = E_ + dEdT * T)
00078   // Switched ON by sending Temperature field in p
00079 
00080   if ( p.isType<std::string>("QP Temperature Name") ) {
00081     Teuchos::RCP<PHX::DataLayout> scalar_dl =
00082       p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00083     PHX::MDField<ScalarT,Cell,QuadPoint>
00084       tmp(p.get<std::string>("QP Temperature Name"), scalar_dl);
00085     Temperature = tmp;
00086     this->addDependentField(Temperature);
00087     isThermoElastic = true;
00088     dEdT_value = elmd_list->get("dEdT Value", 0.0);
00089     refTemp = p.get<RealType>("Reference Temperature", 0.0);
00090     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00091                                 "dEdT Value", this, paramLib);
00092   }
00093   else {
00094     isThermoElastic=false;
00095     dEdT_value=0.0;
00096   }
00097 
00098   // Optional dependence on porosity (E = E_ + dEdT * T)
00099     // Switched ON by sending Temperature field in p
00100 
00101   if ( p.isType<std::string>("Porosity Name") ) {
00102       Teuchos::RCP<PHX::DataLayout> scalar_dl =
00103         p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00104       PHX::MDField<ScalarT,Cell,QuadPoint>
00105         tporo(p.get<std::string>("Porosity Name"), scalar_dl);
00106       porosity = tporo;
00107       this->addDependentField(porosity);
00108       isPoroElastic = true;
00109 
00110     }
00111     else {
00112       isPoroElastic= false;
00113     }
00114 
00115 
00116   this->addEvaluatedField(elasticModulus);
00117   this->setName("Elastic Modulus"+PHX::TypeString<EvalT>::value);
00118 }
00119 
00120 // **********************************************************************
00121 template<typename EvalT, typename Traits>
00122 void ElasticModulus<EvalT, Traits>::
00123 postRegistrationSetup(typename Traits::SetupData d,
00124                       PHX::FieldManager<Traits>& fm)
00125 {
00126   this->utils.setFieldData(elasticModulus,fm);
00127   if (!is_constant)    this->utils.setFieldData(coordVec,fm);
00128   if (isThermoElastic) this->utils.setFieldData(Temperature,fm);
00129   if (isPoroElastic)   this->utils.setFieldData(porosity,fm);
00130 }
00131 
00132 // **********************************************************************
00133 template<typename EvalT, typename Traits>
00134 void ElasticModulus<EvalT, Traits>::
00135 evaluateFields(typename Traits::EvalData workset)
00136 {
00137   std::size_t numCells = workset.numCells;
00138 
00139   if (is_constant) {
00140     for (std::size_t cell=0; cell < numCells; ++cell) {
00141       for (std::size_t qp=0; qp < numQPs; ++qp) {
00142   elasticModulus(cell,qp) = constant_value;
00143       }
00144     }
00145   }
00146   else {
00147     for (std::size_t cell=0; cell < numCells; ++cell) {
00148       for (std::size_t qp=0; qp < numQPs; ++qp) {
00149   Teuchos::Array<MeshScalarT> point(numDims);
00150   for (std::size_t i=0; i<numDims; i++)
00151     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00152   elasticModulus(cell,qp) = exp_rf_kl->evaluate(point, rv);
00153       }
00154     }
00155   }
00156   if (isThermoElastic) {
00157     for (std::size_t cell=0; cell < numCells; ++cell) {
00158       for (std::size_t qp=0; qp < numQPs; ++qp) {
00159   elasticModulus(cell,qp) += dEdT_value * (Temperature(cell,qp) - refTemp);
00160       }
00161     }
00162   }
00163   if (isPoroElastic) {
00164       for (std::size_t cell=0; cell < numCells; ++cell) {
00165         for (std::size_t qp=0; qp < numQPs; ++qp) {
00166     // porosity dependent Young's Modulus. It will be replaced by
00167     // the hyperelasticity model in (Borja, Tamagnini and Amorosi, ASCE JGGE 1997).
00168     elasticModulus(cell,qp) = constant_value;
00169  //       *sqrt(2.0 - porosity(cell,qp));
00170         }
00171       }
00172     }
00173 }
00174 
00175 // **********************************************************************
00176 template<typename EvalT,typename Traits>
00177 typename ElasticModulus<EvalT,Traits>::ScalarT& 
00178 ElasticModulus<EvalT,Traits>::getValue(const std::string &n)
00179 {
00180   if (n == "Elastic Modulus")
00181     return constant_value;
00182   else if (n == "dEdT Value")
00183     return dEdT_value;
00184   for (int i=0; i<rv.size(); i++) {
00185     if (n == Albany::strint("Elastic Modulus KL Random Variable",i))
00186       return rv[i];
00187   }
00188   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00189            std::endl <<
00190            "Error! Logic error in getting paramter " << n
00191            << " in ElasticModulus::getValue()" << std::endl);
00192   return constant_value;
00193 }
00194 
00195 // **********************************************************************
00196 // **********************************************************************
00197 }
00198 

Generated on Wed Mar 26 2014 18:36:38 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1