GRINS-0.7.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 ()
 

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, 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 parse_var_sections (const GetPot &input, std::set< std::string > &sections)
 
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
 
- 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 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 46 of file default_bc_builder.h.

Constructor & Destructor Documentation

GRINS::DefaultBCBuilder::DefaultBCBuilder ( )
inline

Definition at line 50 of file default_bc_builder.h.

51  : BCBuilder()
52  {};
GRINS::DefaultBCBuilder::~DefaultBCBuilder ( )
inline

Definition at line 54 of file default_bc_builder.h.

54 {};

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 130 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().

137  {
138  // Axisymmetric depends on the variable type.
139  // So, we loop over all the Variable names and prepend the type
140  // with the variable name. So, the factories for axisymmetric
141  // need to be instantiated with, e.g., Velocity_axisymmetric
142  for( std::set<std::string>::const_iterator vars = var_sections.begin();
143  vars != var_sections.end(); ++vars )
144  {
145  const std::string& var_section = *vars;
146 
147  // Grab FEVariable
148  const FEVariablesBase& fe_var =
150 
151  // We don't need to do anything for constraint variables
152  if( fe_var.is_constraint_var() )
153  continue;
154 
155  std::string full_bc_type = var_section+"_"+bc_type;
156 
157  // Axisymmetric is Dirichlet or Neumann, depending on the Variable
158  if( this->is_dirichlet_bc_type(full_bc_type) )
159  {
160  this->construct_dbc_core( input,
161  system,
162  bc_ids,
163  fe_var,
164  std::string("dummy"), // No section needed for axisymmetric
165  full_bc_type,
166  dof_map );
167  }
168  else if( this->is_neumann_bc_type(full_bc_type) )
169  {
170  this->construct_nbc_core( input,
171  system,
172  bc_ids,
173  fe_var,
174  std::string("dummy"), // No section needed for axisymmetric
175  full_bc_type,
176  neumann_bcs );
177  }
178  else
179  libmesh_error_msg("ERROR! Invalid axisymmetric type: "+full_bc_type);
180  }
181  }
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 const 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_bcs ( const GetPot &  input,
MultiphysicsSystem system,
std::vector< SharedPtr< NeumannBCContainer > > &  neumann_bcs 
)
protectedvirtual

Implements GRINS::BCBuilder.

Definition at line 45 of file default_bc_builder.C.

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

47  {
48  libMesh::DofMap& dof_map = system.get_dof_map();
49 
50  // Setup map between boundary section name and the boundary ids
51  std::map<std::string,std::set<BoundaryID> > bc_id_map;
52  this->parse_and_build_bc_id_map( input, bc_id_map );
53 
54  // The user must have mapped all boundary ids present in the mesh.
55  // Additionally, they can't have specified an boundary id that doesn't
56  // exist.
57  this->verify_bc_ids_with_mesh( system, bc_id_map );
58 
59  // Parse variable sections. These correspond to the names of FEVariables
60  // We will be looking for all variable sections in all boundary sections
61  // that are not "constraint" variables (i.e. FEVariablesBase::_is_constraint)
62  std::set<std::string> var_sections;
63  this->parse_var_sections( input, var_sections );
64 
65  for( std::map<std::string,std::set<BoundaryID> >::const_iterator bc_it = bc_id_map.begin();
66  bc_it != bc_id_map.end(); ++bc_it )
67  {
68  const std::string& bc_name = bc_it->first;
69 
70  const std::set<BoundaryID>& bc_ids = bc_it->second;
71 
72  // First check for special types of boundary conditions that can be
73  // specified as a single type that applies to all boundaries
74  std::string type_input_section = BoundaryConditionNames::bc_section()+"/"+bc_name;
75  std::string type_input = type_input_section+"/type";
76 
77  if( input.have_variable(type_input) )
78  {
79  this->build_type_based_bcs(input,system,bc_ids,dof_map,
80  type_input_section,var_sections,
81  neumann_bcs);
82  }
83  // Otherwise, we look for each Variable section and the
84  // boundary conditions specified therein
85  else
86  {
87  this->build_bcs_by_var_section(input,system,bc_name,bc_ids,dof_map,
88  var_sections,neumann_bcs);
89  }
90  }
91  }
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 parse_var_sections(const GetPot &input, std::set< std::string > &sections)
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, std::vector< SharedPtr< NeumannBCContainer > > &neumann_bcs)
Helper function to build boundary conditions using Variable sections.
static std::string bc_section()
Outer section name for boundary conditions section in input file.
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,
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 183 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(), and GRINS::BCBuilder::is_neumann_bc_type().

