https://mooseframework.inl.gov
ADIntegratedBC.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 "ADIntegratedBC.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "SubProblem.h"
15 #include "SystemBase.h"
16 #include "MooseVariableFE.h"
17 #include "MooseVariableScalar.h"
18 #include "ADUtils.h"
19 
20 #include "libmesh/quadrature.h"
21 
22 template <typename T>
25 {
28  return params;
29 }
30 
31 template <typename T>
33  : IntegratedBCBase(parameters),
34  MooseVariableInterface<T>(this,
35  false,
36  "variable",
40  ADFunctorInterface(this),
41  _var(*this->mooseVariable()),
42  _normals(_assembly.adNormals()),
43  _ad_q_points(_assembly.adQPointsFace()),
44  _test(_var.phiFace()),
45  _grad_test(_var.adGradPhiFace()),
46  _u(_var.adSln()),
47  _grad_u(_var.adGradSln()),
48  _ad_JxW(_assembly.adJxWFace()),
49  _ad_coord(_assembly.adCoordTransformation()),
50  _phi(_assembly.phi(_var)),
51  _use_displaced_mesh(getParam<bool>("use_displaced_mesh"))
52 {
54 
56 
57  _save_in.resize(_save_in_strings.size());
58  _diag_save_in.resize(_diag_save_in_strings.size());
59 
60  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
61  {
63 
64  if (var->feType() != _var.feType())
65  paramError(
66  "save_in",
67  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
69  _save_in[i] = var;
72  }
73 
74  _has_save_in = _save_in.size() > 0;
75 
76  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
77  {
79 
80  if (var->feType() != _var.feType())
81  paramError(
82  "diag_save_in",
83  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
85 
86  _diag_save_in[i] = var;
89  }
90 
91  _has_diag_save_in = _diag_save_in.size() > 0;
92 }
93 
94 template <typename T>
95 void
97 {
98  _residuals.resize(_test.size(), 0);
99  for (auto & r : _residuals)
100  r = 0;
101 
102  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
103  {
104  const auto jxw_c = _JxW[_qp] * _coord[_qp];
105  for (_i = 0; _i < _test.size(); _i++)
106  _residuals[_i] += jxw_c * raw_value(computeQpResidual());
107  }
108 
109  addResiduals(_assembly, _residuals, _var.dofIndices(), _var.scalingFactor());
110 
111  if (_has_save_in)
112  for (unsigned int i = 0; i < _save_in.size(); i++)
113  _save_in[i]->sys().solution().add_vector(_residuals.data(), _save_in[i]->dofIndices());
114 }
115 
116 template <typename T>
117 void
119 {
120  if (_residuals_and_jacobians.size() != _test.size())
121  _residuals_and_jacobians.resize(_test.size(), 0);
122  for (auto & r : _residuals_and_jacobians)
123  r = 0;
124 
125  if (_use_displaced_mesh)
126  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
127  {
128  _r = _ad_JxW[_qp];
129  _r *= _ad_coord[_qp];
130  for (_i = 0; _i < _test.size(); _i++)
131  _residuals_and_jacobians[_i] += _r * computeQpResidual();
132  }
133  else
134  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
135  {
136  const auto jxw_c = _JxW[_qp] * _coord[_qp];
137  for (_i = 0; _i < _test.size(); _i++)
138  _residuals_and_jacobians[_i] += jxw_c * computeQpResidual();
139  }
140 }
141 
142 template <typename T>
143 void
145 {
146  computeResidualsForJacobian();
147  addResidualsAndJacobian(
148  _assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor());
149 }
150 
151 template <typename T>
152 void
154 {
155  computeADJacobian();
156 
157  if (_has_diag_save_in && !_sys.computingScalingJacobian())
158  mooseError("_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use "
159  "the tagging system instead.");
160 }
161 
162 template <typename T>
163 void
165 {
166  computeResidualsForJacobian();
167  addJacobian(_assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor());
168 }
169 
170 template <typename T>
171 void
173 {
174  // Only need to do this once because AD does all the derivatives at once
175  if (jvar == _var.number())
176  computeADJacobian();
177 }
178 
179 template <typename T>
180 void
182 {
183 }
184 
185 template class ADIntegratedBCTempl<Real>;
VarFieldType
Definition: MooseTypes.h:721
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
std::vector< MooseVariableFEBase * > _diag_save_in
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
std::vector< MooseVariableFEBase * > _save_in
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
void computeResidual() override
Compute this object&#39;s contribution to the residual.
static InputParameters validParams()
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...
THREAD_ID _tid
The thread ID for this kernel.
static InputParameters validParams()
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
ADIntegratedBCTempl(const InputParameters &parameters)
MooseVariableFE< T > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
void computeADJacobian()
compute all the Jacobian entries
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
std::vector< AuxVariableName > _save_in_strings
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:714
Base class for deriving any boundary condition of a integrated type.
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.
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 ...
MooseVariableFE< T > & _var
The variable that this IntegratedBC operates on.
bool _has_save_in
The aux variables to save the residual contributions to.
void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
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
std::vector< AuxVariableName > _diag_save_in_strings
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:173
static InputParameters validParams()
Interface for objects that need to get values of MooseVariables.
Base class for deriving any boundary condition of a integrated type.
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:179
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
SystemBase & sys()
Get the system this variable is part of.
virtual void computeResidualsForJacobian()
compute the _residuals member for filling the Jacobian.
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.