www.mooseframework.org
LStableDirk2.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 "LStableDirk2.h"
11 #include "NonlinearSystem.h"
12 #include "FEProblem.h"
13 #include "PetscSupport.h"
14 
16 
19 {
21  params.addClassDescription(
22  "Second order diagonally implicit Runge Kutta method (Dirk) with two stages.");
23  return params;
24 }
25 
27  : TimeIntegrator(parameters),
28  _stage(1),
29  _residual_stage1(_nl.addVector("residual_stage1", false, GHOSTED)),
30  _residual_stage2(_nl.addVector("residual_stage2", false, GHOSTED)),
31  _alpha(1. - 0.5 * std::sqrt(2))
32 {
33  mooseInfo("LStableDirk2 and other multistage TimeIntegrators are known not to work with "
34  "Materials/AuxKernels that accumulate 'state' and should be used with caution.");
35 }
36 
37 void
39 {
40  // We are multiplying by the method coefficients in postResidual(), so
41  // the time derivatives are of the same form at every stage although
42  // the current solution varies depending on the stage.
43  if (!_sys.solutionUDot())
44  mooseError("LStableDirk2: Time derivative of solution (`u_dot`) is not stored. Please set "
45  "uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
46 
48  u_dot = *_solution;
50  u_dot.close();
51  _du_dot_du = 1. / _dt;
52 }
53 
54 void
56  const dof_id_type & dof,
57  DualReal & /*ad_u_dotdot*/) const
58 {
60 }
61 
62 void
64 {
65  // Time at end of step
66  Real time_new = _fe_problem.time();
67 
68  // Time at beginning of step
69  Real time_old = _fe_problem.timeOld();
70 
71  // Time at stage 1
72  Real time_stage1 = time_old + _alpha * _dt;
73 
74  // Reset iteration counts
77 
78  // Compute first stage
80  _console << "1st stage" << std::endl;
81  _stage = 1;
82  _fe_problem.time() = time_stage1;
83  _nl.system().solve();
86 
87  // Abort time step immediately on stage failure - see TimeIntegrator doc page
89  return;
90 
91  // Compute second stage
93  _console << "2nd stage" << std::endl;
94  _stage = 2;
95  _fe_problem.timeOld() = time_stage1;
96  _fe_problem.time() = time_new;
97  _nl.system().solve();
100 
101  // Reset time at beginning of step to its original value
102  _fe_problem.timeOld() = time_old;
103 }
104 
105 void
107 {
108  if (_stage == 1)
109  {
110  // In the standard RK notation, the stage 1 residual is given by:
111  //
112  // R := (Y_1 - y_n)/dt - alpha*f(t_n + alpha*dt, Y_1) = 0
113  //
114  // where:
115  // .) f(t_n + alpha*dt, Y_1) corresponds to the residuals of the
116  // non-time kernels. We'll save this as "_residual_stage1" to use
117  // later.
118  // .) (Y_1 - y_n)/dt corresponds to the residual of the time kernels.
119  // .) The minus sign in front of alpha is already "baked in" to
120  // the non-time residuals, so it does not appear here.
123 
124  residual.add(1., _Re_time);
125  residual.add(_alpha, _residual_stage1);
126  residual.close();
127  }
128  else if (_stage == 2)
129  {
130  // In the standard RK notation, the stage 2 residual is given by:
131  //
132  // R := (Y_2 - y_n)/dt - (1-alpha)*f(t_n + alpha*dt, Y_1) - alpha*f(t_n + dt, Y_2) = 0
133  //
134  // where:
135  // .) f(t_n + alpha*dt, Y_1) has already been saved as _residual_stage1.
136  // .) f(t_n + dt, Y_2) will now be saved as "_residual_stage2".
137  // .) (Y_2 - y_n)/dt corresponds to the residual of the time kernels.
138  // .) The minus signs are once again "baked in" to the non-time
139  // residuals, so they do not appear here.
140  //
141  // The solution at the end of stage 2, i.e. Y_2, is also the final solution.
144 
145  residual.add(1., _Re_time);
146  residual.add(1. - _alpha, _residual_stage1);
147  residual.add(_alpha, _residual_stage2);
148  residual.close();
149  }
150  else
151  mooseError("LStableDirk2::postResidual(): Member variable _stage can only have values 1 or 2.");
152 }
NonlinearSystemBase & _nl
virtual void initPetscOutputAndSomeSolverSettings()
Reinitialize PETSc output for proper linear/nonlinear iteration display.
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: LStableDirk2.C:55
virtual NumericVector< Number > * solutionUDot()=0
virtual Real & time() const
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: LStableDirk2.C:38
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
FEProblemBase & _fe_problem
void mooseInfo(Args &&... args) const
LStableDirk2(const InputParameters &parameters)
Definition: LStableDirk2.C:26
DualNumber< Real, DNDerivativeType, true > DualReal
Definition: DualReal.h:49
NumericVector< Number > & _Re_non_time
residual vector for non-time contributions
const Real _alpha
Definition: LStableDirk2.h:70
static InputParameters validParams()
Definition: LStableDirk2.C:18
SystemBase & _sys
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: LStableDirk2.C:63
MetaPhysicL::DualNumber< V, D, asd > sqrt(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: ADReal.h:35
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
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
Definition: LStableDirk2.C:106
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
Definition: LStableDirk2.h:75
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.
unsigned int _stage
Indicates the current stage (1 or 2).
Definition: LStableDirk2.h:61
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
NumericVector< Number > & _Re_time
residual vector for time contributions
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
unsigned int _n_nonlinear_iterations
Total number of nonlinear iterations over all stages of the time step.
NumericVector< Number > & _residual_stage2
Buffer to store non-time residual from second stage solve.
Definition: LStableDirk2.h:67
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.
Second order diagonally implicit Runge Kutta method (Dirk) with two stages.
Definition: LStableDirk2.h:38
virtual void add(const numeric_index_type i, const Number value)=0
registerMooseObject("MooseApp", LStableDirk2)
NumericVector< Number > & _residual_stage1
Buffer to store non-time residual from first stage solve.
Definition: LStableDirk2.h:64
static InputParameters validParams()
uint8_t dof_id_type