00001
00002
00003
00004
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 NSForchheimerTerm<EvalT, Traits>::
00017 NSForchheimerTerm(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 rho (p.get<std::string> ("Density 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 F (p.get<std::string> ("Forchheimer QP Variable Name"),
00027 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00028 ForchTerm (p.get<std::string> ("Forchheimer Term"),
00029 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") )
00030
00031 {
00032 if (p.isType<bool>("Disable Transient"))
00033 enableTransient = !p.get<bool>("Disable Transient");
00034 else enableTransient = true;
00035
00036 this->addDependentField(V);
00037 this->addDependentField(rho);
00038 this->addDependentField(phi);
00039 this->addDependentField(K);
00040 this->addDependentField(F);
00041
00042 this->addEvaluatedField(ForchTerm);
00043
00044 Teuchos::RCP<PHX::DataLayout> vector_dl =
00045 p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout");
00046 std::vector<PHX::DataLayout::size_type> dims;
00047 vector_dl->dimensions(dims);
00048 numNodes = dims[1];
00049 numQPs = dims[2];
00050 numDims = dims[3];
00051
00052
00053 normV.resize(dims[0], numQPs);
00054
00055 this->setName("NSForchheimerTerm"+PHX::TypeString<EvalT>::value);
00056 }
00057
00058
00059 template<typename EvalT, typename Traits>
00060 void NSForchheimerTerm<EvalT, Traits>::
00061 postRegistrationSetup(typename Traits::SetupData d,
00062 PHX::FieldManager<Traits>& fm)
00063 {
00064 this->utils.setFieldData(V,fm);
00065 this->utils.setFieldData(rho,fm);
00066 this->utils.setFieldData(phi,fm);
00067 this->utils.setFieldData(K,fm);
00068 this->utils.setFieldData(F,fm);
00069
00070 this->utils.setFieldData(ForchTerm,fm);
00071 }
00072
00073
00074 template<typename EvalT, typename Traits>
00075 void NSForchheimerTerm<EvalT, Traits>::
00076 evaluateFields(typename Traits::EvalData workset)
00077 {
00078 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00079 for (std::size_t qp=0; qp < numQPs; ++qp) {
00080 normV(cell,qp) = 0.0;
00081 for (std::size_t i=0; i < numDims; ++i) {
00082 normV(cell,qp) += V(cell,qp,i)*V(cell,qp,i);
00083 }
00084 if (normV(cell,qp) > 0)
00085 normV(cell,qp) = std::sqrt(normV(cell,qp));
00086 else
00087 normV(cell,qp) = 0.0;
00088 for (std::size_t i=0; i < numDims; ++i) {
00089 ForchTerm(cell,qp,i) = phi(cell,qp)*rho(cell,qp)*F(cell,qp)*normV(cell,qp)*V(cell,qp,i)/std::sqrt(K(cell,qp));
00090 }
00091 }
00092 }
00093 }
00094
00095 }
00096