https://mooseframework.inl.gov
ADKernel.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 "ADKernel.h"
11 #include "Assembly.h"
12 #include "MooseVariable.h"
13 #include "Problem.h"
14 #include "SubProblem.h"
15 #include "NonlinearSystemBase.h"
16 #include "ADUtils.h"
17 #include "FEProblemBase.h"
18 
19 // libmesh includes
20 #include "libmesh/threads.h"
21 
22 template <typename T>
25 {
26  auto params = KernelBase::validParams();
27  if (std::is_same<T, Real>::value)
28  params.registerBase("Kernel");
29  else if (std::is_same<T, RealVectorValue>::value)
30  params.registerBase("VectorKernel");
31  else
32  ::mooseError("unsupported ADKernelTempl specialization");
34  return params;
35 }
36 
37 template <typename T>
39  : KernelBase(parameters),
40  MooseVariableInterface<T>(this,
41  false,
42  "variable",
46  ADFunctorInterface(this),
47  _var(*this->mooseVariable()),
48  _test(_var.phi()),
49  _grad_test(_var.adGradPhi()),
50  _regular_grad_test(_var.gradPhi()),
51  _u(_var.adSln()),
52  _grad_u(_var.adGradSln()),
53  _ad_JxW(_assembly.adJxW()),
54  _ad_coord(_assembly.adCoordTransformation()),
55  _ad_q_point(_assembly.adQPoints()),
56  _phi(_assembly.phi(_var)),
57  _grad_phi(_assembly.template adGradPhi<T>(_var)),
58  _regular_grad_phi(_assembly.gradPhi(_var)),
59  _my_elem(nullptr)
60 {
62 
64  _save_in.resize(_save_in_strings.size());
65  _diag_save_in.resize(_diag_save_in_strings.size());
66 
67  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
68  {
70 
72  paramError("save_in", "cannot use solution variable as save-in variable");
73 
74  if (var->feType() != _var.feType())
75  paramError(
76  "save_in",
77  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
79 
80  _save_in[i] = var;
83  }
84 
85  _has_save_in = _save_in.size() > 0;
86 
87  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
88  {
90 
92  paramError("diag_save_in", "cannot use solution variable as diag save-in variable");
93 
94  if (var->feType() != _var.feType())
95  paramError(
96  "diag_save_in",
97  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
99 
100  _diag_save_in[i] = var;
103  }
104 
105  _has_diag_save_in = _diag_save_in.size() > 0;
106 }
107 
108 template <typename T>
109 void
111 {
112  precalculateResidual();
113 
114  if (_residuals_nonad.size() != _test.size())
115  _residuals_nonad.resize(_test.size(), 0);
116  for (auto & r : _residuals_nonad)
117  r = 0;
118 
119  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
120  {
121  const auto jxw_c = _JxW[_qp] * _coord[_qp];
122  for (_i = 0; _i < _test.size(); _i++)
123  _residuals_nonad[_i] += jxw_c * raw_value(computeQpResidual());
124  }
125 
126  addResiduals(_assembly, _residuals_nonad, _var.dofIndices(), _var.scalingFactor());
127 
128  if (_has_save_in)
129  for (unsigned int i = 0; i < _save_in.size(); i++)
130  _save_in[i]->sys().solution().add_vector(_residuals_nonad.data(), _save_in[i]->dofIndices());
131 }
132 
133 template <typename T>
134 void
136 {
137  if (_residuals.size() != _test.size())
138  _residuals.resize(_test.size(), 0);
139  for (auto & r : _residuals)
140  r = 0;
141 
142  precalculateResidual();
143  if (_use_displaced_mesh)
144  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
145  {
146  _r = _ad_JxW[_qp];
147  _r *= _ad_coord[_qp];
148  for (_i = 0; _i < _test.size(); _i++)
149  _residuals[_i] += _r * computeQpResidual();
150  }
151  else
152  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
153  {
154  const auto jxw_c = _JxW[_qp] * _coord[_qp];
155  for (_i = 0; _i < _test.size(); _i++)
156  _residuals[_i] += jxw_c * computeQpResidual();
157  }
158 }
159 
160 template <typename T>
161 void
163 {
164  computeADJacobian();
165 
166  if (_has_diag_save_in && !_sys.computingScalingJacobian())
167  mooseError("_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use "
168  "the tagging system instead.");
169 }
170 
171 template <typename T>
172 void
174 {
175  computeResidualsForJacobian();
176  addJacobian(_assembly, _residuals, dofIndices(), _var.scalingFactor());
177 }
178 
179 template <typename T>
180 void
182 {
183  _my_elem = nullptr;
184 }
185 
186 template <typename T>
187 void
189 {
190  if (_my_elem != _current_elem)
191  {
192  computeADJacobian();
193  _my_elem = _current_elem;
194  }
195 }
196 
197 template <typename T>
198 void
200 {
201 }
202 
203 template <typename T>
204 void
206 {
207  computeResidualsForJacobian();
208  addResidualsAndJacobian(_assembly, _residuals, _var.dofIndices(), _var.scalingFactor());
209 }
210 
211 template class ADKernelTempl<Real>;
212 template class ADKernelTempl<RealVectorValue>;
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ADKernel.C:199
VarFieldType
Definition: MooseTypes.h:722
void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADKernel.C:110
const libMesh::FEType & feType() const
Get the type of finite element object.
virtual void haveADObjects(bool have_ad_objects)
Method for setting whether we have any ad objects.
Definition: SubProblem.h:767
std::vector< MooseVariableFEBase * > _save_in
Definition: KernelBase.h:65
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
std::vector< AuxVariableName > _diag_save_in_strings
Definition: KernelBase.h:71
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:333
static InputParameters validParams()
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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
THREAD_ID _tid
The thread ID for this kernel.
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
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:69
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
SystemBase & _sys
Reference to the EquationSystem object.
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: ADKernel.C:205
MooseVariableFE< T > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: KernelBase.h:70
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:715
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: KernelBase.h:64
FEProblemBase & _fe_problem
Reference to this kernel&#39;s FEProblemBase.
virtual MooseVariable & getStandardVariable(const THREAD_ID tid, const std::string &var_name)=0
Returns the variable reference for requested MooseVariable which may be in any system.
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
static InputParameters validParams()
Definition: ADKernel.C:24
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1149
std::vector< AuxVariableName > _save_in_strings
Definition: KernelBase.h:66
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void addVariableToZeroOnResidual(std::string var_name)
Adds this variable to the list of variables to be zeroed during each residual evaluation.
Definition: SystemBase.C:174
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.
virtual void addVariableToZeroOnJacobian(std::string var_name)
Adds this variable to the list of variables to be zeroed during each Jacobian evaluation.
Definition: SystemBase.C:180
MooseVariableFE< T > & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: ADKernel.h:80
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
static InputParameters validParams()
Definition: KernelBase.C:21
SystemBase & sys()
Get the system this variable is part of.
void computeOffDiagJacobian(unsigned int) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
Definition: ADKernel.C:188