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