19 #include "libmesh/edge_edge3.h" 20 #include "libmesh/face_tri7.h" 21 #include "libmesh/enum_io_package.h" 22 #include "libmesh/enum_order.h" 24 #ifdef LIBMESH_ENABLE_AMR 49 #ifdef LIBMESH_ENABLE_AMR 56 { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0},
57 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0},
58 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0},
59 { .375, -.125, 0.0, .75, 0.0, 0.0, 0.0},
60 {.09375,-.03125,-.03125, .125, -.125, .125,.84375},
61 { .375, 0.0, -.125, 0.0, 0.0, .75, 0.0},
62 { 5/r18,-1/r18,-1/r18,4/r18,-2/r18,4/r18, 0.5}
68 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0},
69 { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0},
70 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0},
71 { -.125, .375, 0.0, .75, 0.0, 0.0, 0.0},
72 { 0.0, .375, -.125, 0.0, .75, 0.0, 0.0},
73 {-.03125,.09375,-.03125, .125, .125, -.125,.84375},
74 {-1/r18, 5/r18,-1/r18,4/r18,4/r18,-2/r18, 0.5}
80 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0},
81 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0},
82 { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0},
83 {-.03125,-.03125,.09375, -.125, .125, .125,.84375},
84 { 0.0, -.125, .375, 0.0, .75, 0.0, 0.0},
85 { -.125, 0.0, .375, 0.0, 0.0, .75, 0.0},
86 {-1/r18,-1/r18, 5/r18,-2/r18,4/r18,4/r18, 0.5}
92 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0},
93 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0},
94 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0},
95 {-.03125,.09375,-.03125, .125, .125,-.125,.84375},
96 {-.03125,-.03125,.09375,-.125, .125, .125,.84375},
97 {.09375,-.03125,-.03125, .125,-.125, .125,.84375},
98 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}
102 const std::vector<std::pair<unsigned char, unsigned char>>
106 { {},{{0,1}},{{0,2}},{{0,3}},{{3,5}},{{0,5}},{{0,6}} },
108 { {{0,1}}, {},{{1,2}},{{1,3}},{{1,4}},{{3,4}},{{1,6}} },
110 { {{0,2}},{{1,2}}, {},{{4,5}},{{2,4}},{{2,5}},{{2,6}} },
112 { {{0,1}},{{1,2}},{{0,2}},{{3,4}},{{4,5}},{{3,5}},{{0,4},{1,5},{2,3}} }
143 const unsigned int s)
const 145 libmesh_assert_less (s,
n_sides());
151 std::vector<unsigned>
154 libmesh_assert_less(s,
n_sides());
158 std::vector<unsigned>
169 ((this->point(4) - this->
point(1))*2))
173 ((this->point(3) - this->
point(0))*2))
177 ((this->point(5) - this->
point(0))*2))
183 ((this->point(6) - this->
point(0))*3))
207 libmesh_assert_less (s, this->
n_sides());
227 libmesh_error_msg(
"Invalid side s = " << s);
234 unsigned int side_node)
const 236 libmesh_assert_less (side, this->
n_sides());
246 return this->simple_build_side_ptr<Edge3, Tri7>(i);
252 const unsigned int i)
254 this->simple_build_side_ptr<Tri7>(side, i,
EDGE3);
261 std::vector<dof_id_type> & conn)
const 310 libmesh_error_msg(
"Invalid sf = " << sf);
324 libmesh_error_msg(
"Unsupported IO package " << iop);
344 for (
unsigned d=0; d<LIBMESH_DIM; ++d)
347 for (
unsigned int p=1; p != 6; ++p)
348 center += this->
point(p)(d);
351 Real hd = std::abs(center - this->
point(0)(d));
352 for (
unsigned int p=1; p != 6; ++p)
353 hd = std::max(hd, std::abs(center - this->
point(p)(d)));
355 pmin(d) = center - hd;
356 pmax(d) = center + hd;
377 libmesh_error_msg(
"Invalid n = " << n);
384 const unsigned int v)
const 386 libmesh_assert_greater_equal (n, this->
n_vertices());
387 libmesh_assert_less (n, this->
n_nodes());
393 libmesh_assert_less (v, 3);
394 return static_cast<unsigned short int>(v);
399 libmesh_assert_less (v, 2);
416 std::pair<unsigned short int, unsigned short int>
419 libmesh_assert_greater_equal (n, this->
n_vertices());
420 libmesh_assert_less (n, this->
n_nodes());
421 return std::pair<unsigned short int, unsigned short int>
447 libmesh_assert_less (perm_num, 3);
449 for (
unsigned int i = 0; i != perm_num; ++i)
480 libmesh_assert_less (s, 3);
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
ElemType
Defines an enum for geometric element types.
void swap2boundaryedges(unsigned short e1, unsigned short e2, BoundaryInfo *boundary_info) const
Swaps two edges in boundary_info, if it is non-null.
Order
defines an enum for polynomial orders.
virtual bool has_affine_map() const override
static const unsigned short int _second_order_adjacent_vertices[num_sides][2]
Matrix that tells which vertices define the location of mid-side (or second-order) nodes...
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int) const override
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
virtual Order default_order() const override
static const int num_sides
Geometric constants for all Tris.
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual unsigned int n_nodes() const override
virtual bool is_face(const unsigned int i) const override
virtual unsigned int n_vertices() const override final
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
static const unsigned short int _second_order_vertex_child_number[num_nodes]
Vector that names a child sharing each second order node.
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
void swap3neighbors(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three neighbor_ptrs, "rotating" them.
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
static const unsigned short int _second_order_vertex_child_index[num_nodes]
Vector that names the child vertex index for each second order node.
virtual unsigned int n_sides() const override final
void swap3nodes(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three node_ptrs, "rotating" them.
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
virtual dof_id_type key() const override final
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
virtual Order default_side_order() const override
static const int num_nodes
Geometric constants for Tri7.
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
virtual bool is_vertex(const unsigned int i) const override
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Defines a Cartesian bounding box by the two corner extremum.
ElemType side_type(const unsigned int s) const override final
static const std::vector< std::pair< unsigned char, unsigned char > > _parent_bracketing_nodes[num_children][num_nodes]
Pairs of nodes that bracket child nodes when doing mesh refinement.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_sub_elem() const override
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
unsigned int center_node_on_side(const unsigned short side) const override final
virtual BoundingBox loose_bounding_box() const override
static dof_id_type compute_key(dof_id_type n0)
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
A Point defines a location in LIBMESH_DIM dimensional Real space.
dof_id_type node_id(const unsigned int i) const
const Point & point(const unsigned int i) const
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
static const int nodes_per_side
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
virtual bool is_edge(const unsigned int i) const override