Referenced by build_bcs().

190  {
191  for( std::set<std::string>::const_iterator vars = var_sections.begin();
192  vars != var_sections.end(); ++vars )
193  {
194  const std::string& var_section = *vars;
195 
196  // Setup the input section we're in. We use this, but
197  // this is also the section where the BCFactory will try
198  // to parse the values, if it needs to, so let's set up
199  // once and reuse it.
200  std::string input_section = std::string(BoundaryConditionNames::bc_section()+"/"+bc_name+"/"+var_section);
201 
202  // Make sure this section is there, unless that variable is a constraint variable
203  if( !input.have_section(input_section) )
204  {
205  const FEVariablesBase& var =
207  if( var.is_constraint_var() )
208  continue;
209  else
210  libmesh_error_msg("ERROR: Could not find boundary condition specification for "+input_section+"!");
211  }
212 
213  // Grab the FEVariable
214  const FEVariablesBase& fe_var =
216 
217  // Grab the type of the boundary condition
218  // There may be more than one type (e.g. pin displacement in two directions and
219  // and load in the third direction).
220  std::string bc_type_section = input_section+"/type";
221  unsigned int n_bc_types = input.vector_variable_size(bc_type_section);
222 
223  for( unsigned int bc_type_idx = 0; bc_type_idx < n_bc_types; bc_type_idx++ )
224  {
225  std::string bc_type = input(bc_type_section, std::string("DIE!"), bc_type_idx);
226 
227  if( this->is_dirichlet_bc_type(bc_type) )
228  {
229  this->construct_dbc_core( input,system, bc_ids,
230  fe_var, input_section,
231  bc_type, dof_map );
232  }
233  else if( this->is_neumann_bc_type(bc_type) )
234  {
235  this->construct_nbc_core( input,system, bc_ids,
236  fe_var, input_section,
237  bc_type, neumann_bcs );
238  }
239  else
240  libmesh_error_msg("ERROR: Invalid bc_type "+bc_type+"!");
241 
242  } // end loop over bc_types
243 
244  } // end loop over variable sections
245  }
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 const 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 376 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().

380  {
381  // Make sure we're not using a ParallelMesh
382  // https://github.com/libMesh/libmesh/issues/977
383  const libMesh::MeshBase& mesh = system.get_mesh();
384  const libMesh::ParallelMesh* pmesh = dynamic_cast<const libMesh::ParallelMesh*>(&mesh);
385  if(pmesh)
386  {
387  std::stringstream error_msg;
388  error_msg << "ERROR: Cannot use ParallelMesh with periodic boundary conditions!"
389  << std::endl
390  << " See https://github.com/libMesh/libmesh/issues/977 for discussion."
391  << std::endl;
392  libmesh_error_msg(error_msg.str());
393  }
394 
395  libMesh::boundary_id_type invalid_bid =
396  std::numeric_limits<libMesh::boundary_id_type>::max();
397 
398  libMesh::boundary_id_type slave_id = invalid_bid;
399  libMesh::boundary_id_type master_id = invalid_bid;
400 
401  this->parse_periodic_master_slave_ids(input,section,master_id,slave_id);
402 
403  if( bc_ids.find(slave_id) == bc_ids.end() ||
404  bc_ids.find(master_id) == bc_ids.end() )
405  libmesh_error_msg("ERROR: Mismatch between bc_ids and master/slave ids for perioid bcs!");
406 
407  libMesh::RealVectorValue offset_vector =
408  this->parse_periodic_offset(input,section);
409 
410  libMesh::DofMap& dof_map = system.get_dof_map();
411 
412  this->add_periodic_bc_to_dofmap( master_id, slave_id, offset_vector, dof_map );
413  }
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 93 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().

