https://mooseframework.inl.gov
ADArrayKernel.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 "KernelBase.h"
13 #include "MooseVariableInterface.h"
14 #include "ADFunctorInterface.h"
15 
19 class ADArrayKernel : public KernelBase,
20  public MooseVariableInterface<RealEigenVector>,
21  public ADFunctorInterface
22 {
23 public:
25 
27 
28  virtual void computeResidual() override;
29  virtual void computeJacobian() override;
30  virtual void computeOffDiagJacobian(unsigned int jvar) override;
31  virtual void jacobianSetup() override;
32 
33  virtual const ArrayMooseVariable & variable() const override { return _var; }
34 
35 protected:
40  virtual void computeQpResidual(ADRealEigenVector & residual) = 0;
41 
45  virtual void initQpResidual() {}
46 
49 
52 
56 
59 
62 
65 
68 
70  const unsigned int _count;
71 
72  // Pre-allocated vector for storing the AD array residual of the kernel
73  std::vector<ADReal> _local_ad_re;
74 
75 private:
78 
80  const Elem * _my_elem = nullptr;
81 };
Eigen::Matrix< ADReal, Eigen::Dynamic, 1 > ADRealEigenVector
Definition: MooseTypes.h:148
ADArrayKernel(const InputParameters &parameters)
Definition: ADArrayKernel.C:29
virtual void initQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
Definition: ADArrayKernel.h:45
std::vector< ADReal > _local_ad_re
Definition: ADArrayKernel.h:73
static InputParameters validParams()
Definition: ADArrayKernel.C:22
Class for stuff related to variables.
Definition: Adaptivity.h:31
OutputTools< RealEigenVector >::VariablePhiValue ArrayVariablePhiValue
Definition: MooseTypes.h:384
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeQpResidual(ADRealEigenVector &residual)=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point, to be filled in residual.
const ArrayVariablePhiValue & _phi
the current shape functions
Definition: ADArrayKernel.h:58
const ArrayVariableTestValue & _test
the current test function
Definition: ADArrayKernel.h:51
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
const Elem * _my_elem
Cache variable to prevent multiple invocations of Jacobian computation for one element (recall that A...
Definition: ADArrayKernel.h:80
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
ADRealEigenVector _work_vector
Work vector for residual and diag jacobian.
Definition: ADArrayKernel.h:77
Base class for array variable (equation) kernels using automatic differentiation. ...
Definition: ADArrayKernel.h:19
std::vector< std::vector< Eigen::Map< RealDIMValue > > > MappedArrayVariablePhiGradient
Definition: MooseTypes.h:386
const ArrayVariablePhiGradient & _grad_phi
gradient of the shape function
Definition: ADArrayKernel.h:61
const ADArrayVariableGradient & _grad_u
Holds the solution gradient at current quadrature points.
Definition: ADArrayKernel.h:67
ArrayMooseVariable & _var
This is an array kernel so we cast to a ArrayMooseVariable.
Definition: ADArrayKernel.h:48
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: ADArrayKernel.C:76
const unsigned int _count
Number of components of the array variable.
Definition: ADArrayKernel.h:70
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADArrayKernel.C:53
forward declarations
const ADArrayVariableValue & _u
Holds the solution at current quadrature points.
Definition: ADArrayKernel.h:64
virtual const ArrayMooseVariable & variable() const override
Returns the variable that this object operates on.
Definition: ADArrayKernel.h:33
Interface for objects that need to get values of MooseVariables.
const ArrayVariableTestGradient & _grad_test
gradient of the test function
Definition: ADArrayKernel.h:54
virtual void jacobianSetup() override
Gets called just before the Jacobian is computed and before this object is asked to do its job...
Definition: ADArrayKernel.C:98
OutputTools< RealEigenVector >::VariableTestGradient ArrayVariableTestGradient
Definition: MooseTypes.h:391
const MappedArrayVariablePhiGradient & _array_grad_test
Definition: ADArrayKernel.h:55