www.mooseframework.org
FaceFaceConstraint.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 FACEFACECONSTRAINT_H
16 #define FACEFACECONSTRAINT_H
17 
18 // MOOSE includes
19 #include "Constraint.h"
21 #include "MooseMesh.h"
22 
23 // Forward Declarations
24 class FaceFaceConstraint;
25 class FEProblemBase;
26 
27 template <>
29 
47 {
48 public:
50 
54  virtual void reinit();
55  virtual void reinitSide(Moose::ConstraintType res_type);
56 
60  virtual void computeResidual();
65  virtual void computeResidualSide(Moose::ConstraintType side);
66 
70  virtual void computeJacobian();
75  virtual void computeJacobianSide(Moose::ConstraintType side);
76 
77 protected:
78  virtual Real computeQpResidual() = 0;
79  virtual Real computeQpResidualSide(Moose::ConstraintType res_type) = 0;
80  virtual Real computeQpJacobian();
82 
84  unsigned int _dim;
85 
90 
92  QBase *& _qrule;
95 
96  std::vector<Real> _JxW_lm;
97 
101  const Elem *& _current_elem;
102 
103  std::vector<std::vector<Real>> _test;
104  std::vector<std::vector<Real>> _phi;
105 
108 
113 
117 
121  std::vector<Real> _u_master;
122  std::vector<RealGradient> _grad_u_master;
123 
127  std::vector<Point> _phys_points_master;
131  const Elem * _elem_master;
137 
142 
146  std::vector<Real> _u_slave;
147  std::vector<RealGradient> _grad_u_slave;
148 
152  std::vector<Point> _phys_points_slave;
156  const Elem * _elem_slave;
162 
167 };
168 
169 #endif /* FACEFACECONSTRAINT_H */
const VariablePhiValue & _phi_slave
Values of shape function on the slave side.
const VariableTestValue & _test_slave
Values of test functions on the slave side.
ConstraintType
Definition: MooseTypes.h:198
const VariablePhiValue & _phi_master
Values of shape function on the master side.
BoundaryID _slave
Boundary ID for the slave surface.
const MooseArray< Real > & _coord
virtual Real computeQpJacobianSide(Moose::ConstraintJacobianType jac_type)
Class for stuff related to variables.
Definition: MooseVariable.h:43
const VariableValue & _lambda
The values of Lagrange multipliers in quadrature points.
virtual Real computeQpJacobian()
User for mortar methods.
BoundaryID _master
Boundary ID for the master surface.
Base class for all Constraint types.
Definition: Constraint.h:42
const Elem * _elem_master
Element on the master side.
std::vector< RealGradient > _grad_u_slave
std::vector< Real > _JxW_lm
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
PenetrationLocator & _slave_penetration_locator
std::vector< std::vector< Real > > _phi
std::vector< Point > _phys_points_slave
Physical points on the slave side.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
std::vector< Real > _u_master
Values of the constrained variable on the master side.
virtual Real computeQpResidualSide(Moose::ConstraintType res_type)=0
MooseVariable & _master_var
virtual void computeResidual()
Computes the residual for the current element.
std::vector< Real > _u_slave
Values of the constrained variable on the slave side.
virtual void reinitSide(Moose::ConstraintType res_type)
InputParameters validParams< FaceFaceConstraint >()
const MooseArray< Point > & _q_point
FEProblemBase & _fe_problem
MooseMesh::MortarInterface & _iface
const VariableTestGradient & _grad_test_slave
const MooseArray< Real > & _JxW
virtual Real computeQpResidual()=0
const VariableTestGradient & _grad_test_master
virtual void computeJacobian()
Computes the Jacobian for the current element (i.e element of the Mortar interface).
const Elem *& _current_elem
Current element on the interface (i.e in the mortar space)
const Elem * _elem_slave
Element on the master side.
const VariableTestValue & _test_master
Values of test functions on the master side.
std::vector< RealGradient > _grad_u_master
virtual void reinit()
Evaluate variables, compute q-points, etc.
ConstraintJacobianType
Definition: MooseTypes.h:204
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:53
FaceFaceConstraint(const InputParameters &parameters)
virtual void computeResidualSide(Moose::ConstraintType side)
Computes residual contributions from master or slave side.
std::vector< Point > _phys_points_master
Physical points on the master side.
PenetrationLocator & _master_penetration_locator
MooseVariable & _slave_var
std::vector< std::vector< Real > > _test
boundary_id_type BoundaryID
Definition: MooseTypes.h:75
virtual void computeJacobianSide(Moose::ConstraintType side)
Computes Jacobian contributions from master or slave side.