Reference documentation for deal.II version 9.0.0
fe_abf.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/utilities.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/quadrature_lib.h>
20 #include <deal.II/base/qprojector.h>
21 #include <deal.II/base/table.h>
22 #include <deal.II/grid/tria.h>
23 #include <deal.II/grid/tria_iterator.h>
24 #include <deal.II/dofs/dof_accessor.h>
25 #include <deal.II/fe/fe.h>
26 #include <deal.II/fe/mapping.h>
27 #include <deal.II/fe/fe_abf.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/fe/fe_tools.h>
30 
31 #include <sstream>
32 #include <iostream>
33 #include <deal.II/base/std_cxx14/memory.h>
34 
35 
36 //TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
37 //adjust_line_dof_index_for_line_orientation_table fields, and write tests
38 //similar to bits/face_orientation_and_fe_q_*
39 
40 
41 DEAL_II_NAMESPACE_OPEN
42 
43 
44 template <int dim>
45 FE_ABF<dim>::FE_ABF (const unsigned int deg)
46  :
47  FE_PolyTensor<PolynomialsABF<dim>, dim> (
48  deg,
49  FiniteElementData<dim>(get_dpo_vector(deg),
50  dim,
51  deg+2,
52  FiniteElementData<dim>::Hdiv),
53  std::vector<bool>(PolynomialsABF<dim>::compute_n_pols(deg), true),
54  std::vector<ComponentMask>(PolynomialsABF<dim>::compute_n_pols(deg),
55  std::vector<bool>(dim,true))),
56  rt_order(deg)
57 {
58  Assert (dim >= 2, ExcImpossibleInDim(dim));
59  const unsigned int n_dofs = this->dofs_per_cell;
60 
62  // First, initialize the
63  // generalized support points and
64  // quadrature weights, since they
65  // are required for interpolation.
67 
68  // Now compute the inverse node matrix, generating the correct
69  // basis functions from the raw ones. For a discussion of what
70  // exactly happens here, see FETools::compute_node_matrix.
72  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
73  this->inverse_node_matrix.invert(M);
74  // From now on, the shape functions provided by FiniteElement::shape_value
75  // and similar functions will be the correct ones, not
76  // the raw shape functions from the polynomial space anymore.
77 
78  // Reinit the vectors of
79  // restriction and prolongation
80  // matrices to the right sizes.
81  // Restriction only for isotropic
82  // refinement
84  // Fill prolongation matrices with embedding operators
85  FETools::compute_embedding_matrices (*this, this->prolongation, false, 1.e-10);
86 
88 
89  // TODO[TL]: for anisotropic refinement we will probably need a table of submatrices with an array for each refine case
90  std::vector<FullMatrix<double> >
91  face_embeddings(1<<(dim-1), FullMatrix<double>(this->dofs_per_face,
92  this->dofs_per_face));
93  // TODO: Something goes wrong there. The error of the least squares fit
94  // is to large ...
95  // FETools::compute_face_embedding_matrices(*this, face_embeddings.data(), 0, 0);
96  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
97  this->dofs_per_face);
98  unsigned int target_row=0;
99  for (unsigned int d=0; d<face_embeddings.size(); ++d)
100  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
101  {
102  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
103  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
104  ++target_row;
105  }
106 }
107 
108 
109 
110 template <int dim>
111 std::string
113 {
114  // note that the
115  // FETools::get_fe_by_name
116  // function depends on the
117  // particular format of the string
118  // this function returns, so they
119  // have to be kept in synch
120 
121  std::ostringstream namebuf;
122 
123  namebuf << "FE_ABF<" << dim << ">(" << rt_order << ")";
124 
125  return namebuf.str();
126 }
127 
128 
129 
130 template <int dim>
131 std::unique_ptr<FiniteElement<dim,dim> >
133 {
134  return std_cxx14::make_unique<FE_ABF<dim>>(rt_order);
135 }
136 
137 
138 //---------------------------------------------------------------------------
139 // Auxiliary and internal functions
140 //---------------------------------------------------------------------------
141 
142 
143 
144 // Version for 2d and higher. See above for 1d version
145 template <int dim>
146 void
148 {
149  QGauss<dim> cell_quadrature(deg+2);
150  const unsigned int n_interior_points = cell_quadrature.size();
151 
152  unsigned int n_face_points = (dim>1) ? 1 : 0;
153  // compute (deg+1)^(dim-1)
154  for (unsigned int d=1; d<dim; ++d)
155  n_face_points *= deg+1;
156 
157  this->generalized_support_points.resize (GeometryInfo<dim>::faces_per_cell*n_face_points
158  + n_interior_points);
159  this->generalized_face_support_points.resize (n_face_points);
160 
161 
162  // These might be required when the faces contribution is computed
163  // Therefore they will be initialized at this point.
164  std::vector<AnisotropicPolynomials<dim>* > polynomials_abf(dim);
165 
166  // Generate x_1^{i} x_2^{r+1} ...
167  for (unsigned int dd=0; dd<dim; ++dd)
168  {
169  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
170  for (unsigned int d=0; d<dim; ++d)
171  poly[d].push_back (Polynomials::Monomial<double> (deg+1));
173 
174  polynomials_abf[dd] = new AnisotropicPolynomials<dim>(poly);
175  }
176 
177  // Number of the point being entered
178  unsigned int current = 0;
179 
180  if (dim>1)
181  {
182  QGauss<dim-1> face_points (deg+1);
183  TensorProductPolynomials<dim-1> legendre =
185 
186  boundary_weights.reinit(n_face_points, legendre.n());
187 
188 // Assert (face_points.size() == this->dofs_per_face,
189 // ExcInternalError());
190 
191  for (unsigned int k=0; k<n_face_points; ++k)
192  {
193  this->generalized_face_support_points[k] = face_points.point(k);
194  // Compute its quadrature
195  // contribution for each
196  // moment.
197  for (unsigned int i=0; i<legendre.n(); ++i)
198  {
199  boundary_weights(k, i)
200  = face_points.weight(k)
201  * legendre.compute_value(i, face_points.point(k));
202  }
203  }
204 
206  for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
207  ++current)
208  {
209  // Enter the support point
210  // into the vector
211  this->generalized_support_points[current] = faces.point(current);
212  }
213 
214 
215  // Now initialize edge interior weights for the ABF elements.
216  // These are completely independent from the usual edge moments. They
217  // stem from applying the Gauss theorem to the nodal values, which
218  // was necessary to cast the ABF elements into the deal.II framework
219  // for vector valued elements.
220  boundary_weights_abf.reinit(faces.size(), polynomials_abf[0]->n() * dim);
221  for (unsigned int k=0; k < faces.size(); ++k)
222  {
223  for (unsigned int i=0; i<polynomials_abf[0]->n() * dim; ++i)
224  {
225  boundary_weights_abf(k,i) = polynomials_abf[i%dim]->
226  compute_value(i / dim, faces.point(k)) * faces.weight(k);
227  }
228  }
229  }
230 
231  // Create Legendre basis for the
232  // space D_xi Q_k
233  if (deg>0)
234  {
235  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
236 
237  for (unsigned int dd=0; dd<dim; ++dd)
238  {
239  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
240  for (unsigned int d=0; d<dim; ++d)
243 
244  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
245  }
246 
247  interior_weights.reinit(TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
248 
249  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
250  {
251  for (unsigned int i=0; i<polynomials[0]->n(); ++i)
252  for (unsigned int d=0; d<dim; ++d)
253  interior_weights(k,i,d) = cell_quadrature.weight(k)
254  * polynomials[d]->compute_value(i,cell_quadrature.point(k));
255  }
256 
257  for (unsigned int d=0; d<dim; ++d)
258  delete polynomials[d];
259  }
260 
261 
262  // Decouple the creation of the generalized support points
263  // from computation of interior weights.
264  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
265  this->generalized_support_points[current++] = cell_quadrature.point(k);
266 
267  // Additional functionality for the ABF elements
268  // TODO: Here the canonical extension of the principle
269  // behind the ABF elements is implemented. It is unclear,
270  // if this really leads to the ABF spaces in 3D!
271  interior_weights_abf.reinit(TableIndices<3>(cell_quadrature.size(),
272  polynomials_abf[0]->n() * dim, dim));
273  Tensor<1, dim> poly_grad;
274 
275  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
276  {
277  for (unsigned int i=0; i<polynomials_abf[0]->n() * dim; ++i)
278  {
279  poly_grad = polynomials_abf[i%dim]->compute_grad(i / dim,cell_quadrature.point(k))
280  * cell_quadrature.weight(k);
281  // The minus sign comes from the use of the Gauss theorem to replace the divergence.
282  for (unsigned int d=0; d<dim; ++d)
283  interior_weights_abf(k,i,d) = -poly_grad[d];
284  }
285  }
286 
287  for (unsigned int d=0; d<dim; ++d)
288  delete polynomials_abf[d];
289 
290  Assert (current == this->generalized_support_points.size(),
291  ExcInternalError());
292 }
293 
294 
295 
296 // This function is the same Raviart-Thomas interpolation performed by
297 // interpolate. Still, we cannot use interpolate, since it was written
298 // for smooth functions. The functions interpolated here are not
299 // smooth, maybe even not continuous. Therefore, we must double the
300 // number of quadrature points in each direction in order to integrate
301 // only smooth functions.
302 
303 // Then again, the interpolated function is chosen such that the
304 // moments coincide with the function to be interpolated.
305 
306 template <int dim>
307 void
309 {
310  if (dim==1)
311  {
313  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_cell; ++i)
314  this->restriction[iso][i].reinit(0,0);
315  return;
316  }
318  QGauss<dim-1> q_base (rt_order+1);
319  const unsigned int n_face_points = q_base.size();
320  // First, compute interpolation on
321  // subfaces
322  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
323  {
324  // The shape functions of the
325  // child cell are evaluated
326  // in the quadrature points
327  // of a full face.
328  Quadrature<dim> q_face
329  = QProjector<dim>::project_to_face(q_base, face);
330  // Store shape values, since the
331  // evaluation suffers if not
332  // ordered by point
333  Table<2,double> cached_values_face(this->dofs_per_cell, q_face.size());
334  for (unsigned int k=0; k<q_face.size(); ++k)
335  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
336  cached_values_face(i,k)
337  = this->shape_value_component(i, q_face.point(k),
339 
340  for (unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
341  {
342  // The weight functions for
343  // the coarse face are
344  // evaluated on the subface
345  // only.
346  Quadrature<dim> q_sub
347  = QProjector<dim>::project_to_subface(q_base, face, sub);
348  const unsigned int child
351 
352  // On a certain face, we must
353  // compute the moments of ALL
354  // fine level functions with
355  // the coarse level weight
356  // functions belonging to
357  // that face. Due to the
358  // orthogonalization process
359  // when building the shape
360  // functions, these weights
361  // are equal to the
362  // corresponding shape
363  // functions.
364  for (unsigned int k=0; k<n_face_points; ++k)
365  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
366  for (unsigned int i_face = 0; i_face < this->dofs_per_face; ++i_face)
367  {
368  // The quadrature
369  // weights on the
370  // subcell are NOT
371  // transformed, so we
372  // have to do it here.
373  this->restriction[iso][child](face*this->dofs_per_face+i_face,
374  i_child)
375  += Utilities::fixed_power<dim-1>(.5) * q_sub.weight(k)
376  * cached_values_face(i_child, k)
377  * this->shape_value_component(face*this->dofs_per_face+i_face,
378  q_sub.point(k),
380  }
381  }
382  }
383 
384  if (rt_order==0) return;
385 
386  // Create Legendre basis for the
387  // space D_xi Q_k. Here, we cannot
388  // use the shape functions
389  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
390  for (unsigned int dd=0; dd<dim; ++dd)
391  {
392  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
393  for (unsigned int d=0; d<dim; ++d)
395  poly[dd] = Polynomials::Legendre::generate_complete_basis(rt_order-1);
396 
397  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
398  }
399 
400  QGauss<dim> q_cell(rt_order+1);
401  const unsigned int start_cell_dofs
402  = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
403 
404  // Store shape values, since the
405  // evaluation suffers if not
406  // ordered by point
407  Table<3,double> cached_values_cell(this->dofs_per_cell, q_cell.size(), dim);
408  for (unsigned int k=0; k<q_cell.size(); ++k)
409  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
410  for (unsigned int d=0; d<dim; ++d)
411  cached_values_cell(i,k,d) = this->shape_value_component(i, q_cell.point(k), d);
412 
413  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
414  {
415  Quadrature<dim> q_sub = QProjector<dim>::project_to_child(q_cell, child);
416 
417  for (unsigned int k=0; k<q_sub.size(); ++k)
418  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
419  for (unsigned int d=0; d<dim; ++d)
420  for (unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
421  {
422  this->restriction[iso][child](start_cell_dofs+i_weight*dim+d,
423  i_child)
424  += q_sub.weight(k)
425  * cached_values_cell(i_child, k, d)
426  * polynomials[d]->compute_value(i_weight, q_sub.point(k));
427  }
428  }
429 
430  for (unsigned int d=0; d<dim; ++d)
431  delete polynomials[d];
432 }
433 
434 
435 
436 template <int dim>
437 std::vector<unsigned int>
438 FE_ABF<dim>::get_dpo_vector (const unsigned int rt_order)
439 {
440  if (dim == 1)
441  {
442  Assert (false, ExcImpossibleInDim(1));
443  return std::vector<unsigned int>();
444  }
445 
446  // the element is face-based (not
447  // to be confused with George
448  // W. Bush's Faith Based
449  // Initiative...), and we have
450  // (rt_order+1)^(dim-1) DoFs per face
451  unsigned int dofs_per_face = 1;
452  for (int d=0; d<dim-1; ++d)
453  dofs_per_face *= rt_order+1;
454 
455  // and then there are interior dofs
456  const unsigned int
457  interior_dofs = dim*(rt_order+1)*dofs_per_face;
458 
459  std::vector<unsigned int> dpo(dim+1);
460  dpo[dim-1] = dofs_per_face;
461  dpo[dim] = interior_dofs;
462 
463  return dpo;
464 }
465 
466 //---------------------------------------------------------------------------
467 // Data field initialization
468 //---------------------------------------------------------------------------
469 
470 template <int dim>
471 bool
472 FE_ABF<dim>::has_support_on_face (const unsigned int shape_index,
473  const unsigned int face_index) const
474 {
475  Assert (shape_index < this->dofs_per_cell,
476  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
479 
480  // Return computed values if we
481  // know them easily. Otherwise, it
482  // is always safe to return true.
483  switch (rt_order)
484  {
485  case 0:
486  {
487  switch (dim)
488  {
489  case 2:
490  {
491  // only on the one
492  // non-adjacent face
493  // are the values
494  // actually zero. list
495  // these in a table
496  return (face_index != GeometryInfo<dim>::opposite_face[shape_index]);
497  }
498 
499  default:
500  return true;
501  };
502  };
503 
504  default: // other rt_order
505  return true;
506  };
507 
508  return true;
509 }
510 
511 
512 
513 template <int dim>
514 void
517  std::vector<double> &nodal_values) const
518 {
519  Assert (support_point_values.size() == this->generalized_support_points.size(),
520  ExcDimensionMismatch(support_point_values.size(), this->generalized_support_points.size()));
521  Assert (support_point_values[0].size() == this->n_components(),
522  ExcDimensionMismatch(support_point_values[0].size(), this->n_components()));
523  Assert (nodal_values.size() == this->dofs_per_cell,
524  ExcDimensionMismatch(nodal_values.size(),this->dofs_per_cell));
525 
526  std::fill(nodal_values.begin(), nodal_values.end(), 0.);
527 
528  const unsigned int n_face_points = boundary_weights.size(0);
529  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
530  for (unsigned int k=0; k<n_face_points; ++k)
531  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
532  {
533  nodal_values[i+face*this->dofs_per_face] += boundary_weights(k,i)
534  * support_point_values[face*n_face_points+k][GeometryInfo<dim>::unit_normal_direction[face]];
535  }
536 
537  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
538  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
539 
540  for (unsigned int k=0; k<interior_weights.size(0); ++k)
541  for (unsigned int i=0; i<interior_weights.size(1); ++i)
542  for (unsigned int d=0; d<dim; ++d)
543  nodal_values[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * support_point_values[k+start_cell_points][d];
544 
545  const unsigned int start_abf_dofs = start_cell_dofs + interior_weights.size(1) * dim;
546 
547  // Cell integral of ABF terms
548  for (unsigned int k=0; k<interior_weights_abf.size(0); ++k)
549  for (unsigned int i=0; i<interior_weights_abf.size(1); ++i)
550  for (unsigned int d=0; d<dim; ++d)
551  nodal_values[start_abf_dofs+i] += interior_weights_abf(k,i,d) * support_point_values[k+start_cell_points][d];
552 
553  // Face integral of ABF terms
554  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
555  {
556  double n_orient = (double) GeometryInfo<dim>::unit_normal_orientation[face];
557  for (unsigned int fp=0; fp < n_face_points; ++fp)
558  {
559  // TODO: Check what the face_orientation, face_flip and face_rotation have to be in 3D
560  unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
561  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
562  nodal_values[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
563  * support_point_values[face*n_face_points+fp][GeometryInfo<dim>::unit_normal_direction[face]];
564  }
565  }
566 
567  // TODO: Check if this "correction" can be removed.
568  for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
569  if (std::fabs (nodal_values[start_abf_dofs+i]) < 1.0e-16)
570  nodal_values[start_abf_dofs+i] = 0.0;
571 }
572 
573 
574 
575 template <int dim>
576 std::size_t
578 {
579  Assert (false, ExcNotImplemented ());
580  return 0;
581 }
582 
583 
584 
585 /*-------------- Explicit Instantiations -------------------------------*/
586 #include "fe_abf.inst"
587 
588 DEAL_II_NAMESPACE_CLOSE
void initialize_restriction()
Definition: fe_abf.cc:308
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1096
FullMatrix< double > interface_constraints
Definition: fe.h:2401
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:874
virtual std::string get_name() const
Definition: fe_abf.cc:112
const Point< dim > & point(const unsigned int i) const
virtual std::size_t memory_consumption() const
Definition: fe_abf.cc:577
friend class FE_ABF
Definition: fe_abf.h:244
STL namespace.
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
FullMatrix< double > compute_node_matrix(const FiniteElement< dim, spacedim > &fe)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_abf.cc:516
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_abf.cc:472
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
void invert(const FullMatrix< number2 > &M)
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:692
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const unsigned int dofs_per_face
Definition: fe_base.h:288
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_abf.cc:438
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize_support_points(const unsigned int rt_degree)
Definition: fe_abf.cc:147
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
Definition: fe_abf.cc:132
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1174