www.mooseframework.org
ImplicitMidpoint.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 "ImplicitMidpoint.h"
11 #include "NonlinearSystem.h"
12 #include "FEProblem.h"
13 #include "PetscSupport.h"
14 
16 
19 {
21  params.addClassDescription("Second-order Runge-Kutta (implicit midpoint) time integration.");
22  return params;
23 }
24 
26  : TimeIntegrator(parameters),
27  _stage(1),
28  _residual_stage1(_nl.addVector("residual_stage1", false, GHOSTED))
29 {
30  mooseInfo("ImplicitMidpoint and other multistage TimeIntegrators are known not to work with "
31  "Materials/AuxKernels that accumulate 'state' and should be used with caution.");
32 }
33 
34 void
36 {
37  // We are multiplying by the method coefficients in postResidual(), so
38  // the time derivatives are of the same form at every stage although
39  // the current solution varies depending on the stage.
40  if (!_sys.solutionUDot())
41  mooseError("ImplicitMidpoint: Time derivative of solution (`u_dot`) is not stored. Please set "
42  "uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
43 
45  u_dot = *_solution;
47  u_dot.close();
48  _du_dot_du = 1. / _dt;
49 }
50 
51 void
53  const dof_id_type & dof,
54  DualReal & /*ad_u_dotdot*/) const
55 {
57 }
58 
59 void
61 {
62  Real time_new = _fe_problem.time();
63  Real time_old = _fe_problem.timeOld();
64  Real time_half = (time_new + time_old) / 2.;
65 
66  // Reset iteration counts
69 
70  // Compute first stage
72  _console << "1st stage" << std::endl;
73  _stage = 1;
74  _fe_problem.time() = time_half;
75  _nl.system().solve();
78 
79  // Abort time step immediately on stage failure - see TimeIntegrator doc page
81  return;
82 
83  // Compute second stage
85  _console << "2nd stage" << std::endl;
86  _stage = 2;
87  _fe_problem.time() = time_new;
88  _nl.system().solve();
91 }
92 
93 void
95 {
96  if (_stage == 1)
97  {
98  // In the standard RK notation, the stage 1 residual is given by:
99  //
100  // R := M*(Y_1 - y_n)/dt - (1/2)*f(t_n + dt/2, Y_1) = 0
101  //
102  // where:
103  // .) M is the mass matrix
104  // .) f(t_n + dt/2, Y_1) is saved in _residual_stage1
105  // .) The minus sign is baked in to the non-time residuals, so it does not appear here.
108 
109  residual.add(1., _Re_time);
110  residual.add(0.5, _Re_non_time);
111  residual.close();
112  }
113  else if (_stage == 2)
114  {
115  // The update step. In the standard RK notation, the update
116  // residual is given by:
117  //
118  // R := M*(y_{n+1} - y_n)/dt - f(t_n + dt/2, Y_1) = 0
119  //
120  // where:
121  // .) M is the mass matrix.
122  // .) f(t_n + dt/2, Y_1) is the residual from stage 1, it has already
123  // been saved as _residual_stage1.
124  // .) The minus signs are "baked in" to the non-time residuals, so
125  // they do not appear here.
126  residual.add(1., _Re_time);
127  residual.add(1., _residual_stage1);
128  residual.close();
129  }
130  else
131  mooseError(
132  "ImplicitMidpoint::postResidual(): _stage = ", _stage, ", only _stage = 1, 2 is allowed.");
133 }
NonlinearSystemBase & _nl
virtual void initPetscOutputAndSomeSolverSettings()
Reinitialize PETSc output for proper linear/nonlinear iteration display.
NumericVector< Number > & _residual_stage1
Buffer to store non-time residual from the first stage.
virtual NumericVector< Number > * solutionUDot()=0
virtual Real & time() const
Second-order Runge-Kutta (implicit midpoint) time integration.
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
static InputParameters validParams()
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...
ImplicitMidpoint(const InputParameters &parameters)
registerMooseObject("MooseApp", ImplicitMidpoint)
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
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 solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
virtual void close()=0
Real & _du_dot_du
Derivative of time derivative with respect to current solution: .
const NumericVector< Number > *const & _solution
solution vectors
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
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.
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
NumericVector< Number > & _Re_time
residual vector for time contributions
unsigned int _stage
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.
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
static InputParameters validParams()
uint8_t dof_id_type