Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
ArrayNodalBC.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 "ArrayNodalBC.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  return params;
23 }
24 
26  : NodalBCBase(parameters),
28  true,
29  "variable",
32  _var(*mooseVariable()),
33  _current_node(_var.node()),
34  _u(_var.nodalValue()),
35  _count(_var.count()),
36  _work_vector(_count)
37 {
39 }
40 
41 void
43 {
44  if (_var.isNodalDefined())
45  {
46  _work_vector.setZero();
48  mooseAssert(_work_vector.size() == _count,
49  "Size of local residual is not equal to the number of array variable components");
50 
52  }
53 }
54 
55 void
57 {
58  if (_var.isNodalDefined())
59  {
60  const RealEigenVector cached_val = computeQpJacobian();
61  const dof_id_type cached_row = _var.nodalDofIndex();
62 
63  for (const auto i : make_range(_var.count()))
65  cached_val(i),
66  cached_row + i,
67  cached_row + i,
68  /*scaling_factor=*/1);
69  }
70 }
71 
72 void
73 ArrayNodalBC::computeOffDiagJacobian(const unsigned int jvar_num)
74 {
75  if (!_var.isNodalDefined())
76  return;
77 
78  const auto & jvar = getVariable(jvar_num);
79 
80  const RealEigenMatrix cached_val =
81  computeQpOffDiagJacobian(const_cast<MooseVariableFieldBase &>(jvar));
82  const dof_id_type cached_row = _var.nodalDofIndex();
83  // Note: this only works for Lagrange variables...
84  const dof_id_type cached_col = _current_node->dof_number(_sys.number(), jvar_num, 0);
85 
86  // Cache the user's computeQpJacobian() value for later use.
87  for (const auto i : make_range(_var.count()))
88  for (const auto j : make_range(jvar.count()))
90  cached_val(i, j),
91  cached_row + i,
92  cached_col + j,
93  /*scaling_factor=*/1);
94 }
95 
98 {
99  return RealEigenVector::Ones(_var.count());
100  ;
101 }
102 
105 {
106  if (jvar.number() == _var.number())
107  return RealEigenMatrix::Identity(_var.count(), jvar.count());
108  else
109  return RealEigenMatrix::Zero(_var.count(), jvar.count());
110 }
VarFieldType
Definition: MooseTypes.h:715
unsigned int number() const
Get variable number coming from libMesh.
const unsigned int _count
Number of components of the array variable.
Definition: ArrayNodalBC.h:63
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: ArrayNodalBC.C:56
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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.
virtual Assembly & assembly(const THREAD_ID tid, const unsigned int sys_num) override
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< RealEigenVector > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:708
FEProblemBase & _fe_problem
Reference to this kernel&#39;s FEProblemBase.
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1130
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:145
Base class for deriving any boundary condition that works at nodes.
Definition: NodalBCBase.h:26
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
Definition: ArrayNodalBC.C:73
ArrayMooseVariable & _var
Definition: ArrayNodalBC.h:35
ArrayNodalBC(const InputParameters &parameters)
Definition: ArrayNodalBC.C:25
const Node *const & _current_node
current node being processed
Definition: ArrayNodalBC.h:38
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
virtual RealEigenMatrix computeQpOffDiagJacobian(MooseVariableFEBase &jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
Definition: ArrayNodalBC.C:104
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
const dof_id_type & nodalDofIndex() const override
RealEigenVector _work_vector
Work vector for residual.
Definition: ArrayNodalBC.h:67
Interface for objects that need to get values of MooseVariables.
IntRange< T > make_range(T beg, T end)
virtual void computeQpResidual(RealEigenVector &residual)=0
Compute this BC&#39;s contribution to the residual at the current quadrature point, to be filled in resid...
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:142
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
static InputParameters validParams()
Definition: ArrayNodalBC.C:19
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ArrayNodalBC.C:42
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.
uint8_t dof_id_type
virtual RealEigenVector computeQpJacobian()
The user can override this function to compute the "on-diagonal" Jacobian contribution for this Vecto...
Definition: ArrayNodalBC.C:97
static InputParameters validParams()
Definition: NodalBCBase.C:13