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

Albany_PNPProblem.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 Albany_PNPPROBLEM_HPP
00008 #define Albany_PNPPROBLEM_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 PNPProblem : public Albany::AbstractProblem {
00026   public:
00027   
00029     PNPProblem(const Teuchos::RCP<Teuchos::ParameterList>& params,
00030      const Teuchos::RCP<ParamLib>& paramLib,
00031      const int numDim_);
00032 
00034     ~PNPProblem();
00035 
00037     virtual void buildProblem(
00038       Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00039       Albany::StateManager& stateMgr);
00040 
00042     virtual int spatialDimension() const { return numDim; }
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> getValidProblemParameters() const;
00055 
00056   private:
00057 
00059     PNPProblem(const PNPProblem&);
00060     
00062     PNPProblem& operator=(const PNPProblem&);
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     void constructNeumannEvaluators(const Teuchos::RCP<Albany::MeshSpecsStruct>& meshSpecs);
00077 
00078 
00079   protected:
00080 
00081    int numDim;        
00082    int numSpecies;        
00083 
00084    Teuchos::RCP<Albany::Layouts> dl; 
00085     
00086   };
00087 
00088 }
00089 
00090 #include "Intrepid_FieldContainer.hpp"
00091 #include "Intrepid_DefaultCubatureFactory.hpp"
00092 #include "Shards_CellTopology.hpp"
00093 
00094 #include "Albany_Utils.hpp"
00095 #include "Albany_ProblemUtils.hpp"
00096 #include "Albany_EvaluatorUtils.hpp"
00097 #include "Albany_ResponseUtilities.hpp"
00098 
00099 #include "PNP_PotentialResid.hpp"
00100 #include "PNP_ConcentrationResid.hpp"
00101 #include "PHAL_Neumann.hpp"
00102 
00103 template <typename EvalT>
00104 Teuchos::RCP<const PHX::FieldTag>
00105 Albany::PNPProblem::constructEvaluators(
00106   PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00107   const Albany::MeshSpecsStruct& meshSpecs,
00108   Albany::StateManager& stateMgr,
00109   Albany::FieldManagerChoice fieldManagerChoice,
00110   const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00111 {
00112   using Teuchos::RCP;
00113   using Teuchos::rcp;
00114   using Teuchos::ParameterList;
00115   using PHX::DataLayout;
00116   using PHX::MDALayout;
00117   using std::vector;
00118   using std::string;
00119   using std::map;
00120   using PHAL::AlbanyTraits;
00121   
00122   RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00123     intrepidBasis = Albany::getIntrepidBasis(meshSpecs.ctd);
00124   RCP<shards::CellTopology> cellType = rcp(new shards::CellTopology (&meshSpecs.ctd));
00125   
00126   const int numNodes = intrepidBasis->getCardinality();
00127   const int worksetSize = meshSpecs.worksetSize;
00128   
00129   Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00130   RCP <Intrepid::Cubature<RealType> > cubature = cubFactory.create(*cellType, meshSpecs.cubatureDegree);
00131   
00132   const int numQPts = cubature->getNumPoints();
00133   const int numVertices = cellType->getNodeCount();
00134   
00135   *out << "Field Dimensions: Workset=" << worksetSize 
00136        << ", Vertices= " << numVertices
00137        << ", Nodes= " << numNodes
00138        << ", QuadPts= " << numQPts
00139        << ", Dim= " << numDim << std::endl;
00140   
00141 
00142    dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim, numSpecies));
00143    Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00144    bool supportsTransient=true;
00145    int offset=0;
00146 
00147    // Temporary variable used numerous times below
00148    Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00149 
00150    // Define Field Names
00151 
00152    
00153    {
00154      Teuchos::ArrayRCP<string> dof_names(1);
00155      Teuchos::ArrayRCP<string> dof_names_dot(1);
00156      Teuchos::ArrayRCP<string> resid_names(1);
00157      dof_names[0] = "Potential";
00158      resid_names[0] = "Potential Residual";
00159      fm0.template registerEvaluator<EvalT>
00160        (evalUtils.constructGatherSolutionEvaluator(false, dof_names, dof_names_dot, offset));
00161 
00162      fm0.template registerEvaluator<EvalT>
00163        (evalUtils.constructDOFInterpolationEvaluator(dof_names[0], offset));
00164 
00165      fm0.template registerEvaluator<EvalT>
00166        (evalUtils.constructDOFGradInterpolationEvaluator(dof_names[0], offset));
00167 
00168      fm0.template registerEvaluator<EvalT>
00169        (evalUtils.constructScatterResidualEvaluator(false, resid_names,offset, "Scatter Potential"));
00170      offset ++;
00171    }
00172 
00173    {
00174      Teuchos::ArrayRCP<string> dof_names(1);
00175      Teuchos::ArrayRCP<string> dof_names_dot(1);
00176      Teuchos::ArrayRCP<string> resid_names(1);
00177      dof_names[0] = "Concentration";
00178      dof_names_dot[0] = dof_names[0]+"_dot";
00179      resid_names[0] = "Concentration Residual";
00180      fm0.template registerEvaluator<EvalT>
00181        (evalUtils.constructGatherSolutionEvaluator(true, dof_names, dof_names_dot, offset));
00182 
00183      fm0.template registerEvaluator<EvalT>
00184        (evalUtils.constructDOFVecInterpolationEvaluator(dof_names[0], offset));
00185 
00186      fm0.template registerEvaluator<EvalT>
00187        (evalUtils.constructDOFVecInterpolationEvaluator(dof_names_dot[0], offset));
00188 
00189      fm0.template registerEvaluator<EvalT>
00190        (evalUtils.constructDOFVecGradInterpolationEvaluator(dof_names[0], offset));
00191 
00192      fm0.template registerEvaluator<EvalT>
00193        (evalUtils.constructScatterResidualEvaluator(true, resid_names,offset, "Scatter Concentration"));
00194      offset += numSpecies;
00195    }
00196 
00197    fm0.template registerEvaluator<EvalT>
00198      (evalUtils.constructGatherCoordinateVectorEvaluator());
00199 
00200    fm0.template registerEvaluator<EvalT>
00201      (evalUtils.constructMapToPhysicalFrameEvaluator(cellType, cubature));
00202 
00203    fm0.template registerEvaluator<EvalT>
00204      (evalUtils.constructComputeBasisFunctionsEvaluator(cellType, intrepidBasis, cubature));
00205 
00206   { // COncentration Resid
00207     RCP<ParameterList> p = rcp(new ParameterList("COncentration Resid"));
00208 
00209     //Input
00210     p->set<string>("Weighted BF Name", "wBF");
00211     p->set<string>("Weighted Gradient BF Name", "wGrad BF");
00212 
00213     // Variable names hardwired to Concentration and Potential in evaluators
00214 
00215     p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00216   
00217     //Output
00218     p->set<string>("Residual Name", "Concentration Residual");
00219 
00220     ev = rcp(new PNP::ConcentrationResid<EvalT,AlbanyTraits>(*p,dl));
00221     fm0.template registerEvaluator<EvalT>(ev);
00222   }
00223   
00224 
00225   { // Potential Resid
00226     RCP<ParameterList> p = rcp(new ParameterList("Potential Resid"));
00227 
00228     //Input
00229     p->set<string>("Weighted BF Name", "wBF");
00230     p->set<string>("Weighted Gradient BF Name", "wGrad BF");
00231     // Variable names hardwired to Concentration and Potential in evaluators
00232 
00233     //Output
00234     p->set<string>("Residual Name", "Potential Residual");
00235 
00236     ev = rcp(new PNP::PotentialResid<EvalT,AlbanyTraits>(*p,dl));
00237     fm0.template registerEvaluator<EvalT>(ev);
00238   }
00239 
00240 
00241   if (fieldManagerChoice == Albany::BUILD_RESID_FM)  {
00242     Teuchos::RCP<const PHX::FieldTag> ret_tag;
00243     PHX::Tag<typename EvalT::ScalarT> con_tag("Scatter Potential", dl->dummy);
00244     fm0.requireField<EvalT>(con_tag);
00245     PHX::Tag<typename EvalT::ScalarT> mom_tag("Scatter Concentration", dl->dummy);
00246     fm0.requireField<EvalT>(mom_tag);
00247     ret_tag = mom_tag.clone();
00248     return ret_tag;
00249   }
00250   else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
00251     Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00252     return respUtils.constructResponses(fm0, *responseList, stateMgr);
00253   }
00254 
00255   return Teuchos::null;
00256 }
00257 #endif // Albany_PNP_HPP

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