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 
18 
21 {
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 
61  NumericVector<Number> & u_dot = *_sys.solutionUDot();
62  u_dot = *_solution;
64  u_dot.close();
65  _du_dot_du = 1. / _dt;
66 }
67 
68 void
69 LStableDirk4::computeADTimeDerivatives(DualReal & ad_u_dot, const dof_id_type & dof) const
70 {
72 }
73 
74 void
76 {
77  // Time at end of step
78  Real time_old = _fe_problem.timeOld();
79 
80  // Reset iteration counts
83 
84  // A for-loop would increment _stage too far, so we use an extra
85  // loop counter.
86  for (unsigned int current_stage = 1; current_stage <= _n_stages; ++current_stage)
87  {
88  // Set the current stage value
89  _stage = current_stage;
90 
91  // This ensures that all the Output objects in the OutputWarehouse
92  // have had solveSetup() called, and sets the default solver
93  // parameters for PETSc.
95 
96  _console << "Stage " << _stage << "\n";
97 
98  // Set the time for this stage
99  _fe_problem.time() = time_old + _c[_stage - 1] * _dt;
100 
101  // Do the solve
103 
104  // Update the iteration counts
107 
108  // Abort time step immediately on stage failure - see TimeIntegrator doc page
109  if (!_fe_problem.converged())
110  return;
111  }
112 }
113 
114 void
115 LStableDirk4::postResidual(NumericVector<Number> & residual)
116 {
117  // Error if _stage got messed up somehow.
118  if (_stage > _n_stages)
119  // the explicit cast prevents strange compiler weirdness with the static
120  // const variable and the variadic mooseError function
121  mooseError("LStableDirk4::postResidual(): Member variable _stage can only have values 1-",
122  (unsigned int)_n_stages,
123  ".");
124 
125  // In the standard RK notation, the residual of stage 1 of s is given by:
126  //
127  // R := M*(Y_i - y_n)/dt - \sum_{j=1}^s a_{ij} * f(t_n + c_j*dt, Y_j) = 0
128  //
129  // where:
130  // .) M is the mass matrix
131  // .) Y_i is the stage solution
132  // .) dt is the timestep, and is accounted for in the _Re_time residual.
133  // .) f are the "non-time" residuals evaluated for a given stage solution.
134  // .) The minus signs are already "baked in" to the residuals and so do not appear below.
135 
136  // Store this stage's non-time residual. We are calling operator=
137  // here, and that calls close().
139 
140  // Build up the residual for this stage.
141  residual.add(1., _Re_time);
142  for (unsigned int j = 0; j < _stage; ++j)
143  residual.add(_a[_stage - 1][j], *_stage_residuals[j]);
144  residual.close();
145 }
TimeIntegrator::_dt
Real & _dt
Definition: TimeIntegrator.h:156
TimeIntegrator::_n_linear_iterations
unsigned int _n_linear_iterations
Total number of linear iterations over all stages of the time step.
Definition: TimeIntegrator.h:167
SystemBase::addVector
virtual 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:543
LStableDirk4::solve
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: LStableDirk4.C:75
FEProblemBase::timeOld
virtual Real & timeOld() const
Definition: FEProblemBase.h:441
MooseObject::mooseInfo
void mooseInfo(Args &&... args) const
Definition: MooseObject.h:162
NonlinearSystemBase.h
FEProblem.h
MooseObject::mooseError
void mooseError(Args &&... args) const
Definition: MooseObject.h:141
LStableDirk4
Fourth-order diagonally implicit Runge Kutta method (Dirk) with five stages.
Definition: LStableDirk4.h:56
LStableDirk4::validParams
static InputParameters validParams()
Definition: LStableDirk4.C:20
TimeIntegrator::_Re_non_time
NumericVector< Number > & _Re_non_time
residual vector for non-time contributions
Definition: TimeIntegrator.h:162
TimeIntegrator
Base class for time integrators.
Definition: TimeIntegrator.h:43
FEProblemBase::time
virtual Real & time() const
Definition: FEProblemBase.h:440
TimeIntegrator::getNumLinearIterationsLastSolve
unsigned int getNumLinearIterationsLastSolve() const
Gets the number of linear iterations in the most recent solve.
Definition: TimeIntegrator.C:66
TimeIntegrator::_solution
const NumericVector< Number > *const & _solution
solution vectors
Definition: TimeIntegrator.h:150
DualReal
DualNumber< Real, DNDerivativeType > DualReal
Definition: DualReal.h:52
LStableDirk4::_a
static const Real _a[_n_stages][_n_stages]
Definition: LStableDirk4.h:94
TimeIntegrator::_sys
SystemBase & _sys
Definition: TimeIntegrator.h:140
LStableDirk4::computeTimeDerivativeHelper
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:99
defineLegacyParams
defineLegacyParams(LStableDirk4)
ConsoleStreamInterface::_console
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
Definition: ConsoleStreamInterface.h:31
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
LStableDirk4::_stage
unsigned int _stage
Definition: LStableDirk4.h:78
TimeIntegrator::validParams
static InputParameters validParams()
Definition: TimeIntegrator.C:21
LStableDirk4::computeADTimeDerivatives
virtual void computeADTimeDerivatives(DualReal &ad_u_dot, const dof_id_type &dof) const override
method for computing local automatic differentiation time derivatives
Definition: LStableDirk4.C:69
TimeIntegrator::_fe_problem
FEProblemBase & _fe_problem
Definition: TimeIntegrator.h:139
TimeIntegrator::_Re_time
NumericVector< Number > & _Re_time
residual vector for time contributions
Definition: TimeIntegrator.h:160
LStableDirk4::_stage_residuals
NumericVector< Number > * _stage_residuals[_n_stages]
Definition: LStableDirk4.h:87
TimeIntegrator::_n_nonlinear_iterations
unsigned int _n_nonlinear_iterations
Total number of nonlinear iterations over all stages of the time step.
Definition: TimeIntegrator.h:165
TimeIntegrator::getNumNonlinearIterationsLastSolve
unsigned int getNumNonlinearIterationsLastSolve() const
Gets the number of nonlinear iterations in the most recent solve.
Definition: TimeIntegrator.C:60
LStableDirk4::_n_stages
static const unsigned int _n_stages
Definition: LStableDirk4.h:84
NonlinearSystemBase::system
virtual System & system() override
Get the reference to the libMesh system.
Definition: NonlinearSystemBase.h:618
registerMooseObject
registerMooseObject("MooseApp", LStableDirk4)
TimeIntegrator::_du_dot_du
Real & _du_dot_du
Derivative of time derivative with respect to current solution: .
Definition: TimeIntegrator.h:148
LStableDirk4::computeTimeDerivatives
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: LStableDirk4.C:52
TimeIntegrator::_solution_old
const NumericVector< Number > & _solution_old
Definition: TimeIntegrator.h:151
FEProblemBase::converged
virtual bool converged() override
Definition: FEProblemBase.C:4704
LStableDirk4::_c
static const Real _c[_n_stages]
Definition: LStableDirk4.h:90
LStableDirk4::postResidual
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
Definition: LStableDirk4.C:115
TimeIntegrator::_nl
NonlinearSystemBase & _nl
Definition: TimeIntegrator.h:141
FEProblemBase::getNonlinearSystemBase
NonlinearSystemBase & getNonlinearSystemBase()
Definition: FEProblemBase.h:560
LStableDirk4.h
FEProblemBase::initPetscOutput
void initPetscOutput()
Reinitialize petsc output for proper linear/nonlinear iteration display.
Definition: FEProblemBase.C:4825
SystemBase::solutionUDot
virtual NumericVector< Number > * solutionUDot()=0
LStableDirk4::LStableDirk4
LStableDirk4(const InputParameters &parameters)
Definition: LStableDirk4.C:36
PetscSupport.h