Go to the documentation of this file.
32 {.125, .125, .25, 0, 0},
33 {-1.5, .75, 1.5, .25, 0},
34 {0, 1. / 6, 2. / 3, -1. / 12, .25}};
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.");
43 for (
unsigned int stage = 0; stage <
_n_stages; ++stage)
45 std::ostringstream oss;
46 oss <<
"residual_stage" << stage + 1;
58 mooseError(
"LStableDirk4: Time derivative of solution (`u_dot`) is not stored. Please set "
59 "uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
86 for (
unsigned int current_stage = 1; current_stage <=
_n_stages; ++current_stage)
121 mooseError(
"LStableDirk4::postResidual(): Member variable _stage can only have values 1-",
142 for (
unsigned int j = 0; j <
_stage; ++j)
unsigned int _n_linear_iterations
Total number of linear iterations over all stages of the time step.
virtual NumericVector< Number > & addVector(const std::string &vector_name, const bool project, const ParallelType type)
Adds a solution length vector to the system.
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
virtual Real & timeOld() const
void mooseInfo(Args &&... args) const
void mooseError(Args &&... args) const
Fourth-order diagonally implicit Runge Kutta method (Dirk) with five stages.
static InputParameters validParams()
NumericVector< Number > & _Re_non_time
residual vector for non-time contributions
Base class for time integrators.
virtual Real & time() const
unsigned int getNumLinearIterationsLastSolve() const
Gets the number of linear iterations in the most recent solve.
const NumericVector< Number > *const & _solution
solution vectors
DualNumber< Real, DNDerivativeType > DualReal
static const Real _a[_n_stages][_n_stages]
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
defineLegacyParams(LStableDirk4)
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
static InputParameters validParams()
virtual void computeADTimeDerivatives(DualReal &ad_u_dot, const dof_id_type &dof) const override
method for computing local automatic differentiation time derivatives
FEProblemBase & _fe_problem
NumericVector< Number > & _Re_time
residual vector for time contributions
NumericVector< Number > * _stage_residuals[_n_stages]
unsigned int _n_nonlinear_iterations
Total number of nonlinear iterations over all stages of the time step.
unsigned int getNumNonlinearIterationsLastSolve() const
Gets the number of nonlinear iterations in the most recent solve.
static const unsigned int _n_stages
virtual System & system() override
Get the reference to the libMesh system.
registerMooseObject("MooseApp", LStableDirk4)
Real & _du_dot_du
Derivative of time derivative with respect to current solution: .
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
const NumericVector< Number > & _solution_old
virtual bool converged() override
static const Real _c[_n_stages]
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
NonlinearSystemBase & _nl
NonlinearSystemBase & getNonlinearSystemBase()
void initPetscOutput()
Reinitialize petsc output for proper linear/nonlinear iteration display.
virtual NumericVector< Number > * solutionUDot()=0
LStableDirk4(const InputParameters ¶meters)