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 #include <optional>
19 
20 class TimeStepper;
21 class FEProblemBase;
22 
27 class TransientBase : public Executioner
28 {
29 public:
32 
34 
35  virtual void init() override;
36 
37  virtual void execute() override;
38 
42  virtual void takeStep(Real input_dt = -1.0);
43 
47  virtual Real computeConstrainedDT();
48  virtual void estimateTimeError();
49 
53  virtual Real getDT();
54 
58  virtual bool keepGoing();
59 
63  virtual bool lastSolveConverged() const override;
64 
65  virtual void preExecute() override;
66 
67  virtual void postExecute() override;
68 
69  void computeDT();
70 
71  virtual void preStep();
72 
73  virtual void postStep();
74 
78  virtual void incrementStepOrReject();
79 
80  virtual void endStep(Real input_time = -1.0);
81 
86  virtual void setTargetTime(Real target_time);
87 
91  virtual Real getTime() const { return _time; };
92 
97  virtual Real getTargetTime() { return _target_time; }
98 
102  virtual void setTime(Real t) { _time = t; };
103 
107  virtual void setTimeOld(Real t) { _time_old = t; };
108 
115  Real computeSolutionChangeNorm(bool check_aux, bool normalize_by_dt) const;
116 
122  virtual Real relativeSolutionDifferenceNorm(bool check_aux) const = 0;
123 
129  const TimeStepper * getTimeStepper() const { return _time_stepper; }
130 
136  void setTimeStepper(TimeStepper & ts);
137 
141  virtual std::string getTimeStepperName() const override;
142 
149  virtual std::set<TimeIntegrator *> getTimeIntegrators() const = 0;
150 
155  virtual std::vector<std::string> getTimeIntegratorNames() const override;
156 
162 
167  std::set<Real> & syncTimes() { return _sync_times; }
168 
173  Real & dtMax() { return _dtmax; }
174 
179  Real & dtMin() { return _dtmin; }
180 
185  Real getStartTime() const { return _start_time; }
186 
191  Real getEndTime() const { return _end_time; }
192 
197  Real & endTime() { return _end_time; }
198 
204 
209  virtual void setTimestepTolerance(const Real & tolerance) { _timestep_tolerance = tolerance; }
210 
215  bool atSyncPoint() { return _at_sync_point; }
216 
222 
223  void parentOutputPositionChanged() override;
224 
229  virtual void forceNumSteps(const unsigned int num_steps) { _num_steps = num_steps; }
230 
232  virtual SolveObject * timeStepSolveObject() { return _fixed_point_solve.get(); }
233 
235  bool convergedToSteadyState() const;
236 
237 protected:
240 
243 
246 
248  int & _t_step;
256 
259 
262 
265 
269  unsigned int _num_steps;
271 
277 
278  std::set<Real> & _sync_times;
279 
280  bool _abort;
284  const bool _error_on_dtmin;
285 
290 
295 
296  void setupTimeIntegrator();
297 
298 private:
300  void constrainDTFromMultiApp(Real & dt_cur,
301  std::ostringstream & diag,
302  const ExecFlagType & execute_on) const;
303 
305  std::optional<int> _test_restep_step;
307  std::optional<Real> _test_restep_time;
310 };
void setTimeStepper(TimeStepper &ts)
Set the timestepper to use.
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.
const Real _steady_state_start_time
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:69
virtual void endStep(Real input_time=-1.0)
Real computeSolutionChangeNorm(bool check_aux, bool normalize_by_dt) const
Compute the relative L2 norm of the change in the solution.
virtual void setTargetTime(Real target_time)
Can be used to set the next "target time" which is a time to nail perfectly.
std::optional< Real > _test_restep_time
If the time is greater than this then we fail and repeat if –test-restep is enabled.
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()
std::optional< int > _test_restep_step
The timestep we fail and repeat if –test-restep is enabled.
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...
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.
const bool _steady_state_detection
Steady state detection variables:
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:91
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:27
bool _xfem_repeat_step
Whether step should be repeated due to xfem modifying the mesh.
Real _timestep_tolerance
bool _testing_restep
Whether or not the last timestep we solved is being repeated with –test-restep.
virtual Real relativeSolutionDifferenceNorm(bool check_aux) const =0
The relative L2 norm of the difference between solution and old solution vector.
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
virtual Real getTargetTime()
Get the current target time.
Definition: TransientBase.h:97
static InputParameters defaultSteadyStateConvergenceParams()
Definition: TransientBase.C:42
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 InputParameters & parameters() const
Get the parameters of the object.
Moose::TimeIntegratorType getTimeScheme() const
Get the time scheme used.
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...