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 #include "Sacado_ParameterRegistration.hpp"
00013
00014 namespace LCM {
00015
00016
00017 template<typename EvalT, typename Traits>
00018 ThermoMechanicalMomentumResidual<EvalT, Traits>::
00019 ThermoMechanicalMomentumResidual(const Teuchos::ParameterList& p) :
00020 stress (p.get<std::string> ("Stress Name"),
00021 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00022 J (p.get<std::string> ("DetDefGrad Name"),
00023 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00024 defgrad (p.get<std::string> ("DefGrad Name"),
00025 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00026 wGradBF (p.get<std::string> ("Weighted Gradient BF Name"),
00027 p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00028 Residual (p.get<std::string> ("Residual Name"),
00029 p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") )
00030 {
00031 this->addDependentField(stress);
00032 this->addDependentField(J);
00033 this->addDependentField(defgrad);
00034 this->addDependentField(wGradBF);
00035
00036 this->addEvaluatedField(Residual);
00037
00038 this->setName("ThermoMechanicalMomentumResidual"+PHX::TypeString<EvalT>::value);
00039
00040 std::vector<PHX::DataLayout::size_type> dims;
00041 wGradBF.fieldTag().dataLayout().dimensions(dims);
00042 numNodes = dims[1];
00043 numQPs = dims[2];
00044 numDims = dims[3];
00045 int worksetSize = dims[0];
00046
00047
00048 F_inv.resize(worksetSize, numQPs, numDims, numDims);
00049 F_invT.resize(worksetSize, numQPs, numDims, numDims);
00050 JF_invT.resize(worksetSize, numQPs, numDims, numDims);
00051 P.resize(worksetSize, numQPs, numDims, numDims);
00052
00053 Teuchos::RCP<ParamLib> paramLib = p.get< Teuchos::RCP<ParamLib> >("Parameter Library");
00054
00055 matModel = p.get<std::string>("Stress Name");
00056
00057 zGrav=0.0;
00058 new Sacado::ParameterRegistration<EvalT, SPL_Traits>("zGrav", this, paramLib);
00059
00060 }
00061
00062
00063 template<typename EvalT, typename Traits>
00064 void ThermoMechanicalMomentumResidual<EvalT, Traits>::
00065 postRegistrationSetup(typename Traits::SetupData d,
00066 PHX::FieldManager<Traits>& fm)
00067 {
00068 this->utils.setFieldData(stress,fm);
00069 this->utils.setFieldData(J,fm);
00070 this->utils.setFieldData(defgrad,fm);
00071 this->utils.setFieldData(wGradBF,fm);
00072
00073 this->utils.setFieldData(Residual,fm);
00074 }
00075
00076
00077 template<typename EvalT, typename Traits>
00078 void ThermoMechanicalMomentumResidual<EvalT, Traits>::
00079 evaluateFields(typename Traits::EvalData workset)
00080 {
00081 std::cout.precision(15);
00082 typedef Intrepid::FunctionSpaceTools FST;
00083 typedef Intrepid::RealSpaceTools<ScalarT> RST;
00084
00085 RST::inverse(F_inv, defgrad);
00086 RST::transpose(F_invT, F_inv);
00087 FST::scalarMultiplyDataData<ScalarT>(JF_invT, J, F_invT);
00088 FST::tensorMultiplyDataData<ScalarT>(P, stress, JF_invT);
00089 for (std::size_t cell=0; cell < workset.numCells; ++cell)
00090 {
00091 for (std::size_t node=0; node < numNodes; ++node)
00092 {
00093 for (std::size_t dim=0; dim<numDims; dim++) Residual(cell,node,dim)=0.0;
00094 for (std::size_t qp=0; qp < numQPs; ++qp)
00095 {
00096 for (std::size_t i=0; i<numDims; i++)
00097 {
00098 for (std::size_t j=0; j<numDims; j++)
00099 {
00100 Residual(cell,node,i) += P(cell, qp, i, j) * wGradBF(cell, node, qp, j);
00101 }
00102 }
00103 }
00104 }
00105 }
00106 }
00107
00108 template<typename EvalT,typename Traits>
00109 typename ThermoMechanicalMomentumResidual<EvalT,Traits>::ScalarT&
00110 ThermoMechanicalMomentumResidual<EvalT,Traits>::getValue(const std::string &n)
00111 {
00112 return zGrav;
00113 }
00114
00115
00116
00117 }
00118