GRINS-0.7.0
simulation.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
27 #include "grins/simulation.h"
28 
29 // GRINS
30 #include "grins/grins_enums.h"
32 #include "grins/multiphysics_sys.h"
33 #include "grins/solver_context.h"
36 #include "grins/physics_builder.h"
38 
39 // libMesh
40 #include "libmesh/dof_map.h"
41 #include "libmesh/parameter_vector.h"
42 #include "libmesh/qoi_set.h"
43 #include "libmesh/sensitivity_data.h"
44 
45 namespace GRINS
46 {
47 
48  Simulation::Simulation( const GetPot& input,
49  SimulationBuilder& sim_builder,
50  const libMesh::Parallel::Communicator &comm )
51  : _mesh( sim_builder.build_mesh(input, comm) ),
52  _equation_system( new libMesh::EquationSystems( *_mesh ) ),
53  _solver( sim_builder.build_solver(input) ),
54  _system_name( input("screen-options/system_name", "GRINS" ) ),
55  _multiphysics_system( &(_equation_system->add_system<MultiphysicsSystem>( _system_name )) ),
56  _vis( sim_builder.build_vis(input, comm) ),
57  _postprocessing( sim_builder.build_postprocessing(input) ),
58  _print_mesh_info( input("screen-options/print_mesh_info", false ) ),
59  _print_log_info( input("screen-options/print_log_info", false ) ),
60  _print_equation_system_info( input("screen-options/print_equation_system_info", false ) ),
61  _print_qoi( input("screen-options/print_qoi", false ) ),
62  _print_scalars( input("screen-options/print_scalars", false ) ),
63  _output_vis( input("vis-options/output_vis", false ) ),
64  _output_adjoint( input("vis-options/output_adjoint", false ) ),
65  _output_residual( input( "vis-options/output_residual", false ) ),
66  _output_residual_sensitivities( input( "vis-options/output_residual_sensitivities", false ) ),
67  _output_solution_sensitivities( input( "vis-options/output_solution_sensitivities", false ) ),
68  _timesteps_per_vis( input("vis-options/timesteps_per_vis", 1 ) ),
69  _timesteps_per_perflog( input("screen-options/timesteps_per_perflog", 0 ) ),
70  _error_estimator_options(input),
71  _error_estimator(), // effectively NULL
72  _do_adjoint_solve(false), // Helper function will set final value
73  _have_restart(false)
74  {
75  libmesh_deprecated();
76 
77  this->init_multiphysics_system(input);
78 
79  this->init_qois(input,sim_builder);
80 
81  this->init_params(input,sim_builder);
82 
83  this->init_adjoint_solve(input,_output_adjoint);
84 
85  // Must be called after setting QoI on the MultiphysicsSystem
86  this->build_error_estimator(input);
87 
89  this->init_restart(input,sim_builder,comm);
90 
91  this->check_for_unused_vars(input, false /*warning only*/);
92 
93  }
94 
95  Simulation::Simulation( const GetPot& input,
96  GetPot& command_line,
97  SimulationBuilder& sim_builder,
98  const libMesh::Parallel::Communicator &comm )
99  : _mesh( sim_builder.build_mesh(input, comm) ),
100  _equation_system( new libMesh::EquationSystems( *_mesh ) ),
101  _solver( sim_builder.build_solver(input) ),
102  _system_name( input("screen-options/system_name", "GRINS" ) ),
103  _multiphysics_system( &(_equation_system->add_system<MultiphysicsSystem>( _system_name )) ),
104  _vis( sim_builder.build_vis(input, comm) ),
105  _postprocessing( sim_builder.build_postprocessing(input) ),
106  _print_mesh_info( input("screen-options/print_mesh_info", false ) ),
107  _print_log_info( input("screen-options/print_log_info", false ) ),
108  _print_equation_system_info( input("screen-options/print_equation_system_info", false ) ),
109  _print_qoi( input("screen-options/print_qoi", false ) ),
110  _print_scalars( input("screen-options/print_scalars", false ) ),
111  _output_vis( input("vis-options/output_vis", false ) ),
112  _output_adjoint( input("vis-options/output_adjoint", false ) ),
113  _output_residual( input( "vis-options/output_residual", false ) ),
114  _output_residual_sensitivities( input( "vis-options/output_residual_sensitivities", false ) ),
115  _output_solution_sensitivities( input( "vis-options/output_solution_sensitivities", false ) ),
116  _timesteps_per_vis( input("vis-options/timesteps_per_vis", 1 ) ),
117  _timesteps_per_perflog( input("screen-options/timesteps_per_perflog", 0 ) ),
118  _error_estimator_options(input),
119  _error_estimator(), // effectively NULL
120  _do_adjoint_solve(false), // Helper function will set final value
121  _have_restart(false)
122  {
123  this->init_multiphysics_system(input);
124 
125  this->init_qois(input,sim_builder);
126 
127  this->init_params(input,sim_builder);
128 
129  this->init_adjoint_solve(input,_output_adjoint);
130 
131  // Must be called after setting QoI on the MultiphysicsSystem
132  this->build_error_estimator(input);
133 
135  this->init_restart(input,sim_builder,comm);
136 
137  bool warning_only = command_line.search("--warn-only-unused-var");
138  this->check_for_unused_vars(input, warning_only );
139 
140  }
141 
142  void Simulation::init_multiphysics_system( const GetPot& input )
143  {
144  // Only print libMesh logging info if the user requests it
145  libMesh::perflog.disable_logging();
146  if( this->_print_log_info ) libMesh::perflog.enable_logging();
147 
148  PhysicsList physics_list = PhysicsBuilder::build_physics_map(input);
149 
150  _multiphysics_system->attach_physics_list( physics_list );
151 
153 
155 
156  /* Postprocessing needs to be initialized before the solver since that's
157  where equation_system gets init'ed */
159 
160  _solver->initialize( input, _equation_system, _multiphysics_system );
161 
162  // Useful for debugging
163  if( input("screen-options/print_dof_constraints", false ) )
164  {
165  _multiphysics_system->get_dof_map().print_dof_constraints();
166  }
167 
168  // This *must* be done after equation_system->init in order to get variable indices
169  // Set any extra quadrature order the user requested. By default, is 0.
170  _multiphysics_system->extra_quadrature_order = StrategiesParsing::extra_quadrature_order(input);
171  }
172 
173  void Simulation::init_qois( const GetPot& input, SimulationBuilder& sim_builder )
174  {
175  // If the user actually asks for a QoI, then we add it.
176  SharedPtr<CompositeQoI> qois = sim_builder.build_qoi( input );
177  if( qois->n_qois() > 0 )
178  {
179  // This *must* be done after equation_system->init in order to get variable indices
180  qois->init(input, *_multiphysics_system );
181 
182  /* Note that we are effectively transfering ownership of the qoi pointer because
183  it will be cloned in _multiphysics_system and all the calculations are done there. */
184  _multiphysics_system->attach_qoi( qois.get() );
185  }
186  else if (_print_qoi)
187  {
188  std::cout << "Error: print_qoi is specified but\n" <<
189  "no QoIs have been specified.\n" << std::endl;
190  libmesh_error();
191  }
192 
193  return;
194  }
195 
196  void Simulation::init_params( const GetPot& input,
197  SimulationBuilder& /*sim_builder*/ )
198  {
199  unsigned int n_adjoint_parameters =
200  input.vector_variable_size("QoI/adjoint_sensitivity_parameters");
201 
202  unsigned int n_forward_parameters =
203  input.vector_variable_size("QoI/forward_sensitivity_parameters");
204 
205  // If the user actually asks for parameter sensitivities, then we
206  // set up the parameter vectors to use.
207  if ( n_adjoint_parameters )
208  {
209  // If we're doing adjoint sensitivities, dq/dp only makes
210  // sense if we have q
211  CompositeQoI* qoi =
212  libMesh::cast_ptr<CompositeQoI*>
213  (this->_multiphysics_system->get_qoi());
214 
215  if (!qoi)
216  {
217  std::cout <<
218  "Error: adjoint_sensitivity_parameters are specified but\n"
219  << "no QoIs have been specified.\n" << std::endl;
220  libmesh_error();
221  }
222 
224  (input, "QoI/adjoint_sensitivity_parameters",
225  *this->_multiphysics_system, qoi);
226  }
227 
228  if ( n_forward_parameters )
229  {
230  // If we're doing forward sensitivities, du/dp can make
231  // sense even with no q defined
232  CompositeQoI* qoi =
233  dynamic_cast<CompositeQoI*>
234  (this->_multiphysics_system->get_qoi());
235 
236  // dynamic_cast returns NULL if our QoI isn't a CompositeQoI;
237  // i.e. if there were no QoIs that made us bother setting up
238  // the CompositeQoI object. Passing NULL tells
239  // ParameterManager not to bother asking for qoi registration
240  // of parameters.
241 
243  (input, "QoI/forward_sensitivity_parameters",
244  *this->_multiphysics_system, qoi);
245  }
246  }
247 
248 
249  void Simulation::init_restart( const GetPot& input, SimulationBuilder& sim_builder,
250  const libMesh::Parallel::Communicator &comm )
251  {
252  this->read_restart( input );
253 
254  /* We do this here only if there's a restart file. Otherwise, this was done
255  at mesh construction time */
256  sim_builder.mesh_builder().do_mesh_refinement_from_input( input, comm, *_mesh );
257 
258  /* \todo Any way to tell if the mesh got refined so we don't unnecessarily
259  call reinit()? */
260  _equation_system->reinit();
261 
262  return;
263  }
264 
265  void Simulation::check_for_unused_vars( const GetPot& input, bool warning_only )
266  {
267  /* Everything should be set up now, so check if there's any unused variables
268  in the input file. If so, then tell the user what they were and error out. */
269  std::vector<std::string> unused_vars = input.unidentified_variables();
270 
271  bool unused_vars_detected = false;
272 
273  // String we use to check if there is a variable in the Materials section
274  std::string mat_string = "Materials/";
275 
276  // If we do have unused variables, make sure they are in the Materials
277  // section since we're allowing those to be present and not used.
278  if( !unused_vars.empty() )
279  {
280  int n_total = unused_vars.size();
281 
282  int n_mats = 0;
283 
284  for( int v = 0; v < n_total; v++ )
285  if( (unused_vars[v]).find(mat_string) != std::string::npos )
286  n_mats += 1;
287 
288  libmesh_assert_greater_equal( n_total, n_mats );
289 
290  if( n_mats < n_total )
291  unused_vars_detected = true;
292  }
293 
294  if( unused_vars_detected )
295  {
296  libMesh::err << "==========================================================" << std::endl;
297  if( warning_only )
298  libMesh::err << "Warning: ";
299  else
300  libMesh::err << "Error: ";
301 
302  libMesh::err << "Found unused variables!" << std::endl;
303 
304  for( std::vector<std::string>::const_iterator it = unused_vars.begin();
305  it != unused_vars.end(); ++it )
306  {
307  // Don't print out any unused Material variables, since we allow these
308  if( (*it).find(mat_string) != std::string::npos )
309  continue;
310 
311  libMesh::err << *it << std::endl;
312  }
313  libMesh::err << "==========================================================" << std::endl;
314 
315  if( !warning_only )
316  libmesh_error();
317  }
318 
319  return;
320  }
321 
323  {
324  this->print_sim_info();
325 
326  SolverContext context;
327  context.system = _multiphysics_system;
329  context.vis = _vis;
333  context.output_vis = _output_vis;
337  context.print_scalars = _print_scalars;
338  context.print_perflog = _print_log_info;
341  context.print_qoi = _print_qoi;
343  context.have_restart = _have_restart;
344 
347  {
348  std::cout <<
349  "Error: output_residual_sensitivities is specified but\n" <<
350  "no forward sensitivity parameters have been specified.\n" <<
351  std::endl;
352  libmesh_error();
353  }
354 
357  {
358  std::cout <<
359  "Error: output_solution_sensitivities is specified but\n" <<
360  "no forward sensitivity parameters have been specified.\n" <<
361  std::endl;
362  libmesh_error();
363  }
364 
365  _solver->solve( context );
366 
367  if ( this->_print_qoi )
368  {
369  _multiphysics_system->assemble_qoi();
370  const CompositeQoI* my_qoi = libMesh::libmesh_cast_ptr<const CompositeQoI*>(this->_multiphysics_system->get_qoi());
371  my_qoi->output_qoi( std::cout );
372  }
373 
375  {
376  // Default: "calculate sensitivities for all QoIs"
377  libMesh::QoISet qois;
378 
379  const libMesh::ParameterVector & params =
381 
382  libMesh::SensitivityData sensitivities
383  (qois, *this->_multiphysics_system, params);
384 
385  _solver->adjoint_qoi_parameter_sensitivity
386  (context, qois, params, sensitivities);
387 
388  std::cout << "Adjoint sensitivities:" << std::endl;
389 
390  for (unsigned int q=0;
391  q != this->_multiphysics_system->qoi.size(); ++q)
392  {
393  for (unsigned int p=0; p != params.size(); ++p)
394  {
395  std::cout << "dq" << q << "/dp" << p << " = " <<
396  sensitivities[q][p] << std::endl;
397  }
398  }
399  }
400 
402  {
403  // Default: "calculate sensitivities for all QoIs"
404  libMesh::QoISet qois;
405 
406  const libMesh::ParameterVector & params =
408 
409  libMesh::SensitivityData sensitivities
410  (qois, *this->_multiphysics_system, params);
411 
412  _solver->forward_qoi_parameter_sensitivity
413  (context, qois, params, sensitivities);
414 
415  std::cout << "Forward sensitivities:" << std::endl;
416 
417  for (unsigned int q=0;
418  q != this->_multiphysics_system->qoi.size(); ++q)
419  {
420  for (unsigned int p=0; p != params.size(); ++p)
421  {
422  std::cout << "dq" << q << "/dp" << p << " = " <<
423  sensitivities[q][p] << std::endl;
424  }
425  }
426  }
427 
428  return;
429  }
430 
432  {
433  // Print mesh info if the user wants it
434  if( this->_print_mesh_info ) this->_mesh->print_info();
435 
436  // Print info if requested
437  if( this->_print_equation_system_info ) this->_equation_system->print_info();
438 
439  return;
440  }
441 
442  SharedPtr<libMesh::EquationSystems> Simulation::get_equation_system()
443  {
444  return _equation_system;
445  }
446 
447  libMesh::Number Simulation::get_qoi_value( unsigned int qoi_index ) const
448  {
449  const CompositeQoI* qoi = libMesh::libmesh_cast_ptr<const CompositeQoI*>(this->_multiphysics_system->get_qoi());
450  return qoi->get_qoi_value(qoi_index);
451  }
452 
453  void Simulation::read_restart( const GetPot& input )
454  {
455  const std::string restart_file = SimulationParsing::restart_file(input);
456 
457  // Most of this was pulled from FIN-S
458  if (restart_file != "none")
459  {
460  std::cout << " ====== Restarting from " << restart_file << std::endl;
461 
462  // Must have correct file type to restart
463  if (restart_file.rfind(".xdr") < restart_file.size())
464  _equation_system->read(restart_file,GRINSEnums::DECODE,
465  //libMesh::EquationSystems::READ_HEADER | // Allow for thermochemistry upgrades
466  libMesh::EquationSystems::READ_DATA |
467  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
468 
469  else if (restart_file.rfind(".xda") < restart_file.size())
470  _equation_system->read(restart_file,GRINSEnums::READ,
471  //libMesh::EquationSystems::READ_HEADER | // Allow for thermochemistry upgrades
472  libMesh::EquationSystems::READ_DATA |
473  libMesh::EquationSystems::READ_ADDITIONAL_DATA);
474 
475  else
476  {
477  std::cerr << "Error: Restart filename must have .xdr or .xda extension!" << std::endl;
478  libmesh_error();
479  }
480 
481  this->_have_restart = true;
482 
483  const std::string system_name = input("screen-options/system_name", "GRINS" );
484 
485  MultiphysicsSystem& system =
486  _equation_system->get_system<MultiphysicsSystem>(system_name);
487 
488  // Update the old data
489  system.update();
490  }
491 
492  return;
493  }
494 
495  void Simulation::init_adjoint_solve( const GetPot& input, bool output_adjoint )
496  {
497  // Check if we're doing an adjoint solve
499 
500  const libMesh::DifferentiableQoI* raw_qoi = _multiphysics_system->get_qoi();
501  const CompositeQoI* qoi = dynamic_cast<const CompositeQoI*>( raw_qoi );
502 
503  // If we are trying to do an adjoint solve without a QoI, that's an error
504  // If there are no QoIs, the CompositeQoI never gets built and qoi will be NULL
505  if( _do_adjoint_solve && !qoi )
506  {
507  libmesh_error_msg("Error: Adjoint solve requested, but no QoIs detected.");
508  }
509 
510  /* If we are not doing an adjoint solve, but adjoint output was requested:
511  that's an error. */
512  if( !_do_adjoint_solve && output_adjoint )
513  {
514  libmesh_error_msg("Error: Adjoint output requested, but no adjoint solve requested.");
515  }
516  }
517 
518  bool Simulation::check_for_adjoint_solve( const GetPot& input ) const
519  {
520  std::string error_estimator = _error_estimator_options.estimator_type();
521 
522  bool do_adjoint_solve = false;
523 
524  // First check if the error estimator requires an adjoint solve
526  do_adjoint_solve = true;
527 
528  // Next check if parameter sensitivies require an adjoint solve
530  do_adjoint_solve = true;
531 
532  // Now check if the user requested to do an adjoint solve regardless
535  do_adjoint_solve = true;
536 
537  return do_adjoint_solve;
538  }
539 
540  void Simulation::build_error_estimator(const GetPot& input)
541  {
542  std::string estimator_type = _error_estimator_options.estimator_type();
543  if( estimator_type != std::string("none") )
544  {
548  _error_estimator.reset( (ErrorEstimatorFactoryBase::build(estimator_type)).release() );
549  }
550  }
551 
552 #ifdef GRINS_USE_GRVY_TIMERS
553  void Simulation::attach_grvy_timer( GRVY::GRVY_Timer_Class* grvy_timer )
554  {
555  this->_multiphysics_system->attach_grvy_timer( grvy_timer );
556  return;
557  }
558 #endif
559 
560 } // namespace GRINS
static PhysicsList build_physics_map(const GetPot &input)
Returns container of constructed Physics objects.
static bool do_adjoint_solve(const GetPot &input)
Option to let user manually trigger adjoint solve.
void attach_physics_list(PhysicsList physics_list)
PhysicsList gets built by GRINS::PhysicsFactory and attached here.
void read_restart(const GetPot &input)
Definition: simulation.C:453
static void set_estimator_options(const ErrorEstimatorOptions &estimator_options)
SharedPtr< libMesh::EquationSystems > equation_system
SharedPtr< libMesh::EquationSystems > get_equation_system()
Definition: simulation.C:442
const std::string & estimator_type() const
void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Each Physics will register their postprocessed quantities with this call.
void init_multiphysics_system(const GetPot &input)
Helper function.
Definition: simulation.C:142
ParameterManager _forward_parameters
Definition: simulation.h:163
SharedPtr< PostProcessedQuantities< libMesh::Real > > postprocessing
void init_params(const GetPot &input, SimulationBuilder &sim_builder)
Helper function.
Definition: simulation.C:196
unsigned int _timesteps_per_vis
Definition: simulation.h:155
static bool have_restart(const GetPot &input)
SharedPtr< PostProcessedQuantities< libMesh::Real > > _postprocessing
Definition: simulation.h:138
SharedPtr< libMesh::ErrorEstimator > error_estimator
SharedPtr< CompositeQoI > build_qoi(const GetPot &input)
GRINS namespace.
unsigned int timesteps_per_perflog
void build_error_estimator(const GetPot &input)
Definition: simulation.C:540
static int extra_quadrature_order(const GetPot &input)
SharedPtr< libMesh::ErrorEstimator > _error_estimator
Definition: simulation.h:159
bool _print_equation_system_info
Definition: simulation.h:143
libMesh::Number get_qoi_value(unsigned int qoi_index) const
Definition: simulation.C:447
const MeshBuilder & mesh_builder() const
unsigned int timesteps_per_vis
SharedPtr< GRINS::Solver > _solver
Definition: simulation.h:128
static void set_system(MultiphysicsSystem &system)
libMesh::ParameterVector parameter_vector
void init_adjoint_solve(const GetPot &input, bool output_adjoint)
Helper function.
Definition: simulation.C:495
GRINS::MultiphysicsSystem * system
SharedPtr< GRINS::Visualization > vis
void init_qois(const GetPot &input, SimulationBuilder &sim_builder)
Helper function.
Definition: simulation.C:173
Interface with libMesh for solving Multiphysics problems.
ParameterManager _adjoint_parameters
Definition: simulation.h:161
std::map< std::string, SharedPtr< GRINS::Physics > > PhysicsList
Container for GRINS::Physics object pointers.
Definition: var_typedefs.h:59
bool _output_residual_sensitivities
Definition: simulation.h:152
bool _output_solution_sensitivities
Definition: simulation.h:153
void output_qoi(std::ostream &out) const
Basic output for computed QoI's.
ErrorEstimatorOptions _error_estimator_options
Definition: simulation.h:158
void check_for_unused_vars(const GetPot &input, bool warning_only)
Helper function.
Definition: simulation.C:265
static std::string restart_file(const GetPot &input)
void init_restart(const GetPot &input, SimulationBuilder &sim_builder, const libMesh::Parallel::Communicator &comm)
Helper function.
Definition: simulation.C:249
Simple class to hold objects passed to Solver::solve.
SharedPtr< libMesh::UnstructuredMesh > _mesh
Definition: simulation.h:124
GRINS::MultiphysicsSystem * _multiphysics_system
Definition: simulation.h:134
unsigned int _timesteps_per_perflog
Definition: simulation.h:156
void do_mesh_refinement_from_input(const GetPot &input, const libMesh::Parallel::Communicator &comm, libMesh::UnstructuredMesh &mesh) const
Refine the mesh based on user input parameters.
Definition: mesh_builder.C:355
SharedPtr< GRINS::Visualization > _vis
Definition: simulation.h:136
virtual void read_input_options(const GetPot &input)
Reads input options for this class and all physics that are enabled.
libMesh::Number get_qoi_value(unsigned int qoi_index) const
Accessor for value of QoI for given qoi_index.
void print_sim_info()
Definition: simulation.C:431
virtual void initialize(const GetPot &input, const std::string &parameters_varname, GRINS::MultiphysicsSystem &system, GRINS::CompositeQoI *qoi)
bool check_for_adjoint_solve(const GetPot &input) const
Helper function.
Definition: simulation.C:518
SharedPtr< libMesh::EquationSystems > _equation_system
Definition: simulation.h:126

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