Reference documentation for deal.II version 8.4.2
quadrature_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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 protected:
92  std::vector<long double>
93  compute_quadrature_points (const unsigned int q,
94  const int alpha,
95  const int beta) const;
96 
104  std::vector<long double>
105  compute_quadrature_weights (const std::vector<long double> &x,
106  const int alpha,
107  const int beta) const;
108 
115  long double JacobiP(const long double x,
116  const int alpha,
117  const int beta,
118  const unsigned int n) const;
119 
124  long double gamma(const unsigned int n) const;
125 };
126 
127 
128 
133 template <int dim>
134 class QMidpoint : public Quadrature<dim>
135 {
136 public:
137  QMidpoint ();
138 };
139 
140 
145 template <int dim>
146 class QSimpson : public Quadrature<dim>
147 {
148 public:
149  QSimpson ();
150 };
151 
152 
153 
166 template <int dim>
167 class QTrapez : public Quadrature<dim>
168 {
169 public:
170  QTrapez ();
171 };
172 
173 
174 
181 template <int dim>
182 class QMilne : public Quadrature<dim>
183 {
184 public:
185  QMilne ();
186 };
187 
188 
195 template <int dim>
196 class QWeddle : public Quadrature<dim>
197 {
198 public:
199  QWeddle ();
200 };
201 
202 
203 
217 template <int dim>
218 class QGaussLog : public Quadrature<dim>
219 {
220 public:
224  QGaussLog(const unsigned int n,
225  const bool revert=false);
226 
227 protected:
231  std::vector<double>
232  set_quadrature_points(const unsigned int n) const;
233 
237  std::vector<double>
238  set_quadrature_weights(const unsigned int n) const;
239 
240 };
241 
242 
243 
244 
281 template <int dim>
282 class QGaussLogR : public Quadrature<dim>
283 {
284 public:
292  QGaussLogR(const unsigned int n,
293  const Point<dim> x0 = Point<dim>(),
294  const double alpha = 1,
295  const bool factor_out_singular_weight=false);
296 
297 protected:
302  const double fraction;
303 };
304 
305 
329 template <int dim>
330 class QGaussOneOverR : public Quadrature<dim>
331 {
332 public:
365  QGaussOneOverR(const unsigned int n,
366  const Point<dim> singularity,
367  const bool factor_out_singular_weight=false);
402  QGaussOneOverR(const unsigned int n,
403  const unsigned int vertex_index,
404  const bool factor_out_singular_weight=false);
405 private:
411  static unsigned int quad_size(const Point<dim> singularity,
412  const unsigned int n);
413 };
414 
415 
416 
426 template <int dim>
427 class QSorted : public Quadrature<dim>
428 {
429 public:
433  QSorted (const Quadrature<dim>);
434 
438  bool operator()(const std::pair<double, Point<dim> > &a,
439  const std::pair<double, Point<dim> > &b);
440 };
441 
498 template <int dim>
499 class QTelles: public Quadrature<dim>
500 {
501 public:
508  QTelles (const Quadrature<1> &base_quad, const Point<dim> &singularity);
514  QTelles (const unsigned int n, const Point<dim> &singularity);
515 
516 };
517 
533 template <int dim>
534 class QGaussChebyshev : public Quadrature<dim>
535 {
536 public:
538  QGaussChebyshev(const unsigned int n);
539 
540 private:
542  static std::vector<double>
543  get_quadrature_points(const unsigned int n);
544 
546  static std::vector<double>
547  get_quadrature_weights(const unsigned int n);
548 
549 };
550 
551 
568 template <int dim>
569 class QGaussRadauChebyshev : public Quadrature<dim>
570 {
571 public:
572  /* EndPoint is used to specify which of the two endpoints of the unit interval
573  * is used also as quadrature point
574  */
575  enum EndPoint { left,right };
577  QGaussRadauChebyshev(const unsigned int n,
578  EndPoint ep=QGaussRadauChebyshev::left);
579 
580 private:
581  const EndPoint ep;
583  static std::vector<double>
584  get_quadrature_points(const unsigned int n, EndPoint ep);
585 
587  static std::vector<double>
588  get_quadrature_weights(const unsigned int n, EndPoint ep);
589 
590 };
591 
607 template <int dim>
609 {
610 public:
612  QGaussLobattoChebyshev(const unsigned int n);
613 
614 private:
616  static std::vector<double>
617  get_quadrature_points(const unsigned int n);
618 
620  static std::vector<double>
621  get_quadrature_weights(const unsigned int n);
622 
623 };
624 
625 /* -------------- declaration of explicit specializations ------------- */
626 
627 template <> QGauss<1>::QGauss (const unsigned int n);
628 template <> QGaussLobatto<1>::QGaussLobatto (const unsigned int n);
629 template <>
630 std::vector<long double> QGaussLobatto<1>::
631 compute_quadrature_points(const unsigned int, const int, const int) const;
632 template <>
633 std::vector<long double> QGaussLobatto<1>::
634 compute_quadrature_weights(const std::vector<long double> &, const int, const int) const;
635 template <>
636 long double QGaussLobatto<1>::
637 JacobiP(const long double, const int, const int, const unsigned int) const;
638 template <>
639 long double
640 QGaussLobatto<1>::gamma(const unsigned int n) const;
641 
642 template <> std::vector<double> QGaussLog<1>::set_quadrature_points(const unsigned int) const;
643 template <> std::vector<double> QGaussLog<1>::set_quadrature_weights(const unsigned int) const;
644 
645 template <> QMidpoint<1>::QMidpoint ();
646 template <> QTrapez<1>::QTrapez ();
647 template <> QSimpson<1>::QSimpson ();
648 template <> QMilne<1>::QMilne ();
649 template <> QWeddle<1>::QWeddle ();
650 template <> QGaussLog<1>::QGaussLog (const unsigned int n, const bool revert);
651 template <> QGaussLogR<1>::QGaussLogR (const unsigned int n, const Point<1> x0, const double alpha, const bool flag);
652 template <> QGaussOneOverR<2>::QGaussOneOverR (const unsigned int n, const unsigned int index, const bool flag);
653 template <> QTelles<1>::QTelles(const Quadrature<1> &base_quad, const Point<1> &singularity);
654 
655 
656 
657 DEAL_II_NAMESPACE_CLOSE
658 
659 #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)
long double gamma(const unsigned int n) const
std::vector< long double > compute_quadrature_points(const unsigned int q, const int alpha, const int beta) const
const double fraction
QGauss(const unsigned int n)
std::vector< double > set_quadrature_points(const unsigned int n) const
std::vector< long double > compute_quadrature_weights(const std::vector< long double > &x, const int alpha, const int beta) const
QGaussLobatto(const unsigned int n)
std::vector< double > set_quadrature_weights(const unsigned int n) const
long double JacobiP(const long double x, const int alpha, const int beta, const unsigned int n) const
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
QGaussLogR(const unsigned int n, const Point< dim > x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)