43 #ifndef PANZER_EQUATION_SET_DEFAULTIMPL_DECL_HPP 
   44 #define PANZER_EQUATION_SET_DEFAULTIMPL_DECL_HPP 
   57   template<
typename Traits> 
class FieldManager;
 
   62   template <
typename EvalT>
 
   69                             const int& default_integration_order,
 
   72                             const bool build_transient_support);
 
  107                                                         const std::string & model_name,
 
  114                                                             const std::string& model_name,
 
  121     virtual const std::vector<std::pair<std::string,Teuchos::RCP<panzer::PureBasis> > > & 
getProvidedDOFs() 
const;
 
  123     virtual const std::vector<std::vector<std::string> > & 
getCoordinateDOFs() 
const;
 
  125     virtual const std::map<int,Teuchos::RCP<panzer::IntegrationRule> > & 
getIntegrationRules() 
const;
 
  131     virtual std::string 
getType() 
const;
 
  144     void getAddedDOFs(std::vector<std::string> & dofNames) 
const;
 
  160     void updateDOF(
const std::string & dofName,
 
  162                    int integrationOrder = -1);
 
  211     void addDOF(
const std::string & dofName,
 
  212                 const std::string & basisType,
 
  213                 const int & basisOrder,
 
  214                 const int integrationOrder = -1,
 
  215                 const std::string residualName = 
"",
 
  216                 const std::string scatterName = 
"");
 
  227                     const std::string & gradName = 
"");
 
  238                     const std::string & curlName = 
"");
 
  248     void addDOFDiv(
const std::string & dofName,
 
  249                    const std::string & divName = 
"");
 
  260                               const std::string & dotName = 
"");
 
  309                                                     const std::string dof_name,
 
  310                                                     const std::vector<std::string>& residual_contributions,
 
  311                                                     const std::string residual_field_name = 
"") 
const;
 
  329                                                     const std::string dof_name,
 
  330                                                     const std::vector<std::string>& residual_contributions,
 
  331                                                     const std::vector<double>& scale_contributions,
 
  332                                                     const std::string residual_field_name = 
"") 
const;
 
  348         , 
grad(std::make_pair(false,
""))
 
  349         , 
curl(std::make_pair(false,
""))
 
  350         , 
div(std::make_pair(false,
""))
 
  361       std::pair<bool,std::string> 
grad;
 
  362       std::pair<bool,std::string> 
curl;
 
  363       std::pair<bool,std::string> 
div;
 
  366       void print(std::ostream & os)
 const {
 
  367         os << 
"DOF Desc = \"" << 
dofName << 
"\": " 
  370            << 
"Grad = (" << 
grad.first << 
", \"" << 
grad.second << 
"\"), " 
  371            << 
"Curl = (" << 
curl.first << 
", \"" << 
curl.second << 
"\"), " 
  372            << 
"Div = (" << 
div.first << 
", \"" << 
div.second << 
"\"), " 
  410     std::vector<std::pair<std::string,Teuchos::RCP<panzer::PureBasis> > >  
m_provided_dofs;
 
void buildAndRegisterResidualSummationEvaluator(PHX::FieldManager< panzer::Traits > &fm, const std::string dof_name, const std::vector< std::string > &residual_contributions, const std::string residual_field_name="") const 
 
Teuchos::RCP< panzer::IntegrationRule > getIntRuleForDOF(const std::string &dof_name) const 
Returns the integration rule associated with the residual contributions for the dof_name. 
 
std::map< std::string, DOFDescriptor >::const_iterator DescriptorIterator
For convenience, declare the DOFDescriptor iterator. 
 
void setupDeprecatedDOFsSupport()
 
std::pair< bool, std::string > residualName
 
virtual std::string getElementBlockId() const 
 
void addDOFGrad(const std::string &dofName, const std::string &gradName="")
 
virtual void setTangentParamNames(const std::vector< std::string > &tangent_param_names)
Set the list of tangent parameter names. 
 
virtual std::string getType() const 
Returns the type of the equation set object. Corresponds to the keyword used by the equation set fact...
 
int getIntegrationOrder(const std::string &dofName) const 
Get the integration order for an existing degree of freedom. 
 
Default implementation for accessing the GlobalData object. 
 
virtual const std::map< int, Teuchos::RCP< panzer::IntegrationRule > > & getIntegrationRules() const 
Return a map of unique integration rules for the equation set, key is the integration order...
 
const bool m_build_transient_support
 
std::map< std::string, Teuchos::RCP< panzer::PureBasis > > m_unique_bases
Key is the basis name from panzer::PureBasis::name() and value is the corresponding PureBasis...
 
virtual void buildAndRegisterInitialConditionEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const std::string &model_name, const Teuchos::ParameterList &models, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const 
 
std::vector< std::string > m_closure_model_ids
 
Teuchos::RCP< Teuchos::ParameterList > m_eval_plist
 
virtual void buildAndRegisterDOFProjectionsToIPEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLayoutLibrary &fl, const Teuchos::RCP< panzer::IntegrationRule > &ir, const Teuchos::Ptr< const panzer::LinearObjFactory< panzer::Traits > > &lof, const Teuchos::ParameterList &user_data) const 
 
int getBasisOrder(const std::string &dofName) const 
Get the basis order for an existing degree of freedom. 
 
std::map< std::string, std::pair< Teuchos::RCP< panzer::PureBasis >, Teuchos::RCP< std::vector< std::string > > > > m_basis_to_dofs
Map that links a common basis to a vector of dof names. Key is the unique basis name, the value is a pair that contains an RCP to a basis and an RCP to a vector of dof names that share the basis. 
 
std::map< int, Teuchos::RCP< panzer::IntegrationRule > > m_int_rules
Key is the integration rule order and the value is the corresponding integration rule. 
 
Teuchos::RCP< panzer::PureBasis > getBasisForDOF(const std::string &dof_name) const 
Returns the PureBasis associated with the residual contributions for the dof_name. 
 
void setDefaultValidParameters(Teuchos::ParameterList &valid_parameters)
 
Teuchos::RCP< Teuchos::ParameterList > getEquationSetParameterList() const 
Returns the parameter list used to build this equation set. 
 
Teuchos::RCP< panzer::PureBasis > basis
 
void setCoordinateDOFs(const std::vector< std::string > &dofNames)
 
Teuchos::RCP< panzer::IntegrationRule > intRule
 
Data for determining cell topology and dimensionality. 
 
void addDOFDiv(const std::string &dofName, const std::string &divName="")
 
virtual const Teuchos::RCP< Teuchos::ParameterList > getEvaluatorParameterList() const 
Returns the parameter list that will be passed off from the equaiton set to the closure model evaluat...
 
virtual void setupDOFs()
Builds the integration rule, basis, DOFs, and default parameter list. This MUST be called in the cons...
 
std::map< std::string, std::pair< Teuchos::RCP< panzer::PureBasis >, Teuchos::RCP< std::vector< std::string > > > >::const_iterator BasisIterator
For convenience, declare a basis iterator. 
 
virtual ~EquationSet_DefaultImpl()
 
std::pair< bool, std::string > div
 
std::map< std::string, DOFDescriptor > m_provided_dofs_desc
Maps the dof name into a DOFDescriptor. Should be private, but is protected so that the aux equaiton ...
 
std::pair< bool, std::string > curl
 
void addDOF(const std::string &dofName, const std::string &basisType, const int &basisOrder, const int integrationOrder=-1, const std::string residualName="", const std::string scatterName="")
 
void updateDOF(const std::string &dofName, int basisOrder, int integrationOrder=-1)
Modifying an existing DOF's basis function and integration rule. 
 
virtual void buildAndRegisterClosureModelEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLayoutLibrary &fl, const Teuchos::RCP< panzer::IntegrationRule > &ir, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const Teuchos::ParameterList &models, const Teuchos::ParameterList &user_data) const 
Register closure model evaluators with the model name internally specified by the equation set...
 
std::vector< std::vector< std::string > > m_coordinate_dofs
 
std::vector< std::string > m_tangent_param_names
 
Teuchos::RCP< panzer::BasisIRLayout > getBasisIRLayoutForDOF(const std::string &dof_name) const 
Returns the BasisIRLayout for the dof_name. 
 
virtual void buildAndRegisterGatherAndOrientationEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const 
 
std::vector< std::pair< std::string, Teuchos::RCP< panzer::PureBasis > > > m_provided_coord_prefixes
Key is the coordinate prefix name and the value is the corresponding basis. 
 
const panzer::CellData m_cell_data
 
void print(std::ostream &os) const 
 
int m_default_integration_order
 
const Teuchos::RCP< Teuchos::ParameterList > m_input_params
 
void getAddedDOFs(std::vector< std::string > &dofNames) const 
 
virtual const std::vector< std::pair< std::string, Teuchos::RCP< panzer::PureBasis > > > & getProvidedDOFs() const 
Return the Basis for the equation set, key is the DOF name (note coordinate DOFs are NOT included) ...
 
virtual const std::vector< std::vector< std::string > > & getCoordinateDOFs() const 
Return a vector of vectors that correspond to DOFs set as coordinate fields. 
 
void addDOFTimeDerivative(const std::string &dofName, const std::string &dotName="")
 
bool buildTransientSupport() const 
Returns true if transient support should be enabled in the equation set. 
 
virtual void buildAndRegisterScatterEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const 
 
void addDOFCurl(const std::string &dofName, const std::string &curlName="")
 
EquationSet_DefaultImpl(const Teuchos::RCP< Teuchos::ParameterList > ¶ms, const int &default_integration_order, const panzer::CellData &cell_data, const Teuchos::RCP< panzer::GlobalData > &global_data, const bool build_transient_support)
 
std::pair< bool, std::string > timeDerivative
 
virtual void setElementBlockId(const std::string &blockId)
 
std::vector< std::pair< std::string, Teuchos::RCP< panzer::PureBasis > > > m_provided_dofs
Key is the dof name and the value is the corresponding basis. 
 
virtual void buildAndRegisterEquationSetEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::FieldLibrary &fl, const Teuchos::ParameterList &user_data) const =0
 
std::pair< bool, std::string > grad
 
void addClosureModel(const std::string &closure_model)