https://mooseframework.inl.gov
NodalKernel.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 "NodalKernel.h"
11 #include "Problem.h"
12 #include "SubProblem.h"
13 #include "SystemBase.h"
14 #include "MooseVariableFE.h"
15 #include "Assembly.h"
16 
19 {
20  auto params = NodalKernelBase::validParams();
21  params.addParam<std::vector<AuxVariableName>>(
22  "save_in",
23  {},
24  "The name of auxiliary variables to save this BC's residual contributions to. "
25  "Everything about that variable must match everything about this variable (the "
26  "type, what blocks it's on, etc.)");
27  params.addParam<std::vector<AuxVariableName>>(
28  "diag_save_in",
29  {},
30  "The name of auxiliary variables to save this BC's diagonal jacobian "
31  "contributions to. Everything about that variable must match everything "
32  "about this variable (the type, what blocks it's on, etc.)");
33  params.addParamNamesToGroup("diag_save_in save_in", "Advanced");
34  return params;
35 }
36 
38  : NodalKernelBase(parameters),
40  true,
41  "variable",
44  _var(*mooseVariable()),
45  _u(_var.dofValues()),
46  _save_in_strings(parameters.get<std::vector<AuxVariableName>>("save_in")),
47  _diag_save_in_strings(parameters.get<std::vector<AuxVariableName>>("diag_save_in"))
48 
49 {
51 
52  _save_in.resize(_save_in_strings.size());
53  _diag_save_in.resize(_diag_save_in_strings.size());
54 
55  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
56  {
58 
59  if (var->feType() != _var.feType())
60  paramError(
61  "save_in",
62  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
64 
65  _save_in[i] = var;
68  }
69 
70  _has_save_in = _save_in.size() > 0;
71 
72  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
73  {
75 
76  if (var->feType() != _var.feType())
77  paramError(
78  "diag_save_in",
79  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
81 
82  _diag_save_in[i] = var;
85  }
86 
87  _has_diag_save_in = _diag_save_in.size() > 0;
88 }
89 
90 void
92 {
93  if (_var.isNodalDefined())
94  {
95  const dof_id_type & dof_idx = _var.nodalDofIndex();
96  _qp = 0;
97  const Real res = computeQpResidual();
99  std::array<Real, 1>{{res}},
100  std::array<dof_id_type, 1>{{dof_idx}},
101  _var.scalingFactor());
102 
103  if (_has_save_in)
104  for (const auto & var : _save_in)
105  var->sys().solution().add(var->nodalDofIndex(), res);
106  }
107 }
108 
109 void
111 {
112  if (_var.isNodalDefined())
113  {
114  _qp = 0;
115  const Real cached_val = computeQpJacobian();
116  const dof_id_type cached_row = _var.nodalDofIndex();
117 
118  addJacobianElement(_assembly, cached_val, cached_row, cached_row, _var.scalingFactor());
119 
120  if (_has_diag_save_in)
121  {
122  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
123  for (const auto & var : _diag_save_in)
124  var->sys().solution().add(var->nodalDofIndex(), cached_val);
125  }
126  }
127 }
128 
129 void
130 NodalKernel::computeOffDiagJacobian(const unsigned int jvar_num)
131 {
132  if (_var.isNodalDefined())
133  {
134  if (jvar_num == _var.number())
135  computeJacobian();
136  else
137  {
138  _qp = 0;
139  const Real cached_val = computeQpOffDiagJacobian(jvar_num);
140  const dof_id_type cached_row = _var.nodalDofIndex();
141 
142  // Note: this only works for equal order Lagrange variables...
143  const dof_id_type cached_col = _current_node->dof_number(_sys.number(), jvar_num, 0);
144 
145  addJacobianElement(_assembly, cached_val, cached_row, cached_col, _var.scalingFactor());
146  }
147  }
148 }
149 
150 Real
152 {
153  return 0.;
154 }
155 
156 Real
158 {
159  return 0.;
160 }
VarFieldType
Definition: MooseTypes.h:723
std::vector< AuxVariableName > _save_in_strings
Definition: NodalKernel.h:87
const libMesh::FEType & feType() const
Get the type of finite element object.
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:439
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
Definition: NodalKernel.C:157
void addResiduals(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals corresponding to the provided dof indices.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeJacobian() override
Compute the Jacobian at one node.
Definition: NodalKernel.C:110
THREAD_ID _tid
The thread ID for this kernel.
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< Real > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
MooseVariable & _var
variable this works on
Definition: NodalKernel.h:79
std::vector< AuxVariableName > _diag_save_in_strings
Definition: NodalKernel.h:92
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:716
virtual void computeResidual() override
Compute the residual at the current node.
Definition: NodalKernel.C:91
virtual Real computeQpJacobian()
The user can override this function to compute the "on-diagonal" Jacobian contribution.
Definition: NodalKernel.C:151
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:1157
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: NodalKernel.h:90
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: NodalKernel.h:85
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
static InputParameters validParams()
Class constructor.
Base class for creating new types of nodal kernels.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 * > _save_in
Definition: NodalKernel.h:86
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const Node *const & _current_node
current node being processed
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: NodalKernel.h:91
static InputParameters validParams()
Class constructor.
Definition: NodalKernel.C:18
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.
SystemBase & sys()
Get the system this variable is part of.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Compute the off-diagonal Jacobian at one node.
Definition: NodalKernel.C:130
const Elem & get(const ElemType type_in)
unsigned int _qp
Quadrature point index.
virtual Real computeQpResidual()=0
The user can override this function to compute the residual at a node.
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
NodalKernel(const InputParameters &parameters)
Definition: NodalKernel.C:37
uint8_t dof_id_type