libMesh
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/elem.h"
22 #include "libmesh/fe.h"
23 #include "libmesh/fe_interface.h"
24 #include "libmesh/fe_macro.h"
25 #include "libmesh/libmesh_logging.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 (auto n : elem->node_index_range())
195  cached_nodes[n] = elem->point(n);
196  }
197  }
198  else
199  {
200  // libmesh_assert_greater (elem->n_nodes(), 1);
201 
202  cached_nodes_still_fit = true;
203  if (cached_nodes.size() != elem->n_nodes())
204  cached_nodes_still_fit = false;
205  else
206  for (auto n : IntRange<unsigned int>(1, elem->n_nodes()))
207  {
208  if (!(elem->point(n) - elem->point(0)).relative_fuzzy_equals
209  ((cached_nodes[n] - cached_nodes[0]), 1e-13))
210  {
211  cached_nodes_still_fit = false;
212  break;
213  }
214  }
215 
216  if (this->shapes_need_reinit() && !cached_nodes_still_fit)
217  {
218  this->_fe_map->template init_reference_to_physical_map<Dim>
219  (this->qrule->get_points(), elem);
220  this->init_shape_functions (this->qrule->get_points(), elem);
221  cached_nodes.resize(elem->n_nodes());
222  for (auto n : elem->node_index_range())
223  cached_nodes[n] = elem->point(n);
224  }
225  }
226 
227  // The shape functions correspond to the qrule
228  this->shapes_on_quadrature = true;
229  }
230  }
231  else // With no defined elem, so mapping or caching to
232  // be done, and our "quadrature rule" is one point for nonlocal
233  // (SCALAR) variables and zero points for local variables.
234  {
235  this->elem_type = INVALID_ELEM;
236  this->_p_level = 0;
237 
238  if (!pts)
239  {
240  if (T == SCALAR)
241  {
242  this->qrule->get_points() =
243  std::vector<Point>(1,Point(0));
244 
245  this->qrule->get_weights() =
246  std::vector<Real>(1,1);
247  }
248  else
249  {
250  this->qrule->get_points().clear();
251  this->qrule->get_weights().clear();
252  }
253 
254  this->init_shape_functions (this->qrule->get_points(), elem);
255  }
256  else
257  this->init_shape_functions (*pts, elem);
258  }
259 
260  // Compute the map for this element.
261  if (pts != nullptr)
262  {
263  if (weights != nullptr)
264  {
265  this->_fe_map->compute_map (this->dim, *weights, elem, this->calculate_d2phi);
266  }
267  else
268  {
269  std::vector<Real> dummy_weights (pts->size(), 1.);
270  this->_fe_map->compute_map (this->dim, dummy_weights, elem, this->calculate_d2phi);
271  }
272  }
273  else
274  {
275  this->_fe_map->compute_map (this->dim, this->qrule->get_weights(), elem, this->calculate_d2phi);
276  }
277 
278  // Compute the shape functions and the derivatives at all of the
279  // quadrature points.
280  if (!cached_nodes_still_fit)
281  {
282  if (pts != nullptr)
283  this->compute_shape_functions (elem,*pts);
284  else
285  this->compute_shape_functions(elem,this->qrule->get_points());
286  }
287 }
288 
289 
290 
291 template <unsigned int Dim, FEFamily T>
292 void FE<Dim,T>::init_shape_functions(const std::vector<Point> & qp,
293  const Elem * elem)
294 {
295  // Start logging the shape function initialization
296  LOG_SCOPE("init_shape_functions()", "FE");
297 
298  // The number of quadrature points.
299  const unsigned int n_qp = cast_int<unsigned int>(qp.size());
300 
301  // Number of shape functions in the finite element approximation
302  // space.
303  const unsigned int n_approx_shape_functions =
304  this->n_shape_functions(this->get_type(),
305  this->get_order());
306 
307  // resize the vectors to hold current data
308  // Phi are the shape functions used for the FE approximation
309  // Phi_map are the shape functions used for the FE mapping
310  if (this->calculate_phi)
311  this->phi.resize (n_approx_shape_functions);
312 
313  if (this->calculate_dphi)
314  {
315  this->dphi.resize (n_approx_shape_functions);
316  this->dphidx.resize (n_approx_shape_functions);
317  this->dphidy.resize (n_approx_shape_functions);
318  this->dphidz.resize (n_approx_shape_functions);
319  }
320 
321  if (this->calculate_dphiref)
322  {
323  if (Dim > 0)
324  this->dphidxi.resize (n_approx_shape_functions);
325 
326  if (Dim > 1)
327  this->dphideta.resize (n_approx_shape_functions);
328 
329  if (Dim > 2)
330  this->dphidzeta.resize (n_approx_shape_functions);
331  }
332 
333  if (this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
334  this->curl_phi.resize(n_approx_shape_functions);
335 
336  if (this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
337  this->div_phi.resize(n_approx_shape_functions);
338 
339 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
340  if (this->calculate_d2phi)
341  {
342  this->d2phi.resize (n_approx_shape_functions);
343  this->d2phidx2.resize (n_approx_shape_functions);
344  this->d2phidxdy.resize (n_approx_shape_functions);
345  this->d2phidxdz.resize (n_approx_shape_functions);
346  this->d2phidy2.resize (n_approx_shape_functions);
347  this->d2phidydz.resize (n_approx_shape_functions);
348  this->d2phidz2.resize (n_approx_shape_functions);
349 
350  if (Dim > 0)
351  this->d2phidxi2.resize (n_approx_shape_functions);
352 
353  if (Dim > 1)
354  {
355  this->d2phidxideta.resize (n_approx_shape_functions);
356  this->d2phideta2.resize (n_approx_shape_functions);
357  }
358  if (Dim > 2)
359  {
360  this->d2phidxidzeta.resize (n_approx_shape_functions);
361  this->d2phidetadzeta.resize (n_approx_shape_functions);
362  this->d2phidzeta2.resize (n_approx_shape_functions);
363  }
364  }
365 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
366 
367  for (unsigned int i=0; i<n_approx_shape_functions; i++)
368  {
369  if (this->calculate_phi)
370  this->phi[i].resize (n_qp);
371  if (this->calculate_dphi)
372  {
373  this->dphi[i].resize (n_qp);
374  this->dphidx[i].resize (n_qp);
375  this->dphidy[i].resize (n_qp);
376  this->dphidz[i].resize (n_qp);
377  }
378 
379  if (this->calculate_dphiref)
380  {
381  if (Dim > 0)
382  this->dphidxi[i].resize(n_qp);
383 
384  if (Dim > 1)
385  this->dphideta[i].resize(n_qp);
386 
387  if (Dim > 2)
388  this->dphidzeta[i].resize(n_qp);
389  }
390 
391  if (this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
392  this->curl_phi[i].resize(n_qp);
393 
394  if (this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR))
395  this->div_phi[i].resize(n_qp);
396 
397 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
398  if (this->calculate_d2phi)
399  {
400  this->d2phi[i].resize (n_qp);
401  this->d2phidx2[i].resize (n_qp);
402  this->d2phidxdy[i].resize (n_qp);
403  this->d2phidxdz[i].resize (n_qp);
404  this->d2phidy2[i].resize (n_qp);
405  this->d2phidydz[i].resize (n_qp);
406  this->d2phidz2[i].resize (n_qp);
407  if (Dim > 0)
408  this->d2phidxi2[i].resize (n_qp);
409  if (Dim > 1)
410  {
411  this->d2phidxideta[i].resize (n_qp);
412  this->d2phideta2[i].resize (n_qp);
413  }
414  if (Dim > 2)
415  {
416  this->d2phidxidzeta[i].resize (n_qp);
417  this->d2phidetadzeta[i].resize (n_qp);
418  this->d2phidzeta2[i].resize (n_qp);
419  }
420  }
421 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
422  }
423 
424 
425 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
426  //------------------------------------------------------------
427  // Initialize the data fields, which should only be used for infinite
428  // elements, to some sensible values, so that using a FE with the
429  // variational formulation of an InfFE, correct element matrices are
430  // returned
431 
432  {
433  this->weight.resize (n_qp);
434  this->dweight.resize (n_qp);
435  this->dphase.resize (n_qp);
436 
437  for (unsigned int p=0; p<n_qp; p++)
438  {
439  this->weight[p] = 1.;
440  this->dweight[p].zero();
441  this->dphase[p].zero();
442  }
443 
444  }
445 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
446 
447  switch (Dim)
448  {
449 
450  //------------------------------------------------------------
451  // 0D
452  case 0:
453  {
454  break;
455  }
456 
457  //------------------------------------------------------------
458  // 1D
459  case 1:
460  {
461  // Compute the value of the approximation shape function i at quadrature point p
462  if (this->calculate_dphiref)
463  for (unsigned int i=0; i<n_approx_shape_functions; i++)
464  for (unsigned int p=0; p<n_qp; p++)
465  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
466 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
467  if (this->calculate_d2phi)
468  for (unsigned int i=0; i<n_approx_shape_functions; i++)
469  for (unsigned int p=0; p<n_qp; p++)
470  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
471 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
472 
473  break;
474  }
475 
476 
477 
478  //------------------------------------------------------------
479  // 2D
480  case 2:
481  {
482  // Compute the value of the approximation shape function i at quadrature point p
483  if (this->calculate_dphiref)
484  for (unsigned int i=0; i<n_approx_shape_functions; i++)
485  for (unsigned int p=0; p<n_qp; p++)
486  {
487  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
488  this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]);
489  }
490 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
491  if (this->calculate_d2phi)
492  for (unsigned int i=0; i<n_approx_shape_functions; i++)
493  for (unsigned int p=0; p<n_qp; p++)
494  {
495  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
496  this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]);
497  this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]);
498  }
499 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
500 
501 
502  break;
503  }
504 
505 
506 
507  //------------------------------------------------------------
508  // 3D
509  case 3:
510  {
511  // Compute the value of the approximation shape function i at quadrature point p
512  if (this->calculate_dphiref)
513  for (unsigned int i=0; i<n_approx_shape_functions; i++)
514  for (unsigned int p=0; p<n_qp; p++)
515  {
516  this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]);
517  this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]);
518  this->dphidzeta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 2, qp[p]);
519  }
520 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
521  if (this->calculate_d2phi)
522  for (unsigned int i=0; i<n_approx_shape_functions; i++)
523  for (unsigned int p=0; p<n_qp; p++)
524  {
525  this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]);
526  this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]);
527  this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]);
528  this->d2phidxidzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 3, qp[p]);
529  this->d2phidetadzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 4, qp[p]);
530  this->d2phidzeta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 5, qp[p]);
531  }
532 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
533 
534  break;
535  }
536 
537 
538  default:
539  libmesh_error_msg("Invalid dimension Dim = " << Dim);
540  }
541 }
542 
543 
544 
545 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
546 
547 template <unsigned int Dim, FEFamily T>
548 void FE<Dim,T>::init_base_shape_functions(const std::vector<Point> & qp,
549  const Elem * e)
550 {
551  this->elem_type = e->type();
552  this->_fe_map->template init_reference_to_physical_map<Dim>(qp, e);
553  init_shape_functions(qp, e);
554 }
555 
556 #endif // LIBMESH_ENABLE_INFINITE_ELEMENTS
557 
558 
559 
560 //--------------------------------------------------------------
561 // Explicit instantiations using macro from fe_macro.h
562 
563 INSTANTIATE_FE(0);
564 
565 INSTANTIATE_FE(1);
566 
567 INSTANTIATE_FE(2);
568 
569 INSTANTIATE_FE(3);
570 
572 
573 } // namespace libMesh
libMesh::QBase
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Definition: quadrature.h:61
libMesh::Elem::is_node_on_side
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
libMesh::Elem::n_edges
virtual unsigned int n_edges() const =0
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::INSTANTIATE_FE
INSTANTIATE_FE(3)
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::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::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::node_index_range
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2170
libMesh::INSTANTIATE_SUBDIVISION_FE
INSTANTIATE_SUBDIVISION_FE
Definition: fe.C:571
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
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
libMesh::FE
A specific instantiation of the FEBase class.
Definition: fe.h:95
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FE::FE
FE(const FEType &fet)
Constructor.
Definition: fe.C:37
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Elem::is_node_on_edge
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const =0
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::TYPE_VECTOR
Definition: enum_fe_family.h:94
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEOutputType
Most finite element types in libMesh are scalar-valued.
Definition: fe.h:49
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Elem::n_sides
virtual unsigned int n_sides() const =0
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