libMesh
fe_nedelec_one_shape_2D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 // C++ includes
20 
21 // Local includes
22 #include "libmesh/fe.h"
23 #include "libmesh/elem.h"
24 
25 namespace libMesh
26 {
27 
28 template <>
30  const Order,
31  const unsigned int,
32  const Point &)
33 {
34  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
35  return RealGradient();
36 }
37 
38 
39 // An excellent discussion of Nedelec shape functions is given in
40 // http://www.dealii.org/developer/reports/nedelec/nedelec.pdf
41 template <>
43  const Order order,
44  const unsigned int i,
45  const Point & p,
46  const bool add_p_level)
47 {
48 #if LIBMESH_DIM > 1
49  libmesh_assert(elem);
50 
51  const Order total_order = static_cast<Order>(order + add_p_level * elem->p_level());
52 
53  switch (total_order)
54  {
55  case FIRST:
56  {
57  switch (elem->type())
58  {
59  case QUAD8:
60  case QUAD9:
61  {
62  libmesh_assert_less (i, 4);
63 
64  const Real xi = p(0);
65  const Real eta = p(1);
66 
67  // Even with a loose inverse_map tolerance we ought to
68  // be nearly on the element interior in master
69  // coordinates
70  libmesh_assert_less_equal ( std::fabs(xi), 1.0+10*TOLERANCE );
71  libmesh_assert_less_equal ( std::fabs(eta), 1.0+10*TOLERANCE );
72 
73  switch(i)
74  {
75  case 0:
76  {
77  if (elem->point(0) > elem->point(1))
78  return RealGradient( -0.25*(1.0-eta), 0.0 );
79  else
80  return RealGradient( 0.25*(1.0-eta), 0.0 );
81  }
82  case 1:
83  {
84  if (elem->point(1) > elem->point(2))
85  return RealGradient( 0.0, -0.25*(1.0+xi) );
86  else
87  return RealGradient( 0.0, 0.25*(1.0+xi) );
88  }
89 
90  case 2:
91  {
92  if (elem->point(2) > elem->point(3))
93  return RealGradient( 0.25*(1.0+eta), 0.0 );
94  else
95  return RealGradient( -0.25*(1.0+eta), 0.0 );
96  }
97  case 3:
98  {
99  if (elem->point(3) > elem->point(0))
100  return RealGradient( 0.0, -0.25*(xi-1.0) );
101  else
102  return RealGradient( 0.0, 0.25*(xi-1.0) );
103  }
104 
105  default:
106  libmesh_error_msg("Invalid i = " << i);
107  }
108 
109  return RealGradient();
110  }
111 
112  case TRI6:
113  {
114  const Real xi = p(0);
115  const Real eta = p(1);
116 
117  libmesh_assert_less (i, 3);
118 
119  switch(i)
120  {
121  case 0:
122  {
123  if (elem->point(0) > elem->point(1))
124  return RealGradient( -1.0+eta, -xi );
125  else
126  return RealGradient( 1.0-eta, xi );
127  }
128  case 1:
129  {
130  if (elem->point(1) > elem->point(2))
131  return RealGradient( eta, -xi );
132  else
133  return RealGradient( -eta, xi );
134  }
135 
136  case 2:
137  {
138  if (elem->point(2) > elem->point(0))
139  return RealGradient( eta, -xi+1.0 );
140  else
141  return RealGradient( -eta, xi-1.0 );
142  }
143 
144  default:
145  libmesh_error_msg("Invalid i = " << i);
146  }
147  }
148 
149  default:
150  libmesh_error_msg("ERROR: Unsupported 2D element type!: " << elem->type());
151  }
152  }
153 
154  // unsupported order
155  default:
156  libmesh_error_msg("ERROR: Unsupported 2D FE order!: " << total_order);
157  }
158 #else // LIBMESH_DIM > 1
159  libmesh_ignore(elem, order, i, p, add_p_level);
160  libmesh_not_implemented();
161 #endif
162 }
163 
164 
165 
166 template <>
168  const Order,
169  const unsigned int,
170  const unsigned int,
171  const Point &)
172 {
173  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
174  return RealGradient();
175 }
176 
177 
178 
179 template <>
181  const Order order,
182  const unsigned int i,
183  const unsigned int j,
184  const Point &,
185  const bool add_p_level)
186 {
187 #if LIBMESH_DIM > 1
188  libmesh_assert(elem);
189  libmesh_assert_less (j, 2);
190 
191  const Order total_order = static_cast<Order>(order + add_p_level * elem->p_level());
192 
193  switch (total_order)
194  {
195  // linear Lagrange shape functions
196  case FIRST:
197  {
198  switch (elem->type())
199  {
200  case QUAD8:
201  case QUAD9:
202  {
203  libmesh_assert_less (i, 4);
204 
205  switch (j)
206  {
207  // d()/dxi
208  case 0:
209  {
210  switch(i)
211  {
212  case 0:
213  case 2:
214  return RealGradient();
215  case 1:
216  {
217  if (elem->point(1) > elem->point(2))
218  return RealGradient( 0.0, -0.25 );
219  else
220  return RealGradient( 0.0, 0.25 );
221  }
222  case 3:
223  {
224  if (elem->point(3) > elem->point(0))
225  return RealGradient( 0.0, -0.25 );
226  else
227  return RealGradient( 0.0, 0.25 );
228  }
229  default:
230  libmesh_error_msg("Invalid i = " << i);
231  }
232  } // j=0
233 
234  // d()/deta
235  case 1:
236  {
237  switch(i)
238  {
239  case 1:
240  case 3:
241  return RealGradient();
242  case 0:
243  {
244  if (elem->point(0) > elem->point(1))
245  return RealGradient( 0.25 );
246  else
247  return RealGradient( -0.25 );
248  }
249  case 2:
250  {
251  if (elem->point(2) > elem->point(3))
252  return RealGradient( 0.25 );
253  else
254  return RealGradient( -0.25 );
255  }
256  default:
257  libmesh_error_msg("Invalid i = " << i);
258  }
259  } // j=1
260 
261  default:
262  libmesh_error_msg("Invalid j = " << j);
263  }
264 
265  return RealGradient();
266  }
267 
268  case TRI6:
269  {
270  libmesh_assert_less (i, 3);
271 
272  // Account for edge flipping
273  Real f = 1.0;
274 
275  switch(i)
276  {
277  case 0:
278  {
279  if (elem->point(0) > elem->point(1))
280  f = -1.0;
281  break;
282  }
283  case 1:
284  {
285  if (elem->point(1) > elem->point(2))
286  f = -1.0;
287  break;
288  }
289  case 2:
290  {
291  if (elem->point(2) > elem->point(0))
292  f = -1.0;
293  break;
294  }
295  default:
296  libmesh_error_msg("Invalid i = " << i);
297  }
298 
299  switch (j)
300  {
301  // d()/dxi
302  case 0:
303  {
304  return RealGradient( 0.0, f*1.0);
305  }
306  // d()/deta
307  case 1:
308  {
309  return RealGradient( f*(-1.0) );
310  }
311  default:
312  libmesh_error_msg("Invalid j = " << j);
313  }
314  }
315 
316  default:
317  libmesh_error_msg("ERROR: Unsupported 2D element type!: " << elem->type());
318  }
319  }
320  // unsupported order
321  default:
322  libmesh_error_msg("ERROR: Unsupported 2D FE order!: " << total_order);
323  }
324 #else // LIBMESH_DIM > 1
325  libmesh_ignore(elem, order, i, j, add_p_level);
326  libmesh_not_implemented();
327 #endif
328 }
329 
330 
331 
332 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
333 
334 template <>
336  const Order,
337  const unsigned int,
338  const unsigned int,
339  const Point &)
340 {
341  libmesh_error_msg("Nedelec elements require the element type \nbecause edge orientation is needed.");
342  return RealGradient();
343 }
344 
345 
346 
347 template <>
349  const Order order,
350  const unsigned int libmesh_dbg_var(i),
351  const unsigned int libmesh_dbg_var(j),
352  const Point &,
353  const bool add_p_level)
354 {
355 #if LIBMESH_DIM > 1
356  libmesh_assert(elem);
357 
358  // j = 0 ==> d^2 phi / dxi^2
359  // j = 1 ==> d^2 phi / dxi deta
360  // j = 2 ==> d^2 phi / deta^2
361  libmesh_assert_less (j, 3);
362 
363  const Order total_order = static_cast<Order>(order + add_p_level * elem->p_level());
364 
365  switch (total_order)
366  {
367  // linear Lagrange shape functions
368  case FIRST:
369  {
370  switch (elem->type())
371  {
372  case QUAD8:
373  case QUAD9:
374  {
375  libmesh_assert_less (i, 4);
376  // All second derivatives for linear quads are zero.
377  return RealGradient();
378  }
379 
380  case TRI6:
381  {
382  libmesh_assert_less (i, 3);
383  // All second derivatives for linear triangles are zero.
384  return RealGradient();
385  }
386 
387  default:
388  libmesh_error_msg("ERROR: Unsupported 2D element type!: " << elem->type());
389 
390  } // end switch (type)
391  } // end case FIRST
392 
393  // unsupported order
394  default:
395  libmesh_error_msg("ERROR: Unsupported 2D FE order!: " << total_order);
396 
397  } // end switch (order)
398 
399 #else // LIBMESH_DIM > 1
400  libmesh_assert(true || i || j);
401  libmesh_ignore(elem, order, add_p_level);
402  libmesh_not_implemented();
403 #endif
404 }
405 
406 #endif
407 
408 } // namespace libMesh
libMesh::FE::shape_second_deriv
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::RealGradient
RealVectorValue RealGradient
Definition: hp_coarsentest.h:49
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::VectorValue< Real >
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::FE::shape_deriv
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::FIRST
Definition: enum_order.h:42
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33