www.mooseframework.org
LStableDirk4.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 #include "LStableDirk4.h"
11 #include "NonlinearSystemBase.h"
12 #include "FEProblem.h"
13 #include "PetscSupport.h"
14 
16 
19 {
21  params.addClassDescription(
22  "Fourth-order diagonally implicit Runge Kutta method (Dirk) with five stages.");
23  return params;
24 }
25 
26 // Initialize static data
27 const Real LStableDirk4::_c[LStableDirk4::_n_stages] = {.25, 0., .5, 1., 1.};
28 
30  {.25, 0, 0, 0, 0},
31  {-.25, .25, 0, 0, 0},
32  {.125, .125, .25, 0, 0},
33  {-1.5, .75, 1.5, .25, 0},
34  {0, 1. / 6, 2. / 3, -1. / 12, .25}};
35 
37  : TimeIntegrator(parameters), _stage(1)
38 {
39  mooseInfo("LStableDirk4 and other multistage TimeIntegrators are known not to work with "
40  "Materials/AuxKernels that accumulate 'state' and should be used with caution.");
41 
42  // Name the stage residuals "residual_stage1", "residual_stage2", etc.
43  for (unsigned int stage = 0; stage < _n_stages; ++stage)
44  {
45  std::ostringstream oss;
46  oss << "residual_stage" << stage + 1;
47  _stage_residuals[stage] = &(_nl.addVector(oss.str(), false, GHOSTED));
48  }
49 }
50 
51 void
53 {
54  // We are multiplying by the method coefficients in postResidual(), so
55  // the time derivatives are of the same form at every stage although
56  // the current solution varies depending on the stage.
57  if (!_sys.solutionUDot())
58  mooseError("LStableDirk4: Time derivative of solution (`u_dot`) is not stored. Please set "
59  "uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
60 
62  u_dot = *_solution;
64  u_dot.close();
65  _du_dot_du = 1. / _dt;
66 }
67 
68 void
70  const dof_id_type & dof,
71  DualReal & /*ad_u_dotdot*/) const
72 {
74 }
75 
76 void
78 {
79  // Time at end of step
80  Real time_old = _fe_problem.timeOld();
81 
82  // Reset iteration counts
85 
86  // A for-loop would increment _stage too far, so we use an extra
87  // loop counter.
88  for (unsigned int current_stage = 1; current_stage <= _n_stages; ++current_stage)
89  {
90  // Set the current stage value
91  _stage = current_stage;
92 
93  // This ensures that all the Output objects in the OutputWarehouse
94  // have had solveSetup() called, and sets the default solver
95  // parameters for PETSc.
97 
98  _console << "Stage " << _stage << std::endl;
99 
100  // Set the time for this stage
101  _fe_problem.time() = time_old + _c[_stage - 1] * _dt;
102 
103  // Do the solve
104  _nl.system().solve();
105 
106  // Update the iteration counts
109 
110  // Abort time step immediately on stage failure - see TimeIntegrator doc page
112  return;
113  }
114 }
115 
116 void
118 {
119  // Error if _stage got messed up somehow.
120  if (_stage > _n_stages)
121  // the explicit cast prevents strange compiler weirdness with the static
122  // const variable and the variadic mooseError function
123  mooseError("LStableDirk4::postResidual(): Member variable _stage can only have values 1-",
124  (unsigned int)_n_stages,
125  ".");
126 
127  // In the standard RK notation, the residual of stage 1 of s is given by:
128  //
129  // R := M*(Y_i - y_n)/dt - \sum_{j=1}^s a_{ij} * f(t_n + c_j*dt, Y_j) = 0
130  //
131  // where:
132  // .) M is the mass matrix
133  // .) Y_i is the stage solution
134  // .) dt is the timestep, and is accounted for in the _Re_time residual.
135  // .) f are the "non-time" residuals evaluated for a given stage solution.
136  // .) The minus signs are already "baked in" to the residuals and so do not appear below.
137 
138  // Store this stage's non-time residual. We are calling operator=
139  // here, and that calls close().
141 
142  // Build up the residual for this stage.
143  residual.add(1., _Re_time);
144  for (unsigned int j = 0; j < _stage; ++j)
145  residual.add(_a[_stage - 1][j], *_stage_residuals[j]);
146  residual.close();
147 }
NonlinearSystemBase & _nl
Fourth-order diagonally implicit Runge Kutta method (Dirk) with five stages.
Definition: LStableDirk4.h:51
virtual void initPetscOutputAndSomeSolverSettings()
Reinitialize PETSc output for proper linear/nonlinear iteration display.
virtual NumericVector< Number > * solutionUDot()=0
virtual Real & time() const
static const Real _a[_n_stages][_n_stages]
Definition: LStableDirk4.h:90
virtual bool converged(const unsigned int nl_sys_num)
Eventually we want to convert this virtual over to taking a nonlinear system number argument...
Definition: SubProblem.h:101
NumericVector< Number > & addVector(const std::string &vector_name, const bool project, const ParallelType type)
Adds a solution length vector to the system.
Definition: SystemBase.C:597
FEProblemBase & _fe_problem
void mooseInfo(Args &&... args) const
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: LStableDirk4.C:52
LStableDirk4(const InputParameters &parameters)
Definition: LStableDirk4.C:36
DualNumber< Real, DNDerivativeType, true > DualReal
Definition: DualReal.h:49
NumericVector< Number > & _Re_non_time
residual vector for non-time contributions
SystemBase & _sys
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeADTimeDerivatives(DualReal &ad_u_dot, const dof_id_type &dof, DualReal &ad_u_dotdot) const override
method for computing local automatic differentiation time derivatives
Definition: LStableDirk4.C:69
unsigned int _stage
Definition: LStableDirk4.h:74
registerMooseObject("MooseApp", LStableDirk4)
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
Definition: LStableDirk4.C:117
GHOSTED
unsigned int _n_linear_iterations
Total number of linear iterations over all stages of the time step.
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1125
virtual void close()=0
static const Real _c[_n_stages]
Definition: LStableDirk4.h:86
Real & _du_dot_du
Derivative of time derivative with respect to current solution: .
const NumericVector< Number > *const & _solution
solution vectors
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
unsigned int getNumLinearIterationsLastSolve() const
Gets the number of linear iterations in the most recent solve.
virtual System & system() override
Get the reference to the libMesh system.
NumericVector< Number > * _stage_residuals[_n_stages]
Definition: LStableDirk4.h:83
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
NumericVector< Number > & _Re_time
residual vector for time contributions
static InputParameters validParams()
Definition: LStableDirk4.C:18
static const unsigned int _n_stages
Definition: LStableDirk4.h:80
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
virtual Real & timeOld() const
const NumericVector< Number > & _solution_old
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: LStableDirk4.C:77
unsigned int _n_nonlinear_iterations
Total number of nonlinear iterations over all stages of the time step.
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
unsigned int getNumNonlinearIterationsLastSolve() const
Gets the number of nonlinear iterations in the most recent solve.
virtual void add(const numeric_index_type i, const Number value)=0
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
Definition: LStableDirk4.h:95
static InputParameters validParams()
uint8_t dof_id_type