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

PHAL_Absorption_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 namespace PHAL {
00013 
00014 template<typename EvalT, typename Traits>
00015 Absorption<EvalT, Traits>::
00016 Absorption(Teuchos::ParameterList& p) :
00017   absorption(p.get<std::string>("QP Variable Name"),
00018         p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout"))
00019 {
00020   Teuchos::ParameterList* cond_list = 
00021     p.get<Teuchos::ParameterList*>("Parameter List");
00022 
00023   Teuchos::RCP<PHX::DataLayout> vector_dl =
00024     p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00025   std::vector<PHX::DataLayout::size_type> dims;
00026   vector_dl->dimensions(dims);
00027   numQPs  = dims[1];
00028   numDims = dims[2];
00029 
00030   std::string type = cond_list->get("Absorption Type", "Constant");
00031   if (type == "Constant") {
00032     is_constant = true;
00033     constant_value = cond_list->get("Value", 1.0);
00034 
00035     // Add absorption as a Sacado-ized parameter
00036     Teuchos::RCP<ParamLib> paramLib = 
00037       p.get< Teuchos::RCP<ParamLib> >("Parameter Library", Teuchos::null);
00038       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00039   "Absorption", this, paramLib);
00040   }
00041   else if (type == "Truncated KL Expansion") {
00042     is_constant = false;
00043     Teuchos::RCP<PHX::DataLayout> scalar_dl =
00044       p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
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>(*cond_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     Teuchos::RCP<ParamLib> paramLib = 
00057       p.get< Teuchos::RCP<ParamLib> >("Parameter Library", Teuchos::null);
00058     for (int i=0; i<num_KL; i++) {
00059       std::string ss = Albany::strint("Absorption KL Random Variable",i);
00060       new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00061       rv[i] = cond_list->get(ss, 0.0);
00062     }
00063   }
00064   else {
00065     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00066            "Invalid absorption type " << type);
00067   } 
00068 
00069   this->addEvaluatedField(absorption);
00070   this->setName("Absorption"+PHX::TypeString<EvalT>::value);
00071 }
00072 
00073 // **********************************************************************
00074 template<typename EvalT, typename Traits>
00075 void Absorption<EvalT, Traits>::
00076 postRegistrationSetup(typename Traits::SetupData d,
00077                       PHX::FieldManager<Traits>& fm)
00078 {
00079   this->utils.setFieldData(absorption,fm);
00080   if (!is_constant) this->utils.setFieldData(coordVec,fm);
00081 }
00082 
00083 // **********************************************************************
00084 template<typename EvalT, typename Traits>
00085 void Absorption<EvalT, Traits>::
00086 evaluateFields(typename Traits::EvalData workset)
00087 {
00088   if (is_constant) {
00089     for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00090       for (std::size_t qp=0; qp < numQPs; ++qp) {
00091   absorption(cell,qp) = constant_value;
00092       }
00093     }
00094   }
00095   /*else {
00096     for (std::size_t cell=0; cell < workset.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   thermalCond(cell,qp) = exp_rf_kl->evaluate(point, rv);
00102       }
00103     }
00104   }*/
00105 }
00106 
00107 // **********************************************************************
00108 template<typename EvalT,typename Traits>
00109 typename Absorption<EvalT,Traits>::ScalarT& 
00110 Absorption<EvalT,Traits>::getValue(const std::string &n)
00111 {
00112   if (is_constant)
00113     return constant_value;
00114   /*for (int i=0; i<rv.size(); i++) {
00115     if (n == Albany::strint("Thermal Conductivity 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 Absorption::getValue()" << std::endl);
00122   return constant_value;
00123 }
00124 
00125 // **********************************************************************
00126 // **********************************************************************
00127 }
00128 

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