https://mooseframework.inl.gov
VectorNodalBC.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 "VectorNodalBC.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 {
36  if (_var.feType().family != LAGRANGE_VEC)
37  mooseError("Vector nodal boundary conditions only make sense for LAGRANGE_VEC variables");
39 }
40 
41 void
43 {
44  const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
45  if (dof_indices.empty())
46  return;
47 
48  const RealVectorValue res = computeQpResidual();
49 
50  for (const auto i : index_range(dof_indices))
51  setResidual(_sys, res(i), dof_indices[i]);
52 }
53 
54 void
56 {
57  const std::vector<dof_id_type> & cached_rows = _var.dofIndices();
58  if (cached_rows.empty())
59  return;
60 
61  const RealVectorValue cached_val = computeQpJacobian();
62 
63  // Cache the user's computeQpJacobian() value for later use.
64  for (const auto i : index_range(cached_rows))
66  cached_val(i),
67  cached_rows[i],
68  cached_rows[i],
69  /*scaling_factor=*/1);
70 }
71 
72 void
73 VectorNodalBC::computeOffDiagJacobian(const unsigned int jvar_num)
74 {
75  if (jvar_num == _var.number())
77  else
78  {
79  const std::vector<dof_id_type> & cached_rows = _var.dofIndices();
80  if (cached_rows.empty())
81  return;
82 
83  const Real cached_val = computeQpOffDiagJacobian(jvar_num);
84  // Note: this only works for Lagrange variables...
85  const dof_id_type cached_col = _current_node->dof_number(_sys.number(), jvar_num, 0);
86 
87  // Cache the user's computeQpJacobian() value for later use.
88  for (const auto i : index_range(cached_rows))
90  cached_val,
91  cached_rows[i],
92  cached_col,
93  /*scaling_factor=*/1);
94  }
95 }
96 
99 {
100  return RealVectorValue(1., 1., 1.);
101 }
102 
103 Real
105 {
106  return 0.;
107 }
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: VectorNodalBC.C:73
const libMesh::FEType & feType() const
Get the type of finite element object.
LAGRANGE_VEC
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: VectorNodalBC.C:55
unsigned int number() const
Get variable number coming from libMesh.
virtual RealVectorValue computeQpJacobian()
The user can override this function to compute the "on-diagonal" Jacobian contribution for this Vecto...
Definition: VectorNodalBC.C:98
static InputParameters validParams()
Definition: VectorNodalBC.C:19
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual RealVectorValue computeQpResidual()=0
virtual Assembly & assembly(const THREAD_ID tid, const unsigned int sys_num) override
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< RealVectorValue > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
VectorNodalBC(const InputParameters &parameters)
Definition: VectorNodalBC.C:25
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:715
FEProblemBase & _fe_problem
Reference to this kernel&#39;s FEProblemBase.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: VectorNodalBC.C:42
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
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
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
Interface for objects that need to get values of MooseVariables.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:267
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
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.
auto index_range(const T &sizable)
const VectorMooseVariable & _var
Definition: VectorNodalBC.h:35
uint8_t dof_id_type
const Node *const & _current_node
current node being processed
Definition: VectorNodalBC.h:38
static InputParameters validParams()
Definition: NodalBCBase.C:13