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

Albany_GPAMProblem.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_GPAMPROBLEM_HPP
00008 #define ALBANY_GPAMPROBLEM_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 GPAMProblem : public AbstractProblem {
00026   public:
00027   
00029     GPAMProblem(const Teuchos::RCP<Teuchos::ParameterList>& params,
00030      const Teuchos::RCP<ParamLib>& paramLib,
00031      const int numDim_);
00032 
00034     ~GPAMProblem();
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     // 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     GPAMProblem(const GPAMProblem&);
00060     
00062     GPAMProblem& operator=(const GPAMProblem&);
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_GPAMResid.hpp"
00096 
00097 #include "PHAL_Source.hpp"
00098 
00099 template <typename EvalT>
00100 Teuchos::RCP<const PHX::FieldTag>
00101 Albany::GPAMProblem::constructEvaluators(
00102   PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00103   const Albany::MeshSpecsStruct& meshSpecs,
00104   Albany::StateManager& stateMgr,
00105   Albany::FieldManagerChoice fieldManagerChoice,
00106   const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00107 {
00108   using Teuchos::RCP;
00109   using Teuchos::rcp;
00110   using Teuchos::ParameterList;
00111   using PHX::DataLayout;
00112   using PHX::MDALayout;
00113   using std::vector;
00114   using std::string;
00115   using std::map;
00116   using PHAL::AlbanyTraits;
00117   
00118   RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00119     intrepidBasis = Albany::getIntrepidBasis(meshSpecs.ctd);
00120   RCP<shards::CellTopology> cellType = rcp(new shards::CellTopology (&meshSpecs.ctd));
00121   
00122   const int numNodes = intrepidBasis->getCardinality();
00123   const int worksetSize = meshSpecs.worksetSize;
00124   
00125   Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00126   RCP <Intrepid::Cubature<RealType> > cubature = cubFactory.create(*cellType, meshSpecs.cubatureDegree);
00127   
00128   const int numQPts = cubature->getNumPoints();
00129   const int numVertices = cellType->getNodeCount();
00130   
00131   *out << "Field Dimensions: Workset=" << worksetSize 
00132        << ", Vertices= " << numVertices
00133        << ", Nodes= " << numNodes
00134        << ", QuadPts= " << numQPts
00135        << ", Dim= " << numDim << std::endl;
00136   
00137    int vecDim = neq;
00138 
00139    RCP<Albany::Layouts> dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim, vecDim));
00140    Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00141    bool supportsTransient=true;
00142    int offset=0;
00143 
00144    // Temporary variable used numerous times below
00145    Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00146 
00147    // Define Field Names
00148 
00149      Teuchos::ArrayRCP<string> dof_names(1);
00150      Teuchos::ArrayRCP<string> dof_names_dot(1);
00151      Teuchos::ArrayRCP<string> resid_names(1);
00152      dof_names[0] = "Concentration";
00153      dof_names_dot[0] = dof_names[0]+"_dot";
00154      resid_names[0] = "GPAM Residual";
00155      fm0.template registerEvaluator<EvalT>
00156        (evalUtils.constructGatherSolutionEvaluator(true, dof_names, dof_names_dot, offset));
00157 
00158      fm0.template registerEvaluator<EvalT>
00159        (evalUtils.constructDOFVecInterpolationEvaluator(dof_names[0], offset));
00160 
00161      fm0.template registerEvaluator<EvalT>
00162        (evalUtils.constructDOFVecInterpolationEvaluator(dof_names_dot[0], offset));
00163 
00164 // Uncommented this
00165      fm0.template registerEvaluator<EvalT>
00166        (evalUtils.constructDOFVecGradInterpolationEvaluator(dof_names[0], offset));
00167 
00168      fm0.template registerEvaluator<EvalT>
00169        (evalUtils.constructScatterResidualEvaluator(true, resid_names,offset, "Scatter GPAM"));
00170      offset += numDim;
00171 
00172    fm0.template registerEvaluator<EvalT>
00173      (evalUtils.constructGatherCoordinateVectorEvaluator());
00174 
00175    fm0.template registerEvaluator<EvalT>
00176      (evalUtils.constructMapToPhysicalFrameEvaluator(cellType, cubature));
00177 
00178    fm0.template registerEvaluator<EvalT>
00179      (evalUtils.constructComputeBasisFunctionsEvaluator(cellType, intrepidBasis, cubature));
00180 
00181 /*
00182    { // Specialized DofVecGrad Interpolation for this problem
00183     
00184      RCP<ParameterList> p = rcp(new ParameterList("DOFVecGrad Interpolation "+dof_names[0]));
00185      // Input
00186      p->set<string>("Variable Name", dof_names[0]);
00187      
00188      p->set<string>("Gradient BF Name", "Grad BF");
00189      
00190      // Output (assumes same Name as input)
00191      p->set<string>("Gradient Variable Name", dof_names[0]+" Gradient");
00192      
00193      ev = rcp(new PHAL::DOFVecGradInterpolation<EvalT,AlbanyTraits>(*p,dl));
00194      fm0.template registerEvaluator<EvalT>(ev);
00195    }
00196 
00197 */
00198   { // GPAM Resid
00199     RCP<ParameterList> p = rcp(new ParameterList("GPAM Resid"));
00200 
00201     //Input
00202     p->set<string>("Weighted BF Name", "wBF");
00203     p->set<string>("Weighted Gradient BF Name", "wGrad BF");
00204     p->set<string>("QP Variable Name", "Concentration");
00205     p->set<string>("QP Time Derivative Variable Name", "Concentration_dot");
00206     p->set<string>("Gradient QP Variable Name", "Concentration Gradient");
00207     
00208     //Output
00209     p->set<string>("Residual Name", "GPAM Residual");
00210 
00211     ev = rcp(new PHAL::GPAMResid<EvalT,AlbanyTraits>(*p,dl));
00212     fm0.template registerEvaluator<EvalT>(ev);
00213   }
00214 
00215   if (fieldManagerChoice == Albany::BUILD_RESID_FM)  {
00216     PHX::Tag<typename EvalT::ScalarT> res_tag("Scatter GPAM", dl->dummy);
00217     fm0.requireField<EvalT>(res_tag);
00218   }
00219   else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
00220     Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00221     return respUtils.constructResponses(fm0, *responseList, Teuchos::null, stateMgr);
00222   }
00223 
00224   return Teuchos::null;
00225 }
00226 #endif // ALBANY_GPAMPROBLEM_HPP

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