https://mooseframework.inl.gov
TransientBase.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "Executioner.h"
13 #include "TimeIntegrator.h"
14 
15 // System includes
16 #include <string>
17 #include <fstream>
18 
19 class TimeStepper;
20 class FEProblemBase;
21 
26 class TransientBase : public Executioner
27 {
28 public:
30 
32 
33  virtual void init() override;
34 
35  virtual void execute() override;
36 
40  virtual void takeStep(Real input_dt = -1.0);
41 
45  virtual Real computeConstrainedDT();
46  virtual void estimateTimeError();
47 
51  virtual Real getDT();
52 
56  virtual bool keepGoing();
57 
61  virtual bool lastSolveConverged() const override;
62 
63  virtual void preExecute() override;
64 
65  virtual void postExecute() override;
66 
67  virtual void computeDT();
68 
69  virtual void preStep();
70 
71  virtual void postStep();
72 
76  virtual void incrementStepOrReject();
77 
78  virtual void endStep(Real input_time = -1.0);
79 
84  virtual void setTargetTime(Real target_time);
85 
89  virtual Real getTime() const { return _time; };
90 
95  virtual Real getTargetTime() { return _target_time; }
96 
100  virtual void setTime(Real t) { _time = t; };
101 
105  virtual void setTimeOld(Real t) { _time_old = t; };
106 
111 
117  const TimeStepper * getTimeStepper() const { return _time_stepper; }
118 
124  void setTimeStepper(TimeStepper & ts);
125 
129  virtual std::string getTimeStepperName() const override;
130 
137  virtual std::set<TimeIntegrator *> getTimeIntegrators() const = 0;
138 
143  virtual std::vector<std::string> getTimeIntegratorNames() const override;
144 
150 
155  std::set<Real> & syncTimes() { return _sync_times; }
156 
161  Real & dtMax() { return _dtmax; }
162 
167  Real & dtMin() { return _dtmin; }
168 
173  Real getStartTime() const { return _start_time; }
174 
179  Real getEndTime() const { return _end_time; }
180 
185  Real & endTime() { return _end_time; }
186 
192 
197  virtual void setTimestepTolerance(const Real & tolerance) { _timestep_tolerance = tolerance; }
198 
203  bool atSyncPoint() { return _at_sync_point; }
204 
210 
211  void parentOutputPositionChanged() override;
212 
216  virtual Real relativeSolutionDifferenceNorm() = 0;
217 
222  virtual void forceNumSteps(const unsigned int num_steps) { _num_steps = num_steps; }
223 
225  virtual SolveObject * timeStepSolveObject() { return _fixed_point_solve.get(); }
226 
227 protected:
230 
233 
235  const bool _check_aux;
236 
239 
241  int & _t_step;
249 
252 
255 
258 
262  unsigned int _num_steps;
264 
271 
272  std::set<Real> & _sync_times;
273 
274  bool _abort;
278  const bool _error_on_dtmin;
279 
284 
289 
291 
292  void setupTimeIntegrator();
293 
298 
300  bool convergedToSteadyState() const;
301 
302 private:
304  void constrainDTFromMultiApp(Real & dt_cur,
305  std::ostringstream & diag,
306  const ExecFlagType & execute_on) const;
307 };
void setTimeStepper(TimeStepper &ts)
Set the timestepper to use.
virtual Real relativeSolutionDifferenceNorm()=0
The relative L2 norm of the difference between solution and old solution vector.
virtual void preExecute() override
Override this for actions that should take place before execution.
Real & dtMax()
Get the maximum dt.
virtual void postExecute() override
Override this for actions that should take place after execution.
bool _steady_state_detection
Steady state detection variables:
virtual SolveObject * timeStepSolveObject()
Return the solve object wrapped by time stepper.
virtual void estimateTimeError()
Real _time_interval_output_interval
virtual std::string getTimeStepperName() const override
Get the name of the timestepper.
bool atSyncPoint()
Is the current step at a sync point (sync times, time interval, target time, etc)?
static InputParameters validParams()
Definition: TransientBase.C:40
virtual void endStep(Real input_time=-1.0)
virtual void setTargetTime(Real target_time)
Can be used to set the next "target time" which is a time to nail perfectly.
bool & _last_solve_converged
Whether or not the last solve converged.
std::set< Real > & _sync_times
virtual void setTimeOld(Real t)
Set the old time.
const TimeStepper * getTimeStepper() const
void setupTimeIntegrator()
virtual void setTime(Real t)
Set the current time.
const bool _error_on_dtmin
This parameter controls how the system will deal with _dt <= _dtmin If true, the time stepper is expe...
virtual void computeDT()
Base class for time stepping.
Definition: TimeStepper.h:22
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Real getStartTime() const
Return the start time.
virtual Real getDT()
Real & _time
Current time.
virtual void execute() override
Pure virtual execute function MUST be overridden by children classes.
unsigned int _num_steps
virtual std::set< TimeIntegrator * > getTimeIntegrators() const =0
Get the time integrators (time integration scheme) used Note that because some systems might be stead...
bool convergedToSteadyState() const
Determines whether the problem has converged to steady state.
std::set< Real > & syncTimes()
Get the set of sync times.
Real & _unconstrained_dt
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
AuxiliarySystem & _aux
Reference to auxiliary system base for faster access.
TimeStepper * _time_stepper
Real & _target_time
virtual Real computeConstrainedDT()
FEProblemBase & _problem
Here for backward compatibility.
virtual void setTimestepTolerance(const Real &tolerance)
Set the timestep tolerance.
virtual std::vector< std::string > getTimeIntegratorNames() const override
Get the name of the time integrator (time integration scheme) used.
TransientBase(const InputParameters &parameters)
virtual bool keepGoing()
Transient loop will continue as long as this keeps returning true.
virtual void preStep()
TimeStepper * getTimeStepper()
Pointer to the TimeStepper.
Real _steady_state_tolerance
bool & _time_interval
if to use time interval output
virtual void takeStep(Real input_dt=-1.0)
Do whatever is necessary to advance one step.
virtual void forceNumSteps(const unsigned int num_steps)
Set the number of time steps.
Real & _time_old
Previous time.
Real unconstrainedDT()
Get the unconstrained dt.
virtual Real getTime() const
Get the current time.
Definition: TransientBase.h:89
Real & _dt
Current delta t... or timestep size.
virtual void postStep()
Executioners are objects that do the actual work of solving your problem.
Definition: Executioner.h:30
Base class for transient executioners that use a FixedPointSolve solve object for multiapp-main app i...
Definition: TransientBase.h:26
bool _xfem_repeat_step
Whether step should be repeated due to xfem modifying the mesh.
Real _steady_state_start_time
Real _timestep_tolerance
virtual void init() override
Initialize the executioner.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real & endTime()
Get a modifiable reference to the end time.
Class for containing MooseEnum item information.
Definition: MooseEnumItem.h:18
Real & timestepTol()
Get the timestep tolerance.
Moose::TimeIntegratorType _time_scheme
Real getSolutionChangeNorm()
Get the Relative L2 norm of the change in the solution.
virtual Real getTargetTime()
Get the current target time.
Definition: TransientBase.h:95
Real & _solution_change_norm
virtual void incrementStepOrReject()
This is where the solve step is actually incremented.
TimeIntegratorType
Time integrators.
Definition: MooseTypes.h:902
Real & dtMin()
Get the minimal dt.
const bool _check_aux
Whether to use the auxiliary system solution to determine steady-states.
const InputParameters & parameters() const
Get the parameters of the object.
Moose::TimeIntegratorType getTimeScheme() const
Get the time scheme used.
const bool _normalize_solution_diff_norm_by_dt
Whether to divide the solution difference norm by dt.
virtual bool lastSolveConverged() const override
Whether or not the last solve converged.
Real getEndTime() const
Get the end time.
void constrainDTFromMultiApp(Real &dt_cur, std::ostringstream &diag, const ExecFlagType &execute_on) const
Constrain the timestep dt_cur by looking at the timesteps for the MultiApps on execute_on.
bool & _at_sync_point
A system that holds auxiliary variables.
std::unique_ptr< FixedPointSolve > _fixed_point_solve
Definition: Executioner.h:151
int & _t_step
Current timestep.
Real _next_interval_output_time
void parentOutputPositionChanged() override
Can be used by subclasses to call parentOutputPositionChanged() on the underlying FEProblemBase...