libMesh
systems_of_equations_ex8.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1> Systems Example 8 - "Small-sliding" contact. </h1>
21 // \author David Knezevic
22 // \date 2015
23 //
24 // In this example, we consider a linear elastic model with contact. We restrict ourselves
25 // to considering "small sliding", which means that we set the contact "load transfer" between
26 // surfaces up front, based on the undeformed mesh, and retain that load transfer throughout
27 // the solve. Even though we consider linear elasticity here, this is a nonlinear problem due
28 // to the contact condition.
29 //
30 // The contact condition is imposed using the augmented Lagrangian method, e.g. see
31 // Simo & Laursen (1992). For the sake of simplicity, in this example we assume that contact
32 // nodes are perfectly aligned (this assumption can be eliminated relatively easily).
33 //
34 // The mesh in this example consists of two disconnected cylinders. We add edge elements into
35 // the mesh in order to ensure correct parallel communication of data on contact surfaces,
36 // and also so that we do not have to manually augment the sparsity pattern.
37 
38 // We impose a displacement boundary condition of -1 in the z-direction on the top cylinder
39 // in order to impose contact.
40 
41 // C++ include files that we need
42 #include <iostream>
43 
44 // libMesh includes
45 #include "libmesh/libmesh.h"
46 #include "libmesh/replicated_mesh.h"
47 #include "libmesh/exodusII_io.h"
48 #include "libmesh/equation_systems.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/numeric_vector.h"
51 #include "libmesh/getpot.h"
52 #include "libmesh/dirichlet_boundaries.h"
53 #include "libmesh/string_to_enum.h"
54 #include "libmesh/zero_function.h"
55 #include "libmesh/nonlinear_solver.h"
56 #include "libmesh/nonlinear_implicit_system.h"
57 #include "libmesh/petsc_macro.h"
58 #include "libmesh/enum_solver_package.h"
59 #include "libmesh/mesh_refinement.h"
60 
61 // Local includes
63 
64 using namespace libMesh;
65 
66 int main (int argc, char ** argv)
67 {
68  LibMeshInit init (argc, argv);
69 
70  // This example uses an ExodusII input file
71 #ifndef LIBMESH_HAVE_EXODUS_API
72  libmesh_example_requires(false, "--enable-exodus");
73 #endif
74 
75  // We use a 3D domain.
76  libmesh_example_requires(LIBMESH_DIM > 2, "--disable-1D-only --disable-2D-only");
77 
78  // We use Dirichlet boundary conditions here
79 #ifndef LIBMESH_ENABLE_DIRICHLET
80  libmesh_example_requires(false, "--enable-dirichlet");
81 #endif
82 
83  // This example requires the PETSc nonlinear solvers
84  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
85 
86  GetPot infile("systems_of_equations_ex8.in");
87  infile.parse_command_line(argc,argv);
88 
89  const std::string approx_order = infile("approx_order", "FIRST");
90  const std::string fe_family = infile("fe_family", "LAGRANGE");
91 
92  const Real young_modulus = infile("Young_modulus", 1.0);
93  const Real poisson_ratio = infile("poisson_ratio", 0.3);
94 
95  const Real nonlinear_abs_tol = infile("nonlinear_abs_tol", 1.e-8);
96  const Real nonlinear_rel_tol = infile("nonlinear_rel_tol", 1.e-8);
97  const unsigned int nonlinear_max_its = infile("nonlinear_max_its", 50);
98  const Real contact_penalty = infile("contact_penalty", 1.e2);
99  const Real gap_function_tol = infile("gap_function_tol", 1.e-8);
100 
101  // This example code has not been written to cope with a distributed mesh
102  ReplicatedMesh mesh(init.comm());
103  mesh.read("systems_of_equations_ex8.exo");
104 
105  const unsigned int n_refinements = infile("n_refinements", 0);
106  // Skip adaptive runs on a non-adaptive libMesh build
107 #ifndef LIBMESH_ENABLE_AMR
108  libmesh_example_requires(n_refinements==0, "--enable-amr");
109 #else
110  MeshRefinement mesh_refinement(mesh);
111  mesh_refinement.uniformly_refine(n_refinements);
112 #endif
113 
114  mesh.print_info();
115 
116  EquationSystems equation_systems (mesh);
117 
118  NonlinearImplicitSystem & system =
119  equation_systems.add_system<NonlinearImplicitSystem> ("NonlinearElasticity");
120 
121  LinearElasticityWithContact le(system, contact_penalty);
122 
123  unsigned int u_var =
124  system.add_variable("u",
125  Utility::string_to_enum<Order> (approx_order),
126  Utility::string_to_enum<FEFamily>(fe_family));
127 
128  unsigned int v_var =
129  system.add_variable("v",
130  Utility::string_to_enum<Order> (approx_order),
131  Utility::string_to_enum<FEFamily>(fe_family));
132 
133  unsigned int w_var =
134  system.add_variable("w",
135  Utility::string_to_enum<Order> (approx_order),
136  Utility::string_to_enum<FEFamily>(fe_family));
137 
138  // Also, initialize an ExplicitSystem to store stresses
139  ExplicitSystem & stress_system =
140  equation_systems.add_system<ExplicitSystem> ("StressSystem");
141 
142  stress_system.add_variable("sigma_00", CONSTANT, MONOMIAL);
143  stress_system.add_variable("sigma_01", CONSTANT, MONOMIAL);
144  stress_system.add_variable("sigma_02", CONSTANT, MONOMIAL);
145  stress_system.add_variable("sigma_11", CONSTANT, MONOMIAL);
146  stress_system.add_variable("sigma_12", CONSTANT, MONOMIAL);
147  stress_system.add_variable("sigma_22", CONSTANT, MONOMIAL);
148  stress_system.add_variable("vonMises", CONSTANT, MONOMIAL);
149 
150  equation_systems.parameters.set<Real> ("nonlinear solver absolute residual tolerance") = nonlinear_abs_tol;
151  equation_systems.parameters.set<Real> ("nonlinear solver relative residual tolerance") = nonlinear_rel_tol;
152  equation_systems.parameters.set<unsigned int> ("nonlinear solver maximum iterations") = nonlinear_max_its;
153 
154  system.nonlinear_solver->residual_and_jacobian_object = &le;
155 
156  equation_systems.parameters.set<Real>("young_modulus") = young_modulus;
157  equation_systems.parameters.set<Real>("poisson_ratio") = poisson_ratio;
158 
159 #ifdef LIBMESH_ENABLE_DIRICHLET
160  // Attach Dirichlet (Clamped) boundary conditions
162 
163  // Most DirichletBoundary users will want to supply a "locally
164  // indexed" functor
166  (DirichletBoundary ({MIN_Z_BOUNDARY}, {u_var, v_var, w_var},
168 
170  (DirichletBoundary ({MAX_Z_BOUNDARY}, {u_var, v_var}, zero,
172 
173  ConstFunction<Number> neg_one(-1.);
174 
176  (DirichletBoundary ({MAX_Z_BOUNDARY}, {w_var}, neg_one,
178 #endif // LIBMESH_ENABLE_DIRICHLET
179  libmesh_ignore(u_var, v_var, w_var);
180 
181  le.initialize_contact_load_paths();
182 
183  libMesh::out << "Mesh before adding edge connectors" << std::endl;
184  mesh.print_info();
185  le.add_contact_edge_elements();
186 
187  libMesh::out << "Mesh after adding edge connectors" << std::endl;
188  mesh.print_info();
189 
190  equation_systems.init();
191  equation_systems.print_info();
192 
193  libMesh::out << "Contact penalty: " << contact_penalty << std::endl << std::endl;
194 
195  Real current_max_gap_function = std::numeric_limits<Real>::max();
196 
197  const unsigned int max_outer_iterations = 10;
198  for (unsigned int outer_iteration = 0;
199  outer_iteration != max_outer_iterations; ++outer_iteration)
200  {
201  if (current_max_gap_function <= gap_function_tol)
202  {
203  libMesh::out << "Outer loop converged" << std::endl;
204  break;
205  }
206 
207  libMesh::out << "Starting outer iteration " << outer_iteration << std::endl;
208 
209  // Perform inner iteration (i.e. Newton's method loop)
210  system.solve();
211  system.nonlinear_solver->print_converged_reason();
212 
213  // Perform augmented Lagrangian update
214  le.update_lambdas();
215 
216  std::pair<Real, Real> least_and_max_gap_function = le.get_least_and_max_gap_function();
217  Real least_gap_fn = least_and_max_gap_function.first;
218  Real max_gap_fn = least_and_max_gap_function.second;
219 
220  libMesh::out << "Finished outer iteration, least gap function: "
221  << least_gap_fn
222  << ", max gap function: "
223  << max_gap_fn
224  << std::endl
225  << std::endl;
226 
227  current_max_gap_function = std::max(std::abs(least_gap_fn), std::abs(max_gap_fn));
228  }
229 
230  // Enforcing constraints imposes the non-zero Dirichlet constraints on the solution.
231  system.get_dof_map().enforce_constraints_exactly(system);
232  system.update();
233 
234  libMesh::out << "Computing stresses..." << std::endl;
235 
236  le.compute_stresses();
237 
238  std::stringstream filename;
239  filename << "solution.exo";
240  ExodusII_IO (mesh).write_equation_systems(filename.str(),
241  equation_systems);
242 
243  return 0;
244 }
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
std::unique_ptr< NonlinearSolver< Number > > nonlinear_solver
The NonlinearSolver defines the default interface used to solve the nonlinear_implicit system...
ConstFunction that simply returns 0.
Definition: zero_function.h:38
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
int main(int argc, char **argv)
This class encapsulate all functionality required for assembling and solving a linear elastic model w...
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
The libMesh namespace provides an interface to certain functionality in the library.
const Number zero
.
Definition: libmesh.h:304
SolverPackage default_solver_package()
Definition: libmesh.C:1117
void libmesh_ignore(const Args &...)
Implements (adaptive) mesh refinement algorithms for a MeshBase.
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1357
virtual void solve() override
Assembles & solves the nonlinear system R(x) = 0.
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:510
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
Parameters parameters
Data structure holding arbitrary parameters.
Function that returns a single value that never changes.
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
const DofMap & get_dof_map() const
Definition: system.h:2374
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems...
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const
Constrains the numeric vector v, which represents a solution defined on the mesh. ...
Definition: dof_map.h:2350