www.mooseframework.org
ADKernel.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 "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 
18 // libmesh includes
19 #include "libmesh/threads.h"
20 
21 template <typename T>
24 {
25  auto params = KernelBase::validParams();
26  if (std::is_same<T, Real>::value)
27  params.registerBase("Kernel");
28  else if (std::is_same<T, RealVectorValue>::value)
29  params.registerBase("VectorKernel");
30  else
31  ::mooseError("unsupported ADKernelTempl specialization");
33  return params;
34 }
35 
36 template <typename T>
38  : KernelBase(parameters),
39  MooseVariableInterface<T>(this,
40  false,
41  "variable",
45  ADFunctorInterface(this),
46  _var(*this->mooseVariable()),
47  _test(_var.phi()),
48  _grad_test(_var.adGradPhi()),
49  _regular_grad_test(_var.gradPhi()),
50  _u(_var.adSln()),
51  _grad_u(_var.adGradSln()),
52  _ad_JxW(_assembly.adJxW()),
53  _ad_coord(_assembly.adCoordTransformation()),
54  _ad_q_point(_assembly.adQPoints()),
55  _phi(_assembly.phi(_var)),
56  _grad_phi(_assembly.template adGradPhi<T>(_var)),
57  _regular_grad_phi(_assembly.gradPhi(_var)),
58  _my_elem(nullptr)
59 {
61 
63  _save_in.resize(_save_in_strings.size());
64  _diag_save_in.resize(_diag_save_in_strings.size());
65 
66  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
67  {
69 
71  paramError("save_in", "cannot use solution variable as save-in variable");
72 
73  if (var->feType() != _var.feType())
74  paramError(
75  "save_in",
76  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
78 
79  _save_in[i] = var;
82  }
83 
84  _has_save_in = _save_in.size() > 0;
85 
86  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
87  {
89 
91  paramError("diag_save_in", "cannot use solution variable as diag save-in variable");
92 
93  if (var->feType() != _var.feType())
94  paramError(
95  "diag_save_in",
96  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
98 
99  _diag_save_in[i] = var;
102  }
103 
104  _has_diag_save_in = _diag_save_in.size() > 0;
105 }
106 
107 template <typename T>
108 void
110 {
111  precalculateResidual();
112 
113  std::vector<Real> residuals(_test.size(), 0);
114 
115  if (_use_displaced_mesh)
116  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
117  for (_i = 0; _i < _test.size(); _i++)
118  residuals[_i] += raw_value(_ad_JxW[_qp] * _ad_coord[_qp] * computeQpResidual());
119  else
120  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
121  for (_i = 0; _i < _test.size(); _i++)
122  residuals[_i] += raw_value(_JxW[_qp] * _coord[_qp] * computeQpResidual());
123 
124  addResiduals(_assembly, residuals, _var.dofIndices(), _var.scalingFactor());
125 
126  if (_has_save_in)
127  for (unsigned int i = 0; i < _save_in.size(); i++)
128  _save_in[i]->sys().solution().add_vector(residuals.data(), _save_in[i]->dofIndices());
129 }
130 
131 template <typename T>
132 void
134 {
135  if (_residuals.size() != _test.size())
136  _residuals.resize(_test.size(), 0);
137  for (auto & r : _residuals)
138  r = 0;
139 
140  precalculateResidual();
141  if (_use_displaced_mesh)
142  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
143  {
144  _r = _ad_JxW[_qp];
145  _r *= _ad_coord[_qp];
146  for (_i = 0; _i < _test.size(); _i++)
147  _residuals[_i] += _r * computeQpResidual();
148  }
149  else
150  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
151  for (_i = 0; _i < _test.size(); _i++)
152  _residuals[_i] += _JxW[_qp] * _coord[_qp] * computeQpResidual();
153 }
154 
155 template <typename T>
156 void
158 {
159  computeADJacobian();
160 
161  if (_has_diag_save_in && !_sys.computingScalingJacobian())
162  mooseError("_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use "
163  "the tagging system instead.");
164 }
165 
166 template <typename T>
167 void
169 {
170  computeResidualsForJacobian();
171  addJacobian(_assembly, _residuals, dofIndices(), _var.scalingFactor());
172 }
173 
174 template <typename T>
175 void
177 {
178  _my_elem = nullptr;
179 }
180 
181 template <typename T>
182 void
184 {
185  if (_my_elem != _current_elem)
186  {
187  computeADJacobian();
188  _my_elem = _current_elem;
189  }
190 }
191 
192 template <typename T>
193 void
195 {
196 }
197 
198 template <typename T>
199 void
201 {
202  computeResidualsForJacobian();
203  addResidualsAndJacobian(_assembly, _residuals, _var.dofIndices(), _var.scalingFactor());
204 }
205 
206 template class ADKernelTempl<Real>;
207 template class ADKernelTempl<RealVectorValue>;
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ADKernel.C:194
VarFieldType
Definition: MooseTypes.h:634
void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADKernel.C:109
virtual void haveADObjects(bool have_ad_objects)
Method for setting whether we have any ad objects.
Definition: SubProblem.h:737
std::vector< MooseVariableFEBase * > _save_in
Definition: KernelBase.h:64
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:24
std::vector< AuxVariableName > _diag_save_in_strings
Definition: KernelBase.h:70
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
static InputParameters validParams()
auto raw_value(const Eigen::Map< T > &in)
Definition: ADReal.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:157
ADKernelTempl(const InputParameters &parameters)
Definition: ADKernel.C:37
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:176
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:68
const FEType & feType() const
Get the type of finite element object.
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:168
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:200
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:22
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: KernelBase.h:69
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:627
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: KernelBase.h:63
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:23
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 ...
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1132
std::vector< AuxVariableName > _save_in_strings
Definition: KernelBase.h:65
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:170
virtual void computeResidualsForJacobian()
compute the _residuals member for filling the Jacobian.
Definition: ADKernel.C:133
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:176
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:183