36 #include "libmesh/libmesh.h" 
   37 #include "libmesh/mesh.h" 
   38 #include "libmesh/linear_implicit_system.h" 
   39 #include "libmesh/equation_systems.h" 
   40 #include "libmesh/fe.h" 
   41 #include "libmesh/quadrature.h" 
   42 #include "libmesh/node.h" 
   43 #include "libmesh/elem.h" 
   44 #include "libmesh/dof_map.h" 
   45 #include "libmesh/quadrature_gauss.h" 
   46 #include "libmesh/vector_value.h" 
   47 #include "libmesh/tensor_value.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/sparse_matrix.h" 
   53 #include "libmesh/numeric_vector.h" 
   54 #include "libmesh/exodusII_io.h" 
   55 #include "libmesh/dirichlet_boundaries.h" 
   56 #include "libmesh/zero_function.h" 
   57 #include "libmesh/linear_solver.h" 
   58 #include "libmesh/getpot.h" 
   59 #include "libmesh/enum_solver_package.h" 
   60 #include "libmesh/enum_solver_type.h" 
   61 #include "libmesh/parallel.h" 
   64 #ifdef LIBMESH_HAVE_EIGEN 
   65 #include "libmesh/ignore_warnings.h" 
   66 # include <Eigen/Dense> 
   67 #include "libmesh/restore_warnings.h" 
   74 typedef Eigen::Matrix<libMesh::Real, Eigen::Dynamic, Eigen::Dynamic> 
MyMatrixXd;
 
   84                      const std::string & system_name);
 
   87 int main (
int argc, 
char ** argv)
 
   93   libmesh_example_requires (3 == LIBMESH_DIM, 
"3D support");
 
   96 #ifndef LIBMESH_ENABLE_DIRICHLET 
   97   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  101 #ifndef LIBMESH_HAVE_EXODUS_API 
  102   libmesh_example_requires (
false, 
"ExodusII support");
 
  105 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  106   libmesh_example_requires (
false, 
"second derivatives enabled");
 
  111 #ifndef LIBMESH_HAVE_EIGEN 
  112   libmesh_example_requires(
false, 
"--enable-eigen");
 
  117 #ifndef LIBMESH_HAVE_XDR 
  118   libmesh_example_requires (
false, 
"XDR support");
 
  123 #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION 
  124   libmesh_example_requires (
false, 
"--disable-quadruple-precision");
 
  128   GetPot command_line (argc, argv);
 
  129   int distributed_load = 0;
 
  130   if (command_line.search(1, 
"-distributed_load"))
 
  131     distributed_load = command_line.next(distributed_load);
 
  181   equation_systems.
parameters.
set<
bool>(
"distributed load")  = (distributed_load != 0);
 
  191 #ifdef LIBMESH_ENABLE_DIRICHLET 
  194     std::set<boundary_id_type> boundary_ids;
 
  195     boundary_ids.insert(7);
 
  196     unsigned int variables[] = {2, 3, 4};
 
  203        std::vector<unsigned int>(variables, variables+3), zf,
 
  209     std::set<boundary_id_type> boundary_ids;
 
  210     boundary_ids.insert(8);
 
  211     unsigned int variables[] = {1, 3, 5};
 
  216        std::vector<unsigned int>(variables, variables+3), zf,
 
  222     std::set<boundary_id_type> boundary_ids;
 
  223     boundary_ids.insert(9);
 
  224     unsigned int variables[] = {0, 4, 5};
 
  229        std::vector<unsigned int>(variables, variables+3), zf,
 
  235     std::set<boundary_id_type> boundary_ids;
 
  236     boundary_ids.insert(10);
 
  237     unsigned int variables[] = {0, 2, 4};
 
  242        std::vector<unsigned int>(variables, variables+3), zf,
 
  246 #endif // LIBMESH_ENABLE_DIRICHLET 
  249   equation_systems.
init();
 
  269   if (distributed_load==0)
 
  272       Node * node_C = 
