https://mooseframework.inl.gov
NodeElemConstraint.h
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 #pragma once
11 
12 // MOOSE includes
13 #include "NodeElemConstraintBase.h"
14 
21 {
22 public:
24 
26 
28  virtual void computeResidual() override;
29 
31  virtual void computeJacobian() override;
32 
34  virtual void computeOffDiagJacobian(unsigned int jvar) override;
35 
37  void getConnectedDofIndices(unsigned int var_num);
38 
39 protected:
41  virtual void prepareSecondaryToPrimaryMap() = 0;
42 
45 
48 
51  unsigned int /*jvar*/)
52  {
53  return 0;
54  }
55 
58 
61 
64 
67 
70 };
OutputTools< Real >::VariableGradient VariableGradient
Definition: MooseTypes.h:316
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)
ConstraintType
Definition: MooseTypes.h:758
virtual Real computeQpResidual(Moose::ConstraintType type)=0
This is the virtual method that derived classes should override for computing the residual...
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:127
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...
const VariableValue & _u_secondary
Value of the unknown variable on the secondary node.
void getConnectedDofIndices(unsigned int var_num)
Gets the indices for all dofs connected to the constraint.
const VariableValue & _u_primary
Holds the current solution at the current quadrature point.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
const VariableTestGradient & _grad_test_primary
Gradient of side shape function.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:315
forward declarations
virtual Real computeQpJacobian(Moose::ConstraintJacobianType type)
This is the virtual method that derived classes should override for computing the Jacobian...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableGradient & _grad_u_primary
Holds the current solution gradient at the current quadrature point.
virtual void computeResidual() override
Computes the residual Nodal residual.
ConstraintJacobianType
Definition: MooseTypes.h:797
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:326
A NodeElemConstraintBase is a base class for constraints between a subdomain with secondary nodes and...
A NodeElemConstraintBase is used when you need to create constraints between a secondary node and a p...
virtual void prepareSecondaryToPrimaryMap()=0
prepare the _secondary_to_primary_map (see NodeElemConstraintBase)
const VariablePhiGradient & _grad_phi_primary
Gradient of side shape function.