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

PHAL_NSPermeabilityTerm_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 NSPermeabilityTerm<EvalT, Traits>::
00017 NSPermeabilityTerm(const Teuchos::ParameterList& p) :
00018   V           (p.get<std::string>                   ("Velocity QP Variable Name"),
00019          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") ),
00020   mu         (p.get<std::string>                   ("Viscosity QP Variable Name"),
00021          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00022   phi         (p.get<std::string>                   ("Porosity QP Variable Name"),
00023          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00024   K         (p.get<std::string>                   ("Permeability QP Variable Name"),
00025          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00026   permTerm   (p.get<std::string>                ("Permeability Term"),
00027          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") )
00028  
00029 {
00030   if (p.isType<bool>("Disable Transient"))
00031     enableTransient = !p.get<bool>("Disable Transient");
00032   else enableTransient = true;
00033 
00034   this->addDependentField(V);
00035   this->addDependentField(mu);
00036   this->addDependentField(phi);
00037   this->addDependentField(K);
00038 
00039   this->addEvaluatedField(permTerm);
00040 
00041   Teuchos::RCP<PHX::DataLayout> vector_dl =
00042     p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout");
00043   std::vector<PHX::DataLayout::size_type> dims;
00044   vector_dl->dimensions(dims);
00045   numNodes = dims[1];
00046   numQPs  = dims[2];
00047   numDims = dims[3];
00048 
00049   this->setName("NSPermeabilityTerm"+PHX::TypeString<EvalT>::value);
00050 }
00051 
00052 //**********************************************************************
00053 template<typename EvalT, typename Traits>
00054 void NSPermeabilityTerm<EvalT, Traits>::
00055 postRegistrationSetup(typename Traits::SetupData d,
00056                       PHX::FieldManager<Traits>& fm)
00057 {
00058   this->utils.setFieldData(V,fm);
00059   this->utils.setFieldData(mu,fm);
00060   this->utils.setFieldData(phi,fm);
00061   this->utils.setFieldData(K,fm);
00062 
00063   this->utils.setFieldData(permTerm,fm); 
00064 }
00065 
00066 //**********************************************************************
00067 template<typename EvalT, typename Traits>
00068 void NSPermeabilityTerm<EvalT, Traits>::
00069 evaluateFields(typename Traits::EvalData workset)
00070 {
00071   for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00072     for (std::size_t qp=0; qp < numQPs; ++qp) {      
00073       for (std::size_t i=0; i < numDims; ++i) {
00074           permTerm(cell,qp,i) = phi(cell,qp)*mu(cell,qp)*V(cell,qp,i)/K(cell,qp);
00075       } 
00076     }
00077   }
00078 }
00079 
00080 }
00081 

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