fe.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 
20 // Local includes
21 #include "libmesh/elem.h"
22 #include "libmesh/fe.h"
23 #include "libmesh/fe_interface.h"
24 #include "libmesh/fe_macro.h"
26 #include "libmesh/quadrature.h"
27 #include "libmesh/tensor_value.h"
28 #include "libmesh/enum_elem_type.h"
29 
30 namespace libMesh
31 {
32 
33 
34 // ------------------------------------------------------------
35 // FE class members
36 template <unsigned int Dim, FEFamily T>
37 FE<Dim,T>::FE (const FEType & fet) :
38  FEGenericBase<typename FEOutputType<T>::type> (Dim,fet),
39  last_side(INVALID_ELEM),
40  last_edge(libMesh::invalid_uint)
41 {
42  // Sanity check. Make sure the
43  // Family specified in the template instantiation
44  // matches the one in the FEType object
45  libmesh_assert_equal_to (T, this->get_family());
46 }
47 
48 
49 template <unsigned int Dim, FEFamily T>
50 unsigned int FE<Dim,T>::n_shape_functions () const
51 {
52  return FE<Dim,T>::n_dofs (this->elem_type,
53  static_cast<Order>(this->fe_type.order + this->_p_level));
54 }
55 
56 
57 template <unsigned int Dim, FEFamily T>
59 {
60  libmesh_assert(q);
61  this->qrule = q;
62  // make sure we don't cache results from a previous quadrature rule
63  this->elem_type = INVALID_ELEM;
64  return;
65 }
66 
67 
68 template <unsigned int Dim, FEFamily T>
69 unsigned int FE<Dim,T>::n_quadrature_points () const
70 {
71  libmesh_assert(this->qrule);
72  return this->qrule->n_points();
73 }
74 
75 
76 template <unsigned int Dim, FEFamily T>
77 void FE<Dim,T>::dofs_on_side(const Elem * const elem,
78  const Order o,
79  unsigned int s,
80  std::vector<unsigned int> & di)
81 {
82  libmesh_assert(elem);
83  libmesh_assert_less (s, elem->n_sides());
84 
85  di.clear();
86  unsigned int nodenum = 0;
87  const unsigned int n_nodes = elem->n_nodes();
88  for (unsigned int n = 0; n != n_nodes; ++n)
89  {
90  const unsigned int n_dofs = n_dofs_at_node(elem->type(),
91  static_cast<Order>(o + elem->p_level()), n);
92  if (elem->is_node_on_side(n, s))
93  for (unsigned int i = 0; i != n_dofs; ++i)
94  di.push_back(nodenum++);
95  else
96  nodenum += n_dofs;
97  }
98 }
99 
100 
101 
102 template <unsigned int Dim, FEFamily T>
103 void FE<Dim,T>::dofs_on_edge(const Elem * const elem,
104  const Order o,
105  unsigned int e,
106  std::vector<unsigned int> & di)
107 {
108  libmesh_assert(elem);
109  libmesh_assert_less (e, elem->n_edges());
110 
111  di.clear();
112  unsigned int nodenum = 0;
113  const unsigned int n_nodes = elem->n_nodes();
114  for (unsigned int n = 0; n != n_nodes; ++n)
115  {
116  const unsigned int n_dofs = n_dofs_at_node(elem->type(),
117  static_cast<Order>(o + elem->p_level()), n);
118  if (elem->is_node_on_edge(n, e))
119  for (unsigned int i = 0; i != n_dofs; ++i)
120  di.push_back(nodenum++);
121  else
122  nodenum += n_dofs;
123  }
124 }
125 
126 
127 
128 template <unsigned int Dim, FEFamily T>
129 void FE<Dim,T>::reinit(const Elem * elem,
130  const std::vector<Point> * const pts,
131  const std::vector<Real> * const weights)
132 {
133  // We can be called with no element. If we're evaluating SCALAR
134  // dofs we'll still have work to do.
135  // libmesh_assert(elem);
136 
137  // We're calculating now! Time to determine what.
138  this->determine_calculations();
139 
140  // Try to avoid calling init_shape_functions
141  // even when shapes_need_reinit
142  bool cached_nodes_still_fit = false;
143 
144  // Most of the hard work happens when we have an actual element
145  if (elem)
146  {
147  // Initialize the shape functions at the user-specified
148  // points
149  if (pts != nullptr)
150  {
151  // Set the type and p level for this element
152  this->elem_type = elem->type();
153  this->_p_level = elem->p_level();
154 
155  // Initialize the shape functions
156  this->_fe_map->template init_reference_to_physical_map<Dim>
157  (*pts, elem);
158  this->init_shape_functions (*pts, elem);
159 
160  // The shape functions do not correspond to the qrule
161  this->shapes_on_quadrature = false;
162  }
163 
164  // If there are no user specified points, we use the
165  // quadrature rule
166 
167  // update the type in accordance to the current cell
168  // and reinit if the cell type has changed or (as in
169  // the case of the hierarchics) the shape functions need
170  // reinit, since they depend on the particular element shape
171  else
172  {
173  libmesh_assert(this->qrule);
174  this->qrule->init(elem->type(), elem->p_level());
175 
176  if (this->qrule->shapes_need_reinit())
177  this->shapes_on_quadrature = false;
178 
179  if (this->elem_type != elem->type() ||
180  this->_p_level != elem->p_level() ||
181  !this->shapes_on_quadrature)
182  {
183  // Set the type and p level for this element
184  this->elem_type = elem->type();
185  this->_p_level = elem->p_level();
186  // Initialize the shape functions
187  this->_fe_map->template init_reference_to_physical_map<Dim>
188  (this->qrule->get_points(), elem);
189  this->init_shape_functions (this->qrule->get_points(), elem);
190 
191  if (this->shapes_need_reinit())
192  {
193  cached_nodes.resize(elem->n_nodes());
194  for (unsigned int n = 0; n != elem->n_nodes(); ++n)
195  {
196  cached_nodes[n] = elem->point(n);
197  }
198  }
199  }
200  else
201  {
202  // libmesh_assert_greater (elem->n_nodes(), 1);
203 
204  cached_nodes_still_fit = true;
205  if (cached_nodes.size() != elem->n_nodes())
206  cached_nodes_still_fit = false;
207  else
208  for (unsigned int n = 1; n < elem->n_nodes(); ++n)
209  {
210  if (!(elem->point(n) - elem->point(0)).relative_fuzzy_equals
211  ((cached_nodes[n] - cached_nodes[0]), 1e-13))
212  {
213  cached_nodes_still_fit = false;
214  break;
215  }
216  }
217 
218  if (this->shapes_need_reinit() && !cached_nodes_still_fit)
219  {
220  this->_fe_map->template init_reference_to_physical_map<Dim>
221  (this->qrule->get_points(), elem);
222  this->init_shape_functions (this->qrule->get_points(), elem);
223  cached_nodes.resize(elem->n_nodes());
224  for (unsigned int n = 0; n != elem->n_nodes(); ++n)
225  cached_nodes[n] = elem->point(n);
226  }
227  }
228 
229  // The shape functions correspond to the qrule
230  this->shapes_on_quadrature = true;
231  }
232  }
233  else // With no defined elem, so mapping or caching to
234  // be done, and our "quadrature rule" is one point for nonlocal
235  // (SCALAR) variables and zero points for local variables.
236  {
237  this->elem_type = INVALID_ELEM;
238  this->_p_level = 0;
239 
240  if (!pts)
241  {
242  if (T == SCALAR)
243  {
244  this->qrule->get_points() =
245  std::vector<Point>(1,Point(0));
246 
247  this->qrule->get_weights() =
248  std::vector<Real>(1,1);
249  }
250  else
251  {
252  this->qrule->get_points().clear();
253  this->qrule->get_weights().clear();
254  }
255 
256  this->init_shape_functions (this->qrule->get_points(), elem);
257  }
258  else
259  this->init_shape_functions (*pts, elem);
260  }
261 
262  // Compute the map for this element.
263  if (pts != nullptr)
264  {
265  if (weights != nullptr)
266  {
267  this->_fe_map->compute_map (this->dim, *weights, elem, this->calculate_d2phi);
268  }
269  else
270  {
271  std::vector<Real> dummy_weights (pts->size(), 1.);
272  this->_fe_map->compute_map (this->dim, dummy_weights, elem, this->calculate_d2phi);
273  }
274  }
275  else
276  {
277  this->_fe_map->compute_map (this->dim, this->qrule->get_weights(), elem, this->calculate_d2phi);
278  }
279 
280  // Compute the shape functions and the derivatives at all of the
281  // quadrature points.
282  if (!cached_nodes_still_fit)
283  {
284  if (pts != nullptr)
285  this->compute_shape_functions (elem,*pts);
286  else
287  this->compute_shape_functions(elem,this->qrule->get_points());
288  }
289 }
290 
291 
292 
293 template <unsigned int Dim, FEFamily T>
294 void FE<Dim,T>::init_shape_functions(const std::vector<Point> & qp,
295  const Elem * elem)
296 {
297  // Start logging the shape function initialization
298  LOG_SCOPE("init_shape_functions()", "FE");
299 
300  // The number of quadrature points.
301  const unsigned int n_qp = cast_int<unsigned int>(qp.size());
302 
303  // Number of shape functions in the finite element approximation
304  // space.
305  const unsigned int n_approx_shape_functions =
306  this->n_shape_functions(this->get_type(),
307  this->get_order());
308 
309  // resize the vectors to hold current data
310  // Phi are the shape functions used for the FE approximation
311  // Phi_map are the shape functions used for the FE mapping
312  if (this->calculate_phi)
313  this->phi.resize (n_approx_shape_functions);
314 
315  if (this->calculate_dphi)
316  {
317  this->dphi.resize (n_approx_shape_functions);
318  this->dphidx.resize (n_approx_shape_functions);
319  this->dphidy.resize (n_approx_shape_functions);
320  this->dphidz.resize (n_approx_shape_functions);
321  }
322 
323  if (this->calculate_dphiref)
324  {
325  if (Dim > 0)
326  this->dphidxi.resize (n_approx_shape_functions);
327 
328  if (Dim > 1)
329  this->dphideta.resize (n_approx_shape_functions);
330 
331  if (Dim > 2)
332  this->dphidzeta.resize (n_approx_shape_functions);
333  }
334 
335  if (this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
336  this->curl_phi.resize(n_approx_shape_functions);
337 
338  if (this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
339  this->div_phi.resize(n_approx_shape_functions);
340 
341 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
342  if (this->calculate_d2phi)
343  {
344  this->d2phi.resize (n_approx_shape_functions);
345  this->d2phidx2.resize (n_approx_shape_functions);
346  this->d2phidxdy.resize (n_approx_shape_functions);
347  this->d2phidxdz.resize (n_approx_shape_functions);
348  this->d2phidy2.resize (n_approx_shape_functions);
349  this->d2phidydz.resize (n_approx_shape_functions);
350  this->d2phidz2.resize (n_approx_shape_functions);
351 
352  if (Dim > 0)
353  this->d2phidxi2.resize (n_approx_shape_functions);
354 
355  if (Dim > 1)
356  {
357  this->d2phidxideta.resize (n_approx_shape_functions);
358  this->d2phideta2.resize (n_approx_shape_functions);
359  }
360  if (Dim > 2)
361  {
362  this->d2phidxidzeta.resize (n_approx_shape_functions);
363  this->d2phidetadzeta.resize (n_approx_shape_functions);
364  this->d2phidzeta2.resize (n_approx_shape_functions);
365  }
366  }
367 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
368 
369  for (unsigned int i=0; i<n_approx_shape_functions; i++)
370  {
371  if (this->calculate_phi)
372  this->phi[i].resize (n_qp);
373  if (this->calculate_dphi)
374  {
375  this->dphi[i].resize (n_qp);
376  this->dphidx[i].resize (n_qp);
377  this->dphidy[i].resize (n_qp);
378  this->dphidz[i].resize (n_qp);
379  }
380 
381  if (this->calculate_dphiref)
382  {
383  if (Dim > 0)
384  this->dphidxi[i].resize(n_qp);
385 
386  if (Dim > 1)
387  this->dphideta[i].resize(n_qp);
388 
389  if (Dim > 2)
390  this->dphidzeta[i].resize(n_qp);
391  }
392 
393  if (this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
394  this->curl_phi[i].resize(n_qp);
395 
396  if (this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
397  this->div_phi[i].resize(n_qp);
398 
399 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
400  if (this->calculate_d2phi)
401  {
402  this->d2phi[i].resize (n_qp);
403  this->d2phidx2[i].resize (n_qp);
404  this->d2phidxdy[i].resize (n_qp);
405  this->d2phidxdz[i].resize (n_qp);
406  this->d2phidy2[i].resize (n_qp);
407  this->d2phidydz[i].resize (n_qp);
408  this->d2phidz2[i].resize (n_qp);
409  if (Dim > 0)
410  this->d2phidxi2[i].resize (n_qp);
411  if (Dim > 1)
412  {
413  this->d2phidxideta[i].resize (n_qp);
414  this->d2phideta2[i].resize (n_qp);
415  }
416  if (Dim > 2)
417  {
418  this->d2phidxidzeta[i].resize (n_qp);
419  this->d2phidetadzeta[i].resize (n_qp);
420  this->d2phidzeta2[i].resize (n_qp);
421  }
422  }
423 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
424  }
425 
426 
427 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
428  //------------------------------------------------------------
429  // Initialize the data fields, which should only be used for infinite
430  // elements, to some sensible values, so that using a FE with the
431  // variational formulation of an InfFE, correct element matrices are
432  // returned
433 
434  {
435  this->weight.resize (n_qp);
436  this->dweight.resize (n_qp);
437  this->dphase.resize (n_qp);
438 
439  for (unsigned int p=0; p<n_qp; p++)
440  {
441  this->weight[p] = 1.;
442  this->dweight[p].zero();
443  this->dphase[p].zero();
444  }
445 
446  }
447 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
448 
449  switch (Dim)
450  {
451 
452  //------------------------------------------------------------
453  // 0D
454  case 0:
455  {
456  break;
457  }
458 
459  //------------------------------------------------------------
460  // 1D
461  case 1:
462  {
463  // Compute the value of the approximation shape function i at quadrature point p
464  if (this->calculate_dphiref)
465  for (unsigned int i=0; i<n_approx_shape_functions; i++)
466  for (unsigned int p=0; p<n_qp; p++)
467  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
468 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
469  if (this->calculate_d2phi)
470  for (unsigned int i=0; i<n_approx_shape_functions; i++)
471  for (unsigned int p=0; p<n_qp; p++)
472  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
473 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
474 
475  break;
476  }
477 
478 
479 
480  //------------------------------------------------------------
481  // 2D
482  case 2:
483  {
484  // Compute the value of the approximation shape function i at quadrature point p
485  if (this->calculate_dphiref)
486  for (unsigned int i=0; i<n_approx_shape_functions; i++)
487  for (unsigned int p=0; p<n_qp; p++)
488  {
489  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
490  this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]);
491  }
492 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
493  if (this->calculate_d2phi)
494  for (unsigned int i=0; i<n_approx_shape_functions; i++)
495  for (unsigned int p=0; p<n_qp; p++)
496  {
497  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
498  this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]);
499  this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]);
500  }
501 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
502 
503 
504  break;
505  }
506 
507 
508 
509  //------------------------------------------------------------
510  // 3D
511  case 3:
512  {
513  // Compute the value of the approximation shape function i at quadrature point p
514  if (this->calculate_dphiref)
515  for (unsigned int i=0; i<n_approx_shape_functions; i++)
516  for (unsigned int p=0; p<n_qp; p++)
517  {
518  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
519  this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]);
520  this->dphidzeta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 2, qp[p]);
521  }
522 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
523  if (this->calculate_d2phi)
524  for (unsigned int i=0; i<n_approx_shape_functions; i++)
525  for (unsigned int p=0; p<n_qp; p++)
526  {
527  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
528  this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]);
529  this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]);
530  this->d2phidxidzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 3, qp[p]);
531  this->d2phidetadzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 4, qp[p]);
532  this->d2phidzeta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 5, qp[p]);
533  }
534 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
535 
536  break;
537  }
538 
539 
540  default:
541  libmesh_error_msg("Invalid dimension Dim = " << Dim);
542  }
543 }
544 
545 
546 
547 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
548 
549 template <unsigned int Dim, FEFamily T>
550 void FE<Dim,T>::init_base_shape_functions(const std::vector<Point> & qp,
551  const Elem * e)
552 {
553  this->elem_type = e->type();
554  this->_fe_map->template init_reference_to_physical_map<Dim>(qp, e);
555  init_shape_functions(qp, e);
556 }
557 
558 #endif // LIBMESH_ENABLE_INFINITE_ELEMENTS
559 
560 
561 
562 //--------------------------------------------------------------
563 // Explicit instantiations using macro from fe_macro.h
564 
565 INSTANTIATE_FE(0);
566 
567 INSTANTIATE_FE(1);
568 
569 INSTANTIATE_FE(2);
570 
571 INSTANTIATE_FE(3);
572 
574 
575 } // namespace libMesh
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
const unsigned int invalid_uint
Definition: libmesh.h:245
INSTANTIATE_SUBDIVISION_FE
Definition: fe.C:573
INSTANTIATE_FE(3)
The base class for all geometric element types.
Definition: elem.h:100
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
FE(const FEType &fet)
Definition: fe.C:37
unsigned int p_level() const
Definition: elem.h:2555
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const =0
Template class which generates the different FE families and orders.
Definition: fe.h:89
const dof_id_type n_nodes
Definition: tecplot_io.C:68
virtual unsigned int n_nodes() const =0
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:233
virtual unsigned int n_edges() const =0
virtual unsigned int n_sides() const =0
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
Base class for all quadrature families and orders.
Definition: quadrature.h:62