Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_FieldManagerBuilder.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #include <vector>
44 #include <string>
45 #include <sstream>
46 #include <fstream>
47 
49 
50 #include "Phalanx_DataLayout_MDALayout.hpp"
51 #include "Phalanx_FieldManager.hpp"
52 
53 #include "Teuchos_FancyOStream.hpp"
54 
55 #include "Shards_CellTopology.hpp"
56 
57 #include "Panzer_Traits.hpp"
58 #include "Panzer_Workset.hpp"
59 #include "Panzer_Workset_Builder.hpp"
60 #include "Panzer_PhysicsBlock.hpp"
64 #include "Panzer_CellData.hpp"
68 #include "Panzer_GlobalIndexer.hpp"
69 
70 //#include "EpetraExt_BlockMapOut.h"
71 
72 //=======================================================================
73 //=======================================================================
74 void panzer::FieldManagerBuilder::print(std::ostream& os) const
75 {
76  os << "panzer::FieldManagerBuilder output: Not implemented yet!";
77 }
78 
79 //=======================================================================
81 FieldManagerBuilder(bool disablePhysicsBlockScatter,
82  bool disablePhysicsBlockGather)
83  : disablePhysicsBlockScatter_(disablePhysicsBlockScatter)
84  , disablePhysicsBlockGather_(disablePhysicsBlockGather)
85  , active_evaluation_types_(Sacado::mpl::size<panzer::Traits::EvalTypes>::value, true)
86 {}
87 
88 //=======================================================================
89 namespace {
90  struct PostRegistrationFunctor {
91 
92  const std::vector<bool>& active_;
94  panzer::Traits::SD& setup_data_;
95 
96  PostRegistrationFunctor(const std::vector<bool>& active,
98  panzer::Traits::SD& setup_data)
99  : active_(active),fm_(fm),setup_data_(setup_data) {}
100 
101  template<typename T>
102  void operator()(T) const {
103  auto index = Sacado::mpl::find<panzer::Traits::EvalTypes,T>::value;
104  if (active_[index])
105  fm_.postRegistrationSetupForType<T>(setup_data_);
106  }
107  };
108 }
109 
110 //=======================================================================
112  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
113  const std::vector<WorksetDescriptor> & wkstDesc,
115  const Teuchos::ParameterList& closure_models,
116  const panzer::LinearObjFactory<panzer::Traits> & lo_factory,
117  const Teuchos::ParameterList& user_data,
118  const GenericEvaluatorFactory & gEvalFact,
119  bool closureModelByEBlock)
120 {
121  using Teuchos::RCP;
122  using Teuchos::rcp;
123 
124  TEUCHOS_TEST_FOR_EXCEPTION(getWorksetContainer()==Teuchos::null,std::logic_error,
125  "panzer::FMB::setupVolumeFieldManagers: method function getWorksetContainer() returns null. "
126  "Plase call setWorksetContainer() before calling this method");
127  TEUCHOS_TEST_FOR_EXCEPTION(physicsBlocks.size()!=wkstDesc.size(),std::runtime_error,
128  "panzer::FMB::setupVolumeFieldManagers: physics block count must match workset descriptor count.");
129 
130  phx_volume_field_managers_.clear();
131 
133 
134  for (std::size_t blkInd=0;blkInd<physicsBlocks.size();++blkInd) {
135  RCP<panzer::PhysicsBlock> pb = physicsBlocks[blkInd];
136  const WorksetDescriptor wd = wkstDesc[blkInd];
137 
138  Traits::SD setupData;
139  setupData.worksets_ = getWorksetContainer()->getWorksets(wd);
140  setupData.orientations_ = getWorksetContainer()->getOrientations();
141  if(setupData.worksets_->size()==0)
142  continue;
143 
144  // sanity check
146 
147  // build a field manager object
150 
151  // use the physics block to register active evaluators
152  pb->setActiveEvaluationTypes(active_evaluation_types_);
153  pb->buildAndRegisterEquationSetEvaluators(*fm, user_data);
154  if(!physicsBlockGatherDisabled())
155  pb->buildAndRegisterGatherAndOrientationEvaluators(*fm,lo_factory,user_data);
156  pb->buildAndRegisterDOFProjectionsToIPEvaluators(*fm,Teuchos::ptrFromRef(lo_factory),user_data);
157  if(!physicsBlockScatterDisabled())
158  pb->buildAndRegisterScatterEvaluators(*fm,lo_factory,user_data);
159 
160  if(closureModelByEBlock)
161  pb->buildAndRegisterClosureModelEvaluators(*fm,cm_factory,pb->elementBlockID(),closure_models,user_data);
162  else
163  pb->buildAndRegisterClosureModelEvaluators(*fm,cm_factory,closure_models,user_data);
164 
165  // Reset active evaluation types
167 
168  // register additional model evaluator from the generic evaluator factory
169  gEvalFact.registerEvaluators(*fm,wd,*pb);
170 
171  // setup derivative information
172  setKokkosExtendedDataTypeDimensions(wd.getElementBlock(),*globalIndexer,user_data,*fm);
173 
174  // call postRegistrationSetup() for each active type
175  Sacado::mpl::for_each_no_kokkos<panzer::Traits::EvalTypes>(PostRegistrationFunctor(active_evaluation_types_,*fm,setupData));
176 
177  // make sure to add the field manager & workset to the list
178  volume_workset_desc_.push_back(wd);
179  phx_volume_field_managers_.push_back(fm);
180  }
181 }
182 
183 //=======================================================================
185  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
187  const Teuchos::ParameterList& closure_models,
188  const panzer::LinearObjFactory<panzer::Traits> & lo_factory,
189  const Teuchos::ParameterList& user_data)
190 {
191  std::vector<WorksetDescriptor> wkstDesc;
192  for(std::size_t i=0;i<physicsBlocks.size();i++)
193  wkstDesc.push_back(blockDescriptor(physicsBlocks[i]->elementBlockID()));
194 
196  setupVolumeFieldManagers(physicsBlocks,wkstDesc,cm_factory,closure_models,lo_factory,user_data,eef);
197 }
198 
199 //=======================================================================
200 //=======================================================================
202 setupBCFieldManagers(const std::vector<panzer::BC> & bcs,
203  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks,
204  const Teuchos::Ptr<const panzer::EquationSetFactory>& /* eqset_factory */,
206  const panzer::BCStrategyFactory& bc_factory,
207  const Teuchos::ParameterList& closure_models,
208  const panzer::LinearObjFactory<panzer::Traits> & lo_factory,
209  const Teuchos::ParameterList& user_data)
210 {
211  TEUCHOS_TEST_FOR_EXCEPTION(getWorksetContainer()==Teuchos::null,std::logic_error,
212  "panzer::FMB::setupBCFieldManagers: method function getWorksetContainer() returns null. "
213  "Plase call setWorksetContainer() before calling this method");
214 
216 
217  // for convenience build a map (element block id => physics block)
218  std::map<std::string,Teuchos::RCP<panzer::PhysicsBlock> > physicsBlocks_map;
219  {
220  std::vector<Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator blkItr;
221  for(blkItr=physicsBlocks.begin();blkItr!=physicsBlocks.end();++blkItr) {
223  std::string blockId = pb->elementBlockID();
224 
225  // add block id, physics block pair to the map
226  physicsBlocks_map.insert(std::make_pair(blockId,pb));
227  }
228  }
229 
230  // ***************************
231  // BCs
232  // ***************************
233  std::vector<panzer::BC>::const_iterator bc;
234  for (bc=bcs.begin(); bc != bcs.end(); ++bc) {
237  currentWkst = getWorksetContainer()->getSideWorksets(wd);
238  if (currentWkst.is_null()) continue;
239 
240  BCType bc_type = bc->bcType();
241 
242  if (bc_type == BCT_Interface) {
243  // Loop over local face indices and setup each field manager
244  for (std::map<unsigned,panzer::Workset>::const_iterator wkst = currentWkst->begin();
245  wkst != currentWkst->end(); ++wkst) {
246  // Build one FieldManager for each local side workset for each bc
247  std::map<unsigned,PHX::FieldManager<panzer::Traits> >& field_managers =
248  bc_field_managers_[*bc];
249 
250  PHX::FieldManager<panzer::Traits>& fm = field_managers[wkst->first];
251 
252  int gid_count = 0;
253  for (int block_id_index = 0; block_id_index < 2; ++block_id_index) {
254  const std::string element_block_id = block_id_index == 0 ? bc->elementBlockID() : bc->elementBlockID2();
255 
256  std::map<std::string,Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator
257  volume_pb_itr = physicsBlocks_map.find(element_block_id);
258 
259  TEUCHOS_TEST_FOR_EXCEPTION(volume_pb_itr == physicsBlocks_map.end(), std::logic_error,
260  "panzer::FMB::setupBCFieldManagers: Cannot find physics block corresponding to element block \""
261  << element_block_id << "\"");
262 
263  const Teuchos::RCP<const panzer::PhysicsBlock> volume_pb = physicsBlocks_map.find(element_block_id)->second;
264  const Teuchos::RCP<const shards::CellTopology> volume_cell_topology = volume_pb->cellData().getCellTopology();
265 
266  // register evaluators from strategy
267  const panzer::CellData side_cell_data(wkst->second.num_cells,
268  wkst->second.details(block_id_index).subcell_index,
269  volume_cell_topology);
270 
271  // Copy the physics block for side integrations
272  Teuchos::RCP<panzer::PhysicsBlock> side_pb = volume_pb->copyWithCellData(side_cell_data);
273 
275  bcstm = bc_factory.buildBCStrategy(*bc, side_pb->globalData());
276 
277  // Iterate over evaluation types
278  int i=0;
280  bcs_type = bcstm->begin(); bcs_type != bcstm->end(); ++bcs_type,++i) {
281  if (active_evaluation_types_[i]) {
282  bcs_type->setDetailsIndex(block_id_index);
283  side_pb->setDetailsIndex(block_id_index);
284  bcs_type->setup(*side_pb, user_data);
285  bcs_type->buildAndRegisterEvaluators(fm, *side_pb, cm_factory, closure_models, user_data);
286  bcs_type->buildAndRegisterGatherAndOrientationEvaluators(fm, *side_pb, lo_factory, user_data);
287  if ( ! physicsBlockScatterDisabled())
288  bcs_type->buildAndRegisterScatterEvaluators(fm, *side_pb, lo_factory, user_data);
289  }
290  }
291 
292  gid_count += globalIndexer->getElementBlockGIDCount(element_block_id);
293  }
294 
295  { // Use gid_count to set up the derivative information.
296  std::vector<PHX::index_size_type> derivative_dimensions;
297  derivative_dimensions.push_back(gid_count);
298  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Jacobian>(derivative_dimensions);
299 
300  #ifdef Panzer_BUILD_HESSIAN_SUPPORT
301  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Hessian>(derivative_dimensions);
302  #endif
303 
304  derivative_dimensions[0] = 1;
305  if (user_data.isType<int>("Tangent Dimension"))
306  derivative_dimensions[0] = user_data.get<int>("Tangent Dimension");
307  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Tangent>(derivative_dimensions);
308  }
309 
310  // Set up the field manager
311  Traits::SD setupData;
312  Teuchos::RCP<std::vector<panzer::Workset> > worksets = Teuchos::rcp(new std::vector<panzer::Workset>);
313  worksets->push_back(wkst->second);
314  setupData.worksets_ = worksets;
315  setupData.orientations_ = getWorksetContainer()->getOrientations();
316 
317  Sacado::mpl::for_each_no_kokkos<panzer::Traits::EvalTypes>(PostRegistrationFunctor(active_evaluation_types_,fm,setupData));
318 
319  }
320  } else {
321  const std::string element_block_id = bc->elementBlockID();
322 
323  std::map<std::string,Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator volume_pb_itr
324  = physicsBlocks_map.find(element_block_id);
325 
326  TEUCHOS_TEST_FOR_EXCEPTION(volume_pb_itr==physicsBlocks_map.end(),std::logic_error,
327  "panzer::FMB::setupBCFieldManagers: Cannot find physics block corresponding to element block \"" << element_block_id << "\"");
328 
329  Teuchos::RCP<const panzer::PhysicsBlock> volume_pb = physicsBlocks_map.find(element_block_id)->second;
330  Teuchos::RCP<const shards::CellTopology> volume_cell_topology = volume_pb->cellData().getCellTopology();
331 
332  // Build one FieldManager for each local side workset for each dirichlet bc
333  std::map<unsigned,PHX::FieldManager<panzer::Traits> >& field_managers =
334  bc_field_managers_[*bc];
335 
336  // Loop over local face indices and setup each field manager
337  for (std::map<unsigned,panzer::Workset>::const_iterator wkst =
338  currentWkst->begin(); wkst != currentWkst->end();
339  ++wkst) {
340 
341  PHX::FieldManager<panzer::Traits>& fm = field_managers[wkst->first];
342 
343  // register evaluators from strategy
344  const panzer::CellData side_cell_data(wkst->second.num_cells,
345  wkst->first,volume_cell_topology);
346 
347  // Copy the physics block for side integrations
348  Teuchos::RCP<panzer::PhysicsBlock> side_pb = volume_pb->copyWithCellData(side_cell_data);
349 
351  bc_factory.buildBCStrategy(*bc,side_pb->globalData());
352 
353  // Iterate over evaluation types
354  int i=0;
356  bcs_type = bcstm->begin(); bcs_type != bcstm->end(); ++bcs_type,++i) {
357  if (active_evaluation_types_[i]) {
358  bcs_type->setup(*side_pb,user_data);
359  bcs_type->buildAndRegisterEvaluators(fm,*side_pb,cm_factory,closure_models,user_data);
360  bcs_type->buildAndRegisterGatherAndOrientationEvaluators(fm,*side_pb,lo_factory,user_data);
361  if(!physicsBlockScatterDisabled())
362  bcs_type->buildAndRegisterScatterEvaluators(fm,*side_pb,lo_factory,user_data);
363  }
364  }
365 
366  // Setup the fieldmanager
367  Traits::SD setupData;
369  Teuchos::rcp(new(std::vector<panzer::Workset>));
370  worksets->push_back(wkst->second);
371  setupData.worksets_ = worksets;
372  setupData.orientations_ = getWorksetContainer()->getOrientations();
373 
374  // setup derivative information
375  setKokkosExtendedDataTypeDimensions(element_block_id,*globalIndexer,user_data,fm);
376 
377  Sacado::mpl::for_each_no_kokkos<panzer::Traits::EvalTypes>(PostRegistrationFunctor(active_evaluation_types_,fm,setupData));
378  }
379  }
380  }
381 }
382 
383 //=======================================================================
384 //=======================================================================
386 writeVolumeGraphvizDependencyFiles(std::string filename_prefix,
387  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks) const
388 {
389  if(phx_volume_field_managers_.size()<1)
390  return; // nothing to see here folks
391 
392  TEUCHOS_ASSERT(phx_volume_field_managers_.size()==physicsBlocks.size());
393 
394  std::vector<Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator blkItr;
395  int index = 0;
396  for (blkItr=physicsBlocks.begin();blkItr!=physicsBlocks.end();++blkItr,++index) {
397  std::string blockId = (*blkItr)->elementBlockID();
398  phx_volume_field_managers_[index]->writeGraphvizFile(filename_prefix+"_VOLUME_"+blockId);
399  }
400 
401 }
402 
403 //=======================================================================
404 //=======================================================================
406 writeBCGraphvizDependencyFiles(std::string filename_prefix) const
407 {
408  typedef std::map<panzer::BC,std::map<unsigned,PHX::FieldManager<panzer::Traits> >,panzer::LessBC> FMMap;
409 
410  FMMap::const_iterator blkItr;
411  int bc_index = 0;
412  for (blkItr=bc_field_managers_.begin();blkItr!=bc_field_managers_.end();++blkItr,++bc_index) {
413  panzer::BC bc = blkItr->first;
414  const PHX::FieldManager<panzer::Traits> & fm = blkItr->second.begin()->second; // get the first field manager
415 
416  BCType bc_type = bc.bcType();
417  std::string type;
418  if (bc_type == BCT_Dirichlet)
419  type = "_Dirichlet_";
420  else if (bc_type == BCT_Neumann)
421  type = "_Neumann_";
422  else if (bc_type == BCT_Interface)
423  type = "_Interface_";
424  else
425  TEUCHOS_ASSERT(false);
426 
427  std::string blockId = bc.elementBlockID();
428  std::string sideId = bc.sidesetID();
429  fm.writeGraphvizFile(filename_prefix+"_BC_"+std::to_string(bc_index)+type+sideId+"_"+blockId);
430  }
431 
432 }
433 
434 //=======================================================================
435 //=======================================================================
437 writeVolumeTextDependencyFiles(std::string filename_prefix,
438  const std::vector<Teuchos::RCP<panzer::PhysicsBlock> >& physicsBlocks) const
439 {
440  if(phx_volume_field_managers_.size()<1)
441  return; // nothing to see here folks
442 
443  TEUCHOS_ASSERT(phx_volume_field_managers_.size()==physicsBlocks.size());
444 
445  std::vector<Teuchos::RCP<panzer::PhysicsBlock> >::const_iterator blkItr;
446  int index = 0;
447  for (blkItr=physicsBlocks.begin();blkItr!=physicsBlocks.end();++blkItr,++index) {
448 
449  std::string blockId = (*blkItr)->elementBlockID();
450 
451  std::string filename = filename_prefix+"_VOLUME_"+blockId+".txt";
452  std::ofstream ofs;
453  ofs.open(filename.c_str());
454 
455  ofs << *(phx_volume_field_managers_[index]) << std::endl;
456 
457  ofs.close();
458  }
459 
460 }
461 
462 //=======================================================================
463 //=======================================================================
465 writeBCTextDependencyFiles(std::string filename_prefix) const
466 {
467  typedef std::map<panzer::BC,std::map<unsigned,PHX::FieldManager<panzer::Traits> >,panzer::LessBC> FMMap;
468 
469  FMMap::const_iterator blkItr;
470  int bc_index = 0;
471  for (blkItr=bc_field_managers_.begin();blkItr!=bc_field_managers_.end();++blkItr,++bc_index) {
472  panzer::BC bc = blkItr->first;
473  const PHX::FieldManager<panzer::Traits> & fm = blkItr->second.begin()->second; // get the first field manager
474 
475  BCType bc_type = bc.bcType();
476  std::string type;
477  if (bc_type == BCT_Dirichlet)
478  type = "_Dirichlet_";
479  else if (bc_type == BCT_Neumann)
480  type = "_Neumann_";
481  else if (bc_type == BCT_Interface)
482  type = "_Interface_";
483  else
484  TEUCHOS_ASSERT(false);
485 
486  std::string blockId = bc.elementBlockID();
487  std::string sideId = bc.sidesetID();
488 
489  std::string filename = filename_prefix+"_BC_"+std::to_string(bc_index)+type+sideId+"_"+blockId+".txt";
490  std::ofstream ofs;
491  ofs.open(filename.c_str());
492 
493  ofs << fm << std::endl;
494 
495  ofs.close();
496  }
497 
498 }
499 
500 //=======================================================================
501 //=======================================================================
503 setKokkosExtendedDataTypeDimensions(const std::string & eblock,
504  const panzer::GlobalIndexer & globalIndexer,
505  const Teuchos::ParameterList& user_data,
507 {
508  // setup Jacobian derivative terms
509  {
510  std::vector<PHX::index_size_type> derivative_dimensions;
511  derivative_dimensions.push_back(globalIndexer.getElementBlockGIDCount(eblock));
512 
513  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Jacobian>(derivative_dimensions);
514 
515  }
516 
517  #ifdef Panzer_BUILD_HESSIAN_SUPPORT
518  {
519  std::vector<PHX::index_size_type> derivative_dimensions;
520  derivative_dimensions.push_back(globalIndexer.getElementBlockGIDCount(eblock));
521 
522  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Hessian>(derivative_dimensions);
523  }
524  #endif
525 
526  {
527  std::vector<PHX::index_size_type> derivative_dimensions;
528  derivative_dimensions.push_back(1);
529  if (user_data.isType<int>("Tangent Dimension"))
530  derivative_dimensions[0] = user_data.get<int>("Tangent Dimension");
531  fm.setKokkosExtendedDataTypeDimensions<panzer::Traits::Tangent>(derivative_dimensions);
532  }
533 }
534 
536 {active_evaluation_types_ = aet;}
537 
538 //=======================================================================
539 //=======================================================================
541 {
542  phx_volume_field_managers_.clear();
543  volume_workset_desc_.clear();
544  if (clearVolumeWorksets)
545  worksetContainer_->clearVolumeWorksets();
546 }
547 
548 //=======================================================================
549 //=======================================================================
550 std::ostream& panzer::operator<<(std::ostream& os, const panzer::FieldManagerBuilder& rfd)
551 {
552  rfd.print(os);
553  return os;
554 }
Interface for constructing a BCStrategy_TemplateManager.
void setupVolumeFieldManagers(const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const Teuchos::ParameterList &closure_models, const LinearObjFactory< panzer::Traits > &lo_factory, const Teuchos::ParameterList &user_data)
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
BCType
Type of boundary condition.
Definition: Panzer_BC.hpp:74
std::string elementBlockID() const
Returns the element block id associated with this sideset.
Definition: Panzer_BC.cpp:200
virtual Teuchos::RCP< const panzer::GlobalIndexer > getRangeGlobalIndexer() const =0
Get the range global indexer object associated with this factory.
T & get(const std::string &name, T def_value)
void writeVolumeGraphvizDependencyFiles(std::string filename_prefix, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks) const
virtual int getElementBlockGIDCount(const std::size_t &blockIndex) const =0
How any GIDs are associate with a particular element block.
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
const std::string & getElementBlock(const int block=0) const
Get element block name.
void writeBCGraphvizDependencyFiles(std::string filename_prefix) const
void clearVolumeFieldManagers(bool clearVolumeWorksets=true)
Delete all volume field managers, retaining the BC ones.
void setKokkosExtendedDataTypeDimensions(const std::string &eblock, const panzer::GlobalIndexer &globalIndexer, const Teuchos::ParameterList &user_data, PHX::FieldManager< panzer::Traits > &fm) const
void buildAndRegisterGatherAndOrientationEvaluators(PHX::FieldManager< panzer::Traits > &fm, const panzer::LinearObjFactory< panzer::Traits > &lof, const Teuchos::ParameterList &user_data) const
virtual Teuchos::RCP< panzer::BCStrategy_TemplateManager< panzer::Traits > > buildBCStrategy(const panzer::BC &bc, const Teuchos::RCP< panzer::GlobalData > &global_data) const =0
void writeBCTextDependencyFiles(std::string filename_prefix) const
Teuchos::RCP< panzer::GlobalData > globalData() const
void setActiveEvaluationTypes(const std::vector< bool > &aet)
Used to save memory by disabling unneeded evaluation types.
void buildAndRegisterDOFProjectionsToIPEvaluators(PHX::FieldManager< panzer::Traits > &fm, const Teuchos::Ptr< const panzer::LinearObjFactory< panzer::Traits > > &lof, const Teuchos::ParameterList &user_data) const
void print(std::ostream &os) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Data for determining cell topology and dimensionality.
void writeVolumeTextDependencyFiles(std::string filename_prefix, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks) const
std::string elementBlockID() const
BCType bcType() const
Returns the boundary condition type (Dirichlet or Neumann or Interface).
Definition: Panzer_BC.cpp:186
virtual bool registerEvaluators(PHX::FieldManager< panzer::Traits > &fm, const WorksetDescriptor &wd, const PhysicsBlock &pb) const =0
void activateAllEvaluationTypes()
Used to reactivate all evaluation types if some were temporarily disabled with a call to setActiveEva...
std::string sidesetID() const
Returns the set id.
Definition: Panzer_BC.cpp:193
std::ostream & operator<<(std::ostream &os, const AssemblyEngineInArgs &in)
bool isType(const std::string &name) const
FieldManagerBuilder(bool disablePhysicsBlockScatter=false, bool disablePhysicsBlockGather=false)
WorksetDescriptor blockDescriptor(const std::string &eBlock)
WorksetDescriptor bcDescriptor(const panzer::BC &bc)
Definition: Panzer_BC.cpp:331
Stores input information for a boundary condition.
Definition: Panzer_BC.hpp:81
void setupBCFieldManagers(const std::vector< panzer::BC > &bcs, const std::vector< Teuchos::RCP< panzer::PhysicsBlock > > &physicsBlocks, const panzer::EquationSetFactory &eqset_factory, const panzer::ClosureModelFactory_TemplateManager< panzer::Traits > &cm_factory, const panzer::BCStrategyFactory &bc_factory, const Teuchos::ParameterList &closure_models, const LinearObjFactory< panzer::Traits > &lo_factory, const Teuchos::ParameterList &user_data)
#define TEUCHOS_ASSERT(assertion_test)
void setActiveEvaluationTypes(const std::vector< bool > &aet)
Set a vector of active evaluation types to allocate.
void buildAndRegisterEquationSetEvaluators(PHX::FieldManager< panzer::Traits > &fm, const Teuchos::ParameterList &user_data) const
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
bool is_null() const