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