libMesh
fe_monomial_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 
29 
30 
31 template <>
33  const Order libmesh_dbg_var(order),
34  const unsigned int i,
35  const Point & p)
36 {
37 #if LIBMESH_DIM > 1
38 
39  libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
40  (static_cast<unsigned int>(order)+2)/2);
41 
42  const Real xi = p(0);
43  const Real eta = p(1);
44 
45  switch (i)
46  {
47  // constant
48  case 0:
49  return 1.;
50 
51  // linear
52  case 1:
53  return xi;
54 
55  case 2:
56  return eta;
57 
58  // quadratics
59  case 3:
60  return xi*xi;
61 
62  case 4:
63  return xi*eta;
64 
65  case 5:
66  return eta*eta;
67 
68  // cubics
69  case 6:
70  return xi*xi*xi;
71 
72  case 7:
73  return xi*xi*eta;
74 
75  case 8:
76  return xi*eta*eta;
77 
78  case 9:
79  return eta*eta*eta;
80 
81  // quartics
82  case 10:
83  return xi*xi*xi*xi;
84 
85  case 11:
86  return xi*xi*xi*eta;
87 
88  case 12:
89  return xi*xi*eta*eta;
90 
91  case 13:
92  return xi*eta*eta*eta;
93 
94  case 14:
95  return eta*eta*eta*eta;
96 
97  default:
98  unsigned int o = 0;
99  for (; i >= (o+1)*(o+2)/2; o++) { }
100  unsigned int ny = i - (o*(o+1)/2);
101  unsigned int nx = o - ny;
102  Real val = 1.;
103  for (unsigned int index=0; index != nx; index++)
104  val *= xi;
105  for (unsigned int index=0; index != ny; index++)
106  val *= eta;
107  return val;
108  }
109 
110 #else // LIBMESH_DIM == 1
111  libmesh_ignore(i, p);
112  libmesh_assert(order);
113  libmesh_not_implemented();
114 #endif
115 }
116 
117 
118 
119 template <>
121  const Order order,
122  const unsigned int i,
123  const Point & p,
124  const bool add_p_level)
125 {
126  libmesh_assert(elem);
127 
128  // by default call the orientation-independent shape functions
129  return FE<2,MONOMIAL>::shape(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, p);
130 }
131 
132 
133 
134 template <>
136  const Order libmesh_dbg_var(order),
137  const unsigned int i,
138  const unsigned int j,
139  const Point & p)
140 {
141 #if LIBMESH_DIM > 1
142 
143 
144  libmesh_assert_less (j, 2);
145 
146  libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
147  (static_cast<unsigned int>(order)+2)/2);
148 
149  const Real xi = p(0);
150  const Real eta = p(1);
151 
152  // monomials. since they are hierarchic we only need one case block.
153 
154  switch (j)
155  {
156  // d()/dxi
157  case 0:
158  {
159  switch (i)
160  {
161  // constants
162  case 0:
163  return 0.;
164 
165  // linears
166  case 1:
167  return 1.;
168 
169  case 2:
170  return 0.;
171 
172  // quadratics
173  case 3:
174  return 2.*xi;
175 
176  case 4:
177  return eta;
178 
179  case 5:
180  return 0.;
181 
182  // cubics
183  case 6:
184  return 3.*xi*xi;
185 
186  case 7:
187  return 2.*xi*eta;
188 
189  case 8:
190  return eta*eta;
191 
192  case 9:
193  return 0.;
194 
195  // quartics
196  case 10:
197  return 4.*xi*xi*xi;
198 
199  case 11:
200  return 3.*xi*xi*eta;
201 
202  case 12:
203  return 2.*xi*eta*eta;
204 
205  case 13:
206  return eta*eta*eta;
207 
208  case 14:
209  return 0.;
210 
211  default:
212  unsigned int o = 0;
213  for (; i >= (o+1)*(o+2)/2; o++) { }
214  unsigned int ny = i - (o*(o+1)/2);
215  unsigned int nx = o - ny;
216  Real val = nx;
217  for (unsigned int index=1; index < nx; index++)
218  val *= xi;
219  for (unsigned int index=0; index != ny; index++)
220  val *= eta;
221  return val;
222  }
223  }
224 
225 
226  // d()/deta
227  case 1:
228  {
229  switch (i)
230  {
231  // constants
232  case 0:
233  return 0.;
234 
235  // linears
236  case 1:
237  return 0.;
238 
239  case 2:
240  return 1.;
241 
242  // quadratics
243  case 3:
244  return 0.;
245 
246  case 4:
247  return xi;
248 
249  case 5:
250  return 2.*eta;
251 
252  // cubics
253  case 6:
254  return 0.;
255 
256  case 7:
257  return xi*xi;
258 
259  case 8:
260  return 2.*xi*eta;
261 
262  case 9:
263  return 3.*eta*eta;
264 
265  // quartics
266  case 10:
267  return 0.;
268 
269  case 11:
270  return xi*xi*xi;
271 
272  case 12:
273  return 2.*xi*xi*eta;
274 
275  case 13:
276  return 3.*xi*eta*eta;
277 
278  case 14:
279  return 4.*eta*eta*eta;
280 
281  default:
282  unsigned int o = 0;
283  for (; i >= (o+1)*(o+2)/2; o++) { }
284  unsigned int ny = i - (o*(o+1)/2);
285  unsigned int nx = o - ny;
286  Real val = ny;
287  for (unsigned int index=0; index != nx; index++)
288  val *= xi;
289  for (unsigned int index=1; index < ny; index++)
290  val *= eta;
291  return val;
292  }
293  }
294 
295  default:
296  libmesh_error_msg("Invalid shape function derivative j = " << j);
297  }
298 
299 #else // LIBMESH_DIM == 1
300  libmesh_ignore(i, j, p);
301  libmesh_assert(order);
302  libmesh_not_implemented();
303 #endif
304 }
305 
306 
307 
308 template <>
310  const Order order,
311  const unsigned int i,
312  const unsigned int j,
313  const Point & p,
314  const bool add_p_level)
315 {
316  libmesh_assert(elem);
317 
318  // by default call the orientation-independent shape functions
319  return FE<2,MONOMIAL>::shape_deriv(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
320 }
321 
322 
323 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
324 
325 template <>
327  const Order libmesh_dbg_var(order),
328  const unsigned int i,
329  const unsigned int j,
330  const Point & p)
331 {
332 #if LIBMESH_DIM > 1
333 
334 
335  libmesh_assert_less_equal (j, 2);
336 
337  libmesh_assert_less (i, (static_cast<unsigned int>(order)+1)*
338  (static_cast<unsigned int>(order)+2)/2);
339 
340  const Real xi = p(0);
341  const Real eta = p(1);
342 
343  // monomials. since they are hierarchic we only need one case block.
344 
345  switch (j)
346  {
347  // d^2()/dxi^2
348  case 0:
349  {
350  switch (i)
351  {
352  // constants
353  case 0:
354  // linears
355  case 1:
356  case 2:
357  return 0.;
358 
359  // quadratics
360  case 3:
361  return 2.;
362 
363  case 4:
364  case 5:
365  return 0.;
366 
367  // cubics
368  case 6:
369  return 6.*xi;
370 
371  case 7:
372  return 2.*eta;
373 
374  case 8:
375  case 9:
376  return 0.;
377 
378  // quartics
379  case 10:
380  return 12.*xi*xi;
381 
382  case 11:
383  return 6.*xi*eta;
384 
385  case 12:
386  return 2.*eta*eta;
387 
388  case 13:
389  case 14:
390  return 0.;
391 
392  default:
393  unsigned int o = 0;
394  for (; i >= (o+1)*(o+2)/2; o++) { }
395  unsigned int ny = i - (o*(o+1)/2);
396  unsigned int nx = o - ny;
397  Real val = nx * (nx - 1);
398  for (unsigned int index=2; index < nx; index++)
399  val *= xi;
400  for (unsigned int index=0; index != ny; index++)
401  val *= eta;
402  return val;
403  }
404  }
405 
406  // d^2()/dxideta
407  case 1:
408  {
409  switch (i)
410  {
411  // constants
412  case 0:
413 
414  // linears
415  case 1:
416  case 2:
417  return 0.;
418 
419  // quadratics
420  case 3:
421  return 0.;
422 
423  case 4:
424  return 1.;
425 
426  case 5:
427  return 0.;
428 
429  // cubics
430  case 6:
431  return 0.;
432  case 7:
433  return 2.*xi;
434 
435  case 8:
436  return 2.*eta;
437 
438  case 9:
439  return 0.;
440 
441  // quartics
442  case 10:
443  return 0.;
444 
445  case 11:
446  return 3.*xi*xi;
447 
448  case 12:
449  return 4.*xi*eta;
450 
451  case 13:
452  return 3.*eta*eta;
453 
454  case 14:
455  return 0.;
456 
457  default:
458  unsigned int o = 0;
459  for (; i >= (o+1)*(o+2)/2; o++) { }
460  unsigned int ny = i - (o*(o+1)/2);
461  unsigned int nx = o - ny;
462  Real val = nx * ny;
463  for (unsigned int index=1; index < nx; index++)
464  val *= xi;
465  for (unsigned int index=1; index < ny; index++)
466  val *= eta;
467  return val;
468  }
469  }
470 
471  // d^2()/deta^2
472  case 2:
473  {
474  switch (i)
475  {
476  // constants
477  case 0:
478 
479  // linears
480  case 1:
481  case 2:
482  return 0.;
483 
484  // quadratics
485  case 3:
486  case 4:
487  return 0.;
488 
489  case 5:
490  return 2.;
491 
492  // cubics
493  case 6:
494  return 0.;
495 
496  case 7:
497  return 0.;
498 
499  case 8:
500  return 2.*xi;
501 
502  case 9:
503  return 6.*eta;
504 
505  // quartics
506  case 10:
507  case 11:
508  return 0.;
509 
510  case 12:
511  return 2.*xi*xi;
512 
513  case 13:
514  return 6.*xi*eta;
515 
516  case 14:
517  return 12.*eta*eta;
518 
519  default:
520  unsigned int o = 0;
521  for (; i >= (o+1)*(o+2)/2; o++) { }
522  unsigned int ny = i - (o*(o+1)/2);
523  unsigned int nx = o - ny;
524  Real val = ny * (ny - 1);
525  for (unsigned int index=0; index != nx; index++)
526  val *= xi;
527  for (unsigned int index=2; index < ny; index++)
528  val *= eta;
529  return val;
530  }
531  }
532 
533  default:
534  libmesh_error_msg("Invalid shape function derivative j = " << j);
535  }
536 
537 #else // LIBMESH_DIM == 1
538  libmesh_assert(order);
539  libmesh_ignore(i, j, p);
540  libmesh_not_implemented();
541 #endif
542 }
543 
544 
545 
546 template <>
548  const Order order,
549  const unsigned int i,
550  const unsigned int j,
551  const Point & p,
552  const bool add_p_level)
553 {
554  libmesh_assert(elem);
555 
556  // by default call the orientation-independent shape functions
557  return FE<2,MONOMIAL>::shape_second_deriv(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
558 }
559 
560 #endif
561 
562 } // 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
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::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
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::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
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::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33