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

Aeras_XZScalarAdvectionProblem.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 AERAS_XZSCALARADVECTIONPROBLEM_HPP
00008 #define AERAS_XZSCALARADVECTIONPROBLEM_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 Aeras {
00020 
00025   class XZScalarAdvectionProblem : public Albany::AbstractProblem {
00026   public:
00027   
00029     XZScalarAdvectionProblem(const Teuchos::RCP<Teuchos::ParameterList>& params,
00030      const Teuchos::RCP<ParamLib>& paramLib,
00031      const int numDim_);
00032 
00034     ~XZScalarAdvectionProblem();
00035 
00037     virtual int spatialDimension() const { return numDim; }
00038 
00040     virtual void buildProblem(
00041       Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00042       Albany::StateManager& stateMgr);
00043 
00044     // Build evaluators
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>
00055     getValidProblemParameters() const;
00056 
00057   private:
00058 
00060     XZScalarAdvectionProblem(const XZScalarAdvectionProblem&);
00061     
00063     XZScalarAdvectionProblem& operator=(const XZScalarAdvectionProblem&);
00064 
00065   public:
00066 
00069     template <typename EvalT> Teuchos::RCP<const PHX::FieldTag>
00070     constructEvaluators(
00071       PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00072       const Albany::MeshSpecsStruct& meshSpecs,
00073       Albany::StateManager& stateMgr,
00074       Albany::FieldManagerChoice fmchoice,
00075       const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00076 
00077     void constructDirichletEvaluators(const Albany::MeshSpecsStruct& meshSpecs);
00078     void constructNeumannEvaluators(const Teuchos::RCP<Albany::MeshSpecsStruct>& meshSpecs);
00079 
00080   protected:
00081     int numDim;
00082     Teuchos::RCP<Albany::Layouts> dl;
00083 
00084   };
00085 
00086 }
00087 
00088 #include "Intrepid_FieldContainer.hpp"
00089 #include "Intrepid_DefaultCubatureFactory.hpp"
00090 #include "Shards_CellTopology.hpp"
00091 
00092 #include "Albany_Utils.hpp"
00093 #include "Albany_ProblemUtils.hpp"
00094 #include "Albany_EvaluatorUtils.hpp"
00095 #include "Albany_ResponseUtilities.hpp"
00096 #include "PHAL_Neumann.hpp"
00097 
00098 #include "Aeras_XZScalarAdvectionResid.hpp"
00099 
00100 template <typename EvalT>
00101 Teuchos::RCP<const PHX::FieldTag>
00102 Aeras::XZScalarAdvectionProblem::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   int vecDim = neq;
00132   
00133   *out << "Field Dimensions: Workset=" << worksetSize 
00134        << ", Vertices= " << numVertices
00135        << ", Nodes= " << numNodes
00136        << ", QuadPts= " << numQPts
00137        << ", Dim= " << numDim 
00138        << ", vecDim= " << vecDim << std::endl;
00139   
00140    dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim, vecDim));
00141    Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00142 
00143    // Temporary variable used numerous times below
00144    Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00145 
00146    // Define Field Names
00147 
00148   Teuchos::ArrayRCP<std::string> dof_names(1);
00149   Teuchos::ArrayRCP<std::string> dof_names_dot(1);
00150   Teuchos::ArrayRCP<std::string> resid_names(1);
00151   dof_names[0] = "rho";
00152   dof_names_dot[0] = dof_names[0]+"_dot";
00153   resid_names[0] = "XZScalarAdvection Residual";
00154 
00155   // Construct Standard FEM evaluators for Vector equation
00156   fm0.template registerEvaluator<EvalT>
00157     (evalUtils.constructGatherSolutionEvaluator(false, dof_names, dof_names_dot));
00158 
00159   fm0.template registerEvaluator<EvalT>
00160     (evalUtils.constructDOFInterpolationEvaluator(dof_names[0]));
00161 
00162   fm0.template registerEvaluator<EvalT>
00163     (evalUtils.constructDOFInterpolationEvaluator(dof_names_dot[0]));
00164 
00165   fm0.template registerEvaluator<EvalT>
00166     (evalUtils.constructDOFGradInterpolationEvaluator(dof_names[0]));
00167 
00168   fm0.template registerEvaluator<EvalT>
00169     (evalUtils.constructScatterResidualEvaluator(false, resid_names, 0, "Scatter XZScalarAdvection"));
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   { // XZScalarAdvection Resid
00181     RCP<ParameterList> p = rcp(new ParameterList("XZScalarAdvection Resid"));
00182    
00183     //Input
00184     p->set<std::string>("Weighted BF Name", "wBF");
00185     p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
00186     p->set<std::string>("QP Variable Name", "rho");
00187     p->set<std::string>("QP Time Derivative Variable Name", "rho_dot");
00188     p->set<std::string>("Gradient QP Variable Name", "rho Gradient");
00189     p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00190     
00191     p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00192 
00193     Teuchos::ParameterList& paramList = params->sublist("XZScalarAdvection Problem");
00194     p->set<Teuchos::ParameterList*>("XZScalarAdvection Problem", &paramList);
00195 
00196     //Output
00197     p->set<std::string>("Residual Name", "XZScalarAdvection Residual");
00198 
00199     ev = rcp(new Aeras::XZScalarAdvectionResid<EvalT,AlbanyTraits>(*p,dl));
00200     fm0.template registerEvaluator<EvalT>(ev);
00201   }
00202 
00203   if (fieldManagerChoice == Albany::BUILD_RESID_FM)  {
00204     PHX::Tag<typename EvalT::ScalarT> res_tag("Scatter XZScalarAdvection", dl->dummy);
00205     fm0.requireField<EvalT>(res_tag);
00206   }
00207   else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
00208     Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00209     return respUtils.constructResponses(fm0, *responseList, Teuchos::null, stateMgr);
00210   }
00211 
00212   return Teuchos::null;
00213 }
00214 #endif // AERAS_XZSCALARADVECTIONPROBLEM_HPP

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