00001
00002
00003
00004
00005
00006
00007 #ifndef PHAL_NEUMANN_HPP
00008 #define PHAL_NEUMANN_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 #include "Albany_ProblemUtils.hpp"
00019 #include "Sacado_ParameterAccessor.hpp"
00020 #include "PHAL_AlbanyTraits.hpp"
00021
00022 #include "QCAD_MaterialDatabase.hpp"
00023
00024
00025 namespace PHAL {
00026
00032 template<typename EvalT, typename Traits>
00033 class NeumannBase :
00034 public PHX::EvaluatorWithBaseImpl<Traits>,
00035 public PHX::EvaluatorDerived<EvalT, Traits>,
00036 public Sacado::ParameterAccessor<EvalT, SPL_Traits> {
00037
00038 public:
00039
00040 enum NEU_TYPE {COORD, NORMAL, INTJUMP, PRESS, ROBIN, BASAL, TRACTION, LATERAL};
00041 enum SIDE_TYPE {OTHER, LINE, TRI, QUAD};
00042
00043 typedef typename EvalT::ScalarT ScalarT;
00044 typedef typename EvalT::MeshScalarT MeshScalarT;
00045
00046 NeumannBase(const Teuchos::ParameterList& p);
00047
00048 void postRegistrationSetup(typename Traits::SetupData d,
00049 PHX::FieldManager<Traits>& vm);
00050
00051 void evaluateFields(typename Traits::EvalData d) = 0;
00052
00053 ScalarT& getValue(const std::string &n);
00054
00055 protected:
00056
00057 const Teuchos::RCP<Albany::Layouts>& dl;
00058 const Teuchos::RCP<Albany::MeshSpecsStruct>& meshSpecs;
00059
00060 int cellDims, sideDims, numQPs, numQPsSide, numNodes;
00061 Teuchos::Array<int> offset;
00062 int numDOFsSet;
00063
00064
00065 std::string betaName;
00066 double L;
00067 MeshScalarT betaXY;
00068 enum BETAXY_NAME {CONSTANT, EXPTRIG, ISMIP_HOM_TEST_C, ISMIP_HOM_TEST_D, CONFINEDSHELF, CIRCULARSHELF, DOMEUQ, SCALAR_FIELD, LATERAL_BACKPRESSURE};
00069 BETAXY_NAME beta_type;
00070
00071
00072
00073
00074
00075 void calc_dudn_const(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00076 const Intrepid::FieldContainer<MeshScalarT>& phys_side_cub_points,
00077 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00078 const shards::CellTopology & celltopo,
00079 const int cellDims,
00080 int local_side_id,
00081 ScalarT scale = 1.0);
00082
00083
00084 void calc_dudn_robin(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00085 const Intrepid::FieldContainer<MeshScalarT>& phys_side_cub_points,
00086 const Intrepid::FieldContainer<ScalarT>& dof_side,
00087 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00088 const shards::CellTopology & celltopo,
00089 const int cellDims,
00090 int local_side_id,
00091 ScalarT scale,
00092 const ScalarT* robin_param_values);
00093
00094
00095 void calc_gradu_dotn_const(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00096 const Intrepid::FieldContainer<MeshScalarT>& phys_side_cub_points,
00097 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00098 const shards::CellTopology & celltopo,
00099 const int cellDims,
00100 int local_side_id);
00101
00102
00103 void calc_traction_components(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00104 const Intrepid::FieldContainer<MeshScalarT>& phys_side_cub_points,
00105 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00106 const shards::CellTopology & celltopo,
00107 const int cellDims,
00108 int local_side_id);
00109
00110
00111 void calc_press(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00112 const Intrepid::FieldContainer<MeshScalarT>& phys_side_cub_points,
00113 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00114 const shards::CellTopology & celltopo,
00115 const int cellDims,
00116 int local_side_id);
00117
00118
00119 #ifdef ALBANY_FELIX
00120 void calc_dudn_basal(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00121 const Intrepid::FieldContainer<ScalarT>& basalFriction_side,
00122 const Intrepid::FieldContainer<ScalarT>& dof_side,
00123 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00124 const shards::CellTopology & celltopo,
00125 const int cellDims,
00126 int local_side_id);
00127
00128 void calc_dudn_lateral(Intrepid::FieldContainer<ScalarT> & qp_data_returned,
00129 const Intrepid::FieldContainer<ScalarT>& thickness_side,
00130 const Intrepid::FieldContainer<ScalarT>& elevation_side,
00131 const Intrepid::FieldContainer<ScalarT>& dof_side,
00132 const Intrepid::FieldContainer<MeshScalarT>& jacobian_side_refcell,
00133 const shards::CellTopology & celltopo,
00134 const int cellDims,
00135 int local_side_id);
00136 #endif
00137
00138 void evaluateNeumannContribution(typename Traits::EvalData d);
00139
00140
00142 PHX::MDField<MeshScalarT,Cell,Vertex,Dim> coordVec;
00143 PHX::MDField<ScalarT,Cell,Node> dof;
00144 PHX::MDField<ScalarT,Cell,Node,VecDim> dofVec;
00145 PHX::MDField<ScalarT,Cell,Node> beta_field;
00146 PHX::MDField<ScalarT,Cell,Node> thickness_field;
00147 PHX::MDField<RealType,Cell,Node> elevation_field;
00148 Teuchos::RCP<shards::CellTopology> cellType;
00149 Teuchos::RCP<shards::CellTopology> sideType;
00150 Teuchos::RCP<Intrepid::Cubature<RealType> > cubatureCell;
00151 Teuchos::RCP<Intrepid::Cubature<RealType> > cubatureSide;
00152
00153
00154 Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > intrepidBasis;
00155
00156
00157 Intrepid::FieldContainer<RealType> cubPointsSide;
00158 Intrepid::FieldContainer<RealType> refPointsSide;
00159 Intrepid::FieldContainer<RealType> cubWeightsSide;
00160 Intrepid::FieldContainer<MeshScalarT> physPointsSide;
00161 Intrepid::FieldContainer<MeshScalarT> jacobianSide;
00162 Intrepid::FieldContainer<MeshScalarT> jacobianSide_det;
00163
00164 Intrepid::FieldContainer<MeshScalarT> physPointsCell;
00165
00166 Intrepid::FieldContainer<MeshScalarT> weighted_measure;
00167 Intrepid::FieldContainer<RealType> basis_refPointsSide;
00168 Intrepid::FieldContainer<MeshScalarT> trans_basis_refPointsSide;
00169 Intrepid::FieldContainer<MeshScalarT> weighted_trans_basis_refPointsSide;
00170
00171 Intrepid::FieldContainer<ScalarT> dofCell;
00172 Intrepid::FieldContainer<ScalarT> dofSide;
00173
00174 Intrepid::FieldContainer<ScalarT> dofCellVec;
00175 Intrepid::FieldContainer<ScalarT> dofSideVec;
00176
00177 Intrepid::FieldContainer<ScalarT> data;
00178
00179
00180 Intrepid::FieldContainer<ScalarT> neumann;
00181
00182 std::string sideSetID;
00183 Teuchos::Array<RealType> inputValues;
00184 std::string inputConditions;
00185 std::string name;
00186
00187 NEU_TYPE bc_type;
00188 SIDE_TYPE side_type;
00189 ScalarT const_val;
00190 ScalarT robin_vals[5];
00191 std::vector<ScalarT> dudx;
00192
00193 std::vector<ScalarT> matScaling;
00194
00195 };
00196
00197 template<typename EvalT, typename Traits> class Neumann;
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209 template<typename Traits>
00210 class Neumann<PHAL::AlbanyTraits::Residual,Traits>
00211 : public NeumannBase<PHAL::AlbanyTraits::Residual, Traits> {
00212 public:
00213 Neumann(Teuchos::ParameterList& p);
00214 void evaluateFields(typename Traits::EvalData d);
00215 private:
00216 typedef typename PHAL::AlbanyTraits::Residual::ScalarT ScalarT;
00217 };
00218
00219
00220
00221
00222 template<typename Traits>
00223 class Neumann<PHAL::AlbanyTraits::Jacobian,Traits>
00224 : public NeumannBase<PHAL::AlbanyTraits::Jacobian, Traits> {
00225 public:
00226 Neumann(Teuchos::ParameterList& p);
00227 void evaluateFields(typename Traits::EvalData d);
00228 private:
00229 typedef typename PHAL::AlbanyTraits::Jacobian::ScalarT ScalarT;
00230 };
00231
00232
00233
00234
00235 template<typename Traits>
00236 class Neumann<PHAL::AlbanyTraits::Tangent,Traits>
00237 : public NeumannBase<PHAL::AlbanyTraits::Tangent, Traits> {
00238 public:
00239 Neumann(Teuchos::ParameterList& p);
00240 void evaluateFields(typename Traits::EvalData d);
00241 private:
00242 typedef typename PHAL::AlbanyTraits::Tangent::ScalarT ScalarT;
00243 };
00244
00245
00246
00247
00248 #ifdef ALBANY_SG_MP
00249 template<typename Traits>
00250 class Neumann<PHAL::AlbanyTraits::SGResidual,Traits>
00251 : public NeumannBase<PHAL::AlbanyTraits::SGResidual, Traits> {
00252 public:
00253 Neumann(Teuchos::ParameterList& p);
00254 void evaluateFields(typename Traits::EvalData d);
00255 private:
00256 typedef typename PHAL::AlbanyTraits::SGResidual::ScalarT ScalarT;
00257 };
00258
00259
00260
00261
00262 template<typename Traits>
00263 class Neumann<PHAL::AlbanyTraits::SGJacobian,Traits>
00264 : public NeumannBase<PHAL::AlbanyTraits::SGJacobian, Traits> {
00265 public:
00266 Neumann(Teuchos::ParameterList& p);
00267 void evaluateFields(typename Traits::EvalData d);
00268 private:
00269 typedef typename PHAL::AlbanyTraits::SGJacobian::ScalarT ScalarT;
00270 };
00271
00272
00273
00274
00275 template<typename Traits>
00276 class Neumann<PHAL::AlbanyTraits::SGTangent,Traits>
00277 : public NeumannBase<PHAL::AlbanyTraits::SGTangent, Traits> {
00278 public:
00279 Neumann(Teuchos::ParameterList& p);
00280 void evaluateFields(typename Traits::EvalData d);
00281 private:
00282 typedef typename PHAL::AlbanyTraits::SGTangent::ScalarT ScalarT;
00283 };
00284
00285
00286
00287
00288 template<typename Traits>
00289 class Neumann<PHAL::AlbanyTraits::MPResidual,Traits>
00290 : public NeumannBase<PHAL::AlbanyTraits::MPResidual, Traits> {
00291 public:
00292 Neumann(Teuchos::ParameterList& p);
00293 void evaluateFields(typename Traits::EvalData d);
00294 private:
00295 typedef typename PHAL::AlbanyTraits::MPResidual::ScalarT ScalarT;
00296 };
00297
00298
00299
00300
00301 template<typename Traits>
00302 class Neumann<PHAL::AlbanyTraits::MPJacobian,Traits>
00303 : public NeumannBase<PHAL::AlbanyTraits::MPJacobian, Traits> {
00304 public:
00305 Neumann(Teuchos::ParameterList& p);
00306 void evaluateFields(typename Traits::EvalData d);
00307 private:
00308 typedef typename PHAL::AlbanyTraits::MPJacobian::ScalarT ScalarT;
00309 };
00310
00311
00312
00313
00314 template<typename Traits>
00315 class Neumann<PHAL::AlbanyTraits::MPTangent,Traits>
00316 : public NeumannBase<PHAL::AlbanyTraits::MPTangent, Traits> {
00317 public:
00318 Neumann(Teuchos::ParameterList& p);
00319 void evaluateFields(typename Traits::EvalData d);
00320 private:
00321 typedef typename PHAL::AlbanyTraits::MPTangent::ScalarT ScalarT;
00322 };
00323 #endif //ALBANY_SG_MP
00324
00325
00326
00327
00328
00329
00330 template<typename EvalT, typename Traits>
00331 class NeumannAggregator
00332 : public PHX::EvaluatorWithBaseImpl<Traits>,
00333 public PHX::EvaluatorDerived<EvalT, Traits>
00334 {
00335 private:
00336
00337 typedef typename EvalT::ScalarT ScalarT;
00338
00339 public:
00340
00341 NeumannAggregator(const Teuchos::ParameterList& p);
00342
00343 void postRegistrationSetup(typename Traits::SetupData d,
00344 PHX::FieldManager<Traits>& vm) {};
00345
00346 void evaluateFields(typename Traits::EvalData d) {};
00347
00348 };
00349
00350 }
00351
00352 #endif