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

RecoveryModulus_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 RecoveryModulus<EvalT, Traits>::
00017 RecoveryModulus(Teuchos::ParameterList& p) :
00018   recoveryModulus(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("Recovery Modulus Type", "Constant");
00035   if (type == "Constant") {
00036     is_constant = true;
00037     constant_value = elmd_list->get("Value", 1.0);
00038 
00039     // Add Recovery Modulus as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Recovery Modulus", this, paramLib);
00042   }
00043   else if (type == "Truncated KL Expansion") {
00044     is_constant = false;
00045     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim>
00046       fx(p.get<std::string>("QP Coordinate Vector Name"), vector_dl);
00047     coordVec = fx;
00048     this->addDependentField(coordVec);
00049 
00050     exp_rf_kl = 
00051       Teuchos::rcp(new Stokhos::KL::ExponentialRandomField<MeshScalarT>(*elmd_list));
00052     int num_KL = exp_rf_kl->stochasticDimension();
00053 
00054     // Add KL random variables as Sacado-ized parameters
00055     rv.resize(num_KL);
00056     for (int i=0; i<num_KL; i++) {
00057       std::string ss = Albany::strint("Recovery Modulus KL Random Variable",i);
00058       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00059       rv[i] = elmd_list->get(ss, 0.0);
00060     }
00061   }
00062   else {
00063     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00064              "Invalid recovery modulus type " << type);
00065   } 
00066 
00067   // Optional dependence on Temperature
00068   // Switched ON by sending Temperature field in p
00069   if ( p.isType<std::string>("QP Temperature Name") ) {
00070     Teuchos::RCP<PHX::DataLayout> scalar_dl =
00071       p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00072     PHX::MDField<ScalarT,Cell,QuadPoint>
00073       tmp(p.get<std::string>("QP Temperature Name"), scalar_dl);
00074     Temperature = tmp;
00075     this->addDependentField(Temperature);
00076     isThermoElastic = true;
00077     c1 = elmd_list->get("c1 Value", 0.0);
00078     c2 = elmd_list->get("c2 Value", 0.0);
00079     refTemp = p.get<RealType>("Reference Temperature", 0.0);
00080     new Sacado::ParameterRegistration<EvalT, SPL_Traits>("c1 Value", this, paramLib);
00081     new Sacado::ParameterRegistration<EvalT, SPL_Traits>("c2 Value", this, paramLib);
00082   }
00083   else {
00084     isThermoElastic=false;
00085     c1=0.0;
00086     c2=0.0;
00087   }
00088 
00089   this->addEvaluatedField(recoveryModulus);
00090   this->setName("Recovery Modulus"+PHX::TypeString<EvalT>::value);
00091 }
00092 
00093 // **********************************************************************
00094 template<typename EvalT, typename Traits>
00095 void RecoveryModulus<EvalT, Traits>::
00096 postRegistrationSetup(typename Traits::SetupData d,
00097                       PHX::FieldManager<Traits>& fm)
00098 {
00099   this->utils.setFieldData(recoveryModulus,fm);
00100   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00101   if (isThermoElastic) this->utils.setFieldData(Temperature,fm);
00102 }
00103 
00104 // **********************************************************************
00105 template<typename EvalT, typename Traits>
00106 void RecoveryModulus<EvalT, Traits>::
00107 evaluateFields(typename Traits::EvalData workset)
00108 {
00109   std::size_t numCells = workset.numCells;
00110 
00111   if (is_constant) {
00112     for (std::size_t cell=0; cell < numCells; ++cell) {
00113       for (std::size_t qp=0; qp < numQPs; ++qp) {
00114         recoveryModulus(cell,qp) = constant_value;
00115       }
00116     }
00117   }
00118   else {
00119     for (std::size_t cell=0; cell < numCells; ++cell) {
00120       for (std::size_t qp=0; qp < numQPs; ++qp) {
00121         Teuchos::Array<MeshScalarT> point(numDims);
00122         for (std::size_t i=0; i<numDims; i++)
00123           point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00124             recoveryModulus(cell,qp) = exp_rf_kl->evaluate(point, rv);
00125       }
00126     }
00127   }
00128   if (isThermoElastic) {
00129     for (std::size_t cell=0; cell < numCells; ++cell) {
00130       for (std::size_t qp=0; qp < numQPs; ++qp) {
00131         if (Temperature(cell,qp) != 0)
00132           recoveryModulus(cell,qp) = c1 * std::exp(c2 / Temperature(cell,qp));
00133         else
00134           recoveryModulus(cell,qp) = 0.0;
00135       }
00136     }
00137   }
00138 }
00139 
00140 // **********************************************************************
00141 template<typename EvalT,typename Traits>
00142 typename RecoveryModulus<EvalT,Traits>::ScalarT&
00143 RecoveryModulus<EvalT,Traits>::getValue(const std::string &n)
00144 {
00145   if (n == "Recovery Modulus")
00146     return constant_value;
00147   else if (n == "c1 Value")
00148     return c1;
00149   else if (n == "c2 Value")
00150   return c2;
00151 
00152   for (int i=0; i<rv.size(); i++) {
00153     if (n == Albany::strint("Recovery Modulus KL Random Variable",i))
00154       return rv[i];
00155   }
00156   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00157            std::endl <<
00158            "Error! Logic error in getting paramter " << n
00159            << " in RecoveryModulus::getValue()" << std::endl);
00160   return constant_value;
00161 }
00162 
00163 // **********************************************************************
00164 // **********************************************************************
00165 }
00166 

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