www.mooseframework.org
ElemElemConstraint.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "Constraint.h"
15 
16 class ElementPairInfo;
17 class FEProblemBase;
18 
22 {
23 public:
25 
27 
31  virtual void reinit(const ElementPairInfo & element_pair_info);
32 
36  virtual void reinitConstraintQuadrature(const ElementPairInfo & element_pair_info);
37 
42 
46  virtual void computeResidual() override;
47 
52 
56  virtual void computeJacobian() override;
57 
61  unsigned int getInterfaceID() const { return _interface_id; };
62 
66  const MooseVariable & variable() const override { return _var; }
67 
68 protected:
70  unsigned int _dim;
71 
72  unsigned int _interface_id;
73 
75 
76  const Elem * const & _current_elem;
77 
79  const Elem * const & _neighbor_elem;
80 
82  std::vector<Point> _constraint_q_point;
84  std::vector<Real> _constraint_weight;
85 
87  unsigned int _i, _j;
88 
90  const VariableValue & _u;
91 
98 
103 
108 
113 
118 
124 
130 };
OutputTools< Real >::VariableGradient VariableGradient
Definition: MooseTypes.h:303
virtual void computeJacobian() override
Computes the jacobian for the current side.
MooseVariable & _var
std::vector< Real > _constraint_weight
Weights of quadrature points used in integration of constraint.
virtual void computeResidual() override
Computes the residual for the current side.
FEProblemBase & _fe_problem
const VariablePhiGradient & _grad_phi_neighbor
Gradient of neighbor shape function.
Class for stuff related to variables.
Definition: Adaptivity.h:31
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
unsigned int _interface_id
std::vector< Point > _constraint_q_point
Quadrature points used in integration of constraint.
const VariableValue & _u_neighbor
Holds the current solution at the current quadrature point on the neighbor element.
Base class for all Constraint types.
Definition: Constraint.h:19
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const VariableTestGradient & _grad_test
Gradient of test function.
DGResidualType
Definition: MooseTypes.h:656
unsigned int getInterfaceID() const
Get the interface ID.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
ElemElemConstraint(const InputParameters &parameters)
Enhances MooseVariableInterface interface provide values from neighbor elements.
const Elem *const & _current_elem
const VariableTestValue & _test_neighbor
Neighbor test function.
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)
Computes the element/neighbor-element/neighbor Jacobian.
static InputParameters validParams()
virtual Real computeQpJacobian(Moose::DGJacobianType type)=0
Compute the Jacobian for one of the constraint quadrature points.
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:312
virtual void reinitConstraintQuadrature(const ElementPairInfo &element_pair_info)
Set information needed for constraint integration.
const VariableTestValue & _test
Test function.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const VariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point.
const VariableGradient & _grad_u_neighbor
Holds the current solution gradient at the current quadrature point on the neighbor element...
This is the ElementPairInfo class.
const VariablePhiValue & _phi_neighbor
Neighbor shape function.
const VariableTestGradient & _grad_test_neighbor
Gradient of neighbor shape function.
virtual Real computeQpResidual(Moose::DGResidualType type)=0
Compute the residual for one of the constraint quadrature points.
const Elem *const & _neighbor_elem
The neighboring element.
DGJacobianType
Definition: MooseTypes.h:662
OutputTools< Real >::VariablePhiGradient VariablePhiGradient
Definition: MooseTypes.h:308
virtual void computeElemNeighResidual(Moose::DGResidualType type)
Computes the residual for this element or the neighbor.
forward declarations
Definition: MooseArray.h:17
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableValue & _u
Holds the current solution at the current quadrature point.
virtual void reinit(const ElementPairInfo &element_pair_info)
reinit element-element constraint
const VariablePhiGradient & _grad_phi
Shape function gradient.
const VariablePhiValue & _phi
Shape function.
const InputParameters & parameters() const
Get the parameters of the object.
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:313
unsigned int _i
Indices for looping over DOFs.
const MooseVariable & variable() const override
The variable number that this object operates on.