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

BiotModulus_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 BiotModulus<EvalT, Traits>::
00017 BiotModulus(Teuchos::ParameterList& p) :
00018   biotModulus(p.get<std::string>("Biot Modulus 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("Biot Modulus Type", "Constant");
00035   if (type == "Constant") {
00036     is_constant = true;
00037     constant_value = elmd_list->get("Value", 1.0);
00038 
00039     // Add Biot Modulus as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Biot 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("Biot 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 Biot modulus type " << type);
00065   } 
00066 
00067   if ( p.isType<std::string>("Porosity Name") ) {
00068      Teuchos::RCP<PHX::DataLayout> scalar_dl =
00069        p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00070      PHX::MDField<ScalarT,Cell,QuadPoint>
00071        tp(p.get<std::string>("Porosity Name"), scalar_dl);
00072      porosity = tp;
00073      this->addDependentField(porosity);
00074      isPoroElastic = true;
00075      FluidBulkModulus = elmd_list->get("Fluid Bulk Modulus Value", 10.0e9);
00076      new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00077                                  "Skeleton Bulk Modulus Parameter Value", this, paramLib);
00078      GrainBulkModulus = elmd_list->get("Grain Bulk Modulus Value", 10.0e12); // typically Kgrain >> Kskeleton
00079      new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00080                                      "Grain Bulk Modulus Value", this, paramLib);
00081    }
00082    else {
00083      isPoroElastic=false;
00084      FluidBulkModulus=10.0e9; // temp value..need to change
00085      GrainBulkModulus = 10.0e12;  // temp value need to change
00086    }
00087 
00088   if ( p.isType<std::string>("Biot Coefficient Name") ) {
00089      Teuchos::RCP<PHX::DataLayout> scalar_dl =
00090        p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00091      PHX::MDField<ScalarT,Cell,QuadPoint>
00092        btp(p.get<std::string>("Biot Coefficient Name"), scalar_dl);
00093      biotCoefficient = btp;
00094      this->addDependentField(biotCoefficient);
00095   }
00096 
00097 
00098   this->addEvaluatedField(biotModulus);
00099   this->setName("Biot Modulus"+PHX::TypeString<EvalT>::value);
00100 }
00101 
00102 // **********************************************************************
00103 template<typename EvalT, typename Traits>
00104 void BiotModulus<EvalT, Traits>::
00105 postRegistrationSetup(typename Traits::SetupData d,
00106                       PHX::FieldManager<Traits>& fm)
00107 {
00108   this->utils.setFieldData(biotModulus,fm);
00109   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00110 //  if (isThermoElastic) this->utils.setFieldData(Temperature,fm);
00111   if(isPoroElastic) this->utils.setFieldData(porosity,fm);
00112   if(isPoroElastic) this->utils.setFieldData(biotCoefficient,fm);
00113 }
00114 
00115 // **********************************************************************
00116 template<typename EvalT, typename Traits>
00117 void BiotModulus<EvalT, Traits>::
00118 evaluateFields(typename Traits::EvalData workset)
00119 {
00120   std::size_t numCells = workset.numCells;
00121 
00122   if (is_constant) {
00123     for (std::size_t cell=0; cell < numCells; ++cell) {
00124       for (std::size_t qp=0; qp < numQPs; ++qp) {
00125         biotModulus(cell,qp) = constant_value;
00126       }
00127     }
00128   }
00129   else {
00130     for (std::size_t cell=0; cell < numCells; ++cell) {
00131       for (std::size_t qp=0; qp < numQPs; ++qp) {
00132   Teuchos::Array<MeshScalarT> point(numDims);
00133   for (std::size_t i=0; i<numDims; i++)
00134     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00135       biotModulus(cell,qp) = exp_rf_kl->evaluate(point, rv);
00136       }
00137     }
00138   }
00139   if (isPoroElastic) {
00140     for (std::size_t cell=0; cell < numCells; ++cell) {
00141       for (std::size_t qp=0; qp < numQPs; ++qp) {
00142           // 1/M = (B-phi)/Ks + phi/Kf
00143           biotModulus(cell,qp) = 1/((biotCoefficient(cell,qp) - porosity(cell,qp))/GrainBulkModulus
00144                                 + porosity(cell,qp)/FluidBulkModulus);
00145 
00146       }
00147     }
00148   }
00149 }
00150 
00151 // **********************************************************************
00152 template<typename EvalT,typename Traits>
00153 typename BiotModulus<EvalT,Traits>::ScalarT&
00154 BiotModulus<EvalT,Traits>::getValue(const std::string &n)
00155 {
00156   if (n == "Biot Modulus")
00157     return constant_value;
00158   else if (n == "Fluid Bulk Modulus Value")
00159      return FluidBulkModulus;
00160    else if (n == "Grain Bulk Modulus Value")
00161        return GrainBulkModulus;
00162   for (int i=0; i<rv.size(); i++) {
00163     if (n == Albany::strint("Biot Modulus KL Random Variable",i))
00164       return rv[i];
00165   }
00166   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00167          std::endl <<
00168          "Error! Logic error in getting paramter " << n
00169          << " in BiotModulus::getValue()" << std::endl);
00170   return constant_value;
00171 }
00172 
00173 // **********************************************************************
00174 // **********************************************************************
00175 }
00176 

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