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

PHAL_CahnHillRhoResid_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 "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009 
00010 #include "Intrepid_FunctionSpaceTools.hpp"
00011 
00012 namespace PHAL {
00013 
00014 //**********************************************************************
00015 template<typename EvalT, typename Traits>
00016 CahnHillRhoResid<EvalT, Traits>::
00017 CahnHillRhoResid(const Teuchos::ParameterList& p) :
00018   wBF         (p.get<std::string>                   ("Weighted BF Name"),
00019          p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Scalar Data Layout") ),
00020   wGradBF     (p.get<std::string>                   ("Weighted Gradient BF Name"),
00021          p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00022   rhoGrad     (p.get<std::string>                   ("Gradient QP Variable Name"),
00023          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") ),
00024   chemTerm    (p.get<std::string>                   ("Chemical Energy Term"),
00025          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00026   rhoResidual (p.get<std::string>                   ("Residual Name"),
00027          p.get<Teuchos::RCP<PHX::DataLayout> >("Node Scalar Data Layout") )
00028 {
00029 
00030   haveNoise = p.get<bool>("Have Noise");
00031 
00032   this->addDependentField(wBF);
00033   this->addDependentField(rhoGrad);
00034   this->addDependentField(wGradBF);
00035   this->addDependentField(chemTerm);
00036   this->addEvaluatedField(rhoResidual);
00037 
00038   if(haveNoise){
00039     noiseTerm = PHX::MDField<ScalarT, Cell, QuadPoint> (p.get<std::string>("Langevin Noise Term"),
00040          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00041     this->addDependentField(noiseTerm);
00042   }
00043 
00044   gamma = p.get<double>("gamma Value");
00045 
00046   Teuchos::RCP<PHX::DataLayout> vector_dl =
00047     p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout");
00048   std::vector<PHX::DataLayout::size_type> dims;
00049   vector_dl->dimensions(dims);
00050   worksetSize = dims[0];
00051   numNodes = dims[1];
00052   numQPs  = dims[2];
00053   numDims = dims[3];
00054 
00055   gamma_term.resize(worksetSize, numQPs, numDims);
00056 
00057   this->setName("CahnHillRhoResid"+PHX::TypeString<EvalT>::value);
00058 
00059 }
00060 
00061 //**********************************************************************
00062 template<typename EvalT, typename Traits>
00063 void CahnHillRhoResid<EvalT, Traits>::
00064 postRegistrationSetup(typename Traits::SetupData d,
00065                       PHX::FieldManager<Traits>& fm)
00066 {
00067   this->utils.setFieldData(wBF,fm);
00068   this->utils.setFieldData(rhoGrad,fm);
00069   this->utils.setFieldData(wGradBF,fm);
00070   this->utils.setFieldData(chemTerm,fm);
00071   if(haveNoise)
00072     this->utils.setFieldData(noiseTerm,fm);
00073 
00074   this->utils.setFieldData(rhoResidual,fm);
00075 }
00076 
00077 //**********************************************************************
00078 template<typename EvalT, typename Traits>
00079 void CahnHillRhoResid<EvalT, Traits>::
00080 evaluateFields(typename Traits::EvalData workset)
00081 {
00082 
00083 // Form Equation 2.2
00084 
00085   typedef Intrepid::FunctionSpaceTools FST;
00086 
00087   for (std::size_t cell=0; cell < workset.numCells; ++cell) 
00088     for (std::size_t qp=0; qp < numQPs; ++qp) 
00089       for (std::size_t i=0; i < numDims; ++i) 
00090 
00091         gamma_term(cell, qp, i) = rhoGrad(cell,qp,i) * gamma; 
00092 
00093   FST::integrate<ScalarT>(rhoResidual, gamma_term, wGradBF, Intrepid::COMP_CPP, false); // "false" overwrites
00094 
00095   FST::integrate<ScalarT>(rhoResidual, chemTerm, wBF, Intrepid::COMP_CPP, true); // "true" sums into
00096 
00097   if(haveNoise)
00098 
00099     FST::integrate<ScalarT>(rhoResidual, noiseTerm, wBF, Intrepid::COMP_CPP, true); // "true" sums into
00100 
00101 
00102 }
00103 
00104 template<typename EvalT, typename Traits>
00105 typename CahnHillRhoResid<EvalT, Traits>::ScalarT&
00106 CahnHillRhoResid<EvalT, Traits>::getValue(const std::string &n) {
00107 
00108   if (n == "gamma") 
00109 
00110     return gamma;
00111 
00112   else 
00113   {
00114     TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter, std::endl <<
00115         "Error! Logic error in getting parameter " << n <<
00116         " in CahnHillRhoResid::getValue()" << std::endl);
00117     return gamma;
00118   }
00119 
00120 }
00121 
00122 //**********************************************************************
00123 }
00124 

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