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