Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
quadrature_lib.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2022 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
37template <int dim>
38class QGauss : public Quadrature<dim>
39{
40public:
45 QGauss(const unsigned int n);
46};
47
48
71template <int dim>
72class QGaussLobatto : public Quadrature<dim>
73{
74public:
79 QGaussLobatto(const unsigned int n);
80};
81
82
83
93template <int dim>
94class QMidpoint : public Quadrature<dim>
95{
96public:
97 QMidpoint();
98};
99
100
105template <int dim>
106class QSimpson : public Quadrature<dim>
107{
108public:
109 QSimpson();
110};
111
112
113
122template <int dim>
123class QTrapezoid : public Quadrature<dim>
124{
125public:
126 QTrapezoid();
127};
128
129
141template <int dim>
143
144
145
152template <int dim>
153class QMilne : public Quadrature<dim>
154{
155public:
156 QMilne();
157};
158
159
166template <int dim>
167class QWeddle : public Quadrature<dim>
168{
169public:
170 QWeddle();
171};
172
173
174
188template <int dim>
189class QGaussLog : public Quadrature<dim>
190{
191public:
195 QGaussLog(const unsigned int n, const bool revert = false);
196
197private:
201 static std::vector<double>
202 get_quadrature_points(const unsigned int n);
203
207 static std::vector<double>
208 get_quadrature_weights(const unsigned int n);
209};
210
211
212
249template <int dim>
250class QGaussLogR : public Quadrature<dim>
251{
252public:
260 QGaussLogR(const unsigned int n,
261 const Point<dim> & x0 = Point<dim>(),
262 const double alpha = 1,
263 const bool factor_out_singular_weight = false);
264
270 QGaussLogR(QGaussLogR<dim> &&) noexcept = default;
271
272protected:
277 const double fraction;
278};
279
280
304template <int dim>
305class QGaussOneOverR : public Quadrature<dim>
306{
307public:
340 QGaussOneOverR(const unsigned int n,
341 const Point<dim> & singularity,
342 const bool factor_out_singular_weight = false);
377 QGaussOneOverR(const unsigned int n,
378 const unsigned int vertex_index,
379 const bool factor_out_singular_weight = false);
380
381private:
387 static unsigned int
388 quad_size(const Point<dim> &singularity, const unsigned int n);
389};
390
391
392
402template <int dim>
403class QSorted : public Quadrature<dim>
404{
405public:
410 QSorted(const Quadrature<dim> &quad);
411
412private:
418 bool
419 compare_weights(const unsigned int a, const unsigned int b) const;
420};
421
476template <int dim>
477class QTelles : public Quadrature<dim>
478{
479public:
486 QTelles(const Quadrature<1> &base_quad, const Point<dim> &singularity);
492 QTelles(const unsigned int n, const Point<dim> &singularity);
493};
494
506template <int dim>
507class QGaussChebyshev : public Quadrature<dim>
508{
509public:
511 QGaussChebyshev(const unsigned int n);
512};
513
514
529template <int dim>
531{
532public:
533 /* EndPoint is used to specify which of the two endpoints of the unit interval
534 * is used also as quadrature point
535 */
537 {
545 right
546 };
548 QGaussRadauChebyshev(const unsigned int n,
549 EndPoint ep = QGaussRadauChebyshev::left);
550
556
557private:
558 const EndPoint ep;
559};
560
574template <int dim>
576{
577public:
579 QGaussLobattoChebyshev(const unsigned int n);
580};
581
615template <int dim>
616class QSimplex : public Quadrature<dim>
617{
618public:
625 QSimplex(const Quadrature<dim> &quad);
626
651 compute_affine_transformation(
652 const std::array<Point<dim>, dim + 1> &vertices) const;
653};
654
673class QTrianglePolar : public QSimplex<2>
674{
675public:
683 QTrianglePolar(const Quadrature<1> &radial_quadrature,
684 const Quadrature<1> &angular_quadrature);
685
692 QTrianglePolar(const unsigned int n);
693};
694
727class QDuffy : public QSimplex<2>
728{
729public:
744 QDuffy(const Quadrature<1> &radial_quadrature,
745 const Quadrature<1> &angular_quadrature,
746 const double beta = 1.0);
747
755 QDuffy(const unsigned int n, const double beta);
756};
757
762template <int dim>
763class QSplit : public Quadrature<dim>
764{
765public:
800 QSplit(const QSimplex<dim> &base, const Point<dim> &split_point);
801};
802
827template <int dim>
828class QGaussSimplex : public QSimplex<dim>
829{
830public:
835 explicit QGaussSimplex(const unsigned int n_points_1D);
836};
837
868template <int dim>
870{
871public:
878 explicit QWitherdenVincentSimplex(const unsigned int n_points_1D,
879 const bool use_odd_order = true);
880};
881
889template <int dim>
890class QIteratedSimplex : public Quadrature<dim>
891{
892public:
893 QIteratedSimplex(const Quadrature<dim> &base_quadrature,
894 const unsigned int n_copies);
895};
896
900template <int dim>
901class QGaussWedge : public Quadrature<dim>
902{
903public:
909 explicit QGaussWedge(const unsigned int n_points_1D);
910};
911
915template <int dim>
916class QGaussPyramid : public Quadrature<dim>
917{
918public:
924 explicit QGaussPyramid(const unsigned int n_points_1D);
925};
926
929/* -------------- declaration of explicit specializations ------------- */
930
931#ifndef DOXYGEN
932template <>
933QGauss<1>::QGauss(const unsigned int n);
934template <>
935QGaussLobatto<1>::QGaussLobatto(const unsigned int n);
936
937template <>
938std::vector<double>
939QGaussLog<1>::get_quadrature_points(const unsigned int);
940template <>
941std::vector<double>
942QGaussLog<1>::get_quadrature_weights(const unsigned int);
943
944template <>
946template <>
948template <>
950template <>
952template <>
954template <>
955QGaussLog<1>::QGaussLog(const unsigned int n, const bool revert);
956template <>
957QGaussLogR<1>::QGaussLogR(const unsigned int n,
958 const Point<1> & x0,
959 const double alpha,
960 const bool flag);
961template <>
962QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
963 const unsigned int index,
964 const bool flag);
965template <>
966QTelles<1>::QTelles(const Quadrature<1> &base_quad,
967 const Point<1> & singularity);
968#endif // DOXYGEN
969
970
971
973#endif
Definition: point.h:111
QGaussLobatto(const unsigned int n)
const double fraction
QGaussLogR(const unsigned int n, const Point< dim > &x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false)
QGaussLogR(QGaussLogR< dim > &&) noexcept=default
QGaussLog(const unsigned int n, const bool revert=false)
static std::vector< double > get_quadrature_points(const unsigned int n)
static std::vector< double > get_quadrature_weights(const unsigned int n)
QGaussOneOverR(const unsigned int n, const unsigned int vertex_index, const bool factor_out_singular_weight=false)
QGaussOneOverR(const unsigned int n, const Point< dim > &singularity, const bool factor_out_singular_weight=false)
static unsigned int quad_size(const Point< dim > &singularity, const unsigned int n)
QGaussRadauChebyshev(QGaussRadauChebyshev< dim > &&) noexcept=default
QGauss(const unsigned int n)
QTelles(const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
Point< 3 > vertices[4]