https://mooseframework.inl.gov
CoupledTiedValueConstraint.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 
11 
12 // MOOSE includes
13 #include "MooseVariableFE.h"
14 #include "SystemBase.h"
15 
16 #include "libmesh/sparse_matrix.h"
17 
19 
22 {
24  params.addClassDescription(
25  "Requires the value of two variables to be the consistent on both sides of an interface.");
26  params.addParam<Real>("scaling", 1, "scaling factor to be applied to constraint equations");
27  params.set<bool>("use_displaced_mesh") = true;
28  return params;
29 }
30 
32  : NodeFaceConstraint(parameters),
33  _scaling(getParam<Real>("scaling")),
34  _residual_copy(_sys.residualGhosted())
35 {
36 }
37 
38 Real
40 {
41  return _u_primary[_qp];
42 }
43 
44 Real
46 {
47  Real scaling_factor = _var.scalingFactor();
48  Real secondary_resid = 0;
49  Real retVal = 0;
50  switch (type)
51  {
52  case Moose::Secondary:
54  break;
55  case Moose::Primary:
56  secondary_resid =
57  _residual_copy(_current_node->dof_number(0, _var.number(), 0)) / scaling_factor;
58  retVal = secondary_resid * _test_primary[_i][_qp];
59  break;
60  default:
61  break;
62  }
63  return retVal;
64 }
65 
66 Real
68 {
69  Real scaling_factor = _var.scalingFactor();
70  Real secondary_jac = 0;
71  Real retVal = 0;
72  switch (type)
73  {
76  break;
78  retVal = 0;
79  break;
81  secondary_jac =
82  (*_jacobian)(_current_node->dof_number(0, _var.number(), 0), _connected_dof_indices[_j]);
83  retVal = secondary_jac * _test_primary[_i][_qp] / scaling_factor;
84  break;
86  retVal = 0;
87  break;
88  default:
89  mooseError("Unsupported type");
90  break;
91  }
92  return retVal;
93 }
94 
95 Real
97  unsigned int jvar)
98 {
99  Real retVal = 0;
100 
101  if (jvar == _primary_var_num)
102  {
103  switch (type)
104  {
106  retVal = 0;
107  break;
109  retVal = -_phi_primary[_j][_qp] * _test_secondary[_i][_qp] * _scaling;
110  break;
112  retVal = 0;
113  break;
115  retVal = 0;
116  break;
117  default:
118  mooseError("Unsupported type");
119  break;
120  }
121  }
122 
123  return retVal;
124 }
const VariableValue & _u_primary
Holds the current solution at the current quadrature point.
ConstraintType
Definition: MooseTypes.h:758
unsigned int _primary_var_num
Number for the primary variable.
A CoupledTiedValueConstraint forces the value of a variable to be the same on both sides of an interf...
MooseVariable & _var
unsigned int number() const
Get variable number coming from libMesh.
CoupledTiedValueConstraint(const InputParameters &parameters)
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
registerMooseObject("MooseApp", CoupledTiedValueConstraint)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
unsigned int _i
Definition: Constraint.h:35
const VariableTestValue & _test_primary
Side test function.
virtual Real computeQpResidual(Moose::ConstraintType type) override
This is the virtual that derived classes should override for computing the residual on neighboring el...
const Node *const & _current_node
current node being processed
unsigned int _j
Definition: Constraint.h:35
A NodeFaceConstraint is used when you need to create constraints between two surfaces in a mesh...
NumericVector< Number > & _residual_copy
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
const VariableValue & _u_secondary
Value of the unknown variable this BC is action on.
VariableTestValue _test_secondary
Shape function on the secondary side. This will always only have one entry and that entry will always...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
ConstraintJacobianType
Definition: MooseTypes.h:797
virtual Real computeQpJacobian(Moose::ConstraintJacobianType type) override
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
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
std::vector< dof_id_type > _connected_dof_indices
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
virtual Real computeQpSecondaryValue() override
Compute the value the secondary node should have at the beginning of a timestep.
const VariablePhiValue & _phi_primary
Side shape function.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
VariablePhiValue _phi_secondary
Shape function on the secondary side. This will always.
virtual Real computeQpOffDiagJacobian(Moose::ConstraintJacobianType type, unsigned int jvar) override
This is the virtual that derived classes should override for computing the off-diag Jacobian...
unsigned int _qp
Definition: Constraint.h:36
static InputParameters validParams()
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.