API

Implementation of the meta-class BetaEnsemble see \(\beta-\)Ensembles with children:

Such objects have 4 main methods:

  • sample_full_model()

  • sample_banded_model()

  • plot() to display a scatter plot of the last sample and eventually the limiting distribution (after normalization)

  • hist() to display a histogram of the last sample and eventually the limiting distribution (after normalization)

class dppy.beta_ensembles.BetaEnsemble(beta=2)[source]

Bases: object

\(\beta\)-Ensemble object parametrized by

Parameters

beta (int, float, default \(2\)) –

\(\beta >= 0\) inverse temperature parameter.

The default beta\(=2\) corresponds to the DPP case, see Orthogonal Polynomial Ensembles

See also

flush_samples()[source]

Empty the list_of_samples attribute.

abstract hist()[source]

Display histogram of the last realization of the underlying \(\beta\)-Ensemble. For some \(\beta\)-Ensembles, a normalization argument is available to display the limiting (or equilibrium) distribution and scale the points accordingly.

abstract normalize_points()[source]

Normalize points ormalization argument is available to display the limiting (or equilibrium) distribution and scale the points accordingly.

abstract plot()[source]

Display last realization of the underlying \(\beta\)-Ensemble. For some \(\beta\)-Ensembles, a normalization argument is available to display the limiting (or equilibrium) distribution and scale the points accordingly.

abstract sample_banded_model()[source]

Sample from underlying \(\beta\)-Ensemble using the corresponding banded matrix model. Arguments are the associated reference measure’s parameters, or the matrix dimensions used in sample_full_model()

abstract sample_full_model()[source]

Sample from underlying \(\beta\)-Ensemble using the corresponding full matrix model. Arguments are the associated matrix dimensions

class dppy.beta_ensembles.CircularEnsemble(beta=2)[source]

Bases: dppy.beta_ensembles.BetaEnsemble

Circular Ensemble object

See also

flush_samples()

Empty the list_of_samples attribute.

hist(normalization=True)[source]

Display the histogram of the angles \(\theta_{1}, \dots, \theta_{N}\) associated to the last realization \(\left\{ e^{i \theta_{1}}, \dots, e^{i \theta_{N}} \right\}\) object.

Parameters

normalization (bool, default True) – When True, the limiting distribution of the angles, i.e., the uniform distribution in \([0, 2\pi]\) is displayed

See also

normalize_points(points)[source]

No need to renormalize the points

plot(normalization=True)[source]

Display the last realization of the CircularEnsemble object.

Parameters

normalization (bool, default True) – When True, the unit circle is displayed

See also

sample_banded_model(size_N=10, random_state=None)[source]

Sample from Quindiagonal matrix model associated to the Circular Ensemble. Available for beta \(\in\mathbb{N}^*\), and the degenerate case beta \(=0\) corresponding to i.i.d. uniform points on the unit circle

Parameters

size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized

Note

To compare sample_banded_model() with sample_full_model() simply use the size_N parameter.

See also

sample_full_model(size_N=10, haar_mode='Hermite', random_state=None)[source]

Sample from tridiagonal matrix model associated to the Circular ensemble. Only available for beta \(\in\{1, 2, 4\}\) and the degenerate case beta \(=0\) corresponding to i.i.d. uniform points on the unit circle

Parameters
  • size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized

  • haar_mode (str, default 'hermite') – Sample Haar measure i.e. uniformly on the orthogonal/unitary/symplectic group using: - ‘QR’, - ‘Hermite’

See also

class dppy.beta_ensembles.GinibreEnsemble(beta=2)[source]

Bases: dppy.beta_ensembles.BetaEnsemble

Ginibre Ensemble object

See also

flush_samples()

Empty the list_of_samples attribute.

hist(normalization=True)[source]

Display the histogram of the radius of the points the last realization of the GinibreEnsemble object

Parameters

normalization (bool, default True) – When True, the points are first normalized so as to concentrate in the unit disk (see normalize_points()) and the limiting density \(2r 1_{[0, 1]}(r)\) of the radii is displayed

See also

normalize_points(points)[source]

Normalize points to concentrate in the unit disk.

\[x \mapsto \frac{x}{\sqrt{N}}\]

See also

plot(normalization=True)[source]

