edge_edge4.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 
20 // Local includes
21 #include "libmesh/edge_edge4.h"
23 #include "libmesh/enum_order.h"
24 
25 namespace libMesh
26 {
27 
28 // Edge4 class static member initializations
29 const int Edge4::num_nodes;
30 const int Edge4::num_children;
31 
32 #ifdef LIBMESH_ENABLE_AMR
33 
35  {
36  // embedding matrix for child 0
37 
38  {
39  // 0 1 2 3 // Shape function index
40  {1.0, 0.0, 0.0, 0.0}, // left, xi = -1
41  {-0.0625, -0.0625, 0.5625, 0.5625}, // right, xi = 0
42  {0.3125, 0.0625, 0.9375, -0.3125}, // middle left, xi = -2/3
43  {0.0, 0.0, 1.0, 0.0} // middle right, xi = -1/3
44  },
45 
46  // embedding matrix for child 1
47  {
48  // 0 1 2 3 // Shape function index
49  {-0.0625, -0.0625, 0.5625, 0.5625}, // left, xi = 0
50  {0.0, 1.0, 0.0, 0.0}, // right, xi = 1
51  {0.0, 0.0, 0.0, 1.0}, // middle left, xi = 1/3
52  {0.0625, 0.3125, -0.3125, 0.9375} // middle right, xi = 2/3
53  }
54  };
55 
56 #endif
57 
58 bool Edge4::is_vertex(const unsigned int i) const
59 {
60  return (i==0) || (i==1);
61 }
62 
63 bool Edge4::is_edge(const unsigned int i) const
64 {
65  return (i==2) || (i==3);
66 }
67 
68 bool Edge4::is_face(const unsigned int ) const
69 {
70  return false;
71 }
72 
73 bool Edge4::is_node_on_side(const unsigned int n,
74  const unsigned int s) const
75 {
76  libmesh_assert_less (s, 2);
77  libmesh_assert_less (n, Edge4::num_nodes);
78  return (s == n);
79 }
80 
81 bool Edge4::is_node_on_edge(const unsigned int,
82  const unsigned int libmesh_dbg_var(e)) const
83 {
84  libmesh_assert_equal_to (e, 0);
85  return true;
86 }
87 
88 
89 
91 {
92  if (!this->point(2).relative_fuzzy_equals
93  ((this->point(0)*2. + this->point(1))/3.))
94  return false;
95  if (!this->point(3).relative_fuzzy_equals
96  ((this->point(0) + this->point(1)*2.)/3.))
97  return false;
98  return true;
99 }
100 
101 
102 
104 {
105  return THIRD;
106 }
107 
108 
109 
110 void Edge4::connectivity(const unsigned int sc,
111  const IOPackage iop,
112  std::vector<dof_id_type> & conn) const
113 {
114  libmesh_assert_less_equal (sc, 2);
115  libmesh_assert_less (sc, this->n_sub_elem());
116  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
117 
118  // Create storage
119  conn.resize(2);
120 
121  switch(iop)
122  {
123  case TECPLOT:
124  {
125  switch (sc)
126  {
127  case 0:
128  conn[0] = this->node_id(0)+1;
129  conn[1] = this->node_id(2)+1;
130  return;
131 
132  case 1:
133  conn[0] = this->node_id(2)+1;
134  conn[1] = this->node_id(3)+1;
135  return;
136 
137  case 2:
138  conn[0] = this->node_id(3)+1;
139  conn[1] = this->node_id(1)+1;
140  return;
141 
142  default:
143  libmesh_error_msg("Invalid sc = " << sc);
144  }
145 
146  }
147 
148  case VTK:
149  {
150 
151  switch (sc)
152  {
153  case 0:
154  conn[0] = this->node_id(0);
155  conn[1] = this->node_id(2);
156  return;
157 
158  case 1:
159  conn[0] = this->node_id(2);
160  conn[1] = this->node_id(3);
161  return;
162 
163  case 2:
164  conn[0] = this->node_id(3);
165  conn[1] = this->node_id(1);
166  return;
167 
168  default:
169  libmesh_error_msg("Invalid sc = " << sc);
170  }
171  }
172 
173  default:
174  libmesh_error_msg("Unsupported IO package " << iop);
175  }
176 }
177 
178 
179 
181 {
182  // This might be a curved line through 2-space or 3-space, in which
183  // case the full bounding box can be larger than the bounding box of
184  // just the nodes.
185  //
186  // FIXME - I haven't yet proven the formula below to be correct for
187  // cubics - RHS
188  Point pmin, pmax;
189 
190  for (unsigned d=0; d<LIBMESH_DIM; ++d)
191  {
192  Real center = (this->point(2)(d) + this->point(3)(d))/2;
193  Real hd = std::max(std::abs(center - this->point(0)(d)),
194  std::abs(center - this->point(1)(d)));
195 
196  pmin(d) = center - hd;
197  pmax(d) = center + hd;
198  }
199 
200  return BoundingBox(pmin, pmax);
201 }
202 
203 
204 
206 {
207  return this->compute_key(this->node_id(0),
208  this->node_id(1),
209  this->node_id(2),
210  this->node_id(3));
211 }
212 
213 
214 
216 {
217  // Make copies of our points. It makes the subsequent calculations a bit
218  // shorter and avoids dereferencing the same pointer multiple times.
219  Point
220  x0 = point(0),
221  x1 = point(1),
222  x2 = point(2),
223  x3 = point(3);
224 
225  // Construct constant data vectors.
226  // \vec{x}_{\xi} = \vec{a1}*xi**2 + \vec{b1}*xi + \vec{c1}
227  // This is copy-pasted directly from the output of a Python script.
228  Point
229  a1 = -27*x0/16 + 27*x1/16 + 81*x2/16 - 81*x3/16,
230  b1 = 9*x0/8 + 9*x1/8 - 9*x2/8 - 9*x3/8,
231  c1 = x0/16 - x1/16 - 27*x2/16 + 27*x3/16;
232 
233  // 4 point quadrature, 7th-order accurate
234  const unsigned int N = 4;
235  const Real q[N] = {-std::sqrt(525 + 70*std::sqrt(30.)) / 35,
236  -std::sqrt(525 - 70*std::sqrt(30.)) / 35,
237  std::sqrt(525 - 70*std::sqrt(30.)) / 35,
238  std::sqrt(525 + 70*std::sqrt(30.)) / 35};
239  const Real w[N] = {(18 - std::sqrt(30.)) / 36,
240  (18 + std::sqrt(30.)) / 36,
241  (18 + std::sqrt(30.)) / 36,
242  (18 - std::sqrt(30.)) / 36};
243 
244  Real vol=0.;
245  for (unsigned int i=0; i<N; ++i)
246  vol += w[i] * (q[i]*q[i]*a1 + q[i]*b1 + c1).norm();
247 
248  return vol;
249 }
250 
251 } // namespace libMesh
static const int num_children
Definition: edge_edge4.h:182
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: edge_edge4.C:73
double abs(double a)
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: edge_edge4.C:110
virtual BoundingBox loose_bounding_box() const override
Definition: edge_edge4.C:180
virtual bool is_edge(const unsigned int i) const override
Definition: edge_edge4.C:63
static const int num_nodes
Definition: edge_edge4.h:181
virtual dof_id_type key() const override
Definition: edge_edge4.C:205
virtual bool is_face(const unsigned int i) const override
Definition: edge_edge4.C:68
long double max(long double a, double b)
virtual unsigned int n_sub_elem() const override
Definition: edge_edge4.h:81
virtual Real volume() const override
Definition: edge_edge4.C:215
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool is_vertex(const unsigned int i) const override
Definition: edge_edge4.C:58
virtual bool has_affine_map() const override
Definition: edge_edge4.C:90
virtual Order default_order() const override
Definition: edge_edge4.C:103
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2754
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Definition: edge_edge4.h:207
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: edge_edge4.C:81
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
uint8_t dof_id_type
Definition: id_types.h:64