https://mooseframework.inl.gov
NodalBC.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 "NodalBC.h"
11 
12 #include "Assembly.h"
13 #include "MooseVariableFE.h"
14 #include "SystemBase.h"
15 #include "NonlinearSystemBase.h"
16 #include "FEProblemBase.h"
17 
20 {
22 
23  return params;
24 }
25 
26 NodalBC::NodalBC(const InputParameters & parameters)
27  : NodalBCBase(parameters),
29  true,
30  "variable",
33  _var(*mooseVariable()),
34  _current_node(_var.node()),
35  _u(_var.dofValues())
36 {
38 
39  _save_in.resize(_save_in_strings.size());
40  _diag_save_in.resize(_diag_save_in_strings.size());
41 
42  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
43  {
45 
46  if (var->feType() != _var.feType())
47  paramError(
48  "save_in",
49  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
51 
52  _save_in[i] = var;
55  }
56 
57  _has_save_in = _save_in.size() > 0;
58 
59  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
60  {
62 
63  if (var->feType() != _var.feType())
64  paramError(
65  "diag_save_in",
66  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
68 
69  _diag_save_in[i] = var;
72  }
73 
74  _has_diag_save_in = _diag_save_in.size() > 0;
75 }
76 
77 void
79 {
80  if (_var.isNodalDefined())
81  {
82  const Real res = computeQpResidual();
83  setResidual(_sys, res, _var);
84 
85  if (_has_save_in)
86  for (unsigned int i = 0; i < _save_in.size(); i++)
87  _save_in[i]->sys().solution().set(_save_in[i]->nodalDofIndex(), res);
88  }
89 }
90 
91 void
93 {
94  // We call the user's computeQpJacobian() function and store the
95  // results in the _assembly object. We can't store them directly in
96  // the element stiffness matrix, as they will only be inserted after
97  // all the assembly is done.
98  if (_var.isNodalDefined())
99  {
100  const Real cached_val = computeQpJacobian();
101  const dof_id_type cached_row = _var.nodalDofIndex();
102 
103  // Cache the user's computeQpJacobian() value for later use.
105  cached_val,
106  cached_row,
107  cached_row,
108  /*scaling_factor=*/1);
109 
110  if (_has_diag_save_in)
111  for (unsigned int i = 0; i < _diag_save_in.size(); i++)
112  _diag_save_in[i]->sys().solution().set(_diag_save_in[i]->nodalDofIndex(), cached_val);
113  }
114 }
115 
116 void
117 NodalBC::computeOffDiagJacobian(const unsigned int jvar_num)
118 {
119  if (jvar_num == _var.number())
120  computeJacobian();
121  else
122  {
123  if (!_var.isNodalDefined())
124  return;
125 
126  const Real cached_val = computeQpOffDiagJacobian(jvar_num);
127 
128  if (cached_val == 0.)
129  // there's no reason to cache this if it's zero, and it can even lead to new nonzero
130  // allocations
131  return;
132 
133  const dof_id_type cached_row = _var.nodalDofIndex();
134  // Note: this only works for Lagrange variables...
135  const dof_id_type cached_col = _current_node->dof_number(_sys.number(), jvar_num, 0);
136 
137  // Cache the user's computeQpJacobian() value for later use.
139  cached_val,
140  cached_row,
141  cached_col,
142  /*scaling_factor=*/1);
143  }
144 }
145 
146 Real
148 {
149  return 1.;
150 }
151 
152 Real
153 NodalBC::computeQpOffDiagJacobian(unsigned int /*jvar*/)
154 {
155  return 0.;
156 }
157 
158 void
160 {
161  computeResidual();
162 
163  for (const auto & [ivariable, jvariable] : _fe_problem.couplingEntries(_tid, _sys.number()))
164  {
165  const unsigned int ivar = ivariable->number();
166  const unsigned int jvar = jvariable->number();
167 
168  if (ivar != _var.number())
169  continue;
170 
171  if (_is_implicit)
173  }
174 
176 }
VarFieldType
Definition: MooseTypes.h:722
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
Definition: NodalBC.C:117
const libMesh::FEType & feType() const
Get the type of finite element object.
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & couplingEntries(const THREAD_ID tid, const unsigned int nl_sys_num)
Class for stuff related to variables.
Definition: Adaptivity.h:31
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 ...
Definition: MooseBase.h:435
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: NodalBCBase.h:41
virtual void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: NodalBC.C:159
std::vector< AuxVariableName > _save_in_strings
Definition: NodalBCBase.h:43
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpResidual()=0
const Node *const & _current_node
current node being processed
Definition: NodalBC.h:41
NodalBC(const InputParameters &parameters)
Definition: NodalBC.C:26
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
THREAD_ID _tid
The thread ID for this kernel.
virtual Assembly & assembly(const THREAD_ID tid, const unsigned int sys_num) override
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: NodalBC.C:78
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< Real > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
MooseVariable & _var
Definition: NodalBC.h:38
virtual Real computeQpJacobian()
The user can override this function to compute the "on-diagonal" Jacobian contribution for this Nodal...
Definition: NodalBC.C:147
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:715
FEProblemBase & _fe_problem
Reference to this kernel&#39;s FEProblemBase.
std::vector< MooseVariableFEBase * > _save_in
Definition: NodalBCBase.h:42
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.
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1149
Base class for deriving any boundary condition that works at nodes.
Definition: NodalBCBase.h:26
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
Definition: NodalBC.C:153
std::vector< AuxVariableName > _diag_save_in_strings
Definition: NodalBCBase.h:48
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
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: NodalBCBase.h:46
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
const dof_id_type & nodalDofIndex() const override
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:174
Interface for objects that need to get values of MooseVariables.
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:180
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: NodalBCBase.h:47
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
static InputParameters validParams()
Definition: NodalBC.C:19
void addJacobianElement(Assembly &assembly, Real value, dof_id_type row_index, dof_id_type column_index, Real scaling_factor)
Add into a single Jacobian element.
void setResidual(SystemBase &sys, const T &residual, MooseVariableFE< T > &var)
Set residual using the variables&#39; insertion API.
SystemBase & sys()
Get the system this variable is part of.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: NodalBC.C:92
bool _is_implicit
If the object is using implicit or explicit form.
uint8_t dof_id_type
static InputParameters validParams()
Definition: NodalBCBase.C:13