29#include <unordered_set>
34#include <Eigen/Eigenvalues>
57 "num_collect::integration::gauss_legendre_kronrod_integrator");
65template <
typename Signature>
75template <
typename Result, base::concepts::real_scalar Variable>
98 "Degree of Legendre function must be at least one.");
109 "Degree of Legendre function must be at least one.");
123 template <base::concepts::invocable_as<result_type(variable_type)> Function>
126 -> std::pair<result_type, result_type> {
145 sum_gauss *= half_width;
146 sum_kronrod *= half_width;
147 return {sum_gauss, sum_kronrod};
159 template <base::concepts::invocable_as<result_type(variable_type)> Function>
163 "The boundaries of the range to integrate on must satisfy left < "
173 constexpr auto stack_size =
static_cast<std::size_t
>(
174 std::numeric_limits<variable_type>::digits);
180 const auto [val_gauss, val_kronrod] =
186 const auto div_rate = (cur_right - cur_left) * inv_width;
191 if (remaining_right.
empty()) {
194 cur_left = cur_right;
195 cur_right = remaining_right.
top();
196 remaining_right.
pop();
198 remaining_right.
push(cur_right);
215 template <base::concepts::invocable_as<result_type(variable_type)> Function>
231 "Tolerance of absolute error must be positive value.");
246 "Tolerance of relative error must be positive value.");
259 val > std::numeric_limits<variable_type>::epsilon(), this->
logger(),
260 "Minimum rate of division of integration region must be larger "
261 "than the machine epsilon.");
268 using vector_type = Eigen::Matrix<variable_type, Eigen::Dynamic, 1>;
272 Eigen::Matrix<variable_type, Eigen::Dynamic, Eigen::Dynamic>;
284 for (
index_type i = 1; i <= (3 * n + 1) / 2; ++i) {
300 jacobi2gauss(a, b, extended_size, nodes_all, weights_all);
304 std::unordered_set<index_type> additional_nodes_index;
305 additional_nodes_index.reserve(
307 for (
index_type i = 0; i < extended_size; ++i) {
308 additional_nodes_index.insert(i);
313 std::numeric_limits<variable_type>::infinity();
318 if (dist < min_dist) {
324 additional_nodes_index.erase(min_ind);
328 const auto additional_nodes_vec = std::vector<index_type>(
329 additional_nodes_index.begin(), additional_nodes_index.end());
330 for (std::size_t i = 0; i < additional_nodes_vec.size(); ++i) {
332 nodes_all[additional_nodes_vec[i]];
334 weights_all[additional_nodes_vec[i]];
354 for (
index_type k = (m + 1) / 2; k > 0; --k) {
356 u += (a[k + n + 1] - a[l]) * t[k + 1] + b[k + n + 1] * s[k] -
363 u += (a[k + n + 1] - a[l]) * t[k + 1] + b[k + n + 1] * s[k] -
373 for (
index_type m = n - 1; m <= 2 * n - 3; ++m) {
376 for (
index_type k = m + 1 - n; k <= (m - 1) / 2; ++k) {
379 u += -(a[k + n + 1] - a[l]) * t[j + 1] -
380 b[k + n + 1] * s[j + 1] + b[l] * s[j + 2];
386 a[k] + (s[j + 1] - b[k + n + 1] * s[j + 2]) / t[j + 2];
389 b[k + n + 1] = s[j + 1] / s[j + 2];
393 a[2 * n] = a[n - 1] - b[2 * n] * s[1] * t[1];
407 matrix_type jacobi = matrix_type::Zero(size, size);
412 jacobi(i - 1, i) = temp;
413 jacobi(i, i - 1) = temp;
416 Eigen::SelfAdjointEigenSolver<matrix_type> eig(
417 jacobi, Eigen::ComputeEigenvectors);
420 nodes = eig.eigenvalues();
421 weights = b[0] * eig.eigenvectors().row(0).transpose().cwiseAbs2();
444 std::numeric_limits<variable_type>::epsilon() *
455 std::numeric_limits<variable_type>::epsilon() *
Definition of assertion macros.
#define NUM_COLLECT_ASSERT(CONDITION)
Macro to check whether a condition is satisfied.
vector_type weights_kronrod_
Additional weights for Gauss-Kronrod quadrature.
void prepare(index_type degree)
Prepare internal parameters.
std::decay_t< Variable > variable_type
Type of variables.
variable_type tol_abs_error_
Tolerance of absolute error.
auto integrate(const Function &function, variable_type left, variable_type right) const -> result_type
Integrate a function adaptively.
vector_type weights_gauss_
Weights for Gauss quadrature.
Eigen::Matrix< variable_type, Eigen::Dynamic, 1 > vector_type
Type of vectors.
static constexpr index_type default_degree
Default degree.
variable_type min_div_rate_
Minimum rate of division of integration region.
index_type degree_
Degree.
variable_type tol_rel_error_
Tolerance of relative error.
auto tol_rel_error(variable_type val) -> gauss_legendre_kronrod_integrator &
Set tolerance of relative error.
static constexpr variable_type default_tol_error
Default tolerance of error.
std::decay_t< Result > result_type
Type of results.
auto integrate_once(const Function &function, variable_type left, variable_type right) const -> std::pair< result_type, result_type >
Integrate a function and return two estimates.
void generate_jacobi_kronrod_matrix(vector_type &a, vector_type &b)
Generate Jacobi-Kronrod matrix as in a pseudo-code in the appendix of laurie1997.
auto tol_abs_error(variable_type val) -> gauss_legendre_kronrod_integrator &
Set tolerance of absolute error.
auto operator()(const Function &function, variable_type left, variable_type right) const -> result_type
Integrate a function.
gauss_legendre_kronrod_integrator(index_type degree=default_degree)
Constructor.
void compute_parameters()
Compute internal parameters.
static void jacobi2gauss(const vector_type &a, const vector_type &b, index_type size, vector_type &nodes, vector_type &weights)
Calculate Gauss quadrature from Jacobi matrix.
static constexpr variable_type default_min_div_rate
Default minimum rate of division of integration region.
vector_type nodes_kronrod_
Additional nodes for Gauss-Kronrod quadrature.
Eigen::Matrix< variable_type, Eigen::Dynamic, Eigen::Dynamic > matrix_type
Type of matrices.
auto min_div_rate(variable_type val) -> gauss_legendre_kronrod_integrator &
Set minimum rate of division of integration region.
vector_type nodes_gauss_
Nodes for Gauss quadrature.
vector_type weights_gauss_for_kronrod_
Weights of nodes in Gauss quadrature for Gauss-Kronrod quadrature.
Class to perform numerical adaptive integration with Gauss-Legendre-Kronrod formula in laurie1997.
Class of tags of logs without memory management.
Class to incorporate logging in algorithms.
logging_mixin(log_tag_view tag)
Constructor.
auto logger() const noexcept -> const num_collect::logging::logger &
Access to the logger.
Class of stacks using static arrays.
void pop() noexcept
Pop a value.
void push(const T &value)
Push a value.
auto top() const noexcept -> const T &
Get the top value.
auto empty() const noexcept -> bool
Check whether this stack is empty.
Definition of exceptions.
Definition of index_type type.
Definition of invocable_as concept.
Definition of log_tag_view class.
Definition of macros for logging.
Definition of logging_mixin class.
std::ptrdiff_t index_type
Type of indices in this library.
auto norm(const Matrix &matrix)
Calculate norm of a matrix.
constexpr T half
Value 0.5.
Namespace of numerical integration.
constexpr auto gauss_legendre_kronrod_integrator_tag
Tag of gauss_legendre_kronrod_integrator.
auto safe_cast(const From &value) -> To
Cast safely.
Definition of norm class.
Definition of NUM_COLLECT_PRECONDITION macro.
#define NUM_COLLECT_PRECONDITION(CONDITION,...)
Check whether a precondition is satisfied and throw an exception if not.
Definition of real_scalar concept.
Definition of safe_cast function.
Definition of static_stack class.