Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #ifndef ADAPT_ELEMENTSIZE_HPP
00008 #define ADAPT_ELEMENTSIZE_HPP
00009
00010 #include "Phalanx_Evaluator_WithBaseImpl.hpp"
00011 #include "Phalanx_Evaluator_Derived.hpp"
00012 #include "Phalanx_MDField.hpp"
00013 #include "Phalanx_DataLayout.hpp"
00014 #include "Teuchos_ParameterList.hpp"
00015 #include "Albany_ProblemUtils.hpp"
00016
00017 namespace Adapt {
00021 template<typename EvalT, typename Traits>
00022 class ElementSizeFieldBase :
00023 public PHX::EvaluatorWithBaseImpl<Traits>,
00024 public PHX::EvaluatorDerived<EvalT, Traits>
00025 {
00026 public:
00027 typedef typename EvalT::ScalarT ScalarT;
00028 typedef typename EvalT::MeshScalarT MeshScalarT;
00029 ElementSizeFieldBase(Teuchos::ParameterList& p,
00030 const Teuchos::RCP<Albany::Layouts>& dl);
00031
00032 void postRegistrationSetup(typename Traits::SetupData d,
00033 PHX::FieldManager<Traits>& vm);
00034
00035
00036 void preEvaluate(typename Traits::PreEvalData d) = 0;
00037 void postEvaluate(typename Traits::PostEvalData d) = 0;
00038 void evaluateFields(typename Traits::EvalData d) = 0;
00039
00040 Teuchos::RCP<const PHX::FieldTag> getEvaluatedFieldTag() const {
00041 return size_field_tag;
00042 }
00043
00044 Teuchos::RCP<const PHX::FieldTag> getResponseFieldTag() const {
00045 return size_field_tag;
00046 }
00047
00048 protected:
00049
00050 typedef enum {NOTSCALED, SCALAR, VECTOR} ScalingType;
00051
00052 Teuchos::RCP<const Teuchos::ParameterList> getValidSizeFieldParameters() const;
00053
00054 void getCellRadius(const std::size_t cell, MeshScalarT& cellRadius) const;
00055
00056 std::string scalingName;
00057 std::string className;
00058
00059 std::size_t numQPs;
00060 std::size_t numDims;
00061 std::size_t numVertices;
00062
00063 PHX::MDField<MeshScalarT,Cell,QuadPoint> qp_weights;
00064 PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim> coordVec;
00065 PHX::MDField<MeshScalarT,Cell,Node,Dim> coordVec_vertices;
00066
00067 bool outputToExodus;
00068 bool outputCellAverage;
00069 bool outputQPData;
00070 bool outputNodeData;
00071 bool isAnisotropic;
00072 ScalingType scalingType;
00073
00074 Teuchos::RCP< PHX::Tag<ScalarT> > size_field_tag;
00075 Albany::StateManager* pStateMgr;
00076
00077 };
00078
00079 template<typename EvalT, typename Traits>
00080 class ElementSizeField
00081 : public ElementSizeFieldBase<EvalT, Traits> {
00082 public:
00083 ElementSizeField(Teuchos::ParameterList& p,
00084 const Teuchos::RCP<Albany::Layouts>& dl) :
00085 ElementSizeFieldBase<EvalT, Traits>(p, dl){}
00086 void preEvaluate(typename Traits::PreEvalData d){}
00087 void postEvaluate(typename Traits::PostEvalData d){}
00088 void evaluateFields(typename Traits::EvalData d){}
00089 };
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100 template<typename Traits>
00101 class ElementSizeField<PHAL::AlbanyTraits::Residual,Traits>
00102 : public ElementSizeFieldBase<PHAL::AlbanyTraits::Residual, Traits> {
00103 public:
00104 ElementSizeField(Teuchos::ParameterList& p,
00105 const Teuchos::RCP<Albany::Layouts>& dl);
00106 void preEvaluate(typename Traits::PreEvalData d);
00107 void postEvaluate(typename Traits::PostEvalData d);
00108 void evaluateFields(typename Traits::EvalData d);
00109 };
00110
00111 }
00112
00113 #endif // Adapt_ElementSizeField.hpp