www.mooseframework.org
LStableDirk2.h
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 #pragma once
11 
12 #include "TimeIntegrator.h"
13 
39 {
40 public:
42 
44 
45  virtual int order() override { return 2; }
46  virtual void computeTimeDerivatives() override;
47  void computeADTimeDerivatives(DualReal & ad_u_dot,
48  const dof_id_type & dof,
49  DualReal & ad_u_dotdot) const override;
50  virtual void solve() override;
51  virtual void postResidual(NumericVector<Number> & residual) override;
52 
53 protected:
57  template <typename T, typename T2>
58  void computeTimeDerivativeHelper(T & u_dot, const T2 & u_old) const;
59 
61  unsigned int _stage;
62 
65 
68 
69  // The parameter of the method, set at construction time and cannot be changed.
70  const Real _alpha;
71 };
72 
73 template <typename T, typename T2>
74 void
75 LStableDirk2::computeTimeDerivativeHelper(T & u_dot, const T2 & u_old) const
76 {
77  u_dot -= u_old;
78  u_dot *= 1. / _dt;
79 }
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
Definition: LStableDirk2.C:55
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: LStableDirk2.C:38
LStableDirk2(const InputParameters &parameters)
Definition: LStableDirk2.C:26
DualNumber< Real, DNDerivativeType, true > DualReal
Definition: DualReal.h:49
const Real _alpha
Definition: LStableDirk2.h:70
static InputParameters validParams()
Definition: LStableDirk2.C:18
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void solve() override
Solves the time step and sets the number of nonlinear and linear iterations.
Definition: LStableDirk2.C:63
virtual int order() override
Definition: LStableDirk2.h:45
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
Definition: LStableDirk2.C:106
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
Definition: LStableDirk2.h:75
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Base class for time integrators.
unsigned int _stage
Indicates the current stage (1 or 2).
Definition: LStableDirk2.h:61
const InputParameters & parameters() const
Get the parameters of the object.
NumericVector< Number > & _residual_stage2
Buffer to store non-time residual from second stage solve.
Definition: LStableDirk2.h:67
Second order diagonally implicit Runge Kutta method (Dirk) with two stages.
Definition: LStableDirk2.h:38
NumericVector< Number > & _residual_stage1
Buffer to store non-time residual from first stage solve.
Definition: LStableDirk2.h:64
uint8_t dof_id_type