Reference documentation for deal.II version 8.1.0
qprojector.h
1 // ---------------------------------------------------------------------
2 // @f$Id: qprojector.h 30036 2013-07-18 16:55:32Z maier @f$
3 //
4 // Copyright (C) 2005 - 2013 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef __deal2__qprojector_h
18 #define __deal2__qprojector_h
19 
20 
21 #include <deal.II/base/quadrature.h>
22 #include <deal.II/base/geometry_info.h>
23 
27 
28 
77 template <int dim>
79 {
80 public:
88  typedef Quadrature<dim-1> SubQuadrature;
89 
98  static void project_to_face (const SubQuadrature &quadrature,
99  const unsigned int face_no,
100  std::vector<Point<dim> > &q_points);
101 
110  static Quadrature<dim>
111  project_to_face (const SubQuadrature &quadrature,
112  const unsigned int face_no);
113 
128  static void project_to_subface (const SubQuadrature &quadrature,
129  const unsigned int face_no,
130  const unsigned int subface_no,
131  std::vector<Point<dim> > &q_points,
133 
149  static Quadrature<dim>
150  project_to_subface (const SubQuadrature &quadrature,
151  const unsigned int face_no,
152  const unsigned int subface_no,
154 
185  static Quadrature<dim>
186  project_to_all_faces (const SubQuadrature &quadrature);
187 
210  static Quadrature<dim>
211  project_to_all_subfaces (const SubQuadrature &quadrature);
212 
232  static
234  project_to_child (const Quadrature<dim> &quadrature,
235  const unsigned int child_no);
236 
255  static
257  project_to_all_children (const Quadrature<dim> &quadrature);
258 
266  static
268  project_to_line(const Quadrature<1> &quadrature,
269  const Point<dim> &p1,
270  const Point<dim> &p2);
271 
298  {
299  public:
310 
322  static DataSetDescriptor cell ();
323 
344  static
346  face (const unsigned int face_no,
347  const bool face_orientation,
348  const bool face_flip,
349  const bool face_rotation,
350  const unsigned int n_quadrature_points);
351 
376  static
378  subface (const unsigned int face_no,
379  const unsigned int subface_no,
380  const bool face_orientation,
381  const bool face_flip,
382  const bool face_rotation,
383  const unsigned int n_quadrature_points,
385 
400  operator unsigned int () const;
401 
402  private:
408  const unsigned int dataset_offset;
409 
417  DataSetDescriptor (const unsigned int dataset_offset);
418  };
419 
420 private:
434  static Quadrature<2> reflect (const Quadrature<2> &q);
435 
452  static Quadrature<2> rotate (const Quadrature<2> &q,
453  const unsigned int n_times);
454 };
455 
459 // ------------------- inline and template functions ----------------
460 
461 
462 
463 template <int dim>
464 inline
466 DataSetDescriptor (const unsigned int dataset_offset)
467  :
468  dataset_offset (dataset_offset)
469 {}
470 
471 
472 template <int dim>
473 inline
476  :
477  dataset_offset (numbers::invalid_unsigned_int)
478 {}
479 
480 
481 
482 template <int dim>
485 {
486  return 0;
487 }
488 
489 
490 
491 template <int dim>
492 inline
494 {
495  return dataset_offset;
496 }
497 
498 
499 /* -------------- declaration of explicit specializations ------------- */
500 
501 #ifndef DOXYGEN
502 
503 
504 template <>
505 void
507  const unsigned int,
508  std::vector<Point<1> > &);
509 template <>
510 void
512  const unsigned int face_no,
513  std::vector<Point<2> > &q_points);
514 template <>
515 void
517  const unsigned int face_no,
518  std::vector<Point<3> > &q_points);
519 
520 template <>
523 
524 
525 template <>
526 void
528  const unsigned int,
529  const unsigned int,
530  std::vector<Point<1> > &,
531  const RefinementCase<0> &);
532 template <>
533 void
535  const unsigned int face_no,
536  const unsigned int subface_no,
537  std::vector<Point<2> > &q_points,
538  const RefinementCase<1> &);
539 template <>
540 void
542  const unsigned int face_no,
543  const unsigned int subface_no,
544  std::vector<Point<3> > &q_points,
545  const RefinementCase<2> &face_ref_case);
546 
547 template <>
550 
551 
552 template <>
553 bool
554 QIterated<1>::uses_both_endpoints (const Quadrature<1> &base_quadrature);
555 
556 template <>
557 QIterated<1>::QIterated (const Quadrature<1> &base_quadrature,
558  const unsigned int n_copies);
559 
560 
561 #endif // DOXYGEN
562 DEAL_II_NAMESPACE_CLOSE
563 
564 #endif
const unsigned int dataset_offset
Definition: qprojector.h:408
static bool uses_both_endpoints(const Quadrature< 1 > &base_quadrature)
Quadrature< dim-1 > SubQuadrature
Definition: qprojector.h:88
static DataSetDescriptor cell()
Definition: qprojector.h:484
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
static Quadrature< dim > project_to_line(const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
static Quadrature< 2 > rotate(const Quadrature< 2 > &q, const unsigned int n_times)
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
static Quadrature< 2 > reflect(const Quadrature< 2 > &q)
static Quadrature< dim > project_to_all_children(const Quadrature< dim > &quadrature)
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: types.h:173
static DataSetDescriptor subface(const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)