20 #include "libmesh/cell_prism6.h" 21 #include "libmesh/edge_edge2.h" 22 #include "libmesh/face_quad4.h" 23 #include "libmesh/face_tri3.h" 24 #include "libmesh/enum_io_package.h" 25 #include "libmesh/enum_order.h" 26 #include "libmesh/fe_lagrange_shape_1D.h" 39 void cell_integral(
const Elem * elem, T & t);
59 struct NodalVolumeIntegral
65 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1};
66 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2};
69 auto tri3_phi = [](
int i,
Real x,
Real y)
83 libmesh_error_msg(
"Invalid shape function index i = " << i);
87 for (
int i=0; i<6; ++i)
97 std::array<Real, 6> V {};
165 const unsigned int s)
const 167 libmesh_assert_less (s,
n_sides());
173 std::vector<unsigned>
176 libmesh_assert_less(s,
n_sides());
177 auto trim = (s > 0 && s < 4) ? 0 : 1;
181 std::vector<unsigned>
184 libmesh_assert_less(e,
n_edges());
189 const unsigned int e)
const 191 libmesh_assert_less (e,
n_edges());
220 libmesh_assert_less (i, this->
n_sides());
222 std::unique_ptr<Elem> face;
229 face = std::make_unique<Tri3>();
236 face = std::make_unique<Quad4>();
240 libmesh_error_msg(
"Invalid side i = " << i);
244 for (
auto n : face->node_index_range())
247 face->set_interior_parent(
this);
248 face->inherit_data_from(*
this);
256 const unsigned int i)
259 side->set_interior_parent(
this);
260 side->inherit_data_from(*
this);
267 return this->simple_build_edge_ptr<Edge2,Prism6>(i);
274 this->simple_build_edge_ptr<Prism6>(edge, i,
EDGE2);
281 std::vector<dof_id_type> & conn)
const 316 libmesh_error_msg(
"Unsupported IO package " << iop);
322 #ifdef LIBMESH_ENABLE_AMR 329 { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0},
330 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0},
331 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0},
332 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0},
333 { .25, .25, 0.0, .25, .25, 0.0},
334 { .25, 0.0, .25, .25, 0.0, .25}
340 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0},
341 { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0},
342 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0},
343 { .25, .25, 0.0, .25, .25, 0.0},
344 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0},
345 { 0.0, .25, .25, 0.0, .25, .25}
351 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0},
352 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0},
353 { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0},
354 { .25, 0.0, .25, .25, 0.0, .25},
355 { 0.0, .25, .25, 0.0, .25, .25},
356 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}
362 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0},
363 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0},
364 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0},
365 { .25, .25, 0.0, .25, .25, 0.0},
366 { 0.0, .25, .25, 0.0, .25, .25},
367 { .25, 0.0, .25, .25, 0.0, .25}
373 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0},
374 { .25, .25, 0.0, .25, .25, 0.0},
375 { .25, 0.0, .25, .25, 0.0, .25},
376 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0},
377 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0},
378 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}
384 { .25, .25, 0.0, .25, .25, 0.0},
385 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0},
386 { 0.0, .25, .25, 0.0, .25, .25},
387 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0},
388 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0},
389 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}
395 { .25, 0.0, .25, .25, 0.0, .25},
396 { 0.0, .25, .25, 0.0, .25, .25},
397 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5},
398 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5},
399 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5},
400 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}
406 { .25, .25, 0.0, .25, .25, 0.0},
407 { 0.0, .25, .25, 0.0, .25, .25},
408 { .25, 0.0, .25, .25, 0.0, .25},
409 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0},
410 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5},
411 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}
421 NodalVolumeIntegral nvi;
422 cell_integral(
this, nvi);
428 for (
int i=0; i<6; ++i)
430 cp += this->
point(i) * nvi.V[i];
440 cell_integral(
this, vi);
454 libmesh_assert_less (perm_num, 6);
455 const unsigned int side = perm_num % 2;
456 const unsigned int rotate = perm_num / 2;
458 for (
unsigned int i = 0; i !=
rotate; ++i)
500 libmesh_assert_less (s, 5);
501 if (s == 0 || s == 4)
511 template <
typename T>
512 void cell_integral(
const Elem * elem, T & t)
523 -x0/2 + x1/2 - x3/2 + x4/2,
524 x0/2 - x1/2 - x3/2 + x4/2,
531 -x0/2 + x2/2 - x3/2 + x5/2,
532 x0/2 - x2/2 - x3/2 + x5/2,
539 x0/2 - x2/2 - x3/2 + x5/2,
540 x0/2 - x1/2 - x3/2 + x4/2
550 static const Real w2D[N2D] =
552 1.5902069087198858469718450103758e-01_R,
553 9.0979309128011415302815498962418e-02_R,
554 1.5902069087198858469718450103758e-01_R,
555 9.0979309128011415302815498962418e-02_R
558 static const Real xi[N2D] =
560 1.5505102572168219018027159252941e-01_R,
561 6.4494897427831780981972840747059e-01_R,
562 1.5505102572168219018027159252941e-01_R,
563 6.4494897427831780981972840747059e-01_R
566 static const Real eta[N2D] =
568 1.7855872826361642311703513337422e-01_R,
569 7.5031110222608118177475598324603e-02_R,
570 6.6639024601470138670269327409637e-01_R,
571 2.8001991549907407200279599420481e-01_R
579 static const Real zeta[N1D] =
581 -std::sqrt(
Real(3))/3,
582 std::sqrt(
Real(3))/3.
585 for (
int i=0; i<N2D; ++i)
588 Point dx_dzeta_q = dx_dzeta[0] + eta[i]*dx_dzeta[1] + xi[i]*dx_dzeta[2];
590 for (
int j=0; j<N1D; ++j)
594 dx_dxi_q = dx_dxi[0] + zeta[j]*dx_dxi[1],
595 dx_deta_q = dx_deta[0] + zeta[j]*dx_deta[1];
598 t.accumulate_qp(xi[i], eta[i], zeta[j], w2D[i] *
triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q));
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.
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
Node ** _nodes
Pointers to the nodes we are connected to.
static const int num_edges
ElemType side_type(const unsigned int s) const override final
virtual Real volume() const override
Specialized function for computing the element volume.
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
virtual unsigned int n_sides() const override final
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
This is the base class from which all geometric element types are derived.
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
virtual BoundingBox loose_bounding_box() const
virtual bool is_edge(const unsigned int i) const override
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
virtual bool is_face(const unsigned int i) const override
virtual BoundingBox loose_bounding_box() const override
Builds a bounding box out of the nodal positions.
The libMesh namespace provides an interface to certain functionality in the library.
virtual Order default_order() 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...
static const int num_children
virtual bool is_vertex(const unsigned int i) const override
void swap3neighbors(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three neighbor_ptrs, "rotating" them.
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
void swap3nodes(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three node_ptrs, "rotating" them.
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a QUAD4 or TRI3 built coincident with face i.
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
virtual unsigned int n_sub_elem() const override
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
static const int nodes_per_side
Defines a Cartesian bounding box by the two corner extremum.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_edges() const override final
static const int num_sides
Geometric constants for all Prisms.
static const int nodes_per_edge
virtual bool has_affine_map() 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...
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the node of the edge to element node numbers.
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE2 or INFEDGE2 built coincident with face i.
static const unsigned int side_elems_map[num_sides][nodes_per_side]
This maps the child elements with the associated side of the parent element.
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
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
static const int num_nodes
Geometric constants for Prism6.
Real fe_lagrange_1D_linear_shape(const unsigned int i, const Real xi)
virtual Point true_centroid() const override
An Optimized numerical quadrature approach for computing the centroid of the Prism6.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override