36 #include "libmesh/libmesh.h" 
   37 #include "libmesh/mesh.h" 
   38 #include "libmesh/mesh_generation.h" 
   39 #include "libmesh/exodusII_io.h" 
   40 #include "libmesh/gnuplot_io.h" 
   41 #include "libmesh/linear_implicit_system.h" 
   42 #include "libmesh/equation_systems.h" 
   43 #include "libmesh/fe.h" 
   44 #include "libmesh/quadrature_gauss.h" 
   45 #include "libmesh/dof_map.h" 
   46 #include "libmesh/sparse_matrix.h" 
   47 #include "libmesh/numeric_vector.h" 
   48 #include "libmesh/dense_matrix.h" 
   49 #include "libmesh/dense_submatrix.h" 
   50 #include "libmesh/dense_vector.h" 
   51 #include "libmesh/dense_subvector.h" 
   52 #include "libmesh/perf_log.h" 
   53 #include "libmesh/elem.h" 
   54 #include "libmesh/boundary_info.h" 
   55 #include "libmesh/zero_function.h" 
   56 #include "libmesh/dirichlet_boundaries.h" 
   57 #include "libmesh/string_to_enum.h" 
   58 #include "libmesh/getpot.h" 
   59 #include "libmesh/enum_solver_package.h" 
   66                          const std::string & system_name);
 
   76 int main (
int argc, 
char ** argv)
 
   83                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   86   const unsigned int dim = 2;
 
   89   libmesh_example_requires(
dim <= LIBMESH_DIM, 
"2D support");
 
   92 #ifndef LIBMESH_ENABLE_DIRICHLET 
   93   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  118 #ifdef LIBMESH_ENABLE_DIRICHLET 
  126   std::set<boundary_id_type> boundary_ids;
 
  127   boundary_ids.insert(3);
 
  130   std::vector<unsigned int> variables(2);
 
  131   variables[0] = u_var; variables[1] = v_var;
 
  144 #endif // LIBMESH_ENABLE_DIRICHLET 
  147   equation_systems.
init();
 
  156 #ifdef LIBMESH_HAVE_EXODUS_API 
  158 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  166                          const std::string & libmesh_dbg_var(system_name))
 
  168   libmesh_assert_equal_to (system_name, 
"Elasticity");
 
  183   fe->attach_quadrature_rule (&qrule);
 
  187   fe_face->attach_quadrature_rule (&qface);
 
  189   const std::vector<Real> & JxW = fe->get_JxW();
 
  190   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  203   std::vector<dof_id_type> dof_indices;
 
  204   std::vector<dof_id_type> dof_indices_u;
 
  205   std::vector<dof_id_type> dof_indices_v;
 
  213       const unsigned int n_dofs   = dof_indices.size();
 
  214       const unsigned int n_u_dofs = dof_indices_u.size();
 
  215       const unsigned int n_v_dofs = dof_indices_v.size();
 
  219       Ke.
resize (n_dofs, n_dofs);
 
  222       Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
 
  223       Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
 
  225       Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
 
  226       Kvv.
reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
 
  228       Fu.reposition (u_var*n_u_dofs, n_u_dofs);
 
  231       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  233           for (
unsigned int i=0; i<n_u_dofs; i++)
 
  234             for (
unsigned int j=0; j<n_u_dofs; j++)
 
  237                 unsigned int C_i, C_j, C_k, C_l;
 
  253           for (
unsigned int i=0; i<n_u_dofs; i++)
 
  254             for (
unsigned int j=0; j<n_v_dofs; j++)
 
  257                 unsigned int C_i, C_j, C_k, C_l;
 
  273           for (
unsigned int i=0; i<n_v_dofs; i++)
 
  274             for (
unsigned int j=0; j<n_u_dofs; j++)
 
  277                 unsigned int C_i, C_j, C_k, C_l;
 
  293           for (
unsigned int i=0; i<n_v_dofs; i++)
 
  294             for (
unsigned int j=0; j<n_v_dofs; j++)
 
  297                 unsigned int C_i, C_j, C_k, C_l;
 
  315         for (
auto side : elem->side_index_range())
 
  316           if (elem->neighbor_ptr(side) == 
nullptr)
 
  318               const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
 
  319               const std::vector<Real> & JxW_face = fe_face->get_JxW();
 
  321               fe_face->reinit(elem, side);
 
  325                   for (
unsigned int qp=0; qp<qface.n_points(); qp++)
 
  326                     for (
unsigned int i=0; i<n_v_dofs; i++)
 
  327                       Fv(i) += JxW_face[qp] * (-1.) * phi_face[i][qp];
 
  348   const Real lambda_1 = nu / ((1. + nu) * (1. - 2.*nu));
 
  349   const Real lambda_2 = 0.5 / (1 + nu);
 
  352   Real delta_ij = (i == j) ? 1. : 0.;
 
  353   Real delta_il = (i == l) ? 1. : 0.;
 
  354   Real delta_ik = (i == k) ? 1. : 0.;
 
  355   Real delta_jl = (j == l) ? 1. : 0.;
 
  356   Real delta_jk = (j == k) ? 1. : 0.;
 
  357   Real delta_kl = (k == l) ? 1. : 0.;
 
  359   return lambda_1 * delta_ij * delta_kl + lambda_2 * (delta_ik * delta_jl + delta_il * delta_jk);