GRINS-0.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
GRINS::Simulation Class Reference

#include <simulation.h>

Collaboration diagram for GRINS::Simulation:
Collaboration graph
[legend]

Public Member Functions

 Simulation (const GetPot &input, SimulationBuilder &sim_builder, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 Simulation (const GetPot &input, GetPot &command_line, SimulationBuilder &sim_builder, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
virtual ~Simulation ()
 
void run ()
 
void print_sim_info ()
 
std::tr1::shared_ptr< libMesh::EquationSystems > get_equation_system ()
 
libMesh::Number get_qoi_value (unsigned int qoi_index) const
 
const std::string & get_multiphysics_system_name () const
 

Protected Member Functions

void read_restart (const GetPot &input)
 
void attach_neumann_bc_funcs (std::map< GRINS::PhysicsName, GRINS::NBCContainer > neumann_bcs, GRINS::MultiphysicsSystem *system)
 
void attach_dirichlet_bc_funcs (std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > dbc_map, GRINS::MultiphysicsSystem *system)
 
void init_multiphysics_system (const GetPot &input, SimulationBuilder &sim_builder)
 Helper function. More...
 
void init_qois (const GetPot &input, SimulationBuilder &sim_builder)
 Helper function. More...
 
void init_params (const GetPot &input, SimulationBuilder &sim_builder)
 Helper function. More...
 
void init_restart (const GetPot &input, SimulationBuilder &sim_builder, const libMesh::Parallel::Communicator &comm)
 Helper function. More...
 
void check_for_unused_vars (const GetPot &input, bool warning_only)
 Helper function. More...
 
bool check_for_adjoint_solve (const GetPot &input) const
 Helper function. More...
 
void init_adjoint_solve (const GetPot &input, bool output_adjoint)
 Helper function. More...
 

Protected Attributes

std::tr1::shared_ptr< libMesh::UnstructuredMesh > _mesh
 
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
 
std::tr1::shared_ptr< GRINS::Solver_solver
 
std::string _system_name
 GRINS::Multiphysics system name. More...
 
GRINS::MultiphysicsSystem_multiphysics_system
 
std::tr1::shared_ptr< GRINS::Visualization_vis
 
std::tr1::shared_ptr< PostProcessedQuantities< libMesh::Real > > _postprocessing
 
bool _print_mesh_info
 
bool _print_log_info
 
bool _print_equation_system_info
 
bool _print_perflog
 
bool _print_qoi
 
bool _print_scalars
 
bool _output_vis
 
bool _output_adjoint
 
bool _output_residual
 
bool _output_residual_sensitivities
 
bool _output_solution_sensitivities
 
unsigned int _timesteps_per_vis
 
unsigned int _timesteps_per_perflog
 
std::tr1::shared_ptr< libMesh::ErrorEstimator > _error_estimator
 
ParameterManager _adjoint_parameters
 
ParameterManager _forward_parameters
 
bool _do_adjoint_solve
 

Private Member Functions

 Simulation ()
 

Detailed Description

Definition at line 62 of file simulation.h.

Constructor & Destructor Documentation

GRINS::Simulation::Simulation ( const GetPot &  input,
SimulationBuilder sim_builder,
const libMesh::Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
)
GRINS::Simulation::Simulation ( const GetPot &  input,
GetPot &  command_line,
SimulationBuilder sim_builder,
const libMesh::Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
)
GRINS::Simulation::~Simulation ( )
virtual

Definition at line 140 of file simulation.C.

141  {
142  return;
143  }
GRINS::Simulation::Simulation ( )
private

Member Function Documentation

void GRINS::Simulation::attach_dirichlet_bc_funcs ( std::multimap< GRINS::PhysicsName, GRINS::DBCContainer dbc_map,
GRINS::MultiphysicsSystem system 
)
protected

Definition at line 492 of file simulation.C.

References GRINS::MultiphysicsSystem::get_physics().

Referenced by init_multiphysics_system().

494  {
495  for( std::multimap< PhysicsName, DBCContainer >::const_iterator it = dbc_map.begin();
496  it != dbc_map.end();
497  it++ )
498  {
499  std::tr1::shared_ptr<Physics> physics = system->get_physics( it->first );
500 
501  physics->attach_dirichlet_bound_func( it->second );
502  }
503  return;
504  }
std::tr1::shared_ptr< GRINS::Physics > get_physics(const std::string physics_name)
void GRINS::Simulation::attach_neumann_bc_funcs ( std::map< GRINS::PhysicsName, GRINS::NBCContainer neumann_bcs,
GRINS::MultiphysicsSystem system 
)
protected

Definition at line 473 of file simulation.C.

References GRINS::MultiphysicsSystem::get_physics().

Referenced by init_multiphysics_system().

475  {
476  //_neumann_bc_funcs = neumann_bcs;
477 
478  if( neumann_bcs.size() > 0 )
479  {
480  for( std::map< std::string, NBCContainer >::iterator bc = neumann_bcs.begin();
481  bc != neumann_bcs.end();
482  bc++ )
483  {
484  std::tr1::shared_ptr<Physics> physics = system->get_physics( bc->first );
485  physics->attach_neumann_bound_func( bc->second );
486  }
487  }
488 
489  return;
490  }
std::tr1::shared_ptr< GRINS::Physics > get_physics(const std::string physics_name)
bool GRINS::Simulation::check_for_adjoint_solve ( const GetPot &  input) const
protected

Helper function.

Todo:
We need to collect these options into one spot
Todo:
This option name WILL change at some point.

Definition at line 529 of file simulation.C.

Referenced by init_adjoint_solve().

530  {
532  std::string error_estimator = input("MeshAdaptivity/estimator_type", "none");
533 
534  bool do_adjoint_solve = false;
535 
536  // First check if the error estimator requires an adjoint solve
537  if( error_estimator.find("adjoint") != std::string::npos )
538  {
539  do_adjoint_solve = true;
540  }
541 
542  // Now check if the user requested to do an adjoint solve regardless
544  if( input( "linear-nonlinear-solver/do_adjoint_solve", false ) )
545  {
546  do_adjoint_solve = true;
547  }
548 
549  return do_adjoint_solve;
550  }
void GRINS::Simulation::check_for_unused_vars ( const GetPot &  input,
bool  warning_only 
)
protected

Helper function.

Definition at line 273 of file simulation.C.

274  {
275  /* Everything should be set up now, so check if there's any unused variables
276  in the input file. If so, then tell the user what they were and error out. */
277  std::vector<std::string> unused_vars = input.unidentified_variables();
278 
279  if( !unused_vars.empty() )
280  {
281  libMesh::err << "==========================================================" << std::endl;
282  if( warning_only )
283  libMesh::err << "Warning: ";
284  else
285  libMesh::err << "Error: ";
286 
287  libMesh::err << "Found unused variables!" << std::endl;
288 
289  for( std::vector<std::string>::const_iterator it = unused_vars.begin();
290  it != unused_vars.end(); ++it )
291  {
292  libMesh::err << *it << std::endl;
293  }
294  libMesh::err << "==========================================================" << std::endl;
295 
296  if( !warning_only )
297  libmesh_error();
298  }
299 
300  return;
301  }
std::tr1::shared_ptr< libMesh::EquationSystems > GRINS::Simulation::get_equation_system ( )

Definition at line 422 of file simulation.C.

References _equation_system.

Referenced by main().

423  {
424  return _equation_system;
425  }
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
const std::string & GRINS::Simulation::get_multiphysics_system_name ( ) const
inline

Definition at line 176 of file simulation.h.

References _system_name.

177  {
178  return this->_system_name;
179  }
std::string _system_name
GRINS::Multiphysics system name.
Definition: simulation.h:133
libMesh::Number GRINS::Simulation::get_qoi_value ( unsigned int  qoi_index) const

Definition at line 427 of file simulation.C.

References _multiphysics_system, and GRINS::CompositeQoI::get_qoi_value().

428  {
429  const CompositeQoI* qoi = libMesh::libmesh_cast_ptr<const CompositeQoI*>(this->_multiphysics_system->get_qoi());
430  return qoi->get_qoi_value(qoi_index);
431  }
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
void GRINS::Simulation::init_adjoint_solve ( const GetPot &  input,
bool  output_adjoint 
)
protected

Helper function.

Definition at line 506 of file simulation.C.

References _do_adjoint_solve, _multiphysics_system, and check_for_adjoint_solve().

507  {
508  // Check if we're doing an adjoint solve
510 
511  const libMesh::DifferentiableQoI* raw_qoi = _multiphysics_system->get_qoi();
512  const CompositeQoI* qoi = dynamic_cast<const CompositeQoI*>( raw_qoi );
513 
514  // If we are trying to do an adjoint solve without a QoI, that's an error
515  // If there are no QoIs, the CompositeQoI never gets built and qoi will be NULL
516  if( _do_adjoint_solve && !qoi )
517  {
518  libmesh_error_msg("Error: Adjoint solve requested, but no QoIs detected.");
519  }
520 
521  /* If we are not doing an adjoint solve, but adjoint output was requested:
522  that's an error. */
523  if( !_do_adjoint_solve && output_adjoint )
524  {
525  libmesh_error_msg("Error: Adjoint output requested, but no adjoint solve requested.");
526  }
527  }
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
bool check_for_adjoint_solve(const GetPot &input) const
Helper function.
Definition: simulation.C:529
void GRINS::Simulation::init_multiphysics_system ( const GetPot &  input,
SimulationBuilder sim_builder 
)
protected

Helper function.

Definition at line 145 of file simulation.C.

References _equation_system, _multiphysics_system, _postprocessing, _print_log_info, _solver, attach_dirichlet_bc_funcs(), attach_neumann_bc_funcs(), GRINS::MultiphysicsSystem::attach_physics_list(), GRINS::SimulationBuilder::build_dirichlet_bcs(), GRINS::SimulationBuilder::build_neumann_bcs(), GRINS::SimulationBuilder::build_physics(), GRINS::MultiphysicsSystem::read_input_options(), and GRINS::MultiphysicsSystem::register_postprocessing_vars().

147  {
148  // Only print libMesh logging info if the user requests it
149  libMesh::perflog.disable_logging();
150  if( this->_print_log_info ) libMesh::perflog.enable_logging();
151 
152  PhysicsList physics_list = sim_builder.build_physics(input);
153 
154  _multiphysics_system->attach_physics_list( physics_list );
155 
157 
159 
160  // This *must* be done before equation_system->init
161  this->attach_dirichlet_bc_funcs( sim_builder.build_dirichlet_bcs(), _multiphysics_system );
162 
163  /* Postprocessing needs to be initialized before the solver since that's
164  where equation_system gets init'ed */
166 
167  _solver->initialize( input, _equation_system, _multiphysics_system );
168 
169  // Useful for debugging
170  if( input("screen-options/print_dof_constraints", false ) )
171  {
172  _multiphysics_system->get_dof_map().print_dof_constraints();
173  }
174 
175  // This *must* be done after equation_system->init in order to get variable indices
176  this->attach_neumann_bc_funcs( sim_builder.build_neumann_bcs( *_equation_system ), _multiphysics_system );
177 
178  return;
179  }
void attach_neumann_bc_funcs(std::map< GRINS::PhysicsName, GRINS::NBCContainer > neumann_bcs, GRINS::MultiphysicsSystem *system)
Definition: simulation.C:473
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
std::map< std::string, std::tr1::shared_ptr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
Definition: var_typedefs.h:57
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
std::tr1::shared_ptr< PostProcessedQuantities< libMesh::Real > > _postprocessing
Definition: simulation.h:140
void attach_dirichlet_bc_funcs(std::multimap< GRINS::PhysicsName, GRINS::DBCContainer > dbc_map, GRINS::MultiphysicsSystem *system)
Definition: simulation.C:492
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
std::tr1::shared_ptr< GRINS::Solver > _solver
Definition: simulation.h:130
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled.
void GRINS::Simulation::init_params ( const GetPot &  input,
SimulationBuilder sim_builder 
)
protected

Helper function.

Definition at line 204 of file simulation.C.

References _adjoint_parameters, _forward_parameters, _multiphysics_system, and GRINS::ParameterManager::initialize().

206  {
207  unsigned int n_adjoint_parameters =
208  input.vector_variable_size("QoI/adjoint_sensitivity_parameters");
209 
210  unsigned int n_forward_parameters =
211  input.vector_variable_size("QoI/forward_sensitivity_parameters");
212 
213  // If the user actually asks for parameter sensitivities, then we
214  // set up the parameter vectors to use.
215  if ( n_adjoint_parameters )
216  {
217  // If we're doing adjoint sensitivities, dq/dp only makes
218  // sense if we have q
219  CompositeQoI* qoi =
220  libMesh::cast_ptr<CompositeQoI*>
221  (this->_multiphysics_system->get_qoi());
222 
223  if (!qoi)
224  {
225  std::cout <<
226  "Error: adjoint_sensitivity_parameters are specified but\n"
227  << "no QoIs have been specified.\n" << std::endl;
228  libmesh_error();
229  }
230 
232  (input, "QoI/adjoint_sensitivity_parameters",
233  *this->_multiphysics_system, qoi);
234  }
235 
236  if ( n_forward_parameters )
237  {
238  // If we're doing forward sensitivities, du/dp can make
239  // sense even with no q defined
240  CompositeQoI* qoi =
241  dynamic_cast<CompositeQoI*>
242  (this->_multiphysics_system->get_qoi());
243 
244  // dynamic_cast returns NULL if our QoI isn't a CompositeQoI;
245  // i.e. if there were no QoIs that made us bother setting up
246  // the CompositeQoI object. Passing NULL tells
247  // ParameterManager not to bother asking for qoi registration
248  // of parameters.
249 
251  (input, "QoI/forward_sensitivity_parameters",
252  *this->_multiphysics_system, qoi);
253  }
254  }
ParameterManager _forward_parameters
Definition: simulation.h:164
ParameterManager _adjoint_parameters
Definition: simulation.h:162
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
virtual void initialize(const GetPot &input, const std::string &parameters_varname, GRINS::MultiphysicsSystem &system, GRINS::CompositeQoI *qoi)
void GRINS::Simulation::init_qois ( const GetPot &  input,
SimulationBuilder sim_builder 
)
protected

Helper function.

Definition at line 181 of file simulation.C.

References _multiphysics_system, _print_qoi, and GRINS::SimulationBuilder::build_qoi().

182  {
183  // If the user actually asks for a QoI, then we add it.
184  std::tr1::shared_ptr<CompositeQoI> qois = sim_builder.build_qoi( input );
185  if( qois->n_qois() > 0 )
186  {
187  // This *must* be done after equation_system->init in order to get variable indices
188  qois->init(input, *_multiphysics_system );
189 
190  /* Note that we are effectively transfering ownership of the qoi pointer because
191  it will be cloned in _multiphysics_system and all the calculations are done there. */
192  _multiphysics_system->attach_qoi( qois.get() );
193  }
194  else if (_print_qoi)
195  {
196  std::cout << "Error: print_qoi is specified but\n" <<
197  "no QoIs have been specified.\n" << std::endl;
198  libmesh_error();
199  }
200 
201  return;
202  }
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
void GRINS::Simulation::init_restart ( const GetPot &  input,
SimulationBuilder sim_builder,
const libMesh::Parallel::Communicator &  comm 
)
protected

Helper function.

Definition at line 257 of file simulation.C.

References _equation_system, _mesh, GRINS::MeshBuilder::do_mesh_refinement_from_input(), GRINS::SimulationBuilder::mesh_builder(), and read_restart().

259  {
260  this->read_restart( input );
261 
262  /* We do this here only if there's a restart file. Otherwise, this was done
263  at mesh construction time */
264  sim_builder.mesh_builder().do_mesh_refinement_from_input( input, comm, *_mesh );
265 
266  /* \todo Any way to tell if the mesh got refined so we don't unnecessarily
267  call reinit()? */
268  _equation_system->reinit();
269 
270  return;
271  }
void read_restart(const GetPot &input)
Definition: simulation.C:433
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
std::tr1::shared_ptr< libMesh::UnstructuredMesh > _mesh
Definition: simulation.h:126
void GRINS::Simulation::print_sim_info ( )

Definition at line 411 of file simulation.C.

References _equation_system, _mesh, _print_equation_system_info, and _print_mesh_info.

Referenced by run().

412  {
413  // Print mesh info if the user wants it
414  if( this->_print_mesh_info ) this->_mesh->print_info();
415 
416  // Print info if requested
417  if( this->_print_equation_system_info ) this->_equation_system->print_info();
418 
419  return;
420  }
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
bool _print_equation_system_info
Definition: simulation.h:145
std::tr1::shared_ptr< libMesh::UnstructuredMesh > _mesh
Definition: simulation.h:126
void GRINS::Simulation::read_restart ( const GetPot &  input)
protected

Definition at line 433 of file simulation.C.

References _equation_system.

Referenced by init_restart().

434  {
435  const std::string restart_file = input( "restart-options/restart_file", "none" );
436 
437  // Most of this was pulled from FIN-S
438  if (restart_file != "none")
439  {
440  std::cout << " ====== Restarting from " << restart_file << std::endl;
441 
442  // Must have correct file type to restart
443  if (restart_file.rfind(".xdr") < restart_file.size())
444  _equation_system->read(restart_file,GRINSEnums::DECODE,
445  //libMesh::EquationSystems::READ_HEADER | // Allow for thermochemistry upgrades
446  libMesh::EquationSystems::READ_DATA |
447  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
448 
449  else if (restart_file.rfind(".xda") < restart_file.size())
450  _equation_system->read(restart_file,GRINSEnums::READ,
451  //libMesh::EquationSystems::READ_HEADER | // Allow for thermochemistry upgrades
452  libMesh::EquationSystems::READ_DATA |
453  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
454 
455  else
456  {
457  std::cerr << "Error: Restart filename must have .xdr or .xda extension!" << std::endl;
458  libmesh_error();
459  }
460 
461  const std::string system_name = input("screen-options/system_name", "GRINS" );
462 
463  MultiphysicsSystem& system =
464  _equation_system->get_system<MultiphysicsSystem>(system_name);
465 
466  // Update the old data
467  system.update();
468  }
469 
470  return;
471  }
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
void GRINS::Simulation::run ( )

Definition at line 303 of file simulation.C.

References _adjoint_parameters, _do_adjoint_solve, _equation_system, _error_estimator, _forward_parameters, _multiphysics_system, _output_adjoint, _output_residual, _output_residual_sensitivities, _output_solution_sensitivities, _output_vis, _postprocessing, _print_log_info, _print_qoi, _print_scalars, _solver, _timesteps_per_perflog, _timesteps_per_vis, _vis, GRINS::SolverContext::do_adjoint_solve, GRINS::SolverContext::equation_system, GRINS::SolverContext::error_estimator, GRINS::SolverContext::output_adjoint, GRINS::CompositeQoI::output_qoi(), GRINS::SolverContext::output_residual, GRINS::SolverContext::output_residual_sensitivities, GRINS::SolverContext::output_solution_sensitivities, GRINS::SolverContext::output_vis, GRINS::ParameterManager::parameter_vector, GRINS::SolverContext::postprocessing, GRINS::SolverContext::print_perflog, GRINS::SolverContext::print_qoi, GRINS::SolverContext::print_scalars, print_sim_info(), GRINS::SolverContext::system, GRINS::SolverContext::timesteps_per_perflog, GRINS::SolverContext::timesteps_per_vis, and GRINS::SolverContext::vis.

Referenced by main(), and run().

304  {
305  this->print_sim_info();
306 
307  SolverContext context;
308  context.system = _multiphysics_system;
309  context.equation_system = _equation_system;
310  context.vis = _vis;
311  context.output_adjoint = _output_adjoint;
312  context.timesteps_per_vis = _timesteps_per_vis;
313  context.timesteps_per_perflog = _timesteps_per_perflog;
314  context.output_vis = _output_vis;
315  context.output_residual = _output_residual;
316  context.output_residual_sensitivities = _output_residual_sensitivities;
317  context.output_solution_sensitivities = _output_solution_sensitivities;
318  context.print_scalars = _print_scalars;
319  context.print_perflog = _print_log_info;
320  context.postprocessing = _postprocessing;
321  context.error_estimator = _error_estimator;
322  context.print_qoi = _print_qoi;
323  context.do_adjoint_solve = _do_adjoint_solve;
324 
327  {
328  std::cout <<
329  "Error: output_residual_sensitivities is specified but\n" <<
330  "no forward sensitivity parameters have been specified.\n" <<
331  std::endl;
332  libmesh_error();
333  }
334 
337  {
338  std::cout <<
339  "Error: output_solution_sensitivities is specified but\n" <<
340  "no forward sensitivity parameters have been specified.\n" <<
341  std::endl;
342  libmesh_error();
343  }
344 
345  _solver->solve( context );
346 
347  if ( this->_print_qoi )
348  {
349  _multiphysics_system->assemble_qoi();
350  const CompositeQoI* my_qoi = libMesh::libmesh_cast_ptr<const CompositeQoI*>(this->_multiphysics_system->get_qoi());
351  my_qoi->output_qoi( std::cout );
352  }
353 
355  {
356  // Default: "calculate sensitivities for all QoIs"
357  libMesh::QoISet qois;
358 
359  const libMesh::ParameterVector & params =
361 
362  libMesh::SensitivityData sensitivities
363  (qois, *this->_multiphysics_system, params);
364 
365  _solver->adjoint_qoi_parameter_sensitivity
366  (context, qois, params, sensitivities);
367 
368  std::cout << "Adjoint sensitivities:" << std::endl;
369 
370  for (unsigned int q=0;
371  q != this->_multiphysics_system->qoi.size(); ++q)
372  {
373  for (unsigned int p=0; p != params.size(); ++p)
374  {
375  std::cout << "dq" << q << "/dp" << p << " = " <<
376  sensitivities[q][p] << std::endl;
377  }
378  }
379  }
380 
382  {
383  // Default: "calculate sensitivities for all QoIs"
384  libMesh::QoISet qois;
385 
386  const libMesh::ParameterVector & params =
388 
389  libMesh::SensitivityData sensitivities
390  (qois, *this->_multiphysics_system, params);
391 
392  _solver->forward_qoi_parameter_sensitivity
393  (context, qois, params, sensitivities);
394 
395  std::cout << "Forward sensitivities:" << std::endl;
396 
397  for (unsigned int q=0;
398  q != this->_multiphysics_system->qoi.size(); ++q)
399  {
400  for (unsigned int p=0; p != params.size(); ++p)
401  {
402  std::cout << "dq" << q << "/dp" << p << " = " <<
403  sensitivities[q][p] << std::endl;
404  }
405  }
406  }
407 
408  return;
409  }
ParameterManager _forward_parameters
Definition: simulation.h:164
std::tr1::shared_ptr< PostProcessedQuantities< libMesh::Real > > _postprocessing
Definition: simulation.h:140
unsigned int _timesteps_per_vis
Definition: simulation.h:157
std::tr1::shared_ptr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:128
libMesh::ParameterVector parameter_vector
std::tr1::shared_ptr< GRINS::Visualization > _vis
Definition: simulation.h:138
ParameterManager _adjoint_parameters
Definition: simulation.h:162
bool _output_residual_sensitivities
Definition: simulation.h:154
bool _output_solution_sensitivities
Definition: simulation.h:155
std::tr1::shared_ptr< GRINS::Solver > _solver
Definition: simulation.h:130
std::tr1::shared_ptr< libMesh::ErrorEstimator > _error_estimator
Definition: simulation.h:160
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:136
unsigned int _timesteps_per_perflog
Definition: simulation.h:158
void print_sim_info()
Definition: simulation.C:411

Member Data Documentation

ParameterManager GRINS::Simulation::_adjoint_parameters
protected

Definition at line 162 of file simulation.h.

Referenced by init_params(), and run().

bool GRINS::Simulation::_do_adjoint_solve
protected

Definition at line 167 of file simulation.h.

Referenced by init_adjoint_solve(), and run().

std::tr1::shared_ptr<libMesh::EquationSystems> GRINS::Simulation::_equation_system
protected
std::tr1::shared_ptr<libMesh::ErrorEstimator> GRINS::Simulation::_error_estimator
protected

Definition at line 160 of file simulation.h.

Referenced by run().

ParameterManager GRINS::Simulation::_forward_parameters
protected

Definition at line 164 of file simulation.h.

Referenced by init_params(), and run().

std::tr1::shared_ptr<libMesh::UnstructuredMesh> GRINS::Simulation::_mesh
protected

Definition at line 126 of file simulation.h.

Referenced by init_restart(), and print_sim_info().

GRINS::MultiphysicsSystem* GRINS::Simulation::_multiphysics_system
protected
bool GRINS::Simulation::_output_adjoint
protected

Definition at line 152 of file simulation.h.

Referenced by run().

bool GRINS::Simulation::_output_residual
protected

Definition at line 153 of file simulation.h.

Referenced by run().

bool GRINS::Simulation::_output_residual_sensitivities
protected

Definition at line 154 of file simulation.h.

Referenced by run().

bool GRINS::Simulation::_output_solution_sensitivities
protected

Definition at line 155 of file simulation.h.

Referenced by run().

bool GRINS::Simulation::_output_vis
protected

Definition at line 151 of file simulation.h.

Referenced by run().

std::tr1::shared_ptr<PostProcessedQuantities<libMesh::Real> > GRINS::Simulation::_postprocessing
protected

Definition at line 140 of file simulation.h.

Referenced by init_multiphysics_system(), and run().

bool GRINS::Simulation::_print_equation_system_info
protected

Definition at line 145 of file simulation.h.

Referenced by print_sim_info().

bool GRINS::Simulation::_print_log_info
protected

Definition at line 144 of file simulation.h.

Referenced by init_multiphysics_system(), and run().

bool GRINS::Simulation::_print_mesh_info
protected

Definition at line 143 of file simulation.h.

Referenced by print_sim_info().

bool GRINS::Simulation::_print_perflog
protected

Definition at line 146 of file simulation.h.

bool GRINS::Simulation::_print_qoi
protected

Definition at line 147 of file simulation.h.

Referenced by init_qois(), and run().

bool GRINS::Simulation::_print_scalars
protected

Definition at line 148 of file simulation.h.

Referenced by run().

std::tr1::shared_ptr<GRINS::Solver> GRINS::Simulation::_solver
protected

Definition at line 130 of file simulation.h.

Referenced by init_multiphysics_system(), and run().

std::string GRINS::Simulation::_system_name
protected

GRINS::Multiphysics system name.

Definition at line 133 of file simulation.h.

Referenced by get_multiphysics_system_name().

unsigned int GRINS::Simulation::_timesteps_per_perflog
protected

Definition at line 158 of file simulation.h.

Referenced by run().

unsigned int GRINS::Simulation::_timesteps_per_vis
protected

Definition at line 157 of file simulation.h.

Referenced by run().

std::tr1::shared_ptr<GRINS::Visualization> GRINS::Simulation::_vis
protected

Definition at line 138 of file simulation.h.

Referenced by run().


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

Generated on Mon Jun 22 2015 21:32:24 for GRINS-0.6.0 by  doxygen 1.8.9.1