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 NSContinuityResid<EvalT, Traits>::
00017 NSContinuityResid(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 VGrad (p.get<std::string> ("Gradient QP Variable Name"),
00021 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00022 rho (p.get<std::string> ("Density QP Variable Name"),
00023 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00024
00025 CResidual (p.get<std::string> ("Residual Name"),
00026 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Scalar Data Layout") ),
00027 havePSPG(p.get<bool>("Have PSPG"))
00028 {
00029 this->addDependentField(wBF);
00030 this->addDependentField(VGrad);
00031 this->addDependentField(rho);
00032 if (havePSPG) {
00033 wGradBF = PHX::MDField<MeshScalarT,Cell,Node,QuadPoint,Dim>(
00034 p.get<std::string>("Weighted Gradient BF Name"),
00035 p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") );
00036 TauM = PHX::MDField<ScalarT,Cell,QuadPoint>(
00037 p.get<std::string>("Tau M Name"),
00038 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") );
00039 Rm = PHX::MDField<ScalarT,Cell,QuadPoint,Dim>(
00040 p.get<std::string>("Rm Name"),
00041 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout") );
00042 this->addDependentField(wGradBF);
00043 this->addDependentField(TauM);
00044 this->addDependentField(Rm);
00045 }
00046
00047 this->addEvaluatedField(CResidual);
00048
00049 Teuchos::RCP<PHX::DataLayout> vector_dl =
00050 p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout");
00051 std::vector<PHX::DataLayout::size_type> dims;
00052 vector_dl->dimensions(dims);
00053 numNodes = dims[1];
00054 numQPs = dims[2];
00055 numDims = dims[3];
00056
00057
00058 divergence.resize(dims[0], numQPs);
00059
00060 this->setName("NSContinuityResid"+PHX::TypeString<EvalT>::value);
00061 }
00062
00063
00064 template<typename EvalT, typename Traits>
00065 void NSContinuityResid<EvalT, Traits>::
00066 postRegistrationSetup(typename Traits::SetupData d,
00067 PHX::FieldManager<Traits>& fm)
00068 {
00069 this->utils.setFieldData(wBF,fm);
00070 this->utils.setFieldData(VGrad,fm);
00071 this->utils.setFieldData(rho,fm);
00072 if (havePSPG) {
00073 this->utils.setFieldData(wGradBF,fm);
00074 this->utils.setFieldData(TauM,fm);
00075 this->utils.setFieldData(Rm,fm);
00076 }
00077
00078 this->utils.setFieldData(CResidual,fm);
00079 }
00080
00081
00082 template<typename EvalT, typename Traits>
00083 void NSContinuityResid<EvalT, Traits>::
00084 evaluateFields(typename Traits::EvalData workset)
00085 {
00086 typedef Intrepid::FunctionSpaceTools FST;
00087
00088 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00089 for (std::size_t qp=0; qp < numQPs; ++qp) {
00090 divergence(cell,qp) = 0.0;
00091 for (std::size_t i=0; i < numDims; ++i) {
00092 divergence(cell,qp) += rho(cell,qp)*VGrad(cell,qp,i,i);
00093 }
00094 }
00095 }
00096
00097 FST::integrate<ScalarT>(CResidual, divergence, wBF, Intrepid::COMP_CPP,
00098 false);
00099
00100 if (havePSPG) {
00101 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00102 for (std::size_t node=0; node < numNodes; ++node) {
00103 for (std::size_t qp=0; qp < numQPs; ++qp) {
00104 for (std::size_t j=0; j < numDims; ++j) {
00105 CResidual(cell,node) +=
00106 rho(cell,qp)*TauM(cell,qp)*Rm(cell,qp,j)*wGradBF(cell,node,qp,j);
00107 }
00108 }
00109 }
00110 }
00111 }
00112
00113 }
00114
00115
00116 }
00117