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

SaturationExponent_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 SaturationExponent<EvalT, Traits>::
00017 SaturationExponent(Teuchos::ParameterList& p) :
00018   satExp(p.get<std::string>("Saturation Exponent Name"),
00019    p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout"))
00020 {
00021   Teuchos::ParameterList* satExp_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 = satExp_list->get("Saturation Exponent Type", "Constant");
00035   if (type == "Constant") {
00036     is_constant = true;
00037     constant_value = satExp_list->get("Value", 0.0);
00038 
00039     // Add Saturation Exponent as a Sacado-ized parameter
00040     new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041   "Saturation Exponent", 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>(*satExp_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("Saturation Exponent KL Random Variable",i);
00058       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00059       rv[i] = satExp_list->get(ss, 0.0);
00060     }
00061   }
00062   else {
00063     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00064              "Invalid saturation exponent type " << type);
00065   } 
00066 
00067   this->addEvaluatedField(satExp);
00068   this->setName("Saturation Exponent"+PHX::TypeString<EvalT>::value);
00069 }
00070 
00071 // **********************************************************************
00072 template<typename EvalT, typename Traits>
00073 void SaturationExponent<EvalT, Traits>::
00074 postRegistrationSetup(typename Traits::SetupData d,
00075                       PHX::FieldManager<Traits>& fm)
00076 {
00077   this->utils.setFieldData(satExp,fm);
00078   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00079 }
00080 
00081 // **********************************************************************
00082 template<typename EvalT, typename Traits>
00083 void SaturationExponent<EvalT, Traits>::
00084 evaluateFields(typename Traits::EvalData workset)
00085 {
00086   std::size_t numCells = workset.numCells;
00087 
00088   if (is_constant) {
00089     for (std::size_t cell=0; cell < numCells; ++cell) {
00090       for (std::size_t qp=0; qp < numQPs; ++qp) {
00091   satExp(cell,qp) = constant_value;
00092       }
00093     }
00094   }
00095   else {
00096     for (std::size_t cell=0; cell < numCells; ++cell) {
00097       for (std::size_t qp=0; qp < numQPs; ++qp) {
00098   Teuchos::Array<MeshScalarT> point(numDims);
00099   for (std::size_t i=0; i<numDims; i++)
00100     point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00101   satExp(cell,qp) = exp_rf_kl->evaluate(point, rv);
00102       }
00103     }
00104   }
00105 }
00106 
00107 // **********************************************************************
00108 template<typename EvalT,typename Traits>
00109 typename SaturationExponent<EvalT,Traits>::ScalarT& 
00110 SaturationExponent<EvalT,Traits>::getValue(const std::string &n)
00111 {
00112   if (n == "Saturation Exponent")
00113     return constant_value;
00114   for (int i=0; i<rv.size(); i++) {
00115     if (n == Albany::strint("Saturation Exponent KL Random Variable",i))
00116       return rv[i];
00117   }
00118   TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00119            std::endl <<
00120            "Error! Logic error in getting paramter " << n
00121            << " in SaturationExponent::getValue()" << std::endl);
00122   return constant_value;
00123 }
00124 
00125 // **********************************************************************
00126 // **********************************************************************
00127 }
00128 

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