www.mooseframework.org
ADNodalBC.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 "ADNodalBC.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 
19 #include "libmesh/quadrature.h"
20 
23 
24 template <typename T, ComputeStage compute_stage>
27 {
29  return params;
30 }
31 
32 template <typename T, ComputeStage compute_stage>
34  : NodalBCBase(parameters),
35  MooseVariableInterface<T>(this,
36  true,
37  "variable",
39  std::is_same<T, Real>::value ? Moose::VarFieldType::VAR_FIELD_STANDARD
41  _var(*this->mooseVariable()),
42  _current_node(_var.node()),
43  _u(_var.template adNodalValue<compute_stage>())
44 {
46 }
47 
48 template <typename T>
49 T &
50 conversionHelper(T & value, const unsigned int &)
51 {
52  return value;
53 }
54 
55 template <typename T>
56 T &
57 conversionHelper(libMesh::VectorValue<T> & value, const unsigned int & i)
58 {
59  return value(i);
60 }
61 
62 template <typename T, ComputeStage compute_stage>
63 void
65 {
66  const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
67 
68  auto residual = computeQpResidual();
69 
70  for (auto tag_id : _vector_tags)
71  if (_sys.hasVector(tag_id))
72  for (size_t i = 0; i < dof_indices.size(); ++i)
73  _sys.getVector(tag_id).set(dof_indices[i], conversionHelper(residual, i));
74 }
75 
76 template <>
77 void
79 {
80 }
81 
82 template <>
83 void
85 {
86 }
87 
88 template <typename T, ComputeStage compute_stage>
89 void
91 {
92  auto ad_offset = _var.number() * _sys.getMaxVarNDofsPerNode();
93  auto residual = computeQpResidual();
94  const std::vector<dof_id_type> & cached_rows = _var.dofIndices();
95 
96  // Cache the user's computeQpJacobian() value for later use.
97  for (auto tag : _matrix_tags)
98  if (_sys.hasMatrix(tag))
99  for (size_t i = 0; i < cached_rows.size(); ++i)
100  _fe_problem.assembly(0).cacheJacobianContribution(
101  cached_rows[i],
102  cached_rows[i],
103  conversionHelper(residual, i).derivatives()[ad_offset + i],
104  tag);
105 }
106 
107 template <>
108 void
110 {
111 }
112 template <>
113 void
115 {
116 }
117 
118 template <typename T, ComputeStage compute_stage>
119 void
121 {
122  if (jvar == _var.number())
123  computeJacobian();
124  else
125  {
126  auto ad_offset = jvar * _sys.getMaxVarNDofsPerNode();
127  auto residual = computeQpResidual();
128  const std::vector<dof_id_type> & cached_rows = _var.dofIndices();
129 
130  // Note: this only works for Lagrange variables...
131  dof_id_type cached_col = _current_node->dof_number(_sys.number(), jvar, 0);
132 
133  // Cache the user's computeQpJacobian() value for later use.
134  for (auto tag : _matrix_tags)
135  if (_sys.hasMatrix(tag))
136  for (size_t i = 0; i < cached_rows.size(); ++i)
137  _fe_problem.assembly(0).cacheJacobianContribution(
138  cached_rows[i],
139  cached_col,
140  conversionHelper(residual, i).derivatives()[ad_offset + i],
141  tag);
142  }
143 }
144 
145 template <>
146 void
148 {
149 }
150 template <>
151 void
153 {
154 }
155 
156 template <typename T, ComputeStage compute_stage>
157 void
159 {
160  auto ad_offset = jvar * _sys.getMaxVarNDofsPerNode();
161  auto residual = computeQpResidual();
162  const std::vector<dof_id_type> & cached_rows = _var.dofIndices();
163 
164  std::vector<dof_id_type> scalar_dof_indices;
165 
166  _sys.dofMap().SCALAR_dof_indices(scalar_dof_indices, jvar);
167 
168  // Our residuals rely on returning a single scalar and we don't provide any arguments to
169  // computeQpResidual so I think it only makes sense to assume that our SCALAR variable should be
170  // order one
171  mooseAssert(scalar_dof_indices.size() == 1,
172  "ADNodalBC only allows coupling of first order SCALAR variables");
173 
174  // Cache the user's computeQpJacobian() value for later use.
175  for (auto tag : _matrix_tags)
176  if (_sys.hasMatrix(tag))
177  for (size_t i = 0; i < cached_rows.size(); ++i)
178  _fe_problem.assembly(0).cacheJacobianContribution(
179  cached_rows[i],
180  scalar_dof_indices[0],
181  conversionHelper(residual, i).derivatives()[ad_offset + i],
182  tag);
183 }
184 
185 template <>
186 void
188 {
189 }
190 
191 template <>
192 void
194 {
195 }
196 
197 template class ADNodalBCTempl<Real, RESIDUAL>;
198 template class ADNodalBCTempl<Real, JACOBIAN>;
Moose::VarFieldType
VarFieldType
Definition: MooseTypes.h:613
Moose
Definition: Moose.h:116
ADNodalBCTempl::ADNodalBCTempl
ADNodalBCTempl(const InputParameters &parameters)
Definition: ADNodalBC.C:33
MooseVariableInterface::mooseVariable
MooseVariableFE< T > * mooseVariable() const
Get the variable that this object is using.
Definition: MooseVariableInterface.C:70
SystemBase.h
ADNodalBCTempl::validParams
static InputParameters validParams()
Definition: ADNodalBC.C:26
MooseVariableDependencyInterface::addMooseVariableDependency
void addMooseVariableDependency(MooseVariableFEBase *var)
Call this function to add the passed in MooseVariableFEBase as a variable that this object depends on...
Definition: MooseVariableDependencyInterface.h:36
ADNodalBCTempl::computeOffDiagJacobianScalar
void computeOffDiagJacobianScalar(unsigned int jvar) override
Compute the off-diagonal contributions from scalar variables.
Definition: ADNodalBC.C:158
Assembly.h
defineADBaseValidParams
defineADBaseValidParams(ADNodalBC, NodalBCBase,)
libMesh::VectorValue
Definition: Assembly.h:30
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
MooseVariableScalar.h
NodalBCBase
Base class for deriving any boundary condition that works at nodes.
Definition: NodalBCBase.h:32
ADNodalBCTempl
Base class for deriving any automatic differentiation boundary condition of a integrated type.
Definition: ADNodalBC.h:19
Moose::VAR_FIELD_VECTOR
Definition: MooseTypes.h:617
MooseVariableFE.h
SubProblem.h
Moose::VAR_NONLINEAR
Definition: MooseTypes.h:608
Moose::VarKindType
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:606
std
Definition: TheWarehouse.h:80
ADNodalBCTempl::computeJacobian
void computeJacobian() override
Definition: ADNodalBC.C:90
NodalBCBase::validParams
static InputParameters validParams()
Definition: NodalBCBase.C:20
ADNodalBCTempl::computeOffDiagJacobian
void computeOffDiagJacobian(unsigned int jvar) override
Definition: ADNodalBC.C:120
MooseVariableInterface
Interface for objects that need to get values of MooseVariables.
Definition: MooseVariableInterface.h:24
conversionHelper
T & conversionHelper(T &value, const unsigned int &)
Definition: ADNodalBC.C:50
ADNodalBCTempl::computeResidual
void computeResidual() override
Definition: ADNodalBC.C:64
Moose::VAR_FIELD_STANDARD
Definition: MooseTypes.h:615
ADNodalBC.h