nullptr;
 
  273       Point point_C(0, 3, 3);
 
  275         Real nearest_dist_sq = std::numeric_limits<Real>::max();
 
  283             if (dist_sq < nearest_dist_sq)
 
  285                 nearest_dist_sq = dist_sq;
 
  291         unsigned int minrank = 0;
 
  292         system.
comm().minloc(nearest_dist_sq, minrank);
 
  298           nearest_node_id = node_C->
id();
 
  299         system.
comm().broadcast(nearest_node_id, minrank);
 
  316       system.
comm().sum(w);
 
  319       Number w_C_bar = -E*h*w/q;
 
  320       const Real w_C_bar_analytic = 164.24;
 
  324       libMesh::out << 
"z-displacement of the point C: " << w_C_bar << std::endl;
 
  325       libMesh::out << 
"Analytic solution: " << w_C_bar_analytic << std::endl;
 
  329       Point point_D(0, 0, 3);
 
  334       const Real v_D_bar_analytic = 4.114;
 
  338       libMesh::out << 
"y-displacement of the point D: " << v_D_bar << std::endl;
 
  339       libMesh::out << 
"Analytic solution: " << v_D_bar_analytic << std::endl;
 
  351                      const std::string & system_name)
 
  355 #if defined(LIBMESH_HAVE_EIGEN) && defined(LIBMESH_ENABLE_SECOND_DERIVATIVES) 
  358   libmesh_assert_equal_to (system_name, 
"Shell");
 
  372   const bool distributed_load  = es.
parameters.
get<
bool> (
"distributed load");
 
  379     0., 0., 0.5 * (1-nu);
 
  380   Hm *= h * E/(1-nu*nu);
 
  387     0., 0., 0.5 * (1-nu);
 
  388   Hf *= h*h*h/12 * E/(1-nu*nu);
 
  392   Hc0 *= h * 5./6*E/(2*(1+nu));
 
  395   Hc1 *= h*h*h/12 * 5./6*E/(2*(1+nu));
 
  403   fe->attach_quadrature_rule (&qrule);
 
  406   const std::vector<Real> & JxW = fe->get_JxW();
 
  410   const std::vector<RealGradient> & dxyzdxi = fe->get_dxyzdxi();
 
  411   const std::vector<RealGradient> & dxyzdeta = fe->get_dxyzdeta();
 
  412   const std::vector<RealGradient> & d2xyzdxi2 = fe->get_d2xyzdxi2();
 
  413   const std::vector<RealGradient> & d2xyzdeta2 = fe->get_d2xyzdeta2();
 
  414   const std::vector<RealGradient> & d2xyzdxideta = fe->get_d2xyzdxideta();
 
  415   const std::vector<std::vector<Real>> & dphidxi = fe->get_dphidxi();
 
  416   const std::vector<std::vector<Real>> & dphideta = fe->get_dphideta();
 
  417   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  446   std::vector<dof_id_type> dof_indices;
 
  447   std::vector<std::vector<dof_id_type>> dof_indices_var(6);
 
  454       for (
unsigned int var=0; var<6; var++)
 
  455         dof_map.
dof_indices (elem, dof_indices_var[var], var);
 
  457       const unsigned int n_dofs   = dof_indices.size();
 
  458       const unsigned int n_var_dofs = dof_indices_var[0].size();
 
  461       std::vector<Point> nodes;
 
  462       for (
unsigned int i=0; i<elem->n_nodes(); ++i)
 
  463         nodes.push_back(elem->reference_elem()->node_ref(i));
 
  464       fe->reinit (elem, &nodes);
 
  467       std::vector<MyMatrix3d> Qnode;
 
  468       for (
unsigned int i=0; i<elem->n_nodes(); ++i)
 
  471           a1 << dxyzdxi[i](0), dxyzdxi[i](1), dxyzdxi[i](2);
 
  473           a2 << dxyzdeta[i](0), dxyzdeta[i](1), dxyzdeta[i](2);
 
  494                 C+1./(1+C)*ny*ny, -1./(1+C)*nx*ny, nx,
 
  495                 -1./(1+C)*nx*ny, C+1./(1+C)*nx*nx, ny,
 
  501       Ke.
