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

SurfaceTLPoroMassResidual.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 
00008 #ifndef SURFACE_TL_PORO_MASS_RESIDUAL_HPP
00009 #define SURFACE_TL_PORO_MASS_RESIDUAL_HPP
00010 
00011 #include "Phalanx_ConfigDefs.hpp"
00012 #include "Phalanx_Evaluator_WithBaseImpl.hpp"
00013 #include "Phalanx_Evaluator_Derived.hpp"
00014 #include "Phalanx_MDField.hpp"
00015 #include "Intrepid_CellTools.hpp"
00016 #include "Intrepid_Cubature.hpp"
00017 
00018 #include "Albany_Layouts.hpp"
00019 
00020 namespace LCM {
00027 template<typename EvalT, typename Traits>
00028 class SurfaceTLPoroMassResidual : public PHX::EvaluatorWithBaseImpl<Traits>,
00029                               public PHX::EvaluatorDerived<EvalT, Traits>  {
00030 
00031 public:
00032 
00033   SurfaceTLPoroMassResidual(const Teuchos::ParameterList& p,
00034                         const Teuchos::RCP<Albany::Layouts>& dl);
00035 
00036   void postRegistrationSetup(typename Traits::SetupData d,
00037            PHX::FieldManager<Traits>& vm);
00038 
00039   void evaluateFields(typename Traits::EvalData d);
00040 
00041 private:
00042 
00043   typedef typename EvalT::ScalarT ScalarT;
00044   typedef typename EvalT::MeshScalarT MeshScalarT;
00045 
00046 
00047   // Input:
00049   ScalarT thickness;
00051   Teuchos::RCP<Intrepid::Cubature<RealType> > cubature;
00053   Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > intrepidBasis;
00055   PHX::MDField<ScalarT,Cell,QuadPoint,Dim> scalarGrad;
00057   PHX::MDField<MeshScalarT,Cell,Node,QuadPoint,Dim> surface_Grad_BF;
00059    PHX::MDField<ScalarT,Cell,QuadPoint> scalarJump;
00061   PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim, Dim> refDualBasis;
00063   PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim> refNormal;
00065   PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim> refArea;
00067   PHX::MDField<ScalarT,Cell,QuadPoint> J;
00069   PHX::MDField<ScalarT,Cell,QuadPoint> porePressure;
00071   PHX::MDField<ScalarT,Cell,Node> nodalPorePressure;
00073   PHX::MDField<ScalarT,Cell,QuadPoint> biotCoefficient;
00075   PHX::MDField<ScalarT,Cell,QuadPoint> biotModulus;
00077   PHX::MDField<ScalarT,Cell,QuadPoint> kcPermeability;
00079   PHX::MDField<ScalarT,Cell,QuadPoint,Dim, Dim> defGrad;
00080 
00081 //  // weight times basis function value at integration point
00082 //  PHX::MDField<MeshScalarT,Cell,Node,QuadPoint> wBF;
00083 
00084   //Data from previous time step
00085    std::string porePressureName, JName;
00086 
00087    // Time
00088    PHX::MDField<ScalarT,Dummy> deltaTime;
00089 
00091   Intrepid::FieldContainer<RealType> refValues;
00092   Intrepid::FieldContainer<RealType> refGrads;
00093   Intrepid::FieldContainer<RealType> refPoints;
00094   Intrepid::FieldContainer<RealType> refWeights;
00095 
00096   // Work space FCs
00097   Intrepid::FieldContainer<ScalarT> F_inv;
00098   Intrepid::FieldContainer<ScalarT> F_invT;
00099   Intrepid::FieldContainer<ScalarT> C;
00100   Intrepid::FieldContainer<ScalarT> Cinv;
00101   Intrepid::FieldContainer<ScalarT> JF_invT;
00102   Intrepid::FieldContainer<ScalarT> KJF_invT;
00103   Intrepid::FieldContainer<ScalarT> Kref;
00104 
00105   // Temporary FieldContainers
00106   Intrepid::FieldContainer<ScalarT> flux;
00107 
00108   // Output:
00109   PHX::MDField<ScalarT,Cell,Node> poroMassResidual;
00110 
00111   unsigned int worksetSize;
00112   unsigned int numNodes;
00113   unsigned int numQPs;
00114   unsigned int numDims;
00115   unsigned int numPlaneNodes;
00116   unsigned int numPlaneDims;
00117 
00118   bool haveMech;
00119 };
00120 }
00121 
00122 #endif

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