Quadrature

The quadrature API contains the numerical integration kernels and basis machinery used by the solver internals.

Quadrature rules and kernels used by the R-matrix solver.

class jitr.quadrature.Kernel(nbasis, basis='Legendre')[source]

Bases: object

Convenience wrapper around a quadrature rule and its basis functions.

Parameters:
f(n, a, s)[source]

Evaluate the n-th Lagrange basis function.

Return type:

complex

Parameters:
radial_grid(radius)[source]

Return the one-dimensional quadrature grid on [0, radius].

Return type:

ndarray[tuple[Any, ...], dtype[double]]

Parameters:

radius (float)

nonlocal_radial_grids(radius)[source]

Return the tensor-product quadrature grids on [0, radius]^2.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

radius (float)

integrate_local(values, radius)[source]

Integrate local values on [0, radius] using Gauss quadrature.

Return type:

cdouble

Parameters:
  • values (ArrayLike)

  • radius (float)

double_integrate_nonlocal(values, radius, is_symmetric=True)[source]

Integrate nonlocal values on [0, radius] × [0, radius].

Return type:

cdouble

Parameters:
  • values (ArrayLike)

  • radius (float)

  • is_symmetric (bool)

fourier_bessel_transform(l, values, k, radius)[source]

Perform a Fourier-Bessel transform of order l.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • l (int)

  • values (ArrayLike)

  • k (ArrayLike)

  • radius (float)

double_fourier_bessel_transform(l, values, k, radius)[source]

Perform a double Fourier-Bessel transform of order l.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • l (int)

  • values (ArrayLike)

  • k (ArrayLike)

  • radius (float)

dwba_local(bra, ket, values)[source]

Return a DWBA matrix element for a local operator.

Return type:

cdouble

Parameters:
dwba_nonlocal(bra, ket, values, radius)[source]

Return a DWBA matrix element for a nonlocal operator.

Return type:

cdouble

Parameters:
matrix_local(values)[source]

Validate and cast local quadrature values.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:

values (ArrayLike)

matrix_nonlocal(values, radius)[source]

Validate and weight nonlocal quadrature values.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • values (ArrayLike)

  • radius (float)

class jitr.quadrature.LagrangeLaguerreQuadrature(abscissa, weights, overlap=None)[source]

Bases: object

Lagrange Laguerre mesh for the Schrödinger equation following ch. 3.3 of Baye, D. (2015). The Lagrange-mesh method. Physics reports, 565, 1-107, with the only difference being the domain is scaled in each channel; e.g r -> s_i = r * k_i, and each channel’s equation is then divided by it’s asymptotic kinetic energy in the channel T_i = E_inc - E_i

Parameters:
  • abscissa (FloatArray)

  • weights (FloatArray)

  • overlap (FloatArray | None)

kinetic_operator_element(n, m, a, l)[source]

Return the (n,m)th matrix element for the kinetic energy operator at channel radius a = k*r with orbital angular momentum l.

Return type:

float

Parameters:
kinetic_matrix(a, l)[source]

Return the kinetic operator matrix in the Lagrange Laguerre basis.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
class jitr.quadrature.LagrangeLegendreQuadrature(abscissa, weights, overlap=None)[source]

Bases: object

Lagrange Legendre mesh for the Schrödinger equation following ch. 3.4 of Baye, D. (2015). The Lagrange-mesh method. Physics reports, 565, 1-107, with the only difference being the domain is scaled in each channel; e.g. r -> s_i = r * k_i, and each channel’s equation is then divided by it’s asymptotic kinetic energy in the channel T_i = E_inc - E_i

Parameters:
  • abscissa (FloatArray)

  • weights (FloatArray)

  • overlap (FloatArray | None)

kinetic_operator_element(n, m, a, l)[source]

Return the (n,m)th matrix element for the kinetic energy + Bloch operator at channel radius a = k*r with orbital angular momentum l.

Return type:

float

Parameters:
kinetic_matrix(a, l)[source]

Return the kinetic operator matrix in the Lagrange Legendre basis.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
jitr.quadrature.generate_laguerre_quadrature(nbasis)[source]

Return zeros and weights for Gauss quadrature using the Lagrange-Laguerre basis. See Ch. 3.3 of Baye, 2015

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

nbasis (int)

jitr.quadrature.generate_legendre_quadrature(nbasis)[source]

Return zeros and weights for Gauss quadrature using the Lagrange-Legendre basis shifted and scaled onto [0,a]. See Ch. 3.4 of Baye, 2015

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

nbasis (int)

jitr.quadrature.laguerre(n, a, s, quadrature)[source]

nth Lagrange-Laguerre function, scaled by a. Eq. 3.70 in Baye, 2015 with alpha = 0.

Note: n is indexed from 1 (constant function is not part of basis)

Return type:

complex

Parameters:
jitr.quadrature.legendre(n, a, s, quadrature)[source]

nth Lagrange-Legendre polynomial shifted onto [0,a_i] and regularized by s. Eq. 3.122 in Baye, 2015

Note: n is indexed from 1 (constant function is not part of basis)

Return type:

complex

Parameters:

Quadrature kernels and transforms on Lagrange meshes.

class jitr.quadrature.kernel.Kernel(nbasis, basis='Legendre')[source]

Bases: object

Convenience wrapper around a quadrature rule and its basis functions.

Parameters:
f(n, a, s)[source]

Evaluate the n-th Lagrange basis function.

Return type:

complex

Parameters:
radial_grid(radius)[source]

Return the one-dimensional quadrature grid on [0, radius].

Return type:

ndarray[tuple[Any, ...], dtype[double]]

Parameters:

radius (float)

nonlocal_radial_grids(radius)[source]

Return the tensor-product quadrature grids on [0, radius]^2.

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

radius (float)

integrate_local(values, radius)[source]

Integrate local values on [0, radius] using Gauss quadrature.

Return type:

cdouble

Parameters:
  • values (ArrayLike)

  • radius (float)

double_integrate_nonlocal(values, radius, is_symmetric=True)[source]

Integrate nonlocal values on [0, radius] × [0, radius].

Return type:

cdouble

Parameters:
  • values (ArrayLike)

  • radius (float)

  • is_symmetric (bool)

fourier_bessel_transform(l, values, k, radius)[source]

Perform a Fourier-Bessel transform of order l.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • l (int)

  • values (ArrayLike)

  • k (ArrayLike)

  • radius (float)

double_fourier_bessel_transform(l, values, k, radius)[source]

Perform a double Fourier-Bessel transform of order l.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • l (int)

  • values (ArrayLike)

  • k (ArrayLike)

  • radius (float)

dwba_local(bra, ket, values)[source]

Return a DWBA matrix element for a local operator.

Return type:

cdouble

Parameters:
dwba_nonlocal(bra, ket, values, radius)[source]

Return a DWBA matrix element for a nonlocal operator.

Return type:

cdouble

Parameters:
matrix_local(values)[source]

Validate and cast local quadrature values.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:

values (ArrayLike)

matrix_nonlocal(values, radius)[source]

Validate and weight nonlocal quadrature values.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
  • values (ArrayLike)

  • radius (float)

Quadrature rules and Lagrange-mesh basis functions.

jitr.quadrature.quadrature.laguerre(n, a, s, quadrature)[source]

nth Lagrange-Laguerre function, scaled by a. Eq. 3.70 in Baye, 2015 with alpha = 0.

Note: n is indexed from 1 (constant function is not part of basis)

Return type:

complex

Parameters:
jitr.quadrature.quadrature.legendre(n, a, s, quadrature)[source]

nth Lagrange-Legendre polynomial shifted onto [0,a_i] and regularized by s. Eq. 3.122 in Baye, 2015

Note: n is indexed from 1 (constant function is not part of basis)

Return type:

complex

Parameters:
jitr.quadrature.quadrature.generate_laguerre_quadrature(nbasis)[source]

Return zeros and weights for Gauss quadrature using the Lagrange-Laguerre basis. See Ch. 3.3 of Baye, 2015

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

nbasis (int)

jitr.quadrature.quadrature.generate_legendre_quadrature(nbasis)[source]

Return zeros and weights for Gauss quadrature using the Lagrange-Legendre basis shifted and scaled onto [0,a]. See Ch. 3.4 of Baye, 2015

Return type:

tuple[ndarray[tuple[Any, ...], dtype[double]], ndarray[tuple[Any, ...], dtype[double]]]

Parameters:

nbasis (int)

class jitr.quadrature.quadrature.LagrangeLaguerreQuadrature(abscissa, weights, overlap=None)[source]

Bases: object

Lagrange Laguerre mesh for the Schrödinger equation following ch. 3.3 of Baye, D. (2015). The Lagrange-mesh method. Physics reports, 565, 1-107, with the only difference being the domain is scaled in each channel; e.g r -> s_i = r * k_i, and each channel’s equation is then divided by it’s asymptotic kinetic energy in the channel T_i = E_inc - E_i

Parameters:
  • abscissa (FloatArray)

  • weights (FloatArray)

  • overlap (FloatArray | None)

kinetic_operator_element(n, m, a, l)[source]

Return the (n,m)th matrix element for the kinetic energy operator at channel radius a = k*r with orbital angular momentum l.

Return type:

float

Parameters:
kinetic_matrix(a, l)[source]

Return the kinetic operator matrix in the Lagrange Laguerre basis.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters:
class jitr.quadrature.quadrature.LagrangeLegendreQuadrature(abscissa, weights, overlap=None)[source]

Bases: object

Lagrange Legendre mesh for the Schrödinger equation following ch. 3.4 of Baye, D. (2015). The Lagrange-mesh method. Physics reports, 565, 1-107, with the only difference being the domain is scaled in each channel; e.g. r -> s_i = r * k_i, and each channel’s equation is then divided by it’s asymptotic kinetic energy in the channel T_i = E_inc - E_i

Parameters:
  • abscissa (FloatArray)

  • weights (FloatArray)

  • overlap (FloatArray | None)

kinetic_operator_element(n, m, a, l)[source]

Return the (n,m)th matrix element for the kinetic energy + Bloch operator at channel radius a = k*r with orbital angular momentum l.

Return type:

float

Parameters:
kinetic_matrix(a, l)[source]

Return the kinetic operator matrix in the Lagrange Legendre basis.

Return type:

ndarray[tuple[Any, ...], dtype[cdouble]]

Parameters: