https://mooseframework.inl.gov
LStableDirk4.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 "TimeIntegrator.h"
13 
52 {
53 public:
55 
57 
58  virtual int order() override { return 4; }
59  virtual void computeTimeDerivatives() override;
60  virtual void computeADTimeDerivatives(ADReal & ad_u_dot,
61  const dof_id_type & dof,
62  ADReal & ad_u_dotdot) const override;
63  virtual void solve() override;
64  virtual void postResidual(NumericVector<Number> & residual) override;
65  virtual bool overridesSolve() const override { return true; }
66 
67 protected:
71  template <typename T, typename T2>
72  void computeTimeDerivativeHelper(T & u_dot, const T2 & u_old) const;
73 
74  // Indicates the current stage.
75  unsigned int _stage;
76 
77  // The number of stages in the method. According to S9.4.2/4 of the
78  // standard, we can specify a constant initializer like this for
79  // integral types, it does not have to appear outside the class
80  // definition.
81  static const unsigned int _n_stages = 5;
82 
83  // Store pointers to the various stage residuals
85 
86  // Butcher tableau "C" parameters derived from _gamma
87  static const Real _c[_n_stages];
88 
89  // Butcher tableau "A" values derived from _gamma. We only use the
90  // lower triangle of this.
91  static const Real _a[_n_stages][_n_stages];
92 };
93 
94 template <typename T, typename T2>
95 void
96 LStableDirk4::computeTimeDerivativeHelper(T & u_dot, const T2 & u_old) const
97 {
98  u_dot -= u_old;
99  u_dot *= 1. / _dt;
100 }
Fourth-order diagonally implicit Runge Kutta method (Dirk) with five stages.
Definition: LStableDirk4.h:51
static const Real _a[_n_stages][_n_stages]
Definition: LStableDirk4.h:91
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: LStableDirk4.C:54
LStableDirk4(const InputParameters &parameters)
Definition: LStableDirk4.C:38
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _stage
Definition: LStableDirk4.h:75
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
Real & _dt
The current time step size.
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the NonLinearTimeIntegratorInterface called immediately after the residuals are computed ...
Definition: LStableDirk4.C:126
virtual int order() override
Definition: LStableDirk4.h:58
static const Real _c[_n_stages]
Definition: LStableDirk4.h:87
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
NumericVector< Number > * _stage_residuals[_n_stages]
Definition: LStableDirk4.h:84
static InputParameters validParams()
Definition: LStableDirk4.C:20
static const unsigned int _n_stages
Definition: LStableDirk4.h:81
const InputParameters & parameters() const
Get the parameters of the object.
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: LStableDirk4.C:79
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
Definition: LStableDirk4.C:71
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:96
uint8_t dof_id_type
virtual bool overridesSolve() const override
Definition: LStableDirk4.h:65