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> 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();
89 Point base_class_centroid = quad4.Elem::true_centroid();
108 for (
const auto & elem :
mesh.element_ptr_range())
110 Point derived_centroid = elem->true_centroid();
111 Point base_centroid = elem->Elem::true_centroid();
122 Real derived_volume = elem->volume();
123 Real base_volume = elem->Elem::volume();
143 test_true_centroid_and_volume(
PYRAMID5);
156 #ifdef LIBMESH_ENABLE_AMR 170 #endif // LIBMESH_ENABLE_AMR 191 prism6->
point(1) *= -1;
197 CPPUNIT_ASSERT_LESS(
TOLERANCE, std::abs(vol-gold_vol));
206 CPPUNIT_ASSERT_EQUAL(static_cast<dof_id_type>(3),
mesh.
n_nodes());
216 auto & middle_node = edge3->
node_ref(2);
217 auto & right_node = edge3->node_ref(1);
222 middle_node =
Point(0.5 + 1.e-3, 0., 0.);
228 middle_node =
Point(0.5 - 1.e-3, 0., 0.);
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);
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);
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);
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);
271 CPPUNIT_ASSERT(!invertible);
288 test_elem_invertible({
Point(-1, 0, 0),
Point(1, 0, 0),
Point(-0.5, 0, 0),
Point(
Real(1)/3, 0, 0)},
290 CPPUNIT_ASSERT(invertible);
294 test_elem_invertible({
Point(-1, 0, 0),
Point(1, 0, 0),
Point(-0.57, 0, 0),
Point(
Real(1)/3, 0, 0)},
296 CPPUNIT_ASSERT(!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)},
305 CPPUNIT_ASSERT(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)},
311 CPPUNIT_ASSERT(!invertible);
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);
334 for (
auto & pt : pts)
337 invertible = test_elem_invertible(pts,
QUAD4);
338 CPPUNIT_ASSERT(invertible);
345 for (
auto & pt : pts)
348 invertible = test_elem_invertible(pts,
QUAD4);
349 CPPUNIT_ASSERT(invertible);
356 for (
int cnt=0; cnt<3; ++cnt)
357 for (
auto & pt : pts)
360 invertible = test_elem_invertible(pts,
QUAD4);
361 CPPUNIT_ASSERT(invertible);
372 const Real alpha = .5;
375 test_elem_invertible({
Point(0, 0, 0),
Point(1, 0, 0),
Point(alpha, alpha, 0),
Point(0, 1, 0)},
QUAD4);
377 CPPUNIT_ASSERT(!invertible);
383 const Real alpha = -0.25;
386 test_elem_invertible({
Point(0, 0, 0),
Point(1, 0, 0),
Point(alpha, 1, 0),
Point(0, 1, 0)},
QUAD4);
388 CPPUNIT_ASSERT(!invertible);
394 const Real alpha = std::log(2);
397 test_elem_invertible({
Point(alpha, alpha, alpha),
398 Point(alpha, alpha, alpha),
399 Point(alpha, alpha, alpha),
402 CPPUNIT_ASSERT(!invertible);
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);
421 LIBMESH_ASSERT_FP_EQUAL(1.0, aspect_ratio,
TOLERANCE);
430 for (
auto & pt : pts)
434 LIBMESH_ASSERT_FP_EQUAL(1.0, aspect_ratio,
TOLERANCE);
441 for (
auto & pt : pts)
445 LIBMESH_ASSERT_FP_EQUAL(1.0, aspect_ratio,
TOLERANCE);
452 for (
int cnt=0; cnt<3; ++cnt)
453 for (
auto & pt : pts)
457 LIBMESH_ASSERT_FP_EQUAL(1.0, aspect_ratio,
TOLERANCE);
466 auto test_rhombus_quad = [
this](
Real theta)
468 Real ct = std::cos(theta);
469 Real st = std::sin(theta);
470 std::vector<Point> pts = {
473 Point(1. + ct, st, 0),
475 auto [elem, nodes] = this->construct_elem(pts,
QUAD4);
480 LIBMESH_ASSERT_FP_EQUAL(1.0/st, aspect_ratio,
TOLERANCE);
496 auto test_rectangle_quad = [
this](
Real a)
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);
504 LIBMESH_ASSERT_FP_EQUAL(a, aspect_ratio,
TOLERANCE);
508 test_rectangle_quad(2.0);
511 test_rectangle_quad(1e6);
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);
530 LIBMESH_ASSERT_FP_EQUAL(2.5, aspect_ratio,
TOLERANCE);
535 auto test_trapezoid_quad = [
this](
Real a)
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);
543 LIBMESH_ASSERT_FP_EQUAL(1./(1. - a), aspect_ratio,
TOLERANCE);
547 test_trapezoid_quad(1./3);
550 test_trapezoid_quad(0.5);
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);
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);
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);
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);
614 LIBMESH_ASSERT_FP_EQUAL(1.0, warpage,
TOLERANCE);
623 for (
auto & pt : pts)
626 warpage = elem->quality(
WARP);
627 LIBMESH_ASSERT_FP_EQUAL(1.0, warpage,
TOLERANCE);
634 for (
auto & pt : pts)
637 warpage = elem->quality(
WARP);
638 LIBMESH_ASSERT_FP_EQUAL(1.0, warpage,
TOLERANCE);
645 for (
int cnt=0; cnt<3; ++cnt)
646 for (
auto & pt : pts)
649 warpage = elem->quality(
WARP);
650 LIBMESH_ASSERT_FP_EQUAL(1.0, warpage,
TOLERANCE);
655 auto test_warped_quad = [
this](
Real h)
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);
665 LIBMESH_ASSERT_FP_EQUAL(
Real(1) / (h*h + 1), warpage,
TOLERANCE);
670 test_warped_quad(0.1);
675 test_warped_quad(0.3);
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);
696 LIBMESH_ASSERT_FP_EQUAL(90, min_angle,
TOLERANCE);
697 LIBMESH_ASSERT_FP_EQUAL(90, max_angle,
TOLERANCE);
707 for (
auto & pt : pts)
713 LIBMESH_ASSERT_FP_EQUAL(90, min_angle,
TOLERANCE);
714 LIBMESH_ASSERT_FP_EQUAL(90, max_angle,
TOLERANCE);
722 for (
auto & pt : pts)
728 LIBMESH_ASSERT_FP_EQUAL(90, min_angle,
TOLERANCE);
729 LIBMESH_ASSERT_FP_EQUAL(90, max_angle,
TOLERANCE);
737 for (
int cnt=0; cnt<3; ++cnt)
738 for (
auto & pt : pts)
744 LIBMESH_ASSERT_FP_EQUAL(90, min_angle,
TOLERANCE);
745 LIBMESH_ASSERT_FP_EQUAL(90, max_angle,
TOLERANCE);
756 auto test_rhombus_quad = [
this](
Real theta)
758 Real ct = std::cos(theta);
759 Real st = std::sin(theta);
760 std::vector<Point> pts = {
763 Point(1. + ct, st, 0),
765 auto [elem, nodes] = this->construct_elem(pts,
QUAD4);
794 auto test_rhombus_quad = [
this](
Real theta)
796 Real ct = std::cos(theta);
797 Real st = std::sin(theta);
798 std::vector<Point> pts = {
801 Point(1. + ct, st, 0),
803 auto [elem, nodes] = this->construct_elem(pts,
QUAD4);
813 LIBMESH_ASSERT_FP_EQUAL(std::abs(std::sin(theta)), jac,
TOLERANCE);
814 LIBMESH_ASSERT_FP_EQUAL(std::abs(std::sin(theta)), scaled_jac,
TOLERANCE);
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);
849 LIBMESH_ASSERT_FP_EQUAL(44.9985676771277, min_angle,
TOLERANCE);
850 LIBMESH_ASSERT_FP_EQUAL(90, max_angle,
TOLERANCE);
853 CPPUNIT_ASSERT_LESS(
Real(1), min_dihedral_angle);
854 CPPUNIT_ASSERT_LESS(
Real(1), max_dihedral_angle);
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);
877 LIBMESH_ASSERT_FP_EQUAL(h, jac,
TOLERANCE);
878 LIBMESH_ASSERT_FP_EQUAL(h/(h*h+1)/std::sqrt(h*h+2),
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),
895 auto [elem, nodes] = this->construct_elem(pts,
TET4);
905 LIBMESH_ASSERT_FP_EQUAL(2, jac,
TOLERANCE);
913 Real expected_volume)
917 const auto np = points.size();
918 std::unique_ptr<Elem> polygon = std::make_unique<C0Polygon>(np);
923 polygon->set_id() = 0;
927 const Real base_volume = elem->Elem::volume();
931 this->testC0PolygonMethods(
mesh, np);
939 testC0Polygon({{0,0}, {1,0}, {1,1}, {0,1}}, 1);
945 testC0Polygon({{0,0}, {1,0}, {1,2}, {-1,1}}, 2.5);
951 testC0Polygon({{0,0}, {1,0}, {1.5,0.5}, {1,1}, {0,1}}, 1.25);
957 testC0Polygon({{0,0}, {1,0}, {1.5,0.5}, {1,1}, {0,1}, {-0.5, 0.5}}, 1.5);
966 std::pair<std::unique_ptr<Elem>, std::vector<std::unique_ptr<Node>>>
970 const unsigned int n_points = pts.size();
973 std::vector<std::unique_ptr<Node>> nodes(n_points);
974 for (
unsigned int i=0; i<n_points; i++)
978 std::unique_ptr<Elem> elem =
Elem::build(elem_type,
nullptr);
981 libmesh_error_msg_if(elem->n_nodes() != n_points,
982 "Wrong number of points " 984 <<
" provided to build a " 987 for (
unsigned int i=0; i<n_points; i++)
988 elem->set_node(i, nodes[i].get());
991 return std::make_pair(std::move(elem), std::move(nodes));
997 unsigned int n_sides)
1000 bool found_elem = elem;
1002 CPPUNIT_ASSERT(found_elem);
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);
1017 const unsigned int nsover2 = n_sides/2;
1028 CPPUNIT_ASSERT(!elem->
is_edge(i));
1029 CPPUNIT_ASSERT(!elem->
is_face(i));
1035 CPPUNIT_ASSERT_EQUAL(nodes.size(), std::size_t(2));
1036 CPPUNIT_ASSERT(nodes[0] == 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);
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));
1066 bool found_elem = elem;
1068 CPPUNIT_ASSERT(found_elem);
1076 const int E = elem->
n_edges();
1077 const int F = elem->
n_faces();
1079 CPPUNIT_ASSERT_EQUAL(
int(elem->
n_nodes()), V);
1080 CPPUNIT_ASSERT_EQUAL(
int(elem->
n_sides()), F);
1083 CPPUNIT_ASSERT_EQUAL(V-E+F, 2);
1088 std::vector<int> sides_on_vertex(V, 0);
1092 CPPUNIT_ASSERT_EQUAL(side->type(),
C0POLYGON);
1093 FE += side->n_edges();
1094 const int SV = side->n_vertices();
1098 CPPUNIT_ASSERT_EQUAL(nos.size(), std::size_t(SV));
1102 Node * sidenode = side->node_ptr(v);
1105 ++sides_on_vertex[n];
1108 CPPUNIT_ASSERT(side->is_vertex(v));
1109 CPPUNIT_ASSERT(!side->is_edge(v));
1111 CPPUNIT_ASSERT(!elem->
is_edge(n));
1112 CPPUNIT_ASSERT(!elem->
is_face(n));
1114 CPPUNIT_ASSERT(std::find(nos.begin(), nos.end(), n) != nos.end());
1119 CPPUNIT_ASSERT_EQUAL(E*2,
FE);
1122 CPPUNIT_ASSERT_LESSEQUAL(FV, V*3);
1123 for (
auto sov : sides_on_vertex)
1124 CPPUNIT_ASSERT_LESSEQUAL(sov, 3);
1132 Real expected_volume)
1136 std::unique_ptr<Elem> polyhedron = std::make_unique<C0Polyhedron>(sides);
1138 polyhedron->set_id() = 0;
1142 const Real base_volume = elem->Elem::volume();
1146 this->testC0PolyhedronMethods(
mesh);
1165 const std::vector<std::vector<unsigned int>> nodes_on_side =
1174 std::vector<std::shared_ptr<Polygon>> sides(nodes_on_side.size());
1178 const auto & nodes_on_s = nodes_on_side[s];
1179 sides[s] = std::make_shared<C0Polygon>(nodes_on_s.size());
1184 testC0Polyhedron(sides, 1);
1196 auto [elem, nodes] = this->construct_elem(pts, elem_type);
1200 return elem->has_invertible_map();
1224 for (
const auto & elem :
mesh.element_ptr_range())
1226 Point derived_centroid = elem->true_centroid();
1227 Point base_centroid = elem->Elem::true_centroid();
1239 Real derived_volume = elem->volume();
1240 Real base_volume = elem->Elem::volume();
virtual Point true_centroid() const
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()
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
void testQuad4Invertible()
A Node is like a Point, but with more information.
void testEdge4Invertible()
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
unsigned int get_node_index(const Node *node_ptr) const
void testPrism6TrueCentroid()
virtual bool is_face(const unsigned int i) const =0
void testC0PolyhedronMethods(MeshBase &mesh)
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i)=0
libMesh::Parallel::Communicator * TestCommWorld
static constexpr Real TOLERANCE
void testQuad4MinMaxAngle()
void testC0PolygonHexagon()
bool test_elem_invertible(const std::vector< Point > &pts, ElemType elem_type)
This is the base class from which all geometric element types are derived.
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
const Parallel::Communicator & comm() const
void testPyramid5TrueCentroid()
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.
void testC0PolygonSquare()
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.
void testC0Polygon(const std::vector< Point > &points, Real expected_volume)
void testEdge3Invertible()
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const =0
A specific instantiation of the FEBase class.
void libmesh_ignore(const Args &...)
void testC0PolygonPentagon()
const Node & node_ref(const unsigned int i) const
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)
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)
void testTri3TrueCentroid()
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 testTet4DihedralAngle()
virtual unsigned int n_edges() const =0
CPPUNIT_TEST_SUITE_REGISTRATION(VolumeTest)
bool absolute_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
void testTri3AspectRatio()
void test_true_centroid_and_volume(ElemType elem_type)
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()
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int) const =0
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)
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)
virtual bool is_vertex(const unsigned int i) const =0
virtual Real volume() const
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...
void testC0Polyhedron(const std::vector< std::shared_ptr< Polygon >> &sides, Real expected_volume)
virtual bool has_invertible_map(Real tol=TOLERANCE *TOLERANCE) const
virtual unsigned int n_faces() const =0
virtual unsigned int n_sub_elem() const =0
void testQuad4AspectRatio()
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.
virtual unsigned int opposite_side(const unsigned int s) const
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
const Point & point(const unsigned int i) const
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
void testC0PolyhedronCube()
virtual bool is_edge(const unsigned int i) const =0
const Elem & get(const ElemType type_in)
virtual dof_id_type n_nodes() const =0
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
void testHex8TrueCentroid()
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.