face_tri.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/face_tri.h"
20 #include "libmesh/edge_edge2.h"
21 #include "libmesh/face_tri3.h"
23 
24 namespace libMesh
25 {
26 
27 
28 // ------------------------------------------------------------
29 // Tri class static member initializations
30 
31 
32 // We need to require C++11...
33 const Real Tri::_master_points[6][3] =
34  {
35  {0, 0},
36  {1, 0},
37  {0, 1},
38  {0.5, 0},
39  {0.5, 0.5},
40  {0, 0.5}
41  };
42 
43 
44 
45 
46 
47 
48 
49 // ------------------------------------------------------------
50 // Tri class member functions
51 dof_id_type Tri::key (const unsigned int s) const
52 {
53  libmesh_assert_less (s, this->n_sides());
54 
55  return this->compute_key(this->node_id(Tri3::side_nodes_map[s][0]),
56  this->node_id(Tri3::side_nodes_map[s][1]));
57 }
58 
59 
60 
61 unsigned int Tri::which_node_am_i(unsigned int side,
62  unsigned int side_node) const
63 {
64  libmesh_assert_less (side, this->n_sides());
65  libmesh_assert_less (side_node, 2);
66 
67  return Tri3::side_nodes_map[side][side_node];
68 }
69 
70 
71 
73 {
74  return this->compute_key(this->node_id(0),
75  this->node_id(1),
76  this->node_id(2));
77 }
78 
79 
80 
81 std::unique_ptr<Elem> Tri::side_ptr (const unsigned int i)
82 {
83  libmesh_assert_less (i, this->n_sides());
84 
85  std::unique_ptr<Elem> edge = libmesh_make_unique<Edge2>();
86 
87  for (unsigned n=0; n<edge->n_nodes(); ++n)
88  edge->set_node(n) = this->node_ptr(Tri3::side_nodes_map[i][n]);
89 
90  return edge;
91 }
92 
93 
94 
95 void Tri::side_ptr (std::unique_ptr<Elem> & side,
96  const unsigned int i)
97 {
98  this->simple_build_side_ptr<Tri3>(side, i, EDGE2);
99 }
100 
101 
102 
103 bool Tri::is_child_on_side(const unsigned int c,
104  const unsigned int s) const
105 {
106  libmesh_assert_less (c, this->n_children());
107  libmesh_assert_less (s, this->n_sides());
108 
109  return (c == s || c == (s+1)%3);
110 }
111 
112 
113 
115 {
116  switch (q)
117  {
118 
122  case DISTORTION:
123  case STRETCH:
124  {
125  const Point & p1 = this->point(0);
126  const Point & p2 = this->point(1);
127  const Point & p3 = this->point(2);
128 
129  Point v1 = p2 - p1;
130  Point v2 = p3 - p1;
131  Point v3 = p3 - p2;
132  const Real l1 = v1.norm();
133  const Real l2 = v2.norm();
134  const Real l3 = v3.norm();
135 
136  // if one length is 0, quality is quite bad!
137  if ((l1 <=0.) || (l2 <= 0.) || (l3 <= 0.))
138  return 0.;
139 
140  const Real s1 = std::sin(std::acos(v1*v2/l1/l2)/2.);
141  v1 *= -1;
142  const Real s2 = std::sin(std::acos(v1*v3/l1/l3)/2.);
143  const Real s3 = std::sin(std::acos(v2*v3/l2/l3)/2.);
144 
145  return 8. * s1 * s2 * s3;
146 
147  }
148  default:
149  return Elem::quality(q);
150  }
151 
157  return Elem::quality(q);
158 
159 }
160 
161 
162 
163 
164 
165 
166 std::pair<Real, Real> Tri::qual_bounds (const ElemQuality q) const
167 {
168  std::pair<Real, Real> bounds;
169 
170  switch (q)
171  {
172 
173  case MAX_ANGLE:
174  bounds.first = 60.;
175  bounds.second = 90.;
176  break;
177 
178  case MIN_ANGLE:
179  bounds.first = 30.;
180  bounds.second = 60.;
181  break;
182 
183  case CONDITION:
184  bounds.first = 1.;
185  bounds.second = 1.3;
186  break;
187 
188  case JACOBIAN:
189  bounds.first = 0.5;
190  bounds.second = 1.155;
191  break;
192 
193  case SIZE:
194  case SHAPE:
195  bounds.first = 0.25;
196  bounds.second = 1.;
197  break;
198 
199  case DISTORTION:
200  bounds.first = 0.6;
201  bounds.second = 1.;
202  break;
203 
204  default:
205  libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
206  bounds.first = -1;
207  bounds.second = -1;
208  }
209 
210  return bounds;
211 }
212 
213 } // namespace libMesh
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: face_tri.C:81
static const Real _master_points[6][3]
Definition: face_tri.h:176
Real norm() const
Definition: type_vector.h:912
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: face_tri3.h:157
unsigned short int side
Definition: xdr_io.C:50
virtual unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const override
Definition: face_tri.C:61
virtual unsigned int n_sides() const override final
Definition: face_tri.h:92
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override final
Definition: face_tri.C:103
virtual dof_id_type key() const override final
Definition: face_tri.C:72
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1957
virtual Real quality(const ElemQuality q) const
Definition: elem.C:1403
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
Definition: face_tri.C:166
virtual unsigned int n_children() const override final
Definition: face_tri.h:107
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2754
virtual Real quality(const ElemQuality q) const override
Definition: face_tri.C:114
OStreamProxy out(std::cout)
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
std::unique_ptr< Elem > side(const unsigned int i) const
Definition: elem.h:2202
uint8_t dof_id_type
Definition: id_types.h:64