QUGaR 0.0.4
Loading...
Searching...
No Matches
quadrature.hpp
Go to the documentation of this file.
1// --------------------------------------------------------------------------
2//
3// Copyright (C) 2025-present by Pablo Antolin
4//
5// This file is part of the QUGaR library.
6//
7// SPDX-License-Identifier: MIT
8//
9// --------------------------------------------------------------------------
10
11#ifndef QUGAR_LIBRARY_QUADRATURE_HPP
12
13#define QUGAR_LIBRARY_QUADRATURE_HPP
14
21
22
23#include <qugar/bbox.hpp>
24#include <qugar/point.hpp>
25#include <qugar/types.hpp>
26
27#include <array>
28#include <memory>
29#include <type_traits>
30
31
32namespace qugar {
33
34
43class Gauss
44{
45public:
47 static constexpr int n_max = 100;
48
55 static real get_abscissa(const int n_points, const int pt_id);
56
63 static real get_weight(const int n_points, const int pt_id);
64
65
66private:
68 static constexpr int n_entries = n_max * (n_max + 1);
69
81 static const std::array<real, n_entries> &get_quad_data();
82};
83
93{
94public:
96 static constexpr int n_max = 100;
97
104 static real get_abscissa(const int n_points, const int pt_id);
105
112 static real get_weight(const int n_points, const int pt_id);
113
114
115private:
117 static constexpr int n_entries = n_max * (n_max + 1);
118
130 static const std::array<real, n_entries> &get_quad_data();
131};
132
134template<int dim> class Quadrature
135{
136public:
138 using Pt = Point<dim>;
139
141 using Domain = std::conditional_t<dim == 1, Interval, BoundBox<dim>>;
142
143
146 Quadrature() = delete;
147
152 Quadrature(const std::vector<Pt> &points, const std::vector<real> &weights);
153
157 [[nodiscard]] static std::shared_ptr<Quadrature> create_Gauss_01(int n_points);
158
162 [[nodiscard]] static std::shared_ptr<Quadrature> create_Gauss_01(const std::array<int, dim> n_points);
163
168 [[nodiscard]] static std::shared_ptr<Quadrature> create_Gauss(const std::array<int, dim> n_points,
169 const Domain &domain);
170
174 [[nodiscard]] static std::shared_ptr<Quadrature> create_tanh_sinh_01(int n_points);
175
179 [[nodiscard]] static std::shared_ptr<Quadrature> create_tanh_sinh_01(const std::array<int, dim> n_points);
180
185 [[nodiscard]] static std::shared_ptr<Quadrature> create_tanh_sinh(const std::array<int, dim> n_points,
186 const Domain &domain);
187
192 void scale_weights(const real ratio);
193
196 [[nodiscard]] const std::vector<Pt> &points() const;
197
200 [[nodiscard]] const std::vector<real> &weights() const;
201
204 [[nodiscard]] std::vector<Pt> &points();
205
208 [[nodiscard]] std::vector<real> &weights();
209
212 [[nodiscard]] std::size_t get_num_points() const;
213
214private:
216 std::vector<Pt> points_;
218 std::vector<real> weights_;
219
226 template<typename QuadType>
227 [[nodiscard]] static std::shared_ptr<Quadrature<dim>> create_tp_quadrature(const std::array<int, dim> n_points,
228 const Domain &domain);
229};
230
231
232}// namespace qugar
233
234#endif// QUGAR_LIBRARY_QUADRATURE_HPP
Definition of Cartesian bounding box class.
Helper class for computing the abscissae and weights of the Gauss-Legendre quadrature rule referred t...
Definition quadrature.hpp:44
static real get_abscissa(const int n_points, const int pt_id)
Gets a Gauss-Legendre abscissa in the [0,1] domain.
static constexpr int n_entries
Number of entries in data vectors.
Definition quadrature.hpp:68
static real get_weight(const int n_points, const int pt_id)
Gets a Gauss-Legendre weight in the [0,1] domain.
static const std::array< real, n_entries > & get_quad_data()
Gets a reference to a static array containing all the information of the quadrature.
static constexpr int n_max
Maximum precomputed degree of Gauss-Legendre points.
Definition quadrature.hpp:47
Class for storing dim-dimensional quadratures (non-tensor product).
Definition quadrature.hpp:135
static std::shared_ptr< Quadrature > create_Gauss(const std::array< int, dim > n_points, const Domain &domain)
Creates a new quadrature class instance with a Gauss-Legendre quadrature in the domain box.
std::vector< real > & weights()
Returns the vector of weights.
static std::shared_ptr< Quadrature > create_Gauss_01(int n_points)
Creates a new quadrature class instance with a Gauss-Legendre quadrature in [0, 1].
const std::vector< real > & weights() const
Returns the vector of weights.
static std::shared_ptr< Quadrature > create_Gauss_01(const std::array< int, dim > n_points)
Creates a new quadrature class instance with a Gauss-Legendre quadrature in [0, 1].
static std::shared_ptr< Quadrature > create_tanh_sinh_01(int n_points)
Creates a new quadrature class instance with a tanh-sinh quadrature in [0, 1].
static std::shared_ptr< Quadrature > create_tanh_sinh_01(const std::array< int, dim > n_points)
Creates a new quadrature class instance with a tanh-sinh quadrature in [0, 1].
std::vector< Pt > points_
Vector of points.
Definition quadrature.hpp:216
void scale_weights(const real ratio)
Scales the weights multiplying them by the given ratio.
const std::vector< Pt > & points() const
Returns the vector of points.
std::vector< real > weights_
Vector of weights.
Definition quadrature.hpp:218
static std::shared_ptr< Quadrature > create_tanh_sinh(const std::array< int, dim > n_points, const Domain &domain)
Creates a new quadrature class instance with a tanh-sinh quadrature in the domain box.
std::size_t get_num_points() const
Returns the number of points (and weights).
Point< dim > Pt
Point type.
Definition quadrature.hpp:138
Quadrature(const std::vector< Pt > &points, const std::vector< real > &weights)
Constructor.
std::conditional_t< dim==1, Interval, BoundBox< dim > > Domain
Domain type.
Definition quadrature.hpp:141
static std::shared_ptr< Quadrature< dim > > create_tp_quadrature(const std::array< int, dim > n_points, const Domain &domain)
Creates a tensor-product quadrature.
std::vector< Pt > & points()
Returns the vector of points.
Quadrature()=delete
Default constructor.
Helper class for computing the abscissae and weights of the tanh-sinh quadrature rule referred to the...
Definition quadrature.hpp:93
static real get_weight(const int n_points, const int pt_id)
Gets a tanh-sinh weight in the [0,1] domain.
static constexpr int n_max
Maximum precomputed degree of tanh-sinh points.
Definition quadrature.hpp:96
static real get_abscissa(const int n_points, const int pt_id)
Gets a tanh-sinh abscissa in the [0,1] domain.
static constexpr int n_entries
Number of entries in data vectors.
Definition quadrature.hpp:117
static const std::array< real, n_entries > & get_quad_data()
Gets a reference to a static array containing all the information of the quadrature.
QUGaR's main namespace.
Definition affine_transf.hpp:28
double real
Definition types.hpp:18
Vector< T, dim > Point
Class representing a dim-dimensional Point.
Definition point.hpp:34
Definition and implementation of Point class.
Declaration of types.