GRINS-0.8.0
List of all members | Public Member Functions | Protected Member Functions
GRINS::DefaultBCBuilder Class Reference

Manages runtime construction of Dirichlet boundary conditions. More...

#include <default_bc_builder.h>

Inheritance diagram for GRINS::DefaultBCBuilder:
Inheritance graph
[legend]
Collaboration diagram for GRINS::DefaultBCBuilder:
Collaboration graph
[legend]

Public Member Functions

 DefaultBCBuilder ()
 
 ~DefaultBCBuilder ()
 
- Public Member Functions inherited from GRINS::BCBuilder
 BCBuilder ()
 
virtual ~BCBuilder ()
 
- Public Member Functions inherited from GRINS::BuilderHelper
 BuilderHelper ()
 
 ~BuilderHelper ()
 

Protected Member Functions

virtual void build_bcs (const GetPot &input, MultiphysicsSystem &system, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 
void build_type_based_bcs (const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, const std::string &type_input_section, std::set< std::string > &var_sections, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 Helper function to build boundary conditions specified by a single type. More...
 
void build_axisymmetric_bcs (const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, const std::string &bc_type, std::set< std::string > &var_sections, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 
void build_bcs_by_var_section (const GetPot &input, MultiphysicsSystem &system, const std::string &bc_name, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, std::set< std::string > &var_sections, const std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &bc_id_to_subdomain_id_map, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 Helper function to build boundary conditions using Variable sections. More...
 
void parse_and_build_bc_id_map (const GetPot &input, std::map< std::string, std::set< BoundaryID > > &bc_id_map)
 
void verify_bc_ids_with_mesh (const MultiphysicsSystem &system, const std::map< std::string, std::set< BoundaryID > > &bc_id_map) const
 
void build_periodic_bc (const GetPot &input, libMesh::System &system, const std::set< BoundaryID > &bc_ids, const std::string &section)
 
void parse_periodic_master_slave_ids (const GetPot &input, const std::string &section, libMesh::boundary_id_type &master_id, libMesh::boundary_id_type &slave_id) const
 
libMesh::RealVectorValue parse_periodic_offset (const GetPot &input, const std::string &section) const
 
void build_bc_to_subdomain_map_check_with_mesh (const MultiphysicsSystem &system, std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &bc_id_to_subdomain_id_map) const
 Build up bc_id to subdomain_id map. More...
 
bool is_var_active (const FEVariablesBase &var, const std::vector< libMesh::subdomain_id_type > &subdomain_ids) const
 Check if the Variable var is active on the given subdomain_id. More...
 
- Protected Member Functions inherited from GRINS::BCBuilder
void construct_dbc_core (const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, libMesh::DofMap &dof_map)
 
void construct_nbc_core (const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 
bool is_dirichlet_bc_type (const std::string &bc_type)
 
bool is_neumann_bc_type (const std::string &bc_type)
 
void add_periodic_bc_to_dofmap (libMesh::boundary_id_type master_id, libMesh::boundary_id_type slave_id, const libMesh::RealVectorValue &offset_vector, libMesh::DofMap &dof_map)
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::BCBuilder
static void build_boundary_conditions (const GetPot &input, MultiphysicsSystem &system, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
 
- Static Public Member Functions inherited from GRINS::BuilderHelper
static void parse_var_sections (const GetPot &input, std::set< std::string > &sections)
 Parses the input file for [Variables] first-level subsections. More...
 
static void parse_var_sections_vector (const GetPot &input, std::vector< std::string > &sections)
 The same as parse_var_sections, except the result is returned in an std::vector. More...
 
- Static Protected Member Functions inherited from GRINS::BCBuilder
static bool is_new_bc_input_style (const GetPot &input)
 
static libMesh::UniquePtr< BCBuilderbuild_builder (const GetPot &input)
 

Detailed Description

Manages runtime construction of Dirichlet boundary conditions.

This will parse the input for the request Dirichlet boundary conditions and manage their construction. Actual construction of the DirichletBoundary objects is delegated to factory classes. This builder classes merely manages tasks around the factories as needed. To add new Dirichlet boundary conditions, the user should instantiate an appropriate factory sub class.

Definition at line 47 of file default_bc_builder.h.

Constructor & Destructor Documentation

GRINS::DefaultBCBuilder::DefaultBCBuilder ( )
inline

Definition at line 52 of file default_bc_builder.h.

53  : BCBuilder(),
55  {};
GRINS::DefaultBCBuilder::~DefaultBCBuilder ( )
inline

Definition at line 57 of file default_bc_builder.h.

57 {};

Member Function Documentation

void GRINS::DefaultBCBuilder::build_axisymmetric_bcs ( const GetPot &  input,
MultiphysicsSystem system,
const std::set< BoundaryID > &  bc_ids,
libMesh::DofMap &  dof_map,
const std::string &  bc_type,
std::set< std::string > &  var_sections,
std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs 
)
protected

Definition at line 142 of file default_bc_builder.C.

References GRINS::BCBuilder::construct_dbc_core(), GRINS::BCBuilder::construct_nbc_core(), GRINS::GRINSPrivate::VariableWarehouse::get_variable(), GRINS::FEVariablesBase::is_constraint_var(), GRINS::BCBuilder::is_dirichlet_bc_type(), and GRINS::BCBuilder::is_neumann_bc_type().

Referenced by build_type_based_bcs().

149  {
150  // Axisymmetric depends on the variable type.
151  // So, we loop over all the Variable names and prepend the type
152  // with the variable name. So, the factories for axisymmetric
153  // need to be instantiated with, e.g., Velocity_axisymmetric
154  for( std::set<std::string>::const_iterator vars = var_sections.begin();
155  vars != var_sections.end(); ++vars )
156  {
157  const std::string& var_section = *vars;
158 
159  // Grab FEVariable
160  const FEVariablesBase& fe_var =
162 
163  // We don't need to do anything for constraint variables
164  if( fe_var.is_constraint_var() )
165  continue;
166 
167  std::string full_bc_type = var_section+"_"+bc_type;
168 
169  // Axisymmetric is Dirichlet or Neumann, depending on the Variable
170  if( this->is_dirichlet_bc_type(full_bc_type) )
171  {
172  this->construct_dbc_core( input,
173  system,
174  bc_ids,
175  fe_var,
176  std::string("dummy"), // No section needed for axisymmetric
177  full_bc_type,
178  dof_map );
179  }
180  else if( this->is_neumann_bc_type(full_bc_type) )
181  {
182  this->construct_nbc_core( input,
183  system,
184  bc_ids,
185  fe_var,
186  std::string("dummy"), // No section needed for axisymmetric
187  full_bc_type,
188  neumann_bcs );
189  }
190  else
191  libmesh_error_msg("ERROR! Invalid axisymmetric type: "+full_bc_type);
192  }
193  }
void construct_dbc_core(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, libMesh::DofMap &dof_map)
Definition: bc_builder.C:76
void construct_nbc_core(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Definition: bc_builder.C:106
bool is_dirichlet_bc_type(const std::string &bc_type)
Definition: bc_builder.C:139
static FEVariablesBase & get_variable(const std::string &var_name)
bool is_neumann_bc_type(const std::string &bc_type)
Definition: bc_builder.C:145
void GRINS::DefaultBCBuilder::build_bc_to_subdomain_map_check_with_mesh ( const MultiphysicsSystem system,
std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &  bc_id_to_subdomain_id_map 
) const
protected

Build up bc_id to subdomain_id map.

we also check and make sure that there's only one subdomain id per boundary id. If not, we throw an error.

Definition at line 463 of file default_bc_builder.C.

Referenced by build_bcs().

465  {
466  std::vector<libMesh::dof_id_type> elem_ids;
467  std::vector<unsigned short int> side_ids;
468  std::vector<BoundaryID> bc_ids;
469 
470  const libMesh::MeshBase& mesh = system.get_mesh();
471 
472  // Extract the list of elements on the boundary
473  mesh.get_boundary_info().build_active_side_list( elem_ids, side_ids, bc_ids );
474 
475  libmesh_assert_equal_to( elem_ids.size(), side_ids.size() );
476  libmesh_assert_equal_to( elem_ids.size(), bc_ids.size() );
477 
478  for( unsigned int i = 0; i < elem_ids.size(); i++ )
479  {
480  const libMesh::Elem* elem_ptr = mesh.elem_ptr(elem_ids[i]);
481 
482  libMesh::subdomain_id_type subdomain_id = elem_ptr->subdomain_id();
483 
484  std::map<BoundaryID,std::vector<libMesh::subdomain_id_type> >::iterator
485  map_it = bc_id_to_subdomain_id_map.find(bc_ids[i]);
486 
487  if( map_it == bc_id_to_subdomain_id_map.end() )
488  {
489  std::vector<libMesh::subdomain_id_type> sid(1,subdomain_id);
490  bc_id_to_subdomain_id_map.insert(std::make_pair(bc_ids[i],sid));
491  }
492  else
493  {
494  std::vector<libMesh::subdomain_id_type>& sids = map_it->second;
495  bool found_sid = false;
496 
497  for( unsigned int i = 0; i < sids.size(); i++ )
498  {
499  if( sids[i] == subdomain_id )
500  found_sid = true;
501  }
502 
503  if( !found_sid && (sids.size() > 2) )
504  libmesh_error_msg("ERROR: How do you have more than 2 subdomain ids for a boundary id?!");
505  else if( !found_sid && sids.size() < 2 )
506  sids.push_back(subdomain_id);
507  }
508  }
509  }
void GRINS::DefaultBCBuilder::build_bcs ( const GetPot &  input,
MultiphysicsSystem system,
std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs 
)
protectedvirtual

Implements GRINS::BCBuilder.

Definition at line 46 of file default_bc_builder.C.

References GRINS::BoundaryConditionNames::bc_section(), build_bc_to_subdomain_map_check_with_mesh(), build_bcs_by_var_section(), build_type_based_bcs(), parse_and_build_bc_id_map(), GRINS::BuilderHelper::parse_var_sections(), and verify_bc_ids_with_mesh().

48  {
49  libMesh::DofMap& dof_map = system.get_dof_map();
50 
51  // Setup map between boundary section name and the boundary ids
52  std::map<std::string,std::set<BoundaryID> > bc_id_map;
53  this->parse_and_build_bc_id_map( input, bc_id_map );
54 
55  // The user must have mapped all boundary ids present in the mesh.
56  // Additionally, they can't have specified an boundary id that doesn't
57  // exist.
58  this->verify_bc_ids_with_mesh( system, bc_id_map );
59 
60  // Parse variable sections. These correspond to the names of FEVariables
61  // We will be looking for all variable sections in all boundary sections
62  // that are not "constraint" variables (i.e. FEVariablesBase::_is_constraint)
63  std::set<std::string> var_sections;
64  this->parse_var_sections( input, var_sections );
65 
66  // Cache the map between boundary ids and subdomain ids
67  // We map to a vector because it's possible to have an "interface"
68  // boundary id that touches multiple elements with differing
69  // subdomain ids
70  std::map<BoundaryID,std::vector<libMesh::subdomain_id_type> >
71  bc_id_to_subdomain_id_map;
72 
74  bc_id_to_subdomain_id_map );
75 
76  for( std::map<std::string,std::set<BoundaryID> >::const_iterator bc_it = bc_id_map.begin();
77  bc_it != bc_id_map.end(); ++bc_it )
78  {
79  const std::string& bc_name = bc_it->first;
80 
81  const std::set<BoundaryID>& bc_ids = bc_it->second;
82 
83  // Check for special types of boundary conditions that can be
84  // specified as a single type that applies to all boundaries
85  std::string type_input_section = BoundaryConditionNames::bc_section()+"/"+bc_name;
86  std::string type_input = type_input_section+"/type";
87 
88  if( input.have_variable(type_input) )
89  {
90  this->build_type_based_bcs(input,system,bc_ids,dof_map,
91  type_input_section,var_sections,
92  neumann_bcs);
93  }
94  // Otherwise, we look for each Variable section and the
95  // boundary conditions specified therein
96  else
97  {
98  this->build_bcs_by_var_section(input,system,bc_name,bc_ids,dof_map,
99  var_sections,bc_id_to_subdomain_id_map,
100  neumann_bcs);
101  }
102  }
103  }
void parse_and_build_bc_id_map(const GetPot &input, std::map< std::string, std::set< BoundaryID > > &bc_id_map)
void build_type_based_bcs(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, const std::string &type_input_section, std::set< std::string > &var_sections, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Helper function to build boundary conditions specified by a single type.
void verify_bc_ids_with_mesh(const MultiphysicsSystem &system, const std::map< std::string, std::set< BoundaryID > > &bc_id_map) const
void build_bcs_by_var_section(const GetPot &input, MultiphysicsSystem &system, const std::string &bc_name, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, std::set< std::string > &var_sections, const std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &bc_id_to_subdomain_id_map, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Helper function to build boundary conditions using Variable sections.
void build_bc_to_subdomain_map_check_with_mesh(const MultiphysicsSystem &system, std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &bc_id_to_subdomain_id_map) const
Build up bc_id to subdomain_id map.
static std::string bc_section()
Outer section name for boundary conditions section in input file.
static void parse_var_sections(const GetPot &input, std::set< std::string > &sections)
Parses the input file for [Variables] first-level subsections.
void GRINS::DefaultBCBuilder::build_bcs_by_var_section ( const GetPot &  input,
MultiphysicsSystem system,
const std::string &  bc_name,
const std::set< BoundaryID > &  bc_ids,
libMesh::DofMap &  dof_map,
std::set< std::string > &  var_sections,
const std::map< BoundaryID, std::vector< libMesh::subdomain_id_type > > &  bc_id_to_subdomain_id_map,
std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs 
)
protected

Helper function to build boundary conditions using Variable sections.

This is the "standard" part. We parse for each Variable section that should have boundary conditions and then parse the boundary condition type.

Definition at line 195 of file default_bc_builder.C.

References GRINS::BoundaryConditionNames::bc_section(), GRINS::BCBuilder::construct_dbc_core(), GRINS::BCBuilder::construct_nbc_core(), GRINS::GRINSPrivate::VariableWarehouse::get_variable(), GRINS::FEVariablesBase::is_constraint_var(), GRINS::BCBuilder::is_dirichlet_bc_type(), GRINS::BCBuilder::is_neumann_bc_type(), and is_var_active().

Referenced by build_bcs().

203  {
204  for( std::set<std::string>::const_iterator vars = var_sections.begin();
205  vars != var_sections.end(); ++vars )
206  {
207  const std::string& var_section = *vars;
208 
209  // Setup the input section we're in. We use this, but
210  // this is also the section where the BCFactory will try
211  // to parse the values, if it needs to, so let's set up
212  // once and reuse it.
213  std::string input_section = std::string(BoundaryConditionNames::bc_section()+"/"+bc_name+"/"+var_section);
214 
215  // All the boundary ids have the same subdomain id (this is checked
216  // earlier), so just grab the first one
217  std::map<BoundaryID,std::vector<libMesh::subdomain_id_type> >::const_iterator
218  subdomain_ids_iter = bc_id_to_subdomain_id_map.find(*bc_ids.begin());
219 
220  // If we're on a DistributedMesh we might need to learn about
221  // subdomain ids from other processors.
222  std::vector<libMesh::subdomain_id_type> subdomain_ids;
223  if (subdomain_ids_iter != bc_id_to_subdomain_id_map.end())
224  subdomain_ids = subdomain_ids_iter->second;
225  system.comm().allgather(subdomain_ids);
226  std::sort( subdomain_ids.begin(), subdomain_ids.end() );
227  subdomain_ids.erase
228  (std::unique(subdomain_ids.begin(), subdomain_ids.end()),
229  subdomain_ids.end() );
230 
231  // Grab the FEVariable
232  const FEVariablesBase& fe_var =
234 
235  bool var_active = this->is_var_active( fe_var, subdomain_ids );
236 
237  // If the variable is *not* active and the section is there,
238  // that's an error.
239  if( !var_active && input.have_section(input_section) )
240  {
241  std::stringstream error_msg;
242  error_msg << "ERROR: Cannot specify boundary condition for variable "
243  << var_section << " on boundary " << bc_name << std::endl
244  << "since it is inactive on the subdomain associated "
245  << "with this boundary." << std::endl;
246  libmesh_error_msg(error_msg.str());
247  }
248 
249  // Make sure this section is there,
250  // unless that variable is a constraint variable
251  // or it's not enabled on this subdomain
252  if( !input.have_section(input_section) )
253  {
254  if( fe_var.is_constraint_var() || !var_active )
255  continue;
256  else
257  libmesh_error_msg("ERROR: Could not find boundary condition specification for "+input_section+"!");
258 
259  }
260 
261  // Grab the type of the boundary condition
262  // There may be more than one type (e.g. pin displacement in
263  // two directions and load in the third direction).
264  std::string bc_type_section = input_section+"/type";
265  unsigned int n_bc_types = input.vector_variable_size(bc_type_section);
266 
267  for( unsigned int bc_type_idx = 0; bc_type_idx < n_bc_types; bc_type_idx++ )
268  {
269  std::string bc_type = input(bc_type_section, std::string("DIE!"), bc_type_idx);
270 
271  if( this->is_dirichlet_bc_type(bc_type) )
272  {
273  this->construct_dbc_core( input,system, bc_ids,
274  fe_var, input_section,
275  bc_type, dof_map );
276  }
277  else if( this->is_neumann_bc_type(bc_type) )
278  {
279  this->construct_nbc_core( input,system, bc_ids,
280  fe_var, input_section,
281  bc_type, neumann_bcs );
282  }
283  else
284  libmesh_error_msg("ERROR: Invalid bc_type "+bc_type+"!");
285 
286  } // end loop over bc_types
287 
288  } // end loop over variable sections
289  }
void construct_dbc_core(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, libMesh::DofMap &dof_map)
Definition: bc_builder.C:76
void construct_nbc_core(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, const FEVariablesBase &fe_var, const std::string &section, const std::string &bc_type, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Definition: bc_builder.C:106
bool is_dirichlet_bc_type(const std::string &bc_type)
Definition: bc_builder.C:139
bool is_var_active(const FEVariablesBase &var, const std::vector< libMesh::subdomain_id_type > &subdomain_ids) const
Check if the Variable var is active on the given subdomain_id.
static FEVariablesBase & get_variable(const std::string &var_name)
static std::string bc_section()
Outer section name for boundary conditions section in input file.
bool is_neumann_bc_type(const std::string &bc_type)
Definition: bc_builder.C:145
void GRINS::DefaultBCBuilder::build_periodic_bc ( const GetPot &  input,
libMesh::System &  system,
const std::set< BoundaryID > &  bc_ids,
const std::string &  section 
)
protected

Definition at line 384 of file default_bc_builder.C.

References GRINS::BCBuilder::add_periodic_bc_to_dofmap(), parse_periodic_master_slave_ids(), and parse_periodic_offset().

Referenced by build_type_based_bcs().

388  {
389  // Make sure we're not using a ParallelMesh
390  // https://github.com/libMesh/libmesh/issues/977
391  const libMesh::MeshBase& mesh = system.get_mesh();
392  const libMesh::ParallelMesh* pmesh = dynamic_cast<const libMesh::ParallelMesh*>(&mesh);
393  if(pmesh)
394  {
395  std::stringstream error_msg;
396  error_msg << "ERROR: Cannot use ParallelMesh with periodic boundary conditions!"
397  << std::endl
398  << " See https://github.com/libMesh/libmesh/issues/977 for discussion."
399  << std::endl;
400  libmesh_error_msg(error_msg.str());
401  }
402 
403  libMesh::boundary_id_type invalid_bid =
404  std::numeric_limits<libMesh::boundary_id_type>::max();
405 
406  libMesh::boundary_id_type slave_id = invalid_bid;
407  libMesh::boundary_id_type master_id = invalid_bid;
408 
409  this->parse_periodic_master_slave_ids(input,section,master_id,slave_id);
410 
411  if( bc_ids.find(slave_id) == bc_ids.end() ||
412  bc_ids.find(master_id) == bc_ids.end() )
413  libmesh_error_msg("ERROR: Mismatch between bc_ids and master/slave ids for perioid bcs!");
414 
415  libMesh::RealVectorValue offset_vector =
416  this->parse_periodic_offset(input,section);
417 
418  libMesh::DofMap& dof_map = system.get_dof_map();
419 
420  this->add_periodic_bc_to_dofmap( master_id, slave_id, offset_vector, dof_map );
421  }
void add_periodic_bc_to_dofmap(libMesh::boundary_id_type master_id, libMesh::boundary_id_type slave_id, const libMesh::RealVectorValue &offset_vector, libMesh::DofMap &dof_map)
Definition: bc_builder.C:151
void parse_periodic_master_slave_ids(const GetPot &input, const std::string &section, libMesh::boundary_id_type &master_id, libMesh::boundary_id_type &slave_id) const
libMesh::RealVectorValue parse_periodic_offset(const GetPot &input, const std::string &section) const
void GRINS::DefaultBCBuilder::build_type_based_bcs ( const GetPot &  input,
MultiphysicsSystem system,
const std::set< BoundaryID > &  bc_ids,
libMesh::DofMap &  dof_map,
const std::string &  type_input_section,
std::set< std::string > &  var_sections,
std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs 
)
protected

Helper function to build boundary conditions specified by a single type.

Examples include axisymmetric and periodic.

Definition at line 105 of file default_bc_builder.C.

References GRINS::BoundaryConditionNames::axisymmetric(), build_axisymmetric_bcs(), build_periodic_bc(), GRINS::Physics::is_axisymmetric(), and GRINS::BoundaryConditionNames::periodic().

Referenced by build_bcs().

112  {
113  std::string type_input = type_input_section+"/type";
114  std::string bc_type = input( type_input, std::string("DIE!") );
115 
116  if( bc_type == BoundaryConditionNames::axisymmetric() )
117  {
118  // Check and make sure the Physics thinks it's axisymmetric, otherwise error
119  if( !Physics::is_axisymmetric() )
120  libmesh_error_msg("ERROR: Must specify Physics/is_axisymmetric = true for axisymmetric BC!");
121 
122  // Now build the boundary condition
123  this->build_axisymmetric_bcs(input,system,bc_ids,dof_map,
124  bc_type,var_sections,neumann_bcs);
125  }
126  else if( bc_type == BoundaryConditionNames::periodic() )
127  {
128  this->build_periodic_bc(input,system,bc_ids,type_input_section);
129  }
130  else
131  {
132  std::string error_msg = "ERROR: Invalid type '"+bc_type+"' for "+type_input+"!\n";
133  error_msg += " Valid values are: "+BoundaryConditionNames::axisymmetric()+"\n";
134  error_msg += " "+BoundaryConditionNames::periodic()+"\n";
135  error_msg += " If your boundary condition is not one of these types, then \n";
136  error_msg += " you must specify the boundary condition type for each Variable\n";
137  error_msg += " section. Please have a look at the examples.\n";
138  libmesh_error_msg(error_msg);
139  }
140  }
static bool is_axisymmetric()
Definition: physics.h:132
void build_periodic_bc(const GetPot &input, libMesh::System &system, const std::set< BoundaryID > &bc_ids, const std::string &section)
void build_axisymmetric_bcs(const GetPot &input, MultiphysicsSystem &system, const std::set< BoundaryID > &bc_ids, libMesh::DofMap &dof_map, const std::string &bc_type, std::set< std::string > &var_sections, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
bool GRINS::DefaultBCBuilder::is_var_active ( const FEVariablesBase var,
const std::vector< libMesh::subdomain_id_type > &  subdomain_ids 
) const
protected

Check if the Variable var is active on the given subdomain_id.

Definition at line 511 of file default_bc_builder.C.

References GRINS::FEVariablesBase::subdomain_ids().

Referenced by build_bcs_by_var_section().

513  {
514  bool var_active = false;
515 
516  // If the var subdomain_ids are empty, then this var is active
517  // on the whole domain and we don't need to do anything
518  if( var.subdomain_ids().empty() )
519  var_active = true;
520 
521  else
522  {
523  // Now check if this variable is enabled on this subdomain
524  const std::set<libMesh::subdomain_id_type>& var_subdomain_ids =
525  var.subdomain_ids();
526 
527  for(std::vector<libMesh::subdomain_id_type>::const_iterator id = subdomain_ids.begin(); id < subdomain_ids.end(); ++id )
528  {
529  if( var_subdomain_ids.find(*id) != var_subdomain_ids.end() )
530  var_active = true;
531  }
532  }
533 
534  return var_active;
535  }
void GRINS::DefaultBCBuilder::parse_and_build_bc_id_map ( const GetPot &  input,
std::map< std::string, std::set< BoundaryID > > &  bc_id_map 
)
protected

Definition at line 291 of file default_bc_builder.C.

References GRINS::BoundaryConditionNames::bc_id_name_map_var(), GRINS::BoundaryConditionNames::bc_ids_var(), and GRINS::StringUtilities::split_string().

Referenced by build_bcs(), GRINSTesting::DefaultBCBuilderTest::test_parse_and_build_bc_id_map(), and GRINSTesting::DefaultBCBuilderTest::test_verify_bc_ids_with_mesh().

293  {
294  // First, make sure the proper input variables have been set
295  if( !input.have_variable(BoundaryConditionNames::bc_ids_var()) ||
296  !input.have_variable(BoundaryConditionNames::bc_id_name_map_var()) )
297  {
298  std::string error_msg = "ERROR: Must specify both "+BoundaryConditionNames::bc_ids_var();
299  error_msg += " and "+BoundaryConditionNames::bc_id_name_map_var()+"!";
300  libmesh_error_msg(error_msg);
301  }
302 
303  // Make sure the vectors are the same size
304  unsigned int n_ids = input.vector_variable_size(BoundaryConditionNames::bc_ids_var());
305  unsigned int n_names = input.vector_variable_size(BoundaryConditionNames::bc_id_name_map_var());
306 
307  if( n_ids != n_names )
308  libmesh_error_msg("ERROR: Must have matching number of boundary id sets and boundary names!");
309 
310  // We'll build this up along the way and then double check
311  // that we have a one-to-one match with the boundary ids
312  // that the mesh knows about
313  std::set<BoundaryID> all_bc_ids;
314 
315  // Now build up the map
316  for( unsigned int i = 0; i < n_names; i++ )
317  {
318  std::string bc_name = input(BoundaryConditionNames::bc_id_name_map_var(),std::string("DIE!"),i);
319  std::string bc_ids_input = input(BoundaryConditionNames::bc_ids_var(),std::string("DIE!"),i);
320 
321  // Users can group multiple bc_ids together with the ':' delimiter
322  std::vector<std::string> bc_ids_str;
323  StringUtilities::split_string( bc_ids_input, ":", bc_ids_str );
324  std::set<BoundaryID> bc_ids;
325 
326  for(std::vector<std::string>::const_iterator it = bc_ids_str.begin();
327  it < bc_ids_str.end(); ++it )
328  {
329  BoundaryID id = StringUtilities::string_to_T<BoundaryID>(*it);
330 
331  // Can only specify a boundary id once
332  if( (bc_ids.find(id) != bc_ids.end()) ||
333  (all_bc_ids.find(id) != all_bc_ids.end()) )
334  libmesh_error_msg("ERROR: Can only specify a boundary ID once!");
335 
336  bc_ids.insert(id);
337  all_bc_ids.insert(id);
338  }
339 
340  // Insert pair into the bc_id_map
341  bc_id_map.insert( std::make_pair( bc_name, bc_ids ) );
342  }
343  }
libMesh::boundary_id_type BoundaryID
More descriptive name of the type used for boundary ids.
Definition: var_typedefs.h:56
void split_string(const std::string &input, const std::string &delimiter, std::vector< std::string > &results)
Definition: string_utils.C:31
static std::string bc_id_name_map_var()
Names of boundaries to correspond with bc_ids.
static std::string bc_ids_var()
Variable for list boundary ids to correspond with bc_id_name_map.
void GRINS::DefaultBCBuilder::parse_periodic_master_slave_ids ( const GetPot &  input,
const std::string &  section,
libMesh::boundary_id_type &  master_id,
libMesh::boundary_id_type &  slave_id 
) const
protected

Definition at line 424 of file default_bc_builder.C.

Referenced by build_periodic_bc(), and GRINSTesting::DefaultBCBuilderTest::test_parse_periodic_master_slave_ids().

427  {
428  if( !input.have_variable(section+"/master_id") )
429  libmesh_error_msg("ERROR: Could not find master_id in section "+section);
430 
431  if( !input.have_variable(section+"/slave_id") )
432  libmesh_error_msg("ERROR: Could not find slave_id in section "+section);
433 
434  libMesh::boundary_id_type invalid_bid =
435  std::numeric_limits<libMesh::boundary_id_type>::max();
436 
437  master_id = input(section+"/master_id", invalid_bid);
438  slave_id = input(section+"/slave_id", invalid_bid);
439  }
libMesh::RealVectorValue GRINS::DefaultBCBuilder::parse_periodic_offset ( const GetPot &  input,
const std::string &  section 
) const
protected

Definition at line 442 of file default_bc_builder.C.

Referenced by build_periodic_bc(), and GRINSTesting::DefaultBCBuilderTest::test_parse_periodic_offset().

443  {
444  std::string input_section = section+"/boundary_offset";
445  if( !input.have_variable(input_section) )
446  libmesh_error_msg("ERROR: Could not find boundary_offset in section "+section);
447 
448  unsigned int n_comps = input.vector_variable_size(input_section);
449 
450  if( n_comps > 3 )
451  libmesh_error_msg("ERROR: Cannot specify more than 3 components for boundary_offset!");
452 
453  libMesh::RealVectorValue offset;
454  libMesh::Real invalid_real = std::numeric_limits<libMesh::Real>::max();
455 
456  for( unsigned int i = 0; i < n_comps; i++ )
457  offset(i) = input(input_section, invalid_real, i );
458 
459  return offset;
460  }
void GRINS::DefaultBCBuilder::verify_bc_ids_with_mesh ( const MultiphysicsSystem system,
const std::map< std::string, std::set< BoundaryID > > &  bc_id_map 
) const
protected

Definition at line 345 of file default_bc_builder.C.

Referenced by build_bcs(), and GRINSTesting::DefaultBCBuilderTest::test_verify_bc_ids_with_mesh().

347  {
348  const libMesh::MeshBase& mesh = system.get_mesh();
349  const libMesh::BoundaryInfo& boundary_info = mesh.get_boundary_info();
350 
351  std::set<BoundaryID> mesh_ids = boundary_info.get_boundary_ids();
352  mesh.comm().set_union(mesh_ids);
353 
354  // Collect all the bc_ids into one set so we can just compare the sets
355  std::set<BoundaryID> all_ids;
356 
357  for( std::map<std::string,std::set<BoundaryID> >::const_iterator bc_it = bc_id_map.begin();
358  bc_it != bc_id_map.end(); ++bc_it )
359  all_ids.insert( (bc_it->second).begin(), (bc_it->second).end() );
360 
361  if( mesh_ids != all_ids )
362  {
363  std::string err_msg = "ERROR: Mismatch between specified boundary ids and the boundary ids in the mesh!\n";
364  err_msg += "User specified ids: ";
365 
366  for( std::set<BoundaryID>::const_iterator it = all_ids.begin();
367  it != all_ids.end(); ++it )
368  err_msg += StringUtilities::T_to_string<BoundaryID>(*it)+" ";
369 
370  err_msg += "\n";
371 
372  err_msg += "Mesh specified ids: ";
373 
374  for( std::set<BoundaryID>::const_iterator it = mesh_ids.begin();
375  it != mesh_ids.end(); ++it )
376  err_msg += StringUtilities::T_to_string<BoundaryID>(*it)+" ";
377 
378  err_msg += "\n";
379 
380  libmesh_error_msg(err_msg);
381  }
382  }

The documentation for this class was generated from the following files:

Generated on Tue Dec 19 2017 12:47:30 for GRINS-0.8.0 by  doxygen 1.8.9.1