https://mooseframework.inl.gov
ExplicitRK2.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 
55 {
56 public:
58 
60 
61  virtual void preSolve() override;
62  virtual int order() override { return 2; }
63 
64  virtual void computeTimeDerivatives() override;
65  void computeADTimeDerivatives(ADReal & ad_u_dot,
66  const dof_id_type & dof,
67  ADReal & ad_u_dotdot) const override;
68  virtual void solve() override;
69  virtual void postResidual(NumericVector<Number> & residual) override;
70  virtual bool overridesSolve() const override { return true; }
71 
72 protected:
76  template <typename T, typename T2, typename T3>
77  void computeTimeDerivativeHelper(T & u_dot, const T2 & u_old, const T3 & u_older) const;
78 
79  unsigned int _stage;
80 
83 
86 
90  virtual Real a() const = 0;
91  virtual Real b1() const = 0;
92  virtual Real b2() const = 0;
93 };
94 
95 template <typename T, typename T2, typename T3>
96 void
97 ExplicitRK2::computeTimeDerivativeHelper(T & u_dot, const T2 & u_old, const T3 & u_older) const
98 {
99  if (_stage < 3)
100  u_dot -= u_old;
101  else
102  u_dot -= u_older;
103 
104  u_dot *= 1. / _dt;
105 }
virtual bool overridesSolve() const override
Definition: ExplicitRK2.h:70
NumericVector< Number > * _residual_old
Buffer to store non-time residual from the first stage.
Definition: ExplicitRK2.h:82
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Base class for three different explicit second-order Runge-Kutta time integration methods: ...
Definition: ExplicitRK2.h:54
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
Real & _dt
The current time step size.
static InputParameters validParams()
Definition: ExplicitRK2.C:18
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: ExplicitRK2.C:70
virtual Real a() const =0
The method coefficients.
unsigned int _stage
Definition: ExplicitRK2.h:79
virtual int order() override
Definition: ExplicitRK2.h:62
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
virtual Real b1() const =0
virtual void preSolve() override
Definition: ExplicitRK2.C:37
ExplicitRK2(const InputParameters &parameters)
Definition: ExplicitRK2.C:25
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the NonLinearTimeIntegratorInterface called immediately after the residuals are computed ...
Definition: ExplicitRK2.C:117
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: ExplicitRK2.C:46
const InputParameters & parameters() const
Get the parameters of the object.
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: ExplicitRK2.C:62
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old, const T3 &u_older) const
Helper function that actually does the math for computing the time derivative.
Definition: ExplicitRK2.h:97
virtual Real b2() const =0
uint8_t dof_id_type
const NumericVector< Number > & _solution_older
The older solution.
Definition: ExplicitRK2.h:85