https://mooseframework.inl.gov
ExplicitEuler.C
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 #include "ExplicitEuler.h"
11 #include "NonlinearSystem.h"
12 #include "FEProblem.h"
13 
15 
18 {
20  params.addClassDescription("Time integration using the explicit Euler method.");
21  return params;
22 }
23 
24 ExplicitEuler::ExplicitEuler(const InputParameters & parameters) : TimeIntegrator(parameters) {}
25 
26 void
28 {
29  if (_dt == _dt_old)
31  else
33 }
34 
35 void
37 {
38  if (!_sys.solutionUDot())
39  mooseError("ExplicitEuler: Time derivative of solution (`u_dot`) is not stored. Please set "
40  "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
41 
43  u_dot = *_solution;
45  u_dot.close();
47 }
48 
49 void
51  const dof_id_type & dof,
52  ADReal & /*ad_u_dotdot*/) const
53 {
55 }
56 
57 void
59 {
60  residual += *_Re_time;
61  residual += *_Re_non_time;
62  residual.close();
63 }
static InputParameters validParams()
Definition: ExplicitEuler.C:17
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: ExplicitEuler.C:50
FEProblemBase & _fe_problem
Reference to the problem.
ExplicitEuler(const InputParameters &parameters)
Definition: ExplicitEuler.C:24
SystemBase & _sys
Reference to the system this time integrator operates on.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeTimeDerivatives() override
Computes the time derivative and the Jacobian of the time derivative.
Definition: ExplicitEuler.C:36
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the NonLinearTimeIntegratorInterface called immediately after the residuals are computed ...
Definition: ExplicitEuler.C:58
Real & _dt
The current time step size.
Real & _dt_old
The previous time step size.
NumericVector< Number > * _Re_time
residual vector for time contributions
virtual NumericVector< Number > * solutionUDot()
Definition: SystemBase.h:252
virtual void close()=0
void computeTimeDerivativeHelper(T &u_dot, const T2 &u_old) const
Helper function that actually does the math for computing the time derivative.
Definition: ExplicitEuler.h:43
const NumericVector< Number > *const & _solution
Explicit Euler time integrator.
Definition: ExplicitEuler.h:17
Base class for time integrators.
void setConstJacobian(bool state)
Set flag that Jacobian is constant (for optimization purposes)
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
const NumericVector< Number > & _solution_old
virtual void preSolve() override
Definition: ExplicitEuler.C:27
registerMooseObject("MooseApp", ExplicitEuler)
NumericVector< Number > * _Re_non_time
residual vector for non-time contributions
static InputParameters validParams()
uint8_t dof_id_type
void computeDuDotDu()
Compute _du_dot_du.