• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

TLPoroPlasticityResidMomentum_Def.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
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 namespace LCM {
00014 
00015 //**********************************************************************
00016 template<typename EvalT, typename Traits>
00017 TLPoroPlasticityResidMomentum<EvalT, Traits>::
00018 TLPoroPlasticityResidMomentum(const Teuchos::ParameterList& p) :
00019   TotalStress      (p.get<std::string>                   ("Total Stress Name"),
00020          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00021   J           (p.get<std::string>                   ("DetDefGrad Name"),
00022          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout") ),
00023   defgrad     (p.get<std::string>                   ("DefGrad Name"),
00024          p.get<Teuchos::RCP<PHX::DataLayout> >("QP Tensor Data Layout") ),
00025   wGradBF     (p.get<std::string>                   ("Weighted Gradient BF Name"),
00026          p.get<Teuchos::RCP<PHX::DataLayout> >("Node QP Vector Data Layout") ),
00027   ExResidual   (p.get<std::string>                   ("Residual Name"),
00028          p.get<Teuchos::RCP<PHX::DataLayout> >("Node Vector Data Layout") )
00029 {
00030   this->addDependentField(TotalStress);
00031   this->addDependentField(wGradBF);
00032   this->addDependentField(J);
00033   this->addDependentField(defgrad);
00034   this->addEvaluatedField(ExResidual);
00035 
00036   if (p.isType<bool>("Disable Transient"))
00037     enableTransient = !p.get<bool>("Disable Transient");
00038   else enableTransient = true;
00039 
00040   if (enableTransient) {
00041     // Two more fields are required for transient capability
00042     Teuchos::RCP<PHX::DataLayout> node_qp_scalar_dl =
00043        p.get< Teuchos::RCP<PHX::DataLayout> >("Node QP Scalar Data Layout");
00044     Teuchos::RCP<PHX::DataLayout> vector_dl =
00045        p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00046 
00047     PHX::MDField<MeshScalarT,Cell,Node,QuadPoint> wBF_tmp
00048         (p.get<std::string>("Weighted BF Name"), node_qp_scalar_dl);
00049     wBF = wBF_tmp;
00050     PHX::MDField<ScalarT,Cell,QuadPoint,Dim> uDotDot_tmp
00051         (p.get<std::string>("Time Dependent Variable Name"), vector_dl);
00052     uDotDot = uDotDot_tmp;
00053 
00054    this->addDependentField(wBF);
00055    this->addDependentField(uDotDot);
00056 
00057   }
00058 
00059 
00060   this->setName("TLPoroPlasticityResidMomentum"+PHX::TypeString<EvalT>::value);
00061 
00062   std::vector<PHX::DataLayout::size_type> dims;
00063   wGradBF.fieldTag().dataLayout().dimensions(dims);
00064   numNodes = dims[1];
00065   numQPs   = dims[2];
00066   numDims  = dims[3];
00067   int worksetSize = dims[0];
00068 
00069   // Works space FCs
00070   F_inv.resize(worksetSize, numQPs, numDims, numDims);
00071   F_invT.resize(worksetSize, numQPs, numDims, numDims);
00072   JF_invT.resize(worksetSize, numQPs, numDims, numDims);
00073 
00074 }
00075 
00076 //**********************************************************************
00077 template<typename EvalT, typename Traits>
00078 void TLPoroPlasticityResidMomentum<EvalT, Traits>::
00079 postRegistrationSetup(typename Traits::SetupData d,
00080                       PHX::FieldManager<Traits>& fm)
00081 {
00082   this->utils.setFieldData(TotalStress,fm);
00083   this->utils.setFieldData(wGradBF,fm);
00084   this->utils.setFieldData(J,fm);
00085   this->utils.setFieldData(defgrad,fm);
00086   this->utils.setFieldData(ExResidual,fm);
00087 
00088   if (enableTransient) this->utils.setFieldData(uDotDot,fm);
00089   if (enableTransient) this->utils.setFieldData(wBF,fm);
00090 
00091 }
00092 
00093 //**********************************************************************
00094 template<typename EvalT, typename Traits>
00095 void TLPoroPlasticityResidMomentum<EvalT, Traits>::
00096 evaluateFields(typename Traits::EvalData workset)
00097 {
00098   typedef Intrepid::FunctionSpaceTools FST;
00099   typedef Intrepid::RealSpaceTools<ScalarT> RST;
00100 
00101  // RST::inverse(F_inv, defgrad);
00102  // RST::transpose(F_invT, F_inv);
00103  // FST::scalarMultiplyDataData<ScalarT>(JF_invT, J, F_invT);
00104  // FST::tensorMultiplyDataData<ScalarT>(P, TotalStress, JF_invT);
00105 
00106     for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00107       for (std::size_t node=0; node < numNodes; ++node) {
00108               for (std::size_t dim=0; dim<numDims; dim++)  ExResidual(cell,node,dim)=0.0;
00109           for (std::size_t qp=0; qp < numQPs; ++qp) {
00110             for (std::size_t i=0; i<numDims; i++) {
00111               for (std::size_t dim=0; dim<numDims; dim++) {
00112                 ExResidual(cell,node,i) += TotalStress(cell, qp, i, dim) * wGradBF(cell, node, qp, dim);
00113     } } } } }
00114 
00115 
00116   if (workset.transientTerms && enableTransient)
00117     for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00118       for (std::size_t node=0; node < numNodes; ++node) {
00119           for (std::size_t qp=0; qp < numQPs; ++qp) {
00120             for (std::size_t i=0; i<numDims; i++) {
00121                 ExResidual(cell,node,i) += uDotDot(cell, qp, i) * wBF(cell, node, qp);
00122     } } } }
00123 
00124 
00125 //  FST::integrate<ScalarT>(ExResidual, TotalStress, wGradBF, Intrepid::COMP_CPP, false); // "false" overwrites
00126 
00127 }
00128 
00129 //**********************************************************************
00130 }
00131 

Generated on Wed Mar 26 2014 18:36:45 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1