Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
NodeElemConstraint.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 "NodeElemConstraint.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseMesh.h"
15 
18 {
20  return params;
21 }
22 
24  : NodeElemConstraintBase(parameters),
25  _u_primary(_primary_var.slnNeighbor()),
26  _u_secondary(_var.dofValues()),
27  _grad_phi_primary(_assembly.gradPhiNeighbor(_primary_var)),
28  _grad_test_primary(_var.gradPhiNeighbor()),
29  _grad_u_primary(_primary_var.gradSlnNeighbor())
30 {
31 }
32 
33 void
35 {
36  _qp = 0;
37 
39  for (_i = 0; _i < _test_primary.size(); _i++)
42 
44  for (_i = 0; _i < _test_secondary.size(); _i++)
47 }
48 
49 void
51 {
53 
56 
57  for (_i = 0; _i < _test_secondary.size(); _i++)
58  // Loop over the connected dof indices so we can get all the jacobian contributions
59  for (_j = 0; _j < _connected_dof_indices.size(); _j++)
61 
63  if (_local_ke.m() && _local_ke.n())
64  for (_i = 0; _i < _test_secondary.size(); _i++)
65  for (_j = 0; _j < _phi_primary.size(); _j++)
68 
69  for (_i = 0; _i < _test_primary.size(); _i++)
70  // Loop over the connected dof indices so we can get all the jacobian contributions
71  for (_j = 0; _j < _connected_dof_indices.size(); _j++)
73 
76  if (_local_ke.m() && _local_ke.n())
77  for (_i = 0; _i < _test_primary.size(); _i++)
78  for (_j = 0; _j < _phi_primary.size(); _j++)
81 }
82 
83 void
84 NodeElemConstraint::computeOffDiagJacobian(const unsigned int jvar_num)
85 {
86  getConnectedDofIndices(jvar_num);
87 
90 
91  for (_i = 0; _i < _test_secondary.size(); _i++)
92  // Loop over the connected dof indices so we can get all the jacobian contributions
93  for (_j = 0; _j < _connected_dof_indices.size(); _j++)
95 
97  for (_i = 0; _i < _test_secondary.size(); _i++)
98  for (_j = 0; _j < _phi_primary.size(); _j++)
101 
102  if (_Kne.m() && _Kne.n())
103  for (_i = 0; _i < _test_primary.size(); _i++)
104  // Loop over the connected dof indices so we can get all the jacobian contributions
105  for (_j = 0; _j < _connected_dof_indices.size(); _j++)
107 
109  for (_i = 0; _i < _test_primary.size(); _i++)
110  for (_j = 0; _j < _phi_primary.size(); _j++)
113 }
114 
115 void
117 {
118  MooseVariableFEBase & var = _sys.getVariable(0, var_num);
119 
120  _connected_dof_indices.clear();
121  std::set<dof_id_type> unique_dof_indices;
122 
123  auto node_to_elem_pair = _node_to_elem_map.find(_current_node->id());
124  mooseAssert(node_to_elem_pair != _node_to_elem_map.end(), "Missing entry in node to elem map");
125  const std::vector<dof_id_type> & elems = node_to_elem_pair->second;
126 
127  // Get the dof indices from each elem connected to the node
128  for (const auto & cur_elem : elems)
129  {
130  std::vector<dof_id_type> dof_indices;
131 
132  var.getDofIndices(_mesh.elemPtr(cur_elem), dof_indices);
133 
134  for (const auto & dof : dof_indices)
135  unique_dof_indices.insert(dof);
136  }
137 
138  for (const auto & dof : unique_dof_indices)
139  _connected_dof_indices.push_back(dof);
140 
142 
143  const dof_id_type current_node_var_dof_index = _sys.getVariable(0, var_num).nodalDofIndex();
144 
145  // Fill up _phi_secondary so that it is 1 when j corresponds to the dof associated with this node
146  // and 0 for every other dof
147  // This corresponds to evaluating all of the connected shape functions at _this_ node
148  _qp = 0;
149  for (unsigned int j = 0; j < _connected_dof_indices.size(); j++)
150  {
151  _phi_secondary[j].resize(1);
152 
153  if (_connected_dof_indices[j] == current_node_var_dof_index)
154  _phi_secondary[j][_qp] = 1.0;
155  else
156  _phi_secondary[j][_qp] = 0.0;
157  }
158 }
159 
160 Real
162 {
163  mooseError("Derived classes must implement computeQpJacobian.");
164  return 0;
165 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
virtual Real computeQpOffDiagJacobian(Moose::ConstraintJacobianType, unsigned int)
This is the virtual method that derived classes should override for computing the off-diag Jacobian...
static InputParameters validParams()
NodeElemConstraint(const InputParameters &parameters)
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual Elem * elemPtr(const dof_id_type i)
Definition: MooseMesh.C:3082
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpResidual(Moose::ConstraintType type)=0
This is the virtual method that derived classes should override for computing the residual...
virtual void computeJacobian() override
Computes the jacobian for the current element.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void getDofIndices(const Elem *, std::vector< dof_id_type > &) const
unsigned int m() const
This class provides an interface for common operations on field variables of both FE and FV types wit...
unsigned int _i
Definition: Constraint.h:35
virtual const dof_id_type & nodalDofIndex() const =0
VariableTestValue _test_secondary
Shape function on the secondary side. This will always only have one entry and that entry will always...
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariable & _var
secondary node variable
void prepareMatrixTagNeighbor(Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type)
Prepare data for computing element jacobian according to the active tags for DG and interface kernels...
void getConnectedDofIndices(unsigned int var_num)
Gets the indices for all dofs connected to the constraint.
unsigned int _j
Definition: Constraint.h:35
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
VariablePhiValue _phi_secondary
Shape function on the secondary side.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
const Node *const & _current_node
current node being processed
MooseVariable & _primary_var
Primary side variable.
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
const std::map< dof_id_type, std::vector< dof_id_type > > & _node_to_elem_map
MooseMesh map of current nodes to the connected elements.
DenseMatrix< Number > _Kee
stiffness matrix holding secondary-secondary jacobian
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
void prepareVectorTagNeighbor(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags for DG and interface kernels...
virtual Real computeQpJacobian(Moose::ConstraintJacobianType type)
This is the virtual method that derived classes should override for computing the Jacobian...
std::vector< dof_id_type > _connected_dof_indices
dofs connected to the secondary node
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeResidual() override
Computes the residual Nodal residual.
DenseMatrix< Number > _Kne
stiffness matrix holding primary-secondary jacobian
ConstraintJacobianType
Definition: MooseTypes.h:792
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
void resize(const unsigned int new_m, const unsigned int new_n)
static InputParameters validParams()
void resize(unsigned int size)
Change the number of elements the array can store.
Definition: MooseArray.h:216
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
MooseVariableFieldBase & getVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:89
unsigned int n() const
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
A NodeElemConstraintBase is used when you need to create constraints between a secondary node and a p...
const VariablePhiValue & _phi_primary
Side shape function.
unsigned int _qp
Definition: Constraint.h:36
const VariableTestValue & _test_primary
Side test function.
uint8_t dof_id_type