libMesh
inf_fe.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 
20 // Local includes
21 #include "libmesh/libmesh_config.h"
22 
23 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
24 
25 #include "libmesh/inf_fe.h"
26 #include "libmesh/quadrature_gauss.h"
27 #include "libmesh/elem.h"
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/int_range.h"
30 #include "libmesh/auto_ptr.h"
31 
32 namespace libMesh
33 {
34 
35 
36 
37 // Constructor
38 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
40  FEBase (Dim, fet),
41 
42  _n_total_approx_sf (0),
43  _n_total_qp (0),
44 
45  // initialize the current_fe_type to all the same
46  // values as \p fet (since the FE families and coordinate
47  // map type should not change), but use an invalid order
48  // for the radial part (since this is the only order
49  // that may change!).
50  // the data structures like \p phi etc are not initialized
51  // through the constructor, but through reinit()
52  current_fe_type (FEType(fet.order,
53  fet.family,
55  fet.radial_family,
56  fet.inf_map))
57 
58 {
59  // Sanity checks
60  libmesh_assert_equal_to (T_radial, fe_type.radial_family);
61  libmesh_assert_equal_to (T_map, fe_type.inf_map);
62 
63  // build the base_fe object
64  if (Dim != 1)
65  base_fe = FEBase::build(Dim-1, fet);
66 }
67 
68 
69 
70 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
72 {
73  libmesh_assert(q);
74  libmesh_assert(base_fe);
75 
76  const Order base_int_order = q->get_order();
77  const Order radial_int_order = static_cast<Order>(2 * (static_cast<unsigned int>(fe_type.radial_order.get_order()) + 1) +2);
78  const unsigned int qrule_dim = q->get_dim();
79 
80  if (Dim != 1)
81  {
82  // build a Dim-1 quadrature rule of the type that we received
83  base_qrule = QBase::build(q->type(), qrule_dim-1, base_int_order);
84  base_fe->attach_quadrature_rule(base_qrule.get());
85  }
86 
87  // in radial direction, always use Gauss quadrature
88  radial_qrule = libmesh_make_unique<QGauss>(1, radial_int_order);
89 
90  // Maybe helpful to store the QBase *
91  // with which we initialized our own quadrature rules.
92  // Used e.g. in \p InfFE::reinit(elem,side)
93  qrule = q;
94 }
95 
96 
97 
98 
99 
100 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base>
102 {
103  base_elem.reset(InfFEBase::build_elem(inf_elem));
104 }
105 
106 
107 
108 
109 
110 
111 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
113  const std::vector<Point> * const pts,
114  const std::vector<Real> * const weights)
115 {
116  libmesh_assert(base_fe.get());
117  libmesh_assert(inf_elem);
118 
119  // I don't understand infinite elements well enough to risk
120  // calculating too little. :-( RHS
121  this->calculate_phi = this->calculate_dphi = this->calculate_dphiref = true;
122  this->get_xyz();
123  this->determine_calculations();
124  base_fe->calculate_phi = base_fe->calculate_dphi = base_fe->calculate_dphiref = true;
125  base_fe->get_xyz();
126  base_fe->determine_calculations();
127 
128  if (pts == nullptr)
129  {
130  libmesh_assert(base_fe->qrule);
131  libmesh_assert_equal_to (base_fe->qrule, base_qrule.get());
132  libmesh_assert(radial_qrule.get());
133 
134  bool init_shape_functions_required = false;
135 
136  // init the radial data fields only when the radial order changes
137  if (current_fe_type.radial_order != fe_type.radial_order)
138  {
139  current_fe_type.radial_order = fe_type.radial_order;
140 
141  // Watch out: this call to QBase->init() only works for
142  // current_fe_type = const! To allow variable Order,
143  // the init() of QBase has to be modified...
144  radial_qrule->init(EDGE2);
145 
146  // initialize the radial shape functions
147  this->init_radial_shape_functions(inf_elem);
148 
149  init_shape_functions_required=true;
150  }
151 
152 
153  bool update_base_elem_required=true;
154 
155  // update the type in accordance to the current cell
156  // and reinit if the cell type has changed or (as in
157  // the case of the hierarchics) the shape functions
158  // depend on the particular element and need a reinit
159  if ((Dim != 1) &&
160  ((this->get_type() != inf_elem->type()) ||
161  (base_fe->shapes_need_reinit())))
162  {
163  // store the new element type, update base_elem
164  // here. Through \p update_base_elem_required,
165  // remember whether it has to be updated (see below).
166  elem_type = inf_elem->type();
167  this->update_base_elem(inf_elem);
168  update_base_elem_required=false;
169 
170  // initialize the base quadrature rule for the new element
171  base_qrule->init(base_elem->type());
172 
173  // initialize the shape functions in the base
174  base_fe->init_base_shape_functions(base_fe->qrule->get_points(),
175  base_elem.get());
176 
177  // compute the shape functions and map functions of base_fe
178  // before using them later in combine_base_radial.
179  base_fe->_fe_map->compute_map (base_fe->dim, base_fe->qrule->get_weights(),
180  base_elem.get(), base_fe->calculate_d2phi);
181  base_fe->compute_shape_functions(base_elem.get(), base_fe->qrule->get_points());
182 
183  init_shape_functions_required=true;
184  }
185 
186 
187  // when either the radial or base part change,
188  // we have to init the whole fields
189  if (init_shape_functions_required)
190  this->init_shape_functions (radial_qrule->get_points(),
191  base_fe->qrule->get_points(),
192  inf_elem);
193 
194  // computing the distance only works when we have the current
195  // base_elem stored. This happens when fe_type is const,
196  // the inf_elem->type remains the same. Then we have to
197  // update the base elem _here_.
198  if (update_base_elem_required)
199  this->update_base_elem(inf_elem);
200 
201  // compute dist (depends on geometry, therefore has to be updated for
202  // each and every new element), throw radial and base part together
203  this->combine_base_radial (inf_elem);
204 
205  this->_fe_map->compute_map (this->dim, _total_qrule_weights, inf_elem, this->calculate_d2phi);
206 
207  // Compute the shape functions and the derivatives
208  // at all quadrature points.
209  this->compute_shape_functions (inf_elem,base_fe->qrule->get_points());
210  }
211 
212  else // if pts != nullptr
213  {
214  // update the elem_type
215  elem_type = inf_elem->type();
216 
217  // We'll assume that pts is a tensor product mesh of points.
218  // That will handle the pts.size()==1 case that we care about
219  // right now, and it will generalize a bit, and it won't break
220  // the assumptions elsewhere in InfFE.
221  std::vector<Point> radial_pts;
222  if (pts->size() > 0)
223  {
224  Real radius = (*pts)[0](Dim-1);
225  radial_pts.push_back(Point(radius));
226  unsigned int n_radial_pts=1;
227  unsigned int n_angular_pts=1;
228  for (auto p : IntRange<std::size_t>(1, pts->size()))
229  {
230  radius = (*pts)[p](Dim-1);
231  // check for repetition of radius: The max. allowed distance is somewhat arbitrary
232  // but the given value should not produce false positives...
233  if (abs(radial_pts[p-n_radial_pts*n_angular_pts](0) - radius)< 1e-4)
234  {
235  // should the next one be in the next segment?
236  if (p+1 == n_radial_pts*(n_angular_pts+1))
237  n_angular_pts++;
238  }
239  // if none yet repeated:
240  else if (n_angular_pts == 1)
241  radial_pts.push_back(Point(radius));
242  // if there was repetition but this does not, the assumed
243  // format does not work:
244  else
245  {
246  libmesh_error_msg("We assumed that the "<<pts->size()
247  <<" points are of tensor-product type with "
248  <<n_radial_pts<<" radial points and "
249  <<n_angular_pts<< " angular points."<<std::endl
250  <<"But apparently point "<<p
251  <<" does not fit that scheme: Its radius is "
252  <<radius <<"but should have "
253  <<radial_pts[p-n_radial_pts*n_angular_pts]<<".");
254  }
255  }
256  }
257  const std::size_t radial_pts_size = radial_pts.size();
258  const std::size_t base_pts_size = pts->size() / radial_pts_size;
259  // If we're a tensor product we should have no remainder
260  libmesh_assert_equal_to
261  (base_pts_size * radial_pts_size, pts->size());
262 
263 
264  std::vector<Point> base_pts;
265  base_pts.reserve(base_pts_size);
266  for (std::size_t p=0, ps=pts->size(); p != ps; p += radial_pts_size)
267  {
268  Point pt = (*pts)[p];
269  pt(Dim-1) = 0;
270  base_pts.push_back(pt);
271  }
272 
273  // init radial shapes
274  this->init_radial_shape_functions(inf_elem, &radial_pts);
275 
276  // update the base
277  this->update_base_elem(inf_elem);
278 
279  // the finite element on the ifem base
280  base_fe = FEBase::build(Dim-1, this->fe_type);
281 
282  base_fe->calculate_phi = base_fe->calculate_dphi = base_fe->calculate_dphiref = true;
283  base_fe->get_xyz();
284  base_fe->determine_calculations();
285 
286  // init base shapes
287  base_fe->init_base_shape_functions(base_pts,
288  base_elem.get());
289 
290  // compute the shape functions and map functions of base_fe
291  // before using them later in combine_base_radial.
292 
293  if (weights)
294  {
295  base_fe->_fe_map->compute_map (base_fe->dim, *weights,
296  base_elem.get(), base_fe->calculate_d2phi);
297  }
298  else
299  {
300  std::vector<Real> dummy_weights (pts->size(), 1.);
301  base_fe->_fe_map->compute_map (base_fe->dim, dummy_weights,
302  base_elem.get(), base_fe->calculate_d2phi);
303  }
304 
305  base_fe->compute_shape_functions(base_elem.get(), *pts);
306 
307  this->init_shape_functions (radial_pts, base_pts, inf_elem);
308 
309  // combine the base and radial shapes
310  this->combine_base_radial (inf_elem);
311 
312  // weights
313  if (weights != nullptr)
314  {
315  this->_fe_map->compute_map (this->dim, *weights, inf_elem, this->calculate_d2phi);
316  }
317  else
318  {
319  std::vector<Real> dummy_weights (pts->size(), 1.);
320  this->_fe_map->compute_map (this->dim, dummy_weights, inf_elem, this->calculate_d2phi);
321  }
322 
323  // finally compute the ifem shapes
324  this->compute_shape_functions (inf_elem,*pts);
325  }
326 
327 }
328 
329 
330 
331 
332 
333 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
334 void
336 init_radial_shape_functions(const Elem * libmesh_dbg_var(inf_elem),
337  const std::vector<Point> * radial_pts)
338 {
339  libmesh_assert(radial_qrule.get() || radial_pts);
340  libmesh_assert(inf_elem);
341 
342  // Start logging the radial shape function initialization
343  LOG_SCOPE("init_radial_shape_functions()", "InfFE");
344 
345  // initialize most of the things related to mapping
346 
347  // The order to use in the radial map (currently independent of the element type)
348  const Order radial_mapping_order = InfFERadial::mapping_order();
349  const unsigned int n_radial_mapping_shape_functions =
350  InfFERadial::n_dofs(radial_mapping_order);
351 
352  // initialize most of the things related to physical approximation
353  const Order radial_approx_order = fe_type.radial_order;
354  const unsigned int n_radial_approx_shape_functions =
355  InfFERadial::n_dofs(radial_approx_order);
356 
357  const std::size_t n_radial_qp =
358  radial_pts ? radial_pts->size() : radial_qrule->n_points();
359  const std::vector<Point> & radial_qp =
360  radial_pts ? *radial_pts : radial_qrule->get_points();
361 
362  // resize the radial data fields
363 
364  // the radial polynomials (eval)
365  mode.resize (n_radial_approx_shape_functions);
366  dmodedv.resize (n_radial_approx_shape_functions);
367 
368  // the (1-v)/2 weight
369  som.resize (n_radial_qp);
370  dsomdv.resize (n_radial_qp);
371 
372  // the radial map
373  radial_map.resize (n_radial_mapping_shape_functions);
374  dradialdv_map.resize (n_radial_mapping_shape_functions);
375 
376 
377  for (unsigned int i=0; i<n_radial_mapping_shape_functions; i++)
378  {
379  radial_map[i].resize (n_radial_qp);
380  dradialdv_map[i].resize (n_radial_qp);
381  }
382 
383 
384  for (unsigned int i=0; i<n_radial_approx_shape_functions; i++)
385  {
386  mode[i].resize (n_radial_qp);
387  dmodedv[i].resize (n_radial_qp);
388  }
389 
390 
391  // compute scalar values at radial quadrature points
392  for (std::size_t p=0; p<n_radial_qp; p++)
393  {
394  som[p] = InfFERadial::decay (Dim, radial_qp[p](0));
395  dsomdv[p] = InfFERadial::decay_deriv (radial_qp[p](0));
396  }
397 
398 
399  // evaluate the mode shapes in radial direction at radial quadrature points
400  for (unsigned int i=0; i<n_radial_approx_shape_functions; i++)
401  for (std::size_t p=0; p<n_radial_qp; p++)
402  {
403  mode[i][p] = InfFE<Dim,T_radial,T_map>::eval (radial_qp[p](0), radial_approx_order, i);
404  dmodedv[i][p] = InfFE<Dim,T_radial,T_map>::eval_deriv (radial_qp[p](0), radial_approx_order, i);
405  }
406 
407 
408  // evaluate the mapping functions in radial direction at radial quadrature points
409  for (unsigned int i=0; i<n_radial_mapping_shape_functions; i++)
410  for (std::size_t p=0; p<n_radial_qp; p++)
411  {
412  radial_map[i][p] = InfFEMap::eval (radial_qp[p](0), radial_mapping_order, i);
413  dradialdv_map[i][p] = InfFEMap::eval_deriv (radial_qp[p](0), radial_mapping_order, i);
414  }
415 }
416 
417 
418 
419 
420 
421 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
422 void InfFE<Dim,T_radial,T_map>::init_shape_functions(const std::vector<Point> & radial_qp,
423  const std::vector<Point> & base_qp,
424  const Elem * inf_elem)
425 {
426  libmesh_assert(inf_elem);
427 
428  // Start logging the radial shape function initialization
429  LOG_SCOPE("init_shape_functions()", "InfFE");
430 
431  // fast access to some const ints for the radial data
432  const unsigned int n_radial_mapping_sf = cast_int<unsigned int>(radial_map.size());
433  const unsigned int n_radial_approx_sf = cast_int<unsigned int>(mode.size());
434  const unsigned int n_radial_qp = cast_int<unsigned int>(som.size());
435 
436 
437  // initialize most of the things related to mapping
438 
439  // The element type and order to use in the base map
440  const Order base_mapping_order = base_elem->default_order();
441 
442  // the number of base shape functions used to construct the map
443  // (Lagrange shape functions are used for mapping in the base)
444  unsigned int n_base_mapping_shape_functions =
445  InfFEBase::n_base_mapping_sf(*base_elem,
446  base_mapping_order);
447 
448  const unsigned int n_total_mapping_shape_functions =
449  n_radial_mapping_sf * n_base_mapping_shape_functions;
450 
451  // initialize most of the things related to physical approximation
452  unsigned int n_base_approx_shape_functions;
453  if (Dim > 1)
454  n_base_approx_shape_functions = base_fe->n_shape_functions();
455  else
456  n_base_approx_shape_functions = 1;
457 
458 
459  const unsigned int n_total_approx_shape_functions =
460  n_radial_approx_sf * n_base_approx_shape_functions;
461 
462  // update class member field
463  _n_total_approx_sf = n_total_approx_shape_functions;
464 
465 
466  // The number of the base quadrature points.
467  const unsigned int n_base_qp = cast_int<unsigned int>(base_qp.size());
468 
469  // The total number of quadrature points.
470  const unsigned int n_total_qp = n_radial_qp * n_base_qp;
471 
472 
473  // update class member field
474  _n_total_qp = n_total_qp;
475 
476 
477 
478  // initialize the node and shape numbering maps
479  {
480  // these vectors work as follows: the i-th entry stores
481  // the associated base/radial node number
482  _radial_node_index.resize(n_total_mapping_shape_functions);
483  _base_node_index.resize(n_total_mapping_shape_functions);
484 
485  // similar for the shapes: the i-th entry stores
486  // the associated base/radial shape number
487  _radial_shape_index.resize(n_total_approx_shape_functions);
488  _base_shape_index.resize(n_total_approx_shape_functions);
489 
490  const ElemType inf_elem_type = inf_elem->type();
491 
492  // fill the node index map
493  for (unsigned int n=0; n<n_total_mapping_shape_functions; n++)
494  {
495  compute_node_indices (inf_elem_type,
496  n,
497  _base_node_index[n],
498  _radial_node_index[n]);
499  libmesh_assert_less (_base_node_index[n], n_base_mapping_shape_functions);
500  libmesh_assert_less (_radial_node_index[n], n_radial_mapping_sf);
501  }
502 
503  // fill the shape index map
504  for (unsigned int n=0; n<n_total_approx_shape_functions; n++)
505  {
506  compute_shape_indices (this->fe_type,
507  inf_elem_type,
508  n,
509  _base_shape_index[n],
510  _radial_shape_index[n]);
511  libmesh_assert_less (_base_shape_index[n], n_base_approx_shape_functions);
512  libmesh_assert_less (_radial_shape_index[n], n_radial_approx_sf);
513  }
514  }
515 
516  // resize the base data fields
517  dist.resize(n_base_mapping_shape_functions);
518 
519  // resize the total data fields
520 
521  // the phase term varies with xi, eta and zeta(v): store it for _all_ qp
522  //
523  // when computing the phase, we need the base approximations
524  // therefore, initialize the phase here, but evaluate it
525  // in combine_base_radial().
526  //
527  // the weight, though, is only needed at the radial quadrature points, n_radial_qp.
528  // but for a uniform interface to the protected data fields
529  // the weight data field (which are accessible from the outside) are expanded to n_total_qp.
530  weight.resize (n_total_qp);
531  dweightdv.resize (n_total_qp);
532  dweight.resize (n_total_qp);
533 
534  dphase.resize (n_total_qp);
535  dphasedxi.resize (n_total_qp);
536  dphasedeta.resize (n_total_qp);
537  dphasedzeta.resize (n_total_qp);
538 
539  // this vector contains the integration weights for the combined quadrature rule
540  _total_qrule_weights.resize(n_total_qp);
541 
542  // InfFE's data fields phi, dphi, dphidx, phi_map etc hold the _total_
543  // shape and mapping functions, respectively
544  {
545  phi.resize (n_total_approx_shape_functions);
546  dphi.resize (n_total_approx_shape_functions);
547  dphidx.resize (n_total_approx_shape_functions);
548  dphidy.resize (n_total_approx_shape_functions);
549  dphidz.resize (n_total_approx_shape_functions);
550  dphidxi.resize (n_total_approx_shape_functions);
551 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
552  libmesh_warning("Warning: Second derivatives for Infinite elements"
553  << " are not yet implemented!"
554  << std::endl);
555 
556  d2phi.resize (n_total_approx_shape_functions);
557  d2phidx2.resize (n_total_approx_shape_functions);
558  d2phidxdy.resize (n_total_approx_shape_functions);
559  d2phidxdz.resize (n_total_approx_shape_functions);
560  d2phidy2.resize (n_total_approx_shape_functions);
561  d2phidydz.resize (n_total_approx_shape_functions);
562  d2phidz2.resize (n_total_approx_shape_functions);
563  d2phidxi2.resize (n_total_approx_shape_functions);
564 
565  if (Dim > 1)
566  {
567  d2phidxideta.resize(n_total_approx_shape_functions);
568  d2phideta2.resize(n_total_approx_shape_functions);
569  }
570 
571  if (Dim > 2)
572  {
573  d2phidetadzeta.resize(n_total_approx_shape_functions);
574  d2phidxidzeta.resize(n_total_approx_shape_functions);
575  d2phidzeta2.resize(n_total_approx_shape_functions);
576  }
577 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
578 
579  if (Dim > 1)
580  dphideta.resize(n_total_approx_shape_functions);
581 
582  if (Dim == 3)
583  dphidzeta.resize(n_total_approx_shape_functions);
584 
585  std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
586  std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
587 
588  phi_map.resize(n_total_mapping_shape_functions);
589  dphidxi_map.resize(n_total_mapping_shape_functions);
590 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
591  std::vector<std::vector<Real>> & d2phidxi2_map = this->_fe_map->get_d2phidxi2_map();
592  d2phidxi2_map.resize(n_total_mapping_shape_functions);
593 
594  if (Dim > 1)
595  {
596  std::vector<std::vector<Real>> & d2phidxideta_map = this->_fe_map->get_d2phidxideta_map();
597  std::vector<std::vector<Real>> & d2phideta2_map = this->_fe_map->get_d2phideta2_map();
598  d2phidxideta_map.resize(n_total_mapping_shape_functions);
599  d2phideta2_map.resize(n_total_mapping_shape_functions);
600  }
601 
602  if (Dim == 3)
603  {
604  std::vector<std::vector<Real>> & d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map();
605  std::vector<std::vector<Real>> & d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map();
606  std::vector<std::vector<Real>> & d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map();
607  d2phidxidzeta_map.resize(n_total_mapping_shape_functions);
608  d2phidetadzeta_map.resize(n_total_mapping_shape_functions);
609  d2phidzeta2_map.resize(n_total_mapping_shape_functions);
610  }
611 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
612 
613  if (Dim > 1)
614  {
615  std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
616  dphideta_map.resize(n_total_mapping_shape_functions);
617  }
618 
619  if (Dim == 3)
620  {
621  std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
622  dphidzeta_map.resize(n_total_mapping_shape_functions);
623  }
624  }
625 
626  // collect all the for loops, where inner vectors are
627  // resized to the appropriate number of quadrature points
628  {
629  for (unsigned int i=0; i<n_total_approx_shape_functions; i++)
630  {
631  phi[i].resize (n_total_qp);
632  dphi[i].resize (n_total_qp);
633  dphidx[i].resize (n_total_qp);
634  dphidy[i].resize (n_total_qp);
635  dphidz[i].resize (n_total_qp);
636  dphidxi[i].resize (n_total_qp);
637 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
638  d2phi[i].resize (n_total_qp);
639  d2phidx2[i].resize (n_total_qp);
640  d2phidxdy[i].resize (n_total_qp);
641  d2phidxdz[i].resize (n_total_qp);
642  d2phidy2[i].resize (n_total_qp);
643  d2phidydz[i].resize (n_total_qp);
644  d2phidy2[i].resize (n_total_qp);
645  d2phidxi2[i].resize (n_total_qp);
646 
647  if (Dim > 1)
648  {
649  d2phidxideta[i].resize (n_total_qp);
650  d2phideta2[i].resize (n_total_qp);
651  }
652  if (Dim > 2)
653  {
654  d2phidxidzeta[i].resize (n_total_qp);
655  d2phidetadzeta[i].resize (n_total_qp);
656  d2phidzeta2[i].resize (n_total_qp);
657  }
658 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
659 
660  if (Dim > 1)
661  dphideta[i].resize (n_total_qp);
662 
663  if (Dim == 3)
664  dphidzeta[i].resize (n_total_qp);
665 
666  }
667 
668  for (unsigned int i=0; i<n_total_mapping_shape_functions; i++)
669  {
670  std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
671  std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
672  phi_map[i].resize (n_total_qp);
673  dphidxi_map[i].resize (n_total_qp);
674 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
675  std::vector<std::vector<Real>> & d2phidxi2_map = this->_fe_map->get_d2phidxi2_map();
676  d2phidxi2_map[i].resize (n_total_qp);
677  if (Dim > 1)
678  {
679  std::vector<std::vector<Real>> & d2phidxideta_map = this->_fe_map->get_d2phidxideta_map();
680  std::vector<std::vector<Real>> & d2phideta2_map = this->_fe_map->get_d2phideta2_map();
681  d2phidxideta_map[i].resize (n_total_qp);
682  d2phideta2_map[i].resize (n_total_qp);
683  }
684 
685  if (Dim > 2)
686  {
687  std::vector<std::vector<Real>> & d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map();
688  std::vector<std::vector<Real>> & d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map();
689  std::vector<std::vector<Real>> & d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map();
690  d2phidxidzeta_map[i].resize (n_total_qp);
691  d2phidetadzeta_map[i].resize (n_total_qp);
692  d2phidzeta2_map[i].resize (n_total_qp);
693  }
694 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
695 
696  if (Dim > 1)
697  {
698  std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
699  dphideta_map[i].resize (n_total_qp);
700  }
701 
702  if (Dim == 3)
703  {
704  std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
705  dphidzeta_map[i].resize (n_total_qp);
706  }
707  }
708  }
709 
710 
711 
712  {
713  // (a) compute scalar values at _all_ quadrature points -- for uniform
714  // access from the outside to these fields
715  // (b) form a std::vector<Real> which contains the appropriate weights
716  // of the combined quadrature rule!
717  libmesh_assert_equal_to (radial_qp.size(), n_radial_qp);
718 
719  if (radial_qrule && base_qrule)
720  {
721  const std::vector<Real> & radial_qw = radial_qrule->get_weights();
722  const std::vector<Real> & base_qw = base_qrule->get_weights();
723  libmesh_assert_equal_to (radial_qw.size(), n_radial_qp);
724  libmesh_assert_equal_to (base_qw.size(), n_base_qp);
725 
726  for (unsigned int rp=0; rp<n_radial_qp; rp++)
727  for (unsigned int bp=0; bp<n_base_qp; bp++)
728  {
729  weight[bp + rp*n_base_qp] = InfFERadial::D(radial_qp[rp](0));
730  dweightdv[bp + rp*n_base_qp] = InfFERadial::D_deriv(radial_qp[rp](0));
731  _total_qrule_weights[bp + rp*n_base_qp] = radial_qw[rp] * base_qw[bp];
732  }
733  }
734  else
735  {
736  for (unsigned int rp=0; rp<n_radial_qp; rp++)
737  for (unsigned int bp=0; bp<n_base_qp; bp++)
738  {
739  weight[bp + rp*n_base_qp] = InfFERadial::D(radial_qp[rp](0));
740  dweightdv[bp + rp*n_base_qp] = InfFERadial::D_deriv(radial_qp[rp](0));
741  }
742  }
743  }
744 }
745 
746 
747 
748 
749 
750 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
752 {
753  libmesh_assert(inf_elem);
754  // at least check whether the base element type is correct.
755  // otherwise this version of computing dist would give problems
756  libmesh_assert_equal_to (base_elem->type(),
757  InfFEBase::get_elem_type(inf_elem->type()));
758 
759  // Start logging the combination of radial and base parts
760  LOG_SCOPE("combine_base_radial()", "InfFE");
761 
762  // zero the phase, since it is to be summed up
763  std::fill (dphasedxi.begin(), dphasedxi.end(), 0.);
764  std::fill (dphasedeta.begin(), dphasedeta.end(), 0.);
765  std::fill (dphasedzeta.begin(), dphasedzeta.end(), 0.);
766 
767 
768  const unsigned int n_base_mapping_sf = cast_int<unsigned int>(dist.size());
769  const Point origin = inf_elem->origin();
770 
771  // for each new infinite element, compute the radial distances
772  for (unsigned int n=0; n<n_base_mapping_sf; n++)
773  dist[n] = Point(base_elem->point(n) - origin).norm();
774 
775 
776  switch (Dim)
777  {
778  // 1D
779  case 1:
780  {
781  libmesh_not_implemented();
782  break;
783  }
784 
785  // 2D
786  case 2:
787  {
788  libmesh_not_implemented();
789  break;
790  }
791 
792  // 3D
793  case 3:
794  {
795  // fast access to the approximation and mapping shapes of base_fe
796  const std::vector<std::vector<Real>> & S = base_fe->phi;
797  const std::vector<std::vector<Real>> & Ss = base_fe->dphidxi;
798  const std::vector<std::vector<Real>> & St = base_fe->dphideta;
799  const std::vector<std::vector<Real>> & S_map = (base_fe->get_fe_map()).get_phi_map();
800  const std::vector<std::vector<Real>> & Ss_map = (base_fe->get_fe_map()).get_dphidxi_map();
801  const std::vector<std::vector<Real>> & St_map = (base_fe->get_fe_map()).get_dphideta_map();
802 
803  const unsigned int n_radial_qp = cast_int<unsigned int>(som.size());
804  if (radial_qrule)
805  libmesh_assert_equal_to(n_radial_qp, radial_qrule->n_points());
806  const unsigned int n_base_qp = cast_int<unsigned int>(S_map[0].size());
807  if (base_qrule)
808  libmesh_assert_equal_to(n_base_qp, base_qrule->n_points());
809 
810  const unsigned int n_total_mapping_sf =
811  cast_int<unsigned int>(radial_map.size()) * n_base_mapping_sf;
812 
813  const unsigned int n_total_approx_sf =
814  InfFERadial::n_dofs(fe_type.radial_order) *
815  base_fe->n_shape_functions();
816 
817 
818  // compute the phase term derivatives
819  {
820  unsigned int tp=0;
821  for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qps
822  for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qps
823  {
824  // sum over all base shapes, to get the average distance
825  for (unsigned int i=0; i<n_base_mapping_sf; i++)
826  {
827  dphasedxi[tp] += Ss_map[i][bp] * dist[i] * radial_map [1][rp];
828  dphasedeta[tp] += St_map[i][bp] * dist[i] * radial_map [1][rp];
829  dphasedzeta[tp] += S_map [i][bp] * dist[i] * dradialdv_map[1][rp];
830  }
831 
832  tp++;
833 
834  } // loop radial and base qps
835  }
836 
837  libmesh_assert_equal_to (phi.size(), n_total_approx_sf);
838  libmesh_assert_equal_to (dphidxi.size(), n_total_approx_sf);
839  libmesh_assert_equal_to (dphideta.size(), n_total_approx_sf);
840  libmesh_assert_equal_to (dphidzeta.size(), n_total_approx_sf);
841 
842  // compute the overall approximation shape functions,
843  // pick the appropriate radial and base shapes through using
844  // _base_shape_index and _radial_shape_index
845  for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qps
846  for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qps
847  for (unsigned int ti=0; ti<n_total_approx_sf; ti++) // over _all_ approx_sf
848  {
849  // let the index vectors take care of selecting the appropriate base/radial shape
850  const unsigned int bi = _base_shape_index [ti];
851  const unsigned int ri = _radial_shape_index[ti];
852  phi [ti][bp+rp*n_base_qp] = S [bi][bp] * mode[ri][rp] * som[rp];
853  dphidxi [ti][bp+rp*n_base_qp] = Ss[bi][bp] * mode[ri][rp] * som[rp];
854  dphideta [ti][bp+rp*n_base_qp] = St[bi][bp] * mode[ri][rp] * som[rp];
855  dphidzeta[ti][bp+rp*n_base_qp] = S [bi][bp]
856  * (dmodedv[ri][rp] * som[rp] + mode[ri][rp] * dsomdv[rp]);
857  }
858 
859  std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
860  std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
861  std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
862  std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
863 
864  libmesh_assert_equal_to (phi_map.size(), n_total_mapping_sf);
865  libmesh_assert_equal_to (dphidxi_map.size(), n_total_mapping_sf);
866  libmesh_assert_equal_to (dphideta_map.size(), n_total_mapping_sf);
867  libmesh_assert_equal_to (dphidzeta_map.size(), n_total_mapping_sf);
868 
869  // compute the overall mapping functions,
870  // pick the appropriate radial and base entries through using
871  // _base_node_index and _radial_node_index
872  for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qps
873  for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qps
874  for (unsigned int ti=0; ti<n_total_mapping_sf; ti++) // over all mapping shapes
875  {
876  // let the index vectors take care of selecting the appropriate base/radial mapping shape
877  const unsigned int bi = _base_node_index [ti];
878  const unsigned int ri = _radial_node_index[ti];
879  phi_map [ti][bp+rp*n_base_qp] = S_map [bi][bp] * radial_map [ri][rp];
880  dphidxi_map [ti][bp+rp*n_base_qp] = Ss_map[bi][bp] * radial_map [ri][rp];
881  dphideta_map [ti][bp+rp*n_base_qp] = St_map[bi][bp] * radial_map [ri][rp];
882  dphidzeta_map[ti][bp+rp*n_base_qp] = S_map [bi][bp] * dradialdv_map[ri][rp];
883  }
884 
885  break;
886  }
887 
888  default:
889  libmesh_error_msg("Unsupported Dim = " << Dim);
890  }
891 }
892 
893 
894 
895 
896 
897 
898 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
900  const std::vector<Point> &)
901 {
902  // Start logging the overall computation of shape functions
903  LOG_SCOPE("compute_shape_functions()", "InfFE");
904 
905  const unsigned int n_total_qp = _n_total_qp;
906 
907  // Compute the shape function values (and derivatives)
908  // at the Quadrature points. Note that the actual values
909  // have already been computed via init_shape_functions
910 
911  // Compute the value of the derivative shape function i at quadrature point p
912  switch (dim)
913  {
914 
915  case 1:
916  {
917  libmesh_not_implemented();
918  break;
919  }
920 
921  case 2:
922  {
923  libmesh_not_implemented();
924  break;
925  }
926 
927  case 3:
928  {
929  const std::vector<Real> & dxidx_map = this->_fe_map->get_dxidx();
930  const std::vector<Real> & dxidy_map = this->_fe_map->get_dxidy();
931  const std::vector<Real> & dxidz_map = this->_fe_map->get_dxidz();
932 
933  const std::vector<Real> & detadx_map = this->_fe_map->get_detadx();
934  const std::vector<Real> & detady_map = this->_fe_map->get_detady();
935  const std::vector<Real> & detadz_map = this->_fe_map->get_detadz();
936 
937  const std::vector<Real> & dzetadx_map = this->_fe_map->get_dzetadx();
938  const std::vector<Real> & dzetady_map = this->_fe_map->get_dzetady();
939  const std::vector<Real> & dzetadz_map = this->_fe_map->get_dzetadz();
940 
941  // These are _all_ shape functions of this infinite element
942  for (auto i : index_range(phi))
943  for (unsigned int p=0; p<n_total_qp; p++)
944  {
945  // dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) + (dphi/dzeta)*(dzeta/dx);
946  dphi[i][p](0) =
947  dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] +
948  dphideta[i][p]*detadx_map[p] +
949  dphidzeta[i][p]*dzetadx_map[p]);
950 
951  // dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) + (dphi/dzeta)*(dzeta/dy);
952  dphi[i][p](1) =
953  dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] +
954  dphideta[i][p]*detady_map[p] +
955  dphidzeta[i][p]*dzetady_map[p]);
956 
957  // dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) + (dphi/dzeta)*(dzeta/dz);
958  dphi[i][p](2) =
959  dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] +
960  dphideta[i][p]*detadz_map[p] +
961  dphidzeta[i][p]*dzetadz_map[p]);
962  }
963 
964 
965  // This is the derivative of the phase term of this infinite element
966  for (unsigned int p=0; p<n_total_qp; p++)
967  {
968  // the derivative of the phase term
969  dphase[p](0) = (dphasedxi[p] * dxidx_map[p] +
970  dphasedeta[p] * detadx_map[p] +
971  dphasedzeta[p] * dzetadx_map[p]);
972 
973  dphase[p](1) = (dphasedxi[p] * dxidy_map[p] +
974  dphasedeta[p] * detady_map[p] +
975  dphasedzeta[p] * dzetady_map[p]);
976 
977  dphase[p](2) = (dphasedxi[p] * dxidz_map[p] +
978  dphasedeta[p] * detadz_map[p] +
979  dphasedzeta[p] * dzetadz_map[p]);
980 
981  // the derivative of the radial weight - varies only in radial direction,
982  // therefore dweightdxi = dweightdeta = 0.
983  dweight[p](0) = dweightdv[p] * dzetadx_map[p];
984 
985  dweight[p](1) = dweightdv[p] * dzetady_map[p];
986 
987  dweight[p](2) = dweightdv[p] * dzetadz_map[p];
988  }
989 
990  break;
991  }
992 
993  default:
994  libmesh_error_msg("Unsupported dim = " << dim);
995  }
996 }
997 
998 
999 
1000 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
1002 {
1003  return false;
1004 }
1005 
1006 } // namespace libMesh
1007 
1008 
1009 // Explicit instantiations
1010 #include "libmesh/inf_fe_instantiate_1D.h"
1011 #include "libmesh/inf_fe_instantiate_2D.h"
1012 #include "libmesh/inf_fe_instantiate_3D.h"
1013 
1014 
1015 
1016 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
libMesh::FEType::inf_map
InfMapType inf_map
The coordinate mapping type of the infinite element.
Definition: fe_type.h:257
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::InfFERadial::decay
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:850
libMesh::FEType::radial_family
FEFamily radial_family
The type of approximation in radial direction.
Definition: fe_type.h:249
libMesh::InfFERadial::mapping_order
static Order mapping_order()
Definition: inf_fe.h:93
libMesh::InfFE::update_base_elem
void update_base_elem(const Elem *inf_elem)
Updates the protected member base_elem to the appropriate base element for the given inf_elem.
Definition: inf_fe.C:101
libMesh::InfFE::InfFE
friend class InfFE
Make all InfFE<Dim,T_radial,T_map> classes friends of each other, so that the protected eval() may be...
Definition: inf_fe.h:841
libMesh::InfFE::eval_deriv
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
libMesh::FEAbstract::fe_type
FEType fe_type
The finite element type for this object.
Definition: fe_abstract.h:585
libMesh::QBase::get_dim
unsigned int get_dim() const
Definition: quadrature.h:135
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::InfFE::combine_base_radial
void combine_base_radial(const Elem *inf_elem)
Combines the shape functions, which were formed in init_shape_functions(Elem *), with geometric data.
Definition: inf_fe.C:751
libMesh::InfFERadial::D_deriv
static Real D_deriv(const Real v)
Definition: inf_fe.h:87
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::QBase::build
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
Definition: quadrature_build.C:41
libMesh::InfFE::shapes_need_reinit
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:1001
libMesh::InfFE::init_shape_functions
void init_shape_functions(const std::vector< Point > &radial_qp, const std::vector< Point > &base_qp, const Elem *inf_elem)
Initialize all the data fields like weight, mode, phi, dphidxi, dphideta, dphidzeta,...
Definition: inf_fe.C:422
libMesh::InfFERadial::decay_deriv
static Real decay_deriv(const Real)
Definition: inf_fe.h:76
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::InfFEMap::eval
static Real eval(Real v, Order o, unsigned int i)
Definition: inf_fe_map_eval.C:27
libMesh::InfFEMap::eval_deriv
static Real eval_deriv(Real v, Order o, unsigned int i)
Definition: inf_fe_map_eval.C:45
libMesh::InfFE::compute_shape_functions
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &) override
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: inf_fe.C:899
radius
const Real radius
Definition: subdomains_ex3.C:48
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::InfFE::reinit
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
Definition: inf_fe.C:112
libMesh::InfFERadial::n_dofs
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:109
libMesh::INVALID_ORDER
Definition: enum_order.h:86
libMesh::InfFE::base_fe
std::unique_ptr< FEBase > base_fe
Have a FE<Dim-1,T_base> handy for base approximation.
Definition: inf_fe.h:794
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::InfFE::attach_quadrature_rule
virtual void attach_quadrature_rule(QBase *q) override
The use of quadrature rules with the InfFE class is somewhat different from the approach of the FE cl...
Definition: inf_fe.C:71
libMesh::InfFEBase::build_elem
static Elem * build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
Definition: inf_fe_base_radial.C:35
libMesh::QBase::type
virtual QuadratureType type() const =0
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
std::norm
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::Elem::origin
virtual Point origin() const
Definition: elem.h:1593
libMesh::InfFE::eval
static Real eval(Real v, Order o_radial, unsigned int i)
libMesh::InfFE::init_radial_shape_functions
void init_radial_shape_functions(const Elem *inf_elem, const std::vector< Point > *radial_pts=nullptr)
Some of the member data only depend on the radial part of the infinite element.
Definition: inf_fe.C:336
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::QBase::get_order
Order get_order() const
Definition: quadrature.h:211
libMesh::MeshTools::weight
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:236
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::InfFERadial::D
static Real D(const Real v)
Definition: inf_fe.h:82
libMesh::InfFEBase::get_elem_type
static ElemType get_elem_type(const ElemType type)
Definition: inf_fe_base_radial.C:49
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::InfFEBase::n_base_mapping_sf
static unsigned int n_base_mapping_sf(const Elem &base_elem, const Order base_mapping_order)
Definition: inf_fe_base_radial.C:95