Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
quadrature_lib.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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_quadrature_lib_h
17 #define dealii_quadrature_lib_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
25 
28 
39 template <int dim>
40 class QGauss : public Quadrature<dim>
41 {
42 public:
47  QGauss(const unsigned int n);
48 };
49 
50 
75 template <int dim>
76 class QGaussLobatto : public Quadrature<dim>
77 {
78 public:
83  QGaussLobatto(const unsigned int n);
84 };
85 
86 
87 
92 template <int dim>
93 class QMidpoint : public Quadrature<dim>
94 {
95 public:
96  QMidpoint();
97 };
98 
99 
104 template <int dim>
105 class QSimpson : public Quadrature<dim>
106 {
107 public:
108  QSimpson();
109 };
110 
111 
112 
125 template <int dim>
126 class QTrapez : public Quadrature<dim>
127 {
128 public:
129  QTrapez();
130 };
131 
132 
133 
140 template <int dim>
141 class QMilne : public Quadrature<dim>
142 {
143 public:
144  QMilne();
145 };
146 
147 
154 template <int dim>
155 class QWeddle : public Quadrature<dim>
156 {
157 public:
158  QWeddle();
159 };
160 
161 
162 
176 template <int dim>
177 class QGaussLog : public Quadrature<dim>
178 {
179 public:
183  QGaussLog(const unsigned int n, const bool revert = false);
184 
185 private:
189  static std::vector<double>
190  get_quadrature_points(const unsigned int n);
191 
195  static std::vector<double>
196  get_quadrature_weights(const unsigned int n);
197 };
198 
199 
200 
237 template <int dim>
238 class QGaussLogR : public Quadrature<dim>
239 {
240 public:
248  QGaussLogR(const unsigned int n,
249  const Point<dim> x0 = Point<dim>(),
250  const double alpha = 1,
251  const bool factor_out_singular_weight = false);
252 
258  QGaussLogR(QGaussLogR<dim> &&) noexcept = default;
259 
260 protected:
265  const double fraction;
266 };
267 
268 
292 template <int dim>
293 class QGaussOneOverR : public Quadrature<dim>
294 {
295 public:
328  QGaussOneOverR(const unsigned int n,
329  const Point<dim> singularity,
330  const bool factor_out_singular_weight = false);
365  QGaussOneOverR(const unsigned int n,
366  const unsigned int vertex_index,
367  const bool factor_out_singular_weight = false);
368 
369 private:
375  static unsigned int
376  quad_size(const Point<dim> singularity, const unsigned int n);
377 };
378 
379 
380 
390 template <int dim>
391 class QSorted : public Quadrature<dim>
392 {
393 public:
398  QSorted(const Quadrature<dim> &quad);
399 
400 private:
406  bool
407  compare_weights(const unsigned int a, const unsigned int b) const;
408 };
409 
466 template <int dim>
467 class QTelles : public Quadrature<dim>
468 {
469 public:
476  QTelles(const Quadrature<1> &base_quad, const Point<dim> &singularity);
482  QTelles(const unsigned int n, const Point<dim> &singularity);
483 };
484 
498 template <int dim>
499 class QGaussChebyshev : public Quadrature<dim>
500 {
501 public:
503  QGaussChebyshev(const unsigned int n);
504 };
505 
506 
523 template <int dim>
524 class QGaussRadauChebyshev : public Quadrature<dim>
525 {
526 public:
527  /* EndPoint is used to specify which of the two endpoints of the unit interval
528  * is used also as quadrature point
529  */
530  enum EndPoint
531  {
540  };
542  QGaussRadauChebyshev(const unsigned int n,
544 
549  QGaussRadauChebyshev(QGaussRadauChebyshev<dim> &&) noexcept = default;
550 
551 private:
552  const EndPoint ep;
553 };
554 
570 template <int dim>
572 {
573 public:
575  QGaussLobattoChebyshev(const unsigned int n);
576 };
577 
610 template <int dim>
611 class QSimplex : public Quadrature<dim>
612 {
613 public:
620  QSimplex(const Quadrature<dim> &quad);
621 
646  compute_affine_transformation(
647  const std::array<Point<dim>, dim + 1> &vertices) const;
648 };
649 
670 class QTrianglePolar : public QSimplex<2>
671 {
672 public:
680  QTrianglePolar(const Quadrature<1> &radial_quadrature,
681  const Quadrature<1> &angular_quadrature);
682 
689  QTrianglePolar(const unsigned int n);
690 };
691 
726 class QDuffy : public QSimplex<2>
727 {
728 public:
743  QDuffy(const Quadrature<1> &radial_quadrature,
744  const Quadrature<1> &angular_quadrature,
745  const double beta = 1.0);
746 
754  QDuffy(const unsigned int n, const double beta);
755 };
756 
763 template <int dim>
764 class QSplit : public Quadrature<dim>
765 {
766 public:
801  QSplit(const QSimplex<dim> &base, const Point<dim> &split_point);
802 };
803 
806 /* -------------- declaration of explicit specializations ------------- */
807 
808 #ifndef DOXYGEN
809 template <>
810 QGauss<1>::QGauss(const unsigned int n);
811 template <>
812 QGaussLobatto<1>::QGaussLobatto(const unsigned int n);
813 
814 template <>
815 std::vector<double>
816 QGaussLog<1>::get_quadrature_points(const unsigned int);
817 template <>
818 std::vector<double>
819 QGaussLog<1>::get_quadrature_weights(const unsigned int);
820 
821 template <>
823 template <>
825 template <>
827 template <>
829 template <>
831 template <>
832 QGaussLog<1>::QGaussLog(const unsigned int n, const bool revert);
833 template <>
834 QGaussLogR<1>::QGaussLogR(const unsigned int n,
835  const Point<1> x0,
836  const double alpha,
837  const bool flag);
838 template <>
839 QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
840  const unsigned int index,
841  const bool flag);
842 template <>
843 QTelles<1>::QTelles(const Quadrature<1> &base_quad,
844  const Point<1> & singularity);
845 #endif // DOXYGEN
846 
847 
848 
850 #endif
QSimpson
Definition: quadrature_lib.h:105
QGaussRadauChebyshev::left
@ left
Definition: quadrature_lib.h:535
QGaussRadauChebyshev::QGaussRadauChebyshev
QGaussRadauChebyshev(const unsigned int n, EndPoint ep=QGaussRadauChebyshev::left)
Generate a formula with n quadrature points.
Definition: quadrature_lib.cc:1118
QSorted
Definition: quadrature_lib.h:391
QMilne::QMilne
QMilne()
Definition: quadrature_lib.cc:846
QGaussOneOverR
Definition: quadrature_lib.h:293
QGaussOneOverR::quad_size
static unsigned int quad_size(const Point< dim > singularity, const unsigned int n)
QGaussLobatto::QGaussLobatto
QGaussLobatto(const unsigned int n)
Definition: quadrature_lib.cc:818
QGaussLobattoChebyshev
Definition: quadrature_lib.h:571
QGaussLog
Definition: quadrature_lib.h:177
QDuffy
Definition: quadrature_lib.h:726
QGaussLog::get_quadrature_points
static std::vector< double > get_quadrature_points(const unsigned int n)
QTelles::QTelles
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
Definition: quadrature_lib.cc:857
QWeddle::QWeddle
QWeddle()
Definition: quadrature_lib.cc:852
QSorted::QSorted
QSorted(const Quadrature< dim > &quad)
Definition: quadrature_lib.cc:769
QSimplex
Definition: quadrature_lib.h:611
QTrapez
Definition: quadrature_lib.h:126
QGaussLog::get_quadrature_weights
static std::vector< double > get_quadrature_weights(const unsigned int n)
QTrapez::QTrapez
QTrapez()
Definition: quadrature_lib.cc:832
QGaussRadauChebyshev::ep
const EndPoint ep
Definition: quadrature_lib.h:552
QMilne
Definition: quadrature_lib.h:141
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
QSplit
Definition: quadrature_lib.h:764
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
QGaussLogR
Definition: quadrature_lib.h:238
QGaussLog::QGaussLog
QGaussLog(const unsigned int n, const bool revert=false)
QGaussRadauChebyshev::right
@ right
Definition: quadrature_lib.h:539
QSimpson::QSimpson
QSimpson()
Definition: quadrature_lib.cc:839
QGaussLobatto
Definition: quadrature_lib.h:76
QMidpoint::QMidpoint
QMidpoint()
Definition: quadrature_lib.cc:825
QGaussLogR::QGaussLogR
QGaussLogR(const unsigned int n, const Point< dim > x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)
QGauss
Definition: quadrature_lib.h:40
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
QSorted::compare_weights
bool compare_weights(const unsigned int a, const unsigned int b) const
Definition: quadrature_lib.cc:801
QMidpoint
Definition: quadrature_lib.h:93
QGaussChebyshev::QGaussChebyshev
QGaussChebyshev(const unsigned int n)
Generate a formula with n quadrature points.
Definition: quadrature_lib.cc:1019
QGaussOneOverR::QGaussOneOverR
QGaussOneOverR(const unsigned int n, const Point< dim > singularity, const bool factor_out_singular_weight=false)
QGauss::QGauss
QGauss(const unsigned int n)
Definition: quadrature_lib.cc:811
QWeddle
Definition: quadrature_lib.h:155
Point< dim >
QTelles
Definition: quadrature_lib.h:467
quadrature.h
config.h
QGaussChebyshev
Definition: quadrature_lib.h:499
QGaussRadauChebyshev::EndPoint
EndPoint
Definition: quadrature_lib.h:530
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature
Definition: quadrature.h:85
QGaussRadauChebyshev
Definition: quadrature_lib.h:524
QTrianglePolar
Definition: quadrature_lib.h:670
QGaussLogR::fraction
const double fraction
Definition: quadrature_lib.h:265
internal::TensorProductManifoldImplementation::split_point
void split_point(const Point< dim1+dim2 > &source, Point< dim1 > &p1, Point< dim2 > &p2)
Definition: tensor_product_manifold.h:161