face_tri6.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 // C++ includes
19 
20 // Local includes
21 #include "libmesh/side.h"
22 #include "libmesh/edge_edge3.h"
23 #include "libmesh/face_tri6.h"
24 
25 namespace libMesh
26 {
27 
28 
29 
30 
31 // ------------------------------------------------------------
32 // Tri6 class static member initializations
33 const unsigned int Tri6::side_nodes_map[3][3] =
34  {
35  {0, 1, 3}, // Side 0
36  {1, 2, 4}, // Side 1
37  {2, 0, 5} // Side 2
38  };
39 
40 
41 #ifdef LIBMESH_ENABLE_AMR
42 
43 const float Tri6::_embedding_matrix[4][6][6] =
44  {
45  // embedding matrix for child 0
46  {
47  // 0 1 2 3 4 5
48  { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
49  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 1
50  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 2
51  {.375, -.125, 0.0, .75, 0.0, 0.0}, // 3
52  { 0.0, -.125, -.125, 0.5, .25, 0.5}, // 4
53  {.375, 0.0, -.125, 0.0, 0.0, .75} // 5
54  },
55 
56  // embedding matrix for child 1
57  {
58  // 0 1 2 3 4 5
59  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 0
60  { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1
61  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 2
62  {-.125, .375, 0.0, .75, 0.0, 0.0}, // 3
63  { 0.0, .375, -.125, 0.0, .75, 0.0}, // 4
64  {-.125, 0.0, -.125, 0.5, 0.5, .25} // 5
65  },
66 
67  // embedding matrix for child 2
68  {
69  // 0 1 2 3 4 5
70  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 0
71  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 1
72  { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2
73  {-.125, -.125, 0.0, .25, 0.5, 0.5}, // 3
74  { 0.0, -.125, .375, 0.0, .75, 0.0}, // 4
75  {-.125, 0.0, .375, 0.0, 0.0, .75} // 5
76  },
77 
78  // embedding matrix for child 3
79  {
80  // 0 1 2 3 4 5
81  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 0
82  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 1
83  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 2
84  {-.125, 0.0, -.125, 0.5, 0.5, .25}, // 3
85  {-.125, -.125, 0.0, .25, 0.5, 0.5}, // 4
86  { 0.0, -.125, -.125, 0.5, .25, 0.5} // 5
87  }
88  };
89 
90 #endif
91 
92 
93 
94 // ------------------------------------------------------------
95 // Tri6 class member functions
96 
97 bool Tri6::is_vertex(const unsigned int i) const
98 {
99  if (i < 3)
100  return true;
101  return false;
102 }
103 
104 bool Tri6::is_edge(const unsigned int i) const
105 {
106  if (i < 3)
107  return false;
108  return true;
109 }
110 
111 bool Tri6::is_face(const unsigned int) const
112 {
113  return false;
114 }
115 
116 bool Tri6::is_node_on_side(const unsigned int n,
117  const unsigned int s) const
118 {
119  libmesh_assert_less (s, n_sides());
120  for (unsigned int i = 0; i != 3; ++i)
121  if (side_nodes_map[s][i] == n)
122  return true;
123  return false;
124 }
125 
126 
127 
129 {
130  // Make sure edges are straight
131  if (!this->point(3).relative_fuzzy_equals
132  ((this->point(0) + this->point(1))/2.))
133  return false;
134  if (!this->point(4).relative_fuzzy_equals
135  ((this->point(1) + this->point(2))/2.))
136  return false;
137  if (!this->point(5).relative_fuzzy_equals
138  ((this->point(2) + this->point(0))/2.))
139  return false;
140 
141  return true;
142 }
143 
144 
145 
146 dof_id_type Tri6::key (const unsigned int s) const
147 {
148  libmesh_assert_less (s, this->n_sides());
149 
150  switch (s)
151  {
152  case 0:
153 
154  return
155  this->compute_key (this->node_id(3));
156 
157  case 1:
158 
159  return
160  this->compute_key (this->node_id(4));
161 
162  case 2:
163 
164  return
165  this->compute_key (this->node_id(5));
166 
167  default:
168  libmesh_error_msg("Invalid side s = " << s);
169  }
170 }
171 
172 
173 
174 unsigned int Tri6::which_node_am_i(unsigned int side,
175  unsigned int side_node) const
176 {
177  libmesh_assert_less (side, this->n_sides());
178  libmesh_assert_less (side_node, 3);
179 
180  return Tri6::side_nodes_map[side][side_node];
181 }
182 
183 
184 
185 std::unique_ptr<Elem> Tri6::build_side_ptr (const unsigned int i,
186  bool proxy)
187 {
188  libmesh_assert_less (i, this->n_sides());
189 
190  if (proxy)
191  return libmesh_make_unique<Side<Edge3,Tri6>>(this,i);
192 
193  else
194  {
195  std::unique_ptr<Elem> edge = libmesh_make_unique<Edge3>();
196  edge->subdomain_id() = this->subdomain_id();
197 
198  // Set the nodes
199  for (unsigned n=0; n<edge->n_nodes(); ++n)
200  edge->set_node(n) = this->node_ptr(Tri6::side_nodes_map[i][n]);
201 
202  return edge;
203  }
204 }
205 
206 
207 void Tri6::connectivity(const unsigned int sf,
208  const IOPackage iop,
209  std::vector<dof_id_type> & conn) const
210 {
211  libmesh_assert_less (sf, this->n_sub_elem());
212  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
213 
214  switch (iop)
215  {
216  case TECPLOT:
217  {
218  conn.resize(4);
219  switch(sf)
220  {
221  case 0:
222  // linear sub-triangle 0
223  conn[0] = this->node_id(0)+1;
224  conn[1] = this->node_id(3)+1;
225  conn[2] = this->node_id(5)+1;
226  conn[3] = this->node_id(5)+1;
227 
228  return;
229 
230  case 1:
231  // linear sub-triangle 1
232  conn[0] = this->node_id(3)+1;
233  conn[1] = this->node_id(1)+1;
234  conn[2] = this->node_id(4)+1;
235  conn[3] = this->node_id(4)+1;
236 
237  return;
238 
239  case 2:
240  // linear sub-triangle 2
241  conn[0] = this->node_id(5)+1;
242  conn[1] = this->node_id(4)+1;
243  conn[2] = this->node_id(2)+1;
244  conn[3] = this->node_id(2)+1;
245 
246  return;
247 
248  case 3:
249  // linear sub-triangle 3
250  conn[0] = this->node_id(3)+1;
251  conn[1] = this->node_id(4)+1;
252  conn[2] = this->node_id(5)+1;
253  conn[3] = this->node_id(5)+1;
254 
255  return;
256 
257  default:
258  libmesh_error_msg("Invalid sf = " << sf);
259  }
260  }
261 
262  case VTK:
263  {
264  // VTK_QUADRATIC_TRIANGLE has same numbering as libmesh TRI6
265  conn.resize(6);
266  conn[0] = this->node_id(0);
267  conn[1] = this->node_id(1);
268  conn[2] = this->node_id(2);
269  conn[3] = this->node_id(3);
270  conn[4] = this->node_id(4);
271  conn[5] = this->node_id(5);
272  return;
273 
274  // Used to write out linear sub-triangles for VTK...
275  /*
276  conn.resize(3);
277  switch(sf)
278  {
279  case 0:
280  // linear sub-triangle 0
281  conn[0] = this->node_id(0);
282  conn[1] = this->node_id(3);
283  conn[2] = this->node_id(5);
284 
285  return;
286 
287  case 1:
288  // linear sub-triangle 1
289  conn[0] = this->node_id(3);
290  conn[1] = this->node_id(1);
291  conn[2] = this->node_id(4);
292 
293  return;
294 
295  case 2:
296  // linear sub-triangle 2
297  conn[0] = this->node_id(5);
298  conn[1] = this->node_id(4);
299  conn[2] = this->node_id(2);
300 
301  return;
302 
303  case 3:
304  // linear sub-triangle 3
305  conn[0] = this->node_id(3);
306  conn[1] = this->node_id(4);
307  conn[2] = this->node_id(5);
308 
309  return;
310 
311  default:
312  libmesh_error_msg("Invalid sf = " << sf);
313  }
314  */
315  }
316 
317  default:
318  libmesh_error_msg("Unsupported IO package " << iop);
319  }
320 }
321 
322 
323 
325 {
326  // This might have curved edges, or might be a curved surface in
327  // 3-space, in which case the full bounding box can be larger than
328  // the bounding box of just the nodes.
329  //
330  //
331  // FIXME - I haven't yet proven the formula below to be correct for
332  // quadratics in 2D - RHS
333  Point pmin, pmax;
334 
335  for (unsigned d=0; d<LIBMESH_DIM; ++d)
336  {
337  Real center = this->point(0)(d);
338  for (unsigned int p=1; p != 6; ++p)
339  center += this->point(p)(d);
340  center /= 6;
341 
342  Real hd = std::abs(center - this->point(0)(d));
343  for (unsigned int p=1; p != 6; ++p)
344  hd = std::max(hd, std::abs(center - this->point(p)(d)));
345 
346  pmin(d) = center - hd;
347  pmax(d) = center + hd;
348  }
349 
350  return BoundingBox(pmin, pmax);
351 }
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  // Construct constant data vectors.
364  // \vec{x}_{\xi} = \vec{a1}*xi + \vec{b1}*eta + \vec{c1}
365  // \vec{x}_{\eta} = \vec{a2}*xi + \vec{b2}*eta + \vec{c2}
366  Point
367  a1 = 4*x0 + 4*x1 - 8*x3,
368  b1 = 4*x0 - 4*x3 + 4*x4 - 4*x5, /*=a2*/
369  c1 = -3*x0 - 1*x1 + 4*x3,
370  b2 = 4*x0 + 4*x2 - 8*x5,
371  c2 = -3*x0 - 1*x2 + 4*x5;
372 
373  // If a1 == b1 == a2 == b2 == 0, this is a TRI6 with straight sides,
374  // and we can use the TRI3 formula to compute the volume.
375  if (a1.relative_fuzzy_equals(Point(0,0,0)) &&
376  b1.relative_fuzzy_equals(Point(0,0,0)) &&
377  b2.relative_fuzzy_equals(Point(0,0,0)))
378  return 0.5 * cross_norm(c1, c2);
379 
380  // 7-point rule, exact for quintics.
381  const unsigned int N = 7;
382 
383  // Parameters of the quadrature rule
384  const static Real
385  w1 = Real(31)/480 + std::sqrt(15.0L)/2400,
386  w2 = Real(31)/480 - std::sqrt(15.0L)/2400,
387  q1 = Real(2)/7 + std::sqrt(15.0L)/21,
388  q2 = Real(2)/7 - std::sqrt(15.0L)/21;
389 
390  const static Real xi[N] = {Real(1)/3, q1, q1, 1-2*q1, q2, q2, 1-2*q2};
391  const static Real eta[N] = {Real(1)/3, q1, 1-2*q1, q1, q2, 1-2*q2, q2};
392  const static Real wts[N] = {Real(9)/80, w1, w1, w1, w2, w2, w2};
393 
394  // Approximate the area with quadrature
395  Real vol=0.;
396  for (unsigned int q=0; q<N; ++q)
397  vol += wts[q] * cross_norm(xi[q]*a1 + eta[q]*b1 + c1,
398  xi[q]*b1 + eta[q]*b2 + c2);
399 
400  return vol;
401 }
402 
403 
404 
405 unsigned short int Tri6::second_order_adjacent_vertex (const unsigned int n,
406  const unsigned int v) const
407 {
408  libmesh_assert_greater_equal (n, this->n_vertices());
409  libmesh_assert_less (n, this->n_nodes());
410  libmesh_assert_less (v, 2);
411  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
412 }
413 
414 
415 
416 const unsigned short int Tri6::_second_order_adjacent_vertices[3][2] =
417  {
418  {0, 1}, // vertices adjacent to node 3
419  {1, 2}, // vertices adjacent to node 4
420  {0, 2} // vertices adjacent to node 5
421  };
422 
423 
424 
425 std::pair<unsigned short int, unsigned short int>
426 Tri6::second_order_child_vertex (const unsigned int n) const
427 {
428  libmesh_assert_greater_equal (n, this->n_vertices());
429  libmesh_assert_less (n, this->n_nodes());
430  return std::pair<unsigned short int, unsigned short int>
433 }
434 
435 
436 
437 const unsigned short int Tri6::_second_order_vertex_child_number[6] =
438  {
439  99,99,99, // Vertices
440  0,1,0 // Edges
441  };
442 
443 
444 
445 const unsigned short int Tri6::_second_order_vertex_child_index[6] =
446  {
447  99,99,99, // Vertices
448  1,2,2 // Edges
449  };
450 
451 } // namespace libMesh
double abs(double a)
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1092
virtual unsigned int n_sub_elem() const libmesh_override
Definition: face_tri6.h:73
virtual bool is_edge(const unsigned int i) const libmesh_override
Definition: face_tri6.C:104
long double max(long double a, double b)
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const libmesh_override
Definition: face_tri6.C:405
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const libmesh_override
Definition: face_tri6.C:207
virtual bool has_affine_map() const libmesh_override
Definition: face_tri6.C:128
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1875
virtual unsigned int n_nodes() const libmesh_override
Definition: face_tri6.h:68
virtual bool is_face(const unsigned int i) const libmesh_override
Definition: face_tri6.C:111
virtual unsigned int n_vertices() const libmesh_override
Definition: face_tri.h:91
subdomain_id_type subdomain_id() const
Definition: elem.h:1952
virtual unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const libmesh_override
Definition: face_tri6.C:174
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const libmesh_override
Definition: face_tri6.C:116
virtual BoundingBox loose_bounding_box() const libmesh_override
Definition: face_tri6.C:324
static const unsigned short int _second_order_vertex_child_number[6]
Definition: face_tri6.h:226
std::unique_ptr< Elem > side(const unsigned int i) const
Definition: elem.h:2106
virtual unsigned int n_sides() const libmesh_override
Definition: face_tri.h:86
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual dof_id_type key() const libmesh_override
Definition: face_tri.C:73
const Point & point(const unsigned int i) const
Definition: elem.h:1810
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy) libmesh_override
Definition: face_tri6.C:185
static const unsigned short int _second_order_adjacent_vertices[3][2]
Definition: face_tri6.h:221
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const libmesh_override
Definition: face_tri6.C:426
static const unsigned int side_nodes_map[3][3]
Definition: face_tri6.h:172
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:1832
static const unsigned short int _second_order_vertex_child_index[6]
Definition: face_tri6.h:231
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2622
virtual Real volume() const libmesh_override
Definition: face_tri6.C:355
static const float _embedding_matrix[4][6][6]
Definition: face_tri6.h:209
A geometric point in (x,y,z) space.
Definition: point.h:38
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:987
virtual bool is_vertex(const unsigned int i) const libmesh_override
Definition: face_tri6.C:97
uint8_t dof_id_type
Definition: id_types.h:64