www.mooseframework.org
ElemElemConstraint.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #ifndef ELEMELEMCONSTRAINT_H
16 #define ELEMELEMCONSTRAINT_H
17 
18 // MOOSE includes
19 #include "Constraint.h"
21 
22 // Forward Declarations
23 class ElemElemConstraint;
24 class ElementPairInfo;
25 class FEProblemBase;
26 
27 template <>
29 
32 {
33 public:
35 
39  virtual void reinit(const ElementPairInfo & element_pair_info);
40 
44  virtual void reinitConstraintQuadrature(const ElementPairInfo & element_pair_info);
45 
50 
54  virtual void computeResidual();
55 
60 
64  virtual void computeJacobian();
65 
66 protected:
68  unsigned int _dim;
69 
70  const Elem *& _current_elem;
71 
73  const Elem *& _neighbor_elem;
74 
76  std::vector<Point> _constraint_q_point;
78  std::vector<Real> _constraint_weight;
79 
81  unsigned int _i, _j;
82 
84  const VariableValue & _u;
85 
92 
97 
102 
107 
112 
117  virtual Real computeQpResidual(Moose::DGResidualType type) = 0;
118 
123  virtual Real computeQpJacobian(Moose::DGJacobianType type) = 0;
124 };
125 
126 #endif /* ELEMELEMCONSTRAINT_H */
std::vector< Real > _constraint_weight
Weights of quadrature points used in integration of constraint.
FEProblemBase & _fe_problem
const VariablePhiGradient & _grad_phi_neighbor
Gradient of neighbor shape function.
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
const Elem *& _neighbor_elem
The neighboring element.
const Elem *& _current_elem
virtual void computeJacobian()
Computes the jacobian for the current side.
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:42
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:184
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
ElemElemConstraint(const InputParameters &parameters)
const VariableTestValue & _test_neighbor
Neighbor test function.
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)
Computes the element/neighbor-element/neighbor Jacobian.
virtual Real computeQpJacobian(Moose::DGJacobianType type)=0
Compute the Jacobian for one of the constraint quadrature points.
virtual void reinitConstraintQuadrature(const ElementPairInfo &element_pair_info)
Set information needed for constraint integration.
const VariableTestValue & _test
Test function.
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.
virtual void computeResidual()
Computes the residual for the current side.
DGJacobianType
Definition: MooseTypes.h:190
virtual void computeElemNeighResidual(Moose::DGResidualType type)
Computes the residual for this element or the neighbor.
MatType type
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 InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:53
const VariablePhiValue & _phi
Shape function.
InputParameters validParams< ElemElemConstraint >()
unsigned int _i
Indices for looping over DOFs.