00001
00002
00003
00004
00005
00006
00007 #ifndef ALBANY_LINCOMPRNSPROBLEM_HPP
00008 #define ALBANY_LINCOMPRNSPROBLEM_HPP
00009
00010 #include "Teuchos_RCP.hpp"
00011 #include "Teuchos_ParameterList.hpp"
00012
00013 #include "Albany_AbstractProblem.hpp"
00014
00015 #include "Phalanx.hpp"
00016 #include "PHAL_Workset.hpp"
00017 #include "PHAL_Dimension.hpp"
00018
00019 namespace Albany {
00020
00025 class LinComprNSProblem : public AbstractProblem {
00026 public:
00027
00029 LinComprNSProblem(const Teuchos::RCP<Teuchos::ParameterList>& params,
00030 const Teuchos::RCP<ParamLib>& paramLib,
00031 const int numDim_);
00032
00034 ~LinComprNSProblem();
00035
00037 virtual int spatialDimension() const { return numDim; }
00038
00040 virtual void buildProblem(
00041 Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> > meshSpecs,
00042 StateManager& stateMgr);
00043
00044
00045 virtual Teuchos::Array< Teuchos::RCP<const PHX::FieldTag> >
00046 buildEvaluators(
00047 PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00048 const Albany::MeshSpecsStruct& meshSpecs,
00049 Albany::StateManager& stateMgr,
00050 Albany::FieldManagerChoice fmchoice,
00051 const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00052
00054 Teuchos::RCP<const Teuchos::ParameterList> getValidProblemParameters() const;
00055
00056 private:
00057
00059 LinComprNSProblem(const LinComprNSProblem&);
00060
00062 LinComprNSProblem& operator=(const LinComprNSProblem&);
00063
00064 public:
00065
00067 template <typename EvalT> Teuchos::RCP<const PHX::FieldTag>
00068 constructEvaluators(
00069 PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00070 const Albany::MeshSpecsStruct& meshSpecs,
00071 Albany::StateManager& stateMgr,
00072 Albany::FieldManagerChoice fmchoice,
00073 const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00074
00075 void constructDirichletEvaluators(const Albany::MeshSpecsStruct& meshSpecs);
00076
00077 protected:
00078 int numDim;
00079
00080 };
00081
00082 }
00083
00084 #include "Intrepid_FieldContainer.hpp"
00085 #include "Intrepid_DefaultCubatureFactory.hpp"
00086 #include "Shards_CellTopology.hpp"
00087
00088 #include "Albany_Utils.hpp"
00089 #include "Albany_ProblemUtils.hpp"
00090 #include "Albany_EvaluatorUtils.hpp"
00091 #include "Albany_ResponseUtilities.hpp"
00092
00093 #include "PHAL_DOFVecGradInterpolation.hpp"
00094
00095 #include "PHAL_LinComprNSResid.hpp"
00096
00097 #include "PHAL_Source.hpp"
00098 #include "PHAL_LinComprNSBodyForce.hpp"
00099
00100 template <typename EvalT>
00101 Teuchos::RCP<const PHX::FieldTag>
00102 Albany::LinComprNSProblem::constructEvaluators(
00103 PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00104 const Albany::MeshSpecsStruct& meshSpecs,
00105 Albany::StateManager& stateMgr,
00106 Albany::FieldManagerChoice fieldManagerChoice,
00107 const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00108 {
00109 using Teuchos::RCP;
00110 using Teuchos::rcp;
00111 using Teuchos::ParameterList;
00112 using PHX::DataLayout;
00113 using PHX::MDALayout;
00114 using std::vector;
00115 using std::string;
00116 using std::map;
00117 using PHAL::AlbanyTraits;
00118
00119 RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00120 intrepidBasis = Albany::getIntrepidBasis(meshSpecs.ctd);
00121 RCP<shards::CellTopology> cellType = rcp(new shards::CellTopology (&meshSpecs.ctd));
00122
00123 const int numNodes = intrepidBasis->getCardinality();
00124 const int worksetSize = meshSpecs.worksetSize;
00125
00126 Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00127 RCP <Intrepid::Cubature<RealType> > cubature = cubFactory.create(*cellType, meshSpecs.cubatureDegree);
00128
00129 const int numQPts = cubature->getNumPoints();
00130 const int numVertices = cellType->getNodeCount();
00131
00132 *out << "Field Dimensions: Workset=" << worksetSize
00133 << ", Vertices= " << numVertices
00134 << ", Nodes= " << numNodes
00135 << ", QuadPts= " << numQPts
00136 << ", Dim= " << numDim << std::endl;
00137
00138 int vecDim = neq;
00139
00140 RCP<Albany::Layouts> dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim, vecDim));
00141 Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00142 bool supportsTransient=true;
00143 int offset=0;
00144
00145
00146 Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00147
00148
00149
00150 Teuchos::ArrayRCP<string> dof_names(1);
00151 Teuchos::ArrayRCP<string> dof_names_dot(1);
00152 Teuchos::ArrayRCP<string> resid_names(1);
00153 dof_names[0] = "qFluct";
00154 dof_names_dot[0] = dof_names[0]+"_dot";
00155 resid_names[0] = "LinComprNS Residual";
00156 fm0.template registerEvaluator<EvalT>
00157 (evalUtils.constructGatherSolutionEvaluator(true, dof_names, dof_names_dot, offset));
00158
00159 fm0.template registerEvaluator<EvalT>
00160 (evalUtils.constructDOFVecInterpolationEvaluator(dof_names[0], offset));
00161
00162 fm0.template registerEvaluator<EvalT>
00163 (evalUtils.constructDOFVecInterpolationEvaluator(dof_names_dot[0], offset));
00164
00165
00166
00167
00168 fm0.template registerEvaluator<EvalT>
00169 (evalUtils.constructScatterResidualEvaluator(true, resid_names, offset, "Scatter LinComprNS"));
00170
00171 fm0.template registerEvaluator<EvalT>
00172 (evalUtils.constructGatherCoordinateVectorEvaluator());
00173
00174 fm0.template registerEvaluator<EvalT>
00175 (evalUtils.constructMapToPhysicalFrameEvaluator(cellType, cubature));
00176
00177 fm0.template registerEvaluator<EvalT>
00178 (evalUtils.constructComputeBasisFunctionsEvaluator(cellType, intrepidBasis, cubature));
00179
00180 {
00181
00182 RCP<ParameterList> p = rcp(new ParameterList("DOFVecGrad Interpolation "+dof_names[0]));
00183
00184 p->set<string>("Variable Name", dof_names[0]);
00185 p->set<string>("Gradient BF Name", "Grad BF");
00186 p->set<int>("Offset of First DOF", offset);
00187
00188
00189 p->set<string>("Gradient Variable Name", dof_names[0]+" Gradient");
00190
00191 ev = rcp(new PHAL::DOFVecGradInterpolation<EvalT,AlbanyTraits>(*p,dl));
00192 fm0.template registerEvaluator<EvalT>(ev);
00193 }
00194
00195 {
00196 RCP<ParameterList> p = rcp(new ParameterList("LinComprNS Resid"));
00197
00198
00199 p->set<string>("Weighted BF Name", "wBF");
00200 p->set<string>("Weighted Gradient BF Name", "wGrad BF");
00201 p->set<string>("QP Variable Name", "qFluct");
00202 p->set<string>("QP Time Derivative Variable Name", "qFluct_dot");
00203 p->set<string>("Gradient QP Variable Name", "qFluct Gradient");
00204 p->set<string>("Body Force Name", "Body Force");
00205
00206 p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00207 p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_vecgradient);
00208 p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
00209 p->set< RCP<DataLayout> >("Node QP Gradient Data Layout", dl->node_qp_gradient);
00210
00211 p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00212 Teuchos::ParameterList& paramList = params->sublist("Equation Set");
00213 p->set<Teuchos::ParameterList*>("Parameter List", ¶mList);
00214
00215
00216 p->set<string>("Residual Name", "LinComprNS Residual");
00217 p->set< RCP<DataLayout> >("Node Vector Data Layout", dl->node_vector);
00218
00219 ev = rcp(new PHAL::LinComprNSResid<EvalT,AlbanyTraits>(*p));
00220 fm0.template registerEvaluator<EvalT>(ev);
00221 }
00222
00223 {
00224 RCP<ParameterList> p = rcp(new ParameterList("Body Force"));
00225
00226
00227 p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00228 p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00229 p->set< RCP<DataLayout> >("QP Gradient Data Layout", dl->qp_gradient);
00230 p->set<string>("Coordinate Vector Name", "Coord Vec");
00231
00232 Teuchos::ParameterList& paramList = params->sublist("Body Force");
00233 p->set<Teuchos::ParameterList*>("Parameter List", ¶mList);
00234
00235
00236 p->set<string>("Body Force Name", "Body Force");
00237
00238 ev = rcp(new PHAL::LinComprNSBodyForce<EvalT,AlbanyTraits>(*p));
00239 fm0.template registerEvaluator<EvalT>(ev);
00240 }
00241
00242
00243 if (fieldManagerChoice == Albany::BUILD_RESID_FM) {
00244 PHX::Tag<typename EvalT::ScalarT> res_tag("Scatter LinComprNS", dl->dummy);
00245 fm0.requireField<EvalT>(res_tag);
00246 }
00247 else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
00248 Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00249 return respUtils.constructResponses(fm0, *responseList, Teuchos::null, stateMgr);
00250 }
00251
00252 return Teuchos::null;
00253 }
00254 #endif // ALBANY_LinComprNSPROBLEM_HPP