cell_prism6.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // Local includes
20 #include "libmesh/side.h"
21 #include "libmesh/cell_prism6.h"
22 #include "libmesh/edge_edge2.h"
23 #include "libmesh/face_quad4.h"
24 #include "libmesh/face_tri3.h"
26 #include "libmesh/enum_order.h"
27 
28 namespace libMesh
29 {
30 
31 
32 
33 // ------------------------------------------------------------
34 // Prism6 class static member initializations
35 const int Prism6::num_nodes;
36 const int Prism6::num_sides;
37 const int Prism6::num_edges;
38 const int Prism6::num_children;
39 const int Prism6::nodes_per_side;
40 const int Prism6::nodes_per_edge;
41 
43  {
44  {0, 2, 1, 99}, // Side 0
45  {0, 1, 4, 3}, // Side 1
46  {1, 2, 5, 4}, // Side 2
47  {2, 0, 3, 5}, // Side 3
48  {3, 4, 5, 99} // Side 4
49  };
50 
52  {
53  {0, 1, 2, 3}, // Side 0
54  {0, 1, 4, 5}, // Side 1
55  {1, 2, 5, 6}, // Side 2
56  {0, 2, 4, 6}, // Side 3
57  {4, 5, 6, 7} // Side 4
58  };
59 
61  {
62  {0, 1}, // Edge 0
63  {1, 2}, // Edge 1
64  {0, 2}, // Edge 2
65  {0, 3}, // Edge 3
66  {1, 4}, // Edge 4
67  {2, 5}, // Edge 5
68  {3, 4}, // Edge 6
69  {4, 5}, // Edge 7
70  {3, 5} // Edge 8
71  };
72 
73 
74 // ------------------------------------------------------------
75 // Prism6 class member functions
76 
77 bool Prism6::is_vertex(const unsigned int) const
78 {
79  return true;
80 }
81 
82 bool Prism6::is_edge(const unsigned int) const
83 {
84  return false;
85 }
86 
87 bool Prism6::is_face(const unsigned int) const
88 {
89  return false;
90 }
91 
92 bool Prism6::is_node_on_side(const unsigned int n,
93  const unsigned int s) const
94 {
95  libmesh_assert_less (s, n_sides());
96  return std::find(std::begin(side_nodes_map[s]),
98  n) != std::end(side_nodes_map[s]);
99 }
100 
101 std::vector<unsigned>
102 Prism6::nodes_on_side(const unsigned int s) const
103 {
104  libmesh_assert_less(s, n_sides());
105  auto trim = (s > 0 && s < 4) ? 0 : 1;
106  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s]) - trim};
107 }
108 
109 bool Prism6::is_node_on_edge(const unsigned int n,
110  const unsigned int e) const
111 {
112  libmesh_assert_less (e, n_edges());
113  return std::find(std::begin(edge_nodes_map[e]),
115  n) != std::end(edge_nodes_map[e]);
116 }
117 
118 
119 
121 {
122  // Make sure z edges are affine
123  Point v = this->point(3) - this->point(0);
124  if (!v.relative_fuzzy_equals(this->point(4) - this->point(1)) ||
125  !v.relative_fuzzy_equals(this->point(5) - this->point(2)))
126  return false;
127  return true;
128 }
129 
130 
131 
133 {
134  return FIRST;
135 }
136 
137 
138 
139 std::unique_ptr<Elem> Prism6::build_side_ptr (const unsigned int i,
140  bool proxy)
141 {
142  libmesh_assert_less (i, this->n_sides());
143 
144  if (proxy)
145  {
146  switch(i)
147  {
148  case 0:
149  case 4:
150  return libmesh_make_unique<Side<Tri3,Prism6>>(this,i);
151 
152  case 1:
153  case 2:
154  case 3:
155  return libmesh_make_unique<Side<Quad4,Prism6>>(this,i);
156 
157  default:
158  libmesh_error_msg("Invalid side i = " << i);
159  }
160  }
161 
162  else
163  {
164  // Return value
165  std::unique_ptr<Elem> face;
166 
167  switch (i)
168  {
169  case 0: // the triangular face at z=-1
170  case 4: // the triangular face at z=1
171  {
172  face = libmesh_make_unique<Tri3>();
173  break;
174  }
175  case 1: // the quad face at y=0
176  case 2: // the other quad face
177  case 3: // the quad face at x=0
178  {
179  face = libmesh_make_unique<Quad4>();
180  break;
181  }
182  default:
183  libmesh_error_msg("Invalid side i = " << i);
184  }
185 
186  face->subdomain_id() = this->subdomain_id();
187 
188  // Set the nodes
189  for (unsigned n=0; n<face->n_nodes(); ++n)
190  face->set_node(n) = this->node_ptr(Prism6::side_nodes_map[i][n]);
191 
192  return face;
193  }
194 }
195 
196 
197 
198 void Prism6::build_side_ptr (std::unique_ptr<Elem> & side,
199  const unsigned int i)
200 {
201  this->side_ptr(side, i);
202 }
203 
204 
205 
206 std::unique_ptr<Elem> Prism6::build_edge_ptr (const unsigned int i)
207 {
208  libmesh_assert_less (i, this->n_edges());
209 
210  return libmesh_make_unique<SideEdge<Edge2,Prism6>>(this,i);
211 }
212 
213 
214 
215 void Prism6::connectivity(const unsigned int libmesh_dbg_var(sc),
216  const IOPackage iop,
217  std::vector<dof_id_type> & conn) const
218 {
219  libmesh_assert(_nodes);
220  libmesh_assert_less (sc, this->n_sub_elem());
221  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
222 
223  switch (iop)
224  {
225  case TECPLOT:
226  {
227  conn.resize(8);
228  conn[0] = this->node_id(0)+1;
229  conn[1] = this->node_id(1)+1;
230  conn[2] = this->node_id(2)+1;
231  conn[3] = this->node_id(2)+1;
232  conn[4] = this->node_id(3)+1;
233  conn[5] = this->node_id(4)+1;
234  conn[6] = this->node_id(5)+1;
235  conn[7] = this->node_id(5)+1;
236  return;
237  }
238 
239  case VTK:
240  {
241  conn.resize(6);
242  conn[0] = this->node_id(0);
243  conn[1] = this->node_id(2);
244  conn[2] = this->node_id(1);
245  conn[3] = this->node_id(3);
246  conn[4] = this->node_id(5);
247  conn[5] = this->node_id(4);
248  return;
249  }
250 
251  default:
252  libmesh_error_msg("Unsupported IO package " << iop);
253  }
254 }
255 
256 
257 
258 #ifdef LIBMESH_ENABLE_AMR
259 
261  {
262  // embedding matrix for child 0
263  {
264  // 0 1 2 3 4 5
265  { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
266  { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 1
267  { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2
268  { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 3
269  { .25, .25, 0.0, .25, .25, 0.0}, // 4
270  { .25, 0.0, .25, .25, 0.0, .25} // 5
271  },
272 
273  // embedding matrix for child 1
274  {
275  // 0 1 2 3 4 5
276  { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0
277  { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1
278  { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 2
279  { .25, .25, 0.0, .25, .25, 0.0}, // 3
280  { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 4
281  { 0.0, .25, .25, 0.0, .25, .25} // 5
282  },
283 
284  // embedding matrix for child 2
285  {
286  // 0 1 2 3 4 5
287  { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 0
288  { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1
289  { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2
290  { .25, 0.0, .25, .25, 0.0, .25}, // 3
291  { 0.0, .25, .25, 0.0, .25, .25}, // 4
292  { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5} // 5
293  },
294 
295  // embedding matrix for child 3
296  {
297  // 0 1 2 3 4 5
298  { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0
299  { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1
300  { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2
301  { .25, .25, 0.0, .25, .25, 0.0}, // 3
302  { 0.0, .25, .25, 0.0, .25, .25}, // 4
303  { .25, 0.0, .25, .25, 0.0, .25} // 5
304  },
305 
306  // embedding matrix for child 4
307  {
308  // 0 1 2 3 4 5
309  { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 0
310  { .25, .25, 0.0, .25, .25, 0.0}, // 1
311  { .25, 0.0, .25, .25, 0.0, .25}, // 2
312  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3
313  { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 4
314  { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5
315  },
316 
317  // embedding matrix for child 5
318  {
319  // 0 1 2 3 4 5
320  { .25, .25, 0.0, .25, .25, 0.0}, // 0
321  { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 1
322  { 0.0, .25, .25, 0.0, .25, .25}, // 2
323  { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3
324  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4
325  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5} // 5
326  },
327 
328  // embedding matrix for child 6
329  {
330  // 0 1 2 3 4 5
331  { .25, 0.0, .25, .25, 0.0, .25}, // 0
332  { 0.0, .25, .25, 0.0, .25, .25}, // 1
333  { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}, // 2
334  { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}, // 3
335  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4
336  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 5
337  },
338 
339  // embedding matrix for child 7
340  {
341  // 0 1 2 3 4 5
342  { .25, .25, 0.0, .25, .25, 0.0}, // 0
343  { 0.0, .25, .25, 0.0, .25, .25}, // 1
344  { .25, 0.0, .25, .25, 0.0, .25}, // 2
345  { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3
346  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4
347  { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5
348  }
349  };
350 
351 #endif
352 
353 
354 
356 {
357  // Make copies of our points. It makes the subsequent calculations a bit
358  // shorter and avoids dereferencing the same pointer multiple times.
359  Point
360  x0 = point(0), x1 = point(1), x2 = point(2),
361  x3 = point(3), x4 = point(4), x5 = point(5);
362 
363  // constant and zeta terms only. These are copied directly from a
364  // Python script.
365  Point dx_dxi[2] =
366  {
367  -x0/2 + x1/2 - x3/2 + x4/2, // constant
368  x0/2 - x1/2 - x3/2 + x4/2, // zeta
369  };
370 
371  // constant and zeta terms only. These are copied directly from a
372  // Python script.
373  Point dx_deta[2] =
374  {
375  -x0/2 + x2/2 - x3/2 + x5/2, // constant
376  x0/2 - x2/2 - x3/2 + x5/2, // zeta
377  };
378 
379  // Constant, xi, and eta terms
380  Point dx_dzeta[3] =
381  {
382  -x0/2 + x3/2, // constant
383  x0/2 - x2/2 - x3/2 + x5/2, // eta
384  x0/2 - x1/2 - x3/2 + x4/2 // xi
385  };
386 
387  // The quadrature rule the Prism6 is a tensor product between a
388  // four-point TRI3 rule (in xi, eta) and a two-point EDGE2 rule (in
389  // zeta) which is capable of integrating cubics exactly.
390 
391  // Number of points in the 2D quadrature rule.
392  const int N2D = 4;
393 
394  static const Real w2D[N2D] =
395  {
396  Real(1.5902069087198858469718450103758e-01L),
397  Real(9.0979309128011415302815498962418e-02L),
398  Real(1.5902069087198858469718450103758e-01L),
399  Real(9.0979309128011415302815498962418e-02L)
400  };
401 
402  static const Real xi[N2D] =
403  {
404  Real(1.5505102572168219018027159252941e-01L),
405  Real(6.4494897427831780981972840747059e-01L),
406  Real(1.5505102572168219018027159252941e-01L),
407  Real(6.4494897427831780981972840747059e-01L)
408  };
409 
410  static const Real eta[N2D] =
411  {
412  Real(1.7855872826361642311703513337422e-01L),
413  Real(7.5031110222608118177475598324603e-02L),
414  Real(6.6639024601470138670269327409637e-01L),
415  Real(2.8001991549907407200279599420481e-01L)
416  };
417 
418  // Number of points in the 1D quadrature rule. The weights of the
419  // 1D quadrature rule are equal to 1.
420  const int N1D = 2;
421 
422  // Points of the 1D quadrature rule
423  static const Real zeta[N1D] =
424  {
425  -std::sqrt(3.)/3,
426  std::sqrt(3.)/3.
427  };
428 
429  Real vol = 0.;
430  for (int i=0; i<N2D; ++i)
431  {
432  // dx_dzeta depends only on the 2D quadrature rule points.
433  Point dx_dzeta_q = dx_dzeta[0] + eta[i]*dx_dzeta[1] + xi[i]*dx_dzeta[2];
434 
435  for (int j=0; j<N1D; ++j)
436  {
437  // dx_dxi and dx_deta only depend on the 1D quadrature rule points.
438  Point
439  dx_dxi_q = dx_dxi[0] + zeta[j]*dx_dxi[1],
440  dx_deta_q = dx_deta[0] + zeta[j]*dx_deta[1];
441 
442  // Compute scalar triple product, multiply by weight, and accumulate volume.
443  vol += w2D[i] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q);
444  }
445  }
446 
447  return vol;
448 }
449 
452 {
453  return Elem::loose_bounding_box();
454 }
455 
456 } // namespace libMesh
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: cell_prism6.h:159
Node ** _nodes
Definition: elem.h:1695
virtual Real volume() const override
Definition: cell_prism6.C:355
virtual unsigned int n_sides() const override final
Definition: cell_prism.h:79
unsigned short int side
Definition: xdr_io.C:50
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: cell_prism.C:106
virtual BoundingBox loose_bounding_box() const
Definition: elem.C:2857
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Definition: cell_prism6.h:205
virtual bool is_edge(const unsigned int i) const override
Definition: cell_prism6.C:82
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: cell_prism6.C:215
virtual bool is_face(const unsigned int i) const override
Definition: cell_prism6.C:87
static const int num_edges
Definition: cell_prism6.h:150
IterBase * end
virtual BoundingBox loose_bounding_box() const override
Definition: cell_prism6.C:451
virtual Order default_order() const override
Definition: cell_prism6.C:132
virtual bool is_vertex(const unsigned int i) const override
Definition: cell_prism6.C:77
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1054
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: cell_prism6.C:109
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy) override
Definition: cell_prism6.C:139
virtual unsigned int n_sub_elem() const override
Definition: cell_prism6.h:78
static const int nodes_per_side
Definition: cell_prism6.h:152
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
Definition: elem.h:2034
virtual unsigned int n_edges() const override final
Definition: cell_prism.h:89
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1957
static const int nodes_per_edge
Definition: cell_prism6.h:153
static const int num_children
Definition: cell_prism6.h:151
static const int num_sides
Definition: cell_prism6.h:149
virtual bool has_affine_map() const override
Definition: cell_prism6.C:120
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
Definition: cell_prism6.h:170
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Definition: cell_prism6.C:206
static const unsigned int side_elems_map[num_sides][nodes_per_side]
Definition: cell_prism6.h:164
A geometric point in (x,y,z) space.
Definition: point.h:38
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:1914
const Point & point(const unsigned int i) const
Definition: elem.h:1892
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:990
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: cell_prism6.C:92
static const int num_nodes
Definition: cell_prism6.h:148
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: cell_prism6.C:102