Reference documentation for deal.II version 9.0.0
quadrature_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 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 #ifndef dealii_quadrature_lib_h
17 #define dealii_quadrature_lib_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/quadrature.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
27 
38 template <int dim>
39 class QGauss : public Quadrature<dim>
40 {
41 public:
46  QGauss (const unsigned int n);
47 };
48 
49 
74 template <int dim>
75 class QGaussLobatto : public Quadrature<dim>
76 {
77 public:
82  QGaussLobatto(const unsigned int n);
83 };
84 
85 
86 
91 template <int dim>
92 class QMidpoint : public Quadrature<dim>
93 {
94 public:
95  QMidpoint ();
96 };
97 
98 
103 template <int dim>
104 class QSimpson : public Quadrature<dim>
105 {
106 public:
107  QSimpson ();
108 };
109 
110 
111 
124 template <int dim>
125 class QTrapez : public Quadrature<dim>
126 {
127 public:
128  QTrapez ();
129 };
130 
131 
132 
139 template <int dim>
140 class QMilne : public Quadrature<dim>
141 {
142 public:
143  QMilne ();
144 };
145 
146 
153 template <int dim>
154 class QWeddle : public Quadrature<dim>
155 {
156 public:
157  QWeddle ();
158 };
159 
160 
161 
175 template <int dim>
176 class QGaussLog : public Quadrature<dim>
177 {
178 public:
182  QGaussLog(const unsigned int n,
183  const bool revert=false);
184 
185 private:
189  static
190  std::vector<double>
191  get_quadrature_points(const unsigned int n);
192 
196  static
197  std::vector<double>
198  get_quadrature_weights(const unsigned int n);
199 };
200 
201 
202 
203 
240 template <int dim>
241 class QGaussLogR : public Quadrature<dim>
242 {
243 public:
251  QGaussLogR(const unsigned int n,
252  const Point<dim> x0 = Point<dim>(),
253  const double alpha = 1,
254  const bool factor_out_singular_weight=false);
255 
260  QGaussLogR(QGaussLogR<dim> &&) noexcept = default;
261 
262 protected:
267  const double fraction;
268 };
269 
270 
294 template <int dim>
295 class QGaussOneOverR : public Quadrature<dim>
296 {
297 public:
330  QGaussOneOverR(const unsigned int n,
331  const Point<dim> singularity,
332  const bool factor_out_singular_weight=false);
367  QGaussOneOverR(const unsigned int n,
368  const unsigned int vertex_index,
369  const bool factor_out_singular_weight=false);
370 private:
376  static unsigned int quad_size(const Point<dim> singularity,
377  const unsigned int n);
378 };
379 
380 
381 
391 template <int dim>
392 class QSorted : public Quadrature<dim>
393 {
394 public:
399  QSorted (const Quadrature<dim> &quad);
400 
401 private:
407  bool compare_weights(const unsigned int a,
408  const unsigned int b) const;
409 };
410 
467 template <int dim>
468 class QTelles: public Quadrature<dim>
469 {
470 public:
477  QTelles (const Quadrature<1> &base_quad, const Point<dim> &singularity);
483  QTelles (const unsigned int n, const Point<dim> &singularity);
484 
485 };
486 
500 template <int dim>
501 class QGaussChebyshev : public Quadrature<dim>
502 {
503 public:
505  QGaussChebyshev(const unsigned int n);
506 };
507 
508 
525 template <int dim>
526 class QGaussRadauChebyshev : public Quadrature<dim>
527 {
528 public:
529  /* EndPoint is used to specify which of the two endpoints of the unit interval
530  * is used also as quadrature point
531  */
532  enum EndPoint
533  {
542  };
544  QGaussRadauChebyshev(const unsigned int n,
546 
551  QGaussRadauChebyshev(QGaussRadauChebyshev<dim> &&) noexcept = default;
552 
553 private:
554  const EndPoint ep;
555 };
556 
572 template <int dim>
574 {
575 public:
577  QGaussLobattoChebyshev(const unsigned int n);
578 };
579 
612 template <int dim>
613 class QSimplex : public Quadrature<dim>
614 {
615 public:
622  QSimplex(const Quadrature<dim> &quad);
623 
647  compute_affine_transformation(const std::array<Point<dim>, dim+1> &vertices) const;
648 
649 };
650 
671 class QTrianglePolar: public QSimplex<2>
672 {
673 public:
681  QTrianglePolar(const Quadrature<1> &radial_quadrature,
682  const Quadrature<1> &angular_quadrature);
683 
690  QTrianglePolar(const unsigned int &n);
691 };
692 
727 class QDuffy: public QSimplex<2>
728 {
729 public:
743  QDuffy(const Quadrature<1> &radial_quadrature,
744  const Quadrature<1> &angular_quadrature,
745  const double beta = 1.0);
746 
753  QDuffy(const unsigned int n,
754  const double beta);
755 };
756 
763 template<int dim>
764 class QSplit : public Quadrature<dim>
765 {
766 public:
801  QSplit(const QSimplex<dim> &base,
802  const Point<dim> &split_point);
803 };
804 
807 /* -------------- declaration of explicit specializations ------------- */
808 
809 template <> QGauss<1>::QGauss (const unsigned int n);
810 template <> QGaussLobatto<1>::QGaussLobatto (const unsigned int n);
811 
812 template <> std::vector<double> QGaussLog<1>::get_quadrature_points(const unsigned int);
813 template <> std::vector<double> QGaussLog<1>::get_quadrature_weights(const unsigned int);
814 
815 template <> QMidpoint<1>::QMidpoint ();
816 template <> QTrapez<1>::QTrapez ();
817 template <> QSimpson<1>::QSimpson ();
818 template <> QMilne<1>::QMilne ();
819 template <> QWeddle<1>::QWeddle ();
820 template <> QGaussLog<1>::QGaussLog (const unsigned int n, const bool revert);
821 template <> QGaussLogR<1>::QGaussLogR (const unsigned int n, const Point<1> x0, const double alpha, const bool flag);
822 template <> QGaussOneOverR<2>::QGaussOneOverR (const unsigned int n, const unsigned int index, const bool flag);
823 template <> QTelles<1>::QTelles(const Quadrature<1> &base_quad, const Point<1> &singularity);
824 
825 
826 
827 DEAL_II_NAMESPACE_CLOSE
828 
829 #endif
QGaussLog(const unsigned int n, const bool revert=false)
QGaussOneOverR(const unsigned int n, const Point< dim > singularity, const bool factor_out_singular_weight=false)
const double fraction
QGaussChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
static std::vector< double > get_quadrature_points(const unsigned int n)
QGauss(const unsigned int n)
static std::vector< double > get_quadrature_weights(const unsigned int n)
QGaussRadauChebyshev(const unsigned int n, EndPoint ep=QGaussRadauChebyshev::left)
Generate a formula with n quadrature points.
bool compare_weights(const unsigned int a, const unsigned int b) const
QGaussLobatto(const unsigned int n)
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
QSorted(const Quadrature< dim > &quad)
QGaussLogR(const unsigned int n, const Point< dim > x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)
static unsigned int quad_size(const Point< dim > singularity, const unsigned int n)