www.mooseframework.org
ADKernelStabilized.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 "ADKernelStabilized.h"
11 #include "MathUtils.h"
12 #include "Assembly.h"
13 #include "SystemBase.h"
14 
15 // libmesh includes
16 #include "libmesh/threads.h"
17 
20 
21 template <typename T, ComputeStage compute_stage>
24 {
26 }
27 
28 template <typename T, ComputeStage compute_stage>
30  const InputParameters & parameters)
31  : ADKernelTempl<T, compute_stage>(parameters)
32 {
33 }
34 
35 template <typename T, ComputeStage compute_stage>
36 void
38 {
39  prepareVectorTag(_assembly, _var.number());
40 
41  precalculateResidual();
42  const unsigned int n_test = _grad_test.size();
43 
44  if (_use_displaced_mesh)
45  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
46  {
47  const auto value = precomputeQpStrongResidual() * _ad_JxW[_qp] * _ad_coord[_qp];
48  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
49  _local_re(_i) += _grad_test[_i][_qp] * computeQpStabilization() * value;
50  }
51  else
52  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
53  {
54  const auto value = precomputeQpStrongResidual() * _JxW[_qp] * _coord[_qp];
55  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
56  _local_re(_i) += _regular_grad_test[_i][_qp] * computeQpStabilization() * value;
57  }
58 
59  accumulateTaggedLocalResidual();
60 
61  if (_has_save_in)
62  {
63  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
64  for (unsigned int i = 0; i < _save_in.size(); i++)
65  _save_in[i]->sys().solution().add_vector(_local_re, _save_in[i]->dofIndices());
66  }
67 }
68 
69 template <>
70 void
72 {
73 }
74 
75 template <>
76 void
78 {
79 }
80 
81 template <typename T, ComputeStage compute_stage>
82 void
84 {
85  prepareMatrixTag(_assembly, _var.number(), _var.number());
86 
87  size_t ad_offset = _var.number() * _sys.getMaxVarNDofsPerElem();
88 
89  precalculateResidual();
90 
91  if (_use_displaced_mesh)
92  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
93  {
94  // This will also compute the derivative with respect to all dofs
95  const auto value = precomputeQpStrongResidual() * _ad_JxW[_qp] * _ad_coord[_qp];
96  for (_i = 0; _i < _grad_test.size(); _i++)
97  {
98  const auto residual = _grad_test[_i][_qp] * computeQpStabilization() * value;
99  for (_j = 0; _j < _var.phiSize(); _j++)
100  _local_ke(_i, _j) += residual.derivatives()[ad_offset + _j];
101  }
102  }
103  else
104  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
105  {
106  // This will also compute the derivative with respect to all dofs
107  const auto value = precomputeQpStrongResidual() * _JxW[_qp] * _coord[_qp];
108  for (_i = 0; _i < _grad_test.size(); _i++)
109  {
110  const auto residual = _regular_grad_test[_i][_qp] * computeQpStabilization() * value;
111  for (_j = 0; _j < _var.phiSize(); _j++)
112  _local_ke(_i, _j) += residual.derivatives()[ad_offset + _j];
113  }
114  }
115 
116  accumulateTaggedLocalMatrix();
117 
118  if (_has_diag_save_in)
119  {
120  unsigned int rows = _local_ke.m();
121  DenseVector<Number> diag(rows);
122  for (unsigned int i = 0; i < rows; i++)
123  diag(i) = _local_ke(i, i);
124 
125  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
126  for (unsigned int i = 0; i < _diag_save_in.size(); i++)
127  _diag_save_in[i]->sys().solution().add_vector(diag, _diag_save_in[i]->dofIndices());
128  }
129 }
130 
131 template <>
132 void
134 {
135 }
136 
137 template <>
138 void
140 {
141 }
142 
143 template <typename T, ComputeStage compute_stage>
144 void
146 {
147  std::vector<DualReal> residuals(_grad_test.size(), 0);
148 
149  precalculateResidual();
150 
151  if (_use_displaced_mesh)
152  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
153  {
154  const auto value = precomputeQpStrongResidual() * _ad_JxW[_qp] * _ad_coord[_qp];
155  for (_i = 0; _i < _grad_test.size(); _i++)
156  residuals[_i] += _grad_test[_i][_qp] * computeQpStabilization() * value;
157  }
158  else
159  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
160  {
161  const auto value = precomputeQpStrongResidual() * _JxW[_qp] * _coord[_qp];
162  for (_i = 0; _i < _grad_test.size(); _i++)
163  residuals[_i] += _regular_grad_test[_i][_qp] * computeQpStabilization() * value;
164  }
165 
166  auto & ce = _assembly.couplingEntries();
167  for (const auto & it : ce)
168  {
169  MooseVariableFEBase & ivariable = *(it.first);
170  MooseVariableFEBase & jvariable = *(it.second);
171 
172  unsigned int ivar = ivariable.number();
173  unsigned int jvar = jvariable.number();
174 
175  if (ivar != _var.number())
176  continue;
177 
178  size_t ad_offset = jvar * _sys.getMaxVarNDofsPerElem();
179 
180  prepareMatrixTag(_assembly, ivar, jvar);
181 
182  if (_local_ke.m() != _grad_test.size() || _local_ke.n() != jvariable.phiSize())
183  continue;
184 
185  precalculateResidual();
186  for (_i = 0; _i < _grad_test.size(); _i++)
187  for (_j = 0; _j < jvariable.phiSize(); _j++)
188  _local_ke(_i, _j) += residuals[_i].derivatives()[ad_offset + _j];
189 
190  accumulateTaggedLocalMatrix();
191  }
192 }
193 
194 template <>
195 void
197 {
198 }
199 
200 template <>
201 void
203 {
204 }
205 
206 template <typename T, ComputeStage compute_stage>
207 ADReal
209 {
210  mooseError("Override precomputeQpStrongResidual() in your ADKernelStabilized derived class!");
211 }
212 
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
MathUtils.h
SystemBase.h
defineADLegacyParams
defineADLegacyParams(ADKernelStabilized)
MooseVariableFEBase::phiSize
virtual size_t phiSize() const =0
Return phi size.
ADKernelStabilizedTempl::ADKernelStabilizedTempl
ADKernelStabilizedTempl(const InputParameters &parameters)
Definition: ADKernelStabilized.C:29
ADKernelTempl::validParams
static InputParameters validParams()
Definition: ADKernel.C:25
ADKernelStabilized.h
ADKernelStabilizedTempl::validParams
static InputParameters validParams()
Definition: ADKernelStabilized.C:23
Assembly.h
mooseError
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application.
Definition: MooseError.h:210
ADKernelStabilizedTempl::computeQpResidual
virtual ADReal computeQpResidual() override final
Compute this Kernel's contribution to the residual at the current quadrature point.
Definition: ADKernelStabilized.C:208
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
ADKernelStabilizedTempl::computeJacobian
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
Definition: ADKernelStabilized.C:83
ADKernelStabilizedTempl
Definition: ADKernelStabilized.h:19
ADKernelStabilizedTempl::computeResidual
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
Definition: ADKernelStabilized.C:37
ADKernelStabilizedTempl::computeADOffDiagJacobian
virtual void computeADOffDiagJacobian() override
Definition: ADKernelStabilized.C:145
ADKernelTempl
Definition: ADKernel.h:66
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48