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

VanGenuchtenPermeability_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 VanGenuchtenPermeability<EvalT, Traits>::
00017 VanGenuchtenPermeability(Teuchos::ParameterList& p) :
00018   vgPermeability(p.get<std::string>("Van Genuchten Permeability 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("Van Genuchten Permeability Type", "Constant");
00035   if (type == "Constant") {
00036     is_constant = true;
00037     constant_value = elmd_list->get("Value", 1.0e-5); // default value=1, identical to Terzaghi stress
00038 
00039     // Add Van Genuchten Permeability as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Van Genuchten Permeability", 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("Van Genuchten Permeability 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 Van Genuchten Permeability type " << type);
00065   } 
00066 
00067   // Optional dependence on Temperature (E = E_ + dEdT * T)
00068   // Switched ON by sending Temperature field in p
00069 
00070   if ( p.isType<std::string>("Porosity Name") ) {
00071     Teuchos::RCP<PHX::DataLayout> scalar_dl =
00072       p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00073     PHX::MDField<ScalarT,Cell,QuadPoint>
00074       tp(p.get<std::string>("Porosity Name"), scalar_dl);
00075     porosity = tp;
00076     this->addDependentField(porosity);
00077     isPoroElastic = true;
00078 
00079   }
00080   else {
00081     isPoroElastic = false;
00082 
00083   }
00084 
00085   if ( p.isType<std::string>("QP Pore Pressure Name") ) {
00086          Teuchos::RCP<PHX::DataLayout> scalar_dl =
00087            p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00088          PHX::MDField<ScalarT,Cell,QuadPoint>
00089            ppn(p.get<std::string>("QP Pore Pressure Name"), scalar_dl);
00090          porePressure = ppn;
00091          isPoroElastic = true;
00092          this->addDependentField(porePressure);
00093 
00094          waterUnitWeight = elmd_list->get("Water Unit Weight Value", 9810.0); // typically Kgrain >> Kskeleton
00095                     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00096           "Water Unit Weight Value", this, paramLib);
00097   }
00098 
00099 
00100   this->addEvaluatedField(vgPermeability);
00101   this->setName("Van Genuchten Permeability"+PHX::TypeString<EvalT>::value);
00102 }
00103 
00104 // **********************************************************************
00105 template<typename EvalT, typename Traits>
00106 void VanGenuchtenPermeability<EvalT, Traits>::
00107 postRegistrationSetup(typename Traits::SetupData d,
00108                       PHX::FieldManager<Traits>& fm)
00109 {
00110   this->utils.setFieldData(vgPermeability,fm);
00111   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00112   if (isPoroElastic) this->utils.setFieldData(porosity,fm);
00113   if (isPoroElastic) this->utils.setFieldData(porePressure,fm);
00114 }
00115 
00116 // **********************************************************************
00117 template<typename EvalT, typename Traits>
00118 void VanGenuchtenPermeability<EvalT, Traits>::
00119 evaluateFields(typename Traits::EvalData workset)
00120 {
00121   std::size_t numCells = workset.numCells;
00122 
00123   if (is_constant) {
00124     for (std::size_t cell=0; cell < numCells; ++cell) {
00125       for (std::size_t qp=0; qp < numQPs; ++qp) {
00126         vgPermeability(cell,qp) = constant_value;
00127       }
00128     }
00129   }
00130   else {
00131     for (std::size_t cell=0; cell < numCells; ++cell) {
00132       for (std::size_t qp=0; qp < numQPs; ++qp) {
00133   Teuchos::Array<MeshScalarT> point(numDims);
00134   for (std::size_t i=0; i<numDims; i++)
00135     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00136       vgPermeability(cell,qp) = exp_rf_kl->evaluate(point, rv);
00137       }
00138     }
00139   }
00140   if (isPoroElastic) {
00141     for (std::size_t cell=0; cell < numCells; ++cell) {
00142       for (std::size_t qp=0; qp < numQPs; ++qp) {
00143         // van Genuchten permeability equation
00144         vgPermeability(cell,qp) = constant_value*porosity(cell,qp)*porosity(cell,qp)*porosity(cell,qp)/
00145                               (  1.0-porosity(cell,qp)*porosity(cell,qp) )/
00146                               std::pow( 1.0 + std::pow(50.0*porePressure(cell,qp)/waterUnitWeight, 4.0),0.9);
00147       }
00148     }
00149   }
00150 }
00151 
00152 // **********************************************************************
00153 template<typename EvalT,typename Traits>
00154 typename VanGenuchtenPermeability<EvalT,Traits>::ScalarT&
00155 VanGenuchtenPermeability<EvalT,Traits>::getValue(const std::string &n)
00156 {
00157   if (n == "Van Genuchten Permeability")
00158     return constant_value;
00159   else if (n == "Water Unit Weight Value")
00160       return waterUnitWeight;
00161   for (int i=0; i<rv.size(); i++) {
00162     if (n == Albany::strint("Van Genuchten Permeability KL Random Variable",i))
00163       return rv[i];
00164   }
00165   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00166          std::endl <<
00167          "Error! Logic error in getting parameter " << n
00168          << " in VanGenuchtenPermeability::getValue()" << std::endl);
00169   return constant_value;
00170 }
00171 
00172 // **********************************************************************
00173 // **********************************************************************
00174 }
00175 

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