face_tri6.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 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  libmesh_error_msg("We'll never get here!");
172  return 0;
173 }
174 
175 
176 
177 unsigned int Tri6::which_node_am_i(unsigned int side,
178  unsigned int side_node) const
179 {
180  libmesh_assert_less (side, this->n_sides());
181  libmesh_assert_less (side_node, 3);
182 
183  return Tri6::side_nodes_map[side][side_node];
184 }
185 
186 
187 
188 UniquePtr<Elem> Tri6::build_side_ptr (const unsigned int i,
189  bool proxy)
190 {
191  libmesh_assert_less (i, this->n_sides());
192 
193  if (proxy)
194  return UniquePtr<Elem>(new Side<Edge3,Tri6>(this,i));
195 
196  else
197  {
198  Elem * edge = new Edge3;
199  edge->subdomain_id() = this->subdomain_id();
200 
201  // Set the nodes
202  for (unsigned n=0; n<edge->n_nodes(); ++n)
203  edge->set_node(n) = this->node_ptr(Tri6::side_nodes_map[i][n]);
204 
205  return UniquePtr<Elem>(edge);
206  }
207 
208  libmesh_error_msg("We'll never get here!");
209  return UniquePtr<Elem>();
210 }
211 
212 
213 void Tri6::connectivity(const unsigned int sf,
214  const IOPackage iop,
215  std::vector<dof_id_type> & conn) const
216 {
217  libmesh_assert_less (sf, this->n_sub_elem());
218  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
219 
220  switch (iop)
221  {
222  case TECPLOT:
223  {
224  conn.resize(4);
225  switch(sf)
226  {
227  case 0:
228  // linear sub-triangle 0
229  conn[0] = this->node_id(0)+1;
230  conn[1] = this->node_id(3)+1;
231  conn[2] = this->node_id(5)+1;
232  conn[3] = this->node_id(5)+1;
233 
234  return;
235 
236  case 1:
237  // linear sub-triangle 1
238  conn[0] = this->node_id(3)+1;
239  conn[1] = this->node_id(1)+1;
240  conn[2] = this->node_id(4)+1;
241  conn[3] = this->node_id(4)+1;
242 
243  return;
244 
245  case 2:
246  // linear sub-triangle 2
247  conn[0] = this->node_id(5)+1;
248  conn[1] = this->node_id(4)+1;
249  conn[2] = this->node_id(2)+1;
250  conn[3] = this->node_id(2)+1;
251 
252  return;
253 
254  case 3:
255  // linear sub-triangle 3
256  conn[0] = this->node_id(3)+1;
257  conn[1] = this->node_id(4)+1;
258  conn[2] = this->node_id(5)+1;
259  conn[3] = this->node_id(5)+1;
260 
261  return;
262 
263  default:
264  libmesh_error_msg("Invalid sf = " << sf);
265  }
266  }
267 
268  case VTK:
269  {
270  // VTK_QUADRATIC_TRIANGLE has same numbering as libmesh TRI6
271  conn.resize(6);
272  conn[0] = this->node_id(0);
273  conn[1] = this->node_id(1);
274  conn[2] = this->node_id(2);
275  conn[3] = this->node_id(3);
276  conn[4] = this->node_id(4);
277  conn[5] = this->node_id(5);
278  return;
279 
280  // Used to write out linear sub-triangles for VTK...
281  /*
282  conn.resize(3);
283  switch(sf)
284  {
285  case 0:
286  // linear sub-triangle 0
287  conn[0] = this->node_id(0);
288  conn[1] = this->node_id(3);
289  conn[2] = this->node_id(5);
290 
291  return;
292 
293  case 1:
294  // linear sub-triangle 1
295  conn[0] = this->node_id(3);
296  conn[1] = this->node_id(1);
297  conn[2] = this->node_id(4);
298 
299  return;
300 
301  case 2:
302  // linear sub-triangle 2
303  conn[0] = this->node_id(5);
304  conn[1] = this->node_id(4);
305  conn[2] = this->node_id(2);
306 
307  return;
308 
309  case 3:
310  // linear sub-triangle 3
311  conn[0] = this->node_id(3);
312  conn[1] = this->node_id(4);
313  conn[2] = this->node_id(5);
314 
315  return;
316 
317  default:
318  libmesh_error_msg("Invalid sf = " << sf);
319  }
320  */
321  }
322 
323  default:
324  libmesh_error_msg("Unsupported IO package " << iop);
325  }
326 }
327 
328 
329 
331 {
332  // This might have curved edges, or might be a curved surface in
333  // 3-space, in which case the full bounding box can be larger than
334  // the bounding box of just the nodes.
335  //
336  //
337  // FIXME - I haven't yet proven the formula below to be correct for
338  // quadratics in 2D - RHS
339  Point pmin, pmax;
340 
341  for (unsigned d=0; d<LIBMESH_DIM; ++d)
342  {
343  Real center = this->point(0)(d);
344  for (unsigned int p=1; p != 6; ++p)
345  center += this->point(p)(d);
346  center /= 6;
347 
348  Real hd = std::abs(center - this->point(0)(d));
349  for (unsigned int p=1; p != 6; ++p)
350  hd = std::max(hd, std::abs(center - this->point(p)(d)));
351 
352  pmin(d) = center - hd;
353  pmax(d) = center + hd;
354  }
355 
356  return BoundingBox(pmin, pmax);
357 }
358 
359 
360 
362 {
363  // Make copies of our points. It makes the subsequent calculations a bit
364  // shorter and avoids dereferencing the same pointer multiple times.
365  Point
366  x0 = point(0), x1 = point(1), x2 = point(2),
367  x3 = point(3), x4 = point(4), x5 = point(5);
368 
369  // Construct constant data vectors.
370  // \vec{x}_{\xi} = \vec{a1}*xi + \vec{b1}*eta + \vec{c1}
371  // \vec{x}_{\eta} = \vec{a2}*xi + \vec{b2}*eta + \vec{c2}
372  Point
373  a1 = 4*x0 + 4*x1 - 8*x3,
374  b1 = 4*x0 - 4*x3 + 4*x4 - 4*x5, /*=a2*/
375  c1 = -3*x0 - 1*x1 + 4*x3,
376  b2 = 4*x0 + 4*x2 - 8*x5,
377  c2 = -3*x0 - 1*x2 + 4*x5;
378 
379  // If a1 == b1 == a2 == b2 == 0, this is a TRI6 with straight sides,
380  // and we can use the TRI3 formula to compute the volume.
381  if (a1.relative_fuzzy_equals(Point(0,0,0)) &&
382  b1.relative_fuzzy_equals(Point(0,0,0)) &&
383  b2.relative_fuzzy_equals(Point(0,0,0)))
384  return 0.5 * cross_norm(c1, c2);
385 
386  // 7-point rule, exact for quintics.
387  const unsigned int N = 7;
388 
389  // Parameters of the quadrature rule
390  const static Real
391  w1 = Real(31)/480 + std::sqrt(15.0L)/2400,
392  w2 = Real(31)/480 - std::sqrt(15.0L)/2400,
393  q1 = Real(2)/7 + std::sqrt(15.0L)/21,
394  q2 = Real(2)/7 - std::sqrt(15.0L)/21;
395 
396  const static Real xi[N] = {Real(1)/3, q1, q1, 1-2*q1, q2, q2, 1-2*q2};
397  const static Real eta[N] = {Real(1)/3, q1, 1-2*q1, q1, q2, 1-2*q2, q2};
398  const static Real wts[N] = {Real(9)/80, w1, w1, w1, w2, w2, w2};
399 
400  // Approximate the area with quadrature
401  Real vol=0.;
402  for (unsigned int q=0; q<N; ++q)
403  vol += wts[q] * cross_norm(xi[q]*a1 + eta[q]*b1 + c1,
404  xi[q]*b1 + eta[q]*b2 + c2);
405 
406  return vol;
407 }
408 
409 
410 
411 unsigned short int Tri6::second_order_adjacent_vertex (const unsigned int n,
412  const unsigned int v) const
413 {
414  libmesh_assert_greater_equal (n, this->n_vertices());
415  libmesh_assert_less (n, this->n_nodes());
416  libmesh_assert_less (v, 2);
417  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
418 }
419 
420 
421 
422 const unsigned short int Tri6::_second_order_adjacent_vertices[3][2] =
423  {
424  {0, 1}, // vertices adjacent to node 3
425  {1, 2}, // vertices adjacent to node 4
426  {0, 2} // vertices adjacent to node 5
427  };
428 
429 
430 
431 std::pair<unsigned short int, unsigned short int>
432 Tri6::second_order_child_vertex (const unsigned int n) const
433 {
434  libmesh_assert_greater_equal (n, this->n_vertices());
435  libmesh_assert_less (n, this->n_nodes());
436  return std::pair<unsigned short int, unsigned short int>
439 }
440 
441 
442 
443 const unsigned short int Tri6::_second_order_vertex_child_number[6] =
444  {
445  99,99,99, // Vertices
446  0,1,0 // Edges
447  };
448 
449 
450 
451 const unsigned short int Tri6::_second_order_vertex_child_index[6] =
452  {
453  99,99,99, // Vertices
454  1,2,2 // Edges
455  };
456 
457 } // namespace libMesh
double abs(double a)
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:1789
UniquePtr< Elem > side(const unsigned int i) const
Definition: elem.h:1906
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1084
virtual UniquePtr< Elem > build_side_ptr(const unsigned int i, bool proxy) libmesh_override
Definition: face_tri6.C:188
virtual unsigned int n_sub_elem() const libmesh_override
Definition: face_tri6.h:73
The base class for all geometric element types.
Definition: elem.h:86
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:411
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const libmesh_override
Definition: face_tri6.C:213
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
virtual unsigned int n_nodes() const =0
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:1724
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:1799
virtual unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const libmesh_override
Definition: face_tri6.C:177
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const libmesh_override
Definition: face_tri6.C:116
A 1D geometric element with 3 nodes.
Definition: edge_edge3.h:43
virtual BoundingBox loose_bounding_box() const libmesh_override
Definition: face_tri6.C:330
Proxy class for efficiently representing an Elem&#39;s side.
Definition: side.h:48
static const unsigned short int _second_order_vertex_child_number[6]
Definition: face_tri6.h:226
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:1661
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:432
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:1683
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:2405
virtual Real volume() const libmesh_override
Definition: face_tri6.C:361
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:979
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