https://mooseframework.inl.gov
ADArrayKernel.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 "ADArrayKernel.h"
11 
12 #include "Assembly.h"
13 #include "MooseVariableFE.h"
14 #include "MooseVariableScalar.h"
15 #include "SubProblem.h"
16 #include "NonlinearSystem.h"
17 
18 #include "libmesh/threads.h"
19 #include "libmesh/quadrature.h"
20 
23 {
25  params.registerBase("ADArrayKernel");
26  return params;
27 }
28 
30  : KernelBase(parameters),
32  false,
33  "variable",
36  ADFunctorInterface(this),
37  _var(*mooseVariable()),
38  _test(_var.phi()),
39  _grad_test(_var.gradPhi()),
40  _array_grad_test(_var.arrayGradPhi()),
41  _phi(_assembly.phi(_var)),
42  _grad_phi(_assembly.gradPhi(_var)),
43  _u(_var.adSln()),
44  _grad_u(_var.adGradSln()),
45  _count(_var.count()),
46  _work_vector(_count)
47 {
50 }
51 
52 void
54 {
55 
57 
59  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
60  {
62  for (_i = 0; _i < _test.size(); _i++)
63  {
64  _work_vector.setZero();
66  auto raw_work_vector = MetaPhysicL::raw_value(_work_vector);
67  raw_work_vector *= _JxW[_qp] * _coord[_qp];
68  _assembly.saveLocalArrayResidual(_local_re, _i, _test.size(), raw_work_vector);
69  }
70  }
71 
73 }
74 
75 void
77 {
78  _local_ad_re.resize(_test.size() * _count);
79  for (auto & residual : _local_ad_re)
80  residual = 0;
81 
83  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
84  {
86  for (_i = 0; _i < _test.size(); _i++)
87  {
91  }
92  }
93 
95 }
96 
97 void
99 {
100  _my_elem = nullptr;
101 }
102 
103 void
105 {
106  if (_my_elem != _current_elem)
107  {
108  computeJacobian();
110  }
111 }
VarFieldType
Definition: MooseTypes.h:770
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
virtual void haveADObjects(bool have_ad_objects)
Method for setting whether we have any ad objects.
Definition: SubProblem.h:767
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
static InputParameters validParams()
Definition: ADArrayKernel.C:22
void saveLocalADArray(std::vector< ADReal > &re, unsigned int i, unsigned int ntest, const ADRealEigenVector &v) const
Definition: Assembly.h:1825
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:52
const std::vector< Real > & arrayScalingFactor() const
unsigned int number() const
Get variable number coming from libMesh.
virtual void precalculateResidual()
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:100
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void saveLocalArrayResidual(DenseVector< Number > &re, unsigned int i, unsigned int ntest, const RealEigenVector &v) const
Helper function for assembling residual contriubutions on local quadrature points for an array kernel...
Definition: Assembly.h:1816
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 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
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
MooseVariableFE< RealEigenVector > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided residual derivatives into the Jacobian for the provided dof indices.
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
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 void precalculateJacobian()
ADRealEigenVector _work_vector
Work vector for residual and diag jacobian.
Definition: ADArrayKernel.h:77
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:49
ArrayMooseVariable & _var
This is an array kernel so we cast to a ArrayMooseVariable.
Definition: ADArrayKernel.h:48
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
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
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:55
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADArrayKernel.C:53
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
Interface for objects that need to get values of MooseVariables.
const Elem *const & _current_elem
Current element.
Definition: KernelBase.h:37
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
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
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
static InputParameters validParams()
Definition: KernelBase.C:21
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43