libMesh
Functions
miscellaneous_ex11.C File Reference

Go to the source code of this file.

Functions

void assemble_shell (EquationSystems &es, const std::string &system_name)
 
int main (int argc, char **argv)
 
void assemble_shell (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_shell() [1/2]

void assemble_shell ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 259 of file miscellaneous_ex11.C.

261 {
262  // It is a good idea to make sure we are assembling
263  // the proper system.
264  libmesh_assert_equal_to (system_name, "Shell");
265 
266  // Get a constant reference to the mesh object.
267  const MeshBase & mesh = es.get_mesh();
268 
269  // Get a reference to the shell system object.
270  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Shell");
271 
272  // Get the shell parameters that we need during assembly.
273  const Real h = es.parameters.get<Real> ("thickness");
274  const Real E = es.parameters.get<Real> ("young's modulus");
275  const Real nu = es.parameters.get<Real> ("poisson ratio");
276  const Real q = es.parameters.get<Real> ("uniform load");
277 
278  // Compute the membrane stiffness K and the bending
279  // rigidity D from these parameters.
280  const Real K = E * h / (1-nu*nu);
281  const Real D = E * h*h*h / (12*(1-nu*nu));
282 
283  // Numeric ids corresponding to each variable in the system.
284  const unsigned int u_var = system.variable_number ("u");
285  const unsigned int v_var = system.variable_number ("v");
286  const unsigned int w_var = system.variable_number ("w");
287 
288  // Get the Finite Element type for "u". Note this will be
289  // the same as the type for "v" and "w".
290  FEType fe_type = system.variable_type (u_var);
291 
292  // Build a Finite Element object of the specified type.
293  std::unique_ptr<FEBase> fe (FEBase::build(2, fe_type));
294 
295  // A Gauss quadrature rule for numerical integration.
296  // For subdivision shell elements, a single Gauss point per
297  // element is sufficient, hence we use extraorder = 0.
298  const int extraorder = 0;
299  std::unique_ptr<QBase> qrule (fe_type.default_quadrature_rule (2, extraorder));
300 
301  // Tell the finite element object to use our quadrature rule.
302  fe->attach_quadrature_rule (qrule.get());
303 
304  // The element Jacobian * quadrature weight at each integration point.
305  const std::vector<Real> & JxW = fe->get_JxW();
306 
307  // The surface tangents in both directions at the quadrature points.
308  const std::vector<RealGradient> & dxyzdxi = fe->get_dxyzdxi();
309  const std::vector<RealGradient> & dxyzdeta = fe->get_dxyzdeta();
310 
311  // The second partial derivatives at the quadrature points.
312  const std::vector<RealGradient> & d2xyzdxi2 = fe->get_d2xyzdxi2();
313  const std::vector<RealGradient> & d2xyzdeta2 = fe->get_d2xyzdeta2();
314  const std::vector<RealGradient> & d2xyzdxideta = fe->get_d2xyzdxideta();
315 
316  // The element shape function and its derivatives evaluated at the
317  // quadrature points.
318  const std::vector<std::vector<Real>> & phi = fe->get_phi();
319  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
320  const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
321 
322  // A reference to the DofMap object for this system. The DofMap
323  // object handles the index translation from node and element numbers
324  // to degree of freedom numbers.
325  const DofMap & dof_map = system.get_dof_map();
326 
327  // Define data structures to contain the element stiffness matrix
328  // and right-hand-side vector contribution. Following
329  // basic finite element terminology we will denote these
330  // "Ke" and "Fe".
333 
335  Kuu(Ke), Kuv(Ke), Kuw(Ke),
336  Kvu(Ke), Kvv(Ke), Kvw(Ke),
337  Kwu(Ke), Kwv(Ke), Kww(Ke);
338 
340  Fu(Fe),
341  Fv(Fe),
342  Fw(Fe);
343 
344  // This vector will hold the degree of freedom indices for
345  // the element. These define where in the global system
346  // the element degrees of freedom get mapped.
347  std::vector<dof_id_type> dof_indices;
348  std::vector<dof_id_type> dof_indices_u;
349  std::vector<dof_id_type> dof_indices_v;
350  std::vector<dof_id_type> dof_indices_w;
351 
352  // Now we will loop over all the elements in the mesh. We will
353  // compute the element matrix and right-hand-side contribution.
354  for (const auto & elem : mesh.active_local_element_ptr_range())
355  {
356  // The ghost elements at the boundaries need to be excluded
357  // here, as they don't belong to the physical shell,
358  // but serve for a proper boundary treatment only.
359  libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION);
360  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *> (elem);
361  if (sd_elem->is_ghost())
362  continue;
363 
364  // Get the degree of freedom indices for the
365  // current element. These define where in the global
366  // matrix and right-hand-side this element will
367  // contribute to.
368  dof_map.dof_indices (elem, dof_indices);
369  dof_map.dof_indices (elem, dof_indices_u, u_var);
370  dof_map.dof_indices (elem, dof_indices_v, v_var);
371  dof_map.dof_indices (elem, dof_indices_w, w_var);
372 
373  const std::size_t n_dofs = dof_indices.size();
374  const std::size_t n_u_dofs = dof_indices_u.size();
375  const std::size_t n_v_dofs = dof_indices_v.size();
376  const std::size_t n_w_dofs = dof_indices_w.size();
377 
378  // Compute the element-specific data for the current
379  // element. This involves computing the location of the
380  // quadrature points and the shape functions
381  // (phi, dphi, d2phi) for the current element.
382  fe->reinit (elem);
383 
384  // Zero the element matrix and right-hand side before
385  // summing them. We use the resize member here because
386  // the number of degrees of freedom might have changed from
387  // the last element.
388  Ke.resize (n_dofs, n_dofs);
389  Fe.resize (n_dofs);
390 
391  // Reposition the submatrices... The idea is this:
392  //
393  // - - - -
394  // | Kuu Kuv Kuw | | Fu |
395  // Ke = | Kvu Kvv Kvw |; Fe = | Fv |
396  // | Kwu Kwv Kww | | Fw |
397  // - - - -
398  //
399  // The DenseSubMatrix.reposition () member takes the
400  // (row_offset, column_offset, row_size, column_size).
401  //
402  // Similarly, the DenseSubVector.reposition () member
403  // takes the (row_offset, row_size)
404  Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
405  Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
406  Kuw.reposition (u_var*n_u_dofs, w_var*n_u_dofs, n_u_dofs, n_w_dofs);
407 
408  Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
409  Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
410  Kvw.reposition (v_var*n_v_dofs, w_var*n_v_dofs, n_v_dofs, n_w_dofs);
411 
412  Kwu.reposition (w_var*n_w_dofs, u_var*n_w_dofs, n_w_dofs, n_u_dofs);
413  Kwv.reposition (w_var*n_w_dofs, v_var*n_w_dofs, n_w_dofs, n_v_dofs);
414  Kww.reposition (w_var*n_w_dofs, w_var*n_w_dofs, n_w_dofs, n_w_dofs);
415 
416  Fu.reposition (u_var*n_u_dofs, n_u_dofs);
417  Fv.reposition (v_var*n_u_dofs, n_v_dofs);
418  Fw.reposition (w_var*n_u_dofs, n_w_dofs);
419 
420  // Now we will build the element matrix and right-hand-side.
421  for (unsigned int qp=0; qp<qrule->n_points(); ++qp)
422  {
423  // First, we compute the external force resulting
424  // from a load q distributed uniformly across the plate.
425  // Since the load is supposed to be transverse to the plate,
426  // it affects the z-direction, i.e. the "w" variable.
427  for (unsigned int i=0; i<n_u_dofs; ++i)
428  Fw(i) += JxW[qp] * phi[i][qp] * q;
429 
430  // Next, we assemble the stiffness matrix. This is only valid
431  // for the linear theory, i.e., for small deformations, where
432  // reference and deformed surface metrics are indistinguishable.
433 
434  // Get the three surface basis vectors.
435  const RealVectorValue & a1 = dxyzdxi[qp];
436  const RealVectorValue & a2 = dxyzdeta[qp];
437  RealVectorValue a3 = a1.cross(a2);
438  const Real jac = a3.norm(); // the surface Jacobian
439  libmesh_assert_greater (jac, 0);
440  a3 /= jac; // the shell director a3 is normalized to unit length
441 
442  // Get the derivatives of the surface tangents.
443  const RealVectorValue & a11 = d2xyzdxi2[qp];
444  const RealVectorValue & a22 = d2xyzdeta2[qp];
445  const RealVectorValue & a12 = d2xyzdxideta[qp];
446 
447  // Compute the three covariant components of the first
448  // fundamental form of the surface.
449  const RealVectorValue a(a1*a1, a2*a2, a1*a2);
450 
451  // The elastic H matrix in Voigt's notation, computed from the
452  // covariant components of the first fundamental form rather
453  // than the contravariant components, exploiting that the
454  // contravariant first fundamental form is the inverse of the
455  // covariant first fundamental form (hence the determinant etc.).
456  RealTensorValue H;
457  H(0,0) = a(1) * a(1);
458  H(0,1) = H(1,0) = nu * a(1) * a(0) + (1-nu) * a(2) * a(2);
459  H(0,2) = H(2,0) = -a(1) * a(2);
460  H(1,1) = a(0) * a(0);
461  H(1,2) = H(2,1) = -a(0) * a(2);
462  H(2,2) = 0.5 * ((1-nu) * a(1) * a(0) + (1+nu) * a(2) * a(2));
463  const Real det = a(0) * a(1) - a(2) * a(2);
464  libmesh_assert_not_equal_to (det * det, 0);
465  H /= det * det;
466 
467  // Precompute come cross products for the bending part below.
468  const RealVectorValue a11xa2 = a11.cross(a2);
469  const RealVectorValue a22xa2 = a22.cross(a2);
470  const RealVectorValue a12xa2 = a12.cross(a2);
471  const RealVectorValue a1xa11 = a1.cross(a11);
472  const RealVectorValue a1xa22 = a1.cross(a22);
473  const RealVectorValue a1xa12 = a1.cross(a12);
474  const RealVectorValue a2xa3 = a2.cross(a3);
475  const RealVectorValue a3xa1 = a3.cross(a1);
476 
477  // Loop over all pairs of nodes I,J.
478  for (unsigned int i=0; i<n_u_dofs; ++i)
479  {
480  for (unsigned int j=0; j<n_u_dofs; ++j)
481  {
482  // The membrane strain matrices in Voigt's notation.
483  RealTensorValue MI, MJ;
484  for (unsigned int k=0; k<3; ++k)
485  {
486  MI(0,k) = dphi[i][qp](0) * a1(k);
487  MI(1,k) = dphi[i][qp](1) * a2(k);
488  MI(2,k) = dphi[i][qp](1) * a1(k)
489  + dphi[i][qp](0) * a2(k);
490 
491  MJ(0,k) = dphi[j][qp](0) * a1(k);
492  MJ(1,k) = dphi[j][qp](1) * a2(k);
493  MJ(2,k) = dphi[j][qp](1) * a1(k)
494  + dphi[j][qp](0) * a2(k);
495  }
496 
497  // The bending strain matrices in Voigt's notation.
498  RealTensorValue BI, BJ;
499  for (unsigned int k=0; k<3; ++k)
500  {
501  const Real term_ik = dphi[i][qp](0) * a2xa3(k)
502  + dphi[i][qp](1) * a3xa1(k);
503  BI(0,k) = -d2phi[i][qp](0,0) * a3(k)
504  +(dphi[i][qp](0) * a11xa2(k)
505  + dphi[i][qp](1) * a1xa11(k)
506  + (a3*a11) * term_ik) / jac;
507  BI(1,k) = -d2phi[i][qp](1,1) * a3(k)
508  +(dphi[i][qp](0) * a22xa2(k)
509  + dphi[i][qp](1) * a1xa22(k)
510  + (a3*a22) * term_ik) / jac;
511  BI(2,k) = 2 * (-d2phi[i][qp](0,1) * a3(k)
512  +(dphi[i][qp](0) * a12xa2(k)
513  + dphi[i][qp](1) * a1xa12(k)
514  + (a3*a12) * term_ik) / jac);
515 
516  const Real term_jk = dphi[j][qp](0) * a2xa3(k)
517  + dphi[j][qp](1) * a3xa1(k);
518  BJ(0,k) = -d2phi[j][qp](0,0) * a3(k)
519  +(dphi[j][qp](0) * a11xa2(k)
520  + dphi[j][qp](1) * a1xa11(k)
521  + (a3*a11) * term_jk) / jac;
522  BJ(1,k) = -d2phi[j][qp](1,1) * a3(k)
523  +(dphi[j][qp](0) * a22xa2(k)
524  + dphi[j][qp](1) * a1xa22(k)
525  + (a3*a22) * term_jk) / jac;
526  BJ(2,k) = 2 * (-d2phi[j][qp](0,1) * a3(k)
527  +(dphi[j][qp](0) * a12xa2(k)
528  + dphi[j][qp](1) * a1xa12(k)
529  + (a3*a12) * term_jk) / jac);
530  }
531 
532  // The total stiffness matrix coupling the nodes
533  // I and J is a sum of membrane and bending
534  // contributions according to the following formula.
535  const RealTensorValue KIJ = JxW[qp] * K * MI.transpose() * H * MJ
536  + JxW[qp] * D * BI.transpose() * H * BJ;
537 
538  // Insert the components of the coupling stiffness
539  // matrix KIJ into the corresponding directional
540  // submatrices.
541  Kuu(i,j) += KIJ(0,0);
542  Kuv(i,j) += KIJ(0,1);
543  Kuw(i,j) += KIJ(0,2);
544 
545  Kvu(i,j) += KIJ(1,0);
546  Kvv(i,j) += KIJ(1,1);
547  Kvw(i,j) += KIJ(1,2);
548 
549  Kwu(i,j) += KIJ(2,0);
550  Kwv(i,j) += KIJ(2,1);
551  Kww(i,j) += KIJ(2,2);
552  }
553  }
554 
555  } // end of the quadrature point qp-loop
556 
557  // The element matrix and right-hand-side are now built
558  // for this element. Add them to the global matrix and
559  // right-hand-side vector. The NumericMatrix::add_matrix()
560  // and NumericVector::add_vector() members do this for us.
561  system.matrix->add_matrix (Ke, dof_indices);
562  system.rhs->add_vector (Fe, dof_indices);
563  } // end of non-ghost element loop
564 
565  // Next, we apply the boundary conditions. In this case,
566  // all boundaries are clamped by the penalty method, using
567  // the special "ghost" nodes along the boundaries. Note
568  // that there are better ways to implement boundary conditions
569  // for subdivision shells. We use the simplest way here,
570  // which is known to be overly restrictive and will lead to
571  // a slightly too small deformation of the plate.
572  for (const auto & elem : mesh.active_local_element_ptr_range())
573  {
574  // For the boundary conditions, we only need to loop over
575  // the ghost elements.
576  libmesh_assert_equal_to (elem->type(), TRI3SUBDIVISION);
577  const Tri3Subdivision * gh_elem = static_cast<const Tri3Subdivision *> (elem);
578  if (!gh_elem->is_ghost())
579  continue;
580 
581  // Find the side which is part of the physical plate boundary,
582  // that is, the boundary of the original mesh without ghosts.
583  for (auto s : elem->side_index_range())
584  {
585  const Tri3Subdivision * nb_elem = static_cast<const Tri3Subdivision *> (elem->neighbor_ptr(s));
586  if (nb_elem == nullptr || nb_elem->is_ghost())
587  continue;
588 
589  /*
590  * Determine the four nodes involved in the boundary
591  * condition treatment of this side. The MeshTools::Subdiv
592  * namespace provides lookup tables next and prev
593  * for an efficient determination of the next and previous
594  * nodes of an element, respectively.
595  *
596  * n4
597  * / \
598  * / gh \
599  * n2 ---- n3
600  * \ nb /
601  * \ /
602  * n1
603  */
604  const Node * nodes [4]; // n1, n2, n3, n4
605  nodes[1] = gh_elem->node_ptr(s); // n2
606  nodes[2] = gh_elem->node_ptr(MeshTools::Subdivision::next[s]); // n3
607  nodes[3] = gh_elem->node_ptr(MeshTools::Subdivision::prev[s]); // n4
608 
609  // The node in the interior of the domain, n1, is the
610  // hardest to find. Walk along the edges of element nb until
611  // we have identified it.
612  unsigned int n_int = 0;
613  nodes[0] = nb_elem->node_ptr(0);
614  while (nodes[0]->id() == nodes[1]->id() || nodes[0]->id() == nodes[2]->id())
615  nodes[0] = nb_elem->node_ptr(++n_int);
616 
617  // The penalty value. \f$ \frac{1}{\epsilon} \f$
618  const Real penalty = 1.e10;
619 
620  // With this simple method, clamped boundary conditions are
621  // obtained by penalizing the displacements of all four nodes.
622  // This ensures that the displacement field vanishes on the
623  // boundary side s.
624  for (unsigned int n=0; n<4; ++n)
625  {
626  const dof_id_type u_dof = nodes[n]->dof_number (system.number(), u_var, 0);
627  const dof_id_type v_dof = nodes[n]->dof_number (system.number(), v_var, 0);
628  const dof_id_type w_dof = nodes[n]->dof_number (system.number(), w_var, 0);
629  system.matrix->add (u_dof, u_dof, penalty);
630  system.matrix->add (v_dof, v_dof, penalty);
631  system.matrix->add (w_dof, w_dof, penalty);
632  }
633  }
634  } // end of ghost element loop
635 }

References libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::TypeVector< T >::cross(), libMesh::FEType::default_quadrature_rule(), libMesh::DofMap::dof_indices(), libMesh::DofObject::dof_number(), libMesh::Parameters::get(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::Tri3Subdivision::is_ghost(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshTools::Subdivision::next, libMesh::Elem::node_ptr(), libMesh::TypeVector< T >::norm(), libMesh::System::number(), libMesh::EquationSystems::parameters, libMesh::MeshTools::Subdivision::prev, libMesh::Real, libMesh::DenseSubVector< T >::reposition(), libMesh::DenseSubMatrix< T >::reposition(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::TypeTensor< T >::transpose(), libMesh::TRI3SUBDIVISION, libMesh::System::variable_number(), and libMesh::System::variable_type().

◆ assemble_shell() [2/2]

void assemble_shell ( EquationSystems es,
const std::string &  system_name 
)

Definition at line 348 of file miscellaneous_ex12.C.

350 {
351  // This example requires Eigen to actually work, but we should still
352  // let it compile and throw a runtime error if you don't.
353 
354  // The same holds for second derivatives,
355  // since they are class-members only depending on the config.
356 #if defined(LIBMESH_HAVE_EIGEN) && defined(LIBMESH_ENABLE_SECOND_DERIVATIVES)
357  // It is a good idea to make sure we are assembling
358  // the proper system.
359  libmesh_assert_equal_to (system_name, "Shell");
360 
361  // Get a constant reference to the mesh object.
362  const MeshBase & mesh = es.get_mesh();
363  const unsigned int dim = mesh.mesh_dimension();
364 
365  // Get a reference to the shell system object.
366  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> (system_name);
367 
368  // Get the shell parameters that we need during assembly.
369  const Real h = es.parameters.get<Real> ("thickness");
370  const Real E = es.parameters.get<Real> ("young's modulus");
371  const Real nu = es.parameters.get<Real> ("poisson ratio");
372  const Real q = es.parameters.get<Real> ("point load");
373  const bool distributed_load = es.parameters.get<bool> ("distributed load");
374 
375  // The membrane elastic matrix.
376  MyMatrix3d Hm;
377  Hm <<
378  1., nu, 0.,
379  nu, 1., 0.,
380  0., 0., 0.5 * (1-nu);
381  Hm *= h * E/(1-nu*nu);
382 
383  // The bending elastic matrix.
384  MyMatrix3d Hf;
385  Hf <<
386  1., nu, 0.,
387  nu, 1., 0.,
388  0., 0., 0.5 * (1-nu);
389  Hf *= h*h*h/12 * E/(1-nu*nu);
390 
391  // The shear elastic matrices.
392  MyMatrix2d Hc0 = MyMatrix2d::Identity();
393  Hc0 *= h * 5./6*E/(2*(1+nu));
394 
395  MyMatrix2d Hc1 = MyMatrix2d::Identity();
396  Hc1 *= h*h*h/12 * 5./6*E/(2*(1+nu));
397 
398  // Get the Finite Element type, this will be
399  // the same for all variables.
400  FEType fe_type = system.variable_type (0);
401 
402  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
403  QGauss qrule (dim, fe_type.default_quadrature_order());
404  fe->attach_quadrature_rule (&qrule);
405 
406  // The element Jacobian * quadrature weight at each integration point.
407  const std::vector<Real> & JxW = fe->get_JxW();
408 
409  // The element shape function and its derivatives evaluated at the
410  // quadrature points.
411  const std::vector<RealGradient> & dxyzdxi = fe->get_dxyzdxi();
412  const std::vector<RealGradient> & dxyzdeta = fe->get_dxyzdeta();
413 
414  const std::vector<RealGradient> & d2xyzdxi2 = fe->get_d2xyzdxi2();
415  const std::vector<RealGradient> & d2xyzdeta2 = fe->get_d2xyzdeta2();
416  const std::vector<RealGradient> & d2xyzdxideta = fe->get_d2xyzdxideta();
417  const std::vector<std::vector<Real>> & dphidxi = fe->get_dphidxi();
418  const std::vector<std::vector<Real>> & dphideta = fe->get_dphideta();
419  const std::vector<std::vector<Real>> & phi = fe->get_phi();
420 
421  // A reference to the DofMap object for this system. The DofMap
422  // object handles the index translation from node and element numbers
423  // to degree of freedom numbers.
424  const DofMap & dof_map = system.get_dof_map();
425 
426  // Define data structures to contain the element stiffness matrix.
428  DenseSubMatrix<Number> Ke_var[6][6] =
429  {
442  };
443 
444  // Define data structures to contain the element rhs vector.
446  DenseSubVector<Number> Fe_w(Fe);
447 
448  std::vector<dof_id_type> dof_indices;
449  std::vector<std::vector<dof_id_type>> dof_indices_var(6);
450 
451  // Now we will loop over all the elements in the mesh. We will
452  // compute the element matrix and right-hand-side contribution.
453  for (const auto & elem : mesh.active_local_element_ptr_range())
454  {
455  dof_map.dof_indices (elem, dof_indices);
456  for (unsigned int var=0; var<6; var++)
457  dof_map.dof_indices (elem, dof_indices_var[var], var);
458 
459  const unsigned int n_dofs = dof_indices.size();
460  const unsigned int n_var_dofs = dof_indices_var[0].size();
461 
462  // First compute element data at the nodes
463  std::vector<Point> nodes;
464  for (auto i : elem->node_index_range())
465  nodes.push_back(elem->reference_elem()->node_ref(i));
466  fe->reinit (elem, &nodes);
467 
468  // Convenient notation for the element node positions
469  MyVector3d X1(elem->node_ref(0)(0), elem->node_ref(0)(1), elem->node_ref(0)(2));
470  MyVector3d X2(elem->node_ref(1)(0), elem->node_ref(1)(1), elem->node_ref(1)(2));
471  MyVector3d X3(elem->node_ref(2)(0), elem->node_ref(2)(1), elem->node_ref(2)(2));
472  MyVector3d X4(elem->node_ref(3)(0), elem->node_ref(3)(1), elem->node_ref(3)(2));
473 
474  //Store covariant basis and local orthonormal basis at the nodes
475  std::vector<MyMatrix3d> F0node;
476  std::vector<MyMatrix3d> Qnode;
477  for (auto i : elem->node_index_range())
478  {
479  MyVector3d a1;
480  a1 << dxyzdxi[i](0), dxyzdxi[i](1), dxyzdxi[i](2);
481  MyVector3d a2;
482  a2 << dxyzdeta[i](0), dxyzdeta[i](1), dxyzdeta[i](2);
483  MyVector3d n;
484  n = a1.cross(a2);
485  n /= n.norm();
486  MyMatrix3d F0;
487  F0 <<
488  a1(0), a2(0), n(0),
489  a1(1), a2(1), n(1),
490  a1(2), a2(2), n(2);
491  F0node.push_back(F0);
492 
493  Real nx = n(0);
494  Real ny = n(1);
495  Real C = n(2);
496  if (std::abs(1.+C)<1e-6)
497  {
498  MyMatrix3d Q;
499  Q <<
500  1, 0, 0,
501  0, -1, 0,
502  0, 0, -1;
503  Qnode.push_back(Q);
504  }
505  else
506  {
507  MyMatrix3d Q;
508  Q <<
509  C+1./(1+C)*ny*ny, -1./(1+C)*nx*ny, nx,
510  -1./(1+C)*nx*ny, C+1./(1+C)*nx*nx, ny,
511  -nx, -ny, C;
512  Qnode.push_back(Q);
513  }
514  }
515 
516  Ke.resize (n_dofs, n_dofs);
517  for (unsigned int var_i=0; var_i<6; var_i++)
518  for (unsigned int var_j=0; var_j<6; var_j++)
519  Ke_var[var_i][var_j].reposition (var_i*n_var_dofs, var_j*n_var_dofs, n_var_dofs, n_var_dofs);
520 
521  Fe.resize(n_dofs);
522  Fe_w.reposition(2*n_var_dofs,n_var_dofs);
523 
524  // Reinit element data at the regular Gauss quadrature points
525  fe->reinit (elem);
526 
527  // Now we will build the element matrix and right-hand-side.
528  for (unsigned int qp=0; qp<qrule.n_points(); ++qp)
529  {
530 
531  //Covariant basis at the quadrature point
532  MyVector3d a1;
533  a1 << dxyzdxi[qp](0), dxyzdxi[qp](1), dxyzdxi[qp](2);
534  MyVector3d a2;
535  a2 << dxyzdeta[qp](0), dxyzdeta[qp](1), dxyzdeta[qp](2);
536  MyVector3d n;
537  n = a1.cross(a2);
538  n /= n.norm();
539  MyMatrix3d F0;
540  F0 <<
541  a1(0), a2(0), n(0),
542  a1(1), a2(1), n(1),
543  a1(2), a2(2), n(2);
544 
545  //Contravariant basis
546  MyMatrix3d F0it;
547  F0it = F0.inverse().transpose();
548 
549  //Local orthonormal basis at the quadrature point
550  Real nx = n(0);
551  Real ny = n(1);
552  Real C = n(2);
553  MyMatrix3d Q;
554  if (std::abs(1.+C) < 1e-6)
555  {
556  Q <<
557  1, 0, 0,
558  0, -1, 0,
559  0, 0, -1;
560  }
561  else
562  {
563  Q <<
564  C+1./(1+C)*ny*ny, -1./(1+C)*nx*ny, nx,
565  -1./(1+C)*nx*ny, C+1./(1+C)*nx*nx, ny,
566  -nx, -ny, C;
567  }
568 
569  MyMatrix2d C0;
570  C0 = F0it.block<3,2>(0,0).transpose()*Q.block<3,2>(0,0);
571 
572  // Normal derivatives in reference coordinates
573  MyVector3d d2Xdxi2(d2xyzdxi2[qp](0), d2xyzdxi2[qp](1), d2xyzdxi2[qp](2));
574  MyVector3d d2Xdeta2(d2xyzdeta2[qp](0), d2xyzdeta2[qp](1), d2xyzdeta2[qp](2));
575  MyVector3d d2Xdxideta(d2xyzdxideta[qp](0), d2xyzdxideta[qp](1), d2xyzdxideta[qp](2));
576 
577 
578  MyMatrix2d b;
579  b <<
580  n.dot(d2Xdxi2), n.dot(d2Xdxideta),
581  n.dot(d2Xdxideta), n.dot(d2Xdeta2);
582 
583  MyVector3d dndxi = -b(0,0)*F0it.col(0) - b(0,1)*F0it.col(1);
584  MyVector3d dndeta = -b(1,0)*F0it.col(0) - b(1,1)*F0it.col(1);
585 
586  MyMatrix2d bhat;
587  bhat <<
588  F0it.col(1).dot(dndeta), -F0it.col(0).dot(dndeta),
589  -F0it.col(1).dot(dndxi), F0it.col(0).dot(dndxi);
590 
591  MyMatrix2d bc;
592  bc = bhat*C0;
593 
594  // Mean curvature
595  Real H = 0.5*(dndxi.dot(F0it.col(0))+dndeta.dot(F0it.col(1)));
596 
597  // Quadrature point reference coordinates
598  Real xi = qrule.qp(qp)(0);
599  Real eta = qrule.qp(qp)(1);
600 
601  // Preassemble the MITC4 shear strain matrix for all nodes as they involve
602  // cross references to midside nodes.
603  // The QUAD4 element has nodes X1,X2,X3,X4 with coordinates (xi,eta)
604  // in the reference element: (-1,-1),(1,-1),(1,1),(-1,1).
605  // The midside nodes are denoted A1=(X1+X2)/2, B2=(X2+X3)/2, A2=(X3+X4)/2, B1=(X4+X1)/2.
606 
607  // Normals at the midside nodes (average of normals at the edge corners).
608  // Multiplication by the assumed shear strain shape function.
609  MyVector3d nA1 = 0.5*(Qnode[0].col(2)+Qnode[1].col(2));
610  nA1 /= nA1.norm();
611  nA1 *= (1-eta)/4;
612  MyVector3d nB2 = 0.5*(Qnode[1].col(2)+Qnode[2].col(2));
613  nB2 /= nB2.norm();
614  nB2 *= (1+xi)/4;
615  MyVector3d nA2 = 0.5*(Qnode[2].col(2)+Qnode[3].col(2));
616  nA2 /= nA2.norm();
617  nA2 *= (1+eta)/4;
618  MyVector3d nB1 = 0.5*(Qnode[3].col(2)+Qnode[0].col(2));
619  nB1 /= nB1.norm();
620  nB1 *= (1-xi)/4;
621 
622  // Edge tangents
623  MyVector3d aA1 = 0.5*(X2-X1);
624  MyVector3d aA2 = 0.5*(X3-X4);
625  MyVector3d aB1 = 0.5*(X4-X1);
626  MyVector3d aB2 = 0.5*(X3-X2);
627 
628  // Contribution of the rotational dofs to the shear strain
629  MyVector2d AS1A1(-aA1.dot(Qnode[0].col(1)), aA1.dot(Qnode[0].col(0)));
630  MyVector2d AS2A1(-aA1.dot(Qnode[1].col(1)), aA1.dot(Qnode[1].col(0)));
631  AS1A1 *= (1-eta)/4;
632  AS2A1 *= (1-eta)/4;
633 
634  MyVector2d AS1A2(-aA2.dot(Qnode[3].col(1)), aA2.dot(Qnode[3].col(0)));
635  MyVector2d AS2A2(-aA2.dot(Qnode[2].col(1)), aA2.dot(Qnode[2].col(0)));
636  AS1A2 *= (1+eta)/4;
637  AS2A2 *= (1+eta)/4;
638 
639  MyVector2d AS1B1(-aB1.dot(Qnode[0].col(1)), aB1.dot(Qnode[0].col(0)));
640  MyVector2d AS2B1(-aB1.dot(Qnode[3].col(1)), aB1.dot(Qnode[3].col(0)));
641  AS1B1 *= (1-xi)/4;
642  AS2B1 *= (1-xi)/4;
643 
644  MyVector2d AS1B2(-aB2.dot(Qnode[1].col(1)), aB2.dot(Qnode[1].col(0)));
645  MyVector2d AS2B2(-aB2.dot(Qnode[2].col(1)), aB2.dot(Qnode[2].col(0)));
646  AS1B2 *= (1+xi)/4;
647  AS2B2 *= (1+xi)/4;
648 
649  // Store previous quantities in the shear strain matrices for each node
650  std::vector<MyMatrixXd> Bcnode;
651  MyMatrixXd Bc(2, 5);
652  // Node 1
653  Bc.block<1,3>(0,0) = -nA1.transpose();
654  Bc.block<1,2>(0,3) = AS1A1.transpose();
655  Bc.block<1,3>(1,0) = -nB1.transpose();
656  Bc.block<1,2>(1,3) = AS1B1.transpose();
657  Bcnode.push_back(Bc);
658  // Node 2
659  Bc.block<1,3>(0,0) = nA1.transpose();
660  Bc.block<1,2>(0,3) = AS2A1.transpose();
661  Bc.block<1,3>(1,0) = -nB2.transpose();
662  Bc.block<1,2>(1,3) = AS1B2.transpose();
663  Bcnode.push_back(Bc);
664  // Node 3
665  Bc.block<1,3>(0,0) = nA2.transpose();
666  Bc.block<1,2>(0,3) = AS2A2.transpose();
667  Bc.block<1,3>(1,0) = nB2.transpose();
668  Bc.block<1,2>(1,3) = AS2B2.transpose();
669  Bcnode.push_back(Bc);
670  // Node 4
671  Bc.block<1,3>(0,0) = -nA2.transpose();
672  Bc.block<1,2>(0,3) = AS1A2.transpose();
673  Bc.block<1,3>(1,0) = nB1.transpose();
674  Bc.block<1,2>(1,3) = AS2B1.transpose();
675  Bcnode.push_back(Bc);
676 
677  // Loop over all pairs of nodes I,J.
678  for (unsigned int i=0; i<n_var_dofs; ++i)
679  {
680  // Matrix B0, zeroth order (through thickness) membrane-bending strain
681  Real C1i = dphidxi[i][qp]*C0(0,0) + dphideta[i][qp]*C0(1,0);
682  Real C2i = dphidxi[i][qp]*C0(0,1) + dphideta[i][qp]*C0(1,1);
683 
684  MyMatrixXd B0I(3, 5);
685  B0I = MyMatrixXd::Zero(3, 5);
686  B0I.block<1,3>(0,0) = C1i*Q.col(0).transpose();
687  B0I.block<1,3>(1,0) = C2i*Q.col(1).transpose();
688  B0I.block<1,3>(2,0) = C2i*Q.col(0).transpose()+C1i*Q.col(1).transpose();
689 
690  // Matrix B1, first order membrane-bending strain
691  Real bc1i = dphidxi[i][qp]*bc(0,0) + dphideta[i][qp]*bc(1,0);
692  Real bc2i = dphidxi[i][qp]*bc(0,1) + dphideta[i][qp]*bc(1,1);
693 
694  MyVector2d V1i(-Q.col(0).dot(Qnode[i].col(1)),
695  Q.col(0).dot(Qnode[i].col(0)));
696 
697  MyVector2d V2i(-Q.col(1).dot(Qnode[i].col(1)),
698  Q.col(1).dot(Qnode[i].col(0)));
699 
700  MyMatrixXd B1I(3,5);
701  B1I = MyMatrixXd::Zero(3,5);
702  B1I.block<1,3>(0,0) = bc1i*Q.col(0).transpose();
703  B1I.block<1,3>(1,0) = bc2i*Q.col(1).transpose();
704  B1I.block<1,3>(2,0) = bc2i*Q.col(0).transpose()+bc1i*Q.col(1).transpose();
705 
706  B1I.block<1,2>(0,3) = C1i*V1i.transpose();
707  B1I.block<1,2>(1,3) = C2i*V2i.transpose();
708  B1I.block<1,2>(2,3) = C2i*V1i.transpose()+C1i*V2i.transpose();
709 
710  // Matrix B2, second order membrane-bending strain
711  MyMatrixXd B2I(3,5);
712  B2I = MyMatrixXd::Zero(3,5);
713 
714  B2I.block<1,2>(0,3) = bc1i*V1i.transpose();
715  B2I.block<1,2>(1,3) = bc2i*V2i.transpose();
716  B2I.block<1,2>(2,3) = bc2i*V1i.transpose()+bc1i*V2i.transpose();
717 
718  // Matrix Bc0, zeroth order shear strain
719  MyMatrixXd Bc0I(2,5);
720  Bc0I = C0.transpose()*Bcnode[i];
721 
722  // Matrix Bc1, first order shear strain
723  MyMatrixXd Bc1I(2,5);
724  Bc1I = bc.transpose()*Bcnode[i];
725 
726  // Drilling dof (in-plane rotation)
727  MyVector2d BdxiI(dphidxi[i][qp],dphideta[i][qp]);
728  MyVector2d BdI = C0.transpose()*BdxiI;
729 
730  for (unsigned int j=0; j<n_var_dofs; ++j)
731  {
732 
733  // Matrix B0, zeroth order membrane-bending strain
734  Real C1j = dphidxi[j][qp]*C0(0,0) + dphideta[j][qp]*C0(1,0);
735  Real C2j = dphidxi[j][qp]*C0(0,1) + dphideta[j][qp]*C0(1,1);
736 
737  MyMatrixXd B0J(3,5);
738  B0J = MyMatrixXd::Zero(3,5);
739  B0J.block<1,3>(0,0) = C1j*Q.col(0).transpose();
740  B0J.block<1,3>(1,0) = C2j*Q.col(1).transpose();
741  B0J.block<1,3>(2,0) = C2j*Q.col(0).transpose()+C1j*Q.col(1).transpose();
742 
743  // Matrix B1, first order membrane-bending strain
744  Real bc1j = dphidxi[j][qp]*bc(0,0) + dphideta[j][qp]*bc(1,0);
745  Real bc2j = dphidxi[j][qp]*bc(0,1) + dphideta[j][qp]*bc(1,1);
746 
747  MyVector2d V1j(-Q.col(0).dot(Qnode[j].col(1)),
748  Q.col(0).dot(Qnode[j].col(0)));
749 
750  MyVector2d V2j(-Q.col(1).dot(Qnode[j].col(1)),
751  Q.col(1).dot(Qnode[j].col(0)));
752 
753  MyMatrixXd B1J(3,5);
754  B1J = MyMatrixXd::Zero(3,5);
755  B1J.block<1,3>(0,0) = bc1j*Q.col(0).transpose();
756  B1J.block<1,3>(1,0) = bc2j*Q.col(1).transpose();
757  B1J.block<1,3>(2,0) = bc2j*Q.col(0).transpose()+bc1j*Q.col(1).transpose();
758 
759  B1J.block<1,2>(0,3) = C1j*V1j.transpose();
760  B1J.block<1,2>(1,3) = C2j*V2j.transpose();
761  B1J.block<1,2>(2,3) = C2j*V1j.transpose()+C1j*V2j.transpose();
762 
763  // Matrix B2, second order membrane-bending strain
764  MyMatrixXd B2J(3,5);
765  B2J = MyMatrixXd::Zero(3,5);
766 
767  B2J.block<1,2>(0,3) = bc1j*V1j.transpose();
768  B2J.block<1,2>(1,3) = bc2j*V2j.transpose();
769  B2J.block<1,2>(2,3) = bc2j*V1j.transpose()+bc1j*V2j.transpose();
770 
771  // Matrix Bc0, zeroth order shear strain
772  MyMatrixXd Bc0J(2, 5);
773  Bc0J = C0.transpose()*Bcnode[j];
774 
775  // Matrix Bc1, first order shear strain
776  MyMatrixXd Bc1J(2, 5);
777  Bc1J = bc.transpose()*Bcnode[j];
778 
779  // Drilling dof
780  MyVector2d BdxiJ(dphidxi[j][qp], dphideta[j][qp]);
781  MyVector2d BdJ = C0.transpose()*BdxiJ;
782 
783  // The total stiffness matrix coupling the nodes
784  // I and J is a sum of membrane, bending and shear contributions.
785  MyMatrixXd local_KIJ(5, 5);
786  local_KIJ = JxW[qp] * (
787  B0I.transpose() * Hm * B0J
788  + B2I.transpose() * Hf * B0J
789  + B0I.transpose() * Hf * B2J
790  + B1I.transpose() * Hf * B1J
791  + 2*H * B0I.transpose() * Hf * B1J
792  + 2*H * B1I.transpose() * Hf * B0J
793  + Bc0I.transpose() * Hc0 * Bc0J
794  + Bc1I.transpose() * Hc1 * Bc1J
795  + 2*H * Bc0I.transpose() * Hc1 * Bc1J
796  + 2*H * Bc1I.transpose() * Hc1 * Bc0J
797  );
798 
799  // Going from 5 to 6 dofs to add drilling dof
800  MyMatrixXd full_local_KIJ(6, 6);
801  full_local_KIJ = MyMatrixXd::Zero(6, 6);
802  full_local_KIJ.block<5,5>(0,0)=local_KIJ;
803 
804  // Drilling dof stiffness contribution
805  // Note that in the original book, there is a coefficient of
806  // alpha between 1e-4 and 1e-7 to make the fictitious
807  // drilling stiffness small while preventing the stiffness
808  // matrix from being singular. For this problem, we can use
809  // alpha = 1 to also get a good result.
810  //
811  // The explicit conversion to Real here is to work
812  // around an Eigen+boost::float128 incompatibility.
813  full_local_KIJ(5,5) = Real(Hf(0,0)*JxW[qp]*BdI.transpose()*BdJ);
814 
815  // Transform the stiffness matrix to global coordinates
816  MyMatrixXd global_KIJ(6,6);
817  MyMatrixXd TI(6,6);
818  TI = MyMatrixXd::Identity(6,6);
819  TI.block<3,3>(3,3) = Qnode[i].transpose();
820  MyMatrixXd TJ(6,6);
821  TJ = MyMatrixXd::Identity(6,6);
822  TJ.block<3,3>(3,3) = Qnode[j].transpose();
823  global_KIJ = TI.transpose()*full_local_KIJ*TJ;
824 
825  // Insert the components of the coupling stiffness
826  // matrix KIJ into the corresponding directional
827  // submatrices.
828  for (unsigned int k=0;k<6;k++)
829  for (unsigned int l=0;l<6;l++)
830  Ke_var[k][l](i,j) += global_KIJ(k,l);
831  }
832  }
833 
834  } // end of the quadrature point qp-loop
835 
836  if (distributed_load)
837  {
838  // Loop on shell faces
839  for (unsigned int shellface=0; shellface<2; shellface++)
840  {
841  std::vector<boundary_id_type> bids;
842  mesh.get_boundary_info().shellface_boundary_ids(elem, shellface, bids);
843 
844  for (std::size_t k=0; k<bids.size(); k++)
845  if (bids[k]==11) // sideset id for surface load
846  for (unsigned int qp=0; qp<qrule.n_points(); ++qp)
847  for (unsigned int i=0; i<n_var_dofs; ++i)
848  Fe_w(i) -= JxW[qp] * phi[i][qp];
849  }
850  }
851 
852  // The element matrix is now built for this element.
853  // Add it to the global matrix.
854 
855  dof_map.constrain_element_matrix_and_vector (Ke,Fe,dof_indices);
856 
857  system.matrix->add_matrix (Ke, dof_indices);
858  system.rhs->add_vector (Fe, dof_indices);
859 
860  }
861 
862  if (!distributed_load)
863  {
864  //Adding point load to the RHS
865 
866  //Pinch position
867  Point C(0, 3, 3);
868 
869  //Finish assembling rhs so we can set one value
870  system.rhs->close();
871 
872  for (const auto & node : mesh.node_ptr_range())
873  if (((*node) - C).norm() < 1e-3)
874  system.rhs->set(node->dof_number(0, 2, 0), -q/4);
875  }
876 
877 #else
878  // Avoid compiler warnings
879  libmesh_ignore(es, system_name);
880 #endif // defined(LIBMESH_HAVE_EIGEN) && defined(LIBMESH_ENABLE_SECOND_DERIVATIVES)
881 }

References std::abs(), libMesh::MeshBase::active_local_element_ptr_range(), libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), libMesh::NumericVector< T >::close(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::FEType::default_quadrature_order(), dim, libMesh::DofMap::dof_indices(), libMesh::Parameters::get(), libMesh::MeshBase::get_boundary_info(), libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::libmesh_ignore(), libMesh::ImplicitSystem::matrix, mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::MeshBase::node_ptr_range(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::DenseSubVector< T >::reposition(), libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, libMesh::SECOND, libMesh::NumericVector< T >::set(), libMesh::BoundaryInfo::shellface_boundary_ids(), and libMesh::System::variable_type().

Referenced by main().

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 79 of file miscellaneous_ex11.C.

80 {
81  // Initialize libMesh.
82  LibMeshInit init (argc, argv);
83 
84  // This example requires a linear solver package.
85  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
86  "--enable-petsc, --enable-trilinos, or --enable-eigen");
87 
88  // Skip this 3D example if libMesh was compiled as 1D/2D-only.
89  libmesh_example_requires (3 == LIBMESH_DIM, "3D support");
90 
91  // Skip this example without --enable-node-valence
92 #ifndef LIBMESH_ENABLE_NODE_VALENCE
93  libmesh_example_requires (false, "--enable-node-valence");
94 #endif
95 
96  // Skip this example without --enable-amr; requires MeshRefinement
97 #ifndef LIBMESH_ENABLE_AMR
98  libmesh_example_requires(false, "--enable-amr");
99 #else
100 
101  // Skip this example without --enable-second; requires d2phi
102 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES
103  libmesh_example_requires(false, "--enable-second");
104 #elif LIBMESH_DIM > 2
105 
106  // Create a 2D mesh distributed across the default MPI communicator.
107  // Subdivision surfaces do not appear to work with DistributedMesh yet.
108  ReplicatedMesh mesh (init.comm(), 2);
109 
110  // Read the coarse square mesh.
111  mesh.read ("square_mesh.off");
112 
113  // Resize the square plate to edge length L.
114  const Real L = 100.;
116 
117  // Quadrisect the mesh triangles a few times to obtain a
118  // finer mesh. Subdivision surface elements require the
119  // refinement data to be removed afterward.
120  MeshRefinement mesh_refinement (mesh);
121  mesh_refinement.uniformly_refine (3);
123 
124  // Write the mesh before the ghost elements are added.
125 #if defined(LIBMESH_HAVE_VTK)
126  VTKIO(mesh).write ("without_ghosts.pvtu");
127 #endif
128 #if defined(LIBMESH_HAVE_EXODUS_API)
129  ExodusII_IO(mesh).write ("without_ghosts.e");
130 #endif
131 
132  // Print information about the triangulated mesh to the screen.
133  mesh.print_info();
134 
135  // Turn the triangulated mesh into a subdivision mesh
136  // and add an additional row of "ghost" elements around
137  // it in order to complete the extended local support of
138  // the triangles at the boundaries. If the second
139  // argument is set to true, the outermost existing
140  // elements are converted into ghost elements, and the
141  // actual physical mesh is thus getting smaller.
143 
144  // Print information about the subdivision mesh to the screen.
145  mesh.print_info();
146 
147  // Write the mesh with the ghost elements added.
148  // Compare this to the original mesh to see the difference.
149 #if defined(LIBMESH_HAVE_VTK)
150  VTKIO(mesh).write ("with_ghosts.pvtu");
151 #endif
152 #if defined(LIBMESH_HAVE_EXODUS_API)
153  ExodusII_IO(mesh).write ("with_ghosts.e");
154 #endif
155 
156  // Create an equation systems object.
157  EquationSystems equation_systems (mesh);
158 
159  // Declare the system and its variables.
160  // Create a linear implicit system named "Shell".
161  LinearImplicitSystem & system = equation_systems.add_system<LinearImplicitSystem> ("Shell");
162 
163  // Add the three translational deformation variables
164  // "u", "v", "w" to "Shell". Since subdivision shell
165  // elements meet the C1-continuity requirement, no
166  // rotational or other auxiliary variables are needed.
167  // Loop Subdivision Elements are always interpolated
168  // by quartic box splines, hence the order must always
169  // be FOURTH.
170  system.add_variable ("u", FOURTH, SUBDIVISION);
171  system.add_variable ("v", FOURTH, SUBDIVISION);
172  system.add_variable ("w", FOURTH, SUBDIVISION);
173 
174  // Give the system a pointer to the matrix and rhs assembly
175  // function.
177 
178  // Use the parameters of the equation systems object to
179  // tell the shell system about the material properties, the
180  // shell thickness, and the external load.
181  const Real h = 1.;
182  const Real E = 1.e7;
183  const Real nu = 0.;
184  const Real q = 1.;
185  equation_systems.parameters.set<Real> ("thickness") = h;
186  equation_systems.parameters.set<Real> ("young's modulus") = E;
187  equation_systems.parameters.set<Real> ("poisson ratio") = nu;
188  equation_systems.parameters.set<Real> ("uniform load") = q;
189 
190  // Initialize the data structures for the equation system.
191  equation_systems.init();
192 
193  // Print information about the system to the screen.
194  equation_systems.print_info();
195 
196  // Solve the linear system.
197  system.solve();
198 
199  // After solving the system, write the solution to a VTK
200  // or ExodusII output file ready for import in, e.g.,
201  // Paraview.
202 #if defined(LIBMESH_HAVE_VTK)
203  VTKIO(mesh).write_equation_systems ("out.pvtu", equation_systems);
204 #endif
205 #if defined(LIBMESH_HAVE_EXODUS_API)
206  ExodusII_IO(mesh).write_equation_systems ("out.e", equation_systems);
207 #endif
208 
209  // Find the center node to measure the maximum deformation of the plate.
210  Node * center_node = 0;
211  Real nearest_dist_sq = mesh.point(0).norm_sq();
212  for (unsigned int nid=1; nid<mesh.n_nodes(); ++nid)
213  {
214  const Real dist_sq = mesh.point(nid).norm_sq();
215  if (dist_sq < nearest_dist_sq)
216  {
217  nearest_dist_sq = dist_sq;
218  center_node = mesh.node_ptr(nid);
219  }
220  }
221 
222  // Finally, we evaluate the z-displacement "w" at the center node.
223  const unsigned int w_var = system.variable_number ("w");
224  dof_id_type w_dof = center_node->dof_number (system.number(), w_var, 0);
225  Number w = 0;
226  if (w_dof >= system.get_dof_map().first_dof() &&
227  w_dof < system.get_dof_map().end_dof())
228  w = system.current_solution(w_dof);
229  system.comm().sum(w);
230 
231 
232  // The analytic solution for the maximum displacement of
233  // a clamped square plate in pure bending, from Taylor,
234  // Govindjee, Commun. Numer. Meth. Eng. 20, 757-765, 2004.
235  const Real D = E * h*h*h / (12*(1-nu*nu));
236  const Real w_analytic = 0.001265319 * L*L*L*L * q / D;
237 
238  // Print the finite element solution and the analytic
239  // prediction of the maximum displacement of the clamped
240  // square plate to the screen.
241  libMesh::out << "z-displacement of the center point: " << w << std::endl;
242  libMesh::out << "Analytic solution for pure bending: " << w_analytic << std::endl;
243 
244 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES, LIBMESH_DIM > 2
245 
246 #endif // #ifdef LIBMESH_ENABLE_AMR
247 
248  // All done.
249  return 0;
250 }

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_shell(), libMesh::System::attach_assemble_function(), libMesh::ParallelObject::comm(), libMesh::System::current_solution(), libMesh::default_solver_package(), libMesh::DofObject::dof_number(), libMesh::DofMap::end_dof(), libMesh::DofMap::first_dof(), libMesh::MeshTools::Modification::flatten(), libMesh::FOURTH, libMesh::System::get_dof_map(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::MeshBase::n_nodes(), libMesh::MeshBase::node_ptr(), libMesh::TypeVector< T >::norm_sq(), libMesh::System::number(), libMesh::out, libMesh::EquationSystems::parameters, libMesh::MeshBase::point(), libMesh::MeshTools::Subdivision::prepare_subdivision_mesh(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshBase::read(), libMesh::Real, libMesh::MeshTools::Modification::scale(), libMesh::Parameters::set(), libMesh::LinearImplicitSystem::solve(), libMesh::SUBDIVISION, libMesh::MeshRefinement::uniformly_refine(), libMesh::System::variable_number(), libMesh::ExodusII_IO::write(), libMesh::VTKIO::write(), and libMesh::MeshOutput< MT >::write_equation_systems().

libMesh::TypeTensor::transpose
TypeTensor< T > transpose() const
Definition: type_tensor.h:1068
MyVector3d
Eigen::Matrix< libMesh::Real, 3, 1 > MyVector3d
Definition: miscellaneous_ex12.C:71
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::MeshTools::Subdivision::next
static const unsigned int next[3]
A lookup table for the increment modulo 3 operation, for iterating through the three nodes per elemen...
Definition: mesh_subdivision_support.h:102
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::MeshBase::read
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.
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DofMap::dof_indices
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:1967
libMesh::MeshBase::point
virtual const Point & point(const dof_id_type i) const =0
MyMatrixXd
Eigen::Matrix< libMesh::Real, Eigen::Dynamic, Eigen::Dynamic > MyMatrixXd
Definition: miscellaneous_ex12.C:72
libMesh::BoundaryInfo::shellface_boundary_ids
void shellface_boundary_ids(const Elem *const elem, const unsigned short int shellface, std::vector< boundary_id_type > &vec_to_fill) const
Definition: boundary_info.C:1133
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
assemble_shell
void assemble_shell(EquationSystems &es, const std::string &system_name)
Definition: miscellaneous_ex12.C:348
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::DenseMatrix< Number >
libMesh::LinearImplicitSystem::solve
virtual void solve() override
Assembles & solves the linear system A*x=b.
Definition: linear_implicit_system.C:108
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::DofMap::first_dof
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map.h:650
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::VTKIO
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::VectorValue< Real >
libMesh::ReplicatedMesh
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Definition: replicated_mesh.h:47
libMesh::System::current_solution
Number current_solution(const dof_id_type global_dof_number) const
Definition: system.C:194
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::System::attach_assemble_function
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:1755
libMesh::NumericVector::add_vector
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].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
MyMatrix3d
Eigen::Matrix< libMesh::Real, 3, 3 > MyMatrix3d
Definition: miscellaneous_ex12.C:69
libMesh::TypeVector::norm_sq
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_vector.h:986
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::Tri3Subdivision::is_ghost
bool is_ghost() const
Definition: face_tri3_subdivision.h:116
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::System::add_variable
unsigned int add_variable(const std::string &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:1069
libMesh::TensorValue
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:53
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::FEType::default_quadrature_rule
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:31
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::SparseMatrix::add
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DofMap::constrain_element_matrix_and_vector
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:2034
libMesh::DenseSubVector< Number >
libMesh::MeshOutput::write_equation_systems
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
libMesh::FOURTH
Definition: enum_order.h:45
libMesh::System::variable_type
const FEType & variable_type(const unsigned int i) const
Definition: system.h:2233
MyMatrix2d
Eigen::Matrix< libMesh::Real, 2, 2 > MyMatrix2d
Definition: miscellaneous_ex12.C:68
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_matrix
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.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshTools::Subdivision::prepare_subdivision_mesh
void prepare_subdivision_mesh(MeshBase &mesh, bool ghosted=false)
Prepares the mesh for use with subdivision elements.
Definition: mesh_subdivision_support.C:160
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::MeshTools::Modification::flatten
void flatten(MeshBase &mesh)
Removes all the refinement tree structure of Mesh, leaving only the highest-level (most-refined) elem...
Definition: mesh_modification.C:1265
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
F0
Definition: assembly.h:127
libMesh::ExodusII_IO::write
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
Definition: exodusII_io.C:1338
libMesh::MeshTools::Subdivision::prev
static const unsigned int prev[3]
A lookup table for the decrement modulo 3 operation, for iterating through the three nodes per elemen...
Definition: mesh_subdivision_support.h:108
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::SUBDIVISION
Definition: enum_fe_family.h:55
libMesh::TypeVector::cross
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:920
libMesh::VTKIO::write
virtual void write(const std::string &) override
Output the mesh without solutions to a .pvtu file.
libMesh::TRI3SUBDIVISION
Definition: enum_elem_type.h:69
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
libMesh::MeshTools::Modification::scale
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
Scales the mesh.
Definition: mesh_modification.C:243
libMesh::TypeVector::norm
auto norm() const -> decltype(std::norm(T()))
Definition: type_vector.h:955
libMesh::Tri3Subdivision
The Tri3Subdivision element is a three-noded subdivision surface shell element used in mechanics calc...
Definition: face_tri3_subdivision.h:40
libMesh::Parameters::get
const T & get(const std::string &) const
Definition: parameters.h:421
libMesh::out
OStreamProxy out
libMesh::Elem::node_ptr
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2009
libMesh::DenseVector< Number >
MyVector2d
Eigen::Matrix< libMesh::Real, 2, 1 > MyVector2d
Definition: miscellaneous_ex12.C:70
libMesh::DofMap::end_dof
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map.h:692
libMesh::EquationSystems::parameters
Parameters parameters
Data structure holding arbitrary parameters.
Definition: equation_systems.h:557