Go to the documentation of this file.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 LCM {
00013
00014
00015 template<typename EvalT, typename Traits>
00016 PoroElasticityResidMomentum<EvalT, Traits>::
00017 PoroElasticityResidMomentum(const Teuchos::ParameterList& p) :
00018 TotalStress (p.get<std::string> ("Total Stress Name"),
00019 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor 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 ExResidual (p.get<std::string> ("Residual Name"),
00023 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") )
00024 {
00025 this->addDependentField(TotalStress);
00026 this->addDependentField(wGradBF);
00027 this->addEvaluatedField(ExResidual);
00028
00029 if (p.isType<bool>("Disable Transient"))
00030 enableTransient = !p.get<bool>("Disable Transient");
00031 else enableTransient = true;
00032
00033 if (enableTransient) {
00034
00035 Teuchos::RCP<PHX::DataLayout> node_qp_scalar_dl =
00036 p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Scalar Data Layout");
00037 Teuchos::RCP<PHX::DataLayout> vector_dl =
00038 p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00039
00040 PHX::MDField<MeshScalarT,Cell,Node,QuadPoint> wBF_tmp
00041 (p.get<std::string>("Weighted BF Name"), node_qp_scalar_dl);
00042 wBF = wBF_tmp;
00043 PHX::MDField<ScalarT,Cell,QuadPoint,Dim> uDotDot_tmp
00044 (p.get<std::string>("Time Dependent Variable Name"), vector_dl);
00045 uDotDot = uDotDot_tmp;
00046
00047 this->addDependentField(wBF);
00048 this->addDependentField(uDotDot);
00049 }
00050
00051
00052 this->setName("PoroElasticityResidMomentum"+PHX::TypeString<EvalT>::value);
00053
00054 std::vector<PHX::DataLayout::size_type> dims;
00055 wGradBF.fieldTag().dataLayout().dimensions(dims);
00056 numNodes = dims[1];
00057 numQPs = dims[2];
00058 numDims = dims[3];
00059 }
00060
00061
00062 template<typename EvalT, typename Traits>
00063 void PoroElasticityResidMomentum<EvalT, Traits>::
00064 postRegistrationSetup(typename Traits::SetupData d,
00065 PHX::FieldManager<Traits>& fm)
00066 {
00067 this->utils.setFieldData(TotalStress,fm);
00068 this->utils.setFieldData(wGradBF,fm);
00069
00070 this->utils.setFieldData(ExResidual,fm);
00071
00072 if (enableTransient) this->utils.setFieldData(uDotDot,fm);
00073 if (enableTransient) this->utils.setFieldData(wBF,fm);
00074
00075 }
00076
00077
00078 template<typename EvalT, typename Traits>
00079 void PoroElasticityResidMomentum<EvalT, Traits>::
00080 evaluateFields(typename Traits::EvalData workset)
00081 {
00082 typedef Intrepid::FunctionSpaceTools FST;
00083
00084 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00085 for (std::size_t node=0; node < numNodes; ++node) {
00086 for (std::size_t dim=0; dim<numDims; dim++) ExResidual(cell,node,dim)=0.0;
00087 for (std::size_t qp=0; qp < numQPs; ++qp) {
00088 for (std::size_t i=0; i<numDims; i++) {
00089 for (std::size_t dim=0; dim<numDims; dim++) {
00090 ExResidual(cell,node,i) += TotalStress(cell, qp, i, dim) * wGradBF(cell, node, qp, dim);
00091 } } } } }
00092
00093
00094 if (workset.transientTerms && enableTransient)
00095 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00096 for (std::size_t node=0; node < numNodes; ++node) {
00097 for (std::size_t qp=0; qp < numQPs; ++qp) {
00098 for (std::size_t i=0; i<numDims; i++) {
00099 ExResidual(cell,node,i) += uDotDot(cell, qp, i) * wBF(cell, node, qp);
00100 } } } }
00101
00102
00103
00104
00105 }
00106
00107
00108 }
00109