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

PHAL_NSTauT_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 NSTauT<EvalT, Traits>::
00017 NSTauT(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   ThermalCond (p.get<std::string>                   ("Thermal Conductivity Name"),
00021          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00022   Gc            (p.get<std::string>                   ("Contravarient Metric Tensor Name"),  
00023                  p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00024   rho       (p.get<std::string>                   ("Density QP Variable Name"),
00025          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00026   Cp       (p.get<std::string>                  ("Specific Heat QP Variable Name"),
00027          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00028   TauT            (p.get<std::string>                 ("Tau T Name"),
00029                  p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") )
00030   
00031 {
00032   this->addDependentField(V);
00033   this->addDependentField(ThermalCond);
00034   this->addDependentField(Gc);
00035   this->addDependentField(rho);
00036   this->addDependentField(Cp);
00037  
00038   this->addEvaluatedField(TauT);
00039 
00040   Teuchos::RCP<PHX::DataLayout> vector_dl =
00041     p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00042   std::vector<PHX::DataLayout::size_type> dims;
00043   vector_dl->dimensions(dims);
00044   numQPs  = dims[1];
00045   numDims = dims[2];
00046 
00047   // Allocate workspace
00048   normGc.resize(dims[0], numQPs);
00049 
00050   this->setName("NSTauT"+PHX::TypeString<EvalT>::value);
00051 }
00052 
00053 //**********************************************************************
00054 template<typename EvalT, typename Traits>
00055 void NSTauT<EvalT, Traits>::
00056 postRegistrationSetup(typename Traits::SetupData d,
00057                       PHX::FieldManager<Traits>& fm)
00058 {
00059   this->utils.setFieldData(V,fm);
00060   this->utils.setFieldData(ThermalCond,fm);
00061   this->utils.setFieldData(Gc,fm);
00062   this->utils.setFieldData(rho,fm);
00063   this->utils.setFieldData(Cp,fm);
00064   
00065   this->utils.setFieldData(TauT,fm);
00066 }
00067 
00068 //**********************************************************************
00069 template<typename EvalT, typename Traits>
00070 void NSTauT<EvalT, Traits>::
00071 evaluateFields(typename Traits::EvalData workset)
00072 { 
00073     for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00074       for (std::size_t qp=0; qp < numQPs; ++qp) {       
00075         TauT(cell,qp) = 0.0;
00076         normGc(cell,qp) = 0.0;
00077         for (std::size_t i=0; i < numDims; ++i) {
00078           for (std::size_t j=0; j < numDims; ++j) {
00079             TauT(cell,qp) += rho(cell,qp) * Cp(cell,qp) * rho(cell,qp) * Cp(cell,qp)* V(cell,qp,i)*Gc(cell,qp,i,j)*V(cell,qp,j);
00080             normGc(cell,qp) += Gc(cell,qp,i,j)*Gc(cell,qp,i,j);          
00081           }
00082         }
00083         TauT(cell,qp) += 12.*ThermalCond(cell,qp)*ThermalCond(cell,qp)*std::sqrt(normGc(cell,qp));
00084         TauT(cell,qp) = 1./std::sqrt(TauT(cell,qp));
00085       }
00086     }
00087   
00088 
00089 }
00090 
00091 //**********************************************************************
00092 }
00093 

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