www.mooseframework.org
ADIntegratedBC.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 "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  if (_use_displaced_mesh)
103  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
104  for (_i = 0; _i < _test.size(); _i++)
105  _residuals[_i] += raw_value(_ad_JxW[_qp] * _ad_coord[_qp] * computeQpResidual());
106  else
107  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
108  for (_i = 0; _i < _test.size(); _i++)
109  _residuals[_i] += raw_value(_JxW[_qp] * _coord[_qp] * computeQpResidual());
110 
111  addResiduals(_assembly, _residuals, _var.dofIndices(), _var.scalingFactor());
112 
113  if (_has_save_in)
114  for (unsigned int i = 0; i < _save_in.size(); i++)
115  _save_in[i]->sys().solution().add_vector(_residuals.data(), _save_in[i]->dofIndices());
116 }
117 
118 template <typename T>
119 void
121 {
122  if (_residuals_and_jacobians.size() != _test.size())
123  _residuals_and_jacobians.resize(_test.size(), 0);
124  for (auto & r : _residuals_and_jacobians)
125  r = 0;
126 
127  if (_use_displaced_mesh)
128  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
129  {
130  _r = _ad_JxW[_qp];
131  _r *= _ad_coord[_qp];
132  for (_i = 0; _i < _test.size(); _i++)
133  _residuals_and_jacobians[_i] += _r * computeQpResidual();
134  }
135  else
136  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
137  for (_i = 0; _i < _test.size(); _i++)
138  _residuals_and_jacobians[_i] += _JxW[_qp] * _coord[_qp] * computeQpResidual();
139 }
140 
141 template <typename T>
142 void
144 {
145  computeResidualsForJacobian();
146  addResidualsAndJacobian(
147  _assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor());
148 }
149 
150 template <typename T>
151 void
153 {
154  computeADJacobian();
155 
156  if (_has_diag_save_in && !_sys.computingScalingJacobian())
157  mooseError("_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use "
158  "the tagging system instead.");
159 }
160 
161 template <typename T>
162 void
164 {
165  computeResidualsForJacobian();
166  addJacobian(_assembly, _residuals_and_jacobians, _var.dofIndices(), _var.scalingFactor());
167 }
168 
169 template <typename T>
170 void
172 {
173  // Only need to do this once because AD does all the derivatives at once
174  if (jvar == _var.number())
175  computeADJacobian();
176 }
177 
178 template <typename T>
179 void
181 {
182 }
183 
184 template class ADIntegratedBCTempl<Real>;
VarFieldType
Definition: MooseTypes.h:634
virtual void haveADObjects(bool have_ad_objects)
Method for setting whether we have any ad objects.
Definition: SubProblem.h:737
std::vector< MooseVariableFEBase * > _diag_save_in
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:24
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:299
void computeResidual() override
Compute this object&#39;s contribution to the residual.
static InputParameters validParams()
auto raw_value(const Eigen::Map< T > &in)
Definition: ADReal.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()
const FEType & feType() const
Get the type of finite element object.
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:627
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:170
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:176
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.