100  {
101  std::string type_input = type_input_section+"/type";
102  std::string bc_type = input( type_input, std::string("DIE!") );
103 
104  if( bc_type == BoundaryConditionNames::axisymmetric() )
105  {
106  // Check and make sure the Physics thinks it's axisymmetric, otherwise error
107  if( !Physics::is_axisymmetric() )
108  libmesh_error_msg("ERROR: Must specify Physics/is_axisymmetric = true for axisymmetric BC!");
109 
110  // Now build the boundary condition
111  this->build_axisymmetric_bcs(input,system,bc_ids,dof_map,
112  bc_type,var_sections,neumann_bcs);
113  }
114  else if( bc_type == BoundaryConditionNames::periodic() )
115  {
116  this->build_periodic_bc(input,system,bc_ids,type_input_section);
117  }
118  else
119  {
120  std::string error_msg = "ERROR: Invalid type '"+bc_type+"' for "+type_input+"!\n";
121  error_msg += " Valid values are: "+BoundaryConditionNames::axisymmetric()+"\n";
122  error_msg += " "+BoundaryConditionNames::periodic()+"\n";
123  error_msg += " If your boundary condition is not one of these types, then \n";
124  error_msg += " you must specify the boundary condition type for each Variable\n";
125  error_msg += " section. Please have a look at the examples.\n";
126  libmesh_error_msg(error_msg);
127  }
128  }
static bool is_axisymmetric()
Definition: physics.h:133
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)
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 247 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().

249  {
250  // First, make sure the proper input variables have been set
251  if( !input.have_variable(BoundaryConditionNames::bc_ids_var()) ||
252  !input.have_variable(BoundaryConditionNames::bc_id_name_map_var()) )
253  {
254  std::string error_msg = "ERROR: Must specify both "+BoundaryConditionNames::bc_ids_var();
255  error_msg += " and "+BoundaryConditionNames::bc_id_name_map_var()+"!";
256  libmesh_error_msg(error_msg);
257  }
258 
259  // Make sure the vectors are the same size
260  unsigned int n_ids = input.vector_variable_size(BoundaryConditionNames::bc_ids_var());
261  unsigned int n_names = input.vector_variable_size(BoundaryConditionNames::bc_id_name_map_var());
262 
263  if( n_ids != n_names )
264  libmesh_error_msg("ERROR: Must have matching number of boundary id sets and boundary names!");
265 
266  // We'll build this up along the way and then double check
267  // that we have a one-to-one match with the boundary ids
268  // that the mesh knows about
269  std::set<BoundaryID> all_bc_ids;
270 
271  // Now build up the map
272  for( unsigned int i = 0; i < n_names; i++ )
273  {
274  std::string bc_name = input(BoundaryConditionNames::bc_id_name_map_var(),std::string("DIE!"),i);
275  std::string bc_ids_input = input(BoundaryConditionNames::bc_ids_var(),std::string("DIE!"),i);
276 
277  // Users can group multiple bc_ids together with the ':' delimiter
278  std::vector<std::string> bc_ids_str;
279  StringUtilities::split_string( bc_ids_input, ":", bc_ids_str );
280  std::set<BoundaryID> bc_ids;
281 
282  for(std::vector<std::string>::const_iterator it = bc_ids_str.begin();
283  it < bc_ids_str.end(); ++it )
284  {
285  BoundaryID id = StringUtilities::string_to_T<BoundaryID>(*it);
286 
287  // Can only specify a boundary id once
288  if( (bc_ids.find(id) != bc_ids.end()) ||
289  (all_bc_ids.find(id) != all_bc_ids.end()) )
290  libmesh_error_msg("ERROR: Can only specify a boundary ID once!");
291 
292  bc_ids.insert(id);
293  all_bc_ids.insert(id);
294  }
295 
296  // Insert pair into the bc_id_map
297  bc_id_map.insert( std::make_pair( bc_name, bc_ids ) );
298  }
299  }
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 416 of file default_bc_builder.C.

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

419  {
420  if( !input.have_variable(section+"/master_id") )
421  libmesh_error_msg("ERROR: Could not find master_id in section "+section);
422 
423  if( !input.have_variable(section+"/slave_id") )
424  libmesh_error_msg("ERROR: Could not find slave_id in section "+section);
425 
426  libMesh::boundary_id_type invalid_bid =
427  std::numeric_limits<libMesh::boundary_id_type>::max();
428 
429  master_id = input(section+"/master_id", invalid_bid);
430  slave_id = input(section+"/slave_id", invalid_bid);
431  }
libMesh::RealVectorValue GRINS::DefaultBCBuilder::parse_periodic_offset ( const GetPot &  input,
const std::string &  section 
) const
protected

Definition at line 434 of file default_bc_builder.C.

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

435  {
436  std::string input_section = section+"/boundary_offset";
437  if( !input.have_variable(input_section) )
438  libmesh_error_msg("ERROR: Could not find boundary_offset in section "+section);
439 
440  unsigned int n_comps = input.vector_variable_size(input_section);
441 
442  if( n_comps > 3 )
443  libmesh_error_msg("ERROR: Cannot specify more than 3 components for boundary_offset!");
444 
445  libMesh::RealVectorValue offset;
446  libMesh::Real invalid_real = std::numeric_limits<libMesh::Real>::max();
447 
448  for( unsigned int i = 0; i < n_comps; i++ )
449  offset(i) = input(input_section, invalid_real, i );
450 
451  return offset;
452  }
void GRINS::DefaultBCBuilder::parse_var_sections ( const GetPot &  input,
std::set< std::string > &  sections 
)
protected
Todo:
This would probably be a good function to add to libMesh::GetPot

Definition at line 339 of file default_bc_builder.C.

References GRINS::StringUtilities::split_string(), and GRINS::VariablesParsing::variables_section().

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

341  {
342  if( !input.have_section(VariablesParsing::variables_section()) )
343  libmesh_error_msg("ERROR: Could not find "+VariablesParsing::variables_section()+" section!");
344 
345  // We need to extract all the Variable sections from the input file
346  // We'll populate the relevant sections in var_sections
348  std::vector<std::string> all_sections = input.get_section_names();
349  for( std::vector<std::string>::const_iterator s = all_sections.begin();
350  s < all_sections.end(); ++s )
351  {
352  // First check that it contains "Variable" as the first slot
353  if( s->find(VariablesParsing::variables_section()) == 0 )
354  {
355  // Now check it only has 2 elements when we split on "/"
356  std::vector<std::string> split_str;
357  StringUtilities::split_string(*s, "/", split_str );
358 
359  // Our Variable should be the second part of the split
360  if( split_str.size() == 2 )
361  {
362  // Make sure we don't already have that section
363  if( sections.find(split_str[1]) != sections.end() )
364  libmesh_error_msg("ERROR: Found duplicate Variable section "+split_str[1]+"!");
365 
366  sections.insert( split_str[1] );
367  }
368  }
369  }
370 
371  // Make sure we found some variable subsections
372  if( sections.empty() )
373  libmesh_error_msg("ERROR: Did not find any Variable subsections!");
374  }
static std::string variables_section()
Helper function to encapsualte the overall [Variables] section name.
void split_string(const std::string &input, const std::string &delimiter, std::vector< std::string > &results)
Definition: string_utils.C:31
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 301 of file default_bc_builder.C.

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

303  {
304  const libMesh::MeshBase& mesh = system.get_mesh();
305  const libMesh::BoundaryInfo& boundary_info = mesh.get_boundary_info();
306 
307  const std::set<BoundaryID> mesh_ids = boundary_info.get_boundary_ids();
308 
309  // Collect all the bc_ids into one set so we can just compare the sets
310  std::set<BoundaryID> all_ids;
311 
312  for( std::map<std::string,std::set<BoundaryID> >::const_iterator bc_it = bc_id_map.begin();
313  bc_it != bc_id_map.end(); ++bc_it )
314  all_ids.insert( (bc_it->second).begin(), (bc_it->second).end() );
315 
316  if( mesh_ids != all_ids )
317  {
318  std::string err_msg = "ERROR: Mismatch between specified boundary ids and the boundary ids in the mesh!\n";
319  err_msg += "User specified ids: ";
320 
321  for( std::set<BoundaryID>::const_iterator it = all_ids.begin();
322  it != all_ids.end(); ++it )
323  err_msg += StringUtilities::T_to_string<BoundaryID>(*it)+" ";
324 
325  err_msg += "\n";
326 
327  err_msg += "Mesh specified ids: ";
328 
329  for( std::set<BoundaryID>::const_iterator it = mesh_ids.begin();
330  it != mesh_ids.end(); ++it )
331  err_msg += StringUtilities::T_to_string<BoundaryID>(*it)+" ";
332 
333  err_msg += "\n";
334 
335  libmesh_error_msg(err_msg);
336  }
337  }

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

Generated on Thu Jun 2 2016 21:52:30 for GRINS-0.7.0 by  doxygen 1.8.10