resize (n_dofs, n_dofs);
 
  502       for (
unsigned int var_i=0; var_i<6; var_i++)
 
  503         for (
unsigned int var_j=0; var_j<6; var_j++)
 
  504           Ke_var[var_i][var_j].reposition (var_i*n_var_dofs, var_j*n_var_dofs, n_var_dofs, n_var_dofs);
 
  513       for (
unsigned int qp=0; qp<qrule.
n_points(); ++qp)
 
  518           a1 << dxyzdxi[qp](0), dxyzdxi[qp](1), dxyzdxi[qp](2);
 
  520           a2 << dxyzdeta[qp](0), dxyzdeta[qp](1), dxyzdeta[qp](2);
 
  532           F0it = 
F0.inverse().transpose();
 
  549                 C+1./(1+C)*ny*ny, -1./(1+C)*nx*ny, nx,
 
  550                 -1./(1+C)*nx*ny, C+1./(1+C)*nx*nx, ny,
 
  555           C0 = F0it.block<3,2>(0,0).transpose()*Q.block<3,2>(0,0);
 
  558           MyVector3d d2Xdxi2(d2xyzdxi2[qp](0), d2xyzdxi2[qp](1), d2xyzdxi2[qp](2));
 
  559           MyVector3d d2Xdeta2(d2xyzdeta2[qp](0), d2xyzdeta2[qp](1), d2xyzdeta2[qp](2));
 
  560           MyVector3d d2Xdxideta(d2xyzdxideta[qp](0), d2xyzdxideta[qp](1), d2xyzdxideta[qp](2));
 
  564             n.dot(d2Xdxi2), n.dot(d2Xdxideta),
 
  565             n.dot(d2Xdxideta), n.dot(d2Xdeta2);
 
  567           MyVector3d dndxi = -b(0,0)*F0it.col(0) - b(0,1)*F0it.col(1);
 
  568           MyVector3d dndeta = -b(1,0)*F0it.col(0) - b(1,1)*F0it.col(1);
 
  572             F0it.col(1).dot(dndeta), -F0it.col(0).dot(dndeta),
 
  573             -F0it.col(1).dot(dndxi), F0it.col(0).dot(dndxi);
 
  579           Real H = 0.5*(dndxi.dot(F0it.col(0))+dndeta.dot(F0it.col(1)));
 
  582           for (
unsigned int i=0; i<n_var_dofs; ++i)
 
  585               Real C1i = dphidxi[i][qp]*C0(0,0) + dphideta[i][qp]*C0(1,0);
 
  586               Real C2i = dphidxi[i][qp]*C0(0,1) + dphideta[i][qp]*C0(1,1);
 
  589               B0I = MyMatrixXd::Zero(3, 5);
 
  590               B0I.block<1,3>(0,0) = C1i*Q.col(0).transpose();
 
  591               B0I.block<1,3>(1,0) = C2i*Q.col(1).transpose();
 
  592               B0I.block<1,3>(2,0) = C2i*Q.col(0).transpose()+C1i*Q.col(1).transpose();
 
  595               Real bc1i = dphidxi[i][qp]*bc(0,0) + dphideta[i][qp]*bc(1,0);
 
  596               Real bc2i = dphidxi[i][qp]*bc(0,1) + dphideta[i][qp]*bc(1,1);
 
  598               MyVector2d V1i(-Q.col(0).dot(Qnode[i].col(1)),
 
  599                                   Q.col(0).dot(Qnode[i].col(0)));
 
  601               MyVector2d V2i(-Q.col(1).dot(Qnode[i].col(1)),
 
  602                                   Q.col(1).dot(Qnode[i].col(0)));
 
  605               B1I = MyMatrixXd::Zero(3,5);
 
  606               B1I.block<1,3>(0,0) = bc1i*Q.col(0).transpose();
 
  607               B1I.block<1,3>(1,0) = bc2i*Q.col(1).transpose();
 
  608               B1I.block<1,3>(2,0) = bc2i*Q.col(0).transpose()+bc1i*Q.col(1).transpose();
 
  610               B1I.block<1,2>(0,3) = C1i*V1i.transpose();
 
  611               B1I.block<1,2>(1,3) = C2i*V2i.transpose();
 
  612               B1I.block<1,2>(2,3) = C2i*V1i.transpose()+C1i*V2i.transpose();
 
  616               B2I = MyMatrixXd::Zero(3,5);
 
  618               B2I.block<1,2>(0,3) = bc1i*V1i.transpose();
 
  619               B2I.block<1,2>(1,3) = bc2i*V2i.transpose();
 
  620               B2I.block<1,2>(2,3) = bc2i*V1i.transpose()+bc1i*V2i.transpose();
 
  624               Bc0I = MyMatrixXd::Zero(2,5);
 
  625               Bc0I.block<1,3>(0,0) = C1i*Q.col(2).transpose();
 
  626               Bc0I.block<1,3>(1,0) = C2i*Q.col(2).transpose();
 
  627               Bc0I.block<1,2>(0,3) = phi[i][qp]*V1i.transpose();
 
  628               Bc0I.block<1,2>(1,3) = phi[i][qp]*V2i.transpose();
 
  632               Bc1I = MyMatrixXd::Zero(2,5);
 
  633               Bc1I.block<1,3>(0,0) = bc1i*Q.col(2).transpose();
 
  634               Bc1I.block<1,3>(1,0) = bc2i*Q.col(2).transpose();
 
  637               MyVector2d BdxiI(dphidxi[i][qp],dphideta[i][qp]);
 
  640               for (
unsigned int j=0; j<n_var_dofs; ++j)
 
  644                   Real C1j = dphidxi[j][qp]*C0(0,0) + dphideta[j][qp]*C0(1,0);
 
  645                   Real C2j = dphidxi[j][qp]*C0(0,1) + dphideta[j][qp]*C0(1,1);
 
  648                   B0J = MyMatrixXd::Zero(3,5);
 
  649                   B0J.block<1,3>(0,0) = C1j*Q.col(0).transpose();
 
  650                   B0J.block<1,3>(1,0) = C2j*Q.col(1).transpose();
 
  651                   B0J.block<1,3>(2,0) = C2j*Q.col(0).transpose()+C1j*Q.col(1).transpose();
 
  654                   Real bc1j = dphidxi[j][qp]*bc(0,0) + dphideta[j][qp]*bc(1,0);
 
  655                   Real bc2j = dphidxi[j][qp]*bc(0,1) + dphideta[j][qp]*bc(1,1);
 
  657                   MyVector2d V1j(-Q.col(0).dot(Qnode[j].col(1)),
 
  658                                       Q.col(0).dot(Qnode[j].col(0)));
 
  660                   MyVector2d V2j(-Q.col(1).dot(Qnode[j].col(1)),
 
  661                                       Q.col(1).dot(Qnode[j].col(0)));
 
  664                   B1J = MyMatrixXd::Zero(3,5);
 
  665                   B1J.block<1,3>(0,0) = bc1j*Q.col(0).transpose();
 
  666                   B1J.block<1,3>(1,0) = bc2j*Q.col(1).transpose();
 
  667                   B1J.block<1,3>(2,0) = bc2j*Q.col(0).transpose()+bc1j*Q.col(1).transpose();
 
  669                   B1J.block<1,2>(0,3) = C1j*V1j.transpose();
 
  670                   B1J.block<1,2>(1,3) = C2j*V2j.transpose();
 
  671                   B1J.block<1,2>(2,3) = C2j*V1j.transpose()+C1j*V2j.transpose();
 
  675                   B2J = MyMatrixXd::Zero(3,5);
 
  677                   B2J.block<1,2>(0,3) = bc1j*V1j.transpose();
 
  678                   B2J.block<1,2>(1,3) = bc2j*V2j.transpose();
 
  679                   B2J.block<1,2>(2,3) = bc2j*V1j.transpose()+bc1j*V2j.transpose();
 
  683                   Bc0J = MyMatrixXd::Zero(2,5);
 
  684                   Bc0J.block<1,3>(0,0) = C1j*Q.col(2).transpose();
 
  685                   Bc0J.block<1,3>(1,0) = C2j*Q.col(2).transpose();
 
  686                   Bc0J.block<1,2>(0,3) = phi[j][qp]*V1j.transpose();
 
  687                   Bc0J.block<1,2>(1,3) = phi[j][qp]*V2j.transpose();
 
  691                   Bc1J = MyMatrixXd::Zero(2,5);
 
  692                   Bc1J.block<1,3>(0,0) = bc1j*Q.col(2).transpose();
 
  693                   Bc1J.block<1,3>(1,0) = bc2j*Q.col(2).transpose();
 
  696                   MyVector2d BdxiJ(dphidxi[j][qp], dphideta[j][qp]);
 
  702                   local_KIJ = JxW[qp] * (
 
  703                                          B0I.transpose() * Hm * B0J
 
  704                                          +  B2I.transpose() * Hf * B0J
 
  705                                          +  B0I.transpose() * Hf * B2J
 
  706                                          +  B1I.transpose() * Hf * B1J
 
  707                                          +  2*H * B0I.transpose() * Hf * B1J
 
  708                                          +  2*H * B1I.transpose() * Hf * B0J
 
  709                                          +  Bc0I.transpose() * Hc0 * Bc0J
 
  710                                          +  Bc1I.transpose() * Hc1 * Bc1J
 
  711                                          +  2*H * Bc0I.transpose() * Hc1 * Bc1J
 
  712                                          +  2*H * Bc1I.transpose() * Hc1 * Bc0J
 
  717                   full_local_KIJ = MyMatrixXd::Zero(6, 6);
 
  718                   full_local_KIJ.block<5,5>(0,0)=local_KIJ;
 
  724                   full_local_KIJ(5,5) = 
Real(Hf(0,0)*JxW[qp]*BdI.transpose()*BdJ);
 
  729                   TI = MyMatrixXd::Identity(6,6);
 
  730                   TI.block<3,3>(3,3) = Qnode[i].transpose();
 
  732                   TJ = MyMatrixXd::Identity(6,6);
 
  733                   TJ.block<3,3>(3,3) = Qnode[j].transpose();
 
  734                   global_KIJ = TI.transpose()*full_local_KIJ*TJ;
 
  739                   for (
unsigned int k=0;k<6;k++)
 
  740                     for (
unsigned int l=0;l<6;l++)
 
  741                       Ke_var[k][l](i,j) += global_KIJ(k,l);
 
  747       if (distributed_load)
 
  750           for (
unsigned int shellface=0; shellface<2; shellface++)
 
  752               std::vector<boundary_id_type> bids;
 
  755               for (std::size_t k=0; k<bids.size(); k++)
 
  757                   for (
unsigned int qp=0; qp<qrule.
n_points(); ++qp)
 
  758                     for (
unsigned int i=0; i<n_var_dofs; ++i)
 
  759                       Fe_w(i) -= JxW[qp] * phi[i][qp];
 
  773   if (!distributed_load)
 
  784         if (((*node) - C).norm() < 1e-3)
 
  785           system.
rhs->
set(node->dof_number(0, 2, 0), -q/4);
 
  791 #endif // defined(LIBMESH_HAVE_EIGEN) && defined(LIBMESH_ENABLE_SECOND_DERIVATIVES)