libMesh
volume_test.C
Go to the documentation of this file.
1 // libmesh includes
2 #include <libmesh/cell_c0polyhedron.h>
3 #include <libmesh/elem.h>
4 #include <libmesh/enum_elem_type.h>
5 #include <libmesh/face_c0polygon.h>
6 #include <libmesh/mesh_generation.h>
7 #include <libmesh/mesh_modification.h>
8 #include <libmesh/mesh.h>
9 #include <libmesh/reference_elem.h>
10 #include <libmesh/replicated_mesh.h>
11 #include <libmesh/node.h>
12 #include <libmesh/enum_to_string.h>
13 #include <libmesh/tensor_value.h>
14 #include <libmesh/enum_elem_quality.h>
15 
16 // unit test includes
17 #include "test_comm.h"
18 #include "libmesh_cppunit.h"
19 
20 // C++ includes
21 #include <iomanip>
22 
23 using namespace libMesh;
24 
25 class VolumeTest : public CppUnit::TestCase
26 {
27 
28 public:
29  LIBMESH_CPPUNIT_TEST_SUITE( VolumeTest );
30  CPPUNIT_TEST( testTwistedVolume );
31  CPPUNIT_TEST( testEdge3Volume );
32  CPPUNIT_TEST( testEdge3Invertible );
33  CPPUNIT_TEST( testEdge4Invertible );
34  CPPUNIT_TEST( testQuad4Invertible );
35  CPPUNIT_TEST( testTri3TrueCentroid );
36  CPPUNIT_TEST( testQuad4TrueCentroid );
37  CPPUNIT_TEST( testPyramid5TrueCentroid );
38  CPPUNIT_TEST( testHex8TrueCentroid );
39  CPPUNIT_TEST( testPrism6TrueCentroid );
40  CPPUNIT_TEST( testHex20PLevelTrueCentroid );
41  CPPUNIT_TEST( testQuad4AspectRatio );
42  CPPUNIT_TEST( testQuad4Warpage );
43  CPPUNIT_TEST( testQuad4MinMaxAngle );
44  CPPUNIT_TEST( testQuad4Jacobian );
45  CPPUNIT_TEST( testTri3AspectRatio );
46  CPPUNIT_TEST( testTet4DihedralAngle );
47  CPPUNIT_TEST( testTet4Jacobian );
48  CPPUNIT_TEST( testC0PolygonSquare );
49  CPPUNIT_TEST( testC0PolygonQuad );
50  CPPUNIT_TEST( testC0PolygonPentagon );
51  CPPUNIT_TEST( testC0PolygonHexagon );
52  CPPUNIT_TEST( testC0PolyhedronCube );
53  CPPUNIT_TEST_SUITE_END();
54 
55 public:
56  void setUp()
57  {
58  }
59 
60  void tearDown()
61  {
62  }
63 
65  {
66  LOG_UNIT_TEST;
67 
68  // The true_centroid() == vertex_average() == (1/3, 1/3) for reference Tri3
69  {
70  const Elem & tri3 = ReferenceElem::get(TRI3);
71  Point true_centroid = tri3.true_centroid();
72  LIBMESH_ASSERT_FP_EQUAL(Real(1)/3, true_centroid(0), TOLERANCE*TOLERANCE);
73  LIBMESH_ASSERT_FP_EQUAL(Real(1)/3, true_centroid(1), TOLERANCE*TOLERANCE);
74  }
75  }
76 
78  {
79  LOG_UNIT_TEST;
80 
81  // Test Quad4::true_centroid() override
82  {
83  const Elem & quad4 = ReferenceElem::get(QUAD4);
84  Point true_centroid = quad4.true_centroid();
85  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(0), TOLERANCE*TOLERANCE);
86  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(1), TOLERANCE*TOLERANCE);
87 
88  // Compare to centroid computed via generic base class implementation
89  Point base_class_centroid = quad4.Elem::true_centroid();
90  CPPUNIT_ASSERT(true_centroid.absolute_fuzzy_equals(base_class_centroid, TOLERANCE*TOLERANCE));
91  }
92 
93  // Check that "optimized" Quad4::true_centroid() gives same result
94  // as "generic" Elem::true_centroid() on a mesh of 10% distorted elements.
95  {
97 
99  /*nx=*/3, /*ny=*/3,
100  /*xmin=*/0., /*xmax=*/1.,
101  /*ymin=*/0., /*ymax=*/1.,
102  QUAD4);
103 
105  /*factor=*/0.1,
106  /*perturb_boundary=*/false);
107 
108  for (const auto & elem : mesh.element_ptr_range())
109  {
110  Point derived_centroid = elem->true_centroid();
111  Point base_centroid = elem->Elem::true_centroid();
112 
113  // Debugging: check results in detail
114  // auto flags = libMesh::out.flags();
115  // libMesh::out << std::scientific << std::setprecision(16);
116  // libMesh::out << "derived_centroid = " << derived_centroid << std::endl;
117  // libMesh::out << "base_centroid = " << base_centroid << std::endl;
118  // libMesh::out.flags(flags);
119 
120  CPPUNIT_ASSERT(derived_centroid.absolute_fuzzy_equals(base_centroid, TOLERANCE*TOLERANCE));
121 
122  Real derived_volume = elem->volume();
123  Real base_volume = elem->Elem::volume();
124  LIBMESH_ASSERT_FP_EQUAL(base_volume, derived_volume, TOLERANCE*TOLERANCE);
125  }
126  }
127  }
128 
130  {
131  LOG_UNIT_TEST;
132 
133  // Test Pyramid5::true_centroid() gives the correct result for a reference element
134  {
135  const Elem & pyr5 = ReferenceElem::get(PYRAMID5);
136  Point true_centroid = pyr5.true_centroid();
137  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(0), TOLERANCE*TOLERANCE);
138  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(1), TOLERANCE*TOLERANCE);
139  LIBMESH_ASSERT_FP_EQUAL(0.25, true_centroid(2), TOLERANCE*TOLERANCE);
140  }
141 
142  // Currently there is not an optimized Pyramid5::true_centroid() to compare against
143  test_true_centroid_and_volume(PYRAMID5);
144  }
145 
146  void testHex8TrueCentroid() { LOG_UNIT_TEST; test_true_centroid_and_volume(HEX8); }
147  void testPrism6TrueCentroid() { LOG_UNIT_TEST; test_true_centroid_and_volume(PRISM6); }
148 
150  {
151  LOG_UNIT_TEST;
152 
153  // Test that Elem base class true_centroid() implementation works
154  // for an elevated p_level HEX20
155  {
156 #ifdef LIBMESH_ENABLE_AMR
159  /*nelem=*/1, /*nelem=*/1, /*nelem=*/1,
160  /*xmin=*/-1, /*xmax=*/1,
161  /*ymin=*/-1, /*ymax=*/1,
162  /*zmin=*/-1, /*zmax=*/1,
163  HEX20);
164  Elem * hex20 = mesh.elem_ptr(0);
165  hex20->set_p_level(1);
166  Point true_centroid = hex20->true_centroid();
167  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(0), TOLERANCE*TOLERANCE);
168  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(1), TOLERANCE*TOLERANCE);
169  LIBMESH_ASSERT_FP_EQUAL(0, true_centroid(2), TOLERANCE*TOLERANCE);
170 #endif // LIBMESH_ENABLE_AMR
171  }
172  }
173 
175  {
176  LOG_UNIT_TEST;
177 
179 
180  // Build an element type that will fall back on our generic
181  // quadrature-based Elem::volume()
183  /*nelem=*/1, /*nelem=*/1, /*nelem=*/1,
184  /*xmin=*/-1, /*xmax=*/1,
185  /*ymin=*/-1, /*ymax=*/1,
186  /*zmin=*/-1, /*zmax=*/1,
187  PRISM21);
188 
189  // Pick an element and twist it
190  Elem * prism6 = mesh.elem_ptr(0);
191  prism6->point(1) *= -1;
192  prism6->point(1) += 2*prism6->point(0);
193 
194  // The real test here is that volume() doesn't throw
195  const Real vol = prism6->volume();
196  const Real gold_vol = 3+Real(5)/9;
197  CPPUNIT_ASSERT_LESS(TOLERANCE, std::abs(vol-gold_vol));
198  }
199 
201  {
202  LOG_UNIT_TEST;
203 
205  MeshTools::Generation::build_line (mesh, /*nelem=*/1, 0., 1., EDGE3);
206  CPPUNIT_ASSERT_EQUAL(static_cast<dof_id_type>(3), mesh.n_nodes());
207 
208  auto edge3 = mesh.query_elem_ptr(0);
209  if (!edge3) // We may be on a distributed mesh
210  return;
211 
212  // Check unperturbed, straight edge case
213  LIBMESH_ASSERT_FP_EQUAL(1.0, edge3->volume(), TOLERANCE*TOLERANCE);
214 
215  // Get references to the individual Edge3 nodes
216  auto & middle_node = edge3->node_ref(2);
217  auto & right_node = edge3->node_ref(1);
218 
219  // Check middle node perturbed in +x direction case. This should
220  // not change the volume because it's still a straight line
221  // element.
222  middle_node = Point(0.5 + 1.e-3, 0., 0.);
223  LIBMESH_ASSERT_FP_EQUAL(1.0, edge3->volume(), TOLERANCE*TOLERANCE);
224 
225  // Check middle node perturbed in -x direction case. This should
226  // not change the volume because it's still a straight line
227  // element.
228  middle_node = Point(0.5 - 1.e-3, 0., 0.);
229  LIBMESH_ASSERT_FP_EQUAL(1.0, edge3->volume(), TOLERANCE*TOLERANCE);
230 
231  // Check volume of actual curved element against pre-computed value.
232  middle_node = Point(0.5, 0.25, 0.);
233  right_node = Point(1., 1., 0.);
234  LIBMESH_ASSERT_FP_EQUAL(1.4789428575446, edge3->volume(), TOLERANCE);
235 
236  // Compare with volume computed by base class Elem::volume() call
237  // which uses quadrature. We don't expect this to have full
238  // floating point accuracy.
239  middle_node = Point(0.5, 0.1, 0.);
240  right_node = Point(1., 0., 0.);
241  LIBMESH_ASSERT_FP_EQUAL(edge3->Elem::volume(), edge3->volume(), 1e-4);
242  }
243 
245  {
246  LOG_UNIT_TEST;
247 
248  // 1.) This is the original test which started the investigation
249  // of determining invertibility. In this test, the actual
250  // midpoint of nodes 0 and 1 is 0.5*(1.100328e2 + 1.176528e2) =
251  // 113.8428, so we can see that the middle node is closer to the
252  // left endpoint. In this case, it is too close and the element is
253  // not invertible.
254  bool invertible = test_elem_invertible({
255  Point(-3.566160e1, -6.690970e-1, 1.100328e2),
256  Point(-3.566160e1, -6.690970e-1, 1.176528e2),
257  Point(-3.566160e1, -6.690970e-1, 1.115568e2)}, EDGE3);
258  CPPUNIT_ASSERT(!invertible);
259 
260  // 2.) Just like case 1, but now node 2 is at the midpoint, so
261  // this case is invertible.
262  invertible = test_elem_invertible({
263  Point(-3.566160e1, -6.690970e-1, 1.100328e2),
264  Point(-3.566160e1, -6.690970e-1, 1.176528e2),
265  Point(-3.566160e1, -6.690970e-1, 113.8428)}, EDGE3);
266  CPPUNIT_ASSERT(invertible);
267 
268  // 3.) Non-collinear case where the mid-edge node is "above" and "way
269  // past" the right endpoint. This case is not invertible
270  invertible = test_elem_invertible({Point(0, 0, 0), Point(1, 0, 0), Point(3.5, 1.5, 0)}, EDGE3);
271  CPPUNIT_ASSERT(!invertible);
272  }
273 
275  {
276  LOG_UNIT_TEST;
277 
278  // Reference Elem should be invertible
279  {
280  const Elem & edge4 = ReferenceElem::get(EDGE4);
281  CPPUNIT_ASSERT(edge4.has_invertible_map());
282  }
283 
284  // If node 2 goes to the left past -5/9 = -.555, the element becomes non-invertible
285  {
286  // x2 > -5/9, the map is still invertible
287  bool invertible =
288  test_elem_invertible({Point(-1, 0, 0), Point(1, 0, 0), Point(-0.5, 0, 0), Point(Real(1)/3, 0, 0)},
289  EDGE4);
290  CPPUNIT_ASSERT(invertible);
291 
292  // x2 < -5/9, it is too close to x0 now
293  invertible =
294  test_elem_invertible({Point(-1, 0, 0), Point(1, 0, 0), Point(-0.57, 0, 0), Point(Real(1)/3, 0, 0)},
295  EDGE4);
296  CPPUNIT_ASSERT(!invertible);
297  }
298 
299  // If node 2 goes to the right past 5/21 ~ 0.2381, the element becomes non-invertible
300  {
301  // x2 < 5/21, the map should still be invertible
302  bool invertible =
303  test_elem_invertible({Point(-1, 0, 0), Point(1, 0, 0), Point(Real(3)/21, 0, 0), Point(Real(1)/3, 0, 0)},
304  EDGE4);
305  CPPUNIT_ASSERT(invertible);
306 
307  // x2 > 5/21, x2 is too close to x3 now
308  invertible =
309  test_elem_invertible({Point(-1, 0, 0), Point(1, 0, 0), Point(Real(6)/21, 0, 0), Point(Real(1)/3, 0, 0)},
310  EDGE4);
311  CPPUNIT_ASSERT(!invertible);
312  }
313  }
314 
316  {
317  LOG_UNIT_TEST;
318 
319  // Case 1: Test that rigid body rotations have no effect on the
320  // invertibility of the reference element
321  {
322  // 1a) The reference element rotated into various different different planes.
323  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 1, 0), Point(0, 1, 0)};
324  bool invertible = test_elem_invertible(pts, QUAD4);
325  CPPUNIT_ASSERT(invertible);
326 
327  // 1b) Rotate all points about x-axis by 90 degrees
328  Real cost = std::cos(.5*libMesh::pi);
329  Real sint = std::sin(.5*libMesh::pi);
330  RealTensorValue Rx(1, 0, 0,
331  0, cost, sint,
332  0, sint, cost);
333 
334  for (auto & pt : pts)
335  pt = Rx * pt;
336 
337  invertible = test_elem_invertible(pts, QUAD4);
338  CPPUNIT_ASSERT(invertible);
339 
340  // 1c) Rotate all points about z-axis by 90 degrees
341  RealTensorValue Rz(cost, -sint, 0,
342  sint, cost, 0,
343  0, 0, 1);
344 
345  for (auto & pt : pts)
346  pt = Rz * pt;
347 
348  invertible = test_elem_invertible(pts, QUAD4);
349  CPPUNIT_ASSERT(invertible);
350 
351  // 1d) Rotate all points about y-axis by 270 degrees
352  RealTensorValue Ry(cost, 0, sint,
353  0, 1, 0,
354  -sint, 0, cost);
355 
356  for (int cnt=0; cnt<3; ++cnt)
357  for (auto & pt : pts)
358  pt = Ry * pt;
359 
360  invertible = test_elem_invertible(pts, QUAD4);
361  CPPUNIT_ASSERT(invertible);
362  }
363 
364  // Case 2: Planar quad with top right vertex displaced to the position
365  // (alpha, alpha). Some different cases are described below.
366  // .) alpha==1: affine case, always invertible
367  // .) 1/2 < alpha < 1: planar case, invertible
368  // .) alpha<=1/2: planar case but node is now at center of the
369  // element, should give a zero/negative Jacobian on the displaced
370  // Node -> not invertible.
371  {
372  const Real alpha = .5;
373 
374  bool invertible =
375  test_elem_invertible({Point(0, 0, 0), Point(1, 0, 0), Point(alpha, alpha, 0), Point(0, 1, 0)}, QUAD4);
376 
377  CPPUNIT_ASSERT(!invertible);
378  }
379 
380  // Case 3) Top right corner is moved to (alpha, 1, 0). Element
381  // becomes non-invertible when alpha < 0.
382  {
383  const Real alpha = -0.25;
384 
385  bool invertible =
386  test_elem_invertible({Point(0, 0, 0), Point(1, 0, 0), Point(alpha, 1, 0), Point(0, 1, 0)}, QUAD4);
387 
388  CPPUNIT_ASSERT(!invertible);
389  }
390 
391  // Case 4) Degenerate case - all 4 points at same location. This
392  // zero-volume element does not have an invertible map.
393  {
394  const Real alpha = std::log(2);
395 
396  bool invertible =
397  test_elem_invertible({Point(alpha, alpha, alpha),
398  Point(alpha, alpha, alpha),
399  Point(alpha, alpha, alpha),
400  Point(alpha, alpha, alpha)}, QUAD4);
401 
402  CPPUNIT_ASSERT(!invertible);
403  }
404  }
405 
407  {
408  LOG_UNIT_TEST;
409 
410  // Case 1: Test that rigid body rotations of a unit square
411  // quadrilateral that have no effect on the quality of the
412  // element.
413  {
414  // Construct unit square QUAD4
415  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 1, 0), Point(0, 1, 0)};
416  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
417  libmesh_ignore(nodes);
418 
419  // 1a) Unit square aspect ratio should be == 1
420  Real aspect_ratio = elem->quality(ASPECT_RATIO);
421  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/aspect_ratio, TOLERANCE);
422 
423  // 1b) Rotate all points about x-axis by 90 degrees
424  Real cost = std::cos(.5*libMesh::pi);
425  Real sint = std::sin(.5*libMesh::pi);
426  RealTensorValue Rx(1, 0, 0,
427  0, cost, sint,
428  0, sint, cost);
429 
430  for (auto & pt : pts)
431  pt = Rx * pt;
432 
433  aspect_ratio = elem->quality(ASPECT_RATIO);
434  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/aspect_ratio, TOLERANCE);
435 
436  // 1c) Rotate all points about z-axis by 90 degrees
437  RealTensorValue Rz(cost, -sint, 0,
438  sint, cost, 0,
439  0, 0, 1);
440 
441  for (auto & pt : pts)
442  pt = Rz * pt;
443 
444  aspect_ratio = elem->quality(ASPECT_RATIO);
445  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/aspect_ratio, TOLERANCE);
446 
447  // 1d) Rotate all points about y-axis by 270 degrees
448  RealTensorValue Ry(cost, 0, sint,
449  0, 1, 0,
450  -sint, 0, cost);
451 
452  for (int cnt=0; cnt<3; ++cnt)
453  for (auto & pt : pts)
454  pt = Ry * pt;
455 
456  aspect_ratio = elem->quality(ASPECT_RATIO);
457  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/aspect_ratio, TOLERANCE);
458  }
459 
460  // Case 2: Rhombus QUAD4. This case should have an aspect ratio of
461  // 1/sin(theta), where theta is the acute interior angle of the
462  // rhombus.
463  {
464  // Helper lambda function that constructs a rhombus quad with
465  // interior acute angle theta.
466  auto test_rhombus_quad = [this](Real theta)
467  {
468  Real ct = std::cos(theta);
469  Real st = std::sin(theta);
470  std::vector<Point> pts = {
471  Point(0, 0, 0),
472  Point(1, 0, 0),
473  Point(1. + ct, st, 0),
474  Point( ct, st, 0)};
475  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
476  libmesh_ignore(nodes);
477 
478  // The expected aspect ratio for the rhombus is 1/sin(theta)
479  Real aspect_ratio = elem->quality(ASPECT_RATIO);
480  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0/st, /*actual=*/aspect_ratio, TOLERANCE);
481  };
482 
483  // 2a) Rhombus with interior angle theta=pi/6. The expected
484  // aspect ratio in this case is 1/std::sin(pi/6) = 2
485  test_rhombus_quad(libMesh::pi / 6);
486 
487  // 2b) Rhombus with interior angle theta=pi/3. The expected
488  // aspect ratio in this case is 1/std::sin(pi/3) = 2/sqrt(3) = 1.155
489  test_rhombus_quad(libMesh::pi / 3);
490  }
491 
492  // Case 3) Rectangle QUAD4. The "old" and "new" aspect ratio metrics
493  // return the same result for this case, whih is simply the ratio of
494  // the longest to shortest side.
495  {
496  auto test_rectangle_quad = [this](Real a)
497  {
498  std::vector<Point> pts = {Point(0, 0, 0), Point(a, 0, 0), Point(a, 1, 0), Point(0, 1, 0)};
499  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
500  libmesh_ignore(nodes);
501 
502  // The expected aspect ratio for the rectangle is "a"
503  Real aspect_ratio = elem->quality(ASPECT_RATIO);
504  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/a, /*actual=*/aspect_ratio, TOLERANCE);
505  };
506 
507  // 3a) Test case twice as long as it is tall
508  test_rectangle_quad(2.0);
509 
510  // 3b) Test no funny numerical stuff with extremely stretched case
511  test_rectangle_quad(1e6);
512  }
513 
514  // Case 4) Degenerate QUAD with zero length side. In the "old"
515  // aspect ratio metric we'd just get zero (as a stand-in for
516  // infinity) for this case since the minimum side length is zero,
517  // but using the new metric we get a non-zero value of 2.5. Thus,
518  // in the new metric it is possible to compare the aspect ratios
519  // of different degenerate QUAD elements rather than simply
520  // assigning all such elements a quality value of 0. Degenerate
521  // quadrilaterals are sometimes used as a "hack" to avoid having a
522  // separate subdomain of TRIs, and (surprisingly) many finite
523  // element operations just "work" on such elements.
524  {
525  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 0, 0), Point(0, 1, 0)};
526  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
527  libmesh_ignore(nodes);
528 
529  Real aspect_ratio = elem->quality(ASPECT_RATIO);
530  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/2.5, /*actual=*/aspect_ratio, TOLERANCE);
531  }
532 
533  // Case 5) Trapezoid QUAD
534  {
535  auto test_trapezoid_quad = [this](Real a)
536  {
537  // 0 <= a <= 1/2
538  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1-a, 1, 0), Point(a, 1, 0)};
539  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
540  libmesh_ignore(nodes);
541 
542  Real aspect_ratio = elem->quality(ASPECT_RATIO);
543  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1./(1. - a), /*actual=*/aspect_ratio, TOLERANCE);
544  };
545 
546  // 3a) Test "simple" trapezoid with expected aspect ratio 1.5
547  test_trapezoid_quad(1./3);
548 
549  // 3b) Test "degenerate" trapezoid with one zero length base
550  test_trapezoid_quad(0.5);
551  }
552  }
553 
555  {
556  LOG_UNIT_TEST;
557 
558  // Case 1) Reference TRI3
559  {
560  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(0, 1, 0)};
561  auto [elem, nodes] = this->construct_elem(pts, TRI3);
562  libmesh_ignore(nodes);
563 
564  // Compute the aspect ratio for the reference Tri
565  // The expected value is ~ 1.44338
566  Real aspect_ratio = elem->quality(ASPECT_RATIO);
567  // libMesh::out << "Unit TRI3 aspect ratio = " << aspect_ratio << std::endl;
568  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(5)/2/std::sqrt(Real(3)), /*actual=*/aspect_ratio, TOLERANCE);
569  }
570 
571  // Case 2) Equilateral TRI3
572  {
573  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(0.5, 0.5*std::sqrt(3), 0)};
574  auto [elem, nodes] = this->construct_elem(pts, TRI3);
575  libmesh_ignore(nodes);
576 
577  // Compute the aspect ratio for the reference Tri
578  // The expected value is 1.0
579  Real aspect_ratio = elem->quality(ASPECT_RATIO);
580  // libMesh::out << "Equilateral TRI3 aspect ratio = " << aspect_ratio << std::endl;
581  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(1), /*actual=*/aspect_ratio, TOLERANCE);
582  }
583 
584  // Case 3) Reference TRI3 with one leg length = L >> 1
585  {
586  Real L = 10.;
587  std::vector<Point> pts = {Point(0, 0, 0), Point(L, 0, 0), Point(0, 1, 0)};
588  auto [elem, nodes] = this->construct_elem(pts, TRI3);
589  libmesh_ignore(nodes);
590 
591  // Compute the aspect ratio for the reference Tri
592  // The expected value is ~ 11.5759
593  Real aspect_ratio = elem->quality(ASPECT_RATIO);
594  // libMesh::out << "TRI3 with leg length L = " << L << ", aspect ratio = " << aspect_ratio << std::endl;
595  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(1)/std::sqrt(Real(3)) * (Real(2)*L + Real(1)/(2*L)), /*actual=*/aspect_ratio, TOLERANCE);
596  }
597  }
598 
600  {
601  LOG_UNIT_TEST;
602 
603  // Case 1: Test that rigid body rotations of a unit square
604  // quadrilateral that have no effect on the quality of the
605  // element.
606  {
607  // Construct unit square QUAD4
608  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 1, 0), Point(0, 1, 0)};
609  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
610  libmesh_ignore(nodes);
611 
612  // 1a) Any flat element should have warp == 1
613  Real warpage = elem->quality(WARP);
614  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/warpage, TOLERANCE);
615 
616  // 1b) Rotate all points about x-axis by 90 degrees
617  Real cost = std::cos(.5*libMesh::pi);
618  Real sint = std::sin(.5*libMesh::pi);
619  RealTensorValue Rx(1, 0, 0,
620  0, cost, sint,
621  0, sint, cost);
622 
623  for (auto & pt : pts)
624  pt = Rx * pt;
625 
626  warpage = elem->quality(WARP);
627  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/warpage, TOLERANCE);
628 
629  // 1c) Rotate all points about z-axis by 90 degrees
630  RealTensorValue Rz(cost, -sint, 0,
631  sint, cost, 0,
632  0, 0, 1);
633 
634  for (auto & pt : pts)
635  pt = Rz * pt;
636 
637  warpage = elem->quality(WARP);
638  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/warpage, TOLERANCE);
639 
640  // 1d) Rotate all points about y-axis by 270 degrees
641  RealTensorValue Ry(cost, 0, sint,
642  0, 1, 0,
643  -sint, 0, cost);
644 
645  for (int cnt=0; cnt<3; ++cnt)
646  for (auto & pt : pts)
647  pt = Ry * pt;
648 
649  warpage = elem->quality(WARP);
650  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/1.0, /*actual=*/warpage, TOLERANCE);
651  }
652 
653  // Case 2: Unit square quadrilateral with Node 2 displaced by a distance h in the z-direction.
654  {
655  auto test_warped_quad = [this](Real h)
656  {
657  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 1, h), Point(0, 1, 0)};
658  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
659  libmesh_ignore(nodes);
660 
661  Real warpage = elem->quality(WARP);
662  // libMesh::out << "QUAD with node 3 displaced by h = "
663  // << h << " in the z-direction , warpage = " << warpage
664  // << std::endl;
665  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(1) / (h*h + 1), /*actual=*/warpage, TOLERANCE);
666  };
667 
668  // For h = 0.1, the expected warpage value is ~ 0.991, so not
669  // much different than the value of 1.0 for a flat element.
670  test_warped_quad(0.1);
671 
672  // For h = 0.3, the expected warpage value is ~ 0.917431, which
673  // is pretty close to the lower bound (0.9) of "acceptable"
674  // warpage, as suggested in the Cubit manual.
675  test_warped_quad(0.3);
676  }
677  }
678 
680  {
681  LOG_UNIT_TEST;
682 
683  // Case 1: Test that rigid body rotations of a unit square
684  // quadrilateral that have no effect on the quality of the
685  // element.
686  {
687  // Construct unit square QUAD4
688  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(1, 1, 0), Point(0, 1, 0)};
689  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
690  libmesh_ignore(nodes);
691 
692  // 1a) Reference Elem should have min angle == max angle == pi/2
693  {
694  Real min_angle = elem->quality(MIN_ANGLE);
695  Real max_angle = elem->quality(MAX_ANGLE);
696  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/min_angle, TOLERANCE);
697  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/max_angle, TOLERANCE);
698  }
699 
700  // 1b) Rotate all points about x-axis by 90 degrees
701  Real cost = std::cos(.5*libMesh::pi);
702  Real sint = std::sin(.5*libMesh::pi);
703  RealTensorValue Rx(1, 0, 0,
704  0, cost, sint,
705  0, sint, cost);
706 
707  for (auto & pt : pts)
708  pt = Rx * pt;
709 
710  {
711  Real min_angle = elem->quality(MIN_ANGLE);
712  Real max_angle = elem->quality(MAX_ANGLE);
713  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/min_angle, TOLERANCE);
714  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/max_angle, TOLERANCE);
715  }
716 
717  // 1c) Rotate all points about z-axis by 90 degrees
718  RealTensorValue Rz(cost, -sint, 0,
719  sint, cost, 0,
720  0, 0, 1);
721 
722  for (auto & pt : pts)
723  pt = Rz * pt;
724 
725  {
726  Real min_angle = elem->quality(MIN_ANGLE);
727  Real max_angle = elem->quality(MAX_ANGLE);
728  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/min_angle, TOLERANCE);
729  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/max_angle, TOLERANCE);
730  }
731 
732  // 1d) Rotate all points about y-axis by 270 degrees
733  RealTensorValue Ry(cost, 0, sint,
734  0, 1, 0,
735  -sint, 0, cost);
736 
737  for (int cnt=0; cnt<3; ++cnt)
738  for (auto & pt : pts)
739  pt = Ry * pt;
740 
741  {
742  Real min_angle = elem->quality(MIN_ANGLE);
743  Real max_angle = elem->quality(MAX_ANGLE);
744  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/min_angle, TOLERANCE);
745  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/max_angle, TOLERANCE);
746  }
747  }
748 
749  // Case 2: Rhombus QUAD4. This case should have an interior min
750  // angle of theta and an interior max angle of pi-theta, where
751  // "theta" is the specified amount that we "sheared" the element
752  // by on creation
753  {
754  // Helper lambda function that constructs a rhombus quad with
755  // interior acute angle theta.
756  auto test_rhombus_quad = [this](Real theta)
757  {
758  Real ct = std::cos(theta);
759  Real st = std::sin(theta);
760  std::vector<Point> pts = {
761  Point(0, 0, 0),
762  Point(1, 0, 0),
763  Point(1. + ct, st, 0),
764  Point( ct, st, 0)};
765  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
766  libmesh_ignore(nodes);
767 
768  Real min_angle = elem->quality(MIN_ANGLE);
769  Real max_angle = elem->quality(MAX_ANGLE);
770  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(180) / libMesh::pi * theta, /*actual=*/min_angle, TOLERANCE);
771  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(180) / libMesh::pi * (libMesh::pi - theta), /*actual=*/max_angle, TOLERANCE);
772  };
773 
774  // 2a) Rhombus with interior angle theta=pi/6.
775  test_rhombus_quad(libMesh::pi / 6);
776 
777  // 2b) Rhombus with interior angle theta=pi/3.
778  test_rhombus_quad(libMesh::pi / 3);
779  }
780  }
781 
783  {
784  LOG_UNIT_TEST;
785 
786  // Rhombus QUAD4. This case should have a JACOBIAN and
787  // SCALED_JACOBIAN value of sin(theta), where "theta" is the
788  // specified amount that we "sheared" the element by on creation.
789  // The JACOBIAN and SCALED_JACOBIAN are the same for this element
790  // because the edge lengths are all = 1.
791  {
792  // Helper lambda function that constructs a rhombus quad with
793  // interior acute angle theta.
794  auto test_rhombus_quad = [this](Real theta)
795  {
796  Real ct = std::cos(theta);
797  Real st = std::sin(theta);
798  std::vector<Point> pts = {
799  Point(0, 0, 0),
800  Point(1, 0, 0),
801  Point(1. + ct, st, 0),
802  Point( ct, st, 0)};
803  auto [elem, nodes] = this->construct_elem(pts, QUAD4);
804  libmesh_ignore(nodes);
805 
806  Real jac = elem->quality(JACOBIAN);
807  Real scaled_jac = elem->quality(SCALED_JACOBIAN);
808 
809  // Debugging
810  // libMesh::out << "jac = " << jac << std::endl;
811  // libMesh::out << "scaled_jac = " << scaled_jac << std::endl;
812 
813  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/std::abs(std::sin(theta)), /*actual=*/jac, TOLERANCE);
814  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/std::abs(std::sin(theta)), /*actual=*/scaled_jac, TOLERANCE);
815  };
816 
817  // 2a) Rhombus with interior angle theta=pi/6.
818  test_rhombus_quad(libMesh::pi / 6);
819 
820  // 2b) Rhombus with interior angle theta=pi/3.
821  test_rhombus_quad(libMesh::pi / 3);
822  }
823  }
824 
826  {
827  LOG_UNIT_TEST;
828 
829  // Construct a tetrahedron whose projection into the x-y plane looks like the unit square,
830  // and where the apex node is just slightly out of the x-y plane. This element should have
831  // reasonable MIN,MAX_ANGLE values but MIN,MAX_DIHEDRAL angles of nearly 0. This is an
832  // example that demonstrates one should not solely use edge angles to determine Elem quality
833  // in 3D.
834  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(0, 1, 0), Point(1, 1, 0.01)};
835  auto [elem, nodes] = this->construct_elem(pts, TET4);
836  libmesh_ignore(nodes);
837 
838  Real min_angle = elem->quality(MIN_ANGLE);
839  Real max_angle = elem->quality(MAX_ANGLE);
840  Real min_dihedral_angle = elem->quality(MIN_DIHEDRAL_ANGLE);
841  Real max_dihedral_angle = elem->quality(MAX_DIHEDRAL_ANGLE);
842 
843  // Debugging
844  // libMesh::out << "Squashed Tet4 min_angle = " << min_angle << " degrees" << std::endl;
845  // libMesh::out << "Squashed Tet4 max_angle = " << max_angle << " degrees" << std::endl;
846  // libMesh::out << "Squashed Tet4 min_dihedral_angle = " << min_dihedral_angle << " degrees" << std::endl;
847  // libMesh::out << "Squashed Tet4 max_dihedral_angle = " << max_dihedral_angle << " degrees" << std::endl;
848 
849  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/44.9985676771277, /*actual=*/min_angle, TOLERANCE);
850  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/90, /*actual=*/max_angle, TOLERANCE);
851 
852  // Assert that both min and max dihedral angles are less than 1 degree (a very low quality element)
853  CPPUNIT_ASSERT_LESS(Real(1), min_dihedral_angle);
854  CPPUNIT_ASSERT_LESS(Real(1), max_dihedral_angle);
855  }
856 
858  {
859  LOG_UNIT_TEST;
860 
861  {
862  // Same element as in testTet4DihedralAngle(). Here we verify that
863  // the element is low quality since the SCALED_JACOBIAN is < O(h)
864  // as h -> 0, and h can be arbitrarily small in the squashed element.
865  Real h = Real(1)/100;
866  std::vector<Point> pts = {Point(0, 0, 0), Point(1, 0, 0), Point(0, 1, 0), Point(1, 1, h)};
867  auto [elem, nodes] = this->construct_elem(pts, TET4);
868  libmesh_ignore(nodes);
869 
870  Real jac = elem->quality(JACOBIAN);
871  Real scaled_jac = elem->quality(SCALED_JACOBIAN);
872 
873  // Debugging
874  // libMesh::out << "Squashed Tet4 jac = " << jac << std::endl;
875  // libMesh::out << "Squashed Tet4 scaled_jac = " << scaled_jac << std::endl;
876 
877  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/h, /*actual=*/jac, TOLERANCE);
878  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/h/(h*h+1)/std::sqrt(h*h+2),
879  /*actual=*/scaled_jac, TOLERANCE * TOLERANCE);
880  }
881 
882  {
883  // A representative Tet generated by calling build_cube(). This
884  // element has minimum interior edge angle of approximately
885  // 35.26 degrees, so it it is not a simple refinement of the
886  // reference element. It is not located at the origin, but still
887  // has a simple to compute value for the JACOBIAN and
888  // SCALED_JACOBIAN quality metrics.
889  std::vector<Point> pts = {
890  Point(2.5, 2.5, 4.5),
891  Point(3.5, 1.5, 5.5),
892  Point(1.5, 1.5, 5.5),
893  Point(2.5, 2.5, 5.5)
894  };
895  auto [elem, nodes] = this->construct_elem(pts, TET4);
896  libmesh_ignore(nodes);
897 
898  Real jac = elem->quality(JACOBIAN);
899  Real scaled_jac = elem->quality(SCALED_JACOBIAN);
900 
901  // Debugging
902  // libMesh::out << "build_cube() Tet4 jac = " << jac << std::endl;
903  // libMesh::out << "build_cube() Tet4 scaled_jac = " << scaled_jac << std::endl;
904 
905  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/2, /*actual=*/jac, TOLERANCE);
906  LIBMESH_ASSERT_FP_EQUAL(/*expected=*/Real(1)/std::sqrt(Real(6)), /*actual=*/scaled_jac, TOLERANCE);
907  }
908  }
909 
910 
911 
912  void testC0Polygon(const std::vector<Point> & points,
913  Real expected_volume)
914  {
916 
917  const auto np = points.size();
918  std::unique_ptr<Elem> polygon = std::make_unique<C0Polygon>(np);
919 
920  for (auto p : make_range(np))
921  polygon->set_node(p, mesh.add_point(points[p], p));
922 
923  polygon->set_id() = 0;
924  Elem * elem = mesh.add_elem(std::move(polygon));
925 
926  const Real derived_volume = elem->volume();
927  const Real base_volume = elem->Elem::volume();
928  LIBMESH_ASSERT_FP_EQUAL(base_volume, derived_volume, TOLERANCE*TOLERANCE);
929  LIBMESH_ASSERT_FP_EQUAL(derived_volume, expected_volume, TOLERANCE*TOLERANCE);
930 
931  this->testC0PolygonMethods(mesh, np);
932  }
933 
934 
935 
937  {
938  LOG_UNIT_TEST;
939  testC0Polygon({{0,0}, {1,0}, {1,1}, {0,1}}, 1);
940  }
941 
943  {
944  LOG_UNIT_TEST;
945  testC0Polygon({{0,0}, {1,0}, {1,2}, {-1,1}}, 2.5);
946  }
947 
949  {
950  LOG_UNIT_TEST;
951  testC0Polygon({{0,0}, {1,0}, {1.5,0.5}, {1,1}, {0,1}}, 1.25);
952  }
953 
955  {
956  LOG_UNIT_TEST;
957  testC0Polygon({{0,0}, {1,0}, {1.5,0.5}, {1,1}, {0,1}, {-0.5, 0.5}}, 1.5);
958  }
959 
960 protected:
961 
962  // Helper function that is called by test_elem_invertible() to build an Elem
963  // of the requested elem_type from the provided Points. Note: the
964  // Nodes which are constructed in order to construct the Elem are
965  // also returned since
966  std::pair<std::unique_ptr<Elem>, std::vector<std::unique_ptr<Node>>>
967  construct_elem(const std::vector<Point> & pts,
968  ElemType elem_type)
969  {
970  const unsigned int n_points = pts.size();
971 
972  // Create Nodes
973  std::vector<std::unique_ptr<Node>> nodes(n_points);
974  for (unsigned int i=0; i<n_points; i++)
975  nodes[i] = Node::build(pts[i], /*id*/ i);
976 
977  // Create Elem, assign nodes
978  std::unique_ptr<Elem> elem = Elem::build(elem_type, /*parent*/ nullptr);
979 
980  // Make sure we were passed consistent input to build this type of Elem
981  libmesh_error_msg_if(elem->n_nodes() != n_points,
982  "Wrong number of points "
983  << n_points
984  << " provided to build a "
985  << Utility::enum_to_string(elem_type));
986 
987  for (unsigned int i=0; i<n_points; i++)
988  elem->set_node(i, nodes[i].get());
989 
990  // Return Elem and Nodes we created
991  return std::make_pair(std::move(elem), std::move(nodes));
992  }
993 
994 
995  // Helper function to factor out common tests
997  unsigned int n_sides)
998  {
999  Elem * elem = mesh.query_elem_ptr(0);
1000  bool found_elem = elem;
1001  mesh.comm().max(found_elem);
1002  CPPUNIT_ASSERT(found_elem);
1003 
1004  if (!elem)
1005  return;
1006 
1007  CPPUNIT_ASSERT_EQUAL(elem->type(), C0POLYGON);
1008  CPPUNIT_ASSERT_EQUAL(elem->n_nodes(), n_sides);
1009  CPPUNIT_ASSERT_EQUAL(elem->n_sub_elem(), n_sides);
1010  CPPUNIT_ASSERT_EQUAL(elem->n_sides(), n_sides);
1011  CPPUNIT_ASSERT_EQUAL(elem->n_vertices(), n_sides);
1012  CPPUNIT_ASSERT_EQUAL(elem->n_edges(), n_sides);
1013 
1014  // Even number of sides
1015  if (!(n_sides%2))
1016  {
1017  const unsigned int nsover2 = n_sides/2;
1018  for (auto i : make_range(nsover2))
1019  {
1020  CPPUNIT_ASSERT_EQUAL(elem->opposite_side(i), i+nsover2);
1021  CPPUNIT_ASSERT_EQUAL(elem->opposite_side(i+nsover2), i);
1022  }
1023  }
1024 
1025  for (unsigned int i : make_range(n_sides))
1026  {
1027  CPPUNIT_ASSERT(elem->is_vertex(i));
1028  CPPUNIT_ASSERT(!elem->is_edge(i));
1029  CPPUNIT_ASSERT(!elem->is_face(i));
1030  CPPUNIT_ASSERT(elem->is_node_on_side(i,i));
1031  CPPUNIT_ASSERT(elem->is_node_on_edge(i,i));
1032  CPPUNIT_ASSERT(elem->is_node_on_side((i+1)%n_sides,i));
1033  CPPUNIT_ASSERT(elem->is_node_on_edge((i+1)%n_sides,i));
1034  std::vector<unsigned int> nodes = elem->nodes_on_side(i);
1035  CPPUNIT_ASSERT_EQUAL(nodes.size(), std::size_t(2));
1036  CPPUNIT_ASSERT(nodes[0] == i ||
1037  nodes[1] == i);
1038  CPPUNIT_ASSERT(nodes[0] == (i+1)%n_sides ||
1039  nodes[1] == (i+1)%n_sides);
1040  std::vector<unsigned int> edge_nodes = elem->nodes_on_edge(i);
1041  CPPUNIT_ASSERT(nodes == edge_nodes);
1042 
1043 
1044  CPPUNIT_ASSERT_EQUAL(elem->local_side_node(i,0), i);
1045  CPPUNIT_ASSERT_EQUAL(elem->local_side_node(i,1), (i+1)%n_sides);
1046  CPPUNIT_ASSERT_EQUAL(elem->local_edge_node(i,0), i);
1047  CPPUNIT_ASSERT_EQUAL(elem->local_edge_node(i,1), (i+1)%n_sides);
1048 
1049  auto edge = elem->side_ptr(i);
1050  CPPUNIT_ASSERT_EQUAL(edge->type(), EDGE2);
1051  CPPUNIT_ASSERT_EQUAL(edge->node_ptr(0), mesh.node_ptr(i));
1052  CPPUNIT_ASSERT_EQUAL(edge->node_ptr(1), mesh.node_ptr((i+1)%n_sides));
1053  }
1054 
1055  CPPUNIT_ASSERT(!elem->is_flipped());
1056  elem->flip(&mesh.get_boundary_info());
1057  CPPUNIT_ASSERT(elem->is_flipped());
1058  elem->flip(&mesh.get_boundary_info());
1059  }
1060 
1061 
1062  // Helper function to factor out common tests
1064  {
1065  Elem * elem = mesh.query_elem_ptr(0);
1066  bool found_elem = elem;
1067  mesh.comm().max(found_elem);
1068  CPPUNIT_ASSERT(found_elem);
1069 
1070  if (!elem)
1071  return;
1072 
1073  CPPUNIT_ASSERT_EQUAL(elem->type(), C0POLYHEDRON);
1074 
1075  const int V = elem->n_vertices();
1076  const int E = elem->n_edges();
1077  const int F = elem->n_faces();
1078 
1079  CPPUNIT_ASSERT_EQUAL(int(elem->n_nodes()), V);
1080  CPPUNIT_ASSERT_EQUAL(int(elem->n_sides()), F);
1081 
1082  // Euler characteristic for polygons homeomorphic to spheres is 2
1083  CPPUNIT_ASSERT_EQUAL(V-E+F, 2);
1084 
1085  int FE = 0;
1086  int FV = 0;
1087 
1088  std::vector<int> sides_on_vertex(V, 0);
1089  for (auto s : make_range(F))
1090  {
1091  auto side = elem->build_side_ptr(s);
1092  CPPUNIT_ASSERT_EQUAL(side->type(), C0POLYGON);
1093  FE += side->n_edges();
1094  const int SV = side->n_vertices();
1095  FV += SV;
1096 
1097  auto nos = elem->nodes_on_side(s);
1098  CPPUNIT_ASSERT_EQUAL(nos.size(), std::size_t(SV));
1099 
1100  for (auto v : make_range(SV))
1101  {
1102  Node * sidenode = side->node_ptr(v);
1103  const unsigned int n = elem->get_node_index(sidenode);
1104  CPPUNIT_ASSERT(n != invalid_uint);
1105  ++sides_on_vertex[n];
1106 
1107  // We're just looking at first order so far
1108  CPPUNIT_ASSERT(side->is_vertex(v));
1109  CPPUNIT_ASSERT(!side->is_edge(v));
1110  CPPUNIT_ASSERT(elem->is_vertex(n));
1111  CPPUNIT_ASSERT(!elem->is_edge(n));
1112  CPPUNIT_ASSERT(!elem->is_face(n));
1113  CPPUNIT_ASSERT(elem->is_node_on_side(n, s));
1114  CPPUNIT_ASSERT(std::find(nos.begin(), nos.end(), n) != nos.end());
1115  }
1116  }
1117 
1118  // We hit every edge from 2 faces
1119  CPPUNIT_ASSERT_EQUAL(E*2, FE);
1120 
1121  // We hit every vertex from at least 3 faces
1122  CPPUNIT_ASSERT_LESSEQUAL(FV, V*3); // V*3 <= FV
1123  for (auto sov : sides_on_vertex)
1124  CPPUNIT_ASSERT_LESSEQUAL(sov, 3); // 3 <= sov
1125 
1126  CPPUNIT_ASSERT(!elem->is_flipped());
1127  }
1128 
1129 
1130 
1131  void testC0Polyhedron(const std::vector<std::shared_ptr<Polygon>> & sides,
1132  Real expected_volume)
1133  {
1135 
1136  std::unique_ptr<Elem> polyhedron = std::make_unique<C0Polyhedron>(sides);
1137 
1138  polyhedron->set_id() = 0;
1139  Elem * elem = mesh.add_elem(std::move(polyhedron));
1140 
1141  const Real derived_volume = elem->volume();
1142  const Real base_volume = elem->Elem::volume();
1143  LIBMESH_ASSERT_FP_EQUAL(base_volume, derived_volume, TOLERANCE*TOLERANCE);
1144  LIBMESH_ASSERT_FP_EQUAL(derived_volume, expected_volume, TOLERANCE*TOLERANCE);
1145 
1146  this->testC0PolyhedronMethods(mesh);
1147  }
1148 
1149 
1150 
1152  {
1154 
1155  mesh.add_point(Point(0, 0, 0), 0);
1156  mesh.add_point(Point(1, 0, 0), 1);
1157  mesh.add_point(Point(1, 1, 0), 2);
1158  mesh.add_point(Point(0, 1, 0), 3);
1159  mesh.add_point(Point(0, 0, 1), 4);
1160  mesh.add_point(Point(1, 0, 1), 5);
1161  mesh.add_point(Point(1, 1, 1), 6);
1162  mesh.add_point(Point(0, 1, 1), 7);
1163 
1164  // See notes in elem_test.h
1165  const std::vector<std::vector<unsigned int>> nodes_on_side =
1166  { {0, 1, 2, 3}, // min z
1167  {0, 1, 5, 4}, // min y
1168  {2, 6, 5, 1}, // max x
1169  {2, 3, 7, 6}, // max y
1170  {0, 4, 7, 3}, // min x
1171  {5, 6, 7, 4} }; // max z
1172 
1173  // Build all the sides.
1174  std::vector<std::shared_ptr<Polygon>> sides(nodes_on_side.size());
1175 
1176  for (auto s : index_range(nodes_on_side))
1177  {
1178  const auto & nodes_on_s = nodes_on_side[s];
1179  sides[s] = std::make_shared<C0Polygon>(nodes_on_s.size());
1180  for (auto i : index_range(nodes_on_s))
1181  sides[s]->set_node(i, mesh.node_ptr(nodes_on_s[i]));
1182  }
1183 
1184  testC0Polyhedron(sides, 1);
1185  }
1186 
1187 
1188 
1189  // Helper function that builds the specified type of Elem from a
1190  // vector of Points and returns the value of has_invertible_map()
1191  // for that Elem.
1192  bool test_elem_invertible(const std::vector<Point> & pts,
1193  ElemType elem_type)
1194  {
1195  // Construct Elem of desired type
1196  auto [elem, nodes] = this->construct_elem(pts, elem_type);
1197  libmesh_ignore(nodes);
1198 
1199  // Return whether or not this Elem has an invertible map
1200  return elem->has_invertible_map();
1201  }
1202 
1203  // Helper function for testing true_centroid() and volume() implementations for
1204  // 3D elements
1206  {
1207  // Check that derived class true_centroid() gives same result as
1208  // the base class Elem::true_centroid() on a 2x2x2 mesh of 10%
1209  // distorted elements.
1210  {
1212 
1214  /*nelem=*/2, /*nelem=*/2, /*nelem=*/2,
1215  /*xmin=*/-1, /*xmax=*/1,
1216  /*ymin=*/-1, /*ymax=*/1,
1217  /*zmin=*/-1, /*zmax=*/1,
1218  elem_type);
1219 
1221  /*factor=*/0.1,
1222  /*perturb_boundary=*/false);
1223 
1224  for (const auto & elem : mesh.element_ptr_range())
1225  {
1226  Point derived_centroid = elem->true_centroid();
1227  Point base_centroid = elem->Elem::true_centroid();
1228 
1229  // Debugging: check results in detail
1230  // auto flags = libMesh::out.flags();
1231  // libMesh::out << std::scientific << std::setprecision(16);
1232  // libMesh::out << "derived_centroid = " << derived_centroid << std::endl;
1233  // libMesh::out << "base_centroid = " << base_centroid << std::endl;
1234  // libMesh::out.flags(flags);
1235 
1236  CPPUNIT_ASSERT(derived_centroid.absolute_fuzzy_equals(base_centroid, 50*TOLERANCE*TOLERANCE));
1237 
1238  // Make sure that base class and "optimized" routines for computing the cell volume agree
1239  Real derived_volume = elem->volume();
1240  Real base_volume = elem->Elem::volume();
1241  LIBMESH_ASSERT_FP_EQUAL(base_volume, derived_volume, 50*TOLERANCE*TOLERANCE);
1242  }
1243  }
1244  }
1245 };
1246 
virtual Point true_centroid() const
Definition: elem.C:594
void set_p_level(const unsigned int p)
Sets the value of the p-refinement level for the element.
ElemType
Defines an enum for geometric element types.
void testQuad4TrueCentroid()
Definition: volume_test.C:77
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void testQuad4Invertible()
Definition: volume_test.C:315
A Node is like a Point, but with more information.
Definition: node.h:52
void testEdge4Invertible()
Definition: volume_test.C:274
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
unsigned int get_node_index(const Node *node_ptr) const
Definition: elem.h:2545
void testPrism6TrueCentroid()
Definition: volume_test.C:147
virtual bool is_face(const unsigned int i) const =0
void testC0PolyhedronMethods(MeshBase &mesh)
Definition: volume_test.C:1063
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i)=0
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
static constexpr Real TOLERANCE
void testQuad4MinMaxAngle()
Definition: volume_test.C:679
void testC0PolygonHexagon()
Definition: volume_test.C:954
void testEdge3Volume()
Definition: volume_test.C:200
void distort(MeshBase &mesh, const Real factor, const bool perturb_boundary=false)
Randomly perturb the nodal locations.
bool test_elem_invertible(const std::vector< Point > &pts, ElemType elem_type)
Definition: volume_test.C:1192
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
void testQuad4Jacobian()
Definition: volume_test.C:782
MeshBase & mesh
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
const Parallel::Communicator & comm() const
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 testPyramid5TrueCentroid()
Definition: volume_test.C:129
The libMesh namespace provides an interface to certain functionality in the library.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:165
void testC0PolygonSquare()
Definition: volume_test.C:936
virtual bool is_flipped() const =0
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
This is the MeshBase class.
Definition: mesh_base.h:75
void testC0Polygon(const std::vector< Point > &points, Real expected_volume)
Definition: volume_test.C:912
void testEdge3Invertible()
Definition: volume_test.C:244
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const =0
A specific instantiation of the FEBase class.
Definition: fe.h:127
void libmesh_ignore(const Args &...)
void testC0PolygonPentagon()
Definition: volume_test.C:948
const Node & node_ref(const unsigned int i) const
Definition: elem.h:2529
virtual unsigned int n_nodes() const =0
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const =0
void testC0PolygonMethods(MeshBase &mesh, unsigned int n_sides)
Definition: volume_test.C:996
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:444
void testTri3TrueCentroid()
Definition: volume_test.C:64
void testTwistedVolume()
Definition: volume_test.C:174
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const =0
Similar to Elem::local_side_node(), but instead of a side id, takes an edge id and a node id on that ...
void testC0PolygonQuad()
Definition: volume_test.C:942
void testTet4DihedralAngle()
Definition: volume_test.C:825
virtual unsigned int n_edges() const =0
CPPUNIT_TEST_SUITE_REGISTRATION(VolumeTest)
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:972
void testTri3AspectRatio()
Definition: volume_test.C:554
void test_true_centroid_and_volume(ElemType elem_type)
Definition: volume_test.C:1205
virtual void flip(BoundaryInfo *boundary_info)=0
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
std::string enum_to_string(const T e)
void testHex20PLevelTrueCentroid()
Definition: volume_test.C:149
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int) const =0
void testQuad4Warpage()
Definition: volume_test.C:599
virtual const Elem * elem_ptr(const dof_id_type i) const =0
virtual unsigned int n_sides() const =0
static std::unique_ptr< Node > build(const Node &n)
Definition: node.h:315
virtual unsigned int n_vertices() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual std::unique_ptr< Elem > side_ptr(unsigned int i)=0
virtual const Elem * query_elem_ptr(const dof_id_type i) const =0
void max(const T &r, T &o, Request &req) const
std::pair< std::unique_ptr< Elem >, std::vector< std::unique_ptr< Node > > > construct_elem(const std::vector< Point > &pts, ElemType elem_type)
Definition: volume_test.C:967
virtual bool is_vertex(const unsigned int i) const =0
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 testC0Polyhedron(const std::vector< std::shared_ptr< Polygon >> &sides, Real expected_volume)
Definition: volume_test.C:1131
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.
virtual bool has_invertible_map(Real tol=TOLERANCE *TOLERANCE) const
Definition: elem.C:2881
virtual unsigned int n_faces() const =0
void setUp()
Definition: volume_test.C:56
virtual unsigned int n_sub_elem() const =0
void testQuad4AspectRatio()
Definition: volume_test.C:406
virtual const Node * node_ptr(const dof_id_type i) const =0
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
void testTet4Jacobian()
Definition: volume_test.C:857
virtual unsigned int opposite_side(const unsigned int s) const
Definition: elem.C:3511
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Point & point(const unsigned int i) const
Definition: elem.h:2453
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 testC0PolyhedronCube()
Definition: volume_test.C:1151
virtual bool is_edge(const unsigned int i) const =0
const Elem & get(const ElemType type_in)
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.
virtual dof_id_type n_nodes() const =0
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
void testHex8TrueCentroid()
Definition: volume_test.C:146
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
void tearDown()
Definition: volume_test.C:60
const Real pi
.
Definition: libmesh.h:299