Skip to content

Commit

Permalink
ENH: constants: add array api support (#20593)
Browse files Browse the repository at this point in the history
Co-authored-by: Ralf Gommers <ralf.gommers@gmail.com>
  • Loading branch information
j-bowhay and rgommers committed Apr 29, 2024
1 parent 9ae5550 commit 4ee0de3
Show file tree
Hide file tree
Showing 5 changed files with 116 additions and 47 deletions.
1 change: 1 addition & 0 deletions .github/workflows/array_api.yml
Expand Up @@ -91,6 +91,7 @@ jobs:
export OMP_NUM_THREADS=2
# expand as new modules are added
python dev.py --no-build test -b all -s cluster -- --durations 3 --timeout=60
python dev.py --no-build test -b all -s constants -- --durations 3 --timeout=60
python dev.py --no-build test -b all -s fft -- --durations 3 --timeout=60
python dev.py --no-build test -b all -t scipy.special.tests.test_support_alternative_backends -- --durations 3 --timeout=60
python dev.py --no-build test -b all -t scipy._lib.tests.test_array_api
Expand Down
1 change: 1 addition & 0 deletions doc/source/dev/api-dev/array_api.rst
Expand Up @@ -95,6 +95,7 @@ variable is set:

- `scipy.cluster.hierarchy`
- `scipy.cluster.vq`
- `scipy.constants`
- `scipy.fft`

Support is provided in `scipy.special` for the following functions:
Expand Down
13 changes: 10 additions & 3 deletions scipy/_lib/_array_api.py
Expand Up @@ -120,9 +120,11 @@ def array_namespace(*arrays):


def _asarray(
array, dtype=None, order=None, copy=None, *, xp=None, check_finite=False
array, dtype=None, order=None, copy=None, *, xp=None, check_finite=False,
subok=False
):
"""SciPy-specific replacement for `np.asarray` with `order` and `check_finite`.
"""SciPy-specific replacement for `np.asarray` with `order`, `check_finite`, and
`subok`.
Memory layout parameter `order` is not exposed in the Array API standard.
`order` is only enforced if the input array implementation
Expand All @@ -131,13 +133,18 @@ def _asarray(
`check_finite` is also not a keyword in the array API standard; included
here for convenience rather than that having to be a separate function
call inside SciPy functions.
`subok` is included to allow this function to preserve the behaviour of
`np.asanyarray` for NumPy based inputs.
"""
if xp is None:
xp = array_namespace(array)
if xp.__name__ in {"numpy", "scipy._lib.array_api_compat.numpy"}:
# Use NumPy API to support order
if copy is True:
array = np.array(array, order=order, dtype=dtype)
array = np.array(array, order=order, dtype=dtype, subok=subok)
elif subok:
array = np.asanyarray(array, order=order, dtype=dtype)
else:
array = np.asarray(array, order=order, dtype=dtype)

Expand Down
32 changes: 19 additions & 13 deletions scipy/constants/_constants.py
Expand Up @@ -13,11 +13,13 @@
from typing import TYPE_CHECKING, Any

from ._codata import value as _cd
import numpy as np

if TYPE_CHECKING:
import numpy.typing as npt

from scipy._lib._array_api import array_namespace, _asarray


"""
BasSw 2006
physical constants: imported from CODATA
Expand Down Expand Up @@ -269,19 +271,21 @@ def convert_temperature(
array([ 233.15, 313.15])
"""
xp = array_namespace(val)
_val = _asarray(val, xp=xp, subok=True)
# Convert from `old_scale` to Kelvin
if old_scale.lower() in ['celsius', 'c']:
tempo = np.asanyarray(val) + zero_Celsius
tempo = _val + zero_Celsius
elif old_scale.lower() in ['kelvin', 'k']:
tempo = np.asanyarray(val)
tempo = _val
elif old_scale.lower() in ['fahrenheit', 'f']:
tempo = (np.asanyarray(val) - 32) * 5 / 9 + zero_Celsius
tempo = (_val - 32) * 5 / 9 + zero_Celsius
elif old_scale.lower() in ['rankine', 'r']:
tempo = np.asanyarray(val) * 5 / 9
tempo = _val * 5 / 9
else:
raise NotImplementedError("%s scale is unsupported: supported scales "
"are Celsius, Kelvin, Fahrenheit, and "
"Rankine" % old_scale)
raise NotImplementedError(f"{old_scale=} is unsupported: supported scales "
"are Celsius, Kelvin, Fahrenheit, and "
"Rankine")
# and from Kelvin to `new_scale`.
if new_scale.lower() in ['celsius', 'c']:
res = tempo - zero_Celsius
Expand All @@ -292,9 +296,9 @@ def convert_temperature(
elif new_scale.lower() in ['rankine', 'r']:
res = tempo * 9 / 5
else:
raise NotImplementedError("'%s' scale is unsupported: supported "
"scales are 'Celsius', 'Kelvin', "
"'Fahrenheit', and 'Rankine'" % new_scale)
raise NotImplementedError(f"{new_scale=} is unsupported: supported "
"scales are 'Celsius', 'Kelvin', "
"'Fahrenheit', and 'Rankine'")

return res

Expand Down Expand Up @@ -329,7 +333,8 @@ def lambda2nu(lambda_: npt.ArrayLike) -> Any:
array([ 2.99792458e+08, 1.00000000e+00])
"""
return c / np.asanyarray(lambda_)
xp = array_namespace(lambda_)
return c / _asarray(lambda_, xp=xp, subok=True)


def nu2lambda(nu: npt.ArrayLike) -> Any:
Expand Down Expand Up @@ -359,4 +364,5 @@ def nu2lambda(nu: npt.ArrayLike) -> Any:
array([ 2.99792458e+08, 1.00000000e+00])
"""
return c / np.asanyarray(nu)
xp = array_namespace(nu)
return c / _asarray(nu, xp=xp, subok=True)
116 changes: 85 additions & 31 deletions scipy/constants/tests/test_constants.py
@@ -1,35 +1,89 @@
from numpy.testing import assert_equal, assert_allclose
import pytest

import scipy.constants as sc
from scipy.conftest import array_api_compatible
from scipy._lib._array_api import xp_assert_equal, xp_assert_close


pytestmark = [array_api_compatible, pytest.mark.usefixtures("skip_xp_backends")]
skip_xp_backends = pytest.mark.skip_xp_backends


class TestConvertTemperature:
def test_convert_temperature(self, xp):
xp_assert_equal(sc.convert_temperature(xp.asarray(32.), 'f', 'Celsius'),
xp.asarray(0.0))
xp_assert_equal(sc.convert_temperature(xp.asarray([0., 0.]),
'celsius', 'Kelvin'),
xp.asarray([273.15, 273.15]))
xp_assert_equal(sc.convert_temperature(xp.asarray([0., 0.]), 'kelvin', 'c'),
xp.asarray([-273.15, -273.15]))
xp_assert_equal(sc.convert_temperature(xp.asarray([32., 32.]), 'f', 'k'),
xp.asarray([273.15, 273.15]))
xp_assert_equal(sc.convert_temperature(xp.asarray([273.15, 273.15]),
'kelvin', 'F'),
xp.asarray([32., 32.]))
xp_assert_equal(sc.convert_temperature(xp.asarray([0., 0.]), 'C', 'fahrenheit'),
xp.asarray([32., 32.]))
xp_assert_close(sc.convert_temperature(xp.asarray([0., 0.], dtype=xp.float64),
'c', 'r'),
xp.asarray([491.67, 491.67], dtype=xp.float64),
rtol=0., atol=1e-13)
xp_assert_close(sc.convert_temperature(xp.asarray([491.67, 491.67],
dtype=xp.float64),
'Rankine', 'C'),
xp.asarray([0., 0.], dtype=xp.float64), rtol=0., atol=1e-13)
xp_assert_close(sc.convert_temperature(xp.asarray([491.67, 491.67],
dtype=xp.float64),
'r', 'F'),
xp.asarray([32., 32.], dtype=xp.float64), rtol=0., atol=1e-13)
xp_assert_close(sc.convert_temperature(xp.asarray([32., 32.], dtype=xp.float64),
'fahrenheit', 'R'),
xp.asarray([491.67, 491.67], dtype=xp.float64),
rtol=0., atol=1e-13)
xp_assert_close(sc.convert_temperature(xp.asarray([273.15, 273.15],
dtype=xp.float64),
'K', 'R'),
xp.asarray([491.67, 491.67], dtype=xp.float64),
rtol=0., atol=1e-13)
xp_assert_close(sc.convert_temperature(xp.asarray([491.67, 0.],
dtype=xp.float64),
'rankine', 'kelvin'),
xp.asarray([273.15, 0.], dtype=xp.float64), rtol=0., atol=1e-13)

@skip_xp_backends(np_only=True, reasons=['Python list input uses NumPy backend'])
def test_convert_temperature_array_like(self):
xp_assert_close(sc.convert_temperature([491.67, 0.], 'rankine', 'kelvin'),
[273.15, 0.], rtol=0., atol=1e-13)


@skip_xp_backends(np_only=True, reasons=['Python int input uses NumPy backend'])
def test_convert_temperature_errors(self, xp):
with pytest.raises(NotImplementedError, match="old_scale="):
sc.convert_temperature(1, old_scale="cheddar", new_scale="kelvin")
with pytest.raises(NotImplementedError, match="new_scale="):
sc.convert_temperature(1, old_scale="kelvin", new_scale="brie")


class TestLambdaToNu:
def test_lambda_to_nu(self, xp):
xp_assert_equal(sc.lambda2nu(xp.asarray([sc.speed_of_light, 1])),
xp.asarray([1, sc.speed_of_light]))


@skip_xp_backends(np_only=True, reasons=['Python list input uses NumPy backend'])
def test_lambda_to_nu_array_like(self, xp):
xp_assert_equal(sc.lambda2nu([sc.speed_of_light, 1]),
[1, sc.speed_of_light])


class TestNuToLambda:
def test_nu_to_lambda(self, xp):
xp_assert_equal(sc.nu2lambda(xp.asarray([sc.speed_of_light, 1])),
xp.asarray([1, sc.speed_of_light]))

def test_convert_temperature():
assert_equal(sc.convert_temperature(32, 'f', 'Celsius'), 0)
assert_equal(sc.convert_temperature([0, 0], 'celsius', 'Kelvin'),
[273.15, 273.15])
assert_equal(sc.convert_temperature([0, 0], 'kelvin', 'c'),
[-273.15, -273.15])
assert_equal(sc.convert_temperature([32, 32], 'f', 'k'), [273.15, 273.15])
assert_equal(sc.convert_temperature([273.15, 273.15], 'kelvin', 'F'),
[32, 32])
assert_equal(sc.convert_temperature([0, 0], 'C', 'fahrenheit'), [32, 32])
assert_allclose(sc.convert_temperature([0, 0], 'c', 'r'), [491.67, 491.67],
rtol=0., atol=1e-13)
assert_allclose(sc.convert_temperature([491.67, 491.67], 'Rankine', 'C'),
[0., 0.], rtol=0., atol=1e-13)
assert_allclose(sc.convert_temperature([491.67, 491.67], 'r', 'F'),
[32., 32.], rtol=0., atol=1e-13)
assert_allclose(sc.convert_temperature([32, 32], 'fahrenheit', 'R'),
[491.67, 491.67], rtol=0., atol=1e-13)
assert_allclose(sc.convert_temperature([273.15, 273.15], 'K', 'R'),
[491.67, 491.67], rtol=0., atol=1e-13)
assert_allclose(sc.convert_temperature([491.67, 0.], 'rankine', 'kelvin'),
[273.15, 0.], rtol=0., atol=1e-13)


def test_lambda_to_nu():
assert_equal(sc.lambda2nu([sc.speed_of_light, 1]), [1, sc.speed_of_light])


def test_nu_to_lambda():
assert_equal(sc.nu2lambda([sc.speed_of_light, 1]), [1, sc.speed_of_light])
@skip_xp_backends(np_only=True, reasons=['Python list input uses NumPy backend'])
def test_nu_to_lambda_array_like(self, xp):
xp_assert_equal(sc.nu2lambda([sc.speed_of_light, 1]),
[1, sc.speed_of_light])

0 comments on commit 4ee0de3

Please sign in to comment.