Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_PhysicsBlock.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #include "Teuchos_RCP.hpp"
13 #include "Teuchos_Assert.hpp"
14 
15 #include "Phalanx_FieldManager.hpp"
16 #include "Phalanx_Evaluator_Factory.hpp"
17 #include "Panzer_Traits.hpp"
18 #include "Panzer_PhysicsBlock.hpp"
19 #include "Panzer_PureBasis.hpp"
22 #include "Shards_CellTopology.hpp"
23 
24 // *******************************************************************
25 
26 void panzer::buildPhysicsBlocks(const std::map<std::string,std::string>& block_ids_to_physics_ids,
27  const std::map<std::string,Teuchos::RCP<const shards::CellTopology> > & block_ids_to_cell_topo,
28  const Teuchos::RCP<Teuchos::ParameterList>& physics_blocks_plist,
29  const int default_integration_order,
30  const std::size_t workset_size,
32  const Teuchos::RCP<panzer::GlobalData>& global_data,
33  const bool build_transient_support,
34  std::vector<Teuchos::RCP<panzer::PhysicsBlock> > & physicsBlocks,
35  const std::vector<std::string>& tangent_param_names)
36 {
37  using Teuchos::RCP;
38  using Teuchos::rcp;
39  using std::map;
40  using std::string;
41 
42  TEUCHOS_ASSERT(nonnull(physics_blocks_plist));
43 
44  // Create a physics block for each element block
45  map<string,string>::const_iterator itr;
46  for (itr = block_ids_to_physics_ids.begin(); itr!=block_ids_to_physics_ids.end();++itr) {
47  string element_block_id = itr->first;
48  string physics_block_id = itr->second;
49 
50  map<string,RCP<const shards::CellTopology> >::const_iterator ct_itr =
51  block_ids_to_cell_topo.find(element_block_id);
52  TEUCHOS_TEST_FOR_EXCEPTION(ct_itr==block_ids_to_cell_topo.end(),
53  std::runtime_error,
54  "Falied to find CellTopology for element block id: \""
55  << element_block_id << "\"!");
56  RCP<const shards::CellTopology> cellTopo = ct_itr->second;
57 
58  const panzer::CellData volume_cell_data(workset_size,cellTopo);
59 
60  // find physics block parameter sublist
61  TEUCHOS_TEST_FOR_EXCEPTION(!physics_blocks_plist->isSublist(physics_block_id),
62  std::runtime_error,
63  "Failed to find physics id: \""
64  << physics_block_id
65  << "\" requested by element block: \""
66  << element_block_id << "\"!");
67 
68  RCP<panzer::PhysicsBlock> pb = rcp(new panzer::PhysicsBlock(Teuchos::sublist(physics_blocks_plist,physics_block_id,true),
69  element_block_id,
70  default_integration_order,
71  volume_cell_data,
72  eqset_factory,
73  global_data,
74  build_transient_support,
75  tangent_param_names
76  ));
77  physicsBlocks.push_back(pb);
78  }
79 }
80 
81 void panzer::readPhysicsBlocks(const std::map<std::string,std::string>& block_ids_to_physics_ids,
82  const Teuchos::RCP<Teuchos::ParameterList>& physics_blocks_plist,
83  std::vector<Teuchos::RCP<panzer::PhysicsBlock> > & physicsBlocks)
84 {
85  using Teuchos::RCP;
86  using Teuchos::rcp;
87  using std::map;
88  using std::string;
89 
90  TEUCHOS_ASSERT(nonnull(physics_blocks_plist));
91 
92  // Create a physics block for each element block
93  map<string,string>::const_iterator itr;
94  for (itr = block_ids_to_physics_ids.begin(); itr!=block_ids_to_physics_ids.end();++itr) {
95  string element_block_id = itr->first;
96  string physics_block_id = itr->second;
97 
98  // find physics block parameter sublist
99  TEUCHOS_TEST_FOR_EXCEPTION(!physics_blocks_plist->isSublist(physics_block_id),
100  std::runtime_error,
101  "Failed to find physics id: \""
102  << physics_block_id
103  << "\" requested by element block: \""
104  << element_block_id << "\"!");
105 
106  RCP<panzer::PhysicsBlock> pb = rcp(new panzer::PhysicsBlock(Teuchos::sublist(physics_blocks_plist,physics_block_id,true),
107  element_block_id));
108  physicsBlocks.push_back(pb);
109  }
110 }
111 
112 // *******************************************************************
113 Teuchos::RCP<panzer::PhysicsBlock> panzer::findPhysicsBlock(const std::string element_block_id,
114  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> > & physics_blocks,
115  bool throw_on_failure)
116 {
117  std::vector<Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator pb = physics_blocks.begin();
118 
119  while (pb != physics_blocks.end()) {
120  if ((*pb)->elementBlockID() == element_block_id)
121  return *pb;
122 
123  ++pb;
124  }
125 
126  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!");
127 
129  return null_pb;
130 }
131 
132 // *******************************************************************
135  const std::string & element_block_id,
136  const int default_integration_order,
137  const panzer::CellData & cell_data,
139  const Teuchos::RCP<panzer::GlobalData>& global_data,
140  const bool build_transient_support,
141  const std::vector<std::string>& tangent_param_names) :
142  m_element_block_id(element_block_id),
143  m_default_integration_order(default_integration_order),
144  m_cell_data(cell_data),
145  m_input_parameters(physics_block_plist),
146  m_build_transient_support(build_transient_support),
147  m_global_data(global_data),
148  m_eqset_factory(factory),
149  m_active_evaluation_types(Sacado::mpl::size<panzer::Traits::EvalTypes>::value,true)
150 {
151  TEUCHOS_ASSERT(nonnull(physics_block_plist));
152  TEUCHOS_ASSERT(nonnull(factory));
153  TEUCHOS_ASSERT(nonnull(global_data));
154 
155  m_physics_id = physics_block_plist->name();
156 
160  m_cell_data,
161  build_transient_support,
162  tangent_param_names);
163 }
164 
165 // *******************************************************************
168  const std::string & element_block_id) :
169  m_element_block_id(element_block_id),
170  m_default_integration_order(1),
171  m_input_parameters(physics_block_plist),
172  m_build_transient_support(false),
173  m_active_evaluation_types(Sacado::mpl::size<panzer::Traits::EvalTypes>::value,true)
174 {
175 }
176 
177 // *******************************************************************
179 PhysicsBlock(const std::string & element_block_id,
180  const std::string & physics_block_id,
181  const int integration_order,
182  const panzer::CellData & cell_data,
183  const Teuchos::RCP<panzer::GlobalData>& global_data,
184  const Teuchos::RCP<panzer::PureBasis> & basis) :
185  m_physics_id(physics_block_id),
186  m_element_block_id(element_block_id),
187  m_default_integration_order(integration_order),
188  m_cell_data(cell_data),
189  m_input_parameters(Teuchos::null),
190  m_build_transient_support(false),
191  m_global_data(global_data),
192  m_eqset_factory(Teuchos::null),
193  m_active_evaluation_types(Sacado::mpl::size<panzer::Traits::EvalTypes>::value,true)
194 {
195  using Teuchos::RCP;
197 
198  // there will be no equation sets
199  m_equation_sets.clear();
200 
201 
202  // Generate list of dof names
203  m_dof_names.push_back("FAKE");
204 
205  // Generate dof name (string) / basis pairs
206  m_provided_dofs.push_back(std::make_pair("FAKE",basis));
207 
208  // Generate unique list of bases
209  m_bases[basis->name()] = basis;
210 
211  // Get a unique list of point rules. NOTE: This assumes that the
212  // same point rules are used for all evaluation types. In the
213  // future we could easily change this by adding another level here
214  // to differentiate the point rules for each evaluation type.
215  // This would require some refactoring of the phsyics block
216  // registration routines and the workset builder to support all
217  // combinations of bases and point rules for each evaluation type.
218  m_integration_rules[integration_order] = Teuchos::rcp(new panzer::IntegrationRule(integration_order,cell_data));
219 
220  // build up field library
222  for(std::vector<StrPureBasisPair>::const_iterator itr=m_provided_dofs.begin();
223  itr!=m_provided_dofs.end();++itr)
224  m_field_lib->addFieldAndBasis(itr->first,itr->second);
225 }
226 
227 // *******************************************************************
230  const panzer::CellData & cell_data) :
231  m_physics_id(pb.m_physics_id),
232  m_element_block_id(pb.m_element_block_id),
233  m_default_integration_order(pb.m_default_integration_order),
234  m_cell_data(cell_data), // NOT copied from pb
235  m_input_parameters(pb.m_input_parameters),
236  m_build_transient_support(pb.m_build_transient_support),
237  m_global_data(pb.m_global_data),
238  m_eqset_factory(pb.m_eqset_factory),
239  m_active_evaluation_types(Sacado::mpl::size<panzer::Traits::EvalTypes>::value,true)
240 {
244  m_cell_data,
246 }
247 
248 // *******************************************************************
249 void panzer::PhysicsBlock::initialize(const int default_integration_order,
250  const bool build_transient_support,
251  const panzer::CellData & cell_data,
253  const Teuchos::RCP<panzer::GlobalData>& global_data,
254  const std::vector<std::string>& tangent_param_names)
255 {
256  m_default_integration_order = default_integration_order;
257  m_build_transient_support = build_transient_support;
258  m_cell_data = cell_data;
259  m_global_data = global_data;
260  m_eqset_factory = eqset_factory;
261 
262  initialize(m_input_parameters,
263  m_default_integration_order,
264  m_element_block_id,
265  m_cell_data,
266  m_build_transient_support,
267  tangent_param_names);
268 }
269 
270 // *******************************************************************
272  const int& default_integration_order,
273  const std::string & element_block_id,
274  const panzer::CellData & cell_data,
275  const bool build_transient_support,
276  const std::vector<std::string>& tangent_param_names)
277 {
278  using Teuchos::RCP;
280 
281  TEUCHOS_TEST_FOR_EXCEPTION(input_parameters->numParams() < 1, std::runtime_error,
282  "The physics block \"" << input_parameters->name()
283  << "\" required by element block \"" << element_block_id
284  << "\" does not have any equation sets associated with it."
285  << " Please add at least one equation set to this physics block!");
286 
287  m_equation_sets.clear();
288 
289  // Loop over equation sets
290  typedef ParameterList::ConstIterator pl_iter;
291  for (pl_iter eq = input_parameters->begin(); eq != input_parameters->end(); ++eq) {
292 
293  TEUCHOS_TEST_FOR_EXCEPTION( !(eq->second.isList()), std::logic_error,
294  "All entries in the physics block \"" << m_physics_id
295  << "\" must be an equation set sublist!" );
296 
297  RCP<ParameterList> eq_set_pl = Teuchos::sublist(input_parameters,eq->first,true);
298 
300  = m_eqset_factory->buildEquationSet(eq_set_pl, default_integration_order, cell_data, m_global_data, build_transient_support);
301 
302  // Set tangent parameter names for each equation set
303  for (auto eq_it = eq_set->begin(); eq_it != eq_set->end(); ++eq_it) {
304  eq_it->setTangentParamNames(tangent_param_names);
305  }
306 
307  // add this equation set in
308  m_equation_sets.push_back(eq_set);
309 
310  // Interrogate DOFs
311  const std::vector<StrPureBasisPair> & sbNames = eq_set->begin()->getProvidedDOFs();
312  for(std::size_t j=0;j<sbNames.size();j++) {
313 
314  // Generate list of dof names
315  m_dof_names.push_back(sbNames[j].first);
316 
317  // Generate dof name (string) / basis pairs
318  m_provided_dofs.push_back(sbNames[j]);
319 
320  // Generate unique list of bases
321  m_bases[sbNames[j].second->name()] = sbNames[j].second;
322 
323  // Generate tangent field names
324  for (std::size_t k=0; k<tangent_param_names.size(); ++k)
325  m_tangent_fields.push_back( StrPureBasisPair( sbNames[j].first + " SENSITIVITY " + tangent_param_names[k],
326  sbNames[j].second ) );
327 
328  }
329 
330  // add coordinate dofs to physics block
331  const std::vector<std::vector<std::string> > & coord_dofs = eq_set->begin()->getCoordinateDOFs();
332  m_coordinate_dofs.insert(m_coordinate_dofs.begin(),coord_dofs.begin(),coord_dofs.end());
333 
334  // Get a unique list of point rules. NOTE: This assumes that the
335  // same point rules are used for all evaluation types. In the
336  // future we could easily change this by adding another level here
337  // to differentiate the point rules for each evaluation type.
338  // This would require some refactoring of the phsyics block
339  // registration routines and the workset builder to support all
340  // combinations of bases and point rules for each evaluation type.
341  const std::map<int,Teuchos::RCP<panzer::IntegrationRule> > & ir_map = eq_set->begin()->getIntegrationRules();
342  for(std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir = ir_map.begin();
343  ir != ir_map.end(); ++ir)
344  m_integration_rules[ir->second->order()] = ir->second;
345 
346  }
347 
348  // build up field library
349  m_field_lib = Teuchos::rcp(new FieldLibrary);
350  for(std::vector<StrPureBasisPair>::const_iterator itr=m_provided_dofs.begin();
351  itr!=m_provided_dofs.end();++itr)
352  m_field_lib->addFieldAndBasis(itr->first,itr->second);
353 
354  // setup element blocks: loop over each evaluation type
355  for(std::size_t eq_i=0;eq_i<m_equation_sets.size();eq_i++) {
356  RCP<panzer::EquationSet_TemplateManager<panzer::Traits> > eq_set = m_equation_sets[eq_i];
358  itr!=eq_set->end();++itr) {
359  itr->setElementBlockId(element_block_id);
360  }
361  }
362 
363 }
364 
365 // *******************************************************************
366 void panzer::PhysicsBlock::setActiveEvaluationTypes(const std::vector<bool>& aet)
367 {
368  TEUCHOS_ASSERT(aet.size() == std::size_t(Sacado::mpl::size<panzer::Traits::EvalTypes>::value));
369  m_active_evaluation_types = aet;
370 }
371 
372 // *******************************************************************
374 {
375  for (auto&& t : m_active_evaluation_types)
376  t = true;
377 }
378 
379 // *******************************************************************
382  const Teuchos::ParameterList& user_data) const
383 {
384  using namespace std;
385  using namespace panzer;
386  using namespace Teuchos;
387 
388  // Loop over equation set template managers
389  vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
390  eq_set = m_equation_sets.begin();
391  for (;eq_set != m_equation_sets.end(); ++eq_set) {
392 
393  // Loop over evaluation types
396  eqstm.begin();
397  int idx = 0;
398  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
399  if (m_active_evaluation_types[idx]) {
400  // Do not loop over integration rules. Only call this for the
401  // ir that the residual is integrated over. Otherwise the
402  // residual gets contributions from multiple integrations of the
403  // same cell! This ir is only known by equaiton set.
404  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
405  eval_type->buildAndRegisterEquationSetEvaluators(fm, *m_field_lib, user_data);
406  eval_type->setDetailsIndex(di);
407  }
408  }
409  }
410 }
411 
412 // *******************************************************************
416  const Teuchos::ParameterList& user_data) const
417 {
418  using namespace std;
419  using namespace panzer;
420  using namespace Teuchos;
421 
422  // Loop over equation set template managers
423  vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
424  eq_set = m_equation_sets.begin();
425  for (;eq_set != m_equation_sets.end(); ++eq_set) {
426 
427  // Loop over evaluation types
430  eqstm.begin();
431  int idx = 0;
432  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
433  if (m_active_evaluation_types[idx]) {
434  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
435  eval_type->buildAndRegisterGatherAndOrientationEvaluators(fm, *m_field_lib, lof, user_data);
436  eval_type->setDetailsIndex(di);
437  }
438  }
439  }
440 }
441 
442 // *******************************************************************
446  const Teuchos::ParameterList& user_data) const
447 {
448  using namespace std;
449  using namespace panzer;
450  using namespace Teuchos;
451 
452  // Loop over equation set template managers
453  vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
454  eq_set = m_equation_sets.begin();
455  for (;eq_set != m_equation_sets.end(); ++eq_set) {
456 
457  // Loop over evaluation types
460  eqstm.begin();
461  int idx = 0;
462  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
463  if (m_active_evaluation_types[idx]) {
464 
465  // Loop over integration rules
466  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir_iter = m_integration_rules.begin();
467  ir_iter != m_integration_rules.end(); ++ ir_iter) {
468 
469  Teuchos::RCP<panzer::IntegrationRule> ir = ir_iter->second;
470 
471  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
472  eval_type->buildAndRegisterDOFProjectionsToIPEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, lof, user_data);
473  eval_type->setDetailsIndex(di);
474  }
475  }
476  }
477  }
478 }
479 
480 // *******************************************************************
484  const Teuchos::ParameterList& user_data) const
485 {
486  using namespace std;
487  using namespace panzer;
488  using namespace Teuchos;
489 
490  // Loop over equation set template managers
491  vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
492  eq_set = m_equation_sets.begin();
493  for (;eq_set != m_equation_sets.end(); ++eq_set) {
494 
495  // Loop over evaluation types
498  eqstm.begin();
499  int idx = 0;
500  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
501  if (m_active_evaluation_types[idx]) {
502  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
503  eval_type->buildAndRegisterScatterEvaluators(fm, *m_field_lib, lof, user_data);
504  eval_type->setDetailsIndex(di);
505  }
506  }
507  }
508 }
509 
510 // *******************************************************************
514  const Teuchos::ParameterList& models,
515  const Teuchos::ParameterList& user_data) const
516 {
517  using namespace std;
518  using namespace panzer;
519  using namespace Teuchos;
520 
521  // Loop over equation set template managers
522  vector< RCP<EquationSet_TemplateManager<panzer::Traits> > >::const_iterator
523  eq_set = m_equation_sets.begin();
524  for (;eq_set != m_equation_sets.end(); ++eq_set) {
525 
526  // Loop over evaluation types
529  eqstm.begin();
530  int idx = 0;
531  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
532  if (m_active_evaluation_types[idx]) {
533 
534  // Loop over integration rules
535  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir_iter = m_integration_rules.begin();
536  ir_iter != m_integration_rules.end(); ++ ir_iter) {
537 
538  Teuchos::RCP<panzer::IntegrationRule> ir = ir_iter->second;
539 
540  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
541  eval_type->buildAndRegisterClosureModelEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, factory, models, user_data);
542  eval_type->setDetailsIndex(di);
543  }
544  }
545  }
546  }
547 }
548 
549 // *******************************************************************
550 
554  const std::string& model_name,
555  const Teuchos::ParameterList& models,
556  const Teuchos::ParameterList& user_data) const
557 {
558  using namespace std;
559  using namespace panzer;
560  using namespace Teuchos;
561 
562  // Loop over equation set template managers
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) {
566 
567  // Loop over evaluation types
570  eqstm.begin();
571  int idx = 0;
572  for (; eval_type != eqstm.end(); ++eval_type,++idx) {
573  if (m_active_evaluation_types[idx]) {
574 
575  // Loop over integration rules
576  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir_iter = m_integration_rules.begin();
577  ir_iter != m_integration_rules.end(); ++ ir_iter) {
578 
579  Teuchos::RCP<panzer::IntegrationRule> ir = ir_iter->second;
580  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
581  eval_type->buildAndRegisterClosureModelEvaluators(fm, *m_field_lib->buildFieldLayoutLibrary(*ir), ir, factory, model_name, models, user_data);
582  eval_type->setDetailsIndex(di);
583  }
584  }
585  }
586  }
587 
588  // if there are no equation sets call the closure model directly
589  if(m_equation_sets.size()==0) {
591  int idx = 0;
592  for (;eval_type != factory.end(); ++eval_type,++idx) {
593  if (m_active_evaluation_types[idx]) {
594 
595  // setup some place holder parameter list, lets hope no one needs is!
597 
598  // Loop over integration rules
599  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir_iter = m_integration_rules.begin();
600  ir_iter != m_integration_rules.end(); ++ ir_iter) {
601 
602  Teuchos::RCP<panzer::IntegrationRule> ir = ir_iter->second;
603 
604  // call directly to the closure models
606  eval_type->buildClosureModels(model_name,
607  models,
608  *m_field_lib->buildFieldLayoutLibrary(*ir),
609  ir,
610  plist,
611  user_data,
612  this->globalData(),
613  fm);
614 
615  // register the constructed evaluators
616  const int di = eval_type->setDetailsIndex(this->getDetailsIndex());
617  eval_type->registerEvaluators(*evaluators,fm);
618  eval_type->setDetailsIndex(di);
619  }
620  }
621  }
622  }
623 }
624 
625 // *******************************************************************
629  const std::string& model_name,
630  const Teuchos::ParameterList& models,
632  const Teuchos::ParameterList& user_data) const
633 {
634  using namespace std;
635  using namespace panzer;
636  using namespace Teuchos;
637 
638  // Only use the <Residual> evaluation type, so pass through to type specific call
639  this->buildAndRegisterInitialConditionEvaluatorsForType<panzer::Traits::Residual>(fm, factory, model_name, models, lof, user_data);
640 }
641 
642 
643 // *******************************************************************
644 const std::vector<std::string>& panzer::PhysicsBlock::getDOFNames() const
645 {
646  return m_dof_names;
647 }
648 
649 // *******************************************************************
650 const std::vector<panzer::StrPureBasisPair>& panzer::PhysicsBlock::getProvidedDOFs() const
651 {
652  return m_provided_dofs;
653 }
654 
655 // *******************************************************************
656 const std::vector<std::vector<std::string> >& panzer::PhysicsBlock::getCoordinateDOFs() const
657 {
658  return m_coordinate_dofs;
659 }
660 
661 // *******************************************************************
662 const std::vector<panzer::StrPureBasisPair>& panzer::PhysicsBlock::getTangentFields() const
663 {
664  return m_tangent_fields;
665 }
666 
667 // *******************************************************************
669 {
670  panzer::WorksetNeeds needs;
671 
672  needs.cellData = this->cellData();
673  const std::map<int,Teuchos::RCP<panzer::IntegrationRule> >& int_rules = this->getIntegrationRules();
674  for (std::map<int,Teuchos::RCP<panzer::IntegrationRule> >::const_iterator ir_itr = int_rules.begin();
675  ir_itr != int_rules.end(); ++ir_itr)
676  needs.int_rules.push_back(ir_itr->second);
677 
678  const std::map<std::string,Teuchos::RCP<panzer::PureBasis> >& bases= this->getBases();
679  const std::vector<StrPureBasisPair>& fieldToBasis = getProvidedDOFs();
680  for(std::map<std::string,Teuchos::RCP<panzer::PureBasis> >::const_iterator b_itr = bases.begin();
681  b_itr != bases.end(); ++b_itr) {
682 
683  needs.bases.push_back(b_itr->second);
684 
685  bool found = false;
686  for(std::size_t d=0;d<fieldToBasis.size();d++) {
687  if(fieldToBasis[d].second->name()==b_itr->second->name()) {
688  // add representative basis for this field
689  needs.rep_field_name.push_back(fieldToBasis[d].first);
690  found = true;
691 
692  break;
693  }
694  }
695 
696  // this should always work if physics blocks are correctly constructed
697  TEUCHOS_ASSERT(found);
698  }
699 
700  return needs;
701 }
702 
703 // *******************************************************************
704 const std::map<std::string,Teuchos::RCP<panzer::PureBasis> >&
706 {
707  return m_bases;
708 }
709 
710 // *******************************************************************
711 const std::map<int,Teuchos::RCP<panzer::IntegrationRule> >&
713 {
714  return m_integration_rules;
715 }
716 
717 // *******************************************************************
718 const shards::CellTopology panzer::PhysicsBlock::getBaseCellTopology() const
719 {
720  TEUCHOS_TEST_FOR_EXCEPTION(m_bases.size() == 0, std::runtime_error,
721  "Cannot return a basis since none exist in this physics block.");
722  return m_bases.begin()->second->getIntrepid2Basis()->getBaseCellTopology();
723 }
724 
725 // *******************************************************************
727 {
728  return m_physics_id;
729 }
730 
731 // *******************************************************************
733 {
734  return m_element_block_id;
735 }
736 
737 // *******************************************************************
739 {
740  return m_cell_data;
741 }
742 
743 // *******************************************************************
745 {
746  return Teuchos::rcp(new panzer::PhysicsBlock(*this,cell_data));
747 }
748 
749 // *******************************************************************
751 {
752  return m_global_data;
753 }
754 
755 // *******************************************************************
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
PHX::TemplateManager< Traits::EvalTypes, panzer::EquationSetBase, panzer::EquationSet< _ > >::iterator begin()
Teuchos::RCP< panzer::GlobalData > globalData() const
void setActiveEvaluationTypes(const std::vector< bool > &aet)
Used to save memory by disabling unneeded evaluation types.
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.
PHX::TemplateManager< Traits::EvalTypes, panzer::EquationSetBase, panzer::EquationSet< _ > >::iterator end()
WorksetNeeds getWorksetNeeds() const
bool nonnull(const boost::shared_ptr< T > &p)
void activateAllEvaluationTypes()
Used to reactivate all evaluation types if some were temporarily disabled with a call to setActiveEva...
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
Teuchos::RCP< FieldLibrary > m_field_lib
#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