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

Albany_PNPProblem.cpp

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 #include "Albany_PNPProblem.hpp"
00008 
00009 #include "Intrepid_FieldContainer.hpp"
00010 #include "Intrepid_DefaultCubatureFactory.hpp"
00011 #include "Shards_CellTopology.hpp"
00012 #include "PHAL_FactoryTraits.hpp"
00013 #include "Albany_Utils.hpp"
00014 #include "Albany_ProblemUtils.hpp"
00015 
00016 Albany::PNPProblem::
00017 PNPProblem( const Teuchos::RCP<Teuchos::ParameterList>& params_,
00018              const Teuchos::RCP<ParamLib>& paramLib_,
00019              const int numDim_) :
00020   Albany::AbstractProblem(params_, paramLib_),
00021   numDim(numDim_)
00022 {
00023 
00024   // Compute number of equations
00025   // TODO: read in num species
00026   numSpecies = 2;
00027   int num_eq = numSpecies + 1;
00028   this->setNumEquations(num_eq);
00029 
00030   // Print out a summary of the problem
00031   *out << "PNP problem: with numSpecies = " << numSpecies << std::endl;
00032 }
00033 
00034 Albany::PNPProblem::
00035 ~PNPProblem()
00036 {
00037 }
00038 
00039 void
00040 Albany::PNPProblem::
00041 buildProblem(
00042   Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00043   Albany::StateManager& stateMgr)
00044 {
00045   using Teuchos::rcp;
00046   /* Construct All Phalanx Evaluators */
00047   int physSets = meshSpecs.size();
00048   std::cout << "PNP Problem Num MeshSpecs: " << physSets << std::endl;
00049   fm.resize(physSets);
00050 
00051   for (int ps=0; ps<physSets; ps++) {
00052     fm[ps]  = Teuchos::rcp(new PHX::FieldManager<PHAL::AlbanyTraits>);
00053     buildEvaluators(*fm[ps], *meshSpecs[ps], stateMgr, BUILD_RESID_FM,
00054                     Teuchos::null);
00055   }
00056 
00057   if(meshSpecs[0]->nsNames.size() > 0) // Build a nodeset evaluator if nodesets are present
00058 
00059     constructDirichletEvaluators(*meshSpecs[0]);
00060 
00061   if(meshSpecs[0]->ssNames.size() > 0) // Build a sideset evaluator if sidesets are present
00062 
00063     constructNeumannEvaluators(meshSpecs[0]);
00064 }
00065 
00066 Teuchos::Array< Teuchos::RCP<const PHX::FieldTag> >
00067 Albany::PNPProblem::
00068 buildEvaluators(
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   // Call constructeEvaluators<EvalT>(*rfm[0], *meshSpecs[0], stateMgr);
00076   // for each EvalT in PHAL::AlbanyTraits::BEvalTypes
00077   Albany::ConstructEvaluatorsOp<PNPProblem> op(
00078     *this, fm0, meshSpecs, stateMgr, fmchoice, responseList);
00079   boost::mpl::for_each<PHAL::AlbanyTraits::BEvalTypes>(op);
00080   return *op.tags;
00081 }
00082 
00083 void
00084 Albany::PNPProblem::constructDirichletEvaluators(
00085         const Albany::MeshSpecsStruct& meshSpecs)
00086 {
00087    // Construct Dirichlet evaluators for all nodesets and names
00088    std::vector<std::string> dirichletNames(neq);
00089    int index = 0;
00090    dirichletNames[index++] = "C1";
00091    // TODO: arbitraty numSpecies
00092    if (numSpecies>=2) dirichletNames[index++] = "C2";
00093    if (numSpecies==3) dirichletNames[index++] = "C3";
00094    dirichletNames[index++] = "Phi";
00095    Albany::BCUtils<Albany::DirichletTraits> dirUtils;
00096    dfm = dirUtils.constructBCEvaluators(meshSpecs.nsNames, dirichletNames,
00097                                           this->params, this->paramLib);
00098 }
00099 
00100 //Neumann BCs
00101 void
00102 Albany::PNPProblem::constructNeumannEvaluators(const Teuchos::RCP<Albany::MeshSpecsStruct>& meshSpecs)
00103 {
00104 
00105    // Note: we only enter this function if sidesets are defined in the mesh file
00106    // i.e. meshSpecs.ssNames.size() > 0
00107    
00108     Albany::BCUtils<Albany::NeumannTraits> nbcUtils;
00109 
00110    // Check to make sure that Neumann BCs are given in the input file
00111 
00112    if(!nbcUtils.haveBCSpecified(this->params)) {
00113       return;
00114    }
00115 
00116    // Construct BC evaluators for all side sets and names
00117    // Note that the string index sets up the equation offset, so ordering is important
00118    //
00119    // Currently we aren't exactly doing this right.  I think to do this
00120    // correctly we need different neumann evaluators for each DOF (velocity,
00121    // pressure, temperature, flux) since velocity is a vector and the 
00122    // others are scalars.  The dof_names stuff is only used
00123    // for robin conditions, so at this point, as long as we don't enable
00124    // robin conditions, this should work.
00125    
00126    std::vector<std::string> nbcNames;
00127    Teuchos::RCP< Teuchos::Array<std::string> > dof_names =
00128      Teuchos::rcp(new Teuchos::Array<std::string>);
00129    // TODO: arbitraty numSpecies
00130    Teuchos::Array<Teuchos::Array<int> > offsets;
00131    int idx = 0;
00132      nbcNames.push_back("C1");
00133      offsets.push_back(Teuchos::Array<int>(1,idx++));
00134      if (numSpecies>=2) {
00135        nbcNames.push_back("C2");
00136        offsets.push_back(Teuchos::Array<int>(1,idx++));
00137      }
00138      if (numSpecies==3) {
00139        nbcNames.push_back("C3");
00140        offsets.push_back(Teuchos::Array<int>(1,idx++));
00141      }
00142      nbcNames.push_back("Phi");
00143      offsets.push_back(Teuchos::Array<int>(1,idx++));
00144      dof_names->push_back("Concentration");
00145      dof_names->push_back("Potential");
00146 
00147    // Construct BC evaluators for all possible names of conditions
00148    // Should only specify flux vector components (dudx, dudy, dudz), or dudn, not both
00149    std::vector<std::string> condNames; //dudx, dudy, dudz, dudn, basal 
00150 
00151 
00152    nfm[0] = nbcUtils.constructBCEvaluators(meshSpecs, nbcNames,
00153                                            Teuchos::arcp(dof_names),
00154                                            true, 0, condNames, offsets, dl,
00155                                            this->params, this->paramLib);
00156 }
00157 
00158 
00159 
00160 Teuchos::RCP<const Teuchos::ParameterList>
00161 Albany::PNPProblem::getValidProblemParameters() const
00162 {
00163   Teuchos::RCP<Teuchos::ParameterList> validPL =
00164     this->getGenericProblemParams("ValidPNPParams");
00165 
00166   validPL->sublist("Body Force", false, "");
00167 
00168   return validPL;
00169 }
00170 

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