<p align="center">
<a href="https://github.com/nschloe/orthopy"><img alt="orthopy" src="https://nschloe.github.io/orthopy/orthopy-logo-with-text.png" width="30%"></a>
<p align="center">All about orthogonal polynomials.</p>
</p>
[![PyPi Version](https://img.shields.io/pypi/v/orthopy.svg?style=flat-square)](https://pypi.org/project/orthopy)
[![PyPI pyversions](https://img.shields.io/pypi/pyversions/orthopy.svg?style=flat-square)](https://pypi.org/pypi/orthopy/)
[![DOI](https://zenodo.org/badge/DOI/10.5281/zenodo.1173151.svg?style=flat-square)](https://doi.org/10.5281/zenodo.1173151)
[![GitHub stars](https://img.shields.io/github/stars/nschloe/orthopy.svg?style=flat-square&logo=github&label=Stars&logoColor=white)](https://github.com/nschloe/orthopy)
[![Downloads](https://pepy.tech/badge/orthopy/month?style=flat-square)](https://pepy.tech/project/orthopy)
<!--[![PyPi downloads](https://img.shields.io/pypi/dm/orthopy.svg?style=flat-square)](https://pypistats.org/packages/orthopy)-->
[![Discord](https://img.shields.io/static/v1?logo=discord&label=chat&message=on%20discord&color=7289da&style=flat-square)](https://discord.gg/hnTJ5MRX2Y)
[![orthogonal](https://img.shields.io/badge/orthogonal-yes-ff69b4.svg?style=flat-square)](https://github.com/nschloe/orthopy)
[![gh-actions](https://img.shields.io/github/workflow/status/nschloe/orthopy/ci?style=flat-square)](https://github.com/nschloe/orthopy/actions?query=workflow%3Aci)
[![codecov](https://img.shields.io/codecov/c/github/nschloe/orthopy.svg?style=flat-square)](https://codecov.io/gh/nschloe/orthopy)
[![LGTM](https://img.shields.io/lgtm/grade/python/github/nschloe/orthopy.svg?style=flat-square)](https://lgtm.com/projects/g/nschloe/orthopy)
[![Code style: black](https://img.shields.io/badge/code%20style-black-000000.svg?style=flat-square)](https://github.com/psf/black)
orthopy provides various orthogonal polynomial classes for
[lines](#line-segment--1-1-with-weight-function-1-x%CE%B1-1-x%CE%B2),
[triangles](#triangle-42),
[disks](#disk-s2),
[spheres](#sphere-u2),
[n-cubes](#n-cube-cn),
[the nD space with weight function exp(-r<sup>2</sup>)](#nd-space-with-weight-function-exp-r2-enr2)
and more.
All computations are done using numerically stable recurrence schemes. Furthermore, all
functions are fully vectorized and can return results in _exact arithmetic_.
### Basic usage
Install orthopy from [PyPi](https://pypi.org/project/orthopy) via
```
pip install orthopy
```
The main function of all submodules is the iterator `Eval` which evaluates the series of
orthogonal polynomials with increasing degree at given points using a recurrence
relation, e.g.,
```python
import orthopy
x = 0.5
evaluator = orthopy.c1.legendre.Eval(x, "classical")
for _ in range(5):
print(next(evaluator))
```
```python
1.0 # P_0(0.5)
0.5 # P_1(0.5)
-0.125 # P_2(0.5)
-0.4375 # P_3(0.5)
-0.2890625 # P_4(0.5)
```
Other ways of getting the first `n` items are
<!--pytest-codeblocks:skip-->
```python
evaluator = Eval(x, "normal")
vals = [next(evaluator) for _ in range(n)]
import itertools
vals = list(itertools.islice(Eval(x, "normal"), n))
```
Instead of evaluating at only one point, you can provide any array for `x`; the
polynomials will then be evaluated for all points at once. You can also use sympy for
symbolic computation:
```python
import itertools
import orthopy
import sympy
x = sympy.Symbol("x")
evaluator = orthopy.c1.legendre.Eval(x, "classical")
for val in itertools.islice(evaluator, 5):
print(sympy.expand(val))
```
```
1
x
3*x**2/2 - 1/2
5*x**3/2 - 3*x/2
35*x**4/8 - 15*x**2/4 + 3/8
```
All `Eval` methods have a `scaling` argument which can have three values:
* `"monic"`: The leading coefficient is 1.
* `"classical"`: The maximum value is 1 (or (n+alpha over n)).
* `"normal"`: The integral of the squared function over the domain is 1.
For univariate ("one-dimensional") integrals, every new iteration contains one function.
For bivariate ("two-dimensional") domains, every level will contain one function more
than the previous, and similarly for multivariate families. See the tree plots below.
### Line segment (-1, +1) with weight function (1-x)<sup>α</sup> (1+x)<sup>β</sup>
<img src="https://nschloe.github.io/orthopy/legendre.svg" width="100%"> | <img src="https://nschloe.github.io/orthopy/chebyshev1.svg" width="100%"> | <img src="https://nschloe.github.io/orthopy/chebyshev2.svg" width="100%">
:-------------------:|:------------------:|:-------------:|
Legendre | Chebyshev 1 | Chebyshev 2 |
Jacobi, Gegenbauer (α=β), Chebyshev 1 (α=β=-1/2), Chebyshev 2 (α=β=1/2), Legendre
(α=β=0) polynomials.
<!--pytest-codeblocks:skip-->
```python
import orthopy
orthopy.c1.legendre.Eval(x, "normal")
orthopy.c1.chebyshev1.Eval(x, "normal")
orthopy.c1.chebyshev2.Eval(x, "normal")
orthopy.c1.gegenbauer.Eval(x, "normal", lmbda)
orthopy.c1.jacobi.Eval(x, "normal", alpha, beta)
```
The plots above are generated with
```python
import orthopy
orthopy.c1.jacobi.show(5, "normal", 0.0, 0.0)
# plot, savefig also exist
```
Recurrence coefficients can be explicitly retrieved by
```python
import orthopy
rc = orthopy.c1.jacobi.RecurrenceCoefficients(
"monic", # or "classical", "normal"
alpha=0, beta=0, symbolic=True
)
print(rc.p0)
for k in range(5):
print(rc[k])
```
```
1
(1, 0, None)
(1, 0, 1/3)
(1, 0, 4/15)
(1, 0, 9/35)
(1, 0, 16/63)
```
### 1D half-space with weight function x<sup>α</sup> exp(-r)
<img src="https://nschloe.github.io/orthopy/e1r.svg" width="45%">
(Generalized) Laguerre polynomials.
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.e1r.Eval(x, alpha=0, scaling="normal")
```
### 1D space with weight function exp(-r<sup>2</sup>)
<img src="https://nschloe.github.io/orthopy/e1r2.svg" width="45%">
Hermite polynomials come in two standardizations:
* `"physicists"` (against the weight function `exp(-x ** 2)`
* `"probabilists"` (against the weight function `1 / sqrt(2 * pi) * exp(-x ** 2 / 2)`
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.e1r2.Eval(
x,
"probabilists", # or "physicists"
"normal"
)
```
#### Associated Legendre "polynomials"
<img src="https://nschloe.github.io/orthopy/associated-legendre.svg" width="45%">
Not all of those are polynomials, so they should really be called associated Legendre
_functions_. The <i>k</i>th iteration contains _2k+1_ functions, indexed from
_-k_ to _k_. (See the color grouping in the above plot.)
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.c1.associated_legendre.Eval(
x, phi=None, standardization="natural", with_condon_shortley_phase=True
)
```
### Triangle (_T<sub>2</sub>_)
<img src="https://nschloe.github.io/orthopy/triangle-tree.png" width="40%">
orthopy's triangle orthogonal polynomials are evaluated in terms of [barycentric
coordinates](https://en.wikipedia.org/wiki/Barycentric_coordinate_system), so the
`X.shape[0]` has to be 3.
```python
import orthopy
bary = [0.1, 0.7, 0.2]
evaluator = orthopy.t2.Eval(bary, "normal")
```
### Disk (_S<sub>2</sub>_)
<img src="https://nschloe.github.io/orthopy/disk-xu-tree.png" width="70%"> | <img src="https://nschloe.github.io/orthopy/disk-zernike-tree.png" width="70%"> | <img src="https://nschloe.github.io/orthopy/disk-zernike2-tree.png" width="70%">
:------------:|:-----------------:|:-----------:|
Xu | [Zernike](https://en.wikipedia.org/wiki/Zernike_polynomials) | Zernike 2 |
orthopy contains several families of orthogonal polynomials on the unit disk: After
[Xu](https://arxiv.org/abs/1701.02709),
[Zernike](https://en.wikipedia.org/wiki/Zernike_polynomials), and a simplified version
of Zernike polynomials.
```python
import orthopy
x = [0.1, -0.3]
evaluator = orthopy.s2.xu.Eval(x, "normal")
# evaluator = orthopy.s2.zernike.Eval(x, "normal")
# evaluator = orthopy.s2.zernike2.Eval(x, "normal")
```
### Sphere (_U<sub>3</sub>_)
<img src="https://nschloe.github.io/orthopy/sph-tree.png" width="50%">
Complex-valued _spherical harmonics,_ plotted with
[cplot](https://github.com/nschloe/cplot/) coloring (black=zero, green=real positive,
pink=real negative, blue=imaginary positive, yellow=imaginary negative). The functions
in the middle are real-valued. The complex angle takes _n_ turns on the <i>n</i>th
level.
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.u3.EvalCartesian(
x,
scaling="quantum mechanic" # or "acoustic", "geodetic", "schmidt"
)
evaluator = orthopy.u3.EvalSpherical(
theta_phi, # polar, azimuthal angles
scaling="quantum mechanic" # or "acoustic", "geodetic", "schmidt"
)
```
To generate the above plot, write the tree mesh to a file
```python
import orthopy
orthopy.u3.write_tree("u3.vtk", 5, "quantum mechanic")
```
and open it with [ParaView](https://www.paraview.org/). Select the _srgb1_ data set and
turn off _Map Scalars_.
### _n_-Cube (_C<sub>n</sub>_)
<img src="https://nschloe.github.io/orthopy/c1.svg" width="100%"> | <img src="https://nschloe.github.io/orthopy/c2.png" width="100%"> | <img src="https://nschloe.github.io/orthopy/c3.png" width="100%">
:-------------------------:|:------------------:|:---------------:|
C<sub>1</sub> (Legendre) | C<sub>2</sub> | C<sub>3</sub> |
Jacobi product polynomials.
All polynomials are normalized on the n-dimensional cube. The dimensionality is
determined by `X.shape[0]`.
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.cn.Eval(X, alpha=0, beta=0)
values, degrees = next(evaluator)
```
### <i>n</i>D space with weight function exp(-r<sup>2</sup>) (_E<sub>n</sub><sup>r<sup>2</sup></sup>_)
<img src="https://nschloe.github.io/orthopy/e1r2.svg" width="100%"> | <img src="https://nschloe.github.io/orthopy/e2r2.png" width="100%"> | <img src="https://nschloe.github.io/orthopy/e3r2.png" width="100%">
:-------------------------:|:------------------:|:---------------:|
_E<sub>1</sub><sup>r<sup>2</sup></sup>_ | _E<sub>2</sub><sup>r<sup>2</sup></sup>_ | _E<sub>3</sub><sup>r<sup>2</sup></sup>_ |
Hermite product polynomials.
All polynomials are normalized over the measure. The dimensionality is determined by
`X.shape[0]`.
<!--pytest-codeblocks:skip-->
```python
evaluator = orthopy.enr2.Eval(
x,
standardization="probabilists" # or "physicists"
)
values, degrees = next(evaluator)
```
### Other tools
* Generating recurrence coefficients for 1D domains with
[Stieltjes](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#stieltjes),
[Golub-Welsch](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#golub-welsch),
[Chebyshev](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#chebyshev), and
[modified
Chebyshev](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#modified-chebyshev).
* The the sanity of recurrence coefficients with test 3 from [Gautschi's article](https://doi.org/10.1007/BF02218441):
computing the weighted sum of orthogonal polynomials:
<!--pytest-codeblocks:skip-->
```python
orthopy.tools.gautschi_test_3(moments, alpha, beta)
```
* [Clenshaw algorithm](https://en.wikipedia.org/wiki/Clenshaw_algorithm) for
computing the weighted sum of orthogonal polynomials:
<!--pytest-codeblocks:skip-->
```python
vals = orthopy.c1.clenshaw(a, alpha, beta, t)
```
### Installation
orthopy is [available from the Python Package
Index](https://pypi.python.org/pypi/orthopy/), so use
```
pip install orthopy
```
to install.
### Testing
To run the tests, simply check out this repository and run
```
pytest
```
### Relevant publications
* [Robert C. Kirby, Singularity-free evaluation of collapsed-coordinate orthogonal polynomials, ACM Transactions on Mathematical Software (TOMS), Volume 37, Issue 1, January 2010](https://doi.org/10.1145/1644001.1644006)
* [Abedallah Rababah, Recurrence Relations for Orthogonal Polynomials on Triangular Domains, MDPI Mathematics 2016, 4(2)](https://doi.org/10.3390/math4020025)
* [Yuan Xu, Orthogonal polynomials of several variables, arxiv.org, January 2017](https://arxiv.org/abs/1701.02709)
### License
This software is published under the [GPLv3 license](https://www.gnu.org/licenses/gpl-3.0.en.html).
Raw data
{
"_id": null,
"home_page": "https://zenodo.org/record/5751934",
"name": "orthopy-gpl",
"maintainer": "Aleksei Gerasimov",
"docs_url": null,
"requires_python": ">=3.7",
"maintainer_email": "ger-alex@seznam.cz",
"keywords": "",
"author": "Nico Schl\u00f6mer",
"author_email": "nico.schloemer@gmail.com",
"download_url": "https://files.pythonhosted.org/packages/57/8b/43ba088d37a2d4bc212cb02bcfbba7e72a9182a7459730c5cfc304a99fb5/orthopy-gpl-0.9.5.tar.gz",
"platform": null,
"description": "<p align=\"center\">\n <a href=\"https://github.com/nschloe/orthopy\"><img alt=\"orthopy\" src=\"https://nschloe.github.io/orthopy/orthopy-logo-with-text.png\" width=\"30%\"></a>\n <p align=\"center\">All about orthogonal polynomials.</p>\n</p>\n\n[![PyPi Version](https://img.shields.io/pypi/v/orthopy.svg?style=flat-square)](https://pypi.org/project/orthopy)\n[![PyPI pyversions](https://img.shields.io/pypi/pyversions/orthopy.svg?style=flat-square)](https://pypi.org/pypi/orthopy/)\n[![DOI](https://zenodo.org/badge/DOI/10.5281/zenodo.1173151.svg?style=flat-square)](https://doi.org/10.5281/zenodo.1173151)\n[![GitHub stars](https://img.shields.io/github/stars/nschloe/orthopy.svg?style=flat-square&logo=github&label=Stars&logoColor=white)](https://github.com/nschloe/orthopy)\n[![Downloads](https://pepy.tech/badge/orthopy/month?style=flat-square)](https://pepy.tech/project/orthopy)\n<!--[![PyPi downloads](https://img.shields.io/pypi/dm/orthopy.svg?style=flat-square)](https://pypistats.org/packages/orthopy)-->\n\n[![Discord](https://img.shields.io/static/v1?logo=discord&label=chat&message=on%20discord&color=7289da&style=flat-square)](https://discord.gg/hnTJ5MRX2Y)\n[![orthogonal](https://img.shields.io/badge/orthogonal-yes-ff69b4.svg?style=flat-square)](https://github.com/nschloe/orthopy)\n\n[![gh-actions](https://img.shields.io/github/workflow/status/nschloe/orthopy/ci?style=flat-square)](https://github.com/nschloe/orthopy/actions?query=workflow%3Aci)\n[![codecov](https://img.shields.io/codecov/c/github/nschloe/orthopy.svg?style=flat-square)](https://codecov.io/gh/nschloe/orthopy)\n[![LGTM](https://img.shields.io/lgtm/grade/python/github/nschloe/orthopy.svg?style=flat-square)](https://lgtm.com/projects/g/nschloe/orthopy)\n[![Code style: black](https://img.shields.io/badge/code%20style-black-000000.svg?style=flat-square)](https://github.com/psf/black)\n\northopy provides various orthogonal polynomial classes for\n[lines](#line-segment--1-1-with-weight-function-1-x%CE%B1-1-x%CE%B2),\n[triangles](#triangle-42),\n[disks](#disk-s2),\n[spheres](#sphere-u2),\n[n-cubes](#n-cube-cn),\n[the nD space with weight function exp(-r<sup>2</sup>)](#nd-space-with-weight-function-exp-r2-enr2)\nand more.\nAll computations are done using numerically stable recurrence schemes. Furthermore, all\nfunctions are fully vectorized and can return results in _exact arithmetic_.\n\n### Basic usage\n\nInstall orthopy from [PyPi](https://pypi.org/project/orthopy) via\n```\npip install orthopy\n```\nThe main function of all submodules is the iterator `Eval` which evaluates the series of\northogonal polynomials with increasing degree at given points using a recurrence\nrelation, e.g.,\n```python\nimport orthopy\n\nx = 0.5\n\nevaluator = orthopy.c1.legendre.Eval(x, \"classical\")\nfor _ in range(5):\n print(next(evaluator))\n```\n```python\n1.0 # P_0(0.5)\n0.5 # P_1(0.5)\n-0.125 # P_2(0.5)\n-0.4375 # P_3(0.5)\n-0.2890625 # P_4(0.5)\n```\nOther ways of getting the first `n` items are\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = Eval(x, \"normal\")\nvals = [next(evaluator) for _ in range(n)]\n\nimport itertools\nvals = list(itertools.islice(Eval(x, \"normal\"), n))\n```\nInstead of evaluating at only one point, you can provide any array for `x`; the\npolynomials will then be evaluated for all points at once. You can also use sympy for\nsymbolic computation:\n```python\nimport itertools\nimport orthopy\nimport sympy\n\nx = sympy.Symbol(\"x\")\n\nevaluator = orthopy.c1.legendre.Eval(x, \"classical\")\nfor val in itertools.islice(evaluator, 5):\n print(sympy.expand(val))\n```\n```\n1\nx\n3*x**2/2 - 1/2\n5*x**3/2 - 3*x/2\n35*x**4/8 - 15*x**2/4 + 3/8\n```\n\nAll `Eval` methods have a `scaling` argument which can have three values:\n\n * `\"monic\"`: The leading coefficient is 1.\n * `\"classical\"`: The maximum value is 1 (or (n+alpha over n)).\n * `\"normal\"`: The integral of the squared function over the domain is 1.\n\nFor univariate (\"one-dimensional\") integrals, every new iteration contains one function.\nFor bivariate (\"two-dimensional\") domains, every level will contain one function more\nthan the previous, and similarly for multivariate families. See the tree plots below.\n\n\n### Line segment (-1, +1) with weight function (1-x)<sup>\u03b1</sup> (1+x)<sup>\u03b2</sup>\n\n<img src=\"https://nschloe.github.io/orthopy/legendre.svg\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/chebyshev1.svg\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/chebyshev2.svg\" width=\"100%\">\n:-------------------:|:------------------:|:-------------:|\nLegendre | Chebyshev 1 | Chebyshev 2 |\n\nJacobi, Gegenbauer (\u03b1=\u03b2), Chebyshev 1 (\u03b1=\u03b2=-1/2), Chebyshev 2 (\u03b1=\u03b2=1/2), Legendre\n(\u03b1=\u03b2=0) polynomials.\n<!--pytest-codeblocks:skip-->\n```python\nimport orthopy\n\northopy.c1.legendre.Eval(x, \"normal\")\northopy.c1.chebyshev1.Eval(x, \"normal\")\northopy.c1.chebyshev2.Eval(x, \"normal\")\northopy.c1.gegenbauer.Eval(x, \"normal\", lmbda)\northopy.c1.jacobi.Eval(x, \"normal\", alpha, beta)\n```\n\nThe plots above are generated with\n```python\nimport orthopy\n\northopy.c1.jacobi.show(5, \"normal\", 0.0, 0.0)\n# plot, savefig also exist\n```\n\nRecurrence coefficients can be explicitly retrieved by\n```python\nimport orthopy\n\nrc = orthopy.c1.jacobi.RecurrenceCoefficients(\n \"monic\", # or \"classical\", \"normal\"\n alpha=0, beta=0, symbolic=True\n)\nprint(rc.p0)\nfor k in range(5):\n print(rc[k])\n```\n```\n1\n(1, 0, None)\n(1, 0, 1/3)\n(1, 0, 4/15)\n(1, 0, 9/35)\n(1, 0, 16/63)\n```\n\n\n### 1D half-space with weight function x<sup>\u03b1</sup> exp(-r)\n<img src=\"https://nschloe.github.io/orthopy/e1r.svg\" width=\"45%\">\n\n(Generalized) Laguerre polynomials.\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.e1r.Eval(x, alpha=0, scaling=\"normal\")\n```\n\n### 1D space with weight function exp(-r<sup>2</sup>)\n<img src=\"https://nschloe.github.io/orthopy/e1r2.svg\" width=\"45%\">\n\nHermite polynomials come in two standardizations:\n\n * `\"physicists\"` (against the weight function `exp(-x ** 2)`\n * `\"probabilists\"` (against the weight function `1 / sqrt(2 * pi) * exp(-x ** 2 / 2)`\n\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.e1r2.Eval(\n x,\n \"probabilists\", # or \"physicists\"\n \"normal\"\n)\n```\n\n#### Associated Legendre \"polynomials\"\n<img src=\"https://nschloe.github.io/orthopy/associated-legendre.svg\" width=\"45%\">\n\nNot all of those are polynomials, so they should really be called associated Legendre\n_functions_. The <i>k</i>th iteration contains _2k+1_ functions, indexed from\n_-k_ to _k_. (See the color grouping in the above plot.)\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.c1.associated_legendre.Eval(\n x, phi=None, standardization=\"natural\", with_condon_shortley_phase=True\n)\n```\n\n### Triangle (_T<sub>2</sub>_)\n<img src=\"https://nschloe.github.io/orthopy/triangle-tree.png\" width=\"40%\">\n\northopy's triangle orthogonal polynomials are evaluated in terms of [barycentric\ncoordinates](https://en.wikipedia.org/wiki/Barycentric_coordinate_system), so the\n`X.shape[0]` has to be 3.\n\n```python\nimport orthopy\n\nbary = [0.1, 0.7, 0.2]\nevaluator = orthopy.t2.Eval(bary, \"normal\")\n```\n\n\n### Disk (_S<sub>2</sub>_)\n\n<img src=\"https://nschloe.github.io/orthopy/disk-xu-tree.png\" width=\"70%\"> | <img src=\"https://nschloe.github.io/orthopy/disk-zernike-tree.png\" width=\"70%\"> | <img src=\"https://nschloe.github.io/orthopy/disk-zernike2-tree.png\" width=\"70%\">\n:------------:|:-----------------:|:-----------:|\nXu | [Zernike](https://en.wikipedia.org/wiki/Zernike_polynomials) | Zernike 2 |\n\northopy contains several families of orthogonal polynomials on the unit disk: After\n[Xu](https://arxiv.org/abs/1701.02709),\n[Zernike](https://en.wikipedia.org/wiki/Zernike_polynomials), and a simplified version\nof Zernike polynomials.\n\n```python\nimport orthopy\n\nx = [0.1, -0.3]\n\nevaluator = orthopy.s2.xu.Eval(x, \"normal\")\n# evaluator = orthopy.s2.zernike.Eval(x, \"normal\")\n# evaluator = orthopy.s2.zernike2.Eval(x, \"normal\")\n```\n\n\n### Sphere (_U<sub>3</sub>_)\n\n<img src=\"https://nschloe.github.io/orthopy/sph-tree.png\" width=\"50%\">\n\nComplex-valued _spherical harmonics,_ plotted with\n[cplot](https://github.com/nschloe/cplot/) coloring (black=zero, green=real positive,\npink=real negative, blue=imaginary positive, yellow=imaginary negative). The functions\nin the middle are real-valued. The complex angle takes _n_ turns on the <i>n</i>th\nlevel.\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.u3.EvalCartesian(\n x,\n scaling=\"quantum mechanic\" # or \"acoustic\", \"geodetic\", \"schmidt\"\n)\n\nevaluator = orthopy.u3.EvalSpherical(\n theta_phi, # polar, azimuthal angles\n scaling=\"quantum mechanic\" # or \"acoustic\", \"geodetic\", \"schmidt\"\n)\n```\nTo generate the above plot, write the tree mesh to a file\n```python\nimport orthopy\n\northopy.u3.write_tree(\"u3.vtk\", 5, \"quantum mechanic\")\n```\nand open it with [ParaView](https://www.paraview.org/). Select the _srgb1_ data set and\nturn off _Map Scalars_.\n\n### _n_-Cube (_C<sub>n</sub>_)\n\n<img src=\"https://nschloe.github.io/orthopy/c1.svg\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/c2.png\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/c3.png\" width=\"100%\">\n:-------------------------:|:------------------:|:---------------:|\nC<sub>1</sub> (Legendre) | C<sub>2</sub> | C<sub>3</sub> |\n\nJacobi product polynomials.\nAll polynomials are normalized on the n-dimensional cube. The dimensionality is\ndetermined by `X.shape[0]`.\n\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.cn.Eval(X, alpha=0, beta=0)\nvalues, degrees = next(evaluator)\n```\n\n### <i>n</i>D space with weight function exp(-r<sup>2</sup>) (_E<sub>n</sub><sup>r<sup>2</sup></sup>_)\n\n<img src=\"https://nschloe.github.io/orthopy/e1r2.svg\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/e2r2.png\" width=\"100%\"> | <img src=\"https://nschloe.github.io/orthopy/e3r2.png\" width=\"100%\">\n:-------------------------:|:------------------:|:---------------:|\n_E<sub>1</sub><sup>r<sup>2</sup></sup>_ | _E<sub>2</sub><sup>r<sup>2</sup></sup>_ | _E<sub>3</sub><sup>r<sup>2</sup></sup>_ |\n\nHermite product polynomials.\nAll polynomials are normalized over the measure. The dimensionality is determined by\n`X.shape[0]`.\n\n<!--pytest-codeblocks:skip-->\n```python\nevaluator = orthopy.enr2.Eval(\n x,\n standardization=\"probabilists\" # or \"physicists\"\n)\nvalues, degrees = next(evaluator)\n```\n\n\n### Other tools\n\n * Generating recurrence coefficients for 1D domains with\n [Stieltjes](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#stieltjes),\n [Golub-Welsch](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#golub-welsch),\n [Chebyshev](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#chebyshev), and\n [modified\n Chebyshev](https://github.com/nschloe/orthopy/wiki/Generating-1D-recurrence-coefficients-for-a-given-weight#modified-chebyshev).\n\n * The the sanity of recurrence coefficients with test 3 from [Gautschi's article](https://doi.org/10.1007/BF02218441):\n computing the weighted sum of orthogonal polynomials:\n <!--pytest-codeblocks:skip-->\n ```python\n orthopy.tools.gautschi_test_3(moments, alpha, beta)\n ```\n\n * [Clenshaw algorithm](https://en.wikipedia.org/wiki/Clenshaw_algorithm) for\n computing the weighted sum of orthogonal polynomials:\n <!--pytest-codeblocks:skip-->\n ```python\n vals = orthopy.c1.clenshaw(a, alpha, beta, t)\n ```\n\n\n### Installation\n\northopy is [available from the Python Package\nIndex](https://pypi.python.org/pypi/orthopy/), so use\n```\npip install orthopy\n```\nto install.\n\n### Testing\n\nTo run the tests, simply check out this repository and run\n```\npytest\n```\n\n### Relevant publications\n\n* [Robert C. Kirby, Singularity-free evaluation of collapsed-coordinate orthogonal polynomials, ACM Transactions on Mathematical Software (TOMS), Volume 37, Issue 1, January 2010](https://doi.org/10.1145/1644001.1644006)\n* [Abedallah Rababah, Recurrence Relations for Orthogonal Polynomials on Triangular Domains, MDPI Mathematics 2016, 4(2)](https://doi.org/10.3390/math4020025)\n* [Yuan Xu, Orthogonal polynomials of several variables, arxiv.org, January 2017](https://arxiv.org/abs/1701.02709)\n\n### License\nThis software is published under the [GPLv3 license](https://www.gnu.org/licenses/gpl-3.0.en.html).\n\n\n",
"bugtrack_url": null,
"license": "GPL-3.0-or-later",
"summary": "The last GPL version of orthopy, published as is",
"version": "0.9.5",
"project_urls": {
"GitHub": "https://github.com/nschloe/orthopy",
"Homepage": "https://zenodo.org/record/5751934",
"PyPI": "https://pypi.org/project/orthopy/"
},
"split_keywords": [],
"urls": [
{
"comment_text": "",
"digests": {
"blake2b_256": "ee7fcc9fa2170eb120056644d9764eb39ba284dc4a9d5d61943aeabea2410037",
"md5": "4b49e1da5ec0e8c4079f64e96711c7b7",
"sha256": "8715d4860acf67a48a5487446aeaada8e37d0e5345f88dff16a8d05a0e82483c"
},
"downloads": -1,
"filename": "orthopy_gpl-0.9.5-py3-none-any.whl",
"has_sig": false,
"md5_digest": "4b49e1da5ec0e8c4079f64e96711c7b7",
"packagetype": "bdist_wheel",
"python_version": "py3",
"requires_python": ">=3.7",
"size": 51166,
"upload_time": "2023-05-30T13:23:56",
"upload_time_iso_8601": "2023-05-30T13:23:56.929450Z",
"url": "https://files.pythonhosted.org/packages/ee/7f/cc9fa2170eb120056644d9764eb39ba284dc4a9d5d61943aeabea2410037/orthopy_gpl-0.9.5-py3-none-any.whl",
"yanked": false,
"yanked_reason": null
},
{
"comment_text": "",
"digests": {
"blake2b_256": "578b43ba088d37a2d4bc212cb02bcfbba7e72a9182a7459730c5cfc304a99fb5",
"md5": "c21ce97bed088f02c6d28ab1a0a903a3",
"sha256": "c1f902bc4150e773c00b67c57389f78c638ee0f97238d3511071fcea93407200"
},
"downloads": -1,
"filename": "orthopy-gpl-0.9.5.tar.gz",
"has_sig": false,
"md5_digest": "c21ce97bed088f02c6d28ab1a0a903a3",
"packagetype": "sdist",
"python_version": "source",
"requires_python": ">=3.7",
"size": 43256,
"upload_time": "2023-05-30T13:24:00",
"upload_time_iso_8601": "2023-05-30T13:24:00.207103Z",
"url": "https://files.pythonhosted.org/packages/57/8b/43ba088d37a2d4bc212cb02bcfbba7e72a9182a7459730c5cfc304a99fb5/orthopy-gpl-0.9.5.tar.gz",
"yanked": false,
"yanked_reason": null
}
],
"upload_time": "2023-05-30 13:24:00",
"github": true,
"gitlab": false,
"bitbucket": false,
"codeberg": false,
"github_user": "nschloe",
"github_project": "orthopy",
"travis_ci": false,
"coveralls": false,
"github_actions": false,
"lcname": "orthopy-gpl"
}