45 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_FieldManager.hpp"
48 #include "Phalanx_Evaluator_Factory.hpp"
54 #include "Shards_CellTopology.hpp"
58 void panzer::buildPhysicsBlocks(
const std::map<std::string,std::string>& block_ids_to_physics_ids,
61 const int default_integration_order,
62 const std::size_t workset_size,
65 const bool build_transient_support,
67 const std::vector<std::string>& tangent_param_names)
77 map<string,string>::const_iterator itr;
78 for (itr = block_ids_to_physics_ids.begin(); itr!=block_ids_to_physics_ids.end();++itr) {
79 string element_block_id = itr->first;
80 string physics_block_id = itr->second;
82 map<string,RCP<const shards::CellTopology> >::const_iterator ct_itr =
83 block_ids_to_cell_topo.find(element_block_id);
86 "Falied to find CellTopology for element block id: \""
87 << element_block_id <<
"\"!");
88 RCP<const shards::CellTopology> cellTopo = ct_itr->second;
95 "Failed to find physics id: \""
97 <<
"\" requested by element block: \""
98 << element_block_id <<
"\"!");
100 RCP<panzer::PhysicsBlock> pb =
rcp(
new panzer::PhysicsBlock(Teuchos::sublist(physics_blocks_plist,physics_block_id,
true),
102 default_integration_order,
106 build_transient_support,
109 physicsBlocks.push_back(pb);
113 void panzer::readPhysicsBlocks(
const std::map<std::string,std::string>& block_ids_to_physics_ids,
125 map<string,string>::const_iterator itr;
126 for (itr = block_ids_to_physics_ids.begin(); itr!=block_ids_to_physics_ids.end();++itr) {
127 string element_block_id = itr->first;
128 string physics_block_id = itr->second;
133 "Failed to find physics id: \""
135 <<
"\" requested by element block: \""
136 << element_block_id <<
"\"!");
138 RCP<panzer::PhysicsBlock> pb =
rcp(
new panzer::PhysicsBlock(Teuchos::sublist(physics_blocks_plist,physics_block_id,
true),
140 physicsBlocks.push_back(pb);
147 bool throw_on_failure)
149 std::vector<Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator pb = physics_blocks.begin();
151 while (pb != physics_blocks.end()) {
152 if ((*pb)->elementBlockID() == element_block_id)
158 TEUCHOS_TEST_FOR_EXCEPTION(throw_on_failure,std::runtime_error,
"Error: panzer::findPhysicsBlock(): The requested physics block for element block\"" << element_block_id <<
"\" was not found in the vecotr of physics blocks!");
167 const std::string & element_block_id,
168 const int default_integration_order,
172 const bool build_transient_support,
173 const std::vector<std::string>& tangent_param_names) :
174 m_element_block_id(element_block_id),
175 m_default_integration_order(default_integration_order),
176 m_cell_data(cell_data),
177 m_input_parameters(physics_block_plist),
178 m_build_transient_support(build_transient_support),
179 m_global_data(global_data),
180 m_eqset_factory(factory)
192 build_transient_support,
193 tangent_param_names);
199 const std::string & element_block_id) :
200 m_element_block_id(element_block_id),
201 m_default_integration_order(1),
202 m_input_parameters(physics_block_plist),
203 m_build_transient_support(false)
210 const std::string & physics_block_id,
211 const int integration_order,
215 m_physics_id(physics_block_id),
216 m_element_block_id(element_block_id),
217 m_default_integration_order(integration_order),
218 m_cell_data(cell_data),
219 m_input_parameters(Teuchos::null),
220 m_build_transient_support(false),
221 m_global_data(global_data),
222 m_eqset_factory(Teuchos::null)
251 for(std::vector<StrPureBasisPair>::const_iterator itr=
m_provided_dofs.begin();
253 m_field_lib->addFieldAndBasis(itr->first,itr->second);
260 m_physics_id(pb.m_physics_id),
261 m_element_block_id(pb.m_element_block_id),
262 m_default_integration_order(pb.m_default_integration_order),
263 m_cell_data(cell_data),
264 m_input_parameters(pb.m_input_parameters),
265 m_build_transient_support(pb.m_build_transient_support),
266 m_global_data(pb.m_global_data),
267 m_eqset_factory(pb.m_eqset_factory)
278 const bool build_transient_support,
282 const std::vector<std::string>& tangent_param_names)
284 m_default_integration_order = default_integration_order;
285 m_build_transient_support = build_transient_support;
286 m_cell_data = cell_data;
287 m_global_data = global_data;
288 m_eqset_factory = eqset_factory;
290 initialize(m_input_parameters,
291 m_default_integration_order,
294 m_build_transient_support,
295 tangent_param_names);
300 const int& default_integration_order,
301 const std::string & element_block_id,
303 const bool build_transient_support,
304 const std::vector<std::string>& tangent_param_names)
310 "The physics block \"" << input_parameters->
name()
311 <<
"\" required by element block \"" << element_block_id
312 <<
"\" does not have any equation sets associated with it."
313 <<
" Please add at least one equation set to this physics block!");
315 m_equation_sets.clear();
318 typedef ParameterList::ConstIterator pl_iter;
319 for (pl_iter eq = input_parameters->
begin(); eq != input_parameters->
end(); ++eq) {
322 "All entries in the physics block \"" << m_physics_id
323 <<
"\" must be an equation set sublist!" );
328 = m_eqset_factory->buildEquationSet(eq_set_pl, default_integration_order, cell_data, m_global_data, build_transient_support);
331 for (
auto eq_it = eq_set->begin(); eq_it != eq_set->end(); ++eq_it) {
332 eq_it->setTangentParamNames(tangent_param_names);
336 m_equation_sets.push_back(eq_set);
339 const std::vector<StrPureBasisPair> & sbNames = eq_set->begin()->getProvidedDOFs();
340 for(std::size_t j=0;j<sbNames.size();j++) {
343 m_dof_names.push_back(sbNames[j].first);
346 m_provided_dofs.push_back(sbNames[j]);
349 m_bases[sbNames[j].second->name()] = sbNames[j].second;
352 for (std::size_t k=0; k<tangent_param_names.size(); ++k)
353 m_tangent_fields.push_back(
StrPureBasisPair( sbNames[j].first +
" SENSITIVITY " + tangent_param_names[k],
354 sbNames[j].second ) );
359 const std::vector<std::vector<std::string> > & coord_dofs = eq_set->begin()->getCoordinateDOFs();
360 m_coordinate_dofs.insert(m_coordinate_dofs.begin(),coord_dofs.begin(),coord_dofs.end());
369 const std::map<int,Teuchos::RCP<panzer::IntegrationRule> > & ir_map = eq_set->begin()->getIntegrationRules();
371 ir != ir_map.end(); ++ir)
372 m_integration_rules[ir->second->order()] = ir->second;
378 for(std::vector<StrPureBasisPair>::const_iterator itr=m_provided_dofs.begin();
379 itr!=m_provided_dofs.end();++itr)
380 m_field_lib->addFieldAndBasis(itr->first,itr->second);
383 for(std::size_t eq_i=0;eq_i<m_equation_sets.size();eq_i++) {
386 itr!=eq_set->end();++itr) {
387 itr->setElementBlockId(element_block_id);
399 using namespace panzer;
400 using namespace Teuchos;
403 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
404 eq_set = m_equation_sets.begin();
405 for (;eq_set != m_equation_sets.end(); ++eq_set) {
411 for (; eval_type != eqstm.end(); ++eval_type) {
417 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
418 eval_type->buildAndRegisterEquationSetEvaluators(fm, *m_field_lib, user_data);
419 eval_type->setDetailsIndex(di);
431 using namespace panzer;
432 using namespace Teuchos;
435 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
436 eq_set = m_equation_sets.begin();
437 for (;eq_set != m_equation_sets.end(); ++eq_set) {
443 for (; eval_type != eqstm.end(); ++eval_type) {
444 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
445 eval_type->buildAndRegisterGatherAndOrientationEvaluators(fm, *m_field_lib, lof, user_data);
446 eval_type->setDetailsIndex(di);
458 using namespace panzer;
459 using namespace Teuchos;
462 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
463 eq_set = m_equation_sets.begin();
464 for (;eq_set != m_equation_sets.end(); ++eq_set) {
470 for (; eval_type != eqstm.end(); ++eval_type) {
474 ir_iter != m_integration_rules.end(); ++ ir_iter) {
478 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
479 eval_type->buildAndRegisterDOFProjectionsToIPEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, lof, user_data);
480 eval_type->setDetailsIndex(di);
493 using namespace panzer;
494 using namespace Teuchos;
497 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
498 eq_set = m_equation_sets.begin();
499 for (;eq_set != m_equation_sets.end(); ++eq_set) {
505 for (; eval_type != eqstm.end(); ++eval_type) {
506 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
507 eval_type->buildAndRegisterScatterEvaluators(fm, *m_field_lib, lof, user_data);
508 eval_type->setDetailsIndex(di);
521 using namespace panzer;
522 using namespace Teuchos;
525 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
526 eq_set = m_equation_sets.begin();
527 for (;eq_set != m_equation_sets.end(); ++eq_set) {
533 for (; eval_type != eqstm.end(); ++eval_type) {
537 ir_iter != m_integration_rules.end(); ++ ir_iter) {
541 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
542 eval_type->buildAndRegisterClosureModelEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, factory, models, user_data);
543 eval_type->setDetailsIndex(di);
554 const std::string& model_name,
559 using namespace panzer;
560 using namespace Teuchos;
563 vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
564 eq_set = m_equation_sets.begin();
565 for (;eq_set != m_equation_sets.end(); ++eq_set) {
571 for (; eval_type != eqstm.end(); ++eval_type) {
575 ir_iter != m_integration_rules.end(); ++ ir_iter) {
578 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
579 eval_type->buildAndRegisterClosureModelEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, factory, model_name, models, user_data);
580 eval_type->setDetailsIndex(di);
586 if(m_equation_sets.size()==0) {
588 for (;eval_type != factory.end(); ++eval_type) {
594 ir_iter != m_integration_rules.end(); ++ ir_iter) {
600 eval_type->buildClosureModels(model_name,
602 *m_field_lib->buildFieldLayoutLibrary(*ir),
610 const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
611 eval_type->registerEvaluators(*evaluators,fm);
612 eval_type->setDetailsIndex(di);
622 const std::string& model_name,
628 using namespace panzer;
629 using namespace Teuchos;
632 this->buildAndRegisterInitialConditionEvaluatorsForType<panzer::Traits::Residual>(fm, factory, model_name, models, lof, user_data);
646 return m_provided_dofs;
652 return m_coordinate_dofs;
658 return m_tangent_fields;
667 const std::map<int,Teuchos::RCP<panzer::IntegrationRule> >& int_rules = this->getIntegrationRules();
669 ir_itr != int_rules.end(); ++ir_itr)
670 needs.
int_rules.push_back(ir_itr->second);
672 const std::map<std::string,Teuchos::RCP<panzer::PureBasis> >& bases= this->getBases();
673 const std::vector<StrPureBasisPair>& fieldToBasis = getProvidedDOFs();
675 b_itr != bases.end(); ++b_itr) {
677 needs.
bases.push_back(b_itr->second);
680 for(std::size_t d=0;d<fieldToBasis.size();d++) {
681 if(fieldToBasis[d].second->name()==b_itr->second->name()) {
698 const std::map<std::string,Teuchos::RCP<panzer::PureBasis> >&
705 const std::map<int,Teuchos::RCP<panzer::IntegrationRule> >&
708 return m_integration_rules;
715 "Cannot return a basis since none exist in this physics block.");
716 return m_bases.begin()->second->getIntrepid2Basis()->getBaseCellTopology();
728 return m_element_block_id;
746 return m_global_data;
std::string name() const
A unique key that is the combination of the basis type and basis order.
Teuchos::RCP< PhysicsBlock > copyWithCellData(const panzer::CellData &cell_data) const
const std::string & name() const
std::vector< Teuchos::RCP< const PureBasis > > bases
const std::vector< std::string > & getDOFNames() const
ConstIterator end() const
std::vector< StrPureBasisPair > m_provided_dofs
Object that contains information on the physics and discretization of a block of elements with the SA...
const std::map< std::string, Teuchos::RCP< panzer::PureBasis > > & getBases() const
Returns the unique set of bases, key is the unique panzer::PureBasis::name() of the basis...
const panzer::CellData & cellData() const
void buildAndRegisterClosureModelEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const Teuchos::ParameterList &models, const Teuchos::ParameterList &user_data) const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void buildAndRegisterScatterEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
Ordinal numParams() const
std::vector< Teuchos::RCP< const IntegrationRule > > int_rules
std::string physicsBlockID() const
void buildAndRegisterInitialConditionEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &factory, const std::string &model_name, const Teuchos::ParameterList &models, const panzer::LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
std::vector< std::string > rep_field_name
void buildAndRegisterGatherAndOrientationEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
Teuchos::RCP< panzer::GlobalData > globalData() const
int m_default_integration_order
std::vector< Teuchos::RCP< panzer::EquationSet_TemplateManager< panzer::Traits > > > m_equation_sets
void buildAndRegisterDOFProjectionsToIPEvaluators(PHX::FieldManager< panzer::Traits > &fm, const Teuchos::Ptr< const panzer::LinearObjFactory< panzer::Traits > > &lof, const Teuchos::ParameterList &user_data) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Data for determining cell topology and dimensionality.
panzer::CellData m_cell_data
bool isSublist(const std::string &name) const
std::string elementBlockID() const
const std::vector< std::vector< std::string > > & getCoordinateDOFs() const
std::map< std::string, Teuchos::RCP< panzer::PureBasis > > m_bases
map of unique bases, key is the panzer::PureBasis::name() corresponding to its value ...
Teuchos::RCP< Teuchos::ParameterList > m_input_parameters
store the input parameter list for copy ctors
std::map< int, Teuchos::RCP< panzer::IntegrationRule > > m_integration_rules
map of unique integration rules, key is panzer::IntegrationRule::order() corresponding to its value ...
ConstIterator begin() const
const shards::CellTopology getBaseCellTopology() const
const std::vector< StrPureBasisPair > & getTangentFields() const
Returns list of tangent fields from DOFs and tangent param names.
WorksetNeeds getWorksetNeeds() const
bool nonnull(const boost::shared_ptr< T > &p)
const std::map< int, Teuchos::RCP< panzer::IntegrationRule > > & getIntegrationRules() const
Returns the unique set of point rules, key is the unique panzer::PointRule::name() ...
void initialize(const int default_integration_order, const bool build_transient_support, const panzer::CellData &cell_data, const Teuchos::RCP< const panzer::EquationSetFactory > &factory, const Teuchos::RCP< panzer::GlobalData > &global_data, const std::vector< std::string > &tangent_param_names=std::vector< std::string >())
std::pair< std::string, Teuchos::RCP< panzer::PureBasis > > StrPureBasisPair
bool m_build_transient_support
Teuchos::RCP< FieldLibrary > m_field_lib
std::string m_element_block_id
#define TEUCHOS_ASSERT(assertion_test)
std::vector< std::string > m_dof_names
void buildAndRegisterEquationSetEvaluators(PHX::FieldManager< panzer::Traits > &fm, const Teuchos::ParameterList &user_data) const
const std::vector< StrPureBasisPair > & getProvidedDOFs() const