https://mooseframework.inl.gov
ADKernelStabilized.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 "ADKernelStabilized.h"
11 #include "MathUtils.h"
12 #include "Assembly.h"
13 #include "SystemBase.h"
14 #include "ADUtils.h"
15 
16 // libmesh includes
17 #include "libmesh/threads.h"
18 
19 template <typename T>
22 {
24 }
25 
26 template <typename T>
28  : ADKernelTempl<T>(parameters)
29 {
30 }
31 
32 template <typename T>
33 void
35 {
36  prepareVectorTag(_assembly, _var.number());
37 
38  precalculateResidual();
39  const unsigned int n_test = _grad_test.size();
40 
41  if (_use_displaced_mesh)
42  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
43  {
44  const auto value = precomputeQpStrongResidual() * _ad_JxW[_qp] * _ad_coord[_qp];
45  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
46  _local_re(_i) += raw_value(_grad_test[_i][_qp] * computeQpStabilization() * value);
47  }
48  else
49  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
50  {
51  const auto value = precomputeQpStrongResidual() * _JxW[_qp] * _coord[_qp];
52  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
53  _local_re(_i) += raw_value(_regular_grad_test[_i][_qp] * computeQpStabilization() * value);
54  }
55 
56  accumulateTaggedLocalResidual();
57 
58  if (_has_save_in)
59  {
60  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
61  for (unsigned int i = 0; i < _save_in.size(); i++)
62  _save_in[i]->sys().solution().add_vector(_local_re, _save_in[i]->dofIndices());
63  }
64 }
65 
66 template <typename T>
67 void
69 {
70  if (_residuals.size() != _grad_test.size())
71  _residuals.resize(_grad_test.size(), 0);
72  for (auto & r : _residuals)
73  r = 0;
74 
75  precalculateResidual();
76 
77  if (_use_displaced_mesh)
78  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
79  {
80  const auto value = precomputeQpStrongResidual() * _ad_JxW[_qp] * _ad_coord[_qp];
81  for (_i = 0; _i < _grad_test.size(); _i++)
82  _residuals[_i] += _grad_test[_i][_qp] * computeQpStabilization() * value;
83  }
84  else
85  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
86  {
87  const auto value = precomputeQpStrongResidual() * _JxW[_qp] * _coord[_qp];
88  for (_i = 0; _i < _grad_test.size(); _i++)
89  _residuals[_i] += _regular_grad_test[_i][_qp] * computeQpStabilization() * value;
90  }
91 }
92 
93 template <typename T>
94 ADReal
96 {
97  mooseError("Override precomputeQpStrongResidual() in your ADKernelStabilized derived class!");
98 }
99 
100 template class ADKernelStabilizedTempl<Real>;
virtual ADReal computeQpResidual() override final
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
void computeResidualsForJacobian() override
compute the _residuals member for filling the Jacobian.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
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...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
static InputParameters validParams()
Definition: ADKernel.C:24
ADKernelStabilizedTempl(const InputParameters &parameters)
static InputParameters validParams()
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.