https://mooseframework.inl.gov
ADKernel.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 "KernelBase.h"
13 #include "ADFunctorInterface.h"
14 #include "MooseVariableInterface.h"
15 
16 // forward declarations
17 template <typename>
19 
22 
23 template <typename T>
24 class ADKernelTempl : public KernelBase, public MooseVariableInterface<T>, public ADFunctorInterface
25 {
26 public:
28 
30 
31  void jacobianSetup() override;
32 
33  const MooseVariableFE<T> & variable() const override { return _var; }
34 
35 protected:
36  void computeJacobian() override;
37  void computeResidualAndJacobian() override;
38  void computeOffDiagJacobian(unsigned int) override;
39  void computeOffDiagJacobianScalar(unsigned int jvar) override;
40 
50  virtual const std::vector<dof_id_type> & dofIndices() const { return _var.dofIndices(); }
51 
52 protected:
53  // See KernelBase base for documentation of these overridden methods
54  void computeResidual() override;
55 
63  virtual void computeResidualsForJacobian();
64 
66  virtual ADReal computeQpResidual() = 0;
67 
71  virtual Real computeQpJacobian() { return 0; }
72 
77  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
78 
81 
84 
87 
88  // gradient of the test function without possible displacement derivatives
90 
93 
96 
99 
102 
105 
108 
110  std::vector<Real> _residuals_nonad;
111  std::vector<ADReal> _residuals;
112 
115 
118 
119 private:
123  void computeADJacobian();
124 
125  const Elem * _my_elem;
126 };
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ADKernel.C:199
void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADKernel.C:110
const ADTemplateVariableTestGradient< T > & _grad_test
gradient of the test function
Definition: ADKernel.h:86
typename OutputTools< typename Moose::ADType< T >::type >::VariablePhiGradient ADTemplateVariablePhiGradient
Definition: MooseTypes.h:637
std::vector< Real > _residuals_nonad
Definition: ADKernel.h:110
virtual const std::vector< dof_id_type > & dofIndices() const
Just as we allow someone deriving from this to modify the _residuals data member in their computeResi...
Definition: ADKernel.h:50
typename OutputTools< T >::VariableTestValue ADTemplateVariableTestValue
Definition: MooseTypes.h:623
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: ADKernel.h:71
const ADTemplateVariablePhiValue< T > & _phi
The current shape functions.
Definition: ADKernel.h:107
typename OutputTools< T >::VariablePhiValue ADTemplateVariablePhiValue
Definition: MooseTypes.h:625
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:603
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const ADTemplateVariableValue< T > & _u
Holds the solution at current quadrature points.
Definition: ADKernel.h:92
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: ADKernel.C:162
ADKernelTempl(const InputParameters &parameters)
Definition: ADKernel.C:38
const ADTemplateVariableTestValue< T > & _test
the current test function
Definition: ADKernel.h:83
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
void jacobianSetup() override
Gets called just before the Jacobian is computed and before this object is asked to do its job...
Definition: ADKernel.C:181
const MooseVariableFE< T > & variable() const override
Returns the variable that this object operates on.
Definition: ADKernel.h:33
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
void computeADJacobian()
compute all the Jacobian entries
Definition: ADKernel.C:173
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: ADKernel.C:205
std::vector< ADReal > _residuals
Definition: ADKernel.h:111
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
virtual ADReal computeQpResidual()=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
const ADTemplateVariablePhiGradient< T > & _grad_phi
The current gradient of the shape functions.
Definition: ADKernel.h:114
typename OutputTools< typename Moose::ADType< T >::type >::VariableTestGradient ADTemplateVariableTestGradient
Definition: MooseTypes.h:631
const OutputTools< T >::VariableTestGradient & _regular_grad_test
Definition: ADKernel.h:89
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:606
const MooseArray< ADPoint > & _ad_q_point
The ad version of q_point.
Definition: ADKernel.h:104
const MooseArray< ADReal > & _ad_JxW
The ad version of JxW.
Definition: ADKernel.h:98
static InputParameters validParams()
Definition: ADKernel.C:24
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: ADKernel.h:77
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const ADTemplateVariableGradient< T > & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: ADKernel.h:95
const Elem * _my_elem
Definition: ADKernel.h:125
virtual void computeResidualsForJacobian()
compute the _residuals member for filling the Jacobian.
Definition: ADKernel.C:135
Interface for objects that need to get values of MooseVariables.
ADReal _r
Definition: ADKernel.h:109
const OutputTools< T >::VariablePhiGradient & _regular_grad_phi
The current gradient of the shape functions without possible displacement derivatives.
Definition: ADKernel.h:117
MooseVariableFE< T > & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: ADKernel.h:80
const InputParameters & parameters() const
Get the parameters of the object.
const MooseArray< ADReal > & _ad_coord
The ad version of coord.
Definition: ADKernel.h:101
void computeOffDiagJacobian(unsigned int) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
Definition: ADKernel.C:188