cell_tet4.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_tet4.h"
22 #include "libmesh/edge_edge2.h"
23 #include "libmesh/face_tri3.h"
24 #include "libmesh/tensor_value.h"
26 #include "libmesh/enum_order.h"
27 
28 namespace libMesh
29 {
30 
31 
32 
33 // ------------------------------------------------------------
34 // Tet4 class static member initializations
35 const int Tet4::num_nodes;
36 const int Tet4::num_sides;
37 const int Tet4::num_edges;
38 const int Tet4::num_children;
39 const int Tet4::nodes_per_side;
40 const int Tet4::nodes_per_edge;
41 
43  {
44  {0, 2, 1}, // Side 0
45  {0, 1, 3}, // Side 1
46  {1, 2, 3}, // Side 2
47  {2, 0, 3} // Side 3
48  };
49 
51  {
52  {0, 1}, // Edge 0
53  {1, 2}, // Edge 1
54  {0, 2}, // Edge 2
55  {0, 3}, // Edge 3
56  {1, 3}, // Edge 4
57  {2, 3} // Edge 5
58  };
59 
60 
61 // ------------------------------------------------------------
62 // Tet4 class member functions
63 
64 bool Tet4::is_vertex(const unsigned int) const
65 {
66  return true;
67 }
68 
69 bool Tet4::is_edge(const unsigned int) const
70 {
71  return false;
72 }
73 
74 bool Tet4::is_face(const unsigned int) const
75 {
76  return false;
77 }
78 
79 bool Tet4::is_node_on_edge(const unsigned int n,
80  const unsigned int e) const
81 {
82  libmesh_assert_less (e, n_edges());
83  return std::find(std::begin(edge_nodes_map[e]),
85  n) != std::end(edge_nodes_map[e]);
86 }
87 
88 
89 
90 
91 #ifdef LIBMESH_ENABLE_AMR
92 
93 // This function only works if LIBMESH_ENABLE_AMR...
94 bool Tet4::is_child_on_side(const unsigned int c,
95  const unsigned int s) const
96 {
97  // OK, for the Tet4, this is pretty obvious... it is sets of nodes
98  // not equal to the current node. But if we want this algorithm to
99  // be generic and work for Tet10 also it helps to do it this way.
100  const unsigned int nodes_opposite[4][3] =
101  {
102  {1,2,3}, // nodes opposite node 0
103  {0,2,3}, // nodes opposite node 1
104  {0,1,3}, // nodes opposite node 2
105  {0,1,2} // nodes opposite node 3
106  };
107 
108  // Call the base class helper function
109  return Tet::is_child_on_side_helper(c, s, nodes_opposite);
110 }
111 
112 #else
113 
114 bool Tet4::is_child_on_side(const unsigned int /*c*/,
115  const unsigned int /*s*/) const
116 {
117  libmesh_not_implemented();
118  return false;
119 }
120 
121 #endif //LIBMESH_ENABLE_AMR
122 
123 
124 
125 
126 bool Tet4::is_node_on_side(const unsigned int n,
127  const unsigned int s) const
128 {
129  libmesh_assert_less (s, n_sides());
130  return std::find(std::begin(side_nodes_map[s]),
132  n) != std::end(side_nodes_map[s]);
133 }
134 
135 std::vector<unsigned>
136 Tet4::nodes_on_side(const unsigned int s) const
137 {
138  libmesh_assert_less(s, n_sides());
139  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
140 }
141 
143 {
144  return FIRST;
145 }
146 
147 std::unique_ptr<Elem> Tet4::build_side_ptr (const unsigned int i,
148  bool proxy)
149 {
150  libmesh_assert_less (i, this->n_sides());
151 
152  if (proxy)
153  return libmesh_make_unique<Side<Tri3,Tet4>>(this,i);
154 
155  else
156  {
157  std::unique_ptr<Elem> face = libmesh_make_unique<Tri3>();
158  face->subdomain_id() = this->subdomain_id();
159 
160  for (unsigned n=0; n<face->n_nodes(); ++n)
161  face->set_node(n) = this->node_ptr(Tet4::side_nodes_map[i][n]);
162 
163  return face;
164  }
165 }
166 
167 
168 
169 void Tet4::build_side_ptr (std::unique_ptr<Elem> & side,
170  const unsigned int i)
171 {
172  this->simple_build_side_ptr<Tet4>(side, i, TRI3);
173 }
174 
175 
176 
177 std::unique_ptr<Elem> Tet4::build_edge_ptr (const unsigned int i)
178 {
179  libmesh_assert_less (i, this->n_edges());
180 
181  return libmesh_make_unique<SideEdge<Edge2,Tet4>>(this,i);
182 }
183 
184 
185 void Tet4::connectivity(const unsigned int libmesh_dbg_var(sc),
186  const IOPackage iop,
187  std::vector<dof_id_type> & conn) const
188 {
189  libmesh_assert(_nodes);
190  libmesh_assert_less (sc, this->n_sub_elem());
191  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
192 
193 
194  switch (iop)
195  {
196  case TECPLOT:
197  {
198  conn.resize(8);
199  conn[0] = this->node_id(0)+1;
200  conn[1] = this->node_id(1)+1;
201  conn[2] = this->node_id(2)+1;
202  conn[3] = this->node_id(2)+1;
203  conn[4] = this->node_id(3)+1;
204  conn[5] = this->node_id(3)+1;
205  conn[6] = this->node_id(3)+1;
206  conn[7] = this->node_id(3)+1;
207  return;
208  }
209 
210  case VTK:
211  {
212  conn.resize(4);
213  conn[0] = this->node_id(0);
214  conn[1] = this->node_id(1);
215  conn[2] = this->node_id(2);
216  conn[3] = this->node_id(3);
217  return;
218  }
219 
220  default:
221  libmesh_error_msg("Unsupported IO package " << iop);
222  }
223 }
224 
225 
226 
227 #ifdef LIBMESH_ENABLE_AMR
228 
230  {
231  // embedding matrix for child 0
232  {
233  // 0 1 2 3
234  {1.0, 0.0, 0.0, 0.0}, // 0
235  {0.5, 0.5, 0.0, 0.0}, // 1
236  {0.5, 0.0, 0.5, 0.0}, // 2
237  {0.5, 0.0, 0.0, 0.5} // 3
238  },
239 
240  // embedding matrix for child 1
241  {
242  // 0 1 2 3
243  {0.5, 0.5, 0.0, 0.0}, // 0
244  {0.0, 1.0, 0.0, 0.0}, // 1
245  {0.0, 0.5, 0.5, 0.0}, // 2
246  {0.0, 0.5, 0.0, 0.5} // 3
247  },
248 
249  // embedding matrix for child 2
250  {
251  // 0 1 2 3
252  {0.5, 0.0, 0.5, 0.0}, // 0
253  {0.0, 0.5, 0.5, 0.0}, // 1
254  {0.0, 0.0, 1.0, 0.0}, // 2
255  {0.0, 0.0, 0.5, 0.5} // 3
256  },
257 
258  // embedding matrix for child 3
259  {
260  // 0 1 2 3
261  {0.5, 0.0, 0.0, 0.5}, // 0
262  {0.0, 0.5, 0.0, 0.5}, // 1
263  {0.0, 0.0, 0.5, 0.5}, // 2
264  {0.0, 0.0, 0.0, 1.0} // 3
265  },
266 
267  // embedding matrix for child 4
268  {
269  // 0 1 2 3
270  {0.5, 0.5, 0.0, 0.0}, // 0
271  {0.0, 0.5, 0.0, 0.5}, // 1
272  {0.5, 0.0, 0.5, 0.0}, // 2
273  {0.5, 0.0, 0.0, 0.5} // 3
274  },
275 
276  // embedding matrix for child 5
277  {
278  // 0 1 2 3
279  {0.5, 0.5, 0.0, 0.0}, // 0
280  {0.0, 0.5, 0.5, 0.0}, // 1
281  {0.5, 0.0, 0.5, 0.0}, // 2
282  {0.0, 0.5, 0.0, 0.5} // 3
283  },
284 
285  // embedding matrix for child 6
286  {
287  // 0 1 2 3
288  {0.5, 0.0, 0.5, 0.0}, // 0
289  {0.0, 0.5, 0.5, 0.0}, // 1
290  {0.0, 0.0, 0.5, 0.5}, // 2
291  {0.0, 0.5, 0.0, 0.5} // 3
292  },
293 
294  // embedding matrix for child 7
295  {
296  // 0 1 2 3
297  {0.5, 0.0, 0.5, 0.0}, // 0
298  {0.0, 0.5, 0.0, 0.5}, // 1
299  {0.0, 0.0, 0.5, 0.5}, // 2
300  {0.5, 0.0, 0.0, 0.5} // 3
301  }
302  };
303 
304 #endif // #ifdef LIBMESH_ENABLE_AMR
305 
306 
307 
308 
309 
311 {
312  // The volume of a tetrahedron is 1/6 the box product formed
313  // by its base and apex vectors
314  Point a = point(3) - point(0);
315 
316  // b is the vector pointing from 0 to 1
317  Point b = point(1) - point(0);
318 
319  // c is the vector pointing from 0 to 2
320  Point c = point(2) - point(0);
321 
322  return triple_product(a, b, c) / 6.;
323 }
324 
325 
326 
327 
328 std::pair<Real, Real> Tet4::min_and_max_angle() const
329 {
330  Point n[4];
331 
332  // Compute the outward normal vectors on each face
333  n[0] = (this->point(2) - this->point(0)).cross(this->point(1) - this->point(0));
334  n[1] = (this->point(1) - this->point(0)).cross(this->point(3) - this->point(0));
335  n[2] = (this->point(2) - this->point(1)).cross(this->point(3) - this->point(1));
336  n[3] = (this->point(0) - this->point(2)).cross(this->point(3) - this->point(2));
337 
338  Real dihedral_angles[6]; // 01, 02, 03, 12, 13, 23
339 
340  // Compute dihedral angles
341  for (unsigned int k=0,i=0; i<4; ++i)
342  for (unsigned int j=i+1; j<4; ++j,k+=1)
343  dihedral_angles[k] = std::acos(n[i]*n[j] / n[i].norm() / n[j].norm()); // return value is between 0 and PI
344 
345  // Return max/min dihedral angles
346  return std::make_pair(*std::min_element(dihedral_angles, dihedral_angles+6),
347  *std::max_element(dihedral_angles, dihedral_angles+6));
348 
349 }
350 
351 
352 
354 {
355  return this->compute_key(this->node_id(0),
356  this->node_id(1),
357  this->node_id(2),
358  this->node_id(3));
359 }
360 
361 
362 
363 bool Tet4::contains_point (const Point & p, Real tol) const
364 {
365  // See the response by Tony Noe on this thread.
366  // http://bbs.dartmouth.edu/~fangq/MATH/download/source/Point_in_Tetrahedron.htm
367  Point
368  col1 = point(1) - point(0),
369  col2 = point(2) - point(0),
370  col3 = point(3) - point(0);
371 
372  Point r;
373 
374  libmesh_try
375  {
376  RealTensorValue(col1(0), col2(0), col3(0),
377  col1(1), col2(1), col3(1),
378  col1(2), col2(2), col3(2)).solve(p - point(0), r);
379  }
380  libmesh_catch (ConvergenceFailure &)
381  {
382  // The Jacobian was singular, therefore the Tet had
383  // zero volume. In this degenerate case, it is not
384  // possible for the Tet to contain any points.
385  return false;
386  }
387 
388  // The point p is inside the tetrahedron if r1 + r2 + r3 < 1 and
389  // 0 <= ri <= 1 for i = 1,2,3.
390  return
391  r(0) > -tol &&
392  r(1) > -tol &&
393  r(2) > -tol &&
394  r(0) + r(1) + r(2) < 1.0 + tol;
395 }
396 
397 
398 
399 #ifdef LIBMESH_ENABLE_AMR
400 float Tet4::embedding_matrix (const unsigned int i,
401  const unsigned int j,
402  const unsigned int k) const
403 {
404  // Choose an optimal diagonal, if one has not already been selected
405  this->choose_diagonal();
406 
407  // Permuted j and k indices
408  unsigned int
409  jp=j,
410  kp=k;
411 
412  if ((i>3) && (this->_diagonal_selection!=DIAG_02_13))
413  {
414  // Just the enum value cast to an unsigned int...
415  const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection);
416 
417  // Permute j, k:
418  // ds==1 ds==2
419  // 0 -> 1 0 -> 2
420  // 1 -> 2 1 -> 0
421  // 2 -> 0 2 -> 1
422  if (jp != 3)
423  jp = (jp+ds)%3;
424 
425  if (kp != 3)
426  kp = (kp+ds)%3;
427  }
428 
429  // Debugging
430  // libMesh::err << "Selected diagonal " << _diagonal_selection << std::endl;
431  // libMesh::err << "j=" << j << std::endl;
432  // libMesh::err << "k=" << k << std::endl;
433  // libMesh::err << "jp=" << jp << std::endl;
434  // libMesh::err << "kp=" << kp << std::endl;
435 
436  // Call embedding matrix with permuted indices
437  return this->_embedding_matrix[i][jp][kp];
438 }
439 
440 
441 
442 // void Tet4::reselect_diagonal (const Diagonal diag)
443 // {
444 // /* Make sure that the element has just been refined. */
445 // libmesh_assert(_children);
446 // libmesh_assert_equal_to (n_children(), 8);
447 // libmesh_assert_equal_to (_children[0]->refinement_flag(), JUST_REFINED);
448 // libmesh_assert_equal_to (_children[1]->refinement_flag(), JUST_REFINED);
449 // libmesh_assert_equal_to (_children[2]->refinement_flag(), JUST_REFINED);
450 // libmesh_assert_equal_to (_children[3]->refinement_flag(), JUST_REFINED);
451 // libmesh_assert_equal_to (_children[4]->refinement_flag(), JUST_REFINED);
452 // libmesh_assert_equal_to (_children[5]->refinement_flag(), JUST_REFINED);
453 // libmesh_assert_equal_to (_children[6]->refinement_flag(), JUST_REFINED);
454 // libmesh_assert_equal_to (_children[7]->refinement_flag(), JUST_REFINED);
455 //
456 // /* Check whether anything has to be changed. */
457 // if (_diagonal_selection!=diag)
458 // {
459 // /* Set new diagonal selection. */
460 // _diagonal_selection = diag;
461 //
462 // /* The first four children do not have to be changed. For the
463 // others, only the nodes have to be changed. Note that we have
464 // to keep track of the nodes ourselves since there is no \p
465 // MeshRefinement object with a valid \p _new_nodes_map
466 // available. */
467 // for (unsigned int c=4; c<this->n_children(); c++)
468 // {
469 // Elem * child = this->child_ptr(c);
470 // for (unsigned int nc=0; nc<child->n_nodes(); nc++)
471 // {
472 // /* Unassign the current node. */
473 // child->set_node(nc) = nullptr;
474 //
475 // /* We have to find the correct new node now. We know
476 // that it exists somewhere. We make use of the fact
477 // that the embedding matrix for these children consists
478 // of entries 0.0 and 0.5 only. Also, we make use of
479 // the properties of the embedding matrix for the first
480 // (unchanged) four children, which allow us to use a
481 // simple mechanism to find the required node. */
482 //
483 //
484 // unsigned int first_05_in_embedding_matrix = libMesh::invalid_uint;
485 // for (unsigned int n=0; n<this->n_nodes(); n++)
486 // {
487 // if (this->embedding_matrix(c,nc,n) != 0.0)
488 // {
489 // /* It must be 0.5 then. Check whether it's the
490 // first or second time that we get a 0.5
491 // value. */
492 // if (first_05_in_embedding_matrix==libMesh::invalid_uint)
493 // {
494 // /* First time, so just remember this position. */
495 // first_05_in_embedding_matrix = n;
496 // }
497 // else
498 // {
499 // /* Second time, so we know now which node to
500 // use. */
501 // child->set_node(nc) = this->child_ptr(n)->node_ptr(first_05_in_embedding_matrix);
502 // }
503 //
504 // }
505 // }
506 //
507 // /* Make sure that a node has been found. */
508 // libmesh_assert(child->node_ptr(nc));
509 // }
510 // }
511 // }
512 // }
513 
514 
515 
516 // void Tet4::reselect_optimal_diagonal (const Diagonal exclude_this)
517 // {
518 // Real diag_01_23 = (this->point(0)+this->point(1)-this->point(2)-this->point(3)).norm_sq();
519 // Real diag_02_13 = (this->point(0)-this->point(1)+this->point(2)-this->point(3)).norm_sq();
520 // Real diag_03_12 = (this->point(0)-this->point(1)-this->point(2)+this->point(3)).norm_sq();
521 //
522 // Diagonal use_this = INVALID_DIAG;
523 // switch (exclude_this)
524 // {
525 // case DIAG_01_23:
526 // use_this = DIAG_02_13;
527 // if (diag_03_12 < diag_02_13)
528 // {
529 // use_this = DIAG_03_12;
530 // }
531 // break;
532 //
533 // case DIAG_02_13:
534 // use_this = DIAG_03_12;
535 // if (diag_01_23 < diag_03_12)
536 // {
537 // use_this = DIAG_01_23;
538 // }
539 // break;
540 //
541 // case DIAG_03_12:
542 // use_this = DIAG_02_13;
543 // if (diag_01_23 < diag_02_13)
544 // {
545 // use_this = DIAG_01_23;
546 // }
547 // break;
548 //
549 // default:
550 // use_this = DIAG_02_13;
551 // if (diag_01_23 < diag_02_13 || diag_03_12 < diag_02_13)
552 // {
553 // if (diag_01_23 < diag_03_12)
554 // {
555 // use_this = DIAG_01_23;
556 // }
557 // else
558 // {
559 // use_this = DIAG_03_12;
560 // }
561 // }
562 // break;
563 // }
564 //
565 // reselect_diagonal (use_this);
566 // }
567 #endif // #ifdef LIBMESH_ENABLE_AMR
568 
569 } // namespace libMesh
virtual float embedding_matrix(const unsigned int i, const unsigned int j, const unsigned int k) const override
Definition: cell_tet4.C:400
Node ** _nodes
Definition: elem.h:1695
std::pair< Real, Real > min_and_max_angle() const
Definition: cell_tet4.C:328
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: cell_tet4.C:136
Real norm() const
Definition: type_vector.h:912
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: cell_tet4.h:178
virtual unsigned int n_sides() const override final
Definition: cell_tet.h:74
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Definition: cell_tet4.h:243
unsigned short int side
Definition: xdr_io.C:50
virtual bool is_edge(const unsigned int i) const override
Definition: cell_tet4.C:69
IterBase * end
static const int num_edges
Definition: cell_tet4.h:169
Diagonal _diagonal_selection
Definition: cell_tet.h:213
virtual unsigned int n_sub_elem() const override
Definition: cell_tet4.h:84
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: cell_tet4.C:185
virtual bool is_vertex(const unsigned int i) const override
Definition: cell_tet4.C:64
virtual dof_id_type key() const override
Definition: cell_tet4.C:353
TensorValue< Real > RealTensorValue
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1054
static const int nodes_per_side
Definition: cell_tet4.h:171
static const int num_nodes
Definition: cell_tet4.h:167
virtual Order default_order() const override
Definition: cell_tet4.C:142
static const int num_children
Definition: cell_tet4.h:170
virtual unsigned int n_edges() const override final
Definition: cell_tet.h:84
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: cell_tet4.C:79
virtual bool is_face(const unsigned int i) const override
Definition: cell_tet4.C:74
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: cell_tet4.C:126
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
Definition: cell_tet4.C:94
subdomain_id_type subdomain_id() const
Definition: elem.h:2034
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1957
virtual Real volume() const override
Definition: cell_tet4.C:310
virtual bool contains_point(const Point &p, Real tol) const override
Definition: cell_tet4.C:363
void choose_diagonal() const
Definition: cell_tet.C:163
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
Definition: cell_tet4.h:184
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2754
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Definition: cell_tet4.C:177
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy) override
Definition: cell_tet4.C:147
bool is_child_on_side_helper(const unsigned int c, const unsigned int s, const unsigned int checked_nodes[][3]) const
Definition: cell_tet.C:111
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
static const int nodes_per_edge
Definition: cell_tet4.h:172
static const int num_sides
Definition: cell_tet4.h:168
std::unique_ptr< Elem > side(const unsigned int i) const
Definition: elem.h:2202
uint8_t dof_id_type
Definition: id_types.h:64