https://mooseframework.inl.gov
ExplicitSSPRungeKutta.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 #pragma once
11 
12 #include "ExplicitTimeIntegrator.h"
13 
18 {
19 public:
21 
23 
24  virtual void computeTimeDerivatives() override;
25  virtual void computeADTimeDerivatives(ADReal & ad_u_dot,
26  const dof_id_type & dof,
27  ADReal & ad_u_dotdot) const override;
28  virtual void solve() override;
29  virtual void postResidual(NumericVector<Number> & residual) override;
30  virtual int order() override { return _order; }
31  virtual bool overridesSolve() const override { return true; }
32 
33 protected:
39  bool solveStage();
40 
41  virtual Real duDotDuCoeff() const override;
42 
44  const MooseEnum & _order;
45 
47  unsigned int _n_stages;
49  std::vector<std::vector<Real>> _a;
51  std::vector<Real> _b;
53  std::vector<Real> _c;
55  unsigned int _stage;
57  std::vector<const NumericVector<Number> *> _solution_stage;
58 
65 };
virtual void computeADTimeDerivatives(ADReal &ad_u_dot, const dof_id_type &dof, ADReal &ad_u_dotdot) const override
method for computing local automatic differentiation time derivatives
std::vector< Real > _c
Runge-Kutta "c" coefficient vector.
unsigned int _n_stages
Number of stages.
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
virtual int order() override
std::vector< Real > _b
Runge-Kutta "b" coefficient vector.
NumericVector< Number > * _solution_intermediate_stage
Solution vector for intermediate stage.
virtual Real duDotDuCoeff() const override
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
const MooseEnum & _order
Order of time integration.
bool solveStage()
Solves a stage of the time integrator.
std::vector< std::vector< Real > > _a
Runge-Kutta "a" coefficient matrix.
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the NonLinearTimeIntegratorInterface called immediately after the residuals are computed ...
virtual bool overridesSolve() const override
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
ExplicitSSPRungeKutta(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
const InputParameters & parameters() const
Get the parameters of the object.
NumericVector< Number > * _tmp_mass_solution_product
Temporary mass-matrix/solution vector product.
unsigned int _stage
Current stage.
NumericVector< Number > * _tmp_solution
Temporary solution vector.
Explicit strong stability preserving Runge-Kutta methods.
Base class for explicit time integrators that are implemented without using a nonlinear solver...
std::vector< const NumericVector< Number > * > _solution_stage
Pointer to solution vector for each stage.
uint8_t dof_id_type