1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2019 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_qprojector_h
17 #define dealii_qprojector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/base/geometry_info.h>
23 #include <deal.II/base/quadrature.h>
24 
25 #include <deal.II/grid/reference_cell.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 /*!@addtogroup Quadrature */
31 /*@{*/
32 
33 
34 /**
35  * This class is a helper class to facilitate the usage of quadrature formulae
36  * on faces or subfaces of cells. It computes the locations of quadrature
37  * points on the unit cell from a quadrature object for a manifold of one
38  * dimension less than that of the cell and the number of the face. For
39  * example, giving the Simpson rule in one dimension and using the
40  * project_to_face() function with face number 1, the returned points will be
41  * (1,0), (1,0.5) and (1,1). Note that faces have an orientation, so when
42  * projecting to face 3, you will get (0,0), (0,0.5) and (0,1), which is in
43  * clockwise sense, while for face 1 the points were in counterclockwise
44  * sense.
45  *
46  * For the projection to subfaces (i.e. to the children of a face of the unit
47  * cell), the same applies as above. Note the order in which the children of a
48  * face are numbered, which in two dimensions coincides with the orientation
49  * of the face.
50  *
51  * The second set of functions generates a quadrature formula by projecting a
52  * given quadrature rule on <b>all</b> faces and subfaces. This is used in the
53  * FEFaceValues and FESubfaceValues classes. Since we now have the quadrature
54  * points of all faces and subfaces in one array, we need to have a way to
55  * find the starting index of the points and weights corresponding to one face
56  * or subface within this array. This is done through the DataSetDescriptor
57  * member class.
58  *
59  * The different functions are grouped into a common class to avoid putting
60  * them into global namespace. However, since they have no local data, all
61  * functions are declared <tt>static</tt> and can be called without creating
62  * an object of this class.
63  *
64  * For the 3d case, you should note that the orientation of faces is even more
65  * intricate than for two dimensions. Quadrature formulae are projected upon
66  * the faces in their standard orientation, not to the inside or outside of
67  * the hexahedron. To make things more complicated, in 3d we allow faces in
68  * two orientations (which can be identified using
69  * <tt>cell->face_orientation(face)</tt>), so we have to project quadrature
70  * formula onto faces and subfaces in two orientations. (Refer to the
71  * documentation of the Triangulation class for a description of the
72  * orientation of the different faces, as well as to
73  * @ref GlossFaceOrientation "the glossary entry on face orientation"
74  * for more information on this.) The DataSetDescriptor member class is used
75  * to identify where each dataset starts.
76  */
77 template <int dim>
78 class QProjector
79 {
80 public:
81   /**
82    * Define an alias for a quadrature that acts on an object of one dimension
83    * less. For cells, this would then be a face quadrature.
84    */
85   using SubQuadrature = Quadrature<dim - 1>;
86 
87   /**
88    * Compute the quadrature points on the cell if the given quadrature formula
89    * is used on face <tt>face_no</tt>. For further details, see the general
90    * doc for this class.
91    *
92    * @deprecated This function makes an implicit assumption that the cell is
93    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
94    *   this function that takes the reference cell type instead.
95    */
96   DEAL_II_DEPRECATED static void
97   project_to_face(const SubQuadrature &    quadrature,
98                   const unsigned int       face_no,
99                   std::vector<Point<dim>> &q_points);
100 
101   /**
102    * Compute the quadrature points on the cell if the given quadrature formula
103    * is used on face <tt>face_no</tt>. For further details, see the general
104    * doc for this class.
105    */
106   static void
107   project_to_face(const ReferenceCell::Type reference_cell_type,
108                   const SubQuadrature &     quadrature,
109                   const unsigned int        face_no,
110                   std::vector<Point<dim>> & q_points);
111 
112   /**
113    * Compute the cell quadrature formula corresponding to using
114    * <tt>quadrature</tt> on face <tt>face_no</tt>. For further details, see
115    * the general doc for this class.
116    *
117    * @deprecated This function makes an implicit assumption that the cell is
118    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
119    *   this function that takes the reference cell type instead.
120    */
121   DEAL_II_DEPRECATED static Quadrature<dim>
122   project_to_face(const SubQuadrature &quadrature, const unsigned int face_no);
123 
124   /**
125    * Compute the cell quadrature formula corresponding to using
126    * <tt>quadrature</tt> on face <tt>face_no</tt>. For further details, see
127    * the general doc for this class.
128    */
129   static Quadrature<dim>
130   project_to_face(const ReferenceCell::Type reference_cell_type,
131                   const SubQuadrature &     quadrature,
132                   const unsigned int        face_no);
133 
134   /**
135    * Compute the quadrature points on the cell if the given quadrature formula
136    * is used on face <tt>face_no</tt>, subface number <tt>subface_no</tt>
137    * corresponding to RefineCase::Type <tt>ref_case</tt>. The last argument is
138    * only used in 3D.
139    *
140    * @note Only the points are transformed. The quadrature weights are the
141    * same as those of the original rule.
142    *
143    * @deprecated This function makes an implicit assumption that the cell is
144    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
145    *   this function that takes the reference cell type instead.
146    */
147   DEAL_II_DEPRECATED static void
148   project_to_subface(const SubQuadrature &          quadrature,
149                      const unsigned int             face_no,
150                      const unsigned int             subface_no,
151                      std::vector<Point<dim>> &      q_points,
152                      const RefinementCase<dim - 1> &ref_case =
153                        RefinementCase<dim - 1>::isotropic_refinement);
154 
155   /**
156    * Compute the quadrature points on the cell if the given quadrature formula
157    * is used on face <tt>face_no</tt>, subface number <tt>subface_no</tt>
158    * corresponding to RefineCase::Type <tt>ref_case</tt>. The last argument is
159    * only used in 3D.
160    *
161    * @note Only the points are transformed. The quadrature weights are the
162    * same as those of the original rule.
163    */
164   static void
165   project_to_subface(const ReferenceCell::Type      reference_cell_type,
166                      const SubQuadrature &          quadrature,
167                      const unsigned int             face_no,
168                      const unsigned int             subface_no,
169                      std::vector<Point<dim>> &      q_points,
170                      const RefinementCase<dim - 1> &ref_case =
171                        RefinementCase<dim - 1>::isotropic_refinement);
172 
173   /**
174    * Compute the cell quadrature formula corresponding to using
175    * <tt>quadrature</tt> on subface <tt>subface_no</tt> of face
176    * <tt>face_no</tt> with RefinementCase<dim-1> <tt>ref_case</tt>. The last
177    * argument is only used in 3D.
178    *
179    * @note Only the points are transformed. The quadrature weights are the
180    * same as those of the original rule.
181    *
182    * @deprecated This function makes an implicit assumption that the cell is
183    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
184    *   this function that takes the reference cell type instead.
185    */
186   DEAL_II_DEPRECATED static Quadrature<dim>
187   project_to_subface(const SubQuadrature &          quadrature,
188                      const unsigned int             face_no,
189                      const unsigned int             subface_no,
190                      const RefinementCase<dim - 1> &ref_case =
191                        RefinementCase<dim - 1>::isotropic_refinement);
192 
193   /**
194    * Compute the cell quadrature formula corresponding to using
195    * <tt>quadrature</tt> on subface <tt>subface_no</tt> of face
196    * <tt>face_no</tt> with RefinementCase<dim-1> <tt>ref_case</tt>. The last
197    * argument is only used in 3D.
198    *
199    * @note Only the points are transformed. The quadrature weights are the
200    * same as those of the original rule.
201    *
202    * @deprecated This function makes an implicit assumption that the cell is
203    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
204    *   this function that takes the reference cell type instead.
205    */
206   static Quadrature<dim>
207   project_to_subface(const ReferenceCell::Type      reference_cell_type,
208                      const SubQuadrature &          quadrature,
209                      const unsigned int             face_no,
210                      const unsigned int             subface_no,
211                      const RefinementCase<dim - 1> &ref_case =
212                        RefinementCase<dim - 1>::isotropic_refinement);
213 
214   /**
215    * Take a face quadrature formula and generate a cell quadrature formula
216    * from it where the quadrature points of the given argument are projected
217    * on all faces.
218    *
219    * The weights of the new rule are replications of the original weights.
220    * Thus, the sum of the weights is not one, but the number of faces, which
221    * is the surface of the reference cell.
222    *
223    * This in particular allows us to extract a subset of points corresponding
224    * to a single face and use it as a quadrature on this face, as is done in
225    * FEFaceValues.
226    *
227    * @note In 3D, this function produces eight sets of quadrature points for
228    * each face, in order to cope possibly different orientations of the mesh.
229    *
230    * @deprecated This function makes an implicit assumption that the cell is
231    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
232    *   this function that takes the reference cell type instead.
233    */
234   DEAL_II_DEPRECATED static Quadrature<dim>
235   project_to_all_faces(const SubQuadrature &quadrature);
236 
237   /**
238    * Take a face quadrature formula and generate a cell quadrature formula
239    * from it where the quadrature points of the given argument are projected
240    * on all faces.
241    *
242    * The weights of the new rule are replications of the original weights.
243    * Thus, the sum of the weights is not one, but the number of faces, which
244    * is the surface of the reference cell.
245    *
246    * This in particular allows us to extract a subset of points corresponding
247    * to a single face and use it as a quadrature on this face, as is done in
248    * FEFaceValues.
249    *
250    * @note In 3D, this function produces eight sets of quadrature points for
251    * each face, in order to cope possibly different orientations of the mesh.
252    */
253   static Quadrature<dim>
254   project_to_all_faces(const ReferenceCell::Type reference_cell_type,
255                        const SubQuadrature &     quadrature);
256 
257   /**
258    * Take a face quadrature formula and generate a cell quadrature formula
259    * from it where the quadrature points of the given argument are projected
260    * on all subfaces.
261    *
262    * Like in project_to_all_faces(), the weights of the new rule sum up to the
263    * number of faces (not subfaces), which is the surface of the reference
264    * cell.
265    *
266    * This in particular allows us to extract a subset of points corresponding
267    * to a single subface and use it as a quadrature on this face, as is done
268    * in FESubfaceValues.
269    *
270    * @deprecated This function makes an implicit assumption that the cell is
271    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
272    *   this function that takes the reference cell type instead.
273    */
274   DEAL_II_DEPRECATED static Quadrature<dim>
275   project_to_all_subfaces(const SubQuadrature &quadrature);
276 
277   /**
278    * Take a face quadrature formula and generate a cell quadrature formula
279    * from it where the quadrature points of the given argument are projected
280    * on all subfaces.
281    *
282    * Like in project_to_all_faces(), the weights of the new rule sum up to the
283    * number of faces (not subfaces), which is the surface of the reference
284    * cell.
285    *
286    * This in particular allows us to extract a subset of points corresponding
287    * to a single subface and use it as a quadrature on this face, as is done
288    * in FESubfaceValues.
289    */
290   static Quadrature<dim>
291   project_to_all_subfaces(const ReferenceCell::Type reference_cell_type,
292                           const SubQuadrature &     quadrature);
293 
294   /**
295    * Project a given quadrature formula to a child of a cell. You may want to
296    * use this function in case you want to extend an integral only over the
297    * area which a potential child would occupy. The child numbering is the
298    * same as the children would be numbered upon refinement of the cell.
299    *
300    * As integration using this quadrature formula now only extends over a
301    * fraction of the cell, the weights of the resulting object are divided by
302    * GeometryInfo<dim>::children_per_cell.
303    *
304    * @deprecated This function makes an implicit assumption that the cell is
305    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
306    *   this function that takes the reference cell type instead.
307    */
308   DEAL_II_DEPRECATED static Quadrature<dim>
309   project_to_child(const Quadrature<dim> &quadrature,
310                    const unsigned int     child_no);
311 
312   /**
313    * Project a given quadrature formula to a child of a cell. You may want to
314    * use this function in case you want to extend an integral only over the
315    * area which a potential child would occupy. The child numbering is the
316    * same as the children would be numbered upon refinement of the cell.
317    *
318    * As integration using this quadrature formula now only extends over a
319    * fraction of the cell, the weights of the resulting object are divided by
320    * GeometryInfo<dim>::children_per_cell.
321    */
322   static Quadrature<dim>
323   project_to_child(const ReferenceCell::Type reference_cell_type,
324                    const Quadrature<dim> &   quadrature,
325                    const unsigned int        child_no);
326 
327   /**
328    * Project a quadrature rule to all children of a cell. Similarly to
329    * project_to_all_subfaces(), this function replicates the formula generated
330    * by project_to_child() for all children, such that the weights sum up to
331    * one, the volume of the total cell again.
332    *
333    * The child numbering is the same as the children would be numbered upon
334    * refinement of the cell.
335    *
336    * @deprecated This function makes an implicit assumption that the cell is
337    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
338    *   this function that takes the reference cell type instead.
339    */
340   DEAL_II_DEPRECATED static Quadrature<dim>
341   project_to_all_children(const Quadrature<dim> &quadrature);
342 
343   /**
344    * Project a quadrature rule to all children of a cell. Similarly to
345    * project_to_all_subfaces(), this function replicates the formula generated
346    * by project_to_child() for all children, such that the weights sum up to
347    * one, the volume of the total cell again.
348    *
349    * The child numbering is the same as the children would be numbered upon
350    * refinement of the cell.
351    */
352   static Quadrature<dim>
353   project_to_all_children(const ReferenceCell::Type reference_cell_type,
354                           const Quadrature<dim> &   quadrature);
355 
356   /**
357    * Project the one dimensional rule <tt>quadrature</tt> to the straight line
358    * connecting the points <tt>p1</tt> and <tt>p2</tt>.
359    *
360    * @deprecated This function makes an implicit assumption that the cell is
361    *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
362    *   this function that takes the reference cell type instead.
363    */
364   DEAL_II_DEPRECATED static Quadrature<dim>
365   project_to_line(const Quadrature<1> &quadrature,
366                   const Point<dim> &   p1,
367                   const Point<dim> &   p2);
368 
369   /**
370    * Project the one dimensional rule <tt>quadrature</tt> to the straight line
371    * connecting the points <tt>p1</tt> and <tt>p2</tt>.
372    */
373   static Quadrature<dim>
374   project_to_line(const ReferenceCell::Type reference_cell_type,
375                   const Quadrature<1> &     quadrature,
376                   const Point<dim> &        p1,
377                   const Point<dim> &        p2);
378 
379   /**
380    * Since the project_to_all_faces() and project_to_all_subfaces() functions
381    * chain together the quadrature points and weights of all projections of a
382    * face quadrature formula to the faces or subfaces of a cell, we need a way
383    * to identify where the starting index of the points and weights for a
384    * particular face or subface is. This class provides this: there are static
385    * member functions that generate objects of this type, given face or
386    * subface indices, and you can then use the generated object in place of an
387    * integer that denotes the offset of a given dataset.
388    */
389   class DataSetDescriptor
390   {
391   public:
392     /**
393      * Default constructor. This doesn't do much except generating an invalid
394      * index, since you didn't give a valid descriptor of the cell, face, or
395      * subface you wanted.
396      */
397     DataSetDescriptor();
398 
399     /**
400      * Static function to generate the offset of a cell. Since we only have
401      * one cell per quadrature object, this offset is of course zero, but we
402      * carry this function around for consistency with the other static
403      * functions.
404      */
405     static DataSetDescriptor
406     cell();
407 
408     /**
409      * Static function to generate an offset object for a given face of a cell
410      * with the given face orientation, flip and rotation. This function of
411      * course is only allowed if <tt>dim>=2</tt>, and the face orientation,
412      * flip and rotation are ignored if the space dimension equals 2.
413      *
414      * The last argument denotes the number of quadrature points the lower-
415      * dimensional face quadrature formula (the one that has been projected
416      * onto the faces) has.
417      *
418      * @deprecated This function makes an implicit assumption that the cell is
419      *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
420      *   this function that takes the reference cell type instead.
421      */
422     DEAL_II_DEPRECATED static DataSetDescriptor
423     face(const unsigned int face_no,
424          const bool         face_orientation,
425          const bool         face_flip,
426          const bool         face_rotation,
427          const unsigned int n_quadrature_points);
428 
429     /**
430      * Static function to generate an offset object for a given face of a cell
431      * with the given face orientation, flip and rotation. This function of
432      * course is only allowed if <tt>dim>=2</tt>, and the face orientation,
433      * flip and rotation are ignored if the space dimension equals 2.
434      *
435      * The last argument denotes the number of quadrature points the lower-
436      * dimensional face quadrature formula (the one that has been projected
437      * onto the faces) has.
438      */
439     static DataSetDescriptor
440     face(const ReferenceCell::Type reference_cell_type,
441          const unsigned int        face_no,
442          const bool                face_orientation,
443          const bool                face_flip,
444          const bool                face_rotation,
445          const unsigned int        n_quadrature_points);
446 
447     /**
448      * Static function to generate an offset object for a given subface of a
449      * cell with the given face orientation, flip and rotation. This function
450      * of course is only allowed if <tt>dim>=2</tt>, and the face orientation,
451      * flip and rotation are ignored if the space dimension equals 2.
452      *
453      * The last but one argument denotes the number of quadrature points the
454      * lower-dimensional face quadrature formula (the one that has been
455      * projected onto the faces) has.
456      *
457      * Through the last argument anisotropic refinement can be respected.
458      *
459      * @deprecated This function makes an implicit assumption that the cell is
460      *   a line (1D), a quad (2D), or a hex (3D). Use the other version of
461      *   this function that takes the reference cell type instead.
462      */
463     DEAL_II_DEPRECATED static DataSetDescriptor
464     subface(const unsigned int               face_no,
465             const unsigned int               subface_no,
466             const bool                       face_orientation,
467             const bool                       face_flip,
468             const bool                       face_rotation,
469             const unsigned int               n_quadrature_points,
470             const internal::SubfaceCase<dim> ref_case =
471               internal::SubfaceCase<dim>::case_isotropic);
472 
473     /**
474      * Static function to generate an offset object for a given subface of a
475      * cell with the given face orientation, flip and rotation. This function
476      * of course is only allowed if <tt>dim>=2</tt>, and the face orientation,
477      * flip and rotation are ignored if the space dimension equals 2.
478      *
479      * The last but one argument denotes the number of quadrature points the
480      * lower-dimensional face quadrature formula (the one that has been
481      * projected onto the faces) has.
482      *
483      * Through the last argument anisotropic refinement can be respected.
484      */
485     static DataSetDescriptor
486     subface(const ReferenceCell::Type        reference_cell_type,
487             const unsigned int               face_no,
488             const unsigned int               subface_no,
489             const bool                       face_orientation,
490             const bool                       face_flip,
491             const bool                       face_rotation,
492             const unsigned int               n_quadrature_points,
493             const internal::SubfaceCase<dim> ref_case =
494               internal::SubfaceCase<dim>::case_isotropic);
495 
496     /**
497      * Conversion operator to an integer denoting the offset of the first
498      * element of this dataset in the set of quadrature formulas all projected
499      * onto faces and subfaces. This conversion operator allows us to use
500      * offset descriptor objects in place of integer offsets.
501      */
502     operator unsigned int() const;
503 
504   private:
505     /**
506      * Store the integer offset for a given cell, face, or subface.
507      */
508     const unsigned int dataset_offset;
509 
510     /**
511      * This is the real constructor, but it is private and thus only available
512      * to the static member functions above.
513      */
514     DataSetDescriptor(const unsigned int dataset_offset);
515   };
516 
517 private:
518   /**
519    * Given a quadrature object in 2d, reflect all quadrature points at the
520    * main diagonal and return them with their original weights.
521    *
522    * This function is necessary for projecting a 2d quadrature rule onto the
523    * faces of a 3d cube, since there we need both orientations.
524    */
525   static Quadrature<2>
526   reflect(const Quadrature<2> &q);
527 
528   /**
529    * Given a quadrature object in 2d, rotate all quadrature points by @p
530    * n_times * 90 degrees counterclockwise and return them with their original
531    * weights.
532    *
533    * This function is necessary for projecting a 2d quadrature rule onto the
534    * faces of a 3d cube, since there we need all rotations to account for
535    * face_flip and face_rotation of non-standard faces.
536    */
537   static Quadrature<2>
538   rotate(const Quadrature<2> &q, const unsigned int n_times);
539 };
540 
541 /*@}*/
542 
543 
544 // -------------------  inline and template functions ----------------
545 
546 
547 
548 template <int dim>
DataSetDescriptor(const unsigned int dataset_offset)549 inline QProjector<dim>::DataSetDescriptor::DataSetDescriptor(
550   const unsigned int dataset_offset)
551   : dataset_offset(dataset_offset)
552 {}
553 
554 
555 template <int dim>
DataSetDescriptor()556 inline QProjector<dim>::DataSetDescriptor::DataSetDescriptor()
557   : dataset_offset(numbers::invalid_unsigned_int)
558 {}
559 
560 
561 
562 template <int dim>
563 typename QProjector<dim>::DataSetDescriptor
cell()564 QProjector<dim>::DataSetDescriptor::cell()
565 {
566   return 0;
567 }
568 
569 
570 
571 template <int dim>
572 inline QProjector<dim>::DataSetDescriptor::operator unsigned int() const
573 {
574   return dataset_offset;
575 }
576 
577 
578 /* -------------- declaration of explicit specializations ------------- */
579 
580 #ifndef DOXYGEN
581 
582 
583 template <>
584 void
585 QProjector<1>::project_to_face(const Quadrature<0> &,
586                                const unsigned int,
587                                std::vector<Point<1>> &);
588 template <>
589 void
590 QProjector<1>::project_to_face(const ReferenceCell::Type reference_cell_type,
591                                const Quadrature<0> &,
592                                const unsigned int,
593                                std::vector<Point<1>> &);
594 template <>
595 void
596 QProjector<2>::project_to_face(const Quadrature<1> &  quadrature,
597                                const unsigned int     face_no,
598                                std::vector<Point<2>> &q_points);
599 template <>
600 void
601 QProjector<2>::project_to_face(const ReferenceCell::Type reference_cell_type,
602                                const Quadrature<1> &     quadrature,
603                                const unsigned int        face_no,
604                                std::vector<Point<2>> &   q_points);
605 template <>
606 void
607 QProjector<3>::project_to_face(const Quadrature<2> &  quadrature,
608                                const unsigned int     face_no,
609                                std::vector<Point<3>> &q_points);
610 template <>
611 void
612 QProjector<3>::project_to_face(const ReferenceCell::Type reference_cell_type,
613                                const Quadrature<2> &     quadrature,
614                                const unsigned int        face_no,
615                                std::vector<Point<3>> &   q_points);
616 
617 template <>
618 Quadrature<1>
619 QProjector<1>::project_to_all_faces(const Quadrature<0> &quadrature);
620 template <>
621 Quadrature<1>
622 QProjector<1>::project_to_all_faces(
623   const ReferenceCell::Type reference_cell_type,
624   const Quadrature<0> &     quadrature);
625 
626 
627 template <>
628 void
629 QProjector<1>::project_to_subface(const Quadrature<0> &,
630                                   const unsigned int,
631                                   const unsigned int,
632                                   std::vector<Point<1>> &,
633                                   const RefinementCase<0> &);
634 template <>
635 void
636 QProjector<1>::project_to_subface(const ReferenceCell::Type reference_cell_type,
637                                   const Quadrature<0> &,
638                                   const unsigned int,
639                                   const unsigned int,
640                                   std::vector<Point<1>> &,
641                                   const RefinementCase<0> &);
642 template <>
643 void
644 QProjector<2>::project_to_subface(const Quadrature<1> &  quadrature,
645                                   const unsigned int     face_no,
646                                   const unsigned int     subface_no,
647                                   std::vector<Point<2>> &q_points,
648                                   const RefinementCase<1> &);
649 template <>
650 void
651 QProjector<2>::project_to_subface(const ReferenceCell::Type reference_cell_type,
652                                   const Quadrature<1> &     quadrature,
653                                   const unsigned int        face_no,
654                                   const unsigned int        subface_no,
655                                   std::vector<Point<2>> &   q_points,
656                                   const RefinementCase<1> &);
657 template <>
658 void
659 QProjector<3>::project_to_subface(const Quadrature<2> &    quadrature,
660                                   const unsigned int       face_no,
661                                   const unsigned int       subface_no,
662                                   std::vector<Point<3>> &  q_points,
663                                   const RefinementCase<2> &face_ref_case);
664 template <>
665 void
666 QProjector<3>::project_to_subface(const ReferenceCell::Type reference_cell_type,
667                                   const Quadrature<2> &     quadrature,
668                                   const unsigned int        face_no,
669                                   const unsigned int        subface_no,
670                                   std::vector<Point<3>> &   q_points,
671                                   const RefinementCase<2> & face_ref_case);
672 
673 template <>
674 Quadrature<1>
675 QProjector<1>::project_to_all_subfaces(const Quadrature<0> &quadrature);
676 template <>
677 Quadrature<1>
678 QProjector<1>::project_to_all_subfaces(
679   const ReferenceCell::Type reference_cell_type,
680   const Quadrature<0> &     quadrature);
681 
682 
683 #endif // DOXYGEN
684 DEAL_II_NAMESPACE_CLOSE
685 
686 #endif
687