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

KCPermeability_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 KCPermeability<EvalT, Traits>::
00017 KCPermeability(Teuchos::ParameterList& p) :
00018   kcPermeability(p.get<std::string>("Kozeny-Carman 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("Kozeny-Carman 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 Kozeny-Carman Permeability as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Kozeny-Carman 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("Kozeny-Carman 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 Kozeny-Carman 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 
00086   this->addEvaluatedField(kcPermeability);
00087   this->setName("Kozeny-Carman Permeability"+PHX::TypeString<EvalT>::value);
00088 }
00089 
00090 // **********************************************************************
00091 template<typename EvalT, typename Traits>
00092 void KCPermeability<EvalT, Traits>::
00093 postRegistrationSetup(typename Traits::SetupData d,
00094                       PHX::FieldManager<Traits>& fm)
00095 {
00096   this->utils.setFieldData(kcPermeability,fm);
00097   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00098   if (isPoroElastic) this->utils.setFieldData(porosity,fm);
00099 }
00100 
00101 // **********************************************************************
00102 template<typename EvalT, typename Traits>
00103 void KCPermeability<EvalT, Traits>::
00104 evaluateFields(typename Traits::EvalData workset)
00105 {
00106   std::size_t numCells = workset.numCells;
00107 
00108   if (is_constant) {
00109     for (std::size_t cell=0; cell < numCells; ++cell) {
00110       for (std::size_t qp=0; qp < numQPs; ++qp) {
00111         kcPermeability(cell,qp) = constant_value;
00112       }
00113     }
00114   }
00115   else {
00116     for (std::size_t cell=0; cell < numCells; ++cell) {
00117       for (std::size_t qp=0; qp < numQPs; ++qp) {
00118   Teuchos::Array<MeshScalarT> point(numDims);
00119   for (std::size_t i=0; i<numDims; i++)
00120     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00121       kcPermeability(cell,qp) = exp_rf_kl->evaluate(point, rv);
00122       }
00123     }
00124   }
00125   if (isPoroElastic) {
00126     for (std::size_t cell=0; cell < numCells; ++cell) {
00127       for (std::size_t qp=0; qp < numQPs; ++qp) {
00128         // Cozeny Karman permeability equation
00129         kcPermeability(cell,qp) = constant_value*porosity(cell,qp)*porosity(cell,qp)*porosity(cell,qp)/
00130                               (1-porosity(cell,qp)*porosity(cell,qp));
00131       }
00132     }
00133   }
00134 }
00135 
00136 // **********************************************************************
00137 template<typename EvalT,typename Traits>
00138 typename KCPermeability<EvalT,Traits>::ScalarT&
00139 KCPermeability<EvalT,Traits>::getValue(const std::string &n)
00140 {
00141   if (n == "Kozeny-Carman Permeability")
00142     return constant_value;
00143   for (int i=0; i<rv.size(); i++) {
00144     if (n == Albany::strint("Kozeny-Carman Permeability KL Random Variable",i))
00145       return rv[i];
00146   }
00147   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00148          std::endl <<
00149          "Error! Logic error in getting paramter " << n
00150          << " in KCPermeability::getValue()" << std::endl);
00151   return constant_value;
00152 }
00153 
00154 // **********************************************************************
00155 // **********************************************************************
00156 }
00157 

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