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