https://mooseframework.inl.gov
ADNodeElemConstraint.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 "ADNodeElemConstraint.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseEnum.h"
15 #include "MooseMesh.h"
16 #include "MooseTypes.h"
17 
18 #include "libmesh/string_to_enum.h"
19 
22 {
24  return params;
25 }
26 
28  : NodeElemConstraintBase(parameters),
29  _u_primary(_primary_var.adSlnNeighbor()),
30  _u_secondary(_var.adDofValues())
31 {
32 }
33 
34 void
36 {
37  _qp = 0;
38 
39  _residuals.resize(_test_primary.size(), 0);
40  for (auto & r : _residuals)
41  r = 0;
42 
43  for (_i = 0; _i < _test_primary.size(); _i++)
45 
48 
49  _residuals.resize(_test_secondary.size(), 0);
50 
51  for (auto & r : _residuals)
52  r = 0;
53 
54  for (_i = 0; _i < _test_secondary.size(); _i++)
56 
58 }
59 
60 void
62 {
63  _qp = 0;
64 
65  std::vector<ADReal> primary_residual(_test_primary.size(), 0);
66 
67  for (_i = 0; _i < _test_primary.size(); _i++)
68  primary_residual[_i] += computeQpResidual(Moose::Primary);
69 
72 
73  std::vector<ADReal> secondary_residual(_test_secondary.size(), 0);
74 
75  for (_i = 0; _i < _test_secondary.size(); _i++)
76  secondary_residual[_i] += computeQpResidual(Moose::Secondary);
77 
78  addJacobian(_assembly, secondary_residual, _var.dofIndices(), _var.scalingFactor());
79 }
80 
81 Real
83 {
84  mooseAssert(false, "Should not be used");
85  return 0;
86 }
87 
88 Real
90  unsigned int /*jvar*/)
91 {
92  mooseAssert(false, "Should not be used");
93  return 0;
94 }
virtual void computeJacobian() override
Computes the jacobian for the current element.
virtual Real computeQpJacobian(Moose::ConstraintJacobianType)
This is the virtual method that derived classes should override for computing the Jacobian...
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.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
Definition: Constraint.h:35
VariableTestValue _test_secondary
Shape function on the secondary side. This will always only have one entry and that entry will always...
MooseVariable & _var
secondary node variable
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided residual derivatives into the Jacobian for the provided dof indices.
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
virtual ADReal computeQpResidual(Moose::ConstraintType type)=0
This is the virtual method that derived classes should override for computing the residual...
std::vector< Real > _residuals
MooseVariable & _primary_var
Primary side variable.
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
ConstraintJacobianType
Definition: MooseTypes.h:796
static InputParameters validParams()
virtual void computeResidual() override
Computes the Nodal residual.
virtual Real computeQpOffDiagJacobian(Moose::ConstraintJacobianType, unsigned int)
This is the virtual method that derived classes should override for computing the off-diag Jacobian...
ADNodeElemConstraint(const InputParameters &parameters)
static InputParameters validParams()
A NodeElemConstraintBase is used when you need to create constraints between a secondary node and a p...
unsigned int _qp
Definition: Constraint.h:36
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
const VariableTestValue & _test_primary
Side test function.