libMesh
linear_elasticity_with_contact.h
Go to the documentation of this file.
1 #ifndef LINEAR_ELASTICITY_WITH_CONTACT_H
2 #define LINEAR_ELASTICITY_WITH_CONTACT_H
3 
4 #include "libmesh/dof_map.h"
5 #include "libmesh/nonlinear_implicit_system.h"
6 
7 // define the subdomain IDs
8 #define TOP_SUBDOMAIN 2
9 #define BOTTOM_SUBDOMAIN 1
10 
11 // define the boundary IDs in the mesh
12 #define MIN_Z_BOUNDARY 1
13 #define MAX_Z_BOUNDARY 2
14 #define CONTACT_BOUNDARY_LOWER 3
15 #define CONTACT_BOUNDARY_UPPER 4
16 
17 using libMesh::DofMap;
20 using libMesh::Point;
21 using libMesh::Real;
22 using libMesh::Number;
23 using libMesh::MeshBase;
26 
32  public NonlinearImplicitSystem::ComputeResidualandJacobian
33 {
34 private:
35 
40 
45 
50  std::map<dof_id_type, Real> _lambda_plus_penalty_values;
51 
56  std::map<dof_id_type, Real> _lambdas;
57 
61  std::map<dof_id_type, dof_id_type> _contact_node_map;
62 
63 public:
64 
69  Real contact_penalty_in);
70 
74  void set_contact_penalty(Real contact_penalty_in);
75 
79  Real get_contact_penalty() const;
80 
84  Real kronecker_delta(unsigned int i,
85  unsigned int j);
86 
90  Real elasticity_tensor(Real young_modulus,
91  Real poisson_ratio,
92  unsigned int i,
93  unsigned int j,
94  unsigned int k,
95  unsigned int l);
96 
101  void move_mesh(MeshBase & input_mesh,
102  const NumericVector<Number> & input_solution);
103 
108 
117 
121  virtual void residual_and_jacobian (const NumericVector<Number> & soln,
122  NumericVector<Number> * residual,
123  SparseMatrix<Number> * jacobian,
124  NonlinearImplicitSystem & /*sys*/);
125 
129  void compute_stresses();
130 
137  std::pair<Real, Real> update_lambdas();
138 
142  std::pair<Real, Real> get_least_and_max_gap_function();
143 };
144 
145 #endif
std::map< dof_id_type, Real > _lambda_plus_penalty_values
Store the intermediate values of lambda plus penalty.
Real _contact_penalty
Penalize overlapping elements.
This class encapsulate all functionality required for assembling and solving a linear elastic model w...
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:43
void initialize_contact_load_paths()
Set up the load paths on the contact surfaces.
Real get_contact_penalty() const
Get the penalty parameter.
This is the MeshBase class.
Definition: mesh_base.h:74
void move_mesh(MeshBase &input_mesh, const NumericVector< Number > &input_solution)
Move the mesh nodes of input_mesh based on the displacement field in input_solution.
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
NonlinearImplicitSystem & _sys
Keep a reference to the NonlinearImplicitSystem.
void add_contact_edge_elements()
Add edge elements into the mesh based on the contact load paths.
std::map< dof_id_type, Real > _lambdas
Augmented Lagrangian values at each contact node.
virtual void residual_and_jacobian(const NumericVector< Number > &soln, NumericVector< Number > *residual, SparseMatrix< Number > *jacobian, NonlinearImplicitSystem &)
Evaluate the Jacobian of the nonlinear system.
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
std::pair< Real, Real > update_lambdas()
Update the lambda parameters in the augmented Lagrangian method.
Real elasticity_tensor(Real young_modulus, Real poisson_ratio, unsigned int i, unsigned int j, unsigned int k, unsigned int l)
Evaluate the fourth order tensor (C_ijkl) that relates stress to strain.
void set_contact_penalty(Real contact_penalty_in)
Update the penalty parameter.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::map< dof_id_type, dof_id_type > _contact_node_map
This provides a map between contact nodes.
void compute_stresses()
Compute the Cauchy stress for the current solution.
LinearElasticityWithContact(NonlinearImplicitSystem &sys_in, Real contact_penalty_in)
Constructor.
Real kronecker_delta(unsigned int i, unsigned int j)
Kronecker delta function.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
std::pair< Real, Real > get_least_and_max_gap_function()
uint8_t dof_id_type
Definition: id_types.h:67