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