libMesh
systems_of_equations_ex5.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 5 - Linear Elastic Cantilever with Constraint </h1>
21 // \author David Knezevic
22 // \date 2012
23 //
24 // In this example we extend systems_of_equations_ex4 to enforce a constraint.
25 // We apply a uniform load on the top surface of the cantilever, and we
26 // determine the traction on the right boundary in order to obtain zero
27 // average vertical displacement on the right boundary of the domain.
28 //
29 // This constraint is enforced via a Lagrange multiplier (SCALAR variable).
30 // The system we solve, therefore, is of the form:
31 // a(u,v) + \lambda g(v) = f(v)
32 // g(u) = 0
33 // Here \lambda tells us the traction required to satisfy the constraint.
34 
35 // C++ include files that we need
36 #include <iostream>
37 #include <algorithm>
38 #include <math.h>
39 
40 // libMesh includes
41 #include "libmesh/libmesh.h"
42 #include "libmesh/mesh.h"
43 #include "libmesh/mesh_generation.h"
44 #include "libmesh/exodusII_io.h"
45 #include "libmesh/linear_implicit_system.h"
46 #include "libmesh/equation_systems.h"
47 #include "libmesh/fe.h"
48 #include "libmesh/quadrature_gauss.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/sparse_matrix.h"
51 #include "libmesh/numeric_vector.h"
52 #include "libmesh/dense_matrix.h"
53 #include "libmesh/dense_submatrix.h"
54 #include "libmesh/dense_vector.h"
55 #include "libmesh/dense_subvector.h"
56 #include "libmesh/perf_log.h"
57 #include "libmesh/elem.h"
58 #include "libmesh/boundary_info.h"
59 #include "libmesh/zero_function.h"
60 #include "libmesh/dirichlet_boundaries.h"
61 #include "libmesh/string_to_enum.h"
62 #include "libmesh/getpot.h"
63 #include "libmesh/enum_solver_package.h"
64 
65 // Bring in everything from the libMesh namespace
66 using namespace libMesh;
67 
68 // Matrix and right-hand side assemble
70  const std::string & system_name);
71 
72 // Define the elasticity tensor, which is a fourth-order tensor
73 // i.e. it has four indices i, j, k, l
74 Real eval_elasticity_tensor(unsigned int i,
75  unsigned int j,
76  unsigned int k,
77  unsigned int l);
78 
79 // Begin the main program.
80 int main (int argc, char ** argv)
81 {
82  // Initialize libMesh and any dependent libraries
83  LibMeshInit init (argc, argv);
84 
85  // This example requires a linear solver package.
86  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
87  "--enable-petsc, --enable-trilinos, or --enable-eigen");
88 
89  // Initialize the cantilever mesh
90  const unsigned int dim = 2;
91 
92  // Skip this 2D example if libMesh was compiled as 1D-only.
93  libmesh_example_requires(dim <= LIBMESH_DIM, "2D support");
94 
95  // We use Dirichlet boundary conditions here
96 #ifndef LIBMESH_ENABLE_DIRICHLET
97  libmesh_example_requires(false, "--enable-dirichlet");
98 #endif
99 
100  // Get the mesh size from the command line.
101  GetPot command_line (argc, argv);
102 
103  const int nx = libMesh::command_line_next("-nx", 50),
104  ny = libMesh::command_line_next("-ny", 10);
105 
106  // Create a 2D mesh distributed across the default MPI communicator.
107  Mesh mesh(init.comm(), dim);
109  nx, ny,
110  0., 1.,
111  0., 0.2,
112  QUAD9);
113 
114  // Print information about the mesh to the screen.
115  mesh.print_info();
116 
117  // Create an equation systems object.
118  EquationSystems equation_systems (mesh);
119 
120  // Declare the system and its variables.
121  // Create a system named "Elasticity"
122  LinearImplicitSystem & system =
123  equation_systems.add_system<LinearImplicitSystem> ("Elasticity");
124 
125 #ifdef LIBMESH_ENABLE_DIRICHLET
126  // Add two displacement variables, u and v, to the system
127  unsigned int u_var = system.add_variable("u", SECOND, LAGRANGE);
128  unsigned int v_var = system.add_variable("v", SECOND, LAGRANGE);
129 
130  // Add a SCALAR variable for the Lagrange multiplier to enforce our constraint
131  system.add_variable("lambda", FIRST, SCALAR);
132 
134 
135  // Create a ZeroFunction to initialize dirichlet_bc
136  ZeroFunction<> zf;
137 
138  // Construct a Dirichlet boundary condition object
139  // We impose a "clamped" boundary condition on the
140  // "left" boundary, i.e. bc_id = 3
141 
142  // Most DirichletBoundary users will want to supply a "locally
143  // indexed" functor
144  DirichletBoundary dirichlet_bc({3}, {u_var, v_var}, zf,
146 
147  // We must add the Dirichlet boundary condition _before_
148  // we call equation_systems.init()
149  system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
150 #endif // LIBMESH_ENABLE_DIRICHLET
151 
152  // Initialize the data structures for the equation system.
153  equation_systems.init();
154 
155  // Print information about the system to the screen.
156  equation_systems.print_info();
157 
158  // Solve the system
159  system.solve();
160 
161  // Plot the solution
162 #ifdef LIBMESH_HAVE_EXODUS_API
163  ExodusII_IO (mesh).write_equation_systems("displacement.e", equation_systems);
164 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
165 
166  // All done.
167  return 0;
168 }
169 
170 
172  const std::string & libmesh_dbg_var(system_name))
173 {
174  libmesh_assert_equal_to (system_name, "Elasticity");
175 
176  const MeshBase & mesh = es.get_mesh();
177 
178  const unsigned int dim = mesh.mesh_dimension();
179 
180  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Elasticity");
181 
182  const unsigned int u_var = system.variable_number ("u");
183  const unsigned int v_var = system.variable_number ("v");
184  const unsigned int lambda_var = system.variable_number ("lambda");
185 
186  const DofMap & dof_map = system.get_dof_map();
187  FEType fe_type = dof_map.variable_type(0);
188  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
189  QGauss qrule (dim, fe_type.default_quadrature_order());
190  fe->attach_quadrature_rule (&qrule);
191 
192  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
193  QGauss qface(dim-1, fe_type.default_quadrature_order());
194  fe_face->attach_quadrature_rule (&qface);
195 
196  const std::vector<Real> & JxW = fe->get_JxW();
197  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
198 
201 
203  Kuu(Ke), Kuv(Ke),
204  Kvu(Ke), Kvv(Ke);
205  DenseSubMatrix<Number> Klambda_v(Ke), Kv_lambda(Ke);
206 
208  Fu(Fe),
209  Fv(Fe);
210 
211  std::vector<dof_id_type> dof_indices;
212  std::vector<dof_id_type> dof_indices_u;
213  std::vector<dof_id_type> dof_indices_v;
214  std::vector<dof_id_type> dof_indices_lambda;
215 
216  SparseMatrix<Number> & matrix = system.get_system_matrix();
217 
218  for (const auto & elem : mesh.active_local_element_ptr_range())
219  {
220  dof_map.dof_indices (elem, dof_indices);
221  dof_map.dof_indices (elem, dof_indices_u, u_var);
222  dof_map.dof_indices (elem, dof_indices_v, v_var);
223  dof_map.dof_indices (elem, dof_indices_lambda, lambda_var);
224 
225  const unsigned int n_dofs = dof_indices.size();
226  const unsigned int n_u_dofs = dof_indices_u.size();
227  const unsigned int n_v_dofs = dof_indices_v.size();
228  const unsigned int n_lambda_dofs = dof_indices_lambda.size();
229 
230  fe->reinit (elem);
231 
232  Ke.resize (n_dofs, n_dofs);
233  Fe.resize (n_dofs);
234 
235  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
236  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
237 
238  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
239  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
240 
241  // Also, add a row and a column to enforce the constraint
242  Kv_lambda.reposition (v_var*n_u_dofs, v_var*n_u_dofs+n_v_dofs, n_v_dofs, 1);
243  Klambda_v.reposition (v_var*n_v_dofs+n_v_dofs, v_var*n_v_dofs, 1, n_v_dofs);
244 
245  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
246  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
247 
248  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
249  {
250  for (unsigned int i=0; i<n_u_dofs; i++)
251  for (unsigned int j=0; j<n_u_dofs; j++)
252  {
253  // Tensor indices
254  unsigned int C_i, C_j, C_k, C_l;
255  C_i=0, C_k=0;
256 
257  C_j=0, C_l=0;
258  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
259 
260  C_j=1, C_l=0;
261  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
262 
263  C_j=0, C_l=1;
264  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
265 
266  C_j=1, C_l=1;
267  Kuu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
268  }
269 
270  for (unsigned int i=0; i<n_u_dofs; i++)
271  for (unsigned int j=0; j<n_v_dofs; j++)
272  {
273  // Tensor indices
274  unsigned int C_i, C_j, C_k, C_l;
275  C_i=0, C_k=1;
276 
277  C_j=0, C_l=0;
278  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
279 
280  C_j=1, C_l=0;
281  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
282 
283  C_j=0, C_l=1;
284  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
285 
286  C_j=1, C_l=1;
287  Kuv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
288  }
289 
290  for (unsigned int i=0; i<n_v_dofs; i++)
291  for (unsigned int j=0; j<n_u_dofs; j++)
292  {
293  // Tensor indices
294  unsigned int C_i, C_j, C_k, C_l;
295  C_i=1, C_k=0;
296 
297  C_j=0, C_l=0;
298  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
299 
300  C_j=1, C_l=0;
301  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
302 
303  C_j=0, C_l=1;
304  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
305 
306  C_j=1, C_l=1;
307  Kvu(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
308  }
309 
310  for (unsigned int i=0; i<n_v_dofs; i++)
311  for (unsigned int j=0; j<n_v_dofs; j++)
312  {
313  // Tensor indices
314  unsigned int C_i, C_j, C_k, C_l;
315  C_i=1, C_k=1;
316 
317  C_j=0, C_l=0;
318  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
319 
320  C_j=1, C_l=0;
321  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
322 
323  C_j=0, C_l=1;
324  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
325 
326  C_j=1, C_l=1;
327  Kvv(i,j) += JxW[qp]*(eval_elasticity_tensor(C_i, C_j, C_k, C_l) * dphi[i][qp](C_j)*dphi[j][qp](C_l));
328  }
329  }
330 
331  {
332  std::vector<boundary_id_type> bc_ids;
333  for (auto side : elem->side_index_range())
334  if (elem->neighbor_ptr(side) == nullptr)
335  {
336  mesh.get_boundary_info().boundary_ids (elem, side, bc_ids);
337 
338  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
339  const std::vector<Real> & JxW_face = fe_face->get_JxW();
340 
341  fe_face->reinit(elem, side);
342 
343  for (std::vector<boundary_id_type>::const_iterator b =
344  bc_ids.begin(); b != bc_ids.end(); ++b)
345  {
346  const boundary_id_type bc_id = *b;
347  for (unsigned int qp=0; qp<qface.n_points(); qp++)
348  {
349  // Add the loading
350  if (bc_id == 2)
351  for (unsigned int i=0; i<n_v_dofs; i++)
352  Fv(i) += JxW_face[qp] * (-1.) * phi_face[i][qp];
353 
354  // Add the constraint contributions
355  if (bc_id == 1)
356  {
357  for (unsigned int i=0; i<n_v_dofs; i++)
358  for (unsigned int j=0; j<n_lambda_dofs; j++)
359  Kv_lambda(i,j) += JxW_face[qp] * (-1.) * phi_face[i][qp];
360 
361  for (unsigned int i=0; i<n_lambda_dofs; i++)
362  for (unsigned int j=0; j<n_v_dofs; j++)
363  Klambda_v(i,j) += JxW_face[qp] * (-1.) * phi_face[j][qp];
364  }
365  }
366  }
367  }
368  }
369 
370  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
371 
372  matrix.add_matrix (Ke, dof_indices);
373  system.rhs->add_vector (Fe, dof_indices);
374  }
375 }
376 
378  unsigned int j,
379  unsigned int k,
380  unsigned int l)
381 {
382  // Define the Poisson ratio
383  const Real nu = 0.3;
384 
385  // Define the Lame constants (lambda_1 and lambda_2) based on Poisson ratio
386  const Real lambda_1 = nu / ((1. + nu) * (1. - 2.*nu));
387  const Real lambda_2 = 0.5 / (1 + nu);
388 
389  // Define the Kronecker delta functions that we need here
390  Real delta_ij = (i == j) ? 1. : 0.;
391  Real delta_il = (i == l) ? 1. : 0.;
392  Real delta_ik = (i == k) ? 1. : 0.;
393  Real delta_jl = (j == l) ? 1. : 0.;
394  Real delta_jk = (j == k) ? 1. : 0.;
395  Real delta_kl = (k == l) ? 1. : 0.;
396 
397  return lambda_1 * delta_ij * delta_kl + lambda_2 * (delta_ik * delta_jl + delta_il * delta_jk);
398 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
T command_line_next(std::string name, T default_value)
Use GetPot&#39;s search()/next() functions to get following arguments from the command line...
Definition: libmesh.C:1011
This is the EquationSystems class.
ConstFunction that simply returns 0.
Definition: zero_function.h:38
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
Definition: dof_map.h:2254
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1992
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:374
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:2144
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Order default_quadrature_order() const
Definition: fe_type.h:357
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
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 BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:159
virtual void solve() override
Assembles & solves the linear system A*x=b.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:74
unsigned int variable_number(std::string_view var) const
Definition: system.C:1557
Defines a dense subvector for use in finite element computations.
SolverPackage default_solver_package()
Definition: libmesh.C:1050
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
void reposition(const unsigned int ioff, const unsigned int n)
Changes the location of the subvector in the parent vector.
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
int8_t boundary_id_type
Definition: id_types.h:51
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:1489
Defines a dense submatrix for use in Finite Element-type computations.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
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:1305
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:2109
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MeshBase & get_mesh() const
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
Real eval_elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
void reposition(const unsigned int ioff, const unsigned int joff, const unsigned int new_m, const unsigned int new_n)
Changes the location of the submatrix in the parent matrix.
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.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
const DofMap & get_dof_map() const
Definition: system.h:2293
const SparseMatrix< Number > & get_system_matrix() const
void assemble_elasticity(EquationSystems &es, const std::string &system_name)
int main(int argc, char **argv)