Reference documentation for deal.II version 9.3.3
\(\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 - 2021 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
88template <int dim>
89class QMidpoint : public Quadrature<dim>
90{
91public:
92 QMidpoint();
93};
94
95
100template <int dim>
101class QSimpson : public Quadrature<dim>
102{
103public:
104 QSimpson();
105};
106
107
108
117template <int dim>
118class QTrapezoid : public Quadrature<dim>
119{
120public:
121 QTrapezoid();
122};
123
124
136template <int dim>
138
139
140
147template <int dim>
148class QMilne : public Quadrature<dim>
149{
150public:
151 QMilne();
152};
153
154
161template <int dim>
162class QWeddle : public Quadrature<dim>
163{
164public:
165 QWeddle();
166};
167
168
169
183template <int dim>
184class QGaussLog : public Quadrature<dim>
185{
186public:
190 QGaussLog(const unsigned int n, const bool revert = false);
191
192private:
196 static std::vector<double>
197 get_quadrature_points(const unsigned int n);
198
202 static std::vector<double>
203 get_quadrature_weights(const unsigned int n);
204};
205
206
207
244template <int dim>
245class QGaussLogR : public Quadrature<dim>
246{
247public:
255 QGaussLogR(const unsigned int n,
256 const Point<dim> & x0 = Point<dim>(),
257 const double alpha = 1,
258 const bool factor_out_singular_weight = false);
259
265 QGaussLogR(QGaussLogR<dim> &&) noexcept = default;
266
267protected:
272 const double fraction;
273};
274
275
299template <int dim>
300class QGaussOneOverR : public Quadrature<dim>
301{
302public:
335 QGaussOneOverR(const unsigned int n,
336 const Point<dim> & singularity,
337 const bool factor_out_singular_weight = false);
372 QGaussOneOverR(const unsigned int n,
373 const unsigned int vertex_index,
374 const bool factor_out_singular_weight = false);
375
376private:
382 static unsigned int
383 quad_size(const Point<dim> &singularity, const unsigned int n);
384};
385
386
387
397template <int dim>
398class QSorted : public Quadrature<dim>
399{
400public:
405 QSorted(const Quadrature<dim> &quad);
406
407private:
413 bool
414 compare_weights(const unsigned int a, const unsigned int b) const;
415};
416
471template <int dim>
472class QTelles : public Quadrature<dim>
473{
474public:
481 QTelles(const Quadrature<1> &base_quad, const Point<dim> &singularity);
487 QTelles(const unsigned int n, const Point<dim> &singularity);
488};
489
501template <int dim>
502class QGaussChebyshev : public Quadrature<dim>
503{
504public:
506 QGaussChebyshev(const unsigned int n);
507};
508
509
524template <int dim>
526{
527public:
528 /* EndPoint is used to specify which of the two endpoints of the unit interval
529 * is used also as quadrature point
530 */
532 {
540 right
541 };
543 QGaussRadauChebyshev(const unsigned int n,
544 EndPoint ep = QGaussRadauChebyshev::left);
545
551
552private:
553 const EndPoint ep;
554};
555
569template <int dim>
571{
572public:
574 QGaussLobattoChebyshev(const unsigned int n);
575};
576
607template <int dim>
608class QSimplex : public Quadrature<dim>
609{
610public:
617 QSimplex(const Quadrature<dim> &quad);
618
643 compute_affine_transformation(
644 const std::array<Point<dim>, dim + 1> &vertices) const;
645};
646
665class QTrianglePolar : public QSimplex<2>
666{
667public:
675 QTrianglePolar(const Quadrature<1> &radial_quadrature,
676 const Quadrature<1> &angular_quadrature);
677
684 QTrianglePolar(const unsigned int n);
685};
686
719class QDuffy : public QSimplex<2>
720{
721public:
736 QDuffy(const Quadrature<1> &radial_quadrature,
737 const Quadrature<1> &angular_quadrature,
738 const double beta = 1.0);
739
747 QDuffy(const unsigned int n, const double beta);
748};
749
754template <int dim>
755class QSplit : public Quadrature<dim>
756{
757public:
792 QSplit(const QSimplex<dim> &base, const Point<dim> &split_point);
793};
794
811template <int dim>
812class QGaussSimplex : public QSimplex<dim>
813{
814public:
819 explicit QGaussSimplex(const unsigned int n_points_1D);
820};
821
842template <int dim>
844{
845public:
850 explicit QWitherdenVincentSimplex(const unsigned int n_points_1D);
851};
852
856template <int dim>
857class QGaussWedge : public Quadrature<dim>
858{
859public:
865 explicit QGaussWedge(const unsigned int n_points_1D);
866};
867
871template <int dim>
872class QGaussPyramid : public Quadrature<dim>
873{
874public:
880 explicit QGaussPyramid(const unsigned int n_points_1D);
881};
882
885/* -------------- declaration of explicit specializations ------------- */
886
887#ifndef DOXYGEN
888template <>
889QGauss<1>::QGauss(const unsigned int n);
890template <>
891QGaussLobatto<1>::QGaussLobatto(const unsigned int n);
892
893template <>
894std::vector<double>
895QGaussLog<1>::get_quadrature_points(const unsigned int);
896template <>
897std::vector<double>
898QGaussLog<1>::get_quadrature_weights(const unsigned int);
899
900template <>
902template <>
904template <>
906template <>
908template <>
910template <>
911QGaussLog<1>::QGaussLog(const unsigned int n, const bool revert);
912template <>
913QGaussLogR<1>::QGaussLogR(const unsigned int n,
914 const Point<1> & x0,
915 const double alpha,
916 const bool flag);
917template <>
918QGaussOneOverR<2>::QGaussOneOverR(const unsigned int n,
919 const unsigned int index,
920 const bool flag);
921template <>
922QTelles<1>::QTelles(const Quadrature<1> &base_quad,
923 const Point<1> & singularity);
924#endif // DOXYGEN
925
926
927
929#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:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
Point< 3 > vertices[4]
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void split_point(const Point< dim1+dim2 > &source, Point< dim1 > &p1, Point< dim2 > &p2)