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

EquivalentInclusionConductivity_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 #include <typeinfo>
00013 
00014 namespace LCM {
00015 
00016 template<typename EvalT, typename Traits>
00017 EquivalentInclusionConductivity<EvalT, Traits>::
00018 EquivalentInclusionConductivity(Teuchos::ParameterList& p) :
00019   effectiveK(p.get<std::string>("QP Variable Name"),
00020      p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout"))
00021 {
00022   Teuchos::ParameterList* elmd_list = 
00023     p.get<Teuchos::ParameterList*>("Parameter List");
00024 
00025   Teuchos::RCP<PHX::DataLayout> vector_dl =
00026     p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00027   std::vector<PHX::DataLayout::size_type> dims;
00028   vector_dl->dimensions(dims);
00029   numQPs  = dims[1];
00030   numDims = dims[2];
00031 
00032   Teuchos::RCP<ParamLib> paramLib = 
00033     p.get< Teuchos::RCP<ParamLib> >("Parameter Library", Teuchos::null);
00034 
00035   std::string type = elmd_list->get("Effective Thermal Conductivity Type", "Constant");
00036   if (type == "Constant") {
00037     is_constant = true;
00038     constant_value = elmd_list->get("Value", 1.0);
00039 
00040     // Add effective thermal conductivity as a Sacado-ized parameter
00041     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00042   "Effective Thermal COnductivity", this, paramLib);
00043   }
00044   else if (type == "Truncated KL Expansion") {
00045     is_constant = false;
00046     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim>
00047       fx(p.get<std::string>("QP Coordinate Vector Name"), vector_dl);
00048     coordVec = fx;
00049     this->addDependentField(coordVec);
00050 
00051     exp_rf_kl = 
00052       Teuchos::rcp(new Stokhos::KL::ExponentialRandomField<MeshScalarT>(*elmd_list));
00053     int num_KL = exp_rf_kl->stochasticDimension();
00054 
00055     // Add KL random variables as Sacado-ized parameters
00056     rv.resize(num_KL);
00057     for (int i=0; i<num_KL; i++) {
00058       std::string ss = Albany::strint("Effective Thermal Conductivity KL Random Variable",i);
00059       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00060       rv[i] = elmd_list->get(ss, 1.0);
00061     }
00062   }
00063   else {
00064     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00065              "Invalid effective thermal conductivity type " << type);
00066   } 
00067 
00068 
00069   if ( p.isType<std::string>("Porosity 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         tporo(p.get<std::string>("Porosity Name"), scalar_dl);
00074       porosity = tporo;
00075       this->addDependentField(porosity);
00076       isPoroElastic = true;
00077     }
00078     else {
00079       isPoroElastic= false;
00080     }
00081 
00082   if ( p.isType<std::string>("Jacobian Name") ) {
00083       Teuchos::RCP<PHX::DataLayout> scalar_dl =
00084         p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00085       PHX::MDField<ScalarT,Cell,QuadPoint>
00086         detf(p.get<std::string>("Jacobian Name"), scalar_dl);
00087       J = detf;
00088       this->addDependentField(J);
00089       isPoroElastic = true;
00090     }
00091   else {
00092     isPoroElastic= false;
00093   }
00094 
00095   condKs = elmd_list->get("Solid Thermal Conducutivity Value", 1.0);
00096   new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00097                           "Solid Thermal Conductivity Value", this, paramLib);
00098   condKf = elmd_list->get("Fluid Thermal Conductivity Value", 100.0);
00099   new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00100                           "Fluid Thermal Conductivity Value", this, paramLib);
00101 
00102   this->addEvaluatedField(effectiveK);
00103   this->setName("Effective Thermal Conductivity"+PHX::TypeString<EvalT>::value);
00104 }
00105 
00106 // **********************************************************************
00107 template<typename EvalT, typename Traits>
00108 void EquivalentInclusionConductivity<EvalT, Traits>::
00109 postRegistrationSetup(typename Traits::SetupData d,
00110                       PHX::FieldManager<Traits>& fm)
00111 {
00112   this->utils.setFieldData(effectiveK,fm);
00113   if (!is_constant)    this->utils.setFieldData(coordVec,fm);
00114   if (isPoroElastic)   this->utils.setFieldData(porosity,fm);
00115   if (isPoroElastic)   this->utils.setFieldData(J,fm);
00116 }
00117 
00118 // **********************************************************************
00119 template<typename EvalT, typename Traits>
00120 void EquivalentInclusionConductivity<EvalT, Traits>::
00121 evaluateFields(typename Traits::EvalData workset)
00122 {
00123   std::size_t numCells = workset.numCells;
00124 
00125   if (is_constant) {
00126     for (std::size_t cell=0; cell < numCells; ++cell) {
00127       for (std::size_t qp=0; qp < numQPs; ++qp) {
00128         effectiveK(cell,qp) = constant_value;
00129       }
00130     }
00131   }
00132   else {
00133     for (std::size_t cell=0; cell < numCells; ++cell) {
00134       for (std::size_t qp=0; qp < numQPs; ++qp) {
00135   Teuchos::Array<MeshScalarT> point(numDims);
00136   for (std::size_t i=0; i<numDims; i++)
00137     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00138   effectiveK(cell,qp) = exp_rf_kl->evaluate(point, rv);
00139       }
00140     }
00141   }
00142 
00143 
00144   if (isPoroElastic) {
00145       for (std::size_t cell=0; cell < numCells; ++cell) {
00146         for (std::size_t qp=0; qp < numQPs; ++qp) {
00147           effectiveK(cell,qp) = porosity(cell,qp)*condKf +
00148                             (J(cell,qp)-porosity(cell,qp))*condKs;
00149        //   effectiveK(cell,qp) = condKf
00150        //                   + ( J(cell,qp) - porosity(cell,qp))*
00151         //                    (condKs - condKf)*condKf/
00152         //                    ((condKs - condKf)*porosity(cell,qp)
00153         //                    + J(cell,qp)*condKf);
00154         }
00155       }
00156     }
00157 
00158 }
00159 
00160 // **********************************************************************
00161 template<typename EvalT,typename Traits>
00162 typename EquivalentInclusionConductivity<EvalT,Traits>::ScalarT&
00163 EquivalentInclusionConductivity<EvalT,Traits>::getValue(const std::string &n)
00164 {
00165   if (n == "Effective Thermal Conductivity")
00166     return constant_value;
00167   else if (n == "Fluid Thermal Conducutivity Value")
00168     return condKf;
00169   else if (n == "Solid Thermal Conducutivity Value")
00170     return condKs;
00171   for (int i=0; i<rv.size(); i++) {
00172     if (n == Albany::strint("Effective Thermal Conductivity KL Random Variable",i))
00173       return rv[i];
00174   }
00175   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00176            std::endl <<
00177            "Error! Logic error in getting paramter " << n
00178            << " in EquivalentInclusionConductivity::getValue()" << std::endl);
00179   return constant_value;
00180 }
00181 
00182 // **********************************************************************
00183 // **********************************************************************
00184 }
00185 

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