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

SurfaceBasis.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 #ifndef SURFACE_BASIS_HPP
00008 #define SURFACE_BASIS_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 "Albany_Layouts.hpp"
00016 
00017 #include "Intrepid_CellTools.hpp"
00018 #include "Intrepid_Cubature.hpp"
00019 
00020 namespace LCM {
00021   
00028   template<typename EvalT, typename Traits>
00029   class SurfaceBasis : public PHX::EvaluatorWithBaseImpl<Traits>,
00030                        public PHX::EvaluatorDerived<EvalT, Traits>  {
00031 
00032   public:
00033     typedef typename EvalT::ScalarT ScalarT;
00034     typedef typename EvalT::MeshScalarT MeshScalarT;
00035     typedef Intrepid::FieldContainer<ScalarT> SFC;
00036     typedef Intrepid::FieldContainer<MeshScalarT> MFC;
00037 
00043     SurfaceBasis(const Teuchos::ParameterList& p,
00044                  const Teuchos::RCP<Albany::Layouts>& dl);
00045 
00049     void postRegistrationSetup(typename Traits::SetupData d,
00050                                PHX::FieldManager<Traits>& vm);
00051 
00055     void evaluateFields(typename Traits::EvalData d);
00056 
00062     void computeReferenceMidplaneCoords(const PHX::MDField<MeshScalarT,Cell,Vertex,Dim> refCoords,
00063                                         MFC & midplaneCoords);
00064 
00070     void computeCurrentMidplaneCoords(const PHX::MDField<ScalarT,Cell,Vertex,Dim> currentCoords,
00071                                       SFC & midplaneCoords);
00072 
00078     void computeReferenceBaseVectors(const MFC & midplaneCoords, 
00079                                      PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> basis);
00080 
00086     void computeCurrentBaseVectors(const SFC & midplaneCoords, 
00087                             PHX::MDField<ScalarT,Cell,QuadPoint,Dim,Dim> basis);
00088 
00096     void computeDualBaseVectors(const MFC & midplaneCoords, 
00097                                 const PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> basis, 
00098                                 PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim> normal, 
00099                                 PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> dualBasis);
00100 
00107     void computeJacobian(const PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> basis,
00108                          const PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> dualBasis,
00109                          PHX::MDField<MeshScalarT,Cell,QuadPoint> area);
00110 
00111   private:
00112     unsigned int  numDims, numNodes, numQPs, numPlaneNodes, numPlaneDims;
00113 
00114     bool needCurrentBasis;
00115 
00119     PHX::MDField<MeshScalarT,Cell,Vertex,Dim> referenceCoords;
00120 
00124     Teuchos::RCP<Intrepid::Cubature<RealType> > cubature;
00125     
00129     Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > intrepidBasis;
00130 
00134     Intrepid::FieldContainer<MeshScalarT> refMidplaneCoords;
00135 
00139     Intrepid::FieldContainer<ScalarT> currentMidplaneCoords;
00140 
00144     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> refBasis;
00145 
00149     PHX::MDField<MeshScalarT,Cell,QuadPoint> refArea;
00150 
00154     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim,Dim> refDualBasis;
00155 
00159     PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim> refNormal;
00160 
00161     // if we need to compute the current bases (for mechanics)
00165     PHX::MDField<ScalarT,Cell,Vertex,Dim> currentCoords;
00166 
00170     PHX::MDField<ScalarT,Cell,QuadPoint,Dim,Dim> currentBasis;
00171 
00175     Intrepid::FieldContainer<RealType> refValues;
00176 
00180     Intrepid::FieldContainer<RealType> refGrads;
00181 
00185     Intrepid::FieldContainer<RealType> refPoints;
00186 
00190     Intrepid::FieldContainer<RealType> refWeights;
00191   };
00192 }
00193 
00194 #endif

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