39 #include "libmesh/libmesh.h" 
   40 #include "libmesh/replicated_mesh.h" 
   41 #include "libmesh/mesh_generation.h" 
   42 #include "libmesh/linear_implicit_system.h" 
   43 #include "libmesh/equation_systems.h" 
   44 #include "libmesh/exact_solution.h" 
   45 #include "libmesh/exodusII_io.h" 
   46 #include "libmesh/fe.h" 
   47 #include "libmesh/quadrature_gauss.h" 
   48 #include "libmesh/dof_map.h" 
   49 #include "libmesh/sparse_matrix.h" 
   50 #include "libmesh/numeric_vector.h" 
   51 #include "libmesh/dense_matrix.h" 
   52 #include "libmesh/dense_vector.h" 
   53 #include "libmesh/elem.h" 
   54 #include "libmesh/dirichlet_boundaries.h" 
   55 #include "libmesh/zero_function.h" 
   56 #include "libmesh/libmesh_logging.h" 
   57 #include "libmesh/getpot.h" 
   58 #include "libmesh/error_vector.h" 
   59 #include "libmesh/kelly_error_estimator.h" 
   60 #include "libmesh/mesh_refinement.h" 
   61 #include "libmesh/enum_solver_package.h" 
   68                       const std::string & system_name);
 
   72 int main (
int argc, 
char ** argv)
 
   74   START_LOG(
"Initialize and create cubes", 
"main");
 
   79                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   82   GetPot command_line (argc, argv);
 
   86     libmesh_error_msg(
"Usage:\n" << 
"\t " << argv[0] << 
" -n 15");
 
   94       for (
int i=1; i<argc; i++)
 
  101   const unsigned int dim = 3;
 
  104   libmesh_example_requires(
dim <= LIBMESH_DIM, 
"3D support");
 
  107 #ifndef LIBMESH_ENABLE_DIRICHLET 
  108   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  113   if (command_line.search(1, 
"-n"))
 
  114     ps = command_line.next(ps);
 
  125   MeshTools::Generation::build_cube (
mesh, ps, ps, ps, -1,    0,    0,  1,  0, 1, 
HEX8);
 
  126   MeshTools::Generation::build_cube (mesh1, ps, ps, ps,    0,  1,    0,  1,  0, 1, 
HEX8);
 
  127   MeshTools::Generation::build_cube (mesh2, ps, ps, ps, -1,    0, -1,    0,  0, 1, 
HEX8);
 
  128   MeshTools::Generation::build_cube (mesh3, ps, ps, ps,    0,  1, -1,    0,  0, 1, 
HEX8);
 
  129   MeshTools::Generation::build_cube (mesh4, ps, ps, ps, -1,    0,    0,  1, -1, 0, 
HEX8);
 
  130   MeshTools::Generation::build_cube (mesh5, ps, ps, ps,    0,  1,    0,  1, -1, 0, 
HEX8);
 
  131   MeshTools::Generation::build_cube (mesh6, ps, ps, ps, -1,    0, -1,    0, -1, 0, 
HEX8);
 
  132   MeshTools::Generation::build_cube (mesh7, ps, ps, ps,    0,  1, -1,    0, -1, 0, 
HEX8);
 
  136   MeshTools::Generation::build_cube (nostitch_mesh, ps*2, ps*2, ps*2, -1, 1, -1, 1, -1, 1, 
HEX8);
 
  137   STOP_LOG(
"Initialize and create cubes", 
"main");
 
  139   START_LOG(
"Stitching", 
"main");
 
  141   mesh.stitch_meshes(mesh1, 2, 4, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  142   mesh2.stitch_meshes(mesh3, 2, 4, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  143   mesh.stitch_meshes(mesh2, 1, 3, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  144   mesh4.stitch_meshes(mesh5, 2, 4, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  145   mesh6.stitch_meshes(mesh7, 2, 4, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  146   mesh4.stitch_meshes(mesh6, 1, 3, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  147   mesh.stitch_meshes(mesh4, 0, 5, 
TOLERANCE, 
true, 
true, 
false, 
false);
 
  148   STOP_LOG(
"Stitching", 
"main");
 
  150   START_LOG(
"Initialize and solve systems", 
"main");
 
  156   STOP_LOG(
"Initialize and solve systems", 
"main");
 
  158   START_LOG(
"Result comparison", 
"main");
 
  164   libMesh::out << 
"L2 error between stitched and non-stitched cases: " << error << std::endl;
 
  165   STOP_LOG(
"Result comparison", 
"main");
 
  167   START_LOG(
"Output", 
"main");
 
  168 #ifdef LIBMESH_HAVE_EXODUS_API 
  170                                            equation_systems_stitch);
 
  172                                                     equation_systems_nostitch);
 
  173 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  174   STOP_LOG(
"Output", 
"main");
 
  187 #ifdef LIBMESH_ENABLE_DIRICHLET 
  193   std::set<boundary_id_type> boundary_ids;
 
  194   for (
int j = 0; j<6; ++j)
 
  195     boundary_ids.insert(j);
 
  198   std::vector<unsigned int> variables(1);
 
  199   variables[0] = u_var;
 
  208 #endif // LIBMESH_ENABLE_DIRICHLET 
  210   equation_systems.
init();
 
  213 #ifdef LIBMESH_ENABLE_AMR 
  220   const unsigned int max_r_steps = 2;
 
  222   for (
unsigned int r_step=0; r_step<=max_r_steps; r_step++)
 
  225       if (r_step != max_r_steps)
 
  242           equation_systems.
reinit();
 
  251                       const std::string & libmesh_dbg_var(system_name))
 
  253   libmesh_assert_equal_to (system_name, 
"Poisson");
 
  261   FEType fe_type = dof_map.variable_type(0);
 
  264   fe->attach_quadrature_rule (&qrule);
 
  267   fe_face->attach_quadrature_rule (&qface);
 
  269   const std::vector<Real> & JxW = fe->get_JxW();
 
  270   const std::vector<std::vector<Real>> & phi = fe->get_phi();
 
  271   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  276   std::vector<dof_id_type> dof_indices;
 
  280       dof_map.dof_indices (elem, dof_indices);
 
  284       Ke.
resize (dof_indices.size(),
 
  287       Fe.
resize (dof_indices.size());
 
  289       for (
unsigned int qp=0; qp<qrule.
n_points(); qp++)
 
  291           for (std::size_t i=0; i<phi.size(); i++)
 
  293               Fe(i) += JxW[qp]*phi[i][qp];
 
  294               for (std::size_t j=0; j<phi.size(); j++)
 
  295                 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
 
  299       dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);