Reference documentation for deal.II version 8.5.1
fe_q_base.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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/quadrature.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/template_constraints.h>
20 #include <deal.II/base/tensor_product_polynomials.h>
21 #include <deal.II/base/tensor_product_polynomials_const.h>
22 #include <deal.II/base/tensor_product_polynomials_bubbles.h>
23 #include <deal.II/base/polynomials_piecewise.h>
24 #include <deal.II/fe/fe_q_base.h>
25 #include <deal.II/fe/fe_dgq.h>
26 #include <deal.II/fe/fe_dgp.h>
27 #include <deal.II/fe/fe_nothing.h>
28 #include <deal.II/fe/fe_tools.h>
29 #include <deal.II/base/quadrature_lib.h>
30 
31 #include <vector>
32 #include <sstream>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 
37 namespace FE_Q_Helper
38 {
39  namespace
40  {
41  // get the renumbering for faces
42  template <int dim>
43  inline
44  std::vector<unsigned int>
45  face_lexicographic_to_hierarchic_numbering (const unsigned int degree)
46  {
47  std::vector<unsigned int> dpo(dim, 1U);
48  for (unsigned int i=1; i<dpo.size(); ++i)
49  dpo[i]=dpo[i-1]*(degree-1);
50  const ::FiniteElementData<dim-1> face_data(dpo,1,degree);
51  std::vector<unsigned int> face_renumber (face_data.dofs_per_cell);
52  FETools::lexicographic_to_hierarchic_numbering (face_data, face_renumber);
53  return face_renumber;
54  }
55 
56  // dummy specialization for dim == 1 to avoid linker errors
57  template <>
58  inline
59  std::vector<unsigned int>
60  face_lexicographic_to_hierarchic_numbering<1> (const unsigned int)
61  {
62  return std::vector<unsigned int>();
63  }
64 
65 
66 
67  // in get_restriction_matrix() and get_prolongation_matrix(), want to undo
68  // tensorization on inner loops for performance reasons. this clears a
69  // dim-array
70  template <int dim>
71  inline
72  void
73  zero_indices (unsigned int (&indices)[dim])
74  {
75  for (unsigned int d=0; d<dim; ++d)
76  indices[d] = 0;
77  }
78 
79 
80 
81  // in get_restriction_matrix() and get_prolongation_matrix(), want to undo
82  // tensorization on inner loops for performance reasons. this increments
83  // tensor product indices
84  template <int dim>
85  inline
86  void
87  increment_indices (unsigned int (&indices)[dim],
88  const unsigned int dofs1d)
89  {
90  ++indices[0];
91  for (int d=0; d<dim-1; ++d)
92  if (indices[d]==dofs1d)
93  {
94  indices[d] = 0;
95  indices[d+1]++;
96  }
97  }
98  }
99 }
100 
101 
102 
107 template <class PolynomialType, int xdim, int xspacedim>
108 struct FE_Q_Base<PolynomialType,xdim,xspacedim>::Implementation
109 {
114  template <int spacedim>
115  static
116  void initialize_constraints (const std::vector<Point<1> > &,
118  {
119  // no constraints in 1d
120  }
121 
122 
123  template <int spacedim>
124  static
125  void initialize_constraints (const std::vector<Point<1> > &/*points*/,
127  {
128  const unsigned int dim = 2;
129 
130  unsigned int q_deg = fe.degree;
131  if (types_are_equal<PolynomialType, TensorProductPolynomialsBubbles<dim> >::value)
132  q_deg = fe.degree-1;
133 
134  // restricted to each face, the traces of the shape functions is an
135  // element of P_{k} (in 2d), or Q_{k} (in 3d), where k is the degree of
136  // the element. from this, we interpolate between mother and cell face.
137 
138  // the interpolation process works as follows: on each subface, we want
139  // that finite element solutions from both sides coincide. i.e. if a and b
140  // are expansion coefficients for the shape functions from both sides, we
141  // seek a relation between a and b such that
142  // sum_j a_j phi^c_j(x) == sum_j b_j phi_j(x)
143  // for all points x on the interface. here, phi^c_j are the shape
144  // functions on the small cell on one side of the face, and phi_j those on
145  // the big cell on the other side. To get this relation, it suffices to
146  // look at a sufficient number of points for which this has to hold. if
147  // there are n functions, then we need n evaluation points, and we choose
148  // them equidistantly.
149  //
150  // we obtain the matrix system
151  // A a == B b
152  // where
153  // A_ij = phi^c_j(x_i)
154  // B_ij = phi_j(x_i)
155  // and the relation we are looking for is
156  // a = A^-1 B b
157  //
158  // for the special case of Lagrange interpolation polynomials, A_ij
159  // reduces to delta_ij, and
160  // a_i = B_ij b_j
161  // Hence, interface_constraints(i,j)=B_ij.
162  //
163  // for the general case, where we don't have Lagrange interpolation
164  // polynomials, this is a little more complicated. Then we would evaluate
165  // at a number of points and invert the interpolation matrix A.
166  //
167  // Note, that we build up these matrices for all subfaces at once, rather
168  // than considering them separately. the reason is that we finally will
169  // want to have them in this order anyway, as this is the format we need
170  // inside deal.II
171 
172  // In the following the points x_i are constructed in following order
173  // (n=degree-1)
174  // *----------*---------*
175  // 1..n 0 n+1..2n
176  // i.e. first the midpoint of the line, then the support points on subface
177  // 0 and on subface 1
178  std::vector<Point<dim-1> > constraint_points;
179  // Add midpoint
180  constraint_points.push_back (Point<dim-1> (0.5));
181 
182  if (q_deg>1)
183  {
184  const unsigned int n=q_deg-1;
185  const double step=1./q_deg;
186  // subface 0
187  for (unsigned int i=1; i<=n; ++i)
188  constraint_points.push_back (
190  // subface 1
191  for (unsigned int i=1; i<=n; ++i)
192  constraint_points.push_back (
194  }
195 
196  // Now construct relation between destination (child) and source (mother)
197  // dofs.
198 
200  .TableBase<2,double>::reinit (fe.interface_constraints_size());
201 
202  // use that the element evaluates to 1 at index 0 and along the line at
203  // zero
204  const std::vector<unsigned int> &index_map_inverse =
205  fe.poly_space.get_numbering_inverse();
206  const std::vector<unsigned int> face_index_map =
207  FE_Q_Helper::face_lexicographic_to_hierarchic_numbering<dim>(q_deg);
208  Assert(std::abs(fe.poly_space.compute_value(index_map_inverse[0],Point<dim>())
209  - 1.) < 1e-14,
210  ExcInternalError());
211 
212  for (unsigned int i=0; i<constraint_points.size(); ++i)
213  for (unsigned int j=0; j<q_deg+1; ++j)
214  {
215  Point<dim> p;
216  p[0] = constraint_points[i](0);
217  fe.interface_constraints(i,face_index_map[j]) =
218  fe.poly_space.compute_value(index_map_inverse[j], p);
219 
220  // if the value is small up to round-off, then simply set it to zero
221  // to avoid unwanted fill-in of the constraint matrices (which would
222  // then increase the number of other DoFs a constrained DoF would
223  // couple to)
224  if (std::fabs(fe.interface_constraints(i,face_index_map[j])) < 1e-13)
225  fe.interface_constraints(i,face_index_map[j]) = 0;
226  }
227  }
228 
229 
230  template <int spacedim>
231  static
232  void initialize_constraints (const std::vector<Point<1> > &/*points*/,
234  {
235  const unsigned int dim = 3;
236 
237  unsigned int q_deg = fe.degree;
238  if (types_are_equal<PolynomialType,TensorProductPolynomialsBubbles<dim> >::value)
239  q_deg = fe.degree-1;
240 
241  // For a detailed documentation of the interpolation see the
242  // FE_Q_Base<2>::initialize_constraints function.
243 
244  // In the following the points x_i are constructed in the order as
245  // described in the documentation of the FiniteElement class (fe_base.h),
246  // i.e.
247  // *--15--4--16--*
248  // | | |
249  // 10 19 6 20 12
250  // | | |
251  // 1--7---0--8---2
252  // | | |
253  // 9 17 5 18 11
254  // | | |
255  // *--13--3--14--*
256  std::vector<Point<dim-1> > constraint_points;
257 
258  // Add midpoint
259  constraint_points.push_back (Point<dim-1> (0.5, 0.5));
260 
261  // Add midpoints of lines of "mother-face"
262  constraint_points.push_back (Point<dim-1> (0, 0.5));
263  constraint_points.push_back (Point<dim-1> (1, 0.5));
264  constraint_points.push_back (Point<dim-1> (0.5, 0));
265  constraint_points.push_back (Point<dim-1> (0.5, 1));
266 
267  if (q_deg>1)
268  {
269  const unsigned int n=q_deg-1;
270  const double step=1./q_deg;
271  std::vector<Point<dim-2> > line_support_points(n);
272  for (unsigned int i=0; i<n; ++i)
273  line_support_points[i](0)=(i+1)*step;
274  Quadrature<dim-2> qline(line_support_points);
275 
276  // auxiliary points in 2d
277  std::vector<Point<dim-1> > p_line(n);
278 
279  // Add nodes of lines interior in the "mother-face"
280 
281  // line 5: use line 9
282  QProjector<dim-1>::project_to_subface(qline, 0, 0, p_line);
283  for (unsigned int i=0; i<n; ++i)
284  constraint_points.push_back (p_line[i] + Point<dim-1> (0.5, 0));
285  // line 6: use line 10
286  QProjector<dim-1>::project_to_subface(qline, 0, 1, p_line);
287  for (unsigned int i=0; i<n; ++i)
288  constraint_points.push_back (p_line[i] + Point<dim-1> (0.5, 0));
289  // line 7: use line 13
290  QProjector<dim-1>::project_to_subface(qline, 2, 0, p_line);
291  for (unsigned int i=0; i<n; ++i)
292  constraint_points.push_back (p_line[i] + Point<dim-1> (0, 0.5));
293  // line 8: use line 14
294  QProjector<dim-1>::project_to_subface(qline, 2, 1, p_line);
295  for (unsigned int i=0; i<n; ++i)
296  constraint_points.push_back (p_line[i] + Point<dim-1> (0, 0.5));
297 
298  // DoFs on bordering lines lines 9-16
299  for (unsigned int face=0; face<GeometryInfo<dim-1>::faces_per_cell; ++face)
300  for (unsigned int subface=0;
301  subface<GeometryInfo<dim-1>::max_children_per_face; ++subface)
302  {
303  QProjector<dim-1>::project_to_subface(qline, face, subface, p_line);
304  constraint_points.insert(constraint_points.end(),
305  p_line.begin(), p_line.end());
306  }
307 
308  // Create constraints for interior nodes
309  std::vector<Point<dim-1> > inner_points(n*n);
310  for (unsigned int i=0, iy=1; iy<=n; ++iy)
311  for (unsigned int ix=1; ix<=n; ++ix)
312  inner_points[i++] = Point<dim-1> (ix*step, iy*step);
313 
314  // at the moment do this for isotropic face refinement only
315  for (unsigned int child=0;
316  child<GeometryInfo<dim-1>::max_children_per_cell; ++child)
317  for (unsigned int i=0; i<inner_points.size(); ++i)
318  constraint_points.push_back (
319  GeometryInfo<dim-1>::child_to_cell_coordinates(inner_points[i], child));
320  }
321 
322  // Now construct relation between destination (child) and source (mother)
323  // dofs.
324  const unsigned int pnts=(q_deg+1)*(q_deg+1);
325 
326  // use that the element evaluates to 1 at index 0 and along the line at
327  // zero
328  const std::vector<unsigned int> &index_map_inverse =
329  fe.poly_space.get_numbering_inverse();
330  const std::vector<unsigned int> face_index_map =
331  FE_Q_Helper::face_lexicographic_to_hierarchic_numbering<dim>(q_deg);
332  Assert(std::abs(fe.poly_space.compute_value(index_map_inverse[0],Point<dim>())
333  - 1.) < 1e-14,
334  ExcInternalError());
335 
337  .TableBase<2,double>::reinit (fe.interface_constraints_size());
338 
339  for (unsigned int i=0; i<constraint_points.size(); ++i)
340  {
341  const double interval = (double) (q_deg * 2);
342  bool mirror[dim - 1];
343  Point<dim> constraint_point;
344 
345  // Eliminate FP errors in constraint points. Due to their origin, they
346  // must all be fractions of the unit interval. If we have polynomial
347  // degree 4, the refined element has 8 intervals. Hence the
348  // coordinates must be 0, 0.125, 0.25, 0.375 etc. Now the coordinates
349  // of the constraint points will be multiplied by the inverse of the
350  // interval size (in the example by 8). After that the coordinates
351  // must be integral numbers. Hence a normal truncation is performed
352  // and the coordinates will be scaled back. The equal treatment of all
353  // coordinates should eliminate any FP errors.
354  for (unsigned int k=0; k<dim-1; ++k)
355  {
356  const int coord_int =
357  static_cast<int> (constraint_points[i](k) * interval + 0.25);
358  constraint_point(k) = 1.*coord_int / interval;
359 
360  // The following lines of code should eliminate the problems with
361  // the Constraint-Matrix, which appeared for P>=4. The
362  // ConstraintMatrix class complained about different constraints
363  // for the same entry of the Constraint-Matrix. Actually this
364  // difference could be attributed to FP errors, as it was in the
365  // range of 1.0e-16. These errors originate in the loss of
366  // symmetry in the FP approximation of the shape-functions.
367  // Considering a 3rd order shape function in 1D, we have
368  // N0(x)=N3(1-x) and N1(x)=N2(1-x). For higher order polynomials
369  // the FP approximations of the shape functions do not satisfy
370  // these equations any more! Thus in the following code
371  // everything is computed in the interval x \in [0..0.5], which is
372  // sufficient to express all values that could come out from a
373  // computation of any shape function in the full interval
374  // [0..1]. If x > 0.5 the computation is done for 1-x with the
375  // shape function N_{p-n} instead of N_n. Hence symmetry is
376  // preserved and everything works fine...
377  //
378  // For a different explanation of the problem, see the discussion
379  // in the FiniteElement class for constraint matrices in 3d.
380  mirror[k] = (constraint_point(k) > 0.5);
381  if (mirror[k])
382  constraint_point(k) = 1.0 - constraint_point(k);
383  }
384 
385  for (unsigned int j=0; j<pnts; ++j)
386  {
387  unsigned int indices[2] = { j % (q_deg+1), j / (q_deg+1) };
388 
389  for (unsigned int k = 0; k<2; ++k)
390  if (mirror[k])
391  indices[k] = q_deg - indices[k];
392 
393  const unsigned int
394  new_index = indices[1] * (q_deg + 1) + indices[0];
395 
396  fe.interface_constraints(i,face_index_map[j]) =
397  fe.poly_space.compute_value (index_map_inverse[new_index],
398  constraint_point);
399 
400  // if the value is small up to round-off, then simply set it to
401  // zero to avoid unwanted fill-in of the constraint matrices
402  // (which would then increase the number of other DoFs a
403  // constrained DoF would couple to)
404  if (std::fabs(fe.interface_constraints(i,face_index_map[j])) < 1e-13)
405  fe.interface_constraints(i,face_index_map[j]) = 0;
406  }
407  }
408  }
409 };
410 
411 
412 
413 template <class PolynomialType, int dim, int spacedim>
415 (const PolynomialType &poly_space,
416  const FiniteElementData<dim> &fe_data,
417  const std::vector<bool> &restriction_is_additive_flags)
418  :
419  FE_Poly<PolynomialType,dim,spacedim>(poly_space, fe_data, restriction_is_additive_flags,
420  std::vector<ComponentMask>(1, std::vector<bool>(1,true))),
421  q_degree (types_are_equal<PolynomialType, TensorProductPolynomialsBubbles<dim> >::value
422  ?this->degree-1
423  :this->degree)
424 {}
425 
426 
427 
428 template <class PolynomialType, int dim, int spacedim>
429 void
431 {
432  Assert (points[0][0] == 0,
433  ExcMessage ("The first support point has to be zero."));
434  Assert (points.back()[0] == 1,
435  ExcMessage ("The last support point has to be one."));
436 
437  // distinguish q/q_dg0 case: need to be flexible enough to allow more
438  // degrees of freedom than there are FE_Q degrees of freedom for derived
439  // class FE_Q_DG0 that otherwise shares 95% of the code.
440  const unsigned int q_dofs_per_cell = Utilities::fixed_power<dim>(q_degree+1);
441  Assert(q_dofs_per_cell == this->dofs_per_cell ||
442  q_dofs_per_cell+1 == this->dofs_per_cell ||
443  q_dofs_per_cell+dim == this->dofs_per_cell , ExcInternalError());
444 
445  {
446  std::vector<unsigned int> renumber(q_dofs_per_cell);
447  const FiniteElementData<dim> fe(get_dpo_vector(q_degree),1,
448  q_degree);
450  for (unsigned int i= q_dofs_per_cell; i<this->dofs_per_cell; ++i)
451  renumber.push_back(i);
452  this->poly_space.set_numbering(renumber);
453  }
454 
455  // finally fill in support points on cell and face
456  initialize_unit_support_points (points);
457  initialize_unit_face_support_points (points);
458 
459  // reinit constraints
460  initialize_constraints (points);
461 
462  // do not initialize embedding and restriction here. these matrices are
463  // initialized on demand in get_restriction_matrix and
464  // get_prolongation_matrix
465 
466  this->initialize_quad_dof_index_permutation();
467 }
468 
469 
470 
471 template <class PolynomialType, int dim, int spacedim>
472 void
475  FullMatrix<double> &interpolation_matrix) const
476 {
477  // go through the list of elements we can interpolate from
478  if (const FE_Q_Base<PolynomialType,dim,spacedim> *source_fe
479  = dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim>*>(&x_source_fe))
480  {
481  // ok, source is a Q element, so we will be able to do the work
482  Assert (interpolation_matrix.m() == this->dofs_per_cell,
483  ExcDimensionMismatch (interpolation_matrix.m(),
484  this->dofs_per_cell));
485  Assert (interpolation_matrix.n() == x_source_fe.dofs_per_cell,
486  ExcDimensionMismatch (interpolation_matrix.m(),
487  x_source_fe.dofs_per_cell));
488 
489  // only evaluate Q dofs
490  const unsigned int q_dofs_per_cell = Utilities::fixed_power<dim>(q_degree+1);
491  const unsigned int source_q_dofs_per_cell = Utilities::fixed_power<dim>(source_fe->degree+1);
492 
493  // evaluation is simply done by evaluating the other FE's basis functions on
494  // the unit support points (FE_Q has the property that the cell
495  // interpolation matrix is a unit matrix, so no need to evaluate it and
496  // invert it)
497  for (unsigned int j=0; j<q_dofs_per_cell; ++j)
498  {
499  // read in a point on this cell and evaluate the shape functions there
500  const Point<dim> p = this->unit_support_points[j];
501 
502  // FE_Q element evaluates to 1 in unit support point and to zero in all
503  // other points by construction
504  Assert(std::abs(this->poly_space.compute_value (j, p)-1.)<1e-13,
505  ExcInternalError());
506 
507  for (unsigned int i=0; i<source_q_dofs_per_cell; ++i)
508  interpolation_matrix(j,i) = source_fe->poly_space.compute_value (i, p);
509  }
510 
511  // for FE_Q_DG0, add one last row of identity
512  if (q_dofs_per_cell < this->dofs_per_cell)
513  {
514  AssertDimension(source_q_dofs_per_cell+1, source_fe->dofs_per_cell);
515  for (unsigned int i=0; i<source_q_dofs_per_cell; ++i)
516  interpolation_matrix(q_dofs_per_cell, i) = 0.;
517  for (unsigned int j=0; j<q_dofs_per_cell; ++j)
518  interpolation_matrix(j, source_q_dofs_per_cell) = 0.;
519  interpolation_matrix(q_dofs_per_cell, source_q_dofs_per_cell) = 1.;
520  }
521 
522  // cut off very small values
523  const double eps = 2e-13*q_degree*dim;
524  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
525  for (unsigned int j=0; j<source_fe->dofs_per_cell; ++j)
526  if (std::fabs(interpolation_matrix(i,j)) < eps)
527  interpolation_matrix(i,j) = 0.;
528 
529  // make sure that the row sum of each of the matrices is 1 at this
530  // point. this must be so since the shape functions sum up to 1
531  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
532  {
533  double sum = 0.;
534  for (unsigned int j=0; j<source_fe->dofs_per_cell; ++j)
535  sum += interpolation_matrix(i,j);
536 
537  Assert (std::fabs(sum-1) < eps, ExcInternalError());
538  }
539  }
540  else if (dynamic_cast<const FE_Nothing<dim>*>(&x_source_fe))
541  {
542  // the element we want to interpolate from is an FE_Nothing. this
543  // element represents a function that is constant zero and has no
544  // degrees of freedom, so the interpolation is simply a multiplication
545  // with a n_dofs x 0 matrix. there is nothing to do here
546 
547  // we would like to verify that the number of rows and columns of
548  // the matrix equals this->dofs_per_cell and zero. unfortunately,
549  // whenever we do FullMatrix::reinit(m,0), it sets both rows and
550  // columns to zero, instead of m and zero. thus, only test the
551  // number of columns
552  Assert (interpolation_matrix.n() == x_source_fe.dofs_per_cell,
553  ExcDimensionMismatch (interpolation_matrix.m(),
554  x_source_fe.dofs_per_cell));
555 
556  }
557  else
558  AssertThrow (false,
560 
561 }
562 
563 
564 
565 template <class PolynomialType, int dim, int spacedim>
566 void
569  FullMatrix<double> &interpolation_matrix) const
570 {
571  Assert (dim > 1, ExcImpossibleInDim(1));
572  get_subface_interpolation_matrix (source_fe, numbers::invalid_unsigned_int,
573  interpolation_matrix);
574 }
575 
576 
577 
578 template <class PolynomialType, int dim, int spacedim>
579 void
582  const unsigned int subface,
583  FullMatrix<double> &interpolation_matrix) const
584 {
585  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
586  ExcDimensionMismatch (interpolation_matrix.m(),
587  x_source_fe.dofs_per_face));
588 
589  // see if source is a Q element
590  if (const FE_Q_Base<PolynomialType,dim,spacedim> *source_fe
591  = dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim> *>(&x_source_fe))
592  {
593  // have this test in here since a table of size 2x0 reports its size as
594  // 0x0
595  Assert (interpolation_matrix.n() == this->dofs_per_face,
596  ExcDimensionMismatch (interpolation_matrix.n(),
597  this->dofs_per_face));
598 
599  // Make sure that the element for which the DoFs should be constrained
600  // is the one with the higher polynomial degree. Actually the procedure
601  // will work also if this assertion is not satisfied. But the matrices
602  // produced in that case might lead to problems in the hp procedures,
603  // which use this method.
604  Assert (this->dofs_per_face <= source_fe->dofs_per_face,
605  (typename FiniteElement<dim,spacedim>::
606  ExcInterpolationNotImplemented ()));
607 
608  // generate a point on this cell and evaluate the shape functions there
609  const Quadrature<dim-1>
610  quad_face_support (source_fe->get_unit_face_support_points ());
611 
612  // Rule of thumb for FP accuracy, that can be expected for a given
613  // polynomial degree. This value is used to cut off values close to
614  // zero.
615  double eps = 2e-13*q_degree*(dim-1);
616 
617  // compute the interpolation matrix by simply taking the value at the
618  // support points.
619 //TODO: Verify that all faces are the same with respect to
620 // these support points. Furthermore, check if something has to
621 // be done for the face orientation flag in 3D.
622  const Quadrature<dim> subface_quadrature
623  = subface == numbers::invalid_unsigned_int
624  ?
625  QProjector<dim>::project_to_face (quad_face_support, 0)
626  :
627  QProjector<dim>::project_to_subface (quad_face_support, 0, subface);
628  for (unsigned int i=0; i<source_fe->dofs_per_face; ++i)
629  {
630  const Point<dim> &p = subface_quadrature.point (i);
631 
632  for (unsigned int j=0; j<this->dofs_per_face; ++j)
633  {
634  double matrix_entry = this->shape_value (this->face_to_cell_index(j, 0), p);
635 
636  // Correct the interpolated value. I.e. if it is close to 1 or
637  // 0, make it exactly 1 or 0. Unfortunately, this is required to
638  // avoid problems with higher order elements.
639  if (std::fabs (matrix_entry - 1.0) < eps)
640  matrix_entry = 1.0;
641  if (std::fabs (matrix_entry) < eps)
642  matrix_entry = 0.0;
643 
644  interpolation_matrix(i,j) = matrix_entry;
645  }
646  }
647 
648  // make sure that the row sum of each of the matrices is 1 at this
649  // point. this must be so since the shape functions sum up to 1
650  for (unsigned int j=0; j<source_fe->dofs_per_face; ++j)
651  {
652  double sum = 0.;
653 
654  for (unsigned int i=0; i<this->dofs_per_face; ++i)
655  sum += interpolation_matrix(j,i);
656 
657  Assert (std::fabs(sum-1) < eps, ExcInternalError());
658  }
659  }
660  else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe) != 0)
661  {
662  // nothing to do here, the FE_Nothing has no degrees of freedom anyway
663  }
664  else
665  AssertThrow (false,(typename FiniteElement<dim,spacedim>::
666  ExcInterpolationNotImplemented()));
667 }
668 
669 
670 
671 template <class PolynomialType, int dim, int spacedim>
672 bool
674 {
675  return true;
676 }
677 
678 
679 
680 
681 template <class PolynomialType, int dim, int spacedim>
682 std::vector<std::pair<unsigned int, unsigned int> >
685 {
686  // we can presently only compute these identities if both FEs are FE_Qs or
687  // if the other one is an FE_Nothing. in the first case, there should be
688  // exactly one single DoF of each FE at a vertex, and they should have
689  // identical value
690  if (dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim>*>(&fe_other) != 0)
691  {
692  return
693  std::vector<std::pair<unsigned int, unsigned int> >
694  (1, std::make_pair (0U, 0U));
695  }
696  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
697  {
698  // the FE_Nothing has no degrees of freedom, so there are no
699  // equivalencies to be recorded
700  return std::vector<std::pair<unsigned int, unsigned int> > ();
701  }
702  else if (fe_other.dofs_per_face == 0)
703  {
704  // if the other element has no elements on faces at all,
705  // then it would be impossible to enforce any kind of
706  // continuity even if we knew exactly what kind of element
707  // we have -- simply because the other element declares
708  // that it is discontinuous because it has no DoFs on
709  // its faces. in that case, just state that we have no
710  // constraints to declare
711  return std::vector<std::pair<unsigned int, unsigned int> > ();
712  }
713  else
714  {
715  Assert (false, ExcNotImplemented());
716  return std::vector<std::pair<unsigned int, unsigned int> > ();
717  }
718 }
719 
720 
721 
722 template <class PolynomialType, int dim, int spacedim>
723 std::vector<std::pair<unsigned int, unsigned int> >
726 {
727  // we can presently only compute these identities if both FEs are FE_Qs or
728  // if the other one is an FE_Nothing
729  if (const FE_Q_Base<PolynomialType,dim,spacedim> *fe_q_other
730  = dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim>*>(&fe_other))
731  {
732  // dofs are located along lines, so two dofs are identical if they are
733  // located at identical positions. if we had only equidistant points, we
734  // could simply check for similarity like (i+1)*q == (j+1)*p, but we
735  // might have other support points (e.g. Gauss-Lobatto
736  // points). Therefore, read the points in unit_support_points for the
737  // first coordinate direction. We take the lexicographic ordering of the
738  // points in the first direction (i.e., x-direction), which we access
739  // between index 1 and p-1 (index 0 and p are vertex dofs).
740  const unsigned int p = this->degree;
741  const unsigned int q = fe_q_other->degree;
742 
743  std::vector<std::pair<unsigned int, unsigned int> > identities;
744 
745  const std::vector<unsigned int> &index_map_inverse=
746  this->poly_space.get_numbering_inverse();
747  const std::vector<unsigned int> &index_map_inverse_other=
748  fe_q_other->poly_space.get_numbering_inverse();
749 
750  for (unsigned int i=0; i<p-1; ++i)
751  for (unsigned int j=0; j<q-1; ++j)
752  if (std::fabs(this->unit_support_points[index_map_inverse[i+1]][0]-
753  fe_q_other->unit_support_points[index_map_inverse_other[j+1]][0])
754  < 1e-14)
755  identities.push_back (std::make_pair(i,j));
756 
757  return identities;
758  }
759  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
760  {
761  // the FE_Nothing has no degrees of freedom, so there are no
762  // equivalencies to be recorded
763  return std::vector<std::pair<unsigned int, unsigned int> > ();
764  }
765  else if (fe_other.dofs_per_face == 0)
766  {
767  // if the other element has no elements on faces at all,
768  // then it would be impossible to enforce any kind of
769  // continuity even if we knew exactly what kind of element
770  // we have -- simply because the other element declares
771  // that it is discontinuous because it has no DoFs on
772  // its faces. in that case, just state that we have no
773  // constraints to declare
774  return std::vector<std::pair<unsigned int, unsigned int> > ();
775  }
776  else
777  {
778  Assert (false, ExcNotImplemented());
779  return std::vector<std::pair<unsigned int, unsigned int> > ();
780  }
781 }
782 
783 
784 
785 template <class PolynomialType, int dim, int spacedim>
786 std::vector<std::pair<unsigned int, unsigned int> >
789 {
790  // we can presently only compute these identities if both FEs are FE_Qs or
791  // if the other one is an FE_Nothing
792  if (const FE_Q_Base<PolynomialType,dim,spacedim> *fe_q_other
793  = dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim>*>(&fe_other))
794  {
795  // this works exactly like the line case above, except that now we have
796  // to have two indices i1, i2 and j1, j2 to characterize the dofs on the
797  // face of each of the finite elements. since they are ordered
798  // lexicographically along the first line and we have a tensor product,
799  // the rest is rather straightforward
800  const unsigned int p = this->degree;
801  const unsigned int q = fe_q_other->degree;
802 
803  std::vector<std::pair<unsigned int, unsigned int> > identities;
804 
805  const std::vector<unsigned int> &index_map_inverse=
806  this->poly_space.get_numbering_inverse();
807  const std::vector<unsigned int> &index_map_inverse_other=
808  fe_q_other->poly_space.get_numbering_inverse();
809 
810  for (unsigned int i1=0; i1<p-1; ++i1)
811  for (unsigned int i2=0; i2<p-1; ++i2)
812  for (unsigned int j1=0; j1<q-1; ++j1)
813  for (unsigned int j2=0; j2<q-1; ++j2)
814  if ((std::fabs(this->unit_support_points[index_map_inverse[i1+1]][0]-
815  fe_q_other->unit_support_points[index_map_inverse_other[j1+1]][0])
816  < 1e-14)
817  &&
818  (std::fabs(this->unit_support_points[index_map_inverse[i2+1]][0]-
819  fe_q_other->unit_support_points[index_map_inverse_other[j2+1]][0])
820  < 1e-14))
821  identities.push_back (std::make_pair(i1*(p-1)+i2,
822  j1*(q-1)+j2));
823 
824  return identities;
825  }
826  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
827  {
828  // the FE_Nothing has no degrees of freedom, so there are no
829  // equivalencies to be recorded
830  return std::vector<std::pair<unsigned int, unsigned int> > ();
831  }
832  else if (fe_other.dofs_per_face == 0)
833  {
834  // if the other element has no elements on faces at all,
835  // then it would be impossible to enforce any kind of
836  // continuity even if we knew exactly what kind of element
837  // we have -- simply because the other element declares
838  // that it is discontinuous because it has no DoFs on
839  // its faces. in that case, just state that we have no
840  // constraints to declare
841  return std::vector<std::pair<unsigned int, unsigned int> > ();
842  }
843  else
844  {
845  Assert (false, ExcNotImplemented());
846  return std::vector<std::pair<unsigned int, unsigned int> > ();
847  }
848 }
849 
850 
851 
852 template <class PolynomialType, int dim, int spacedim>
856 {
857  if (const FE_Q_Base<PolynomialType,dim,spacedim> *fe_q_other
858  = dynamic_cast<const FE_Q_Base<PolynomialType,dim,spacedim>*>(&fe_other))
859  {
860  if (this->degree < fe_q_other->degree)
862  else if (this->degree == fe_q_other->degree)
864  else
866  }
867  else if (const FE_Nothing<dim,spacedim> *fe_nothing
868  = dynamic_cast<const FE_Nothing<dim,spacedim>*>(&fe_other))
869  {
870  if (fe_nothing->is_dominating())
871  {
873  }
874  else
875  {
876  // the FE_Nothing has no degrees of freedom and it is typically used in
877  // a context where we don't require any continuity along the interface
879  }
880  }
881  else if ((dynamic_cast<const FE_DGQ<dim,spacedim>*>(&fe_other) != 0)
882  ||
883  (dynamic_cast<const FE_DGP<dim,spacedim>*>(&fe_other) != 0))
884  {
885  // there are no requirements between continuous and
886  // discontinuous elements
888  }
889 
890  Assert (false, ExcNotImplemented());
892 }
893 
894 
895 //---------------------------------------------------------------------------
896 // Auxiliary functions
897 //---------------------------------------------------------------------------
898 
899 
900 
901 template <class PolynomialType, int dim, int spacedim>
903 (const std::vector<Point<1> > &points)
904 {
905  const std::vector<unsigned int> &index_map_inverse=
906  this->poly_space.get_numbering_inverse();
907 
908  Quadrature<1> support_1d(points);
909  Quadrature<dim> support_quadrature(support_1d);
910  this->unit_support_points.resize(support_quadrature.size());
911 
912  for (unsigned int k=0; k<support_quadrature.size(); k++)
913  this->unit_support_points[index_map_inverse[k]] = support_quadrature.point(k);
914 }
915 
916 
917 
918 template <class PolynomialType, int dim, int spacedim>
920 (const std::vector<Point<1> > &points)
921 {
922  // no faces in 1d, so nothing to do
923  if (dim == 1)
924  return;
925 
926  const unsigned int codim = dim-1;
927  this->unit_face_support_points.resize(Utilities::fixed_power<codim>(q_degree+1));
928 
929  // find renumbering of faces and assign from values of quadrature
930  std::vector<unsigned int> face_index_map =
931  FE_Q_Helper::face_lexicographic_to_hierarchic_numbering<dim>(q_degree);
932  Quadrature<1> support_1d(points);
933  Quadrature<codim> support_quadrature(support_1d);
934  this->unit_face_support_points.resize(support_quadrature.size());
935 
936  for (unsigned int k=0; k<support_quadrature.size(); k++)
937  this->unit_face_support_points[face_index_map[k]] = support_quadrature.point(k);
938 }
939 
940 
941 
942 template <class PolynomialType, int dim, int spacedim>
943 void
945 {
946  // for 1D and 2D, do nothing
947  if (dim < 3)
948  return;
949 
950  Assert (this->adjust_quad_dof_index_for_face_orientation_table.n_elements()==8*this->dofs_per_quad,
951  ExcInternalError());
952 
953  const unsigned int n=q_degree-1;
954  Assert(n*n==this->dofs_per_quad, ExcInternalError());
955 
956  // alias for the table to fill
957  Table<2,int> &data=this->adjust_quad_dof_index_for_face_orientation_table;
958 
959  // the dofs on a face are connected to a n x n matrix. for example, for
960  // degree==4 we have the following dofs on a quad
961 
962  // ___________
963  // | |
964  // | 6 7 8 |
965  // | |
966  // | 3 4 5 |
967  // | |
968  // | 0 1 2 |
969  // |___________|
970  //
971  // we have dof_no=i+n*j with index i in x-direction and index j in
972  // y-direction running from 0 to n-1. to extract i and j we can use
973  // i=dof_no%n and j=dof_no/n. i and j can then be used to construct the
974  // rotated and mirrored numbers.
975 
976 
977  for (unsigned int local=0; local<this->dofs_per_quad; ++local)
978  // face support points are in lexicographic ordering with x running
979  // fastest. invert that (y running fastest)
980  {
981  unsigned int i=local%n,
982  j=local/n;
983 
984  // face_orientation=false, face_flip=false, face_rotation=false
985  data(local,0)=j + i *n - local;
986  // face_orientation=false, face_flip=false, face_rotation=true
987  data(local,1)=i + (n-1-j)*n - local;
988  // face_orientation=false, face_flip=true, face_rotation=false
989  data(local,2)=(n-1-j) + (n-1-i)*n - local;
990  // face_orientation=false, face_flip=true, face_rotation=true
991  data(local,3)=(n-1-i) + j *n - local;
992  // face_orientation=true, face_flip=false, face_rotation=false
993  data(local,4)=0;
994  // face_orientation=true, face_flip=false, face_rotation=true
995  data(local,5)=j + (n-1-i)*n - local;
996  // face_orientation=true, face_flip=true, face_rotation=false
997  data(local,6)=(n-1-i) + (n-1-j)*n - local;
998  // face_orientation=true, face_flip=true, face_rotation=true
999  data(local,7)=(n-1-j) + i *n - local;
1000  }
1001 
1002  // additionally initialize reordering of line dofs
1003  for (unsigned int i=0; i<this->dofs_per_line; ++i)
1004  this->adjust_line_dof_index_for_line_orientation_table[i]=this->dofs_per_line-1-i - i;
1005 }
1006 
1007 
1008 
1009 template <class PolynomialType, int dim, int spacedim>
1010 unsigned int
1012 face_to_cell_index (const unsigned int face_index,
1013  const unsigned int face,
1014  const bool face_orientation,
1015  const bool face_flip,
1016  const bool face_rotation) const
1017 {
1018  Assert (face_index < this->dofs_per_face,
1019  ExcIndexRange(face_index, 0, this->dofs_per_face));
1022 
1023 //TODO: we could presumably solve the 3d case below using the
1024 // adjust_quad_dof_index_for_face_orientation_table field. for the
1025 // 2d case, we can't use adjust_line_dof_index_for_line_orientation_table
1026 // since that array is empty (presumably because we thought that
1027 // there are no flipped edges in 2d, but these can happen in
1028 // DoFTools::make_periodicity_constraints, for example). so we
1029 // would need to either fill this field, or rely on derived classes
1030 // implementing this function, as we currently do
1031 
1032  // we need to distinguish between DoFs on vertices, lines and in 3d quads.
1033  // do so in a sequence of if-else statements
1034  if (face_index < this->first_face_line_index)
1035  // DoF is on a vertex
1036  {
1037  // get the number of the vertex on the face that corresponds to this DoF,
1038  // along with the number of the DoF on this vertex
1039  const unsigned int face_vertex = face_index / this->dofs_per_vertex;
1040  const unsigned int dof_index_on_vertex = face_index % this->dofs_per_vertex;
1041 
1042  // then get the number of this vertex on the cell and translate
1043  // this to a DoF number on the cell
1044  return (GeometryInfo<dim>::face_to_cell_vertices(face, face_vertex,
1045  face_orientation,
1046  face_flip,
1047  face_rotation)
1048  * this->dofs_per_vertex
1049  +
1050  dof_index_on_vertex);
1051  }
1052  else if (face_index < this->first_face_quad_index)
1053  // DoF is on a face
1054  {
1055  // do the same kind of translation as before. we need to only consider
1056  // DoFs on the lines, i.e., ignoring those on the vertices
1057  const unsigned int index = face_index - this->first_face_line_index;
1058 
1059  const unsigned int face_line = index / this->dofs_per_line;
1060  const unsigned int dof_index_on_line = index % this->dofs_per_line;
1061 
1062  // we now also need to adjust the line index for the case of
1063  // face orientation, face flips, etc
1064  unsigned int adjusted_dof_index_on_line;
1065  switch (dim)
1066  {
1067  case 1:
1068  Assert (false, ExcInternalError());
1069 
1070  case 2:
1071  // in 2d, only face_flip has a meaning. if it is set, consider
1072  // dofs in reverse order
1073  if (face_flip == false)
1074  adjusted_dof_index_on_line = dof_index_on_line;
1075  else
1076  adjusted_dof_index_on_line = this->dofs_per_line - dof_index_on_line - 1;
1077  break;
1078 
1079  case 3:
1080  // in 3d, things are difficult. someone will have to think
1081  // about how this code here should look like, by drawing a bunch
1082  // of pictures of how all the faces can look like with the various
1083  // flips and rotations.
1084  //
1085  // that said, the Q2 case is easy enough to implement, as is the case
1086  // where everything is in standard orientation
1087  Assert ((this->dofs_per_line <= 1) ||
1088  ((face_orientation == true) &&
1089  (face_flip == false) &&
1090  (face_rotation == false)),
1091  ExcNotImplemented());
1092  adjusted_dof_index_on_line = dof_index_on_line;
1093  break;
1094  }
1095 
1096  return (this->first_line_index
1097  + GeometryInfo<dim>::face_to_cell_lines(face, face_line,
1098  face_orientation,
1099  face_flip,
1100  face_rotation)
1101  * this->dofs_per_line
1102  +
1103  adjusted_dof_index_on_line);
1104  }
1105  else
1106  // DoF is on a quad
1107  {
1108  Assert (dim >= 3, ExcInternalError());
1109 
1110  // ignore vertex and line dofs
1111  const unsigned int index = face_index - this->first_face_quad_index;
1112 
1113  // the same is true here as above for the 3d case -- someone will
1114  // just have to draw a bunch of pictures. in the meantime,
1115  // we can implement the Q2 case in which it is simple
1116  Assert ((this->dofs_per_quad <= 1) ||
1117  ((face_orientation == true) &&
1118  (face_flip == false) &&
1119  (face_rotation == false)),
1120  ExcNotImplemented());
1121  return (this->first_quad_index
1122  + face * this->dofs_per_quad
1123  + index);
1124  }
1125 }
1126 
1127 
1128 
1129 
1130 template <class PolynomialType, int dim, int spacedim>
1131 std::vector<unsigned int>
1133 {
1135  AssertThrow(degree>0, typename FEQ::ExcFEQCannotHaveDegree0());
1136  std::vector<unsigned int> dpo(dim+1, 1U);
1137  for (unsigned int i=1; i<dpo.size(); ++i)
1138  dpo[i]=dpo[i-1]*(degree-1);
1139  return dpo;
1140 }
1141 
1142 
1143 
1144 template <class PolynomialType, int dim, int spacedim>
1145 void
1147 (const std::vector<Point<1> > &points)
1148 {
1149  Implementation::initialize_constraints (points, *this);
1150 }
1151 
1152 
1153 
1154 template <class PolynomialType, int dim, int spacedim>
1155 const FullMatrix<double> &
1157 ::get_prolongation_matrix (const unsigned int child,
1158  const RefinementCase<dim> &refinement_case) const
1159 {
1162  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
1163  ExcMessage("Prolongation matrices are only available for refined cells!"));
1164  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
1165  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
1166 
1167  // initialization upon first request
1168  if (this->prolongation[refinement_case-1][child].n() == 0)
1169  {
1170  Threads::Mutex::ScopedLock lock(this->mutex);
1171 
1172  // if matrix got updated while waiting for the lock
1173  if (this->prolongation[refinement_case-1][child].n() ==
1174  this->dofs_per_cell)
1175  return this->prolongation[refinement_case-1][child];
1176 
1177  // distinguish q/q_dg0 case: only treat Q dofs first
1178  const unsigned int q_dofs_per_cell = Utilities::fixed_power<dim>(q_degree+1);
1179 
1180  // compute the interpolation matrices in much the same way as we do for
1181  // the constraints. it's actually simpler here, since we don't have this
1182  // weird renumbering stuff going on. The trick is again that we the
1183  // interpolation matrix is formed by a permutation of the indices of the
1184  // cell matrix. The value eps is used a threshold to decide when certain
1185  // evaluations of the Lagrange polynomials are zero or one.
1186  const double eps = 1e-15*q_degree*dim;
1187 
1188 #ifdef DEBUG
1189  // in DEBUG mode, check that the evaluation of support points in the
1190  // current numbering gives the identity operation
1191  for (unsigned int i=0; i<q_dofs_per_cell; ++i)
1192  {
1193  Assert (std::fabs (1.-this->poly_space.compute_value
1194  (i, this->unit_support_points[i])) < eps,
1195  ExcInternalError("The Lagrange polynomial does not evaluate "
1196  "to one or zero in a nodal point. "
1197  "This typically indicates that the "
1198  "polynomial interpolation is "
1199  "ill-conditioned such that round-off "
1200  "prevents the sum to be one."));
1201  for (unsigned int j=0; j<q_dofs_per_cell; ++j)
1202  if (j!=i)
1203  Assert (std::fabs (this->poly_space.compute_value
1204  (i, this->unit_support_points[j])) < eps,
1205  ExcInternalError("The Lagrange polynomial does not evaluate "
1206  "to one or zero in a nodal point. "
1207  "This typically indicates that the "
1208  "polynomial interpolation is "
1209  "ill-conditioned such that round-off "
1210  "prevents the sum to be one."));
1211  }
1212 #endif
1213 
1214  // to efficiently evaluate the polynomial at the subcell, make use of
1215  // the tensor product structure of this element and only evaluate 1D
1216  // information from the polynomial. This makes the cost of this function
1217  // almost negligible also for high order elements
1218  const unsigned int dofs1d = q_degree+1;
1219  std::vector<Table<2,double> >
1220  subcell_evaluations (dim, Table<2,double>(dofs1d, dofs1d));
1221  const std::vector<unsigned int> &index_map_inverse =
1222  this->poly_space.get_numbering_inverse();
1223 
1224  // helper value: step size how to walk through diagonal and how many
1225  // points we have left apart from the first dimension
1226  unsigned int step_size_diag = 0;
1227  {
1228  unsigned int factor = 1;
1229  for (unsigned int d=0; d<dim; ++d)
1230  {
1231  step_size_diag += factor;
1232  factor *= dofs1d;
1233  }
1234  }
1235 
1236  FullMatrix<double> prolongate (this->dofs_per_cell, this->dofs_per_cell);
1237 
1238  // go through the points in diagonal to capture variation in all
1239  // directions simultaneously
1240  for (unsigned int j=0; j<dofs1d; ++j)
1241  {
1242  const unsigned int diag_comp = index_map_inverse[j*step_size_diag];
1243  const Point<dim> p_subcell = this->unit_support_points[diag_comp];
1244  const Point<dim> p_cell =
1246  refinement_case);
1247  for (unsigned int i=0; i<dofs1d; ++i)
1248  for (unsigned int d=0; d<dim; ++d)
1249  {
1250  // evaluate along line where only x is different from zero
1251  Point<dim> point;
1252  point[0] = p_cell[d];
1253  const double cell_value =
1254  this->poly_space.compute_value(index_map_inverse[i], point);
1255 
1256  // cut off values that are too small. note that we have here
1257  // Lagrange interpolation functions, so they should be zero at
1258  // almost all points, and one at the others, at least on the
1259  // subcells. so set them to their exact values
1260  //
1261  // the actual cut-off value is somewhat fuzzy, but it works
1262  // for 2e-13*degree*dim (see above), which is kind of
1263  // reasonable given that we compute the values of the
1264  // polynomials via an degree-step recursion and then multiply
1265  // the 1d-values. this gives us a linear growth in degree*dim,
1266  // times a small constant.
1267  //
1268  // the embedding matrix is given by applying the inverse of
1269  // the subcell matrix on the cell_interpolation matrix. since
1270  // the subcell matrix is actually only a permutation vector,
1271  // all we need to do is to switch the rows we write the data
1272  // into. moreover, cut off very small values here
1273  if (std::fabs(cell_value) < eps)
1274  subcell_evaluations[d](j,i) = 0;
1275  else
1276  subcell_evaluations[d](j,i) = cell_value;
1277  }
1278  }
1279 
1280  // now expand from 1D info. block innermost dimension (x_0) in order to
1281  // avoid difficult checks at innermost loop
1282  unsigned int j_indices[dim];
1283  FE_Q_Helper::zero_indices<dim> (j_indices);
1284  for (unsigned int j=0; j<q_dofs_per_cell; j+=dofs1d)
1285  {
1286  unsigned int i_indices[dim];
1287  FE_Q_Helper::zero_indices<dim> (i_indices);
1288  for (unsigned int i=0; i<q_dofs_per_cell; i+=dofs1d)
1289  {
1290  double val_extra_dim = 1.;
1291  for (unsigned int d=1; d<dim; ++d)
1292  val_extra_dim *= subcell_evaluations[d](j_indices[d-1],
1293  i_indices[d-1]);
1294 
1295  // innermost sum where we actually compute. the same as
1296  // prolongate(j,i) = this->poly_space.compute_value (i, p_cell)
1297  for (unsigned int jj=0; jj<dofs1d; ++jj)
1298  {
1299  const unsigned int j_ind = index_map_inverse[j+jj];
1300  for (unsigned int ii=0; ii<dofs1d; ++ii)
1301  prolongate(j_ind,index_map_inverse[i+ii])
1302  = val_extra_dim * subcell_evaluations[0](jj,ii);
1303  }
1304 
1305  // update indices that denote the tensor product position. a bit
1306  // fuzzy and therefore not done for innermost x_0 direction
1307  FE_Q_Helper::increment_indices<dim> (i_indices, dofs1d);
1308  }
1309  Assert (i_indices[dim-1] == 1, ExcInternalError());
1310  FE_Q_Helper::increment_indices<dim> (j_indices, dofs1d);
1311  }
1312 
1313  // the discontinuous node is simply mapped on the discontinuous node on
1314  // the child element
1315  if (q_dofs_per_cell < this->dofs_per_cell)
1316  prolongate(q_dofs_per_cell,q_dofs_per_cell) = 1.;
1317 
1318  // and make sure that the row sum is 1. this must be so since for this
1319  // element, the shape functions add up to one
1320 #ifdef DEBUG
1321  for (unsigned int row=0; row<this->dofs_per_cell; ++row)
1322  {
1323  double sum = 0;
1324  for (unsigned int col=0; col<this->dofs_per_cell; ++col)
1325  sum += prolongate(row,col);
1326  Assert (std::fabs(sum-1.) <
1327  std::max(eps, 5e-16*std::sqrt(this->dofs_per_cell)),
1328  ExcInternalError("The entries in a row of the local "
1329  "prolongation matrix do not add to one. "
1330  "This typically indicates that the "
1331  "polynomial interpolation is "
1332  "ill-conditioned such that round-off "
1333  "prevents the sum to be one."));
1334  }
1335 #endif
1336 
1337  // swap matrices
1338  prolongate.swap(const_cast<FullMatrix<double> &>
1339  (this->prolongation[refinement_case-1][child]));
1340  }
1341 
1342  // finally return the matrix
1343  return this->prolongation[refinement_case-1][child];
1344 }
1345 
1346 
1347 
1348 template <class PolynomialType, int dim, int spacedim>
1349 const FullMatrix<double> &
1351 ::get_restriction_matrix (const unsigned int child,
1352  const RefinementCase<dim> &refinement_case) const
1353 {
1356  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
1357  ExcMessage("Restriction matrices are only available for refined cells!"));
1358  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
1359  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
1360 
1361  // initialization upon first request
1362  if (this->restriction[refinement_case-1][child].n() == 0)
1363  {
1364  Threads::Mutex::ScopedLock lock(this->mutex);
1365 
1366  // if matrix got updated while waiting for the lock...
1367  if (this->restriction[refinement_case-1][child].n() ==
1368  this->dofs_per_cell)
1369  return this->restriction[refinement_case-1][child];
1370 
1371  FullMatrix<double> my_restriction(this->dofs_per_cell, this->dofs_per_cell);
1372  // distinguish q/q_dg0 case
1373  const unsigned int q_dofs_per_cell = Utilities::fixed_power<dim>(q_degree+1);
1374 
1375  // for Lagrange interpolation polynomials based on equidistant points,
1376  // construction of the restriction matrices is relatively simple. the
1377  // reason is that in this case the interpolation points on the mother
1378  // cell are always also interpolation points for some shape function on
1379  // one or the other child.
1380  //
1381  // in the general case with non-equidistant points, we need to actually
1382  // do an interpolation. thus, we take the interpolation points on the
1383  // mother cell and evaluate the shape functions of the child cell on
1384  // those points. it does not hurt in the equidistant case because then
1385  // simple one shape function evaluates to one and the others to zero.
1386  //
1387  // this element is non-additive in all its degrees of freedom by
1388  // default, which requires care in downstream use. fortunately, even the
1389  // interpolation on non-equidistant points is invariant under the
1390  // assumption that whenever a row makes a non-zero contribution to the
1391  // mother's residual, the correct value is interpolated.
1392 
1393  const double eps = 1e-15*q_degree*dim;
1394  const std::vector<unsigned int> &index_map_inverse =
1395  this->poly_space.get_numbering_inverse();
1396 
1397  const unsigned int dofs1d = q_degree+1;
1398  std::vector<Tensor<1,dim> > evaluations1d (dofs1d);
1399 
1400  my_restriction.reinit(this->dofs_per_cell, this->dofs_per_cell);
1401 
1402  for (unsigned int i=0; i<q_dofs_per_cell; ++i)
1403  {
1404  unsigned int mother_dof = index_map_inverse[i];
1405  const Point<dim> p_cell = this->unit_support_points[mother_dof];
1406 
1407  // check whether this interpolation point is inside this child cell
1408  const Point<dim> p_subcell
1410  refinement_case);
1412  {
1413  // same logic as in initialize_embedding to evaluate the
1414  // polynomial faster than from the tensor product: since we
1415  // evaluate all polynomials, it is much faster to just compute
1416  // the 1D values for all polynomials before and then get the
1417  // dim-data.
1418  for (unsigned int j=0; j<dofs1d; ++j)
1419  for (unsigned int d=0; d<dim; ++d)
1420  {
1421  Point<dim> point;
1422  point[0] = p_subcell[d];
1423  evaluations1d[j][d] =
1424  this->poly_space.compute_value(index_map_inverse[j], point);
1425  }
1426  unsigned int j_indices[dim];
1427  FE_Q_Helper::zero_indices<dim> (j_indices);
1428  double sum_check = 0;
1429  for (unsigned int j = 0; j<q_dofs_per_cell; j += dofs1d)
1430  {
1431  double val_extra_dim = 1.;
1432  for (unsigned int d=1; d<dim; ++d)
1433  val_extra_dim *= evaluations1d[j_indices[d-1]][d];
1434  for (unsigned int jj=0; jj<dofs1d; ++jj)
1435  {
1436 
1437  // find the child shape function(s) corresponding to
1438  // this point. Usually this is just one function;
1439  // however, when we use FE_Q on arbitrary nodes a parent
1440  // support point might not be a child support point, and
1441  // then we will get more than one nonzero value per
1442  // row. Still, the values should sum up to 1
1443  const double val
1444  = val_extra_dim * evaluations1d[jj][0];
1445  const unsigned int child_dof =
1446  index_map_inverse[j+jj];
1447  if (std::fabs (val-1.) < eps)
1448  my_restriction(mother_dof,child_dof)=1.;
1449  else if (std::fabs(val) > eps)
1450  my_restriction(mother_dof,child_dof)=val;
1451  sum_check += val;
1452  }
1453  FE_Q_Helper::increment_indices<dim> (j_indices, dofs1d);
1454  }
1455  Assert (std::fabs(sum_check-1) <
1456  std::max(eps, 5e-16*std::sqrt(this->dofs_per_cell)),
1457  ExcInternalError("The entries in a row of the local "
1458  "restriction matrix do not add to one. "
1459  "This typically indicates that the "
1460  "polynomial interpolation is "
1461  "ill-conditioned such that round-off "
1462  "prevents the sum to be one."));
1463  }
1464 
1465  // part for FE_Q_DG0
1466  if (q_dofs_per_cell < this->dofs_per_cell)
1467  my_restriction(this->dofs_per_cell-1,this->dofs_per_cell-1) =
1469  }
1470 
1471  // swap the just computed restriction matrix into the
1472  // element of the vector stored in the base class
1473  my_restriction.swap(const_cast<FullMatrix<double> &>
1474  (this->restriction[refinement_case-1][child]));
1475  }
1476 
1477  return this->restriction[refinement_case-1][child];
1478 }
1479 
1480 
1481 
1482 //---------------------------------------------------------------------------
1483 // Data field initialization
1484 //---------------------------------------------------------------------------
1485 
1486 
1487 template <class PolynomialType, int dim, int spacedim>
1488 bool
1490 (const unsigned int shape_index,
1491  const unsigned int face_index) const
1492 {
1493  Assert (shape_index < this->dofs_per_cell,
1494  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
1497 
1498  // in 1d, things are simple. since there is only one degree of freedom per
1499  // vertex in this class, the first is on vertex 0 (==face 0 in some sense),
1500  // the second on face 1:
1501  if (dim == 1)
1502  return (((shape_index == 0) && (face_index == 0)) ||
1503  ((shape_index == 1) && (face_index == 1)));
1504 
1505  // first, special-case interior shape functions, since they have no support
1506  // no-where on the boundary
1507  if (((dim==2) && (shape_index>=this->first_quad_index))
1508  ||
1509  ((dim==3) && (shape_index>=this->first_hex_index)))
1510  return false;
1511 
1512  // let's see whether this is a vertex
1513  if (shape_index < this->first_line_index)
1514  {
1515  // for Q elements, there is one dof per vertex, so
1516  // shape_index==vertex_number. check whether this vertex is on the given
1517  // face. thus, for each face, give a list of vertices
1518  const unsigned int vertex_no = shape_index;
1520  ExcInternalError());
1521 
1522  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_face; ++v)
1523  if (GeometryInfo<dim>::face_to_cell_vertices(face_index, v) == vertex_no)
1524  return true;
1525 
1526  return false;
1527  }
1528  else if (shape_index < this->first_quad_index)
1529  // ok, dof is on a line
1530  {
1531  const unsigned int line_index
1532  = (shape_index - this->first_line_index) / this->dofs_per_line;
1534  ExcInternalError());
1535 
1536  // in 2d, the line is the face, so get the line index
1537  if (dim == 2)
1538  return (line_index == face_index);
1539  else if (dim == 3)
1540  {
1541  // silence compiler warning
1542  const unsigned int lines_per_face =
1543  dim == 3 ? GeometryInfo<dim>::lines_per_face : 1;
1544  // see whether the given line is on the given face.
1545  for (unsigned int l=0; l<lines_per_face; ++l)
1546  if (GeometryInfo<3>::face_to_cell_lines(face_index, l) == line_index)
1547  return true;
1548 
1549  return false;
1550  }
1551  else
1552  Assert (false, ExcNotImplemented());
1553  }
1554  else if (shape_index < this->first_hex_index)
1555  // dof is on a quad
1556  {
1557  const unsigned int quad_index
1558  = (shape_index - this->first_quad_index) / this->dofs_per_quad;
1559  Assert (static_cast<signed int>(quad_index) <
1560  static_cast<signed int>(GeometryInfo<dim>::quads_per_cell),
1561  ExcInternalError());
1562 
1563  // in 2d, cell bubble are zero on all faces. but we have treated this
1564  // case above already
1565  Assert (dim != 2, ExcInternalError());
1566 
1567  // in 3d, quad_index=face_index
1568  if (dim == 3)
1569  return (quad_index == face_index);
1570  else
1571  Assert (false, ExcNotImplemented());
1572  }
1573  else
1574  // dof on hex
1575  {
1576  // can only happen in 3d, but this case has already been covered above
1577  Assert (false, ExcNotImplemented());
1578  return false;
1579  }
1580 
1581  // we should not have gotten here
1582  Assert (false, ExcInternalError());
1583  return false;
1584 }
1585 
1586 
1587 
1588 template <typename PolynomialType, int dim, int spacedim>
1589 std::pair<Table<2,bool>, std::vector<unsigned int> >
1591 {
1592  Table<2,bool> constant_modes(1, this->dofs_per_cell);
1593  // We here just care for the constant mode due to the polynomial space
1594  // without any enrichments
1595  // As there may be more constant modes derived classes may to implement this
1596  // themselves. An example for this is FE_Q_DG0.
1597  for (unsigned int i=0; i<Utilities::fixed_power<dim>(q_degree+1); ++i)
1598  constant_modes(0, i) = true;
1599  return std::pair<Table<2,bool>, std::vector<unsigned int> >
1600  (constant_modes, std::vector<unsigned int>(1, 0));
1601 }
1602 
1603 
1604 // explicit instantiations
1605 #include "fe_q_base.inst"
1606 
1607 DEAL_II_NAMESPACE_CLOSE
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
size_type m() const
FE_Q_Base(const PolynomialType &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
Definition: fe_q_base.cc:415
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
void swap(TableBase< N, T > &v)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_q_base.cc:1132
FullMatrix< double > interface_constraints
Definition: fe.h:2211
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:684
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:568
Definition: fe_dgq.h:105
void hierarchic_to_lexicographic_numbering(unsigned int degree, std::vector< unsigned int > &h2l)
Definition: fe_tools.cc:2928
const unsigned int degree
Definition: fe_base.h:311
const Point< dim > & point(const unsigned int i) const
static Point< dim > cell_to_child_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:788
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Definition: point.h:89
Definition: fe_dgp.h:309
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_q_base.cc:1157
size_type n() const
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:313
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_q_base.cc:1351
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)
virtual bool hp_constraints_are_implemented() const
Definition: fe_q_base.cc:673
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe_q_base.cc:1012
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:725
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_q_base.cc:1490
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_q_base.cc:855
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_q_base.cc:1590
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
void lexicographic_to_hierarchic_numbering(const FiniteElementData< dim > &fe_data, std::vector< unsigned int > &l2h)
Definition: fe_tools.cc:3120
void initialize(const std::vector< Point< 1 > > &support_points_1d)
Definition: fe_q_base.cc:430
void initialize_unit_face_support_points(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:920
void initialize_constraints(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:1147
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:581
static void initialize_constraints(const std::vector< Point< 1 > > &, FE_Q_Base< PolynomialType, 1, spacedim > &)
Definition: fe_q_base.cc:116
void initialize_quad_dof_index_permutation()
Definition: fe_q_base.cc:944
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
void initialize_unit_support_points(const std::vector< Point< 1 > > &points)
Definition: fe_q_base.cc:903
const unsigned int dofs_per_face
Definition: fe_base.h:288
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_q_base.cc:474
static ::ExceptionBase & ExcNotImplemented()
PolynomialType poly_space
Definition: fe_poly.h:444
TableIndices< 2 > interface_constraints_size() const
Definition: fe.cc:823
static ::ExceptionBase & ExcInternalError()