libMesh
mesh_input.C
Go to the documentation of this file.
1 #include <libmesh/distributed_mesh.h>
2 #include <libmesh/dof_map.h>
3 #include <libmesh/equation_systems.h>
4 #include <libmesh/linear_implicit_system.h>
5 #include <libmesh/mesh.h>
6 #include <libmesh/mesh_communication.h>
7 #include <libmesh/mesh_generation.h>
8 #include <libmesh/numeric_vector.h>
9 #include <libmesh/replicated_mesh.h>
10 #include <libmesh/enum_norm_type.h>
11 #include <libmesh/enum_to_string.h>
12 
13 #include <libmesh/abaqus_io.h>
14 #include <libmesh/dyna_io.h>
15 #include <libmesh/exodusII_io.h>
16 #include <libmesh/gmsh_io.h>
17 #include <libmesh/nemesis_io.h>
18 #include <libmesh/stl_io.h>
19 #include <libmesh/vtk_io.h>
20 #include <libmesh/tetgen_io.h>
21 
22 #include "test_comm.h"
23 #include "libmesh_cppunit.h"
24 
25 #include <regex>
26 
27 using namespace libMesh;
28 
29 
31  const Parameters&,
32  const std::string&,
33  const std::string&)
34 {
35  const Real & x = p(0);
36  const Real & y = p(1);
37 
38  return 6*x + 60*y;
39 }
40 
41 
43  const Parameters&,
44  const std::string&,
45  const std::string&)
46 {
47  const Real & x = p(0);
48  const Real & y = p(1);
49 
50  return sin(x) + cos(y);
51 }
52 
53 constexpr int added_sides_nxyz[] = {2,2,2};
54 
56  const Parameters & param,
57  const std::string &,
58  const std::string &)
59 {
60  const Real & x = p(0);
61  const Real & y = p(1);
62  const Real & z = p(2);
63 
64  short facedim = param.have_parameter<short>("face") ?
65  param.get<short>("face") : -1;
66 
67  // What face are we on?
68  auto is_on_face = [facedim](Real r, short rdim) {
69  if (facedim == rdim)
70  return true;
71  if (facedim >= 0)
72  return false;
73  const Real numerator = r * added_sides_nxyz[rdim];
74  return (std::abs(numerator - std::round(numerator)) <
76  };
77 
78 
79 
80  // x/y/z components of a div-free flux,
81  // curl([x^2yz, xy^2z, xyz])
82  if (is_on_face(x, 0))
83  {
84  libmesh_assert(!is_on_face(y, 1));
85  libmesh_assert(!is_on_face(z, 2));
86  return (x*z-x*y*y);
87  }
88  if (is_on_face(y, 1))
89  {
90  libmesh_assert(!is_on_face(z, 2));
91  return (x*x*y-y*z);
92  }
93 
94  libmesh_assert(is_on_face(z, 2));
95  return (y*y*z-x*x*z);
96 }
97 
98 
99 class MeshInputTest : public CppUnit::TestCase {
100 public:
101  LIBMESH_CPPUNIT_TEST_SUITE( MeshInputTest );
102 
103 #if LIBMESH_DIM > 1
104 #ifdef LIBMESH_HAVE_VTK
105  CPPUNIT_TEST( testVTKPreserveElemIds );
106  CPPUNIT_TEST( testVTKPreserveSubdomainIds );
107 #endif
108 
109 #ifdef LIBMESH_HAVE_EXODUS_API
110  CPPUNIT_TEST( testExodusCopyNodalSolutionDistributed );
111  CPPUNIT_TEST( testExodusCopyElementSolutionDistributed );
112  CPPUNIT_TEST( testExodusCopyNodalSolutionReplicated );
113  CPPUNIT_TEST( testExodusCopyElementSolutionReplicated );
114  CPPUNIT_TEST( testExodusReadHeader );
115 #if LIBMESH_DIM > 2
116  CPPUNIT_TEST( testExodusIGASidesets );
117  CPPUNIT_TEST( testLowOrderEdgeBlocks );
118 #endif
119 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
120  CPPUNIT_TEST( testExodusCopyElementVectorDistributed );
121  CPPUNIT_TEST( testExodusCopyElementVectorReplicated );
122 
123  // Eventually this will support complex numbers.
124  CPPUNIT_TEST( testExodusWriteElementDataFromDiscontinuousNodalData );
125 #endif // !LIBMESH_USE_COMPLEX_NUMBERS
126 
127  CPPUNIT_TEST( testExodusWriteAddedSidesEdgeC0 );
128  CPPUNIT_TEST( testExodusDiscWriteAddedSidesEdgeC0 );
129  CPPUNIT_TEST( testExodusWriteAddedSidesMixedEdgeC0 );
130  CPPUNIT_TEST( testExodusDiscWriteAddedSidesMixedEdgeC0 );
131  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesEdgeDisc ); // need is_on_face fixes
132  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesEdgeDisc ); // need is_on_face fixes
133  CPPUNIT_TEST( testExodusWriteAddedSidesTriC0 );
134  CPPUNIT_TEST( testExodusDiscWriteAddedSidesTriC0 );
135  CPPUNIT_TEST( testExodusWriteAddedSidesMixedTriC0 );
136  CPPUNIT_TEST( testExodusDiscWriteAddedSidesMixedTriC0 );
137  // CPPUNIT_TEST( testExodusWriteAddedSidesTriDisc ); // Need aligned faces
138  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesTriDisc ); // Need aligned faces
139  CPPUNIT_TEST( testExodusWriteAddedSidesQuadC0 );
140  CPPUNIT_TEST( testExodusDiscWriteAddedSidesQuadC0 );
141  CPPUNIT_TEST( testExodusWriteAddedSidesMixedQuadC0 );
142  CPPUNIT_TEST( testExodusDiscWriteAddedSidesMixedQuadC0 );
143  // CPPUNIT_TEST( testExodusWriteAddedSidesQuadDisc ); // need is_on_face fixes
144  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesQuadDisc ); // need is_on_face fixes
145  // CPPUNIT_TEST( testExodusWriteAddedSidesTetC0 ); // BROKEN!?! WHY!?!
146  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesTetC0 ); // BROKEN!?! WHY!?!
147  // CPPUNIT_TEST( testExodusWriteAddedSidesTetDisc );
148  // CPPUNIT_TEST( testExodusDiscWriteAddedSidesTetDisc );
149  CPPUNIT_TEST( testExodusWriteAddedSidesHexC0 );
150  CPPUNIT_TEST( testExodusDiscWriteAddedSidesHexC0 );
151  CPPUNIT_TEST( testExodusWriteAddedSidesMixedHexC0 );
152  CPPUNIT_TEST( testExodusDiscWriteAddedSidesMixedHexC0 );
153  CPPUNIT_TEST( testExodusWriteAddedSidesHexDisc );
154  CPPUNIT_TEST( testExodusDiscWriteAddedSidesHexDisc );
155 
156  CPPUNIT_TEST( testExodusFileMappingsPlateWithHole);
157  CPPUNIT_TEST( testExodusFileMappingsTwoBlocks);
158  CPPUNIT_TEST( testExodusFileMappingsTwoElemIGA);
159  CPPUNIT_TEST( testExodusFileMappingsCyl3d);
160 
161  CPPUNIT_TEST( testExodusDiscPlateWithHole);
162  CPPUNIT_TEST( testExodusDiscTwoBlocks);
163  CPPUNIT_TEST( testExodusDiscTwoElemIGA);
164  CPPUNIT_TEST( testExodusDiscCyl3d);
165 #endif // LIBMESH_HAVE_EXODUS_API
166 
167 #if defined(LIBMESH_HAVE_EXODUS_API) && defined(LIBMESH_HAVE_NEMESIS_API)
168  CPPUNIT_TEST( testNemesisReadReplicated );
169  CPPUNIT_TEST( testNemesisReadDistributed );
170 
171  CPPUNIT_TEST( testNemesisCopyNodalSolutionDistributed );
172  CPPUNIT_TEST( testNemesisCopyNodalSolutionReplicated );
173  CPPUNIT_TEST( testNemesisCopyElementSolutionDistributed );
174  CPPUNIT_TEST( testNemesisCopyElementSolutionReplicated );
175 
176  CPPUNIT_TEST( testNemesisSingleElementDistributed );
177  CPPUNIT_TEST( testNemesisSingleElementReplicated );
178 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
179  CPPUNIT_TEST( testNemesisCopyElementVectorDistributed );
180  CPPUNIT_TEST( testNemesisCopyElementVectorReplicated );
181 #endif // !LIBMESH_USE_COMPLEX_NUMBERS
182 #endif // defined(LIBMESH_HAVE_EXODUS_API) && defined(LIBMESH_HAVE_NEMESIS_API)
183 
184 #ifdef LIBMESH_HAVE_GZSTREAM
185  CPPUNIT_TEST( testAbaqusReadFirst );
186  CPPUNIT_TEST( testAbaqusReadSecond );
187  CPPUNIT_TEST( testDynaReadElem );
188  CPPUNIT_TEST( testDynaNoSplines );
189  CPPUNIT_TEST( testDynaReadPatch );
190  CPPUNIT_TEST( testDynaFileMappingsFEMEx5);
191  CPPUNIT_TEST( testDynaFileMappingsBlockWithHole);
192  CPPUNIT_TEST( testDynaFileMappingsPlateWithHole);
193  CPPUNIT_TEST( testDynaFileMappingsCyl3d);
194 #endif // LIBMESH_HAVE_GZSTREAM
195 #endif // LIBMESH_DIM > 1
196  //
197 #if LIBMESH_DIM > 1
198  CPPUNIT_TEST( testBadGmsh );
199  CPPUNIT_TEST( testGoodGmsh );
200 #endif
201 
202 #if LIBMESH_DIM > 2
203  CPPUNIT_TEST( testGoodSTL );
204  CPPUNIT_TEST( testGoodSTLBinary );
205 
206  CPPUNIT_TEST( testGmshBCIDOverlap );
207 
208 #ifdef LIBMESH_HAVE_TETGEN
209  CPPUNIT_TEST( testTetgenIO );
210 #endif
211 #endif
212 
213  CPPUNIT_TEST_SUITE_END();
214 
215 private:
216 
217 public:
218  void setUp()
219  {}
220 
221  void tearDown()
222  {}
223 
224 #ifdef LIBMESH_HAVE_VTK
226  {
227  LOG_UNIT_TEST;
228 
229  // Come up with some crazy numbering. Make all the new ids higher
230  // than the existing ids so we don't have to worry about conflicts
231  // while renumbering.
232  dof_id_type start_id;
233 
234  // first scope: write file
235  {
237  mesh.allow_renumbering(false);
238  MeshTools::Generation::build_square (mesh, 3, 3, 0., 1., 0., 1.);
239 
240  start_id = mesh.max_elem_id();
241 
242  // Use a separate container that won't invalidate iterators when
243  // we renumber
244  std::set<Elem *> elements {mesh.elements_begin(), mesh.elements_end()};
245  for (Elem * elem : elements)
246  {
247  const Point center = elem->vertex_average();
248  const int xn = int(center(0)*3);
249  const int yn = int(center(1)*3);
250  const dof_id_type new_id = start_id + yn*5 + xn;
251  mesh.renumber_elem(elem->id(), new_id);
252  }
253 
254  // Explicit writer object here to be absolutely sure we get VTK
255  VTKIO vtk(mesh);
256  vtk.write("read_elem_ids_test.pvtu");
257  }
258 
259  // Make sure that the writing is done before the reading starts.
261 
262  // second scope: read file
263  {
265  mesh.allow_renumbering(false);
266 
267  mesh.read("read_elem_ids_test.pvtu");
269 
270  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), dof_id_type(16));
271  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), dof_id_type(9));
272 
273  for (const auto & elem : mesh.element_ptr_range())
274  {
275  const Point center = elem->vertex_average();
276  const int xn = int(center(0)*3);
277  const int yn = int(center(1)*3);
278  const dof_id_type expected_id = start_id + yn*5 + xn;
279  CPPUNIT_ASSERT_EQUAL(elem->id(), expected_id);
280  }
281  }
282  }
283 
285  {
286  LOG_UNIT_TEST;
287 
288  // first scope: write file
289  {
291  mesh.allow_renumbering(false);
292  MeshTools::Generation::build_square (mesh, 3, 3, 0., 1., 0., 1.);
293 
294  for (const auto & elem : mesh.element_ptr_range())
295  {
296  const Point center = elem->vertex_average();
297  const int xn = int(center(0)*3);
298  const int yn = int(center(1)*3);
299  const subdomain_id_type new_id = yn*4 + xn;
300  elem->subdomain_id() = new_id;
301  }
302 
303  // Explicit writer object here to be absolutely sure we get VTK
304  VTKIO vtk(mesh);
305  vtk.write("read_sbd_ids_test.pvtu");
306  }
307 
308  // Make sure that the writing is done before the reading starts.
310 
311  // second scope: read file
312  {
314  mesh.allow_renumbering(false);
315 
316  mesh.read("read_sbd_ids_test.pvtu");
318 
319  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), dof_id_type(16));
320  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), dof_id_type(9));
321 
322  for (const auto & elem : mesh.element_ptr_range())
323  {
324  const Point center = elem->vertex_average();
325  const int xn = int(center(0)*3);
326  const int yn = int(center(1)*3);
327  const subdomain_id_type expected_id = yn*4 + xn;
328  CPPUNIT_ASSERT_EQUAL(elem->subdomain_id(), expected_id);
329  }
330  }
331  }
332 #endif // LIBMESH_HAVE_VTK
333 
334 
335 #ifdef LIBMESH_HAVE_EXODUS_API
337  {
338  LOG_UNIT_TEST;
339 
340  // first scope: write file
341  {
343  MeshTools::Generation::build_square (mesh, 3, 3, 0., 1., 0., 1.);
344  ExodusII_IO exii(mesh);
345  mesh.write("read_header_test.e");
346  }
347 
348  // Make sure that the writing is done before the reading starts.
350 
351  // second scope: read header
352  // Note: The header information is read from file on processor 0
353  // and then broadcast to the other procs, so with this test we are
354  // checking both that the header information is read correctly and
355  // that it is correctly communicated to other procs.
356  {
358  ExodusII_IO exii(mesh);
359  ExodusHeaderInfo header_info = exii.read_header("read_header_test.e");
360 
361  // Make sure the header information is as expected.
362  CPPUNIT_ASSERT_EQUAL(std::string(header_info.title.data()), std::string("read_header_test.e"));
363  CPPUNIT_ASSERT_EQUAL(header_info.num_dim, 2);
364  CPPUNIT_ASSERT_EQUAL(header_info.num_elem, 9);
365  CPPUNIT_ASSERT_EQUAL(header_info.num_elem_blk, 1);
366  CPPUNIT_ASSERT_EQUAL(header_info.num_node_sets, 4);
367  CPPUNIT_ASSERT_EQUAL(header_info.num_side_sets, 4);
368  CPPUNIT_ASSERT_EQUAL(header_info.num_edge_blk, 0);
369  CPPUNIT_ASSERT_EQUAL(header_info.num_edge, 0);
370  }
371  }
372 
374  {
375  LOG_UNIT_TEST;
376 
378  ExodusII_IO exii(mesh);
379 
380  if (mesh.processor_id() == 0)
381  exii.read("meshes/mesh_with_low_order_edge_blocks.e");
382 
385 
386  // Check that we see the boundary ids we expect
388 
389  // On a ReplicatedMesh, check that the number of edge boundary
390  // conditions is as expected. The real test is that we can read
391  // this file in at all. Prior to the changes in #3491, the Exodus
392  // reader threw an exception while trying to read this mesh.
393  if (mesh.is_serial())
394  {
395  // Mesh has 26 boundary ids total (including edge and side ids).
396  // ss_prop1 = 200, 201 ;
397  // ed_prop1 = 8000, 8001, 8002, 8003, 8004, 8005, 8006, 8007, 8008, 8009, 8010,
398  // 8011, 9001, 9002, 9003, 9004, 9005, 9006, 9007, 9008, 9009, 9010,
399  // 9011, 9012 ;
400  CPPUNIT_ASSERT_EQUAL(static_cast<std::size_t>(26), bi.n_boundary_ids());
401 
402  // We can binary_search() the build_edge_list() which is sorted
403  // in lexicographical order before it's returned.
404  auto edge_list = bi.build_edge_list();
405 
406  // Search for some tuples we expect to be present
407  CPPUNIT_ASSERT(std::binary_search(edge_list.begin(), edge_list.end(), std::make_tuple(4, 1, 8007)));
408  CPPUNIT_ASSERT(std::binary_search(edge_list.begin(), edge_list.end(), std::make_tuple(10, 6, 8001)));
409 
410  // And make sure we don't have entries we shouldn't have
411  CPPUNIT_ASSERT(!std::binary_search(edge_list.begin(), edge_list.end(), std::make_tuple(1, 8, 8009)));
412  CPPUNIT_ASSERT(!std::binary_search(edge_list.begin(), edge_list.end(), std::make_tuple(2, 10, 9011)));
413  }
414  }
415 
417  {
418  LOG_UNIT_TEST;
419 
421 
422  // Block here so we trigger exii destructor early; I thought I
423  // might have had a bug in there at one point
424  {
425  ExodusII_IO exii(mesh);
426  // IGA Exodus meshes require ExodusII 8 or higher
427  if (exii.get_exodus_version() < 800)
428  return;
429 
430  if (mesh.processor_id() == 0)
431  exii.read("meshes/Cube_With_Sidesets.e");
432 
433  }
434 
437 
438  // 5^3 spline nodes + 7^3 Rational Bezier nodes
439  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(468));
440  // 5^3 spline elements + 3^3 Rational Bezier elements
441  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(152));
442 
443  // Check that we see the boundary ids we expect
445 
446  // On a ReplicatedMesh, we should see all 6 boundary ids on each processor
447  if (mesh.is_serial())
448  CPPUNIT_ASSERT_EQUAL(static_cast<std::size_t>(6), bi.n_boundary_ids());
449 
450  // On any mesh, we should see each id on *some* processor
451  {
452  const std::set<boundary_id_type> & bc_ids = bi.get_boundary_ids();
453  // CoreForm gave me a file with 1-based numbering! (faints)
454  for (boundary_id_type i = 1 ; i != 7; ++i)
455  {
456  bool has_bcid = bc_ids.count(i);
457  mesh.comm().max(has_bcid);
458  CPPUNIT_ASSERT(has_bcid);
459  }
460  }
461 
462  // Indexed by bcid-1, because we count from 0, like God and
463  // Dijkstra intended!
464  std::vector<int> side_counts(6, 0);
465 
466  // Map from our side numbers to the file's BCIDs
467  const boundary_id_type libmesh_side_to_bcid[] = {1, 4, 6, 3, 5, 2};
468 
469  for (const auto & elem : mesh.active_local_element_ptr_range())
470  {
471  if (elem->type() == NODEELEM)
472  continue;
473 
474  for (unsigned short side=0; side<elem->n_sides(); side++)
475  {
476  if (elem->neighbor_ptr(side))
477  CPPUNIT_ASSERT_EQUAL(bi.n_boundary_ids(elem, side), 0u);
478  else
479  {
480  CPPUNIT_ASSERT_EQUAL(bi.n_boundary_ids(elem, side), 1u);
481  std::vector<boundary_id_type> bids;
482  bi.boundary_ids(elem, side, bids);
483  side_counts[bids[0]-1]++;
484  CPPUNIT_ASSERT_EQUAL(libmesh_side_to_bcid[side], bids[0]);
485  }
486  }
487  }
488 
489  for (auto bc_count : side_counts)
490  {
491  // We should have 3^2 sides with each id
492  mesh.comm().sum(bc_count);
493  CPPUNIT_ASSERT_EQUAL(bc_count, 9);
494  }
495 
496  // Test a write when we're done reading; I was getting weirdness
497  // from NetCDF at this point in a Moose output test.
498  {
499  ExodusII_IO exii(mesh);
500 
501  exii.write("Cube_With_Sidesets_out.e");
502  }
503  }
504 
505 
506  template <typename MeshType, typename IOType>
507  void testCopyNodalSolutionImpl (const std::string & filename)
508  {
509  {
510  MeshType mesh(*TestCommWorld);
511 
512  EquationSystems es(mesh);
513  System &sys = es.add_system<System> ("SimpleSystem");
514  sys.add_variable("n", FIRST, LAGRANGE);
515 
517  3, 3,
518  0., 1., 0., 1.);
519 
520  es.init();
522 
523  IOType meshoutput(mesh);
524 
525  meshoutput.write_equation_systems(filename, es);
526  }
527 
528  {
529  MeshType mesh(*TestCommWorld);
530  IOType meshinput(mesh);
531 
532  // Avoid getting Nemesis solution values mixed up
533  if (meshinput.is_parallel_format())
534  {
535  mesh.allow_renumbering(false);
537  }
538 
539  EquationSystems es(mesh);
540  System &sys = es.add_system<System> ("SimpleSystem");
541  sys.add_variable("testn", FIRST, LAGRANGE);
542 
543  if (mesh.processor_id() == 0 || meshinput.is_parallel_format())
544  meshinput.read(filename);
545  if (!meshinput.is_parallel_format())
548 
549  es.init();
550 
551  // Read the solution e into variable teste.
552  //
553  // With complex numbers, we'll only bother reading the real
554  // part.
555 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
556  meshinput.copy_nodal_solution(sys, "testn", "r_n");
557 #else
558  meshinput.copy_nodal_solution(sys, "testn", "n");
559 #endif
560 
561  // Exodus only handles double precision
562  Real exotol = std::max(TOLERANCE*TOLERANCE, Real(1e-12));
563 
564  for (Real x = 0; x < 1 + TOLERANCE; x += Real(1.L/3.L))
565  for (Real y = 0; y < 1 + TOLERANCE; y += Real(1.L/3.L))
566  {
567  Point p(x,y);
568  LIBMESH_ASSERT_NUMBERS_EQUAL
569  (sys.point_value(0,p), 6*x+60*y, exotol);
570  }
571  }
572  }
573 
574 
576  { LOG_UNIT_TEST; testCopyNodalSolutionImpl<ReplicatedMesh,ExodusII_IO>("repl_with_nodal_soln.e"); }
577 
579  { LOG_UNIT_TEST; testCopyNodalSolutionImpl<DistributedMesh,ExodusII_IO>("dist_with_nodal_soln.e"); }
580 
581 #if defined(LIBMESH_HAVE_NEMESIS_API)
583  { LOG_UNIT_TEST; testCopyNodalSolutionImpl<ReplicatedMesh,Nemesis_IO>("repl_with_nodal_soln.nem"); }
584 
586  { LOG_UNIT_TEST; testCopyNodalSolutionImpl<DistributedMesh,Nemesis_IO>("dist_with_nodal_soln.nem"); }
587 #endif
588 
589 
590  template <typename MeshType, typename IOType>
591  void testCopyElementSolutionImpl (const std::string & filename)
592  {
593  {
594  MeshType mesh(*TestCommWorld);
595 
596  EquationSystems es(mesh);
597  System &sys = es.add_system<System> ("SimpleSystem");
598  sys.add_variable("e", CONSTANT, MONOMIAL);
599 
601  3, 3,
602  0., 1., 0., 1.);
603 
604  es.init();
606 
607  IOType meshinput(mesh);
608 
609  // Don't try to write element data as nodal data
610  std::set<std::string> sys_list;
611  meshinput.write_equation_systems(filename, es, &sys_list);
612 
613  // Just write it as element data
614  meshinput.write_element_data(es);
615  }
616 
617  {
618  MeshType mesh(*TestCommWorld);
619  IOType meshinput(mesh);
620 
621  // Avoid getting Nemesis solution values mixed up
622  if (meshinput.is_parallel_format())
623  {
624  mesh.allow_renumbering(false);
626  }
627 
628  EquationSystems es(mesh);
629  System &sys = es.add_system<System> ("SimpleSystem");
630  sys.add_variable("teste", CONSTANT, MONOMIAL);
631 
632  if (mesh.processor_id() == 0 || meshinput.is_parallel_format())
633  meshinput.read(filename);
634  if (!meshinput.is_parallel_format())
637 
638  es.init();
639 
640  // Read the solution e into variable teste.
641  //
642  // With complex numbers, we'll only bother reading the real
643  // part.
644 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
645  meshinput.copy_elemental_solution(sys, "teste", "r_e");
646 #else
647  meshinput.copy_elemental_solution(sys, "teste", "e");
648 #endif
649 
650  // Exodus only handles double precision
651  Real exotol = std::max(TOLERANCE*TOLERANCE, Real(1e-12));
652 
653  for (Real x = Real(1.L/6.L); x < 1; x += Real(1.L/3.L))
654  for (Real y = Real(1.L/6.L); y < 1; y += Real(1.L/3.L))
655  {
656  Point p(x,y);
657  LIBMESH_ASSERT_NUMBERS_EQUAL
658  (sys.point_value(0,p), 6*x+60*y, exotol);
659  }
660  }
661  }
662 
663 
665  { LOG_UNIT_TEST; testCopyElementSolutionImpl<ReplicatedMesh,ExodusII_IO>("repl_with_elem_soln.e"); }
666 
668  { LOG_UNIT_TEST; testCopyElementSolutionImpl<DistributedMesh,ExodusII_IO>("dist_with_elem_soln.e"); }
669 
670 #if defined(LIBMESH_HAVE_NEMESIS_API)
672  { LOG_UNIT_TEST; testCopyElementSolutionImpl<ReplicatedMesh,Nemesis_IO>("repl_with_elem_soln.nem"); }
673 
675  { LOG_UNIT_TEST; testCopyElementSolutionImpl<DistributedMesh,Nemesis_IO>("dist_with_elem_soln.nem"); }
676 
677 
678  // This tests that a single-element mesh solution makes it through all of the API to write out a
679  // set of Nemesis files. When this executes in parallel, there are ((number of processors) - 1)
680  // subdomains with zero elements and therefore nothing to write, but the API should handle this.
681  template <typename MeshType, typename IOType>
682  void testSingleElementImpl(const std::string & filename)
683  {
684  {
685  // Generate a single 1x1 square element mesh
686  MeshType mesh(*TestCommWorld);
688 
689  EquationSystems es(mesh);
690  auto & sys = es.add_system<System>("SimpleSystem");
691  sys.add_variable("e", CONSTANT, MONOMIAL);
692 
693  // Set an arbitrary solution for the single element DOF
694  es.init();
695  sys.project_solution(six_x_plus_sixty_y, nullptr, es.parameters);
696 
697  // Write the solution to Nemesis file(s) - only proc 0 should have anything to write!
698  Nemesis_IO nem_io(mesh);
699  std::set<std::string> sys_list;
700  nem_io.write_equation_systems(filename, es, &sys_list);
701  nem_io.write_element_data(es);
702  }
703 
704  // If we can read and copy the correct element value back into the mesh, we know the Nemesis
705  // file(s) were written properly.
706  {
707  MeshType mesh(*TestCommWorld);
708  EquationSystems es(mesh);
709  auto & sys = es.add_system<System>("SimpleSystem");
710  sys.add_variable("teste", CONSTANT, MONOMIAL);
711 
712  Nemesis_IO nem_io(mesh);
713  if (mesh.processor_id() == 0 || nem_io.is_parallel_format())
714  nem_io.read(filename);
715 
717  es.init();
718 
719 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
720  nem_io.copy_elemental_solution(sys, "teste", "r_e");
721 #else
722  nem_io.copy_elemental_solution(sys, "teste", "e");
723 #endif
724 
725  // The result should be '\frac{6 + 60}{2} = 33' at all points in the element domain
726  CPPUNIT_ASSERT_EQUAL(int(sys.solution->size()), 1);
727  CPPUNIT_ASSERT_EQUAL(libmesh_real(sys.point_value(0, Point(0.5, 0.5))), Real(33));
728  }
729  }
730 
731 
733  { LOG_UNIT_TEST; testSingleElementImpl<ReplicatedMesh,Nemesis_IO>("repl_with_single_elem.nem"); }
734 
736  { LOG_UNIT_TEST; testSingleElementImpl<DistributedMesh,Nemesis_IO>("dist_with_single_elem.nem"); }
737 #endif //defined(LIBMESH_HAVE_NEMESIS_API)
738 
739 
740 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
741  // So this tester runs through pretty much the same process as 'testCopyElementSolutionImpl()'
742  // except for a CONSTANT MONOMIAL_VEC variable. It mainly serves as a test for writing elemental
743  // vector variables to elements in ExodusII files, and is not actually a test for reading and
744  // copying an elemental solution.
745  template <typename MeshType, typename IOType>
746  void testCopyElementVectorImpl(const std::string & filename)
747  {
748  {
749  MeshType mesh(*TestCommWorld);
750 
751  EquationSystems es(mesh);
752  System & sys = es.add_system<System> ("SimpleSystem");
753  auto e_var = sys.add_variable("e", CONSTANT, MONOMIAL_VEC);
754 
756  3, 3,
757  0., 1., 0., 1.);
758 
759  es.init();
760 
761  // Here, we're going to manually set up the solution because the 'project_solution()' and
762  // 'project_vector()' methods don't work so well with CONSTANT MONOMIAL_VEC variables. They
763  // each lead to an error downstream asserting positive-definiteness when Cholesky decomposing.
764  // Interestingly, the error is only invoked for CONSTANT MONOMIAL_VEC, and not, e.g.,
765  // CONSTANT MONOMIAL nor FIRST LAGRANGE_VEC.
766  //
767  // Anyways, the important thing here is that we test the ExodusII and Nemesis writers, how
768  // the solution is set is hardly important, and we're pretty much following the same
769  // philosophy as the 'test2DProjectVectorFE()' unit tester in 'systems_test.C'
770  Parameters params;
771  for (const auto & elem : mesh.active_local_element_ptr_range())
772  {
773  const Point & p = elem->vertex_average();
774 
775  // Set the x-component with the value from 'six_x_plus_sixty_y()' and the y-component
776  // with that from 'sin_x_plus_cos_y()' at the element centroid (vertex average)
777  sys.current_local_solution->set(
778  elem->dof_number(sys.number(), e_var, 0), six_x_plus_sixty_y(p, params, "", ""));
779  sys.current_local_solution->set(
780  elem->dof_number(sys.number(), e_var, 1), sin_x_plus_cos_y(p, params, "", ""));
781  }
782 
783  // After setting values, we need to assemble
784  sys.current_local_solution->close();
785 
786  IOType meshinput(mesh);
787 
788  // Don't try to write element data as nodal data
789  std::set<std::string> sys_list;
790  meshinput.write_equation_systems(filename, es, &sys_list);
791 
792  // Just write it as element data
793  meshinput.write_element_data(es);
794  }
795 
796  {
797  MeshType mesh(*TestCommWorld);
798  IOType meshinput(mesh);
799 
800  // Avoid getting Nemesis solution values mixed up
801  if (meshinput.is_parallel_format())
802  {
803  mesh.allow_renumbering(false);
805  }
806 
807  EquationSystems es(mesh);
808  System & sys = es.add_system<System> ("SimpleSystem");
809 
810  // We have to read the CONSTANT MONOMIAL_VEC var "e" into separate CONSTANT MONOMIAL vars
811  // "e_x" and "e_y" because 'copy_elemental_solution()' currently doesn't support vectors.
812  // Again, this isn't a test for reading/copying an elemental vector solution, only writing.
813  sys.add_variable("teste_x", CONSTANT, MONOMIAL);
814  sys.add_variable("teste_y", CONSTANT, MONOMIAL);
815 
816  if (mesh.processor_id() == 0 || meshinput.is_parallel_format())
817  meshinput.read(filename);
818  if (!meshinput.is_parallel_format())
821 
822  es.init();
823 
824  // Read the solution e_x and e_y into variable teste_x and teste_y, respectively.
825  meshinput.copy_elemental_solution(sys, "teste_x", "e_x");
826  meshinput.copy_elemental_solution(sys, "teste_y", "e_y");
827 
828  // Exodus only handles double precision
829  Real exotol = std::max(TOLERANCE*TOLERANCE, Real(1e-12));
830 
831  for (Real x = Real(1.L/6.L); x < 1; x += Real(1.L/3.L))
832  for (Real y = Real(1.L/6.L); y < 1; y += Real(1.L/3.L))
833  {
834  Point p(x,y);
835  LIBMESH_ASSERT_NUMBERS_EQUAL
836  (sys.point_value(0,p), 6*x+60*y, exotol);
837  LIBMESH_ASSERT_NUMBERS_EQUAL
838  (sys.point_value(1,p), sin(x)+cos(y), exotol);
839  }
840  }
841  }
842 
844  { LOG_UNIT_TEST; testCopyElementVectorImpl<ReplicatedMesh, ExodusII_IO>("repl_with_elem_vec.e"); }
845 
847  { LOG_UNIT_TEST; testCopyElementVectorImpl<DistributedMesh,ExodusII_IO>("dist_with_elem_vec.e"); }
848 
849 #if defined(LIBMESH_HAVE_NEMESIS_API)
851  { LOG_UNIT_TEST; testCopyElementVectorImpl<ReplicatedMesh,Nemesis_IO>("repl_with_elem_vec.nem"); }
852 
854  { LOG_UNIT_TEST; testCopyElementVectorImpl<DistributedMesh,Nemesis_IO>("dist_with_elem_vec.nem"); }
855 #endif
856 
857 
859  {
860  LOG_UNIT_TEST;
861 
862  // first scope: write file
863  {
865 
866  EquationSystems es(mesh);
867  System & sys = es.add_system<System> ("SimpleSystem");
868  sys.add_variable("u", FIRST, L2_LAGRANGE);
869 
871  (mesh, 2, 2, 2, 0., 1., 0., 1., 0., 1., HEX8);
872 
873  es.init();
874 
875  // Set solution u^e_i = i, for the ith vertex of a given element e.
876  const DofMap & dof_map = sys.get_dof_map();
877  std::vector<dof_id_type> dof_indices;
878  for (const auto & elem : mesh.element_ptr_range())
879  {
880  dof_map.dof_indices(elem, dof_indices, /*var_id=*/0);
881  for (unsigned int i=0; i<dof_indices.size(); ++i)
882  sys.solution->set(dof_indices[i], i);
883  }
884  sys.solution->close();
885 
886  // Now write to file.
887  ExodusII_IO exii(mesh);
888 
889  // Don't try to write element data as averaged nodal data.
890  std::set<std::string> sys_list;
891  exii.write_equation_systems("elemental_from_nodal.e", es, &sys_list);
892 
893  // Write one elemental data field per vertex value.
894  sys_list = {"SimpleSystem"};
895 
897  (es, &sys_list, /*var_suffix=*/"_elem_corner_");
898  } // end first scope
899 
900  // second scope: read values back in, verify they are correct.
901  {
902  std::vector<std::string> file_var_names =
903  {"u_elem_corner_0",
904  "u_elem_corner_1",
905  "u_elem_corner_2",
906  "u_elem_corner_3"};
907  std::vector<Real> expected_values = {0., 1., 2., 3.};
908 
909  // copy_elemental_solution currently requires ReplicatedMesh
911 
912  EquationSystems es(mesh);
913  System & sys = es.add_system<System> ("SimpleSystem");
914  for (auto i : index_range(file_var_names))
915  sys.add_variable(file_var_names[i], CONSTANT, MONOMIAL);
916 
917  ExodusII_IO exii(mesh);
918 
919  if (mesh.processor_id() == 0)
920  exii.read("elemental_from_nodal.e");
923 
924  es.init();
925 
926  for (auto i : index_range(file_var_names))
928  (sys, sys.variable_name(i), file_var_names[i]);
929 
930  // Check that the values we read back in are as expected.
931  for (const auto & elem : mesh.active_element_ptr_range())
932  for (auto i : index_range(file_var_names))
933  {
934  Real read_val = sys.point_value(i, elem->vertex_average());
935  LIBMESH_ASSERT_FP_EQUAL
936  (expected_values[i], read_val, TOLERANCE*TOLERANCE);
937  }
938  } // end second scope
939  } // end testExodusWriteElementDataFromDiscontinuousNodalData
940 
941 #endif // !LIBMESH_USE_COMPLEX_NUMBERS
942 
943 
944  void testExodusWriteAddedSides
945  (Number (*exact_sol)(const Point &, const Parameters &, const
946  std::string &, const std::string &),
947  const ElemType elem_type,
948  const Order order,
949  const bool write_discontinuous = false,
950  const std::vector<FEType> earlier_vars = {},
951  const std::vector<FEType> later_vars = {})
952  {
953  constexpr unsigned int nx = added_sides_nxyz[0],
954  ny = added_sides_nxyz[1],
955  nz = added_sides_nxyz[2];
956 
957  const unsigned int dim = Elem::type_to_dim_map[elem_type];
958  const bool is_tensor = (Elem::build(elem_type)->n_sides() == dim * 2);
959 
960  // Figure out how many fake and true elements to expect
961  dof_id_type n_fake_elem = 0;
962  dof_id_type n_true_elem = 0;
963 
964  dof_id_type n_fake_nodes = 0;
965  dof_id_type n_true_nodes = 0;
966 
967  const std::string filename =
968  "side_discontinuous_"+Utility::enum_to_string<ElemType>(elem_type)+(write_discontinuous?"_disc":"")+".e";
969 
970  // first scope: write file
971  {
973 
974  EquationSystems es(mesh);
975  System & sys = es.add_system<System> ("SimpleSystem");
976  int varnum = 1;
977  for (auto vartype : earlier_vars)
978  sys.add_variable("earlier_"+std::to_string(varnum++), vartype);
979 
980  sys.add_variable("u", order, SIDE_HIERARCHIC);
981 
982  varnum = 1;
983  for (auto vartype : later_vars)
984  sys.add_variable("later_"+std::to_string(varnum++), vartype);
985 
986  if (dim == 3)
988  (mesh, nx, ny, nz, 0., 1., 0., 1., 0., 1., elem_type);
989  else if (dim == 2)
991  (mesh, nx, ny, 0., 1., 0., 1., elem_type);
992  else
994  (mesh, nx, 0., 1., elem_type);
995 
996  n_true_elem = mesh.n_elem();
997  n_true_nodes = mesh.n_nodes();
998  CPPUNIT_ASSERT_LESS(n_true_nodes, n_true_elem); // Ne < Nn
999 
1000  const unsigned int our_ny = dim>1 ? ny : 1;
1001  const unsigned int our_nz = dim>2 ? nz : 1;
1002 
1003  dof_id_type min_n_elem = nx * our_ny * our_nz;
1004  CPPUNIT_ASSERT_LESSEQUAL(n_true_elem, min_n_elem); // "backwards" API...
1005 
1006  for (const auto & elem : mesh.active_local_element_ptr_range())
1007  {
1008  for (auto s : make_range(elem->n_sides()))
1009  if (!elem->neighbor_ptr(s) || elem->neighbor_ptr(s)->id() < elem->id())
1010  {
1011  ++n_fake_elem;
1012  auto side = elem->build_side_ptr(s);
1013  n_fake_nodes += side->n_nodes();
1014  }
1015  }
1016  mesh.comm().sum(n_fake_elem);
1017  mesh.comm().sum(n_fake_nodes);
1018 
1019  const dof_id_type expected_fakes = [elem_type]() {
1020  switch (elem_type)
1021  {
1022  case EDGE3:
1023  return nx+1;
1024  case TRI6:
1025  return 3*nx*ny + nx + ny;
1026  case QUAD8:
1027  case QUAD9:
1028  return 2*nx*ny + nx + ny;
1029  case TET14:
1030  return 48*nx*ny*nz + 4*(nx*ny+nx*nz+ny*nz);
1031  case HEX27:
1032  return 3*nx*ny*nz + nx*ny + nx*nz + ny*nz;
1033  default:
1034  libmesh_error();
1035  }
1036  } (); // Invoke anonymous lambda
1037 
1038  CPPUNIT_ASSERT_EQUAL(n_fake_elem, expected_fakes); // "backwards" API...
1039 
1040  es.init();
1041  sys.project_solution(exact_sol, nullptr,
1042  es.parameters);
1043 
1044  // Set solution u^e_i = i, for the ith vertex of a given element e.
1045 
1046  // Now write to file.
1047  ExodusII_IO exii(mesh);
1048  exii.write_added_sides(true);
1049 
1050  if (write_discontinuous)
1051  exii.write_discontinuous_equation_systems(filename, es);
1052  else
1053  exii.write_equation_systems(filename, es);
1054  } // end first scope
1055 
1056  // second scope: read file, verify extra elements exist
1057  {
1059  mesh.allow_renumbering(false);
1060 
1061  ExodusII_IO exii(mesh);
1062 
1063  if (mesh.processor_id() == 0)
1064  exii.read(filename);
1067 
1068  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), n_true_elem + n_fake_elem);
1069  if (write_discontinuous)
1070  {
1071  const dof_id_type nodes_per_elem = Elem::build(elem_type)->n_nodes();
1072  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(),
1073  n_true_elem*nodes_per_elem + n_fake_nodes);
1074  }
1075  else
1076  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), n_true_nodes + n_fake_nodes);
1077 
1078  EquationSystems es(mesh);
1079  System & sys = es.add_system<System> ("SimpleSystem");
1080  // Read back into a LAGRANGE variable for testing; we still
1081  // can't use Exodus for a proper restart.
1082  sys.add_variable("ul", SECOND);
1083  es.init();
1084 
1085  const DofMap & dof_map = sys.get_dof_map();
1086 
1087 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1088  exii.copy_nodal_solution(sys, "ul", "r_u");
1089 #else
1090  exii.copy_nodal_solution(sys, "ul", "u");
1091 #endif
1092 
1093  dof_id_type n_side_nodes = 0;
1094  const std::string nullstr;
1095  const std::string facestr = "face";
1096 
1097  // Debugging this in parallel is tricky. Let's make sure that
1098  // if we have a failure on one rank we see it on all the others
1099  // and we can go on to other tests.
1100 #ifdef LIBMESH_ENABLE_EXCEPTIONS
1101  bool threw_exception = false;
1102  try
1103 #endif // LIBMESH_ENABLE_EXCEPTIONS
1104  {
1105  for (const auto & elem : mesh.active_local_element_ptr_range())
1106  {
1107  // Just look at side elements, not interiors
1108  if (elem->dim() == dim)
1109  continue;
1110 
1111  std::vector<dof_id_type> dof_indices;
1112  dof_map.dof_indices(elem, dof_indices, 0);
1113 
1114  // Find what face direction we're looking at, to
1115  // disambiguate when testing against a discontinuous
1116  // function, since we're evaluating on nodes that overlap
1117  // multiple faces
1118  const Point normal = [elem](){
1119  if (elem->dim() == 2)
1120  return Point((elem->point(1) - elem->point(0)).cross
1121  (elem->point(2) - elem->point(0)));
1122  else if (elem->dim() == 1)
1123  return Point
1124  (elem->point(1)(1)-elem->point(0)(1),
1125  elem->point(0)(0)-elem->point(1)(0));
1126  else
1127  return Point(1);
1128  } (); // Invoke anonymous lambda
1129 
1130  short faceval = -1;
1131  if (is_tensor)
1132  {
1133  if (std::abs(normal(0)) > TOLERANCE)
1134  {
1135  faceval = 0;
1136  libmesh_assert_less(std::abs(normal(1)), TOLERANCE);
1137  libmesh_assert_less(std::abs(normal(2)), TOLERANCE);
1138  }
1139  else if (std::abs(normal(1)) > TOLERANCE)
1140  {
1141  faceval = 1;
1142  libmesh_assert_less(std::abs(normal(2)), TOLERANCE);
1143  }
1144  else
1145  {
1146  faceval = 2;
1147  libmesh_assert_greater(std::abs(normal(2)), TOLERANCE);
1148  }
1149  libmesh_assert_greater_equal(faceval, 0);
1150  es.parameters.set<short>(facestr) = faceval;
1151  }
1152 
1153  for (auto i : index_range(dof_indices))
1154  {
1155  const Point node_pt = elem->point(i);
1156  const Real nodal_coef =
1157  libmesh_real((*sys.current_local_solution)(dof_indices[i]));
1158  const Real exact_val =
1159  libmesh_real(exact_sol
1160  (node_pt, es.parameters, nullstr,
1161  nullstr));
1162  LIBMESH_ASSERT_FP_EQUAL
1163  (nodal_coef, exact_val,
1164  std::max(Real(2),nodal_coef+exact_val)*
1165  TOLERANCE*std::sqrt(TOLERANCE));
1166  ++n_side_nodes;
1167  }
1168  }
1169  }
1170 #ifdef LIBMESH_ENABLE_EXCEPTIONS
1171  catch (...)
1172  {
1173  threw_exception = true;
1174  TestCommWorld->max(threw_exception);
1175  throw;
1176  }
1177  if (!threw_exception)
1178  TestCommWorld->max(threw_exception);
1179  CPPUNIT_ASSERT(!threw_exception);
1180 #endif // LIBMESH_ENABLE_EXCEPTIONS
1181 
1182  TestCommWorld->sum(n_side_nodes);
1183  CPPUNIT_ASSERT_EQUAL(n_side_nodes, n_fake_nodes);
1184  } // end second scope
1185  } // end testExodusWriteAddedSides
1186 
1188  {
1189  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, FIRST);
1190  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, SECOND);
1191  }
1192 
1194  {
1195  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, FIRST, true);
1196  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, SECOND, true);
1197  }
1198 
1200  {
1201  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, FIRST, false, {{FIRST, LAGRANGE}});
1202  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, SECOND, false, {}, {{FIRST, LAGRANGE}});
1203  }
1204 
1206  {
1207  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, FIRST, true, {{FIRST, LAGRANGE}});
1208  testExodusWriteAddedSides(six_x_plus_sixty_y, EDGE3, SECOND, true, {}, {{FIRST, LAGRANGE}});
1209  }
1210 
1212  {
1213  testExodusWriteAddedSides(designed_for_side_elems, EDGE3, SECOND);
1214  }
1215 
1217  {
1218  testExodusWriteAddedSides(designed_for_side_elems, EDGE3, SECOND, true);
1219  }
1220 
1222  {
1223  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, FIRST);
1224  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, SECOND);
1225  }
1226 
1228  {
1229  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, FIRST, true);
1230  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, SECOND, true);
1231  }
1232 
1234  {
1235  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, FIRST, false, {{SECOND, HIERARCHIC}});
1236  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, SECOND, false, {}, {{SECOND, SZABAB}});
1237  }
1238 
1240  {
1241  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, FIRST, true, {{SECOND, HIERARCHIC}});
1242  testExodusWriteAddedSides(six_x_plus_sixty_y, TRI6, SECOND, true, {}, {{SECOND, SZABAB}});
1243  }
1244 
1246  {
1247  testExodusWriteAddedSides(designed_for_side_elems, TRI6, SECOND);
1248  }
1249 
1251  {
1252  testExodusWriteAddedSides(designed_for_side_elems, TRI6, SECOND, true);
1253  }
1254 
1256  {
1257  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, FIRST);
1258  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, SECOND);
1259  }
1260 
1262  {
1263  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, FIRST, true);
1264  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, SECOND, true);
1265  }
1266 
1268  {
1269  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, FIRST, false, {{SECOND, LAGRANGE}});
1270  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, SECOND, false, {}, {{FIRST, LAGRANGE}});
1271  }
1272 
1274  {
1275  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, FIRST, true, {{SECOND, LAGRANGE}});
1276  testExodusWriteAddedSides(six_x_plus_sixty_y, QUAD9, SECOND, true, {}, {{FIRST, LAGRANGE}});
1277  }
1278 
1280  {
1281  testExodusWriteAddedSides(designed_for_side_elems, QUAD9, SECOND);
1282  }
1283 
1285  {
1286  testExodusWriteAddedSides(designed_for_side_elems, QUAD9, SECOND, true);
1287  }
1288 
1290  {
1291  testExodusWriteAddedSides(six_x_plus_sixty_y, TET14, FIRST);
1292  testExodusWriteAddedSides(six_x_plus_sixty_y, TET14, SECOND);
1293  }
1294 
1296  {
1297  testExodusWriteAddedSides(six_x_plus_sixty_y, TET14, FIRST, true);
1298  testExodusWriteAddedSides(six_x_plus_sixty_y, TET14, SECOND, true);
1299  }
1300 
1302  {
1303  testExodusWriteAddedSides(designed_for_side_elems, TET14, SECOND);
1304  }
1305 
1307  {
1308  testExodusWriteAddedSides(designed_for_side_elems, TET14, SECOND, true);
1309  }
1310 
1312  {
1313  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, FIRST);
1314  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, SECOND);
1315  }
1316 
1318  {
1319  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, FIRST, true);
1320  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, SECOND, true);
1321  }
1322 
1324  {
1325  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, FIRST, false, {{FIRST, LAGRANGE}});
1326  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, SECOND, false, {}, {{SECOND, HIERARCHIC}});
1327  }
1328 
1330  {
1331  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, FIRST, true, {{FIRST, LAGRANGE}});
1332  testExodusWriteAddedSides(six_x_plus_sixty_y, HEX27, SECOND, true, {}, {{SECOND, HIERARCHIC}});
1333  }
1334 
1336  {
1337  testExodusWriteAddedSides(designed_for_side_elems, HEX27, SECOND);
1338  }
1339 
1341  {
1342  testExodusWriteAddedSides(designed_for_side_elems, HEX27, SECOND, true);
1343  }
1344 
1345 #endif // LIBMESH_HAVE_EXODUS_API
1346 
1347 
1348 #if defined(LIBMESH_HAVE_EXODUS_API) && defined(LIBMESH_HAVE_NEMESIS_API)
1349  template <typename MeshType>
1351  {
1352  // first scope: write file
1353  {
1354  MeshType mesh(*TestCommWorld);
1355  MeshTools::Generation::build_square (mesh, 3, 3, 0., 1., 0., 1.);
1356  mesh.write("test_nemesis_read.nem");
1357  }
1358 
1359  // Make sure that the writing is done before the reading starts.
1361 
1362  // second scope: read file
1363  {
1364  MeshType mesh(*TestCommWorld);
1365  Nemesis_IO nem(mesh);
1366 
1367  nem.read("test_nemesis_read.nem");
1369  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(9));
1370  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(16));
1371  }
1372  }
1373 
1375  { LOG_UNIT_TEST; testNemesisReadImpl<ReplicatedMesh>(); }
1376 
1378  { LOG_UNIT_TEST; testNemesisReadImpl<DistributedMesh>(); }
1379 #endif
1380 
1381 
1383  {
1384  auto locator = mesh.sub_point_locator();
1385 
1386  for (auto & elem : mesh.element_ptr_range())
1387  {
1388  Point master_pt = {}; // center, for tensor product elements
1389 
1390  // But perturb it to try and trigger any mapping weirdness
1391  if (elem->dim() > 0)
1392  master_pt(0) = 0.25;
1393 
1394  if (elem->dim() > 1)
1395  master_pt(1) = -0.25;
1396 
1397  if (elem->dim() > 2)
1398  master_pt(2) = 0.75;
1399 
1400  Point physical_pt = FEMap::map(elem->dim(), elem, master_pt);
1401 
1402  Point inverse_pt = FEMap::inverse_map(elem->dim(), elem,
1403  physical_pt);
1404 
1405  CPPUNIT_ASSERT((inverse_pt-master_pt).norm() < TOLERANCE);
1406 
1407  CPPUNIT_ASSERT(elem->contains_point(physical_pt));
1408 
1409  // We only want to find elements in the same block
1410  std::set<subdomain_id_type> my_subdomain { elem->subdomain_id() };
1411 
1412  // We can *still* have overlapping NodeElem from a slit mesh
1413  // input file; better check them all
1414  std::set<const Elem * > located_elems;
1415  (*locator)(physical_pt, located_elems, &my_subdomain);
1416 
1417  CPPUNIT_ASSERT(located_elems.count(elem));
1418  }
1419  }
1420 
1421 
1422 
1424  {
1425  CPPUNIT_ASSERT_EQUAL(mesh.default_mapping_type(),
1427 
1428  unsigned char weight_index = mesh.default_mapping_data();
1429 
1430  bool found_the_quad = false;
1431 
1432  for (auto & elem : mesh.element_ptr_range())
1433  {
1434  if (elem->type() == NODEELEM)
1435  continue;
1436 
1437  CPPUNIT_ASSERT_EQUAL(elem->type(), QUAD9);
1438  found_the_quad = true;
1439 
1440  for (unsigned int n=0; n != 9; ++n)
1441  CPPUNIT_ASSERT_EQUAL
1442  (elem->node_ref(n).get_extra_datum<Real>(weight_index),
1443  Real(0.75));
1444 
1445  CPPUNIT_ASSERT_EQUAL(elem->point(0)(0), Real(0.5));
1446  CPPUNIT_ASSERT_EQUAL(elem->point(0)(1), Real(0.5));
1447  CPPUNIT_ASSERT_EQUAL(elem->point(1)(0), Real(1.5));
1448  CPPUNIT_ASSERT_EQUAL(elem->point(1)(1), Real(0.5));
1449  CPPUNIT_ASSERT_EQUAL(elem->point(2)(0), Real(1.5));
1450  CPPUNIT_ASSERT_EQUAL(elem->point(2)(1), Real(1.5));
1451  CPPUNIT_ASSERT_EQUAL(elem->point(3)(0), Real(0.5));
1452  CPPUNIT_ASSERT_EQUAL(elem->point(3)(1), Real(1.5));
1453  CPPUNIT_ASSERT(elem->has_affine_map());
1454 #if LIBMESH_DIM > 2
1455  for (unsigned int v=0; v != 4; ++v)
1456  CPPUNIT_ASSERT_EQUAL(elem->point(v)(2), Real(0));
1457 #endif
1458  }
1459 
1460  TestCommWorld->max(found_the_quad);
1461  CPPUNIT_ASSERT(found_the_quad);
1462 
1463  testMasterCenters(mesh);
1464  }
1465 
1466 
1467  void testAbaqusRead (const std::string & fname,
1470  {
1472 
1473  AbaqusIO abaqus(mesh);
1474 
1475  if (mesh.processor_id() == 0)
1476  abaqus.read(fname);
1478 
1480 
1481  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), n_elem);
1482  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), n_nodes);
1483  }
1484 
1485 
1487  {
1488  LOG_UNIT_TEST;
1489  testAbaqusRead("meshes/tensile_sample_test1.inp.gz", 728, 1166);
1490  }
1491 
1492 
1494  {
1495  LOG_UNIT_TEST;
1496  testAbaqusRead("meshes/poly_sample_test2.inp.gz", 1280, 1625);
1497  }
1498 
1499 
1501  {
1502  LOG_UNIT_TEST;
1503 
1505 
1506  DynaIO dyna(mesh);
1507 
1508  if (mesh.processor_id() == 0)
1509  dyna.read("meshes/1_quad.bxt.gz");
1511 
1513 
1514  // We have 1 QUAD9 finite element, attached via a trivial map to 9
1515  // spline Node+NodeElem objects
1516  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(10));
1517  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(18));
1518 
1519  helperTestingDynaQuad(mesh);
1520  }
1521 
1522  void testBadGmsh ()
1523  {
1524  LOG_UNIT_TEST;
1525 
1527 
1528  GmshIO gmsh_io(mesh);
1529 
1530 #ifdef LIBMESH_ENABLE_EXCEPTIONS
1531  std::string what = "";
1532  try
1533  {
1534  if (mesh.processor_id() == 0)
1535  gmsh_io.read("meshes/block.msh");
1536  }
1537  catch (libMesh::LogicError & e)
1538  {
1539  what = e.what();
1540  }
1541 
1542  TestCommWorld->broadcast(what);
1543  std::regex msg_regex("outside entity physical bounding box");
1544  CPPUNIT_ASSERT(std::regex_search(what, msg_regex));
1545 #endif
1546  }
1547 
1549  {
1550  LOG_UNIT_TEST;
1551 
1553 
1554  GmshIO gmsh_io(mesh);
1555 
1556  if (mesh.processor_id() == 0)
1557  gmsh_io.read("meshes/circle.msh");
1559 
1560  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), dof_id_type(14));
1561  }
1562 
1564  {
1565  LOG_UNIT_TEST;
1566 
1568 
1569  GmshIO gmsh_io(mesh);
1570 
1571  if (mesh.processor_id() == 0)
1572  gmsh_io.read("meshes/bcid_overlap.msh");
1574 
1575  CPPUNIT_ASSERT_EQUAL(mesh.get_boundary_info().get_sideset_name_map().size(),
1576  std::size_t(2));
1577  CPPUNIT_ASSERT_EQUAL(mesh.get_boundary_info().sideset_name(1),
1578  std::string("srfBC4A"));
1579  CPPUNIT_ASSERT_EQUAL(mesh.get_boundary_info().sideset_name(2),
1580  std::string("srfBC4B"));
1581  CPPUNIT_ASSERT_EQUAL(mesh.get_subdomain_name_map().size(),
1582  std::size_t(2));
1583  CPPUNIT_ASSERT_EQUAL(mesh.subdomain_name(1),
1584  std::string("volBC3A"));
1585  CPPUNIT_ASSERT_EQUAL(mesh.subdomain_name(2),
1586  std::string("volBC3B"));
1587  }
1588 
1589  void testGoodSTL ()
1590  {
1591  LOG_UNIT_TEST;
1592 
1594 
1595  STLIO stl_io(mesh);
1596 
1597  if (mesh.processor_id() == 0)
1598  stl_io.read("meshes/Cluster_34.stl");
1600 
1601  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), dof_id_type(40));
1602  }
1603 
1605  {
1606  LOG_UNIT_TEST;
1607 
1609 
1610  STLIO stl_io(mesh);
1611 
1612  if (mesh.processor_id() == 0)
1613  stl_io.read("meshes/engraving.stl");
1615 
1616  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), dof_id_type(426));
1617  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), dof_id_type(215));
1618  }
1619 
1621  {
1622 #ifdef LIBMESH_HAVE_TETGEN
1623  LOG_UNIT_TEST;
1624 
1626 
1627  TetGenIO tetgen_io(mesh);
1628 
1629  if (mesh.processor_id() == 0)
1630  tetgen_io.read("meshes/tetgen_one_tet10.ele");
1632 
1634 
1635  // Mesh should contain 1 TET10 finite element
1636  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(1));
1637  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(10));
1638 
1639  // Element should have TET10 reference element volume
1640  const Elem * const elem = mesh.query_elem_ptr(0);
1641 
1642  // On a serial mesh we have every element everywhere
1643  if (mesh.is_serial())
1644  CPPUNIT_ASSERT(elem);
1645  else
1646  {
1647  bool have_elem = elem;
1648  mesh.comm().max(have_elem);
1649  CPPUNIT_ASSERT(have_elem);
1650  }
1651 
1652  if (elem)
1653  {
1654  const Real vol = elem->volume();
1655  LIBMESH_ASSERT_FP_EQUAL(vol, 1./6, TOLERANCE*TOLERANCE);
1656  }
1657 #endif
1658  }
1659 
1661  {
1662  LOG_UNIT_TEST;
1663 
1665 
1666  DynaIO dyna(mesh, /* keep_spline_nodes = */ false);
1667 
1668  if (mesh.processor_id() == 0)
1669  dyna.read("meshes/1_quad.bxt.gz");
1671 
1673 
1674  // We have 1 QUAD9 finite element
1675  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(1));
1676  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(9));
1677 
1678  helperTestingDynaQuad(mesh);
1679  }
1680 
1681 
1683  {
1684  LOG_UNIT_TEST;
1685 
1687 
1688  DynaIO dyna(mesh);
1689  if (mesh.processor_id() == 0)
1690  dyna.read("meshes/25_quad.bxt.gz");
1692 
1694 
1695  // We have 5^2 QUAD9 elements, with 11^2 nodes,
1696  // tied to 49 Node/NodeElem spline nodes
1697  CPPUNIT_ASSERT_EQUAL(mesh.n_elem(), static_cast<dof_id_type>(25+49));
1698  CPPUNIT_ASSERT_EQUAL(mesh.n_nodes(), static_cast<dof_id_type>(121+49));
1699 
1700  CPPUNIT_ASSERT_EQUAL(mesh.default_mapping_type(),
1702 
1703  unsigned char weight_index = mesh.default_mapping_data();
1704 
1705  for (const auto & elem : mesh.active_element_ptr_range())
1706  {
1707  if (elem->type() == NODEELEM)
1708  continue;
1709  LIBMESH_ASSERT_FP_EQUAL(Real(0.04), elem->volume(), TOLERANCE);
1710 
1711  for (unsigned int n=0; n != 9; ++n)
1712  CPPUNIT_ASSERT_EQUAL
1713  (elem->node_ref(n).get_extra_datum<Real>(weight_index),
1714  Real(1.0));
1715 
1716  unsigned int n_neighbors = 0, n_neighbors_expected = 2;
1717  for (unsigned int side=0; side != 4; ++side)
1718  if (elem->neighbor_ptr(side))
1719  n_neighbors++;
1720  Point c = elem->vertex_average();
1721 
1722  if (c(0) > 0.2 && c(0) < 0.8)
1723  n_neighbors_expected++;
1724  if (c(1) > 0.2 && c(1) < 0.8)
1725  n_neighbors_expected++;
1726 
1727  CPPUNIT_ASSERT_EQUAL(n_neighbors, n_neighbors_expected);
1728  }
1729 
1730  testMasterCenters(mesh);
1731 
1732 #ifdef LIBMESH_HAVE_SOLVER
1733 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1734  // Now test whether we can assign the desired constraint equations
1735  EquationSystems es(mesh);
1736  System & sys = es.add_system<LinearImplicitSystem>("test");
1737  sys.add_variable("u", SECOND); // to match QUAD9
1738  es.init();
1739 
1740  // We should have a constraint on every FE dof
1741  CPPUNIT_ASSERT_EQUAL(sys.get_dof_map().n_constrained_dofs(), static_cast<dof_id_type>(121));
1742 #endif // LIBMESH_ENABLE_CONSTRAINTS
1743 #endif // LIBMESH_HAVE_SOLVER
1744  }
1745 
1746  void testProjectionRegression(MeshBase & mesh, std::array<Real, 4> expected_norms)
1747  {
1748  int order = 0;
1749  for (const auto elem : mesh.element_ptr_range())
1750  order = std::max(order, int(elem->default_order()));
1751  TestCommWorld->max(order);
1752  CPPUNIT_ASSERT (order > 0);
1753 
1754  // Let's test that IGA constraints are preserved (in a relative
1755  // sense) when we clone a mesh.
1756  std::unique_ptr<MeshBase> mesh_clone = mesh.clone();
1757  CPPUNIT_ASSERT(*mesh_clone == mesh);
1758 
1759  EquationSystems es(mesh);
1760  System &sys = es.add_system<System> ("SimpleSystem");
1761  sys.add_variable("n", Order(order), RATIONAL_BERNSTEIN);
1762 
1763  es.init();
1764 
1765  sys.project_solution(sin_x_plus_cos_y, nullptr, es.parameters);
1766 
1767  // Make this easy to tweak in the future
1768  const Real my_tolerance = TOLERANCE;
1769 
1770  // Calculate some norms, skipping the spline points, and compare
1771  // to regression standard values
1772  std::set<unsigned int> skip_dimensions {0};
1773  const Real L2_norm =
1774  sys.calculate_norm(*sys.solution, 0, L2, &skip_dimensions);
1775 // std::cout.precision(16);
1776 // std::cout << "L2_norm = " << L2_norm << std::endl;
1777  LIBMESH_ASSERT_FP_EQUAL(L2_norm, expected_norms[0], my_tolerance);
1778  const Real Linf_norm =
1779  sys.calculate_norm(*sys.solution, 0, L_INF, &skip_dimensions);
1780 // std::cout << "Linf_norm = " << Linf_norm << std::endl;
1781  LIBMESH_ASSERT_FP_EQUAL(Linf_norm, expected_norms[1], my_tolerance);
1782  const Real H1_norm =
1783  sys.calculate_norm(*sys.solution, 0, H1_SEMINORM, &skip_dimensions);
1784 // std::cout << "H1_norm = " << H1_norm << std::endl;
1785  LIBMESH_ASSERT_FP_EQUAL(H1_norm, expected_norms[2], my_tolerance);
1786  const Real W1inf_norm =
1787  sys.calculate_norm(*sys.solution, 0, W1_INF_SEMINORM, &skip_dimensions);
1788 // std::cout << "W1inf_norm = " << W1inf_norm << std::endl;
1789  // W1_inf seems more sensitive to FP error...
1790  LIBMESH_ASSERT_FP_EQUAL(W1inf_norm, expected_norms[3], 10*my_tolerance);
1791  }
1792 
1793  void testDynaFileMappings (const std::string & filename, std::array<Real, 4> expected_norms)
1794  {
1796 
1797  DynaIO dyna(mesh);
1798  if (mesh.processor_id() == 0)
1799  dyna.read(filename);
1801 
1803 
1804  CPPUNIT_ASSERT_EQUAL(mesh.default_mapping_type(),
1806 
1807  // Useful when trying out different projection functions
1808  // std::cout << filename << ":" << std::endl;
1809 
1810  testMasterCenters(mesh);
1811 
1812  testProjectionRegression(mesh, expected_norms);
1813  }
1814 
1816  {
1817  LOG_UNIT_TEST;
1818 
1819  testDynaFileMappings("meshes/PressurizedCyl_Patch6_256Elem.bxt.gz",
1820  // Regression values for sin_x_plus_cos_y
1821  {{0.9639857809698268, 1.839870171669186,
1822  0.7089812562241862, 1.306121188539059}});
1823  }
1824 
1826  {
1827  LOG_UNIT_TEST;
1828 
1829  testDynaFileMappings("meshes/BlockWithHole_Patch9.bxt.gz",
1830  // Regression values for sin_x_plus_cos_y
1831  {{3.22612556930183, 1.97405365384733,
1832  2.53376235803176, 1.41374070517223}});
1833  }
1834 
1836  {
1837  LOG_UNIT_TEST;
1838 
1839  testDynaFileMappings("meshes/PlateWithHole_Patch8.bxt.gz",
1840  // Regression values for sin_x_plus_cos_y
1841  {{2.2812154374012, 1.974049990211937,
1842  1.791640772215248, 1.413679237529376}});
1843  }
1844 
1846  {
1847  LOG_UNIT_TEST;
1848 
1849  testDynaFileMappings("meshes/PressurizedCyl3d_Patch1_8Elem.bxt.gz",
1850  // Regression values for sin_x_plus_cos_y
1851  {{0.963612880188165, 1.82329452603503,
1852  0.707998701597943, 1.31399222566683}});
1853  }
1854 
1855  void testExodusFileMappings (const std::string & filename,
1856  std::array<Real, 4> expected_norms,
1857  bool use_disc_bex = false)
1858  {
1860 
1861  ExodusII_IO exii(mesh);
1862  // IGA Exodus meshes require ExodusII 8 or higher
1863  if (exii.get_exodus_version() < 800)
1864  return;
1865 
1866  // This should default to false
1867  if (use_disc_bex)
1868  exii.set_discontinuous_bex(true);
1869 
1870  if (mesh.processor_id() == 0)
1871  exii.read(filename);
1873 
1875 
1876  CPPUNIT_ASSERT_EQUAL(mesh.default_mapping_type(),
1878 
1879  testMasterCenters(mesh);
1880 
1881  testProjectionRegression(mesh, expected_norms);
1882 
1883  // Test a write when we're done reading; I was getting weirdness
1884  // from NetCDF at this point in a Moose output test.
1885  {
1886  ExodusII_IO exii(mesh);
1887 
1888  exii.write("exodus_file_mapping_out.e");
1889  }
1890 
1891 #ifdef LIBMESH_HAVE_VTK
1892  {
1893  VTKIO vtkout(mesh);
1894 
1895  vtkout.write("vtk_file_mapping_out.pvtu");
1896  }
1897 #endif
1898  }
1899 
1901  {
1902  LOG_UNIT_TEST;
1903 
1904  testExodusFileMappings("meshes/PlateWithHole_Patch8.e",
1905  // Regression values for sin_x_plus_cos_y
1906  {{2.2812154374012, 1.974049990211937,
1907  1.791640772215248, 1.413679237529376}});
1908  }
1909 
1911  {
1912  LOG_UNIT_TEST;
1913 
1914  testExodusFileMappings("meshes/two_quads_two_blocks.e",
1915  // Regression values for sin_x_plus_cos_y
1916  {{2.03496953073072, 1.97996853164955,
1917  1.18462134113435, 1.03085301158959}});
1918  }
1920  {
1921  LOG_UNIT_TEST;
1922 
1923  testExodusFileMappings("meshes/two_element_iga_in.e",
1924  // Regression values for sin_x_plus_cos_y
1925  {{1.26865962862531, 1.42562070158386,
1926  1.54905363492342, 1.29782906548366}});
1927  }
1928 
1930  {
1931  LOG_UNIT_TEST;
1932 
1933  testExodusFileMappings("meshes/PressurizedCyl3d_Patch1_8Elem.e",
1934  {{0.963612880188165, 1.82329452603503,
1935  0.707998701597943, 1.31399222566683}});
1936  }
1937 
1939  {
1940  LOG_UNIT_TEST;
1941 
1942  testExodusFileMappings("meshes/PlateWithHole_Patch8.e",
1943  // Regression values for sin_x_plus_cos_y
1944  //
1945  // These are *not* the same as for the continuous plate, because
1946  // we do the C^TKCx=C^Tf trick to pull back projections to the
1947  // spline nodes, and the pseudoinverse here minimizes the error in
1948  // a discretization-dependent norm, not in a Sobolev norm. For
1949  // these coarse meshes, our Sobolev norms can end up being ~0.1%
1950  // different.
1951  {{2.28234312456534, 1.97439548757586,
1952  1.79290449809266, 1.41075128955985}},
1953  true);
1954  }
1955 
1957  {
1958  LOG_UNIT_TEST;
1959 
1960  testExodusFileMappings("meshes/two_quads_two_blocks.e",
1961  // Regression values for sin_x_plus_cos_y
1962  {{2.03496953073072, 1.97996853164955,
1963  1.18462134113435, 1.03085301158959}},
1964  true);
1965  }
1967  {
1968  LOG_UNIT_TEST;
1969 
1970  testExodusFileMappings("meshes/two_element_iga_in.e",
1971  // Regression values for sin_x_plus_cos_y
1972  {{1.26877626663365, 1.42553698909339,
1973  1.54810114917177, 1.29792704408979}},
1974  true);
1975  }
1976 
1978  {
1979  LOG_UNIT_TEST;
1980 
1981  testExodusFileMappings("meshes/PressurizedCyl3d_Patch1_8Elem.e",
1982  {{0.963855209590556, 1.8234396424318,
1983  0.708286572453382, 1.31468940958327}},
1984  true);
1985  }
1986 };
1987 
virtual void read(const std::string &name) override
This method implements reading a mesh from a specified file.
Definition: abaqus_io.C:227
T libmesh_real(T a)
CPPUNIT_TEST_SUITE_REGISTRATION(MeshInputTest)
ElemType
Defines an enum for geometric element types.
virtual void read(const std::string &name) override
Reads in a mesh in the Dyna format from the ASCII file given by name.
Definition: dyna_io.C:139
This is the EquationSystems class.
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.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
Reading and writing meshes in the Gmsh format.
Definition: gmsh_io.h:51
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
void testExodusDiscWriteAddedSidesTriDisc()
Definition: mesh_input.C:1250
bool have_parameter(std::string_view) const
Definition: parameters.h:395
virtual void read(const std::string &base_filename) override
Implements reading the mesh from several different files.
Definition: nemesis_io.C:214
void testCopyNodalSolutionImpl(const std::string &filename)
Definition: mesh_input.C:507
void testGoodSTL()
Definition: mesh_input.C:1589
Reading and writing meshes in (a subset of) LS-DYNA format.
Definition: dyna_io.h:52
virtual void read(const std::string &) override
This method implements reading a mesh from a specified file in TetGen format.
Definition: tetgen_io.C:34
std::unique_ptr< PointLocatorBase > sub_point_locator() const
Definition: mesh_base.C:1638
This class implements reading and writing triangle meshes in the STL format.
Definition: stl_io.h:45
The AbaqusIO class is a preliminary implementation for reading Abaqus mesh files in ASCII format...
Definition: abaqus_io.h:41
void setUp()
Definition: mesh_input.C:218
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:67
void testExodusFileMappingsCyl3d()
Definition: mesh_input.C:1929
This class is used as both an external data structure for passing around Exodus file header informati...
void testDynaNoSplines()
Definition: mesh_input.C:1660
void allow_renumbering(bool allow)
If false is passed in then this mesh will no longer be renumbered when being prepared for use...
Definition: mesh_base.h:1196
void testExodusWriteAddedSidesHexC0()
Definition: mesh_input.C:1311
void testExodusDiscWriteAddedSidesTetDisc()
Definition: mesh_input.C:1306
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
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:969
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
void testExodusDiscWriteAddedSidesHexDisc()
Definition: mesh_input.C:1340
static constexpr Real TOLERANCE
unsigned int dim
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
Definition: fe_map.C:1628
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void testExodusDiscWriteAddedSidesMixedHexC0()
Definition: mesh_input.C:1329
void testExodusCopyNodalSolutionReplicated()
Definition: mesh_input.C:575
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void testDynaFileMappingsFEMEx5()
Definition: mesh_input.C:1815
void testExodusFileMappingsPlateWithHole()
Definition: mesh_input.C:1900
void testExodusCopyNodalSolutionDistributed()
Definition: mesh_input.C:578
void testExodusCopyElementSolutionReplicated()
Definition: mesh_input.C:664
void prepare_for_use(const bool skip_renumber_nodes_and_elements, const bool skip_find_neighbors)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:759
void testNemesisReadImpl()
Definition: mesh_input.C:1350
static const unsigned int type_to_dim_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the geometric dimension of the ele...
Definition: elem.h:635
Number six_x_plus_sixty_y(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: mesh_input.C:30
void sum(T &r) const
void testExodusWriteAddedSidesMixedEdgeC0()
Definition: mesh_input.C:1199
const std::map< boundary_id_type, std::string > & get_sideset_name_map() const
void barrier() const
virtual void read(const std::string &mesh_file) override
This method implements reading a mesh from a specified file.
Definition: stl_io.C:165
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
MeshBase & mesh
virtual void write(const std::string &) override
Output the mesh without solutions to a .pvtu file.
void testExodusFileMappingsTwoElemIGA()
Definition: mesh_input.C:1919
Number point_value(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
Definition: system.C:2421
void skip_noncritical_partitioning(bool skip)
If true is passed in then the elements on this mesh will no longer be (re)partitioned, and the nodes on this mesh will only be repartitioned if they are found "orphaned" via coarsening or other removal of the last element responsible for their node/element processor id consistency.
Definition: mesh_base.h:1236
const Parallel::Communicator & comm() const
void testLowOrderEdgeBlocks()
Definition: mesh_input.C:373
void testAbaqusReadFirst()
Definition: mesh_input.C:1486
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
void tearDown()
Definition: mesh_input.C:221
void boundary_ids(const Node *node, std::vector< boundary_id_type > &vec_to_fill) const
Fills a user-provided std::vector with the boundary ids associated with Node node.
void testNemesisCopyNodalSolutionReplicated()
Definition: mesh_input.C:582
This class implements reading and writing meshes in the TetGen format.
Definition: tetgen_io.h:47
The libMesh namespace provides an interface to certain functionality in the library.
void testAbaqusRead(const std::string &fname, dof_id_type n_elem, dof_id_type n_nodes)
Definition: mesh_input.C:1467
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
void testNemesisCopyElementSolutionDistributed()
Definition: mesh_input.C:674
virtual std::unique_ptr< MeshBase > clone() const =0
Virtual "copy constructor".
void testNemesisCopyElementVectorReplicated()
Definition: mesh_input.C:850
void testDynaFileMappingsCyl3d()
Definition: mesh_input.C:1845
This is the MeshBase class.
Definition: mesh_base.h:75
std::size_t n_boundary_ids() const
ElemMappingType default_mapping_type() const
Returns the default master space to physical space mapping basis functions to be used on newly added ...
Definition: mesh_base.h:812
void testNemesisReadReplicated()
Definition: mesh_input.C:1374
void testNemesisCopyNodalSolutionDistributed()
Definition: mesh_input.C:585
void testExodusFileMappings(const std::string &filename, std::array< Real, 4 > expected_norms, bool use_disc_bex=false)
Definition: mesh_input.C:1855
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
void testCopyElementSolutionImpl(const std::string &filename)
Definition: mesh_input.C:591
void testNemesisSingleElementReplicated()
Definition: mesh_input.C:732
virtual bool is_serial() const
Definition: mesh_base.h:211
void testExodusCopyElementSolutionDistributed()
Definition: mesh_input.C:667
const dof_id_type n_nodes
Definition: tecplot_io.C:67
unsigned int number() const
Definition: system.h:2350
void testExodusDiscWriteAddedSidesTetC0()
Definition: mesh_input.C:1295
const T & get(std::string_view) const
Definition: parameters.h:426
int8_t boundary_id_type
Definition: id_types.h:51
void testAbaqusReadSecond()
Definition: mesh_input.C:1493
void testExodusDiscTwoBlocks()
Definition: mesh_input.C:1956
This is the MeshCommunication class.
const std::map< subdomain_id_type, std::string > & get_subdomain_name_map() const
Definition: mesh_base.h:1694
void testBadGmsh()
Definition: mesh_input.C:1522
The Nemesis_IO class implements reading parallel meshes in the Nemesis file format from Sandia Nation...
Definition: nemesis_io.h:51
unsigned char default_mapping_data() const
Returns any default data value used by the master space to physical space mapping.
Definition: mesh_base.h:830
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
bool is_parallel_format() const
Returns true iff this mesh file format and input class are parallelized, so that all processors can r...
Definition: mesh_input.h:87
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:444
void testDynaFileMappingsPlateWithHole()
Definition: mesh_input.C:1835
void testExodusWriteAddedSidesEdgeC0()
Definition: mesh_input.C:1187
void testExodusWriteAddedSidesMixedQuadC0()
Definition: mesh_input.C:1267
void testExodusDiscPlateWithHole()
Definition: mesh_input.C:1938
void write_element_data(const EquationSystems &es)
Write out element solution in parallel, without localizing the solution vector.
Definition: nemesis_io.C:1466
virtual dof_id_type max_elem_id() const =0
void copy_elemental_solution(System &system, std::string system_var_name, std::string exodus_var_name, unsigned int timestep=1)
If we read in a elemental solution while reading in a mesh, we can attempt to copy that elemental sol...
Definition: exodusII_io.C:1147
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1593
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
Definition: system.C:1724
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
void testExodusWriteAddedSidesQuadC0()
Definition: mesh_input.C:1255
void testNemesisCopyElementSolutionReplicated()
Definition: mesh_input.C:671
unsigned int add_variable(std::string_view 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:1357
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2478
std::string & subdomain_name(subdomain_id_type id)
Definition: mesh_base.C:1692
void testGoodSTLBinary()
Definition: mesh_input.C:1604
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
std::vector< char > title
virtual void read(const std::string &name) override
This method implements reading a mesh from a specified file.
Definition: exodusII_io.C:244
void testNemesisSingleElementDistributed()
Definition: mesh_input.C:735
void helperTestingDynaQuad(const MeshBase &mesh)
Definition: mesh_input.C:1423
void broadcast(T &data, const unsigned int root_id=0, const bool identical_sizes=false) const
void testSingleElementImpl(const std::string &filename)
Definition: mesh_input.C:682
Number sin_x_plus_cos_y(const Point &p, const Parameters &, const std::string &, const std::string &)
Definition: mesh_input.C:42
void testNemesisReadDistributed()
Definition: mesh_input.C:1377
std::string & sideset_name(boundary_id_type id)
virtual void write(const std::string &name) const =0
void testExodusWriteAddedSidesEdgeDisc()
Definition: mesh_input.C:1211
void testDynaFileMappingsBlockWithHole()
Definition: mesh_input.C:1825
void testCopyElementVectorImpl(const std::string &filename)
Definition: mesh_input.C:746
void testExodusWriteAddedSidesTetDisc()
Definition: mesh_input.C:1301
const std::set< boundary_id_type > & get_boundary_ids() const
void testExodusDiscWriteAddedSidesMixedQuadC0()
Definition: mesh_input.C:1273
void testExodusIGASidesets()
Definition: mesh_input.C:416
void testMasterCenters(const MeshBase &mesh)
Definition: mesh_input.C:1382
void testDynaFileMappings(const std::string &filename, std::array< Real, 4 > expected_norms)
Definition: mesh_input.C:1793
void testProjectionRegression(MeshBase &mesh, std::array< Real, 4 > expected_norms)
Definition: mesh_input.C:1746
static int get_exodus_version()
Definition: exodusII_io.C:168
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void testDynaReadElem()
Definition: mesh_input.C:1500
virtual const Elem * query_elem_ptr(const dof_id_type i) const =0
void testGmshBCIDOverlap()
Definition: mesh_input.C:1563
dof_id_type n_constrained_dofs() const
virtual void renumber_elem(dof_id_type old_id, dof_id_type new_id)=0
Changes the id of element old_id, both by changing elem(old_id)->id() and by moving elem(old_id) in t...
void max(const T &r, T &o, Request &req) const
A class to represent the internal "this should never happen" errors, to be thrown by "libmesh_error()...
virtual void write(const std::string &fname) override
This method implements writing a mesh to a specified file.
Definition: exodusII_io.C:2180
void build_edge_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &edge_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, edges, and boundary ids for those edges.
constexpr int added_sides_nxyz[]
Definition: mesh_input.C:53
void testExodusCopyElementVectorReplicated()
Definition: mesh_input.C:843
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
Definition: fe_map.C:2069
void testExodusDiscWriteAddedSidesTriC0()
Definition: mesh_input.C:1227
void testExodusDiscWriteAddedSidesHexC0()
Definition: mesh_input.C:1317
void broadcast(MeshBase &) const
Finds all the processors that may contain elements that neighbor my elements.
void copy_elemental_solution(System &system, std::string system_var_name, std::string exodus_var_name, unsigned int timestep=1)
If we read in a elemental solution while reading in a mesh, we can attempt to copy that elemental sol...
Definition: nemesis_io.C:1719
void testTetgenIO()
Definition: mesh_input.C:1620
virtual Real volume() const
Definition: elem.C:3429
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
void testExodusWriteAddedSidesMixedHexC0()
Definition: mesh_input.C:1323
Parameters parameters
Data structure holding arbitrary parameters.
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
Definition: system.h:1605
ExodusHeaderInfo read_header(const std::string &name)
Read only the header information, instead of the entire mesh.
Definition: exodusII_io.C:923
void testExodusDiscCyl3d()
Definition: mesh_input.C:1977
void build_line(UnstructuredMesh &mesh, const unsigned int nx, const Real xmin=0., const Real xmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 1D meshes.
void testDynaReadPatch()
Definition: mesh_input.C:1682
void write_element_data_from_discontinuous_nodal_data(const EquationSystems &es, const std::set< std::string > *system_names=nullptr, const std::string &var_suffix="_elem_node_")
Similar to the function above, but instead of only handling (CONSTANT, MONOMIAL) data, writes out a general discontinuous solution field, e.g.
Definition: exodusII_io.C:1446
void testVTKPreserveElemIds()
Definition: mesh_input.C:225
void set_discontinuous_bex(bool disc_bex)
Set to true (false is the default) to generate independent nodes for every Bezier Extraction element ...
Definition: exodusII_io.C:2413
virtual void init()
Initialize all the systems.
void testExodusWriteAddedSidesTriC0()
Definition: mesh_input.C:1221
void testExodusReadHeader()
Definition: mesh_input.C:336
void testExodusDiscWriteAddedSidesEdgeDisc()
Definition: mesh_input.C:1216
void testExodusCopyElementVectorDistributed()
Definition: mesh_input.C:846
void testExodusWriteAddedSidesQuadDisc()
Definition: mesh_input.C:1279
void testExodusWriteAddedSidesMixedTriC0()
Definition: mesh_input.C:1233
virtual void read(const std::string &name) override
Reads in a mesh in the Gmsh *.msh format from the ASCII file given by name.
Definition: gmsh_io.C:149
virtual dof_id_type n_elem() const =0
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
processor_id_type processor_id() const
void testExodusDiscTwoElemIGA()
Definition: mesh_input.C:1966
const DofMap & get_dof_map() const
Definition: system.h:2374
void testExodusDiscWriteAddedSidesMixedEdgeC0()
Definition: mesh_input.C:1205
void testExodusDiscWriteAddedSidesQuadDisc()
Definition: mesh_input.C:1284
void testExodusWriteElementDataFromDiscontinuousNodalData()
Definition: mesh_input.C:858
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
void testExodusFileMappingsTwoBlocks()
Definition: mesh_input.C:1910
void testNemesisCopyElementVectorDistributed()
Definition: mesh_input.C:853
void testExodusWriteAddedSidesHexDisc()
Definition: mesh_input.C:1335
void ErrorVector unsigned int
Definition: adjoints_ex3.C:360
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
void testExodusDiscWriteAddedSidesQuadC0()
Definition: mesh_input.C:1261
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.
void testExodusDiscWriteAddedSidesEdgeC0()
Definition: mesh_input.C:1193
virtual dof_id_type n_nodes() const =0
uint8_t dof_id_type
Definition: id_types.h:67
Number designed_for_side_elems(const Point &p, const Parameters &param, const std::string &, const std::string &)
Definition: mesh_input.C:55
void testExodusDiscWriteAddedSidesMixedTriC0()
Definition: mesh_input.C:1239
void testExodusWriteAddedSidesTetC0()
Definition: mesh_input.C:1289
void testGoodGmsh()
Definition: mesh_input.C:1548
void testExodusWriteAddedSidesTriDisc()
Definition: mesh_input.C:1245
void testVTKPreserveSubdomainIds()
Definition: mesh_input.C:284