41 #include "libmesh/libmesh.h" 
   42 #include "libmesh/mesh.h" 
   43 #include "libmesh/mesh_generation.h" 
   44 #include "libmesh/exodusII_io.h" 
   45 #include "libmesh/linear_implicit_system.h" 
   46 #include "libmesh/equation_systems.h" 
   47 #include "libmesh/fe.h" 
   48 #include "libmesh/quadrature_gauss.h" 
   49 #include "libmesh/dof_map.h" 
   50 #include "libmesh/sparse_matrix.h" 
   51 #include "libmesh/numeric_vector.h" 
   52 #include "libmesh/dense_matrix.h" 
   53 #include "libmesh/dense_submatrix.h" 
   54 #include "libmesh/dense_vector.h" 
   55 #include "libmesh/dense_subvector.h" 
   56 #include "libmesh/perf_log.h" 
   57 #include "libmesh/elem.h" 
   58 #include "libmesh/boundary_info.h" 
   59 #include "libmesh/zero_function.h" 
   60 #include "libmesh/dirichlet_boundaries.h" 
   61 #include "libmesh/string_to_enum.h" 
   62 #include "libmesh/getpot.h" 
   63 #include "libmesh/enum_solver_package.h" 
   70                          const std::string & system_name);
 
   80 int main (
int argc, 
char ** argv)
 
   87                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   93   const unsigned int dim = 2;
 
   96   libmesh_example_requires(
dim <= LIBMESH_DIM, 
"2D support");
 
   99 #ifndef LIBMESH_ENABLE_DIRICHLET 
  100   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  122 #ifdef LIBMESH_ENABLE_DIRICHLET 
  135   std::set<boundary_id_type> boundary_ids;
 
  136   boundary_ids.insert(3);
 
  139   std::vector<unsigned int> variables(2);
 
  140   variables[0] = u_var;
 
  141   variables[1] = v_var;
 
  154 #endif // LIBMESH_ENABLE_DIRICHLET 
  157   equation_systems.
init();
 
  166 #ifdef LIBMESH_HAVE_EXODUS_API 
  168 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  176                          const std::string & libmesh_dbg_var(system_name))
 
  178   libmesh_assert_equal_to (system_name, 
"Elasticity");
 
  194   fe->attach_quadrature_rule (&qrule);
 
  198   fe_face->attach_quadrature_rule (&qface);
 
  200   const std::vector<Real> & JxW = fe->get_JxW();
 
  201   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  215   std::vector<dof_id_type> dof_indices;
 
  216   std::vector<dof_id_type> dof_indices_u;
 
  217   std::vector<dof_id_type> dof_indices_v;
 
  218   std::vector<dof_id_type> dof_indices_lambda;
 
  225       dof_map.
dof_indices (elem, dof_indices_lambda, lambda_var);
 
  227       const unsigned int n_dofs   = dof_indices.size();
 
  228       const unsigned int n_u_dofs = dof_indices_u.size();
 
  229       const unsigned int n_v_dofs = dof_indices_v.size();
 
  230       const unsigned int n_lambda_dofs = dof_indices_lambda.size();
 
  234       Ke.
resize (n_dofs, n_dofs);
 
  237       Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
 
  238       Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
 
  240       Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
 
  241       Kvv.
reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
 
  244       Kv_lambda.
reposition (v_var*n_u_dofs, v_var*n_u_dofs+n_v_dofs, n_v_dofs, 1);
 
  245       Klambda_v.reposition (v_var*n_v_dofs+n_v_dofs, v_var*n_v_dofs, 1, n_v_dofs);
 
  247       Fu.reposition (u_var*n_u_dofs, n_u_dofs);
 
  250       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  252           for (
unsigned int i=0; i<n_u_dofs; i++)
 
  253             for (
unsigned int j=0; j<n_u_dofs; j++)
 
  256                 unsigned int C_i, C_j, C_k, C_l;
 
  272           for (
unsigned int i=0; i<n_u_dofs; i++)
 
  273             for (
unsigned int j=0; j<n_v_dofs; j++)
 
  276                 unsigned int C_i, C_j, C_k, C_l;
 
  292           for (
unsigned int i=0; i<n_v_dofs; i++)
 
  293             for (
unsigned int j=0; j<n_u_dofs; j++)
 
  296                 unsigned int C_i, C_j, C_k, C_l;
 
  312           for (
unsigned int i=0; i<n_v_dofs; i++)
 
  313             for (
unsigned int j=0; j<n_v_dofs; j++)
 
  316                 unsigned int C_i, C_j, C_k, C_l;
 
  334         std::vector<boundary_id_type> bc_ids;
 
  335         for (
auto side : elem->side_index_range())
 
  336           if (elem->neighbor_ptr(side) == 
nullptr)
 
  340               const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
 
  341               const std::vector<Real> & JxW_face = fe_face->get_JxW();
 
  343               fe_face->reinit(elem, side);
 
  345               for (std::vector<boundary_id_type>::const_iterator b =
 
  346                      bc_ids.begin(); b != bc_ids.end(); ++b)
 
  349                   for (
unsigned int qp=0; qp<qface.n_points(); qp++)
 
  353                         for (
unsigned int i=0; i<n_v_dofs; i++)
 
  354                           Fv(i) += JxW_face[qp] * (-1.) * phi_face[i][qp];
 
  359                           for (
unsigned int i=0; i<n_v_dofs; i++)
 
  360                             for (
unsigned int j=0; j<n_lambda_dofs; j++)
 
  361                               Kv_lambda(i,j) += JxW_face[qp] * (-1.) * phi_face[i][qp];
 
  363                           for (
unsigned int i=0; i<n_lambda_dofs; i++)
 
  364                             for (
unsigned int j=0; j<n_v_dofs; j++)
 
  365                               Klambda_v(i,j) += JxW_face[qp] * (-1.) * phi_face[j][qp];
 
  388   const Real lambda_1 = nu / ((1. + nu) * (1. - 2.*nu));
 
  389   const Real lambda_2 = 0.5 / (1 + nu);
 
  392   Real delta_ij = (i == j) ? 1. : 0.;
 
  393   Real delta_il = (i == l) ? 1. : 0.;
 
  394   Real delta_ik = (i == k) ? 1. : 0.;
 
  395   Real delta_jl = (j == l) ? 1. : 0.;
 
  396   Real delta_jk = (j == k) ? 1. : 0.;
 
  397   Real delta_kl = (k == l) ? 1. : 0.;
 
  399   return lambda_1 * delta_ij * delta_kl + lambda_2 * (delta_ik * delta_jl + delta_il * delta_jk);