Display the last realization of the GinibreEnsemble object

Parameters

normalization (bool, default True) – When True, the points are first normalized so as to concentrate in the unit disk (see normalize_points()) and the unit circle is displayed

See also

sample_banded_model(*args, **kwargs)[source]

No banded model is known for Ginibre, use sample_full_model()

sample_full_model(size_N=10, random_state=None)[source]

Sample from full matrix model associated to the Ginibre ensemble. Only available for beta \(=2\)

Parameters

size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized

See also

class dppy.beta_ensembles.HermiteEnsemble(beta=2)[source]

Bases: dppy.beta_ensembles.BetaEnsemble

Hermite Ensemble object

See also

flush_samples()

Empty the list_of_samples attribute.

hist(normalization=True)[source]

Display the histogram of the last realization of the HermiteEnsemble object.

Parameters

normalization (bool, default True) –

When True, the points are first normalized (see normalize_points()) so that they concentrate as

  • If beta \(=0\), the \(\mathcal{N}(0, 2)\) reference measure associated to full full matrix model

  • If beta \(>0\), the limiting distribution, i.e., the semi-circle distribution

in both cases, the corresponding p.d.f. is displayed

See also

normalize_points(points)[source]

Normalize points obtained after sampling to fit the limiting distribution, i.e., the semi-circle

\[f(x) = \frac{1}{2\pi} \sqrt{4-x^2}\]
Parameters

points (array_like) – A sample from Hermite ensemble, accessible through the list_of_samples attribute

  • If sampled using sample_banded_model() with reference measure \(\mathcal{N}(\mu,\sigma^2)\)

    1. Normalize the points to fit the p.d.f. of \(\mathcal{N}(0,2)\) reference measure of the full matrix model

      \[x \mapsto \sqrt{2}\frac{x-\mu}{\sigma}\]
    2. If beta \(>0\), normalize the points to fit the semi-circle distribution

      \[x \mapsto \frac{x}{\beta N}\]

      Otherwise if beta \(=0\) do nothing more

  • If sampled using sample_full_model(), apply 2. above

Note

This method is called in plot() and hist() when normalization=True

plot(normalization=True)[source]

Display the last realization of the HermiteEnsemble object

Parameters

normalization (bool, default True) –

When True, the points are first normalized (see normalize_points()) so that they concentrate as

  • If beta \(=0\), the \(\mathcal{N}(0, 2)\) reference measure associated to full full matrix model

  • If beta \(>0\), the limiting distribution, i.e., the semi-circle distribution

in both cases, the corresponding p.d.f. is displayed

See also

sample_banded_model(loc=0.0, scale=1.4142135623730951, size_N=10, random_state=None)[source]

Sample from tridiagonal matrix model associated to the Hermite Ensemble. Available for beta \(>0\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the Gaussian \(\mathcal{N}(\mu,\sigma^2)\) reference measure

Parameters
  • loc (float, default \(0\)) – Mean \(\mu\) of the Gaussian \(\mathcal{N}(\mu, \sigma^2)\)

  • scale (float, default \(\sqrt{2}\)) – Standard deviation \(\sigma\) of the Gaussian \(\mathcal{N}(\mu, \sigma^2)\)

  • size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized

Note

The reference measure associated with the full matrix model is \(\mathcal{N}(0,2)\). For this reason, in the sampling_params attribute, the default values are set to loc\(=0\) and scale\(=\sqrt{2}\).

To compare sample_banded_model() with sample_full_model() simply use the size_N parameter.

See also

sample_full_model(size_N=10, random_state=None)[source]

Sample from full matrix model associated to the Hermite ensemble. Only available for beta \(\in\{1, 2, 4\}\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the Gaussian \(\mathcal{N}(\mu,\sigma^2)\) reference measure

Parameters

size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized

Note

The reference measure associated with the full matrix model is \(\mathcal{N}(0,2)\). For this reason, in the sampling_params attribute, the values of the parameters are set to loc\(=0\) and scale\(=\sqrt{2}\).

To compare sample_banded_model() with sample_full_model() simply use the size_N parameter.

See also

class dppy.beta_ensembles.JacobiEnsemble(beta=2)[source]

Bases: dppy.beta_ensembles.BetaEnsemble

Jacobi Ensemble object

See also

flush_samples()

Empty the list_of_samples attribute.

hist(normalization=True)[source]

Display the histogram of the last realization of the JacobiEnsemble object.

Parameters

normalization (bool, default True) –

When True

  • If beta \(=0\), display the p.d.f. of the \(\operatorname{Beta}(a, b)\)

  • If beta \(>0\), display the limiting distribution, i.e., the Wachter distribution

See also

normalize_points(points)[source]

No need to renormalize the points

plot(normalization=True)[source]

Display the last realization of the JacobiEnsemble object

Parameters

normalization (bool, default True) –

When True

  • If beta \(=0\), display the p.d.f. of the \(\operatorname{Beta}(a, b)\)

  • If beta \(>0\), display the limiting distribution, i.e., the Wachter distribution

See also

sample_banded_model(a=1.0, b=2.0, size_N=10, size_M1=None, size_M2=None, random_state=None)[source]

Sample from tridiagonal matrix model associated to the Jacobi ensemble. Available for beta \(>0\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the \(\operatorname{Beta}(a,b)\) reference measure

Parameters
  • shape (float, default \(1\)) – Shape parameter \(k\) of \(\Gamma(k, \theta)\) reference measure

  • scale (float, default \(2.0\)) – Scale parameter \(\theta\) of \(\Gamma(k, \theta)\) reference measure

  • size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized. Equivalent to the first dimension \(N\) of the matrices used in the full matrix model.

  • size_M1 (int) – Equivalent to the second dimension \(M_1\) of the first matrix used in the full matrix model.

  • size_M2 (int) – Equivalent to the second dimension \(M_2\) of the second matrix used in the full matrix model.

Note

The reference measure associated with the full matrix model is :

\[\operatorname{Beta}\left(\frac{\beta}{2}(M_1-N+1), \frac{\beta}{2}(M_2-N+1)\right)\]

For this reason, in the sampling_params attribute, the values of the parameters are set to a\(=\frac{\beta}{2}(M_1-N+1)\) and b\(=\frac{\beta}{2}(M_2-N+1)\).

To compare sample_banded_model() with sample_full_model() simply use the size_N, size_M2 and size_M2 parameters.

  • If size_M1 and size_M2 are not provided:

    In the sampling_params attribute, size_M1,2 are set to size_M1\(= \frac{2a}{\beta} + N - 1\) and size_M2\(= \frac{2b}{\beta} + N - 1\), to give an idea of the corresponding second dimensions \(M_{1,2}\).

  • If size_M1 and size_M2 are provided:

    In the sampling_params attribute, a and b are set to: a\(=\frac{\beta}{2}(M_1-N+1)\) and b\(=\frac{\beta}{2}(M_2-N+1)\).

See also

sample_full_model(size_N=100, size_M1=150, size_M2=200, random_state=None)[source]

Sample from full matrix model associated to the Jacobi ensemble. Only available for beta \(\in\{1, 2, 4\}\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the \(\operatorname{Beta}(a,b)\) reference measure

Parameters
  • size_N (int, default \(100\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized. First dimension of the matrix used to form the covariance matrix to be diagonalized, see full matrix model.

  • size_M1 (int, default \(150\)) – Second dimension \(M_1\) of the first matrix used to form the matrix to be diagonalized, see full matrix model.

  • size_M2 (int, default \(200\)) – Second dimension \(M_2\) of the second matrix used to form the matrix to be diagonalized, see full matrix model.

Note

The reference measure associated with the full matrix model is

\[\operatorname{Beta}\left(\frac{\beta}{2}(M_1-N+1), \frac{\beta}{2}(M_2-N+1)\right)\]

For this reason, in the sampling_params attribute, the values of the parameters are set to a\(=\frac{\beta}{2}(M_1-N+1)\) and b\(=\frac{\beta}{2}(M_2-N+1)\).

To compare sample_banded_model() with sample_full_model() simply use the size_N, size_M2 and size_M2 parameters.

See also

class dppy.beta_ensembles.LaguerreEnsemble(beta=2)[source]

Bases: dppy.beta_ensembles.BetaEnsemble

Laguerre Ensemble object

See also

flush_samples()

Empty the list_of_samples attribute.

hist(normalization=True)[source]

Display the histogram of the last realization of the LaguerreEnsemble object.

Parameters

normalization (bool, default True) –

When True, the points are first normalized (see normalize_points()) so that they concentrate as

  • If beta \(=0\), the \(\Gamma(k, 2)\) reference measure associated to full full matrix model

  • If beta \(>0\), the limiting distribution, i.e., the Marcenko-Pastur distribution

in both cases the corresponding p.d.f. is displayed

See also

normalize_points(points)[source]

Normalize points obtained after sampling to fit the limiting distribution, i.e., the Marcenko-Pastur distribution

\[\frac{1}{2\pi} \frac{\sqrt{(\lambda_+-x)(x-\lambda_-)}}{cx} 1_{[\lambda_-,\lambda_+]} dx\]

where \(c = \frac{M}{N}\) and \(\lambda_\pm = (1\pm\sqrt{c})^2\)

Parameters

points (array_like) – A sample from Laguerre ensemble, accessible through the list_of_samples attribute

  • If sampled using sample_banded_model() with reference measure \(\Gamma(k,\theta)\)

    \[x \mapsto \frac{2x}{\theta} \quad \text{and} \quad x \mapsto \frac{x}{\beta M}\]
  • If sampled using sample_full_model()

    \[x \mapsto \frac{x}{\beta M}\]

Note

This method is called in plot() and hist() when normalization=True.

plot(normalization=True)[source]

Display the last realization of the LaguerreEnsemble object

Parameters

normalization (bool, default True) –

When True, the points are first normalized (see normalize_points()) so that they concentrate as

  • If beta \(=0\), the \(\Gamma(k, 2)\) reference measure associated to full full matrix model

  • If beta \(>0\), the limiting distribution, i.e., the Marcenko-Pastur distribution

in both cases the corresponding p.d.f. is displayed

See also

sample_banded_model(shape=1.0, scale=2.0, size_N=10, size_M=None, random_state=None)[source]

Sample from tridiagonal matrix model associated to the Laguerre ensemble. Available for beta \(>0\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the \(\Gamma(k,\theta)\) reference measure

Parameters
  • shape (float, default \(1\)) – Shape parameter \(k\) of \(\Gamma(k, \theta)\) reference measure

  • scale (float, default \(2.0\)) – Scale parameter \(\theta\) of \(\Gamma(k, \theta)\) reference measure

  • size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized. Equivalent to the first dimension \(N\) of the matrix used to form the covariance matrix in the full matrix model.

  • size_M (int, default None) – Equivalent to the second dimension \(M\) of the matrix used to form the covariance matrix in the full matrix model.

  • If size_M is not provided:

    In the sampling_params, size_M is set to size_M\(= \frac{2k}{\beta} + N - 1\), to give an idea of the corresponding second dimension \(M\).

  • If size_M is provided:

    In the sampling_params, shape and scale are set to: shape\(=\frac{1}{2} \beta (M-N+1)\) and scale\(=2\)

Note

The reference measure associated with the full matrix model is \(\Gamma\left(\frac{\beta}{2}(M-N+1), 2\right)\). This explains the role of the size_M parameter.

To compare sample_banded_model() with sample_full_model() simply use the size_N and size_M parameters

See also

sample_full_model(size_N=10, size_M=100, random_state=None)[source]

Sample from full matrix model associated to the Laguerre ensemble. Only available for beta \(\in\{1, 2, 4\}\) and the degenerate case beta \(=0\) corresponding to i.i.d. points from the \(\Gamma(k,\theta)\) reference measure

Parameters
  • size_N (int, default \(10\)) – Number \(N\) of points, i.e., size of the matrix to be diagonalized. First dimension of the matrix used to form the covariance matrix to be diagonalized, see full matrix model.

  • size_M (int, default \(100\)) – Second dimension \(M\) of the matrix used to form the covariance matrix to be diagonalized, see full matrix model.

Note

The reference measure associated with the full matrix model is \(\Gamma\left(\frac{\beta}{2}(M-N+1), 2\right)\). For this reason, in the sampling_params, the values of the parameters are set to shape\(=\frac{\beta}{2}(M-N+1)\) and scale\(=2\).

To compare sample_banded_model() with sample_full_model() simply use the size_N and size_M parameters.

See also