https://mooseframework.inl.gov
ADDGKernel.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 "ADDGKernel.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 
23 {
25  params.addClassDescription(
26  "Base class for all DG kernels making use of automatic differentiation");
27  return params;
28 }
29 
31  : DGKernelBase(parameters),
34  _var(*mooseVariable()),
35  _phi(_assembly.phiFace(_var)),
36  _grad_phi(_assembly.gradPhiFace(_var)),
37 
38  _test(_var.phiFace()),
39  _grad_test(_var.gradPhiFace()),
40 
41  _phi_neighbor(_assembly.phiFaceNeighbor(_var)),
42  _grad_phi_neighbor(_assembly.gradPhiFaceNeighbor(_var)),
43 
44  _test_neighbor(_var.phiFaceNeighbor()),
45  _grad_test_neighbor(_var.gradPhiFaceNeighbor()),
46 
47  _u(_var.adSln()),
48  _grad_u(_var.adGradSln()),
49  _u_neighbor(_var.adSlnNeighbor()),
50  _grad_u_neighbor(_var.adGradSlnNeighbor())
51 {
53 
55 
56  _save_in.resize(_save_in_strings.size());
57  _diag_save_in.resize(_diag_save_in_strings.size());
58 
59  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
60  {
65 
67  mooseError("Trying to use solution variable " + _save_in_strings[i] +
68  " as a save_in variable in " + name());
69 
70  if (var->feType() != _var.feType())
71  paramError(
72  "save_in",
73  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
75 
76  _save_in[i] = var;
79  }
80 
81  _has_save_in = _save_in.size() > 0;
82 
83  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
84  {
89 
91  mooseError("Trying to use solution variable " + _diag_save_in_strings[i] +
92  " as a diag_save_in variable in " + name());
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 void
110 {
111  bool is_elem;
112  if (type == Moose::Element)
113  is_elem = true;
114  else
115  is_elem = false;
116 
117  const VariableTestValue & test_space = is_elem ? _test : _test_neighbor;
118 
119  if (is_elem)
121  else
123 
124  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
125  for (_i = 0; _i < test_space.size(); _i++)
127 
129 
130  if (_has_save_in)
131  for (const auto & var : _save_in)
132  {
133  const std::vector<dof_id_type> & dof_indices =
134  is_elem ? var->dofIndices() : var->dofIndicesNeighbor();
135  var->sys().solution().add_vector(_local_re, dof_indices);
136  }
137 }
138 
139 void
141 {
142  // AD only needs to do one computation for one variable because it does the derivatives all at
143  // once
144  if (!excludeBoundary())
145  {
148  }
149 }
150 
151 void
153 {
155  "With AD you should need one call per side");
156 
158 
159  std::vector<ADReal> residuals(test_space.size(), 0);
160 
161  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
162  for (_i = 0; _i < test_space.size(); _i++)
163  residuals[_i] +=
164  _JxW[_qp] * _coord[_qp] *
166 
168  residuals,
170  _var.scalingFactor());
171 
172  if (_has_diag_save_in)
173  {
174  unsigned int rows = _local_ke.m();
175  DenseVector<Number> diag(rows);
176  for (unsigned int i = 0; i < rows; i++)
177  diag(i) = _local_ke(i, i);
178 
179  for (const auto & var : _diag_save_in)
180  {
182  var->sys().solution().add_vector(diag, var->dofIndices());
183  else
184  var->sys().solution().add_vector(diag, var->dofIndicesNeighbor());
185  }
186  }
187 }
188 
189 void
190 ADDGKernel::computeOffDiagJacobian(const unsigned int jvar_num)
191 {
192  // AD only needs to do one computation for one variable because it does the derivatives all at
193  // once
194  if (!excludeBoundary() && jvar_num == _var.number())
195  {
196  const auto & jvar = getVariable(jvar_num);
199  }
200 }
201 
202 void
204 {
206  "With AD you should need one call per side");
207 
209 
210  std::vector<ADReal> residuals(test_space.size(), 0);
211 
212  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
213  for (_i = 0; _i < test_space.size(); _i++)
214  residuals[_i] +=
215  _JxW[_qp] * _coord[_qp] *
217 
219  residuals,
221  _var.scalingFactor());
222 }
VarFieldType
Definition: MooseTypes.h:770
const QBase *const & _qrule
Quadrature rule.
Definition: DGKernelBase.h:107
ADDGKernel(const InputParameters &parameters)
Definition: ADDGKernel.C:30
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernelBase.C:26
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
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
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:439
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
unsigned int number() const
Get variable number coming from libMesh.
std::vector< AuxVariableName > _save_in_strings
Definition: DGKernelBase.h:118
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...
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: DGKernelBase.h:121
unsigned int m() const
This class provides an interface for common operations on field variables of both FE and FV types wit...
static InputParameters validParams()
Definition: ADDGKernel.C:22
unsigned int _i
Definition: DGKernelBase.h:136
DGResidualType
Definition: MooseTypes.h:792
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: DGKernelBase.h:122
THREAD_ID _tid
The thread ID for this kernel.
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
SystemBase & _sys
Reference to the EquationSystem object.
unsigned int _qp
Definition: DGKernelBase.h:134
std::vector< MooseVariableFEBase * > _save_in
Definition: DGKernelBase.h:117
MooseVariableFE< Real > * 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.
Enhances MooseVariableInterface interface provide values from neighbor elements.
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
void computeElemNeighResidual(Moose::DGResidualType type) override final
Computes the residual for this element or the neighbor.
Definition: ADDGKernel.C:109
const VariableTestValue & _test_neighbor
Side test function.
Definition: ADDGKernel.h:50
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:353
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
std::vector< AuxVariableName > _diag_save_in_strings
Definition: DGKernelBase.h:123
void computeOffDiagJacobian(unsigned int jvar) override final
Computes d-residual / d-jvar...
Definition: ADDGKernel.C:190
bool excludeBoundary() const
Check current element if it contains broken boundary.
Definition: DGKernelBase.C:184
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:93
const MooseArray< Real > & _coord
Coordinate transform mainly for curvilinear coordinates.
Definition: DGKernelBase.h:111
virtual const MooseVariableFieldBase & getVariable(const THREAD_ID tid, const std::string &var_name, Moose::VarKindType expected_var_type=Moose::VarKindType::VAR_ANY, Moose::VarFieldType expected_var_field_type=Moose::VarFieldType::VAR_FIELD_ANY) const =0
Returns the variable reference for requested variable which must be of the expected_var_type (Nonline...
virtual bool hasVariable(const std::string &var_name) const
Query a system for a variable.
Definition: SystemBase.C:851
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
MooseVariable & _var
Variable this kernel operates on.
Definition: ADDGKernel.h:36
DGJacobianType
Definition: MooseTypes.h:798
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
void prepareVectorTagNeighbor(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags for DG and interface kernels...
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
const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
const MooseArray< Real > & _JxW
Jacobian det times quadrature weighting on quadrature points.
Definition: DGKernelBase.h:109
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
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
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void computeElemNeighJacobian(Moose::DGJacobianType type) override final
Computes the element/neighbor-element/neighbor Jacobian.
Definition: ADDGKernel.C:152
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.
virtual ADReal computeQpResidual(Moose::DGResidualType type)=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar) override final
Computes the element-element off-diagonal Jacobian.
Definition: ADDGKernel.C:203
SystemBase & sys()
Get the system this variable is part of.
const VariableTestValue & _test
test functions
Definition: ADDGKernel.h:42
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: DGKernelBase.h:116
void computeJacobian() override final
Computes the jacobian for the current side.
Definition: ADDGKernel.C:140