basix More...
Functions | |
| xt::xtensor< double, 2 > | compute_jacobi_deriv (double a, std::size_t n, std::size_t nderiv, const xtl::span< const double > &x) | 
| std::vector< double > | compute_gauss_jacobi_points (double a, int m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | compute_gauss_jacobi_rule (double a, int m) | 
| Gauss-Jacobi quadrature rule (points and weights)  More... | |
| std::pair< xt::xarray< double >, std::vector< double > > | make_quadrature_line (int m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | make_quadrature_triangle_collapsed (std::size_t m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | make_quadrature_tetrahedron_collapsed (std::size_t m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | make_quadrature (const std::string &rule, cell::type celltype, int m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | make_gll_line (int m) | 
| std::pair< xt::xarray< double >, std::vector< double > > | compute_gll_rule (int m) | 
| GLL quadrature rule (points and weights)  | |
basix
Integration using Gauss-Jacobi quadrature on simplices. Other shapes can be obtained by using a product.
| std::vector< double > basix::quadrature::compute_gauss_jacobi_points | ( | double | a, | 
| int | m | ||
| ) | 
Finds the m roots of \(P_{m}^{a,0}\) on [-1,1] by Newton's method.
| [in] | a | weight in Jacobi (b=0) | 
| [in] | m | order | 
Computes the m roots of \(P_{m}^{a,0}\) on [-1,1] by Newton's method. The initial guesses are the Chebyshev points. Algorithm implemented from the pseudocode given by Karniadakis and Sherwin
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::compute_gauss_jacobi_rule | ( | double | a, | 
| int | m | ||
| ) | 
Gauss-Jacobi quadrature rule (points and weights)
| xt::xtensor< double, 2 > basix::quadrature::compute_jacobi_deriv | ( | double | a, | 
| std::size_t | n, | ||
| std::size_t | nderiv, | ||
| const xtl::span< const double > & | x | ||
| ) | 
Evaluate the nth Jacobi polynomial and derivatives with weight parameters (a, 0) at points x
| [in] | a | Jacobi weight a | 
| [in] | n | Order of polynomial | 
| [in] | nderiv | Number of derivatives (if zero, just compute polynomial itself) | 
| [in] | x | Points at which to evaluate | 
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::make_gll_line | ( | int | m | ) | 
Compute GLL line quadrature rule on [0, 1]
| m | order | 
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::make_quadrature | ( | const std::string & | rule, | 
| cell::type | celltype, | ||
| int | m | ||
| ) | 
Utility for quadrature rule on reference cell
| [in] | rule | Name of quadrature rule (or use "default") | 
| [in] | celltype | |
| [in] | m | Maximum degree of polynomial that this quadrature rule will integrate exactly | 
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::make_quadrature_line | ( | int | m | ) | 
Compute line quadrature rule on [0, 1]
| m | order | 
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::make_quadrature_tetrahedron_collapsed | ( | std::size_t | m | ) | 
Compute tetrahedron quadrature rule on [0, 1]x[0, 1]x[0, 1]
| [in] | m | order | 
| std::pair< xt::xarray< double >, std::vector< double > > basix::quadrature::make_quadrature_triangle_collapsed | ( | std::size_t | m | ) | 
Compute triangle quadrature rule on [0, 1]x[0, 1]
| [in] | m | order | 
 1.8.17