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