libMesh
miscellaneous_ex11.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>Miscellaneous Example 11 - Using Loop Subdivision Shell Elements</h1>
21 // \author Roman Vetter
22 // \author Norbert Stoop
23 // \date 2014
24 //
25 // This example demonstrates how subdivision surface shell elements
26 // are used, and how boundary conditions can be applied to them. To
27 // keep it simple, we solve the static deflection of a clamped, square,
28 // linearly elastic Kirchhoff-Love plate subject to a uniform
29 // load distribution. Refer to Cirak et al., Int. J. Numer. Meth.
30 // Engng. 2000; 47: 2039-2072, for a detailed description of what's
31 // implemented. In fact, this example follows that paper very closely.
32 
33 
34 // C++ include files that we need
35 #include <iostream>
36 
37 // LibMesh includes
38 #include "libmesh/libmesh.h"
39 #include "libmesh/replicated_mesh.h"
40 #include "libmesh/mesh_refinement.h"
41 #include "libmesh/mesh_modification.h"
42 #include "libmesh/mesh_tools.h"
43 #include "libmesh/linear_implicit_system.h"
44 #include "libmesh/equation_systems.h"
45 #include "libmesh/fe.h"
46 #include "libmesh/quadrature.h"
47 #include "libmesh/node.h"
48 #include "libmesh/elem.h"
49 #include "libmesh/dof_map.h"
50 #include "libmesh/vector_value.h"
51 #include "libmesh/tensor_value.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/sparse_matrix.h"
57 #include "libmesh/numeric_vector.h"
58 #include "libmesh/vtk_io.h"
59 #include "libmesh/exodusII_io.h"
60 #include "libmesh/enum_solver_package.h"
61 #include "libmesh/parallel.h"
62 #include "libmesh/getpot.h"
63 
64 // These are the include files typically needed for subdivision elements.
65 #include "libmesh/face_tri3_subdivision.h"
66 #include "libmesh/mesh_subdivision_support.h"
67 
68 // Bring in everything from the libMesh namespace
69 using namespace libMesh;
70 
71 #if defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2
72 // Function prototype. This is the function that will assemble
73 // the stiffness matrix and the right-hand-side vector ready
74 // for solution.
76  const std::string & system_name);
77 #endif
78 
79 // Begin the main program.
80 int main (int argc, char ** argv)
81 {
82  // Initialize libMesh.
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  // Skip this 3D example if libMesh was compiled as 1D/2D-only.
90  libmesh_example_requires (3 == LIBMESH_DIM, "3D support");
91 
92  // Skip this example without --enable-node-valence
93 #ifndef LIBMESH_ENABLE_NODE_VALENCE
94  libmesh_example_requires (false, "--enable-node-valence");
95 #endif
96 
97  // Skip this example without --enable-amr; requires MeshRefinement
98 #ifndef LIBMESH_ENABLE_AMR
99  libmesh_example_requires(false, "--enable-amr");
100 #else
101 
102  // Skip this example without --enable-second; requires d2phi
103 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES
104  libmesh_example_requires(false, "--enable-second");
105 #elif LIBMESH_DIM > 2
106 
107  // Create a 2D mesh distributed across the default MPI communicator.
108  // Subdivision surfaces do not appear to work with DistributedMesh yet.
109  ReplicatedMesh mesh (init.comm(), 2);
110 
111  // Read the coarse square mesh.
112  mesh.read ("square_mesh.off");
113 
114  // Resize the square plate to edge length L.
115  const Real L = 100.;
117 
118  // Get the number of mesh refinements from the command line
119  const int n_refinements =
120  libMesh::command_line_next("-n_refinements", 3);
121 
122  // Quadrisect the mesh triangles a few times to obtain a
123  // finer mesh. Subdivision surface elements require the
124  // refinement data to be removed afterward.
125  MeshRefinement mesh_refinement (mesh);
126  mesh_refinement.uniformly_refine (n_refinements);
128 
129  // Write the mesh before the ghost elements are added.
130 #if defined(LIBMESH_HAVE_VTK)
131  VTKIO(mesh).write ("without_ghosts.pvtu");
132 #endif
133 #if defined(LIBMESH_HAVE_EXODUS_API)
134  ExodusII_IO(mesh).write ("without_ghosts.e");
135 #endif
136 
137  // Print information about the triangulated mesh to the screen.
138  mesh.print_info();
139 
140  // Turn the triangulated mesh into a subdivision mesh
141  // and add an additional row of "ghost" elements around
142  // it in order to complete the extended local support of
143  // the triangles at the boundaries. If the second
144  // argument is set to true, the outermost existing
145  // elements are converted into ghost elements, and the
146  // actual physical mesh is thus getting smaller.
148 
149  // Print information about the subdivision mesh to the screen.
150  mesh.print_info();
151 
152  // Write the mesh with the ghost elements added.
153  // Compare this to the original mesh to see the difference.
154 #if defined(LIBMESH_HAVE_VTK)
155  VTKIO(mesh).write ("with_ghosts.pvtu");
156 #endif
157 #if defined(LIBMESH_HAVE_EXODUS_API)
158  ExodusII_IO(mesh).write ("with_ghosts.e");
159 #endif
160 
161  // Create an equation systems object.
162  EquationSystems equation_systems (mesh);
163 
164  // Declare the system and its variables.
165  // Create a linear implicit system named "Shell".
166  LinearImplicitSystem & system = equation_systems.add_system<LinearImplicitSystem> ("Shell");
167 
168  // Add the three translational deformation variables
169  // "u", "v", "w" to "Shell". Since subdivision shell
170  // elements meet the C1-continuity requirement, no
171  // rotational or other auxiliary variables are needed.
172  // Loop Subdivision Elements are always interpolated
173  // by quartic box splines, hence the order must always
174  // be FOURTH.
175  system.add_variable ("u", FOURTH, SUBDIVISION);
176  system.add_variable ("v", FOURTH, SUBDIVISION);
177  system.add_variable ("w", FOURTH, SUBDIVISION);
178 
179  // Give the system a pointer to the matrix and rhs assembly
180  // function.
182 
183  // Use the parameters of the equation systems object to
184  // tell the shell system about the material properties, the
185  // shell thickness, and the external load.
186  const Real h = 1.;
187  const Real E = 1.e7;
188  const Real nu = 0.;
189  const Real q = 1.;
190  equation_systems.parameters.set<Real> ("thickness") = h;
191  equation_systems.parameters.set<Real> ("young's modulus") = E;
192  equation_systems.parameters.set<Real> ("poisson ratio") = nu;
193  equation_systems.parameters.set<Real> ("uniform load") = q;
194 
195  // Initialize the data structures for the equation system.
196  equation_systems.init();
197 
198  // Print information about the system to the screen.
199  equation_systems.print_info();
200 
201  // Solve the linear system.
202  system.solve();
203 
204  // After solving the system, write the solution to a VTK
205  // or ExodusII output file ready for import in, e.g.,
206  // Paraview.
207 #if defined(LIBMESH_HAVE_VTK)
208  VTKIO(mesh).write_equation_systems ("out.pvtu", equation_systems);
209 #endif
210 #if defined(LIBMESH_HAVE_EXODUS_API)
211  ExodusII_IO(mesh).write_equation_systems ("out.e", equation_systems);
212 #endif
213 
214  // Find the center node to measure the maximum deformation of the plate.
215  Node * center_node = 0;
216  Real nearest_dist_sq = mesh.point(0).norm_sq();
217  for (unsigned int nid=1; nid<mesh.n_nodes(); ++nid)
218  {
219  const Real dist_sq = mesh.point(nid).norm_sq();
220  if (dist_sq < nearest_dist_sq)
221  {
222  nearest_dist_sq = dist_sq;
223  center_node = mesh.node_ptr(nid);
224  }
225  }
226 
227  // Finally, we evaluate the z-displacement "w" at the center node.
228  const unsigned int w_var = system.variable_number ("w");
229  dof_id_type w_dof = center_node->dof_number (system.number(), w_var, 0);
230  Number w = 0;
231  if (w_dof >= system.get_dof_map().first_dof() &&
232  w_dof < system.get_dof_map().end_dof())
233  w = system.current_solution(w_dof);
234  system.comm().sum(w);
235 
236 
237  // The analytic solution for the maximum displacement of
238  // a clamped square plate in pure bending, from Taylor,
239  // Govindjee, Commun. Numer. Meth. Eng. 20, 757-765, 2004.
240  const Real D = E * h*h*h / (12*(1-nu*nu));
241  const Real w_analytic = 0.001265319 * L*L*L*L * q / D;
242 
243  // Print the finite element solution and the analytic
244  // prediction of the maximum displacement of the clamped
245  // square plate to the screen.
246  libMesh::out << "z-displacement of the center point: " << w << std::endl;
247  libMesh::out << "Analytic solution for pure bending: " << w_analytic << std::endl;
248 
249 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES, LIBMESH_DIM > 2
250 
251 #endif // #ifdef LIBMESH_ENABLE_AMR
252 
253  // All done.
254  return 0;
255 }
256 
257 #if defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2
258 
259 // We now define the matrix and rhs vector assembly function
260 // for the shell system. This function implements the
261 // linear Kirchhoff-Love theory for thin shells. At the
262 // end we also take into account the boundary conditions
263 // here, using the penalty method.
265  const std::string & libmesh_dbg_var(system_name))
266 {
267  // It is a good idea to make sure we are assembling
268  // the proper system.
269  libmesh_assert_equal_to (system_name, "Shell");
270 
271  // Get a constant reference to the mesh object.
272  const MeshBase & mesh = es.get_mesh();
273 
274  // Get a reference to the shell system object.
275  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Shell");
276 
277  // Get the shell parameters that we need during assembly.
278  const Real h = es.parameters.get<Real> ("thickness");
279  const Real E = es.parameters.get<Real> ("young's modulus");
280  const Real nu = es.parameters.get<Real> ("poisson ratio");
281  const Real q = es.parameters.get<Real> ("uniform load");
282 
283  // Compute the membrane stiffness K and the bending
284  // rigidity D from these parameters.
285  const Real K = E * h / (1-nu*nu);
286  const Real D = E * h*h*h / (12*(1-nu*nu));
287 
288  // Numeric ids corresponding to each variable in the system.
289  const unsigned int u_var = system.variable_number ("u");
290  const unsigned int v_var = system.variable_number ("v");
291  const unsigned int w_var = system.variable_number ("w");
292 
293  // Get the Finite Element type for "u". Note this will be
294  // the same as the type for "v" and "w".
295  FEType fe_type = system.variable_type (u_var);
296 
297  // Build a Finite Element object of the specified type.
298  std::unique_ptr<FEBase> fe (FEBase::build(2, fe_type));
299 
300  // A Gauss quadrature rule for numerical integration.
301  // For subdivision shell elements, a single Gauss point per
302  // element is sufficient, hence we use extraorder = 0.
303  const int extraorder = 0;
304  std::unique_ptr<QBase> qrule (fe_type.default_quadrature_rule (2, extraorder));
305 
306  // Tell the finite element object to use our quadrature rule.
307  fe->attach_quadrature_rule (qrule.get());
308 
309  // The element Jacobian * quadrature weight at each integration point.
310  const std::vector<Real> & JxW = fe->get_JxW();
311 
312  // The surface tangents in both directions at the quadrature points.
313  const std::vector<RealGradient> & dxyzdxi = fe->get_dxyzdxi();
314  const std::vector<RealGradient> & dxyzdeta = fe->get_dxyzdeta();
315 
316  // The second partial derivatives at the quadrature points.
317  const std::vector<RealGradient> & d2xyzdxi2 = fe->get_d2xyzdxi2();
318  const std::vector<RealGradient> & d2xyzdeta2 = fe->get_d2xyzdeta2();
319  const std::vector<RealGradient> & d2xyzdxideta = fe->get_d2xyzdxideta();
320 
321  // The element shape function and its derivatives evaluated at the
322  // quadrature points.
323  const std::vector<std::vector<Real>> & phi = fe->get_phi();
324  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
325  const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
326 
327  // A reference to the DofMap object for this system. The DofMap
328  // object handles the index translation from node and element numbers
329  // to degree of freedom numbers.
330  const DofMap & dof_map = system.get_dof_map();
331 
332  // The global system matrix
333  SparseMatrix<Number> & matrix = system.get_system_matrix();
334 
335  // Define data structures to contain the element stiffness matrix
336  // and right-hand-side vector contribution. Following
337  // basic finite element terminology we will denote these
338  // "Ke" and "Fe".
341 
343  Kuu(Ke), Kuv(Ke), Kuw(Ke),
344  Kvu(Ke), Kvv(Ke), Kvw(Ke),
345  Kwu(Ke), Kwv(Ke), Kww(Ke);
346 
348  Fu(Fe),
349  Fv(Fe),
350  Fw(Fe);
351 
352  // This vector will hold the degree of freedom indices for
353  // the element. These define where in the global system
354  // the element degrees of freedom get mapped.
355  std::vector<dof_id_type> dof_indices;
356  std::vector<dof_id_type> dof_indices_u;
357  std::vector<dof_id_type> dof_indices_v;
358  std::vector<dof_id_type> dof_indices_w;
359 
360  // Now we will loop over all the elements in the mesh. We will
361  // compute the element matrix and right-hand-side contribution.
362  for (const auto & elem : mesh.active_local_element_ptr_range())
363  {
364  // The ghost elements at the boundaries need to be excluded
365  // here, as they don't belong to the physical shell,
366  // but serve for a proper boundary treatment only.
367  libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION);
368  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *> (elem);
369  if (sd_elem->is_ghost())
370  continue;
371 
372  // Get the degree of freedom indices for the
373  // current element. These define where in the global
374  // matrix and right-hand-side this element will
375  // contribute to.
376  dof_map.dof_indices (elem, dof_indices);
377  dof_map.dof_indices (elem, dof_indices_u, u_var);
378  dof_map.dof_indices (elem, dof_indices_v, v_var);
379  dof_map.dof_indices (elem, dof_indices_w, w_var);
380 
381  const std::size_t n_dofs = dof_indices.size();
382  const std::size_t n_u_dofs = dof_indices_u.size();
383  const std::size_t n_v_dofs = dof_indices_v.size();
384  const std::size_t n_w_dofs = dof_indices_w.size();
385 
386  // Compute the element-specific data for the current
387  // element. This involves computing the location of the
388  // quadrature points and the shape functions
389  // (phi, dphi, d2phi) for the current element.
390  fe->reinit (elem);
391 
392  // Zero the element matrix and right-hand side before
393  // summing them. We use the resize member here because
394  // the number of degrees of freedom might have changed from
395  // the last element.
396  Ke.resize (n_dofs, n_dofs);
397  Fe.resize (n_dofs);
398 
399  // Reposition the submatrices... The idea is this:
400  //
401  // - - - -
402  // | Kuu Kuv Kuw | | Fu |
403  // Ke = | Kvu Kvv Kvw |; Fe = | Fv |
404  // | Kwu Kwv Kww | | Fw |
405  // - - - -
406  //
407  // The DenseSubMatrix.reposition () member takes the
408  // (row_offset, column_offset, row_size, column_size).
409  //
410  // Similarly, the DenseSubVector.reposition () member
411  // takes the (row_offset, row_size)
412  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
413  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
414  Kuw.reposition (u_var*n_u_dofs, w_var*n_u_dofs, n_u_dofs, n_w_dofs);
415 
416  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
417  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
418  Kvw.reposition (v_var*n_v_dofs, w_var*n_v_dofs, n_v_dofs, n_w_dofs);
419 
420  Kwu.reposition (w_var*n_w_dofs, u_var*n_w_dofs, n_w_dofs, n_u_dofs);
421  Kwv.reposition (w_var*n_w_dofs, v_var*n_w_dofs, n_w_dofs, n_v_dofs);
422  Kww.reposition (w_var*n_w_dofs, w_var*n_w_dofs, n_w_dofs, n_w_dofs);
423 
424  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
425  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
426  Fw.reposition (w_var*n_u_dofs, n_w_dofs);
427 
428  // Now we will build the element matrix and right-hand-side.
429  for (unsigned int qp=0; qp<qrule->n_points(); ++qp)
430  {
431  // First, we compute the external force resulting
432  // from a load q distributed uniformly across the plate.
433  // Since the load is supposed to be transverse to the plate,
434  // it affects the z-direction, i.e. the "w" variable.
435  for (unsigned int i=0; i<n_u_dofs; ++i)
436  Fw(i) += JxW[qp] * phi[i][qp] * q;
437 
438  // Next, we assemble the stiffness matrix. This is only valid
439  // for the linear theory, i.e., for small deformations, where
440  // reference and deformed surface metrics are indistinguishable.
441 
442  // Get the three surface basis vectors.
443  const RealVectorValue & a1 = dxyzdxi[qp];
444  const RealVectorValue & a2 = dxyzdeta[qp];
445  RealVectorValue a3 = a1.cross(a2);
446  const Real jac = a3.norm(); // the surface Jacobian
447  libmesh_assert_greater (jac, 0);
448  a3 /= jac; // the shell director a3 is normalized to unit length
449 
450  // Get the derivatives of the surface tangents.
451  const RealVectorValue & a11 = d2xyzdxi2[qp];
452  const RealVectorValue & a22 = d2xyzdeta2[qp];
453  const RealVectorValue & a12 = d2xyzdxideta[qp];
454 
455  // Compute the three covariant components of the first
456  // fundamental form of the surface.
457  const RealVectorValue a(a1*a1, a2*a2, a1*a2);
458 
459  // The elastic H matrix in Voigt's notation, computed from the
460  // covariant components of the first fundamental form rather
461  // than the contravariant components, exploiting that the
462  // contravariant first fundamental form is the inverse of the
463  // covariant first fundamental form (hence the determinant etc.).
464  RealTensorValue H;
465  H(0,0) = a(1) * a(1);
466  H(0,1) = H(1,0) = nu * a(1) * a(0) + (1-nu) * a(2) * a(2);
467  H(0,2) = H(2,0) = -a(1) * a(2);
468  H(1,1) = a(0) * a(0);
469  H(1,2) = H(2,1) = -a(0) * a(2);
470  H(2,2) = 0.5 * ((1-nu) * a(1) * a(0) + (1+nu) * a(2) * a(2));
471  const Real det = a(0) * a(1) - a(2) * a(2);
472  libmesh_assert_not_equal_to (det * det, 0);
473  H /= det * det;
474 
475  // Precompute come cross products for the bending part below.
476  const RealVectorValue a11xa2 = a11.cross(a2);
477  const RealVectorValue a22xa2 = a22.cross(a2);
478  const RealVectorValue a12xa2 = a12.cross(a2);
479  const RealVectorValue a1xa11 = a1.cross(a11);
480  const RealVectorValue a1xa22 = a1.cross(a22);
481  const RealVectorValue a1xa12 = a1.cross(a12);
482  const RealVectorValue a2xa3 = a2.cross(a3);
483  const RealVectorValue a3xa1 = a3.cross(a1);
484 
485  // Loop over all pairs of nodes I,J.
486  for (unsigned int i=0; i<n_u_dofs; ++i)
487  {
488  for (unsigned int j=0; j<n_u_dofs; ++j)
489  {
490  // The membrane strain matrices in Voigt's notation.
491  RealTensorValue MI, MJ;
492  for (unsigned int k=0; k<3; ++k)
493  {
494  MI(0,k) = dphi[i][qp](0) * a1(k);
495  MI(1,k) = dphi[i][qp](1) * a2(k);
496  MI(2,k) = dphi[i][qp](1) * a1(k)
497  + dphi[i][qp](0) * a2(k);
498 
499  MJ(0,k) = dphi[j][qp](0) * a1(k);
500  MJ(1,k) = dphi[j][qp](1) * a2(k);
501  MJ(2,k) = dphi[j][qp](1) * a1(k)
502  + dphi[j][qp](0) * a2(k);
503  }
504 
505  // The bending strain matrices in Voigt's notation.
506  RealTensorValue BI, BJ;
507  for (unsigned int k=0; k<3; ++k)
508  {
509  const Real term_ik = dphi[i][qp](0) * a2xa3(k)
510  + dphi[i][qp](1) * a3xa1(k);
511  BI(0,k) = -d2phi[i][qp](0,0) * a3(k)
512  +(dphi[i][qp](0) * a11xa2(k)
513  + dphi[i][qp](1) * a1xa11(k)
514  + (a3*a11) * term_ik) / jac;
515  BI(1,k) = -d2phi[i][qp](1,1) * a3(k)
516  +(dphi[i][qp](0) * a22xa2(k)
517  + dphi[i][qp](1) * a1xa22(k)
518  + (a3*a22) * term_ik) / jac;
519  BI(2,k) = 2 * (-d2phi[i][qp](0,1) * a3(k)
520  +(dphi[i][qp](0) * a12xa2(k)
521  + dphi[i][qp](1) * a1xa12(k)
522  + (a3*a12) * term_ik) / jac);
523 
524  const Real term_jk = dphi[j][qp](0) * a2xa3(k)
525  + dphi[j][qp](1) * a3xa1(k);
526  BJ(0,k) = -d2phi[j][qp](0,0) * a3(k)
527  +(dphi[j][qp](0) * a11xa2(k)
528  + dphi[j][qp](1) * a1xa11(k)
529  + (a3*a11) * term_jk) / jac;
530  BJ(1,k) = -d2phi[j][qp](1,1) * a3(k)
531  +(dphi[j][qp](0) * a22xa2(k)
532  + dphi[j][qp](1) * a1xa22(k)
533  + (a3*a22) * term_jk) / jac;
534  BJ(2,k) = 2 * (-d2phi[j][qp](0,1) * a3(k)
535  +(dphi[j][qp](0) * a12xa2(k)
536  + dphi[j][qp](1) * a1xa12(k)
537  + (a3*a12) * term_jk) / jac);
538  }
539 
540  // The total stiffness matrix coupling the nodes
541  // I and J is a sum of membrane and bending
542  // contributions according to the following formula.
543  const RealTensorValue KIJ = JxW[qp] * K * MI.transpose() * H * MJ
544  + JxW[qp] * D * BI.transpose() * H * BJ;
545 
546  // Insert the components of the coupling stiffness
547  // matrix KIJ into the corresponding directional
548  // submatrices.
549  Kuu(i,j) += KIJ(0,0);
550  Kuv(i,j) += KIJ(0,1);
551  Kuw(i,j) += KIJ(0,2);
552 
553  Kvu(i,j) += KIJ(1,0);
554  Kvv(i,j) += KIJ(1,1);
555  Kvw(i,j) += KIJ(1,2);
556 
557  Kwu(i,j) += KIJ(2,0);
558  Kwv(i,j) += KIJ(2,1);
559  Kww(i,j) += KIJ(2,2);
560  }
561  }
562 
563  } // end of the quadrature point qp-loop
564 
565  // The element matrix and right-hand-side are now built
566  // for this element. Add them to the global matrix and
567  // right-hand-side vector. The NumericMatrix::add_matrix()
568  // and NumericVector::add_vector() members do this for us.
569  matrix.add_matrix (Ke, dof_indices);
570  system.rhs->add_vector (Fe, dof_indices);
571  } // end of non-ghost element loop
572 
573  // Next, we apply the boundary conditions. In this case,
574  // all boundaries are clamped by the penalty method, using
575  // the special "ghost" nodes along the boundaries. Note
576  // that there are better ways to implement boundary conditions
577  // for subdivision shells. We use the simplest way here,
578  // which is known to be overly restrictive and will lead to
579  // a slightly too small deformation of the plate.
580  for (const auto & elem : mesh.active_local_element_ptr_range())
581  {
582  // For the boundary conditions, we only need to loop over
583  // the ghost elements.
584  libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION);
585  const Tri3Subdivision * gh_elem = static_cast<const Tri3Subdivision *> (elem);
586  if (!gh_elem->is_ghost())
587  continue;
588 
589  // Find the side which is part of the physical plate boundary,
590  // that is, the boundary of the original mesh without ghosts.
591  for (auto s : elem->side_index_range())
592  {
593  const Tri3Subdivision * nb_elem = static_cast<const Tri3Subdivision *> (elem->neighbor_ptr(s));
594  if (nb_elem == nullptr || nb_elem->is_ghost())
595  continue;
596 
597  /*
598  * Determine the four nodes involved in the boundary
599  * condition treatment of this side. The MeshTools::Subdiv
600  * namespace provides lookup tables next and prev
601  * for an efficient determination of the next and previous
602  * nodes of an element, respectively.
603  *
604  * n4
605  * / \
606  * / gh \
607  * n2 ---- n3
608  * \ nb /
609  * \ /
610  * n1
611  */
612  const Node * nodes [4]; // n1, n2, n3, n4
613  nodes[1] = gh_elem->node_ptr(s); // n2
614  nodes[2] = gh_elem->node_ptr(MeshTools::Subdivision::next[s]); // n3
615  nodes[3] = gh_elem->node_ptr(MeshTools::Subdivision::prev[s]); // n4
616 
617  // The node in the interior of the domain, n1, is the
618  // hardest to find. Walk along the edges of element nb until
619  // we have identified it.
620  unsigned int n_int = 0;
621  nodes[0] = nb_elem->node_ptr(0);
622  while (nodes[0]->id() == nodes[1]->id() || nodes[0]->id() == nodes[2]->id())
623  nodes[0] = nb_elem->node_ptr(++n_int);
624 
625  // The penalty value. \f$ \frac{1}{\epsilon} \f$
626  const Real penalty = 1.e10;
627 
628  // With this simple method, clamped boundary conditions are
629  // obtained by penalizing the displacements of all four nodes.
630  // This ensures that the displacement field vanishes on the
631  // boundary side s.
632  for (unsigned int n=0; n<4; ++n)
633  {
634  const dof_id_type u_dof = nodes[n]->dof_number (system.number(), u_var, 0);
635  const dof_id_type v_dof = nodes[n]->dof_number (system.number(), v_var, 0);
636  const dof_id_type w_dof = nodes[n]->dof_number (system.number(), w_var, 0);
637  matrix.add (u_dof, u_dof, penalty);
638  matrix.add (v_dof, v_dof, penalty);
639  matrix.add (w_dof, w_dof, penalty);
640  }
641  }
642  } // end of ghost element loop
643 }
644 
645 #endif // defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) && LIBMESH_DIM > 2
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
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:1078
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map_base.h:191
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.
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:1032
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
A Node is like a Point, but with more information.
Definition: node.h:52
auto norm() const -> decltype(std::norm(T()))
Definition: type_vector.h:907
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
Definition: dof_map.C:2164
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
Scales the mesh.
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 ...
TypeTensor< T > transpose() const
Definition: type_tensor.h:1050
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:396
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]...
void sum(T &r) const
static const unsigned int prev[3]
A lookup table for the decrement modulo 3 operation, for iterating through the three nodes per elemen...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
virtual void write(const std::string &) override
Output the mesh without solutions to a .pvtu file.
void assemble_shell(EquationSystems &es, const std::string &system_name)
NumericVector< Number > * rhs
The system matrix.
const Parallel::Communicator & comm() const
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.
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
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:75
Number current_solution(const dof_id_type global_dof_number) const
Definition: system.C:165
unsigned int variable_number(std::string_view var) const
Definition: system.C:1609
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
Defines a dense subvector for use in finite element computations.
SolverPackage default_solver_package()
Definition: libmesh.C:1117
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:34
void reposition(const unsigned int ioff, const unsigned int n)
Changes the location of the subvector in the parent vector.
unsigned int number() const
Definition: system.h:2350
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.
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_vector.h:926
const T & get(std::string_view) const
Definition: parameters.h:426
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
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.
The Tri3Subdivision element is a three-noded subdivision surface shell element used in mechanics calc...
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
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:2161
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:884
static const unsigned int next[3]
A lookup table for the increment modulo 3 operation, for iterating through the three nodes per elemen...
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2508
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2507
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
Definition: exodusII_io.C:2180
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
const MeshBase & get_mesh() const
void prepare_subdivision_mesh(MeshBase &mesh, bool ghosted=false)
Prepares the mesh for use with subdivision elements.
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
Parameters parameters
Data structure holding arbitrary parameters.
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 const Point & point(const dof_id_type i) const =0
virtual void init()
Initialize all the systems.
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map_base.h:185
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.
virtual const Node * node_ptr(const dof_id_type i) const =0
const DofMap & get_dof_map() const
Definition: system.h:2374
int main(int argc, char **argv)
const SparseMatrix< Number > & get_system_matrix() const
void flatten(MeshBase &mesh)
Removes all the refinement tree structure of Mesh, leaving only the highest-level (most-refined) elem...
virtual dof_id_type n_nodes() const =0
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
uint8_t dof_id_type
Definition: id_types.h:67
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.