Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/basix/v0.9.0/cpp/quadrature_8h_source.html
Home Installation C++ docs Python docs
quadrature.h
1 // Copyright (c) 2020 Chris Richardson
2 // FEniCS Project
3 // SPDX-License-Identifier: MIT
4 
5 #pragma once
6 
7 #include "cell.h"
8 #include <xtl/xspan.hpp>
9 #include <utility>
10 #include <vector>
11 #include <xtensor/xtensor.hpp>
12 
14 
19 {
29 xt::xtensor<double, 2> compute_jacobi_deriv(double a, std::size_t n,
30  std::size_t nderiv,
31  const xtl::span<const double>& x);
32 
33 // Computes Gauss-Jacobi quadrature points
38 std::vector<double> compute_gauss_jacobi_points(double a, int m);
39 
41 std::pair<xt::xarray<double>, std::vector<double>>
42 compute_gauss_jacobi_rule(double a, int m);
43 
47 std::pair<xt::xarray<double>, std::vector<double>> make_quadrature_line(int m);
48 
52 std::pair<xt::xarray<double>, std::vector<double>>
54 
59 std::pair<xt::xarray<double>, std::vector<double>>
61 
69 std::pair<xt::xarray<double>, std::vector<double>>
70 make_quadrature(const std::string& rule, cell::type celltype, int m);
71 
75 std::pair<xt::xarray<double>, std::vector<double>> make_gll_line(int m);
76 
78 std::pair<xt::xarray<double>, std::vector<double>> compute_gll_rule(int m);
79 
80 } // namespace basix::quadrature
basix::quadrature::make_quadrature_triangle_collapsed
std::pair< xt::xarray< double >, std::vector< double > > make_quadrature_triangle_collapsed(std::size_t m)
Definition: quadrature.cpp:1655
basix::quadrature::make_quadrature
std::pair< xt::xarray< double >, std::vector< double > > make_quadrature(const std::string &rule, cell::type celltype, int m)
Definition: quadrature.cpp:1705
basix::cell::type
type
Cell type.
Definition: cell.h:16
basix::quadrature::compute_gll_rule
std::pair< xt::xarray< double >, std::vector< double > > compute_gll_rule(int m)
GLL quadrature rule (points and weights)
Definition: quadrature.cpp:1614
basix::quadrature::compute_jacobi_deriv
xt::xtensor< double, 2 > compute_jacobi_deriv(double a, std::size_t n, std::size_t nderiv, const xtl::span< const double > &x)
Definition: quadrature.cpp:1502
basix::quadrature
basix
Definition: quadrature.h:18
basix::quadrature::make_quadrature_tetrahedron_collapsed
std::pair< xt::xarray< double >, std::vector< double > > make_quadrature_tetrahedron_collapsed(std::size_t m)
Definition: quadrature.cpp:1677
basix::quadrature::compute_gauss_jacobi_rule
std::pair< xt::xarray< double >, std::vector< double > > compute_gauss_jacobi_rule(double a, int m)
Gauss-Jacobi quadrature rule (points and weights)
Definition: quadrature.cpp:1591
basix::quadrature::compute_gauss_jacobi_points
std::vector< double > compute_gauss_jacobi_points(double a, int m)
Definition: quadrature.cpp:1552
basix::quadrature::make_gll_line
std::pair< xt::xarray< double >, std::vector< double > > make_gll_line(int m)
Definition: quadrature.cpp:1646
basix::quadrature::make_quadrature_line
std::pair< xt::xarray< double >, std::vector< double > > make_quadrature_line(int m)
Definition: quadrature.cpp:1637