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