00001
00002
00003
00004
00005
00006
00007 #ifndef TLPOROPLASTICITYRESIDMASS_HPP
00008 #define TLPOROPLASTICITYRESIDMASS_HPP
00009
00010 #include "Phalanx_ConfigDefs.hpp"
00011 #include "Phalanx_Evaluator_WithBaseImpl.hpp"
00012 #include "Phalanx_Evaluator_Derived.hpp"
00013 #include "Phalanx_MDField.hpp"
00014
00015 #include "Intrepid_CellTools.hpp"
00016 #include "Intrepid_Cubature.hpp"
00017
00018 namespace LCM {
00025 template<typename EvalT, typename Traits>
00026 class TLPoroPlasticityResidMass : public PHX::EvaluatorWithBaseImpl<Traits>,
00027 public PHX::EvaluatorDerived<EvalT, Traits> {
00028
00029 public:
00030
00031 TLPoroPlasticityResidMass(Teuchos::ParameterList& p);
00032
00033 void postRegistrationSetup(typename Traits::SetupData d,
00034 PHX::FieldManager<Traits>& vm);
00035
00036 void evaluateFields(typename Traits::EvalData d);
00037
00038 private:
00039
00040 typedef typename EvalT::ScalarT ScalarT;
00041 typedef typename EvalT::MeshScalarT MeshScalarT;
00042
00043
00044 PHX::MDField<MeshScalarT,Cell,Node,QuadPoint> wBF;
00045 PHX::MDField<ScalarT,Cell,QuadPoint> porePressure;
00046 PHX::MDField<ScalarT,Cell,QuadPoint> Tdot;
00047
00048 PHX::MDField<ScalarT,Cell,QuadPoint> ThermalCond;
00049 PHX::MDField<ScalarT,Cell,QuadPoint> kcPermeability;
00050 PHX::MDField<ScalarT,Cell,QuadPoint> porosity;
00051 PHX::MDField<ScalarT,Cell,QuadPoint> biotCoefficient;
00052 PHX::MDField<ScalarT,Cell,QuadPoint> biotModulus;
00053 PHX::MDField<MeshScalarT,Cell,Node,QuadPoint,Dim> wGradBF;
00054 PHX::MDField<ScalarT,Cell,QuadPoint,Dim> TGrad;
00055 PHX::MDField<ScalarT,Cell,QuadPoint> Source;
00056 Teuchos::Array<double> convectionVels;
00057 PHX::MDField<ScalarT,Cell,QuadPoint> rhoCp;
00058 PHX::MDField<ScalarT,Cell,QuadPoint> Absorption;
00059 PHX::MDField<ScalarT,Cell,QuadPoint,Dim,Dim> strain;
00060
00061 PHX::MDField<ScalarT,Cell,QuadPoint,Dim,Dim> defgrad;
00062 PHX::MDField<ScalarT,Cell,QuadPoint> J;
00063 PHX::MDField<ScalarT,Cell,QuadPoint> elementLength;
00064
00065
00066 PHX::MDField<MeshScalarT,Cell,Vertex,Dim> coordVec;
00067 Teuchos::RCP<Intrepid::Cubature<RealType> > cubature;
00068 Teuchos::RCP<shards::CellTopology> cellType;
00069 PHX::MDField<MeshScalarT,Cell,QuadPoint> weights;
00070
00071
00072 PHX::MDField<ScalarT,Dummy> deltaTime;
00073
00074
00075 std::string strainName, porePressureName, porosityName,
00076 JName;
00077
00078
00079
00080 bool haveSource;
00081 bool haveConvection;
00082 bool haveAbsorption;
00083 bool enableTransient;
00084 bool haverhoCp;
00085 bool haveMechanics;
00086 unsigned int numNodes;
00087 unsigned int numQPs;
00088 unsigned int numDims;
00089 unsigned int worksetSize;
00090
00091
00092 Intrepid::FieldContainer<ScalarT> flux;
00093 Intrepid::FieldContainer<ScalarT> fluxdt;
00094 Intrepid::FieldContainer<ScalarT> pterm;
00095 Intrepid::FieldContainer<ScalarT> tpterm;
00096 Intrepid::FieldContainer<ScalarT> aterm;
00097
00098 Intrepid::FieldContainer<RealType> refPoints;
00099 Intrepid::FieldContainer<RealType> refWeights;
00100 Intrepid::FieldContainer<MeshScalarT> jacobian;
00101 Intrepid::FieldContainer<MeshScalarT> jacobian_inv;
00102 Intrepid::FieldContainer<MeshScalarT> Gc;
00103
00104
00105 Intrepid::FieldContainer<ScalarT> F_inv;
00106 Intrepid::FieldContainer<ScalarT> F_invT;
00107 Intrepid::FieldContainer<ScalarT> C;
00108 Intrepid::FieldContainer<ScalarT> Cinv;
00109 Intrepid::FieldContainer<ScalarT> JF_invT;
00110 Intrepid::FieldContainer<ScalarT> KJF_invT;
00111 Intrepid::FieldContainer<ScalarT> Kref;
00112
00113 ScalarT porePbar, vol;
00114 ScalarT trialPbar;
00115
00116
00117
00118
00119
00120 PHX::MDField<ScalarT,Cell,Node> TResidual;
00121
00122 RealType stab_param_;
00123
00124 };
00125 }
00126
00127 #endif