libMesh
edge_edge3.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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 
20 // Local includes
21 #include "libmesh/edge_edge3.h"
22 #include "libmesh/enum_io_package.h"
23 #include "libmesh/enum_order.h"
24 
25 namespace libMesh
26 {
27 
28 // Edge3 class static member initializations
29 const int Edge3::num_nodes;
30 const int Edge3::num_sides;
31 const int Edge3::num_edges;
32 const int Edge3::num_children;
33 const int Edge3::nodes_per_side;
34 const int Edge3::nodes_per_edge;
35 
36 #ifdef LIBMESH_ENABLE_AMR
37 
38 const Real Edge3::_embedding_matrix[2][3][3] =
39  {
40  // embedding matrix for child 0
41  {
42  // 0 1 2
43  {1.0, 0.0, 0.0}, // left
44  {0.0, 0.0, 1.0}, // right
45  {0.375,-0.125,0.75} // middle
46  },
47 
48  // embedding matrix for child 1
49  {
50  // 0 1 2
51  {0.0, 0.0, 1.0}, // left
52  {0.0, 1.0, 0.0}, // right
53  {-0.125,0.375,0.75} // middle
54  }
55  };
56 
57 #endif
58 
59 bool Edge3::is_vertex(const unsigned int i) const
60 {
61  if (i < 2)
62  return true;
63  return false;
64 }
65 
66 bool Edge3::is_edge(const unsigned int i) const
67 {
68  if (i < 2)
69  return false;
70  return true;
71 }
72 
73 bool Edge3::is_face(const unsigned int ) const
74 {
75  return false;
76 }
77 
78 bool Edge3::is_node_on_side(const unsigned int n,
79  const unsigned int s) const
80 {
81  libmesh_assert_less (s, 2);
82  libmesh_assert_less (n, Edge3::num_nodes);
83  return (s == n);
84 }
85 
86 bool Edge3::is_node_on_edge(const unsigned int,
87  const unsigned int libmesh_dbg_var(e)) const
88 {
89  libmesh_assert_equal_to (e, 0);
90  return true;
91 }
92 
93 
94 
96 {
97  Point v = this->point(1) - this->point(0);
98  return (v.relative_fuzzy_equals
99  ((this->point(2) - this->point(0))*2, affine_tol));
100 }
101 
102 
103 
105 {
106  // At the moment this only makes sense for Lagrange elements
107  libmesh_assert_equal_to(this->mapping_type(), LAGRANGE_MAP);
108 
109  // The "Jacobian vector" (dx/dxi, dy/dxi, dz/dxi) is:
110  // j(xi) := a*xi + b, where
111  Point a = this->point(0) + this->point(1) - 2 * this->point(2);
112  Point b = Real(.5) * (this->point(1) - this->point(0));
113 
114  // Now we solve for the point xi_m where j(xi_m) \cdot j(0) = 0.
115  // If this occurs somewhere on the reference element, then the
116  // element is not invertible.
117  // j(xi_m) . j(0) = 0
118  // <=>
119  // (a*xi_m + b) . b = 0
120  // <=>
121  // (a.b)*xi_m + b.b = 0
122  // <=>
123  // xi_m = -(b.b) / (a.b)
124 
125  // 1.) If b.b==0, then the endpoints of the Edge3 are at the same
126  // location, and the map is therefore not invertible.
127  Real b_norm2 = b.norm_sq();
128  if (b_norm2 <= tol*tol)
129  return false;
130 
131  // 2.) If a.b==0, but b != 0 (see above), then the element is
132  // invertible but we don't want to divide by zero in the
133  // formula, so simply return true.
134  Real ab = a * b;
135  if (std::abs(ab) <= tol*tol)
136  return true;
137 
138  Real xi_m = -b_norm2 / ab;
139  return (xi_m < -1.) || (xi_m > 1.);
140 }
141 
142 
143 
145 {
146  return SECOND;
147 }
148 
149 
150 
151 void Edge3::connectivity(const unsigned int sc,
152  const IOPackage iop,
153  std::vector<dof_id_type> & conn) const
154 {
155  libmesh_assert_less_equal (sc, 1);
156  libmesh_assert_less (sc, this->n_sub_elem());
157  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
158 
159  // Create storage
160  conn.resize(2);
161 
162  switch (iop)
163  {
164  case TECPLOT:
165  {
166  switch (sc)
167  {
168  case 0:
169  conn[0] = this->node_id(0)+1;
170  conn[1] = this->node_id(2)+1;
171  return;
172 
173  case 1:
174  conn[0] = this->node_id(2)+1;
175  conn[1] = this->node_id(1)+1;
176  return;
177 
178  default:
179  libmesh_error_msg("Invalid sc = " << sc);
180  }
181  }
182 
183 
184  case VTK:
185  {
186  conn.resize(3);
187  conn[0] = this->node_id(0);
188  conn[1] = this->node_id(1);
189  conn[2] = this->node_id(2);
190  return;
191 
192  /*
193  switch (sc)
194  {
195  case 0:
196  conn[0] = this->node_id(0);
197  conn[1] = this->node_id(2);
198 
199  return;
200 
201  case 1:
202  conn[0] = this->node_id(2);
203  conn[1] = this->node_id(1);
204 
205  return;
206 
207  default:
208  libmesh_error_msg("Invalid sc = " << sc);
209  }
210  */
211  }
212 
213  default:
214  libmesh_error_msg("Unsupported IO package " << iop);
215  }
216 }
217 
218 
219 
220 std::pair<unsigned short int, unsigned short int>
221 Edge3::second_order_child_vertex (const unsigned int) const
222 {
223  return std::pair<unsigned short int, unsigned short int>(0,0);
224 }
225 
226 
227 
229 {
230  // Finding the (exact) length of a general quadratic element
231  // is a surprisingly complicated formula.
232  Point A = this->point(0) + this->point(1) - 2*this->point(2);
233  Point B = (this->point(1) - this->point(0))/2;
234 
235  const Real a = A.norm_sq();
236  const Real c = B.norm_sq();
237 
238  // Degenerate straight line case
239  if (a == 0.)
240  return 2. * std::sqrt(c);
241 
242  const Real b = -2.*std::abs(A*B);
243 
244  const Real sqrt_term1 = a - b + c;
245  const Real sqrt_term2 = a + b + c;
246 
247  // Fall back on straight line case instead of computing nan
248  // Note: b can be positive or negative so we have to check both cases.
249  if (sqrt_term1 < 0. || sqrt_term2 < 0.)
250  return 2. * std::sqrt(c);
251 
252  const Real r1 = std::sqrt(sqrt_term1);
253  const Real r2 = std::sqrt(sqrt_term2);
254  const Real rsum = r1 + r2;
255  const Real root_a = std::sqrt(a);
256  const Real b_over_root_a = b / root_a;
257 
258  // Pre-compute the denominator of the log term. If it's zero, fall
259  // back on the straight line case.
260  const Real log_term_denom = -root_a - 0.5*b_over_root_a + r2;
261  if (log_term_denom == 0. || rsum == 0.)
262  return 2. * std::sqrt(c);
263 
264  Real log1p_arg = 2*(root_a - b/rsum) / log_term_denom;
265 
266  return
267  0.5*(rsum + b_over_root_a*b_over_root_a/rsum +
268  (c-0.25*b_over_root_a*b_over_root_a)*std::log1p(log1p_arg)/root_a);
269 }
270 
271 
272 
274 {
275  // This might be a curved line through 2-space or 3-space, in which
276  // case the full bounding box can be larger than the bounding box of
277  // just the nodes.
278  Point pmin, pmax;
279 
280  for (unsigned d=0; d<LIBMESH_DIM; ++d)
281  {
282  Real center = this->point(2)(d);
283  Real hd = std::max(std::abs(center - this->point(0)(d)),
284  std::abs(center - this->point(1)(d)));
285 
286  pmin(d) = center - hd;
287  pmax(d) = center + hd;
288  }
289 
290  return BoundingBox(pmin, pmax);
291 }
292 
293 
294 
296 {
297  return this->compute_key(this->node_id(2));
298 }
299 
300 
301 void Edge3::flip(BoundaryInfo * boundary_info)
302 {
303  swap2nodes(0,1);
304  swap2neighbors(0,1);
305  swap2boundarysides(0,1,boundary_info);
306 }
307 
308 
309 } // namespace libMesh
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: edge_edge3.C:221
virtual bool is_face(const unsigned int i) const override
Definition: edge_edge3.C:73
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: edge_edge3.C:86
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
Definition: edge_edge3.C:301
virtual BoundingBox loose_bounding_box() const override
Definition: edge_edge3.C:273
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: edge_edge3.C:78
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
Definition: elem.C:3155
virtual Real volume() const override
An optimized method for computing the length of a 3-node edge.
Definition: edge_edge3.C:228
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:53
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
Definition: edge_edge3.h:203
virtual Order default_order() const override
Definition: edge_edge3.C:144
virtual unsigned int n_sub_elem() const override
Definition: edge_edge3.h:83
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:57
Definition: assembly.h:38
static const int nodes_per_edge
Definition: edge_edge3.h:205
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: edge_edge3.C:151
auto norm_sq() const -> decltype(std::norm(T()))
Definition: type_vector.h:948
ElemMappingType mapping_type() const
Definition: elem.h:2957
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:1933
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
static const int nodes_per_side
Definition: edge_edge3.h:204
virtual bool is_edge(const unsigned int i) const override
Definition: edge_edge3.C:66
static const int num_sides
Definition: edge_edge3.h:201
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
Definition: elem.h:1902
static const int num_edges
Definition: edge_edge3.h:202
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:1943
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool has_invertible_map(Real tol) const override
Definition: edge_edge3.C:104
static const int num_nodes
Geometric constants for Edge3.
Definition: edge_edge3.h:200
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:3131
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
Definition: edge_edge3.h:232
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2299
virtual dof_id_type key() const override
Compute a unique key for this element which is suitable for hashing (not necessarily unique...
Definition: edge_edge3.C:295
const Point & point(const unsigned int i) const
Definition: elem.h:2277
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:1004
virtual bool has_affine_map() const override
Definition: edge_edge3.C:95
virtual bool is_vertex(const unsigned int i) const override
Definition: edge_edge3.C:59
uint8_t dof_id_type
Definition: id_types.h:67