https://mooseframework.inl.gov
PenaltyFrictionUserObject.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 
15 #include "TwoVector.h"
16 
23  virtual public WeightedVelocitiesUserObject
24 {
25 public:
27 
29 
30  virtual const ADVariableValue & contactTangentialPressureDirOne() const override;
31  virtual const ADVariableValue & contactTangentialPressureDirTwo() const override;
32 
33  virtual void initialize() override;
34  virtual void finalize() override;
35  virtual void reinit() override;
36  virtual void timestepSetup() override;
37 
38  virtual Real getFrictionalContactPressure(const Node * const node,
39  const unsigned int component) const override;
40  virtual Real getAccumulatedSlip(const Node * const node,
41  const unsigned int component) const override;
42  virtual Real getTangentialVelocity(const Node * const node,
43  const unsigned int component) const override;
44  virtual Real getDeltaTangentialLagrangeMultiplier(const Node * const node,
45  const unsigned int component) const override;
46 
47  virtual bool isAugmentedLagrangianConverged() override;
48  virtual void updateAugmentedLagrangianMultipliers() override;
49 
50 protected:
51  virtual const VariableTestValue & test() const override;
52  virtual bool constrainedByOwner() const override { return false; }
53 
55  const Real _penalty;
56 
59 
62 
65 
67  std::unordered_map<const DofObject *, std::pair<TwoVector, TwoVector>> _dof_to_step_slip;
68 
70  std::unordered_map<const DofObject *, std::pair<TwoVector, TwoVector>> _dof_to_accumulated_slip;
71 
73  std::unordered_map<const DofObject *, std::pair<ADTwoVector, TwoVector>>
75 
78 
81 
83  std::unordered_map<const DofObject *, TwoVector> _dof_to_frictional_lagrange_multipliers;
84 
86  std::unordered_map<const DofObject *, Real> _dof_to_local_penalty_friction;
87 
90 
93 
96 };
ADVariableValue _frictional_contact_traction_one
The first frictional contact pressure on the mortar segment quadrature points.
virtual void finalize() override
virtual bool constrainedByOwner() const override
ADVariableValue _frictional_contact_traction_two
The second frictional contact pressure on the mortar segment quadrature points.
std::unordered_map< const DofObject *, std::pair< TwoVector, TwoVector > > _dof_to_accumulated_slip
Map from degree of freedom to current and old accumulated slip.
virtual Real getAccumulatedSlip(const Node *const node, const unsigned int component) const override
const Real _epsilon_tolerance
Tolerance to avoid NaN/Inf in automatic differentiation operations.
virtual void updateAugmentedLagrangianMultipliers() override
virtual bool isAugmentedLagrangianConverged() override
static const std::string component
Definition: NS.h:153
virtual const ADVariableValue & contactTangentialPressureDirTwo() const override
AdaptivityFrictionalPenalty
The adaptivity method for the penalty factor at augmentations.
PenaltyFrictionUserObject(const InputParameters &parameters)
std::unordered_map< const DofObject *, TwoVector > _dof_to_frictional_lagrange_multipliers
Map from degree of freedom to augmented lagrange multiplier.
std::unordered_map< const DofObject *, std::pair< ADTwoVector, TwoVector > > _dof_to_tangential_traction
Map from degree of freedom to current and old tangential traction.
User object for computing weighted gaps and contact pressure for penalty based mortar constraints...
const Real _penalty_multiplier_friction
Penalty growth factor for augmented Lagrange.
Creates dof object to weighted tangential velocities map.
virtual void timestepSetup() override
std::unordered_map< const DofObject *, std::pair< TwoVector, TwoVector > > _dof_to_step_slip
Map from degree of freedom to current and old step slip.
User object that computes tangential pressures due to friction using a penalty approach, following J.C.
const Real _penalty
The normal penalty factor.
std::unordered_map< const DofObject *, Real > _dof_to_local_penalty_friction
Map from degree of freedom to local friction penalty value.
enum PenaltyFrictionUserObject::AdaptivityFrictionalPenalty _adaptivity_friction
OutputTools< Real >::VariableTestValue VariableTestValue
const Real _friction_coefficient
The friction coefficient.
virtual Real getDeltaTangentialLagrangeMultiplier(const Node *const node, const unsigned int component) const override
virtual const ADVariableValue & contactTangentialPressureDirOne() const override
const Real _slip_tolerance
Acceptable slip distance for augmented Lagrange convergence.
const Real _penalty_friction
The penalty factor for the frictional constraints.
static InputParameters validParams()
virtual Real getFrictionalContactPressure(const Node *const node, const unsigned int component) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual const VariableTestValue & test() const override
const InputParameters & parameters() const
virtual void initialize() override
virtual Real getTangentialVelocity(const Node *const node, const unsigned int component) const override