00001
00002
00003
00004
00005
00006
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009
00010 #include "Intrepid_FunctionSpaceTools.hpp"
00011 #include "Intrepid_RealSpaceTools.hpp"
00012
00013 #include <typeinfo>
00014
00015 namespace LCM {
00016
00017
00018 template<typename EvalT, typename Traits>
00019 ThermoPoroPlasticityResidMomentum<EvalT, Traits>::
00020 ThermoPoroPlasticityResidMomentum(const Teuchos::ParameterList& p) :
00021 TotalStress (p.get<std::string> ("Total Stress Name"),
00022 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00023 J (p.get<std::string> ("DetDefGrad Name"),
00024 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00025 Bulk (p.get<std::string> ("Bulk Modulus Name"),
00026 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00027 Temp (p.get<std::string> ("Temperature Name"),
00028 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00029 TempRef (p.get<std::string> ("Reference Temperature Name"),
00030 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00031 alphaSkeleton (p.get<std::string> ("Skeleton Thermal Expansion Name"),
00032 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00033 defgrad (p.get<std::string> ("DefGrad Name"),
00034 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00035 wGradBF (p.get<std::string> ("Weighted Gradient BF Name"),
00036 p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00037 ExResidual (p.get<std::string> ("Residual Name"),
00038 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") )
00039 {
00040 this->addDependentField(TotalStress);
00041 this->addDependentField(wGradBF);
00042 this->addDependentField(J);
00043 this->addDependentField(Bulk);
00044 this->addDependentField(alphaSkeleton);
00045 this->addDependentField(Temp);
00046 this->addDependentField(TempRef);
00047 this->addDependentField(defgrad);
00048
00049
00050 if (p.isType<bool>("Disable Transient"))
00051 enableTransient = !p.get<bool>("Disable Transient");
00052 else enableTransient = true;
00053
00054 if (enableTransient) {
00055
00056 Teuchos::RCP<PHX::DataLayout> node_qp_scalar_dl =
00057 p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Scalar Data Layout");
00058 Teuchos::RCP<PHX::DataLayout> vector_dl =
00059 p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00060
00061 PHX::MDField<MeshScalarT,Cell,Node,QuadPoint> wBF_tmp
00062 (p.get<std::string>("Weighted BF Name"), node_qp_scalar_dl);
00063 wBF = wBF_tmp;
00064 PHX::MDField<ScalarT,Cell,QuadPoint,Dim> uDotDot_tmp
00065 (p.get<std::string>("Time Dependent Variable Name"), vector_dl);
00066 uDotDot = uDotDot_tmp;
00067
00068 this->addDependentField(wBF);
00069 this->addDependentField(uDotDot);
00070
00071 }
00072
00073 this->addEvaluatedField(ExResidual);
00074
00075
00076 this->setName("ThermoPoroPlasticityResidMomentum"+PHX::TypeString<EvalT>::value);
00077
00078 std::vector<PHX::DataLayout::size_type> dims;
00079 wGradBF.fieldTag().dataLayout().dimensions(dims);
00080 numNodes = dims[1];
00081 numQPs = dims[2];
00082 numDims = dims[3];
00083 int worksetSize = dims[0];
00084
00085
00086 F_inv.resize(worksetSize, numQPs, numDims, numDims);
00087 F_invT.resize(worksetSize, numQPs, numDims, numDims);
00088 JF_invT.resize(worksetSize, numQPs, numDims, numDims);
00089
00090 thermoEPS.resize(worksetSize, numQPs, numDims, numDims);
00091 }
00092
00093
00094 template<typename EvalT, typename Traits>
00095 void ThermoPoroPlasticityResidMomentum<EvalT, Traits>::
00096 postRegistrationSetup(typename Traits::SetupData d,
00097 PHX::FieldManager<Traits>& fm)
00098 {
00099 this->utils.setFieldData(TotalStress,fm);
00100 this->utils.setFieldData(wGradBF,fm);
00101 this->utils.setFieldData(J,fm);
00102 this->utils.setFieldData(Bulk,fm);
00103 this->utils.setFieldData(alphaSkeleton,fm);
00104 this->utils.setFieldData(Temp,fm);
00105 this->utils.setFieldData(TempRef,fm);
00106 this->utils.setFieldData(defgrad,fm);
00107 this->utils.setFieldData(ExResidual,fm);
00108
00109 if (enableTransient) this->utils.setFieldData(uDotDot,fm);
00110 if (enableTransient) this->utils.setFieldData(wBF,fm);
00111
00112 }
00113
00114
00115 template<typename EvalT, typename Traits>
00116 void ThermoPoroPlasticityResidMomentum<EvalT, Traits>::
00117 evaluateFields(typename Traits::EvalData workset)
00118 {
00119 typedef Intrepid::FunctionSpaceTools FST;
00120 typedef Intrepid::RealSpaceTools<ScalarT> RST;
00121
00122 RST::inverse(F_inv, defgrad);
00123 RST::transpose(F_invT, F_inv);
00124 FST::scalarMultiplyDataData<ScalarT>(JF_invT, J, F_invT);
00125 FST::scalarMultiplyDataData<ScalarT>(thermoEPS, alphaSkeleton , JF_invT);
00126
00127 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00128 for (std::size_t node=0; node < numNodes; ++node) {
00129 for (std::size_t dim=0; dim<numDims; dim++) {
00130 ExResidual(cell,node,dim)=0.0;
00131 }
00132 for (std::size_t qp=0; qp < numQPs; ++qp) {
00133
00134 dTemp = Temp(cell,qp) - TempRef(cell,qp);
00135
00136 for (std::size_t i=0; i<numDims; i++) {
00137 for (std::size_t dim=0; dim<numDims; dim++) {
00138 ExResidual(cell,node,i) += (TotalStress(cell, qp, i, dim)
00139 - Bulk(cell,qp)*thermoEPS(cell, qp, i, dim)
00140 *dTemp)
00141 * wGradBF(cell, node, qp, dim);
00142 } } } } }
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156 }
00157
00158
00159 }
00160