www.mooseframework.org
ExplicitEuler.C
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 #include "ExplicitEuler.h"
11 #include "NonlinearSystem.h"
12 #include "FEProblem.h"
13 
15 
16 template <>
19 {
21 
22  return params;
23 }
24 
25 ExplicitEuler::ExplicitEuler(const InputParameters & parameters) : TimeIntegrator(parameters) {}
26 
27 void
29 {
30  if (_dt == _dt_old)
32  else
34 }
35 
36 void
38 {
39  if (!_sys.solutionUDot())
40  mooseError("ExplicitEuler: Time derivative of solution (`u_dot`) is not stored. Please set "
41  "uDotRequested() to true in FEProblemBase befor requesting `u_dot`.");
42 
43  NumericVector<Number> & u_dot = *_sys.solutionUDot();
44  u_dot = *_solution;
46  u_dot.close();
47 
48  _du_dot_du = 1.0 / _dt;
49 }
50 
51 void
52 ExplicitEuler::computeADTimeDerivatives(DualReal & ad_u_dot, const dof_id_type & dof) const
53 {
55 }
56 
57 void
58 ExplicitEuler::postResidual(NumericVector<Number> & residual)
59 {
60  residual += _Re_time;
61  residual += _Re_non_time;
62  residual.close();
63 }
virtual NumericVector< Number > * solutionUDot()=0
FEProblemBase & _fe_problem
ExplicitEuler(const InputParameters &parameters)
Definition: ExplicitEuler.C:25
NumericVector< Number > & _Re_non_time
residual vector for non-time contributions
SystemBase & _sys
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
DualNumber< Real, NumberArray< AD_MAX_DOFS_PER_ELEM, Real > > DualReal
Definition: DualReal.h:29
void mooseError(Args &&... args) const
Definition: MooseObject.h:147
virtual void computeTimeDerivatives() override
Definition: ExplicitEuler.C:37
virtual void postResidual(NumericVector< Number > &residual) override
Callback to the TimeIntegrator called immediately after the residuals are computed in NonlinearSystem...
Definition: ExplicitEuler.C:58
InputParameters validParams< TimeIntegrator >()
void computeADTimeDerivatives(DualReal &ad_u_dot, const dof_id_type &dof) const override
method for computing local automatic differentiation time derivatives
Definition: ExplicitEuler.C:52
Real & _du_dot_du
solution vector for
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
solution vectors
Explicit Euler time integrator.
Definition: ExplicitEuler.h:22
Base class for time integrators.
void setConstJacobian(bool state)
Set flag that Jacobian is constant (for optimization purposes)
NumericVector< Number > & _Re_time
residual vector for time contributions
InputParameters validParams< ExplicitEuler >()
Definition: ExplicitEuler.C:18
const NumericVector< Number > & _solution_old
virtual void preSolve() override
Definition: ExplicitEuler.C:28
registerMooseObject("MooseApp", ExplicitEuler)