Reference documentation for deal.II version 9.0.0
fe_dgq.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2018 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/quadrature_lib.h>
19 #include <deal.II/lac/vector.h>
20 #include <deal.II/fe/fe.h>
21 #include <deal.II/fe/fe_dgq.h>
22 #include <deal.II/fe/fe_tools.h>
23 
24 
25 #include <iostream>
26 #include <sstream>
27 #include <deal.II/base/std_cxx14/memory.h>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 namespace internal
34 {
35  namespace FE_DGQ
36  {
37  namespace
38  {
39  std::vector<Point<1> >
40  get_QGaussLobatto_points (const unsigned int degree)
41  {
42  if (degree > 0)
43  return QGaussLobatto<1>(degree+1).get_points();
44  else
45  return std::vector<Point<1> >(1, Point<1>(0.5));
46  }
47  }
48  }
49 }
50 
51 
52 
53 template <int dim, int spacedim>
54 FE_DGQ<dim, spacedim>::FE_DGQ (const unsigned int degree)
55  :
56  FE_Poly<TensorProductPolynomials<dim>, dim, spacedim>
57  (TensorProductPolynomials<dim>(Polynomials::generate_complete_Lagrange_basis
58  (internal::FE_DGQ::get_QGaussLobatto_points(degree))),
59  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
60  std::vector<bool>(FiniteElementData<dim>(get_dpo_vector(degree),1, degree).dofs_per_cell, true),
61  std::vector<ComponentMask>(FiniteElementData<dim>(get_dpo_vector(degree),1, degree).dofs_per_cell, std::vector<bool>(1,true)))
62 {
63  // Compute support points, which are the tensor product of the Lagrange
64  // interpolation points in the constructor.
65  this->unit_support_points =
66  Quadrature<dim>(internal::FE_DGQ::get_QGaussLobatto_points(degree)).get_points();
67 
68  // do not initialize embedding and restriction here. these matrices are
69  // initialized on demand in get_restriction_matrix and
70  // get_prolongation_matrix
71 
72  // note: no face support points for DG elements
73 }
74 
75 
76 
77 template <int dim, int spacedim>
79  :
80  FE_Poly<TensorProductPolynomials<dim>, dim, spacedim> (
81  TensorProductPolynomials<dim>(polynomials),
82  FiniteElementData<dim>(get_dpo_vector(polynomials.size()-1), 1, polynomials.size()-1, FiniteElementData<dim>::L2),
83  std::vector<bool>(FiniteElementData<dim>(get_dpo_vector(polynomials.size()-1),1, polynomials.size()-1).dofs_per_cell, true),
84  std::vector<ComponentMask>(FiniteElementData<dim>(get_dpo_vector(polynomials.size()-1),1, polynomials.size()-1).dofs_per_cell, std::vector<bool>(1,true)))
85 {
86  // No support points can be defined in general. Derived classes might define
87  // support points like the class FE_DGQArbitraryNodes
88 
89  // do not initialize embedding and restriction here. these matrices are
90  // initialized on demand in get_restriction_matrix and
91  // get_prolongation_matrix
92 }
93 
94 
95 
96 template <int dim, int spacedim>
97 std::string
99 {
100  // note that the FETools::get_fe_by_name function depends on the
101  // particular format of the string this function returns, so they have to be
102  // kept in sync
103 
104  std::ostringstream namebuf;
105  namebuf << "FE_DGQ<"
106  << Utilities::dim_string(dim,spacedim)
107  << ">(" << this->degree << ")";
108  return namebuf.str();
109 }
110 
111 
112 
113 template <int dim, int spacedim>
114 void
117  std::vector<double> &nodal_values) const
118 {
119  AssertDimension (support_point_values.size(),
120  this->get_unit_support_points().size());
121  AssertDimension (support_point_values.size(),
122  nodal_values.size());
123  AssertDimension (this->dofs_per_cell,
124  nodal_values.size());
125 
126  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
127  {
128  AssertDimension (support_point_values[i].size(), 1);
129 
130  nodal_values[i] = support_point_values[i](0);
131  }
132 }
133 
134 
135 
136 template <int dim, int spacedim>
137 std::unique_ptr<FiniteElement<dim,spacedim> >
139 {
140  return std_cxx14::make_unique<FE_DGQ<dim, spacedim>>(*this);
141 }
142 
143 
144 //---------------------------------------------------------------------------
145 // Auxiliary functions
146 //---------------------------------------------------------------------------
147 
148 
149 template <int dim, int spacedim>
150 std::vector<unsigned int>
152 {
153  std::vector<unsigned int> dpo(dim+1, 0U);
154  dpo[dim] = deg+1;
155  for (unsigned int i=1; i<dim; ++i)
156  dpo[dim] *= deg+1;
157  return dpo;
158 }
159 
160 
161 
162 template <int dim, int spacedim>
163 void
164 FE_DGQ<dim, spacedim>::rotate_indices (std::vector<unsigned int> &numbers,
165  const char direction) const
166 {
167  const unsigned int n = this->degree+1;
168  unsigned int s = n;
169  for (unsigned int i=1; i<dim; ++i)
170  s *= n;
171  numbers.resize (s);
172 
173  unsigned int l = 0;
174 
175  if (dim==1)
176  {
177  // Mirror around midpoint
178  for (unsigned int i=n; i>0;)
179  numbers[l++]=--i;
180  }
181  else
182  {
183  switch (direction)
184  {
185  // Rotate xy-plane
186  // counter-clockwise
187  case 'z':
188  for (unsigned int iz=0; iz<((dim>2) ? n:1); ++iz)
189  for (unsigned int j=0; j<n; ++j)
190  for (unsigned int i=0; i<n; ++i)
191  {
192  unsigned int k = n*i-j+n-1 + n*n*iz;
193  numbers[l++] = k;
194  }
195  break;
196  // Rotate xy-plane
197  // clockwise
198  case 'Z':
199  for (unsigned int iz=0; iz<((dim>2) ? n:1); ++iz)
200  for (unsigned int iy=0; iy<n; ++iy)
201  for (unsigned int ix=0; ix<n; ++ix)
202  {
203  unsigned int k = n*ix-iy+n-1 + n*n*iz;
204  numbers[k] = l++;
205  }
206  break;
207  // Rotate yz-plane
208  // counter-clockwise
209  case 'x':
210  Assert (dim>2, ExcDimensionMismatch (dim,3));
211  for (unsigned int iz=0; iz<n; ++iz)
212  for (unsigned int iy=0; iy<n; ++iy)
213  for (unsigned int ix=0; ix<n; ++ix)
214  {
215  unsigned int k = n*(n*iy-iz+n-1) + ix;
216  numbers[l++] = k;
217  }
218  break;
219  // Rotate yz-plane
220  // clockwise
221  case 'X':
222  Assert (dim>2, ExcDimensionMismatch (dim,3));
223  for (unsigned int iz=0; iz<n; ++iz)
224  for (unsigned int iy=0; iy<n; ++iy)
225  for (unsigned int ix=0; ix<n; ++ix)
226  {
227  unsigned int k = n*(n*iy-iz+n-1) + ix;
228  numbers[k] = l++;
229  }
230  break;
231  default:
232  Assert (false, ExcNotImplemented ());
233  }
234  }
235 }
236 
237 
238 
239 template <int dim, int spacedim>
240 void
243  FullMatrix<double> &interpolation_matrix) const
244 {
245  // this is only implemented, if the
246  // source FE is also a
247  // DGQ element
248  typedef FiniteElement<dim, spacedim> FE;
249  AssertThrow ((dynamic_cast<const FE_DGQ<dim, spacedim>*>(&x_source_fe) != nullptr),
250  typename FE::ExcInterpolationNotImplemented() );
251 
252  // ok, source is a Q element, so
253  // we will be able to do the work
254  const FE_DGQ<dim, spacedim> &source_fe
255  = dynamic_cast<const FE_DGQ<dim, spacedim>&>(x_source_fe);
256 
257  Assert (interpolation_matrix.m() == this->dofs_per_cell,
258  ExcDimensionMismatch (interpolation_matrix.m(),
259  this->dofs_per_cell));
260  Assert (interpolation_matrix.n() == source_fe.dofs_per_cell,
261  ExcDimensionMismatch (interpolation_matrix.n(),
262  source_fe.dofs_per_cell));
263 
264 
265  // compute the interpolation
266  // matrices in much the same way as
267  // we do for the embedding matrices
268  // from mother to child.
269  FullMatrix<double> cell_interpolation (this->dofs_per_cell,
270  this->dofs_per_cell);
271  FullMatrix<double> source_interpolation (this->dofs_per_cell,
272  source_fe.dofs_per_cell);
273  FullMatrix<double> tmp (this->dofs_per_cell,
274  source_fe.dofs_per_cell);
275  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
276  {
277  // generate a point on this
278  // cell and evaluate the
279  // shape functions there
280  const Point<dim> p = this->unit_support_points[j];
281  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
282  cell_interpolation(j,i)
283  = this->poly_space.compute_value (i, p);
284 
285  for (unsigned int i=0; i<source_fe.dofs_per_cell; ++i)
286  source_interpolation(j,i)
287  = source_fe.poly_space.compute_value (i, p);
288  }
289 
290  // then compute the
291  // interpolation matrix matrix
292  // for this coordinate
293  // direction
294  cell_interpolation.gauss_jordan ();
295  cell_interpolation.mmult (interpolation_matrix,
296  source_interpolation);
297 
298  // cut off very small values
299  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
300  for (unsigned int j=0; j<source_fe.dofs_per_cell; ++j)
301  if (std::fabs(interpolation_matrix(i,j)) < 1e-15)
302  interpolation_matrix(i,j) = 0.;
303 
304  // make sure that the row sum of
305  // each of the matrices is 1 at
306  // this point. this must be so
307  // since the shape functions sum up
308  // to 1
309  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
310  {
311  double sum = 0.;
312  for (unsigned int j=0; j<source_fe.dofs_per_cell; ++j)
313  sum += interpolation_matrix(i,j);
314 
315  Assert (std::fabs(sum-1) < 5e-14*std::max(this->degree,1U)*dim,
316  ExcInternalError());
317  }
318 }
319 
320 
321 
322 template <int dim, int spacedim>
323 void
326  FullMatrix<double> &interpolation_matrix) const
327 {
328  // this is only implemented, if the source
329  // FE is also a DGQ element. in that case,
330  // both elements have no dofs on their
331  // faces and the face interpolation matrix
332  // is necessarily empty -- i.e. there isn't
333  // much we need to do here.
334  (void)interpolation_matrix;
335  typedef FiniteElement<dim,spacedim> FE;
336  AssertThrow ((dynamic_cast<const FE_DGQ<dim, spacedim>*>(&x_source_fe) != nullptr),
337  typename FE::ExcInterpolationNotImplemented());
338 
339  Assert (interpolation_matrix.m() == 0,
340  ExcDimensionMismatch (interpolation_matrix.m(),
341  0));
342  Assert (interpolation_matrix.n() == 0,
343  ExcDimensionMismatch (interpolation_matrix.m(),
344  0));
345 }
346 
347 
348 
349 template <int dim, int spacedim>
350 void
353  const unsigned int,
354  FullMatrix<double> &interpolation_matrix) const
355 {
356  // this is only implemented, if the source
357  // FE is also a DGQ element. in that case,
358  // both elements have no dofs on their
359  // faces and the face interpolation matrix
360  // is necessarily empty -- i.e. there isn't
361  // much we need to do here.
362  (void)interpolation_matrix;
363  typedef FiniteElement<dim, spacedim> FE;
364  AssertThrow ((dynamic_cast<const FE_DGQ<dim, spacedim>*>(&x_source_fe) != nullptr),
365  typename FE::ExcInterpolationNotImplemented());
366 
367  Assert (interpolation_matrix.m() == 0,
368  ExcDimensionMismatch (interpolation_matrix.m(),
369  0));
370  Assert (interpolation_matrix.n() == 0,
371  ExcDimensionMismatch (interpolation_matrix.m(),
372  0));
373 }
374 
375 
376 
377 template <int dim, int spacedim>
378 const FullMatrix<double> &
380 ::get_prolongation_matrix (const unsigned int child,
381  const RefinementCase<dim> &refinement_case) const
382 {
385  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
386  ExcMessage("Prolongation matrices are only available for refined cells!"));
387  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
388  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
389 
390  // initialization upon first request
391  if (this->prolongation[refinement_case-1][child].n() == 0)
392  {
393  Threads::Mutex::ScopedLock lock(this->mutex);
394 
395  // if matrix got updated while waiting for the lock
396  if (this->prolongation[refinement_case-1][child].n() ==
397  this->dofs_per_cell)
398  return this->prolongation[refinement_case-1][child];
399 
400  // now do the work. need to get a non-const version of data in order to
401  // be able to modify them inside a const function
402  FE_DGQ<dim,spacedim> &this_nonconst = const_cast<FE_DGQ<dim,spacedim>& >(*this);
403  if (refinement_case == RefinementCase<dim>::isotropic_refinement)
404  {
405  std::vector<std::vector<FullMatrix<double> > >
406  isotropic_matrices(RefinementCase<dim>::isotropic_refinement);
407  isotropic_matrices.back().
408  resize(GeometryInfo<dim>::n_children(RefinementCase<dim>(refinement_case)),
409  FullMatrix<double>(this->dofs_per_cell, this->dofs_per_cell));
410  if (dim == spacedim)
411  FETools::compute_embedding_matrices (*this, isotropic_matrices, true);
412  else
414  isotropic_matrices, true);
415  this_nonconst.prolongation[refinement_case-1].swap(isotropic_matrices.back());
416  }
417  else
418  {
419  // must compute both restriction and prolongation matrices because
420  // we only check for their size and the reinit call initializes them
421  // all
423  if (dim == spacedim)
424  {
425  FETools::compute_embedding_matrices (*this, this_nonconst.prolongation);
426  FETools::compute_projection_matrices (*this, this_nonconst.restriction);
427  }
428  else
429  {
430  FE_DGQ<dim> tmp(this->degree);
433  }
434  }
435  }
436 
437  // finally return the matrix
438  return this->prolongation[refinement_case-1][child];
439 }
440 
441 
442 
443 template <int dim, int spacedim>
444 const FullMatrix<double> &
446 ::get_restriction_matrix (const unsigned int child,
447  const RefinementCase<dim> &refinement_case) const
448 {
451  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
452  ExcMessage("Restriction matrices are only available for refined cells!"));
453  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
454  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
455 
456  // initialization upon first request
457  if (this->restriction[refinement_case-1][child].n() == 0)
458  {
459  Threads::Mutex::ScopedLock lock(this->mutex);
460 
461  // if matrix got updated while waiting for the lock...
462  if (this->restriction[refinement_case-1][child].n() ==
463  this->dofs_per_cell)
464  return this->restriction[refinement_case-1][child];
465 
466  // now do the work. need to get a non-const version of data in order to
467  // be able to modify them inside a const function
468  FE_DGQ<dim,spacedim> &this_nonconst = const_cast<FE_DGQ<dim,spacedim>& >(*this);
469  if (refinement_case == RefinementCase<dim>::isotropic_refinement)
470  {
471  std::vector<std::vector<FullMatrix<double> > >
472  isotropic_matrices(RefinementCase<dim>::isotropic_refinement);
473  isotropic_matrices.back().
474  resize(GeometryInfo<dim>::n_children(RefinementCase<dim>(refinement_case)),
475  FullMatrix<double>(this->dofs_per_cell, this->dofs_per_cell));
476  if (dim == spacedim)
477  FETools::compute_projection_matrices (*this, isotropic_matrices, true);
478  else
480  isotropic_matrices, true);
481  this_nonconst.restriction[refinement_case-1].swap(isotropic_matrices.back());
482  }
483  else
484  {
485  // must compute both restriction and prolongation matrices because
486  // we only check for their size and the reinit call initializes them
487  // all
489  if (dim == spacedim)
490  {
491  FETools::compute_embedding_matrices (*this, this_nonconst.prolongation);
492  FETools::compute_projection_matrices (*this, this_nonconst.restriction);
493  }
494  else
495  {
496  FE_DGQ<dim> tmp(this->degree);
499  }
500  }
501  }
502 
503  // finally return the matrix
504  return this->restriction[refinement_case-1][child];
505 }
506 
507 
508 
509 template <int dim, int spacedim>
510 bool
512 {
513  return true;
514 }
515 
516 
517 
518 template <int dim, int spacedim>
519 std::vector<std::pair<unsigned int, unsigned int> >
522 {
523  // this element is discontinuous, so by definition there can
524  // be no identities between its dofs and those of any neighbor
525  // (of whichever type the neighbor may be -- after all, we have
526  // no face dofs on this side to begin with)
527  return
528  std::vector<std::pair<unsigned int, unsigned int> > ();
529 }
530 
531 
532 
533 template <int dim, int spacedim>
534 std::vector<std::pair<unsigned int, unsigned int> >
537 {
538  // this element is discontinuous, so by definition there can
539  // be no identities between its dofs and those of any neighbor
540  // (of whichever type the neighbor may be -- after all, we have
541  // no face dofs on this side to begin with)
542  return
543  std::vector<std::pair<unsigned int, unsigned int> > ();
544 }
545 
546 
547 
548 template <int dim, int spacedim>
549 std::vector<std::pair<unsigned int, unsigned int> >
552 {
553  // this element is discontinuous, so by definition there can
554  // be no identities between its dofs and those of any neighbor
555  // (of whichever type the neighbor may be -- after all, we have
556  // no face dofs on this side to begin with)
557  return
558  std::vector<std::pair<unsigned int, unsigned int> > ();
559 }
560 
561 
562 
563 template <int dim, int spacedim>
566 {
567  // this is a discontinuous element, so by definition there will
568  // be no constraints wherever this element comes together with
569  // any other kind of element
571 }
572 
573 
574 
575 template <int dim, int spacedim>
576 bool
577 FE_DGQ<dim, spacedim>::has_support_on_face (const unsigned int shape_index,
578  const unsigned int face_index) const
579 {
580  Assert (shape_index < this->dofs_per_cell,
581  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
584 
585  unsigned int n = this->degree+1;
586 
587  // For DGQ elements that do not define support points, we cannot define
588  // whether they have support at the boundary easily, so return true in any
589  // case
590  if (this->unit_support_points.empty())
591  return true;
592 
593  // for DGQ(0) elements or arbitrary node DGQ with support points not located
594  // at the element boundary, the single shape functions is constant and
595  // therefore lives on the boundary
596  bool support_points_on_boundary = true;
597  for (unsigned int d=0; d<dim; ++d)
598  if (std::abs(this->unit_support_points[0][d]) > 1e-13)
599  support_points_on_boundary = false;
600  for (unsigned int d=0; d<dim; ++d)
601  if (std::abs(this->unit_support_points.back()[d]-1.) > 1e-13)
602  support_points_on_boundary = false;
603  if (support_points_on_boundary == false)
604  return true;
605 
606  unsigned int n2 = n*n;
607 
608  switch (dim)
609  {
610  case 1:
611  {
612  // in 1d, things are simple. since
613  // there is only one degree of
614  // freedom per vertex in this
615  // class, the first is on vertex 0
616  // (==face 0 in some sense), the
617  // second on face 1:
618  return (((shape_index == 0) && (face_index == 0)) ||
619  ((shape_index == this->degree) && (face_index == 1)));
620  };
621 
622  case 2:
623  {
624  if (face_index==0 && (shape_index % n) == 0)
625  return true;
626  if (face_index==1 && (shape_index % n) == this->degree)
627  return true;
628  if (face_index==2 && shape_index < n)
629  return true;
630  if (face_index==3 && shape_index >= this->dofs_per_cell-n)
631  return true;
632  return false;
633  };
634 
635  case 3:
636  {
637  const unsigned int in2 = shape_index % n2;
638 
639  // x=0
640  if (face_index==0 && (shape_index % n) == 0)
641  return true;
642  // x=1
643  if (face_index==1 && (shape_index % n) == n-1)
644  return true;
645  // y=0
646  if (face_index==2 && in2 < n )
647  return true;
648  // y=1
649  if (face_index==3 && in2 >= n2-n)
650  return true;
651  // z=0
652  if (face_index==4 && shape_index < n2)
653  return true;
654  // z=1
655  if (face_index==5 && shape_index >= this->dofs_per_cell - n2)
656  return true;
657  return false;
658  };
659 
660  default:
661  Assert (false, ExcNotImplemented());
662  }
663  return true;
664 }
665 
666 
667 
668 template <int dim, int spacedim>
669 std::pair<Table<2,bool>, std::vector<unsigned int> >
671 {
672  Table<2,bool> constant_modes(1, this->dofs_per_cell);
673  constant_modes.fill(true);
674  return std::pair<Table<2,bool>, std::vector<unsigned int> >
675  (constant_modes, std::vector<unsigned int>(1, 0));
676 }
677 
678 
679 
680 
681 template <int dim, int spacedim>
682 std::size_t
684 {
685  Assert (false, ExcNotImplemented ());
686  return 0;
687 }
688 
689 
690 
691 // ------------------------------ FE_DGQArbitraryNodes -----------------------
692 
693 template <int dim, int spacedim>
695  : FE_DGQ<dim,spacedim>(Polynomials::generate_complete_Lagrange_basis(points.get_points()))
696 {
697  Assert (points.size() > 0,
699  this->unit_support_points = Quadrature<dim>(points).get_points();
700 }
701 
702 
703 
704 template <int dim, int spacedim>
705 std::string
707 {
708  // note that the FETools::get_fe_by_name function does not work for
709  // FE_DGQArbitraryNodes since there is no initialization by a degree value.
710  std::ostringstream namebuf;
711  bool equidistant = true;
712  std::vector<double> points(this->degree+1);
713 
714  std::vector<unsigned int> lexicographic = this->poly_space.get_numbering_inverse();
715  for (unsigned int j=0; j<=this->degree; j++)
716  points[j] = this->unit_support_points[lexicographic[j]][0];
717 
718  // Check whether the support points are equidistant.
719  for (unsigned int j=0; j<=this->degree; j++)
720  if (std::abs(points[j] - (double)j/this->degree) > 1e-15)
721  {
722  equidistant = false;
723  break;
724  }
725  if (this->degree == 0 && std::abs(points[0]-0.5) < 1e-15)
726  equidistant = true;
727 
728  if (equidistant == true)
729  {
730  if (this->degree > 2)
731  namebuf << "FE_DGQArbitraryNodes<" << Utilities::dim_string(dim,spacedim) << ">(QIterated(QTrapez()," << this->degree << "))";
732  else
733  namebuf << "FE_DGQ<" << Utilities::dim_string(dim,spacedim) << ">(" << this->degree << ")";
734  return namebuf.str();
735  }
736 
737  // Check whether the support points come from QGaussLobatto.
738  const QGaussLobatto<1> points_gl(this->degree+1);
739  bool gauss_lobatto = true;
740  for (unsigned int j=0; j<=this->degree; j++)
741  if (points[j] != points_gl.point(j)(0))
742  {
743  gauss_lobatto = false;
744  break;
745  }
746 
747  if (gauss_lobatto == true)
748  {
749  namebuf << "FE_DGQ<" << Utilities::dim_string(dim,spacedim) << ">(" << this->degree << ")";
750  return namebuf.str();
751  }
752 
753  // Check whether the support points come from QGauss.
754  const QGauss<1> points_g(this->degree+1);
755  bool gauss = true;
756  for (unsigned int j=0; j<=this->degree; j++)
757  if (points[j] != points_g.point(j)(0))
758  {
759  gauss = false;
760  break;
761  }
762 
763  if (gauss == true)
764  {
765  namebuf << "FE_DGQArbitraryNodes<" << Utilities::dim_string(dim,spacedim) << ">(QGauss(" << this->degree+1 << "))";
766  return namebuf.str();
767  }
768 
769  // Check whether the support points come from QGauss.
770  const QGaussLog<1> points_glog(this->degree+1);
771  bool gauss_log = true;
772  for (unsigned int j=0; j<=this->degree; j++)
773  if (points[j] != points_glog.point(j)(0))
774  {
775  gauss_log = false;
776  break;
777  }
778 
779  if (gauss_log == true)
780  {
781  namebuf << "FE_DGQArbitraryNodes<" << Utilities::dim_string(dim,spacedim) << ">(QGaussLog(" << this->degree+1 << "))";
782  return namebuf.str();
783  }
784 
785  // All guesses exhausted
786  namebuf << "FE_DGQArbitraryNodes<" << Utilities::dim_string(dim,spacedim) << ">(QUnknownNodes(" << this->degree+1 << "))";
787  return namebuf.str();
788 }
789 
790 
791 
792 template <int dim, int spacedim>
793 void
796  std::vector<double> &nodal_values) const
797 {
798  AssertDimension (support_point_values.size(),
799  this->get_unit_support_points().size());
800  AssertDimension (support_point_values.size(),
801  nodal_values.size());
802  AssertDimension (this->dofs_per_cell,
803  nodal_values.size());
804 
805  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
806  {
807  AssertDimension (support_point_values[i].size(), 1);
808 
809  nodal_values[i] = support_point_values[i](0);
810  }
811 }
812 
813 
814 
815 
816 template <int dim, int spacedim>
817 std::unique_ptr<FiniteElement<dim,spacedim> >
819 {
820  // Construct a dummy quadrature formula containing the FE's nodes:
821  std::vector<Point<1> > qpoints(this->degree+1);
822  std::vector<unsigned int> lexicographic = this->poly_space.get_numbering_inverse();
823  for (unsigned int i=0; i<=this->degree; ++i)
824  qpoints[i] = Point<1>(this->unit_support_points[lexicographic[i]][0]);
825  Quadrature<1> pquadrature(qpoints);
826 
827  return std_cxx14::make_unique<FE_DGQArbitraryNodes<dim,spacedim>>(pquadrature);
828 }
829 
830 
831 
832 // ---------------------------------- FE_DGQLegendre -------------------------
833 
834 template <int dim, int spacedim>
836  : FE_DGQ<dim,spacedim>(Polynomials::Legendre::generate_complete_basis(degree))
837 {}
838 
839 
840 
841 template <int dim, int spacedim>
842 std::pair<Table<2,bool>, std::vector<unsigned int> >
844 {
845  // Legendre represents a constant function by one in the first basis
846  // function and zero in all others
847  Table<2,bool> constant_modes(1, this->dofs_per_cell);
848  constant_modes(0,0) = true;
849  return std::pair<Table<2,bool>, std::vector<unsigned int> >
850  (constant_modes, std::vector<unsigned int>(1, 0));
851 }
852 
853 
854 
855 template <int dim, int spacedim>
856 std::string
858 {
859  return "FE_DGQLegendre<" + Utilities::dim_string(dim,spacedim) + ">("
860  + Utilities::int_to_string(this->degree) + ")";
861 }
862 
863 
864 
865 template <int dim, int spacedim>
866 std::unique_ptr<FiniteElement<dim,spacedim> >
868 {
869  return std_cxx14::make_unique<FE_DGQLegendre<dim,spacedim>>(this->degree);
870 }
871 
872 
873 
874 // ---------------------------------- FE_DGQHermite --------------------------
875 
876 template <int dim, int spacedim>
878  : FE_DGQ<dim,spacedim>(Polynomials::HermiteLikeInterpolation::generate_complete_basis(degree))
879 {}
880 
881 
882 
883 template <int dim, int spacedim>
884 std::string
886 {
887  return "FE_DGQHermite<" + Utilities::dim_string(dim,spacedim) + ">("
888  + Utilities::int_to_string(this->degree) + ")";
889 }
890 
891 
892 
893 template <int dim, int spacedim>
894 std::unique_ptr<FiniteElement<dim,spacedim> >
896 {
897  return std_cxx14::make_unique<FE_DGQHermite<dim,spacedim>>(this->degree);
898 }
899 
900 
901 
902 // explicit instantiations
903 #include "fe_dgq.inst"
904 
905 
906 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcFEHasNoSupportPoints()
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_dgq.cc:795
size_type m() const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:352
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
void gauss_jordan()
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:818
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:565
FE_DGQLegendre(const unsigned int degree)
Definition: fe_dgq.cc:835
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:867
Definition: fe_dgq.h:105
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_dgq.cc:116
const unsigned int degree
Definition: fe_base.h:311
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:325
const Point< dim > & point(const unsigned int i) const
STL namespace.
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
virtual std::size_t memory_consumption() const
Definition: fe_dgq.cc:683
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Definition: point.h:104
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_dgq.cc:380
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2413
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:138
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgq.cc:895
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_dgq.cc:670
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
virtual std::string get_name() const
Definition: fe_dgq.cc:885
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_dgq.cc:151
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_dgq.cc:577
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_dgq.cc:446
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_dgq.cc:242
void rotate_indices(std::vector< unsigned int > &indices, const char direction) const
Definition: fe_dgq.cc:164
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:98
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:176
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
FE_DGQArbitraryNodes(const Quadrature< 1 > &points)
Definition: fe_dgq.cc:694
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
FE_DGQHermite(const unsigned int degree)
Definition: fe_dgq.cc:877
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:551
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:536
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
static ::ExceptionBase & ExcNotImplemented()
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_dgq.cc:843
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgq.cc:521
double compute_value(const unsigned int i, const Point< dim > &p) const
virtual std::string get_name() const
Definition: fe_dgq.cc:857
void fill(InputIterator entries, const bool C_style_indexing=true)
friend class FE_DGQ
Definition: fe_dgq.h:376
virtual bool hp_constraints_are_implemented() const
Definition: fe_dgq.cc:511
virtual std::string get_name() const
Definition: fe_dgq.cc:706
static ::ExceptionBase & ExcInternalError()
virtual std::string get_name() const
Definition: fe_dgq.cc:98