Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_xyz_shape_2D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 
24 #include "libmesh/elem.h"
25 #include "libmesh/int_range.h"
26 
27 
28 namespace libMesh
29 {
30 
31 
33 
34 
35 template <>
36 Real FE<2,XYZ>::shape(const Elem * elem,
37  const Order libmesh_dbg_var(order),
38  const unsigned int i,
39  const Point & point_in,
40  const bool libmesh_dbg_var(add_p_level))
41 {
42 #if LIBMESH_DIM > 1
43 
44  libmesh_assert(elem);
45 
46  Point avg = elem->vertex_average();
47  Point max_distance = Point(0.,0.,0.);
48  for (const Point & p : elem->node_ref_range())
49  for (unsigned int d = 0; d < 2; d++)
50  {
51  const Real distance = std::abs(avg(d) - p(d));
52  max_distance(d) = std::max(distance, max_distance(d));
53  }
54 
55  const Real x = point_in(0);
56  const Real y = point_in(1);
57  const Real xc = avg(0);
58  const Real yc = avg(1);
59  const Real distx = max_distance(0);
60  const Real disty = max_distance(1);
61  const Real dx = (x - xc)/distx;
62  const Real dy = (y - yc)/disty;
63 
64 #ifndef NDEBUG
65  // totalorder is only used in the assertion below, so
66  // we avoid declaring it when asserts are not active.
67  const unsigned int totalorder = order + add_p_level * elem->p_level();
68 #endif
69  libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2);
70 
71 
72  // monomials. since they are hierarchic we only need one case block.
73  switch (i)
74  {
75  // constant
76  case 0:
77  return 1.;
78 
79  // linear
80  case 1:
81  return dx;
82 
83  case 2:
84  return dy;
85 
86  // quadratics
87  case 3:
88  return dx*dx;
89 
90  case 4:
91  return dx*dy;
92 
93  case 5:
94  return dy*dy;
95 
96  // cubics
97  case 6:
98  return dx*dx*dx;
99 
100  case 7:
101  return dx*dx*dy;
102 
103  case 8:
104  return dx*dy*dy;
105 
106  case 9:
107  return dy*dy*dy;
108 
109  // quartics
110  case 10:
111  return dx*dx*dx*dx;
112 
113  case 11:
114  return dx*dx*dx*dy;
115 
116  case 12:
117  return dx*dx*dy*dy;
118 
119  case 13:
120  return dx*dy*dy*dy;
121 
122  case 14:
123  return dy*dy*dy*dy;
124 
125  default:
126  unsigned int o = 0;
127  for (; i >= (o+1)*(o+2)/2; o++) { }
128  unsigned int i2 = i - (o*(o+1)/2);
129  Real val = 1.;
130  for (unsigned int index=i2; index != o; index++)
131  val *= dx;
132  for (unsigned int index=0; index != i2; index++)
133  val *= dy;
134  return val;
135  }
136 
137 #else // LIBMESH_DIM <= 1
138  libmesh_assert(true || order || add_p_level);
139  libmesh_ignore(elem, i, point_in);
140  libmesh_not_implemented();
141 #endif
142 }
143 
144 
145 
146 
147 template <>
149  const Order,
150  const unsigned int,
151  const Point &)
152 {
153  libmesh_error_msg("XYZ polynomials require the element.");
154  return 0.;
155 }
156 
157 
158 
159 template <>
161  const Elem * elem,
162  const unsigned int i,
163  const Point & p,
164  const bool add_p_level)
165 {
166  return FE<2,XYZ>::shape(elem, fet.order, i, p, add_p_level);
167 }
168 
169 
170 
171 
172 
173 template <>
175  const Order libmesh_dbg_var(order),
176  const unsigned int i,
177  const unsigned int j,
178  const Point & point_in,
179  const bool libmesh_dbg_var(add_p_level))
180 {
181 #if LIBMESH_DIM > 1
182 
183 
184  libmesh_assert_less (j, 2);
185  libmesh_assert(elem);
186 
187  Point avg = elem->vertex_average();
188  Point max_distance = Point(0.,0.,0.);
189  for (const Point & p : elem->node_ref_range())
190  for (unsigned int d = 0; d < 2; d++)
191  {
192  const Real distance = std::abs(avg(d) - p(d));
193  max_distance(d) = std::max(distance, max_distance(d));
194  }
195 
196  const Real x = point_in(0);
197  const Real y = point_in(1);
198  const Real xc = avg(0);
199  const Real yc = avg(1);
200  const Real distx = max_distance(0);
201  const Real disty = max_distance(1);
202  const Real dx = (x - xc)/distx;
203  const Real dy = (y - yc)/disty;
204 
205 #ifndef NDEBUG
206  // totalorder is only used in the assertion below, so
207  // we avoid declaring it when asserts are not active.
208  const unsigned int totalorder = order + add_p_level * elem->p_level();
209 #endif
210  libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2);
211 
212  // monomials. since they are hierarchic we only need one case block.
213 
214  switch (j)
215  {
216  // d()/dx
217  case 0:
218  {
219  switch (i)
220  {
221  // constants
222  case 0:
223  return 0.;
224 
225  // linears
226  case 1:
227  return 1./distx;
228 
229  case 2:
230  return 0.;
231 
232  // quadratics
233  case 3:
234  return 2.*dx/distx;
235 
236  case 4:
237  return dy/distx;
238 
239  case 5:
240  return 0.;
241 
242  // cubics
243  case 6:
244  return 3.*dx*dx/distx;
245 
246  case 7:
247  return 2.*dx*dy/distx;
248 
249  case 8:
250  return dy*dy/distx;
251 
252  case 9:
253  return 0.;
254 
255  // quartics
256  case 10:
257  return 4.*dx*dx*dx/distx;
258 
259  case 11:
260  return 3.*dx*dx*dy/distx;
261 
262  case 12:
263  return 2.*dx*dy*dy/distx;
264 
265  case 13:
266  return dy*dy*dy/distx;
267 
268  case 14:
269  return 0.;
270 
271  default:
272  unsigned int o = 0;
273  for (; i >= (o+1)*(o+2)/2; o++) { }
274  unsigned int i2 = i - (o*(o+1)/2);
275  Real val = o - i2;
276  for (unsigned int index=i2+1; index < o; index++)
277  val *= dx;
278  for (unsigned int index=0; index != i2; index++)
279  val *= dy;
280  return val/distx;
281  }
282  }
283 
284 
285  // d()/dy
286  case 1:
287  {
288  switch (i)
289  {
290  // constants
291  case 0:
292  return 0.;
293 
294  // linears
295  case 1:
296  return 0.;
297 
298  case 2:
299  return 1./disty;
300 
301  // quadratics
302  case 3:
303  return 0.;
304 
305  case 4:
306  return dx/disty;
307 
308  case 5:
309  return 2.*dy/disty;
310 
311  // cubics
312  case 6:
313  return 0.;
314 
315  case 7:
316  return dx*dx/disty;
317 
318  case 8:
319  return 2.*dx*dy/disty;
320 
321  case 9:
322  return 3.*dy*dy/disty;
323 
324  // quartics
325  case 10:
326  return 0.;
327 
328  case 11:
329  return dx*dx*dx/disty;
330 
331  case 12:
332  return 2.*dx*dx*dy/disty;
333 
334  case 13:
335  return 3.*dx*dy*dy/disty;
336 
337  case 14:
338  return 4.*dy*dy*dy/disty;
339 
340  default:
341  unsigned int o = 0;
342  for (; i >= (o+1)*(o+2)/2; o++) { }
343  unsigned int i2 = i - (o*(o+1)/2);
344  Real val = i2;
345  for (unsigned int index=i2; index != o; index++)
346  val *= dx;
347  for (unsigned int index=1; index <= i2; index++)
348  val *= dy;
349  return val/disty;
350  }
351  }
352 
353 
354  default:
355  libmesh_error_msg("Invalid j = " << j);
356  }
357 
358 #else // LIBMESH_DIM <= 1
359  libmesh_assert(true || order || add_p_level);
360  libmesh_ignore(elem, i, j, point_in);
361  libmesh_not_implemented();
362 #endif
363 }
364 
365 
366 template <>
368  const Order,
369  const unsigned int,
370  const unsigned int,
371  const Point &)
372 {
373  libmesh_error_msg("XYZ polynomials require the element.");
374  return 0.;
375 }
376 
377 
378 template <>
380  const Elem * elem,
381  const unsigned int i,
382  const unsigned int j,
383  const Point & p,
384  const bool add_p_level)
385 {
386  return FE<2,XYZ>::shape_deriv(elem, fet.order, i, j, p, add_p_level);
387 }
388 
389 
390 
391 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
392 
393 template <>
395  const Order libmesh_dbg_var(order),
396  const unsigned int i,
397  const unsigned int j,
398  const Point & point_in,
399  const bool libmesh_dbg_var(add_p_level))
400 {
401 #if LIBMESH_DIM > 1
402 
403  libmesh_assert_less_equal (j, 2);
404  libmesh_assert(elem);
405 
406  Point avg = elem->vertex_average();
407  Point max_distance = Point(0.,0.,0.);
408  for (auto p : make_range(elem->n_nodes()))
409  for (unsigned int d = 0; d < 2; d++)
410  {
411  const Real distance = std::abs(avg(d) - elem->point(p)(d));
412  max_distance(d) = std::max(distance, max_distance(d));
413  }
414 
415  const Real x = point_in(0);
416  const Real y = point_in(1);
417  const Real xc = avg(0);
418  const Real yc = avg(1);
419  const Real distx = max_distance(0);
420  const Real disty = max_distance(1);
421  const Real dx = (x - xc)/distx;
422  const Real dy = (y - yc)/disty;
423  const Real dist2x = pow(distx,2.);
424  const Real dist2y = pow(disty,2.);
425  const Real distxy = distx * disty;
426 
427 #ifndef NDEBUG
428  // totalorder is only used in the assertion below, so
429  // we avoid declaring it when asserts are not active.
430  const unsigned int totalorder = order + add_p_level * elem->p_level();
431 #endif
432  libmesh_assert_less (i, (totalorder+1)*(totalorder+2)/2);
433 
434  // monomials. since they are hierarchic we only need one case block.
435 
436  switch (j)
437  {
438  // d^2()/dx^2
439  case 0:
440  {
441  switch (i)
442  {
443  // constants
444  case 0:
445  // linears
446  case 1:
447  case 2:
448  return 0.;
449 
450  // quadratics
451  case 3:
452  return 2./dist2x;
453 
454  case 4:
455  case 5:
456  return 0.;
457 
458  // cubics
459  case 6:
460  return 6.*dx/dist2x;
461 
462  case 7:
463  return 2.*dy/dist2x;
464 
465  case 8:
466  case 9:
467  return 0.;
468 
469  // quartics
470  case 10:
471  return 12.*dx*dx/dist2x;
472 
473  case 11:
474  return 6.*dx*dy/dist2x;
475 
476  case 12:
477  return 2.*dy*dy/dist2x;
478 
479  case 13:
480  case 14:
481  return 0.;
482 
483  default:
484  unsigned int o = 0;
485  for (; i >= (o+1)*(o+2)/2; o++) { }
486  unsigned int i2 = i - (o*(o+1)/2);
487  Real val = (o - i2) * (o - i2 - 1);
488  for (unsigned int index=i2+2; index < o; index++)
489  val *= dx;
490  for (unsigned int index=0; index != i2; index++)
491  val *= dy;
492  return val/dist2x;
493  }
494  }
495 
496  // d^2()/dxdy
497  case 1:
498  {
499  switch (i)
500  {
501  // constants
502  case 0:
503 
504  // linears
505  case 1:
506  case 2:
507  return 0.;
508 
509  // quadratics
510  case 3:
511  return 0.;
512 
513  case 4:
514  return 1./distxy;
515 
516  case 5:
517  return 0.;
518 
519  // cubics
520  case 6:
521  return 0.;
522  case 7:
523  return 2.*dx/distxy;
524 
525  case 8:
526  return 2.*dy/distxy;
527 
528  case 9:
529  return 0.;
530 
531  // quartics
532  case 10:
533  return 0.;
534 
535  case 11:
536  return 3.*dx*dx/distxy;
537 
538  case 12:
539  return 4.*dx*dy/distxy;
540 
541  case 13:
542  return 3.*dy*dy/distxy;
543 
544  case 14:
545  return 0.;
546 
547  default:
548  unsigned int o = 0;
549  for (; i >= (o+1)*(o+2)/2; o++) { }
550  unsigned int i2 = i - (o*(o+1)/2);
551  Real val = (o - i2) * i2;
552  for (unsigned int index=i2+1; index < o; index++)
553  val *= dx;
554  for (unsigned int index=1; index < i2; index++)
555  val *= dy;
556  return val/distxy;
557  }
558  }
559 
560  // d^2()/dy^2
561  case 2:
562  {
563  switch (i)
564  {
565  // constants
566  case 0:
567 
568  // linears
569  case 1:
570  case 2:
571  return 0.;
572 
573  // quadratics
574  case 3:
575  case 4:
576  return 0.;
577 
578  case 5:
579  return 2./dist2y;
580 
581  // cubics
582  case 6:
583  return 0.;
584 
585  case 7:
586  return 0.;
587 
588  case 8:
589  return 2.*dx/dist2y;
590 
591  case 9:
592  return 6.*dy/dist2y;
593 
594  // quartics
595  case 10:
596  case 11:
597  return 0.;
598 
599  case 12:
600  return 2.*dx*dx/dist2y;
601 
602  case 13:
603  return 6.*dx*dy/dist2y;
604 
605  case 14:
606  return 12.*dy*dy/dist2y;
607 
608  default:
609  unsigned int o = 0;
610  for (; i >= (o+1)*(o+2)/2; o++) { }
611  unsigned int i2 = i - (o*(o+1)/2);
612  Real val = i2 * (i2 - 1);
613  for (unsigned int index=i2; index != o; index++)
614  val *= dx;
615  for (unsigned int index=2; index < i2; index++)
616  val *= dy;
617  return val/dist2y;
618  }
619  }
620 
621  default:
622  libmesh_error_msg("Invalid shape function derivative j = " << j);
623  }
624 
625 #else // LIBMESH_DIM <= 1
626  libmesh_assert(true || order || add_p_level);
627  libmesh_ignore(elem, i, j, point_in);
628  libmesh_not_implemented();
629 #endif
630 }
631 
632 
633 template <>
635  const Order,
636  const unsigned int,
637  const unsigned int,
638  const Point &)
639 {
640  libmesh_error_msg("XYZ polynomials require the element.");
641  return 0.;
642 }
643 
644 
645 
646 template <>
648  const Elem * elem,
649  const unsigned int i,
650  const unsigned int j,
651  const Point & p,
652  const bool add_p_level)
653 {
654  return FE<2,XYZ>::shape_second_deriv(elem, fet.order, i, j, p, add_p_level);
655 }
656 
657 
658 #endif
659 
660 } // namespace libMesh
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
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:3105
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:215
The libMesh namespace provides an interface to certain functionality in the library.
Real distance(const Point &p)
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
void libmesh_ignore(const Args &...)
T pow(const T &x)
Definition: utility.h:328
virtual unsigned int n_nodes() const =0
libmesh_assert(ctx)
SimpleRange< NodeRefIter > node_ref_range()
Returns a range with all nodes of an element, usable in range-based for loops.
Definition: elem.h:2634
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Point & point(const unsigned int i) const
Definition: elem.h:2422
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
Point vertex_average() const
Definition: elem.C:691