https://mooseframework.inl.gov
ArrayNodalKernel.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 "ArrayNodalKernel.h"
11 #include "SubProblem.h"
12 #include "SystemBase.h"
13 #include "MooseVariableFE.h"
14 #include "Assembly.h"
15 
18 {
20 }
21 
23  : NodalKernelBase(parameters),
25  true,
26  "variable",
29  _var(*mooseVariable()),
30  _u(_var.dofValues()),
31  _count(_var.count()),
32  _work_vector(_count),
33  _work_dofs(_count)
34 {
36 }
37 
38 void
40 {
41  const dof_id_type root_dof_idx = _var.nodalDofIndex();
42  std::iota(_work_dofs.begin(), _work_dofs.end(), root_dof_idx);
43 }
44 
45 void
47 {
48  if (!_var.isNodalDefined())
49  return;
50  prepareDofs();
51  _qp = 0;
54 }
55 
56 void
58 {
59  if (!_var.isNodalDefined())
60  return;
61  prepareDofs();
62  _qp = 0;
63  const auto jacobian = computeQpJacobian();
64  for (const auto i : make_range(_count))
66  _assembly, jacobian(i), _work_dofs[i], _work_dofs[i], _var.arrayScalingFactor()[i]);
67 }
68 
69 void
70 ArrayNodalKernel::computeOffDiagJacobian(const unsigned int jvar_num)
71 {
72  if (!_var.isNodalDefined())
73  return;
74 
75  if (jvar_num == _var.number())
77  else
78  {
79  prepareDofs();
80 
81  const auto & jvar = getVariable(jvar_num);
82  const auto jacobian = computeQpOffDiagJacobian(jvar);
83  const auto root_jvar_idx = _current_node->dof_number(_sys.number(), jvar_num, 0);
84 
85  for (const auto i : make_range(_var.count()))
86  for (const auto j : make_range(jvar.count()))
88  jacobian(i, j),
89  _work_dofs[i],
90  root_jvar_idx + j,
92  }
93 }
94 
97 {
98  return RealEigenVector::Zero(_count);
99 }
100 
103 {
104  return RealEigenMatrix::Zero(_count, jvar.count());
105 }
VarFieldType
Definition: MooseTypes.h:723
RealEigenVector _work_vector
Work vector for residual.
virtual void computeResidual() override
Compute and assemble the residual at the current node.
const std::vector< Real > & arrayScalingFactor() const
unsigned int number() const
Get variable number coming from libMesh.
virtual void computeJacobian() override
Compute and assemble the Jacobian at one node.
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.
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
void prepareDofs()
Prepare our array variable degrees of freedom.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFieldBase &jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
This class provides an interface for common operations on field variables of both FE and FV types wit...
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< RealEigenVector > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
MooseVariableFE< RealEigenVector > & _var
variable this works on
static InputParameters validParams()
Class constructor.
ArrayNodalKernel(const InputParameters &parameters)
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:716
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1157
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:150
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
virtual void computeQpResidual(RealEigenVector &residual)=0
The user must override this function to compute the residual at a node.
std::vector< dof_id_type > _work_dofs
Work vector for the degree of freedom indices.
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.
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
const dof_id_type & nodalDofIndex() const override
const unsigned int _count
Number of components of the array variable.
Interface for objects that need to get values of MooseVariables.
IntRange< T > make_range(T beg, T end)
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
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
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.
virtual RealEigenVector computeQpJacobian()
The user can override this function to compute the "on-diagonal" Jacobian contribution.
unsigned int _qp
Quadrature point index.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Compute the off-diagonal Jacobian at one node.
uint8_t dof_id_type