Data structures

Contents

Data structures#

class pyxel.data_structure.Scene[source]#

Bases: object

Scene class defining and storing information of all multi-wavelength photons (unit: ph / (cm2 nm s)).

Multi-wavelength photon information are store in form of xarray Datasets within a hierarchical structure.

Examples

Create an empty Scene

>>> from pyxel.detectors import CCD, CCDGeometry, Characteristics, Environment
>>> detector = CCD(
...     geometry=CCDGeometry(row=5, col=5),
...     characteristics=Characteristics(),
...     environment=Environment(),
... )
>>> detector.scene
Scene<no source>

Add a source

>>> import xarray as xr
>>> source = xr.Dataset(...)
>>> source
<xarray.Dataset>
Dimensions:     (ref: 345, wavelength: 343)
Coordinates:
  * ref         (ref) int64 0 1 2 3 4 5 6 7 ... 337 338 339 340 341 342 343 344
  * wavelength  (wavelength) float64 336.0 338.0 340.0 ... 1.018e+03 1.02e+03
Data variables:
    x           (ref) float64 2.057e+05 2.058e+05 ... 2.031e+05 2.03e+05
    y           (ref) float64 8.575e+04 8.58e+04 ... 8.795e+04 8.807e+04
    weight      (ref) float64 11.49 14.13 15.22 14.56 ... 15.21 11.51 8.727
    flux        (ref, wavelength) float64 0.03769 0.04137 ... 1.813 1.896
Attributes:
    right_ascension:  56.75 deg
    declination:      24.1167 deg
    fov_radius:       0.5 deg
>>> detector.scene.add_source(source)
>>> detector.scene
Scene<1 source(s)>

Get the pointing coordinates of the first source >>> detector.scene.get_pointing_coordinates() or >>> detector.scene.get_pointing_coordinates(source_idx=0) SceneCoordinates(right_ascension=<Quantity 56.75 deg>, declination=<Quantity 24.1167 deg>, fov=<Quantity 0.5 deg>)

add_source(source)[source]#

Add a source to the current scene.

Parameters:

source (Dataset)

Raises:
  • TypeError – If ‘source’ is not a Dataset object.

  • ValueError – If ‘source’ has not the expected format.

Examples

>>> from pyxel.detectors import CCD
>>> detector = CCD(...)
>>> source
<xarray.Dataset>
Dimensions:     (ref: 345, wavelength: 343)
Coordinates:
  * ref         (ref) int64 0 1 2 3 4 5 6 7 ... 337 338 339 340 341 342 343 344
  * wavelength  (wavelength) float64 336.0 338.0 340.0 ... 1.018e+03 1.02e+03
Data variables:
    x           (ref) float64 1.334e+03 1.434e+03 ... -1.271e+03 -1.381e+03
    y           (ref) float64 -1.009e+03 -956.1 -797.1 ... 1.195e+03 1.309e+03
    weight      (ref) float64 11.49 14.13 15.22 14.56 ... 15.21 11.51 8.727
    flux        (ref, wavelength) float64 0.003769 0.004137 ... 0.1813 0.1896
Attributes:
    right_ascension:  56.75 deg
    declination:      24.1167 deg
    fov_radius:       0.5 deg
>>> detector.scene.add_source(source)
>>> detector.scene.data
DataTree('scene', parent=None)
└── DataTree('list')
    └── DataTree('0')
            Dimensions:     (ref: 4, wavelength: 4)
            Coordinates:
              * ref         (ref) int64 0 1 2 3
              * wavelength  (wavelength) float64 336.0 338.0 1.018e+03 1.02e+03
            Data variables:
                x           (ref) float64 64.97 11.44 -55.75 -20.66
                y           (ref) float64 89.62 -129.3 -48.16 87.87
                weight      (ref) float64 14.73 12.34 14.63 14.27
                flux        (ref, wavelength) float64 0.003769 0.004137 ... 0.1813 0.1896
            Attributes:
                right_ascension:  56.75 deg
                declination:      24.1167 deg
                fov_radius:       0.5 deg
get_pointing_coordinates(source_idx=0)[source]#

Get the SceneCoordinates from a source in the scene.

property data#

Get a multi-wavelength object.

empty()[source]#

Create a new empty source.

from_scopesim(source)[source]#

Convert a ScopeSim Source object into a Scene object.

Parameters:

source (scopesim.Source) – Object to convert to a Scene object.

Raises:
  • RuntimeError – If package ‘scopesim’ is not installed.

  • TypeError – If input parameter ‘source’ is not a ScopeSim Source object.

Notes

More information about ScopeSim Source objects at this link: https://scopesim.readthedocs.io/en/latest/reference/scopesim.source.source.html

to_scopesim()[source]#

Convert this Scene object into a ScopeSim Source object.

Returns:

Source – A ScopeSim Source object.

Notes

More information about ScopeSim Source objects at this link: https://scopesim.readthedocs.io/en/latest/reference/scopesim.source.source.html

to_dict()[source]#

Convert an instance of Scene to a dict.

classmethod from_dict(dct)[source]#

Create a new instance of a Scene object from a dict.

to_xarray()[source]#

Convert current scene to a xarray Dataset.

Returns:

xr.Dataset

Examples

>>> ds = scene.to_xarray()
>>> ds
<xarray.Dataset>
Dimensions:     (ref: 345, wavelength: 343)
Coordinates:
  * ref         (ref) int64 0 1 2 3 4 5 6 7 ... 337 338 339 340 341 342 343 344
  * wavelength  (wavelength) float64 336.0 338.0 340.0 ... 1.018e+03 1.02e+03
Data variables:
    x           (ref) float64 2.057e+05 2.058e+05 ... 2.031e+05 2.03e+05
    y           (ref) float64 8.575e+04 8.58e+04 ... 8.795e+04 8.807e+04
    weight      (ref) float64 11.49 14.13 15.22 14.56 ... 15.21 11.51 8.727
    flux        (ref, wavelength) float64 0.03769 0.04137 ... 1.813 1.896
Attributes:
    right_ascension:  56.75 deg
    declination:      24.1167 deg
    fov_radius:       0.5 deg
>>> ds["wavelength"]
<xarray.DataArray 'wavelength' (wavelength: 343)>
array([ 336.,  338.,  340., ..., 1016., 1018., 1020.])
Coordinates:
  * wavelength  (wavelength) float64 336.0 338.0 340.0 ... 1.018e+03 1.02e+03
Attributes:
    units:    nm
>>> ds["flux"]
<xarray.DataArray 'flux' (ref: 345, wavelength: 343)>
array([[3.76907117e-02, 4.13740861e-02, ..., 3.98815404e-02, 7.96581117e-01],
       [1.15190254e-02, 1.02210366e-02, ..., 2.00486326e-02, 2.05518196e-02],
       ...,
       [1.01187592e-01, 9.57637374e-02, ..., 2.71410354e-01, 2.85997559e-01],
       [1.80093381e+00, 1.69864354e+00, ..., 1.81295134e+00, 1.89642359e+00]])
Coordinates:
  * ref         (ref) int64 0 1 2 3 4 5 6 7 ... 337 338 339 340 341 342 343 344
  * wavelength  (wavelength) float64 336.0 338.0 340.0 ... 1.018e+03 1.02e+03
Attributes:
    units:    ph / (cm2 nm s)
class pyxel.data_structure.Photon(geo)[source]#

Bases: object

Photon class designed to handle the storage of monochromatic (unit: ph or multi-wavelength photons (unit ph/nm).

Monochromatic photons are stored in a 2D Numpy array and multi-wavelength photons are stored in a 3D Xarray DataArray.

Accepted array types: np.float16, np.float32, np.float64

Examples

Create an empty Photon container

>>> import numpy as np
>>> from pyxel.detectors import CCD, CCDGeometry, Characteristics, Environment
>>> detector = CCD(
...     geometry=CCDGeometry(row=5, col=5),
...     characteristics=Characteristics(),
...     environment=Environment(),
... )
>>> detector.photon
Photon<UNINITIALIZED, shape=(5, 5)>
>>> detector.ndim
0
>>> detector.shape
()

Use monochromatic photons

>>> detector.photon.empty()
>>> detector.photon.array = np.zeros(shape=(5, 5), dtype=float)
>>> detector.photon
Photon<shape=(5, 5), dtype=float64>
>>> detector.photon.array = detector.photon.array + np.ones(shape=(5, 5))
>>> detector.photon += np.ones(shape=(5, 5))
>>> detector.photon.array
array([[2., ...., 2.],
       ...,
       [2., ..., 2.]])
>>> detector.photon.to_xarray()
<xarray.DataArray 'photon' (y: 5, x: 5)>
array([[2., ...., 2.],
       ...,
       [2., ..., 2.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4
  * x        (x) int64 0 1 2 3 4
Attributes:
    units:      Ph
    long_name:  Photon

Use multi-wavelength photons

>>> detector.photon.empty()
>>> new_photon_3d = xr.DataArray(
...     np.ones(shape=(4, 5, 5), dtype=float),
...     dims=["wavelength", "y", "x"],
...     coords={"wavelength": [400.0, 420.0, 440.0, 460.0]},
... )
>>> detector.photon.array_3d = detector.photon.array_3d + new_photon_3d
>>> detector.photon += new_photon_3d
>>> detector.photon.array_3d
<xarray.DataArray 'photon' (wavelength: 4, y: 5, x: 5)>
array([[[2., ...., 2.],
        ...,
        [2., ..., 2.]]])
Coordinates:
  * wavelength  (wavelength) float 400.0 ... 460.0
TYPE_LIST = (dtype('float16'), dtype('float32'), dtype('float64'))#
property shape#

Return the shape of the Photon container.

property ndim#
property dtype#
property array#

Two-dimensional numpy array storing monochromatic photon.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

Examples

>>> from pyxel.detectors import CCD, CCDGeometry
>>> detector = CCD(geometry=CCDGeometry(row=2, col=3))
>>> import numpy as np
>>> detector.photon.array = np.array([[1.0, 2.0, 3.0], [4.0, 5.0, 6.0]])
>>> detector.photon.array
array([[1., 2., 3.],
      [4., 5., 6.]])
property array_2d#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property array_3d#

Three-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

to_dict()[source]#
classmethod from_dict(geometry, data)[source]#
to_xarray(dtype=None)[source]#
plot(robust=True)[source]#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png
empty()[source]#

Empty the data container.

class pyxel.data_structure.Charge(geo)[source]#

Bases: object

Charge class representing charge distribution (unit: e⁻).

This class manipulates charge data in the form of a Numpy array and Pandas dataframe.

EXP_TYPE#

alias of float

TYPE_LIST = (dtype('float16'), dtype('float32'), dtype('float64'))#
static create_charges(*, particle_type, particles_per_cluster, init_energy, init_ver_position, init_hor_position, init_z_position, init_ver_velocity, init_hor_velocity, init_z_velocity)[source]#

Create new charge(s) or group of charge(s) as a DataFrame.

Parameters:
  • particle_type (str) – Type of particle. Valid values: ‘e’ for an electron or ‘h’ for a hole.

  • particles_per_cluster (array-like)

  • init_energy (array)

  • init_ver_position (array)

  • init_hor_position (array)

  • init_z_position (array)

  • init_ver_velocity (array)

  • init_hor_velocity (array)

  • init_z_velocity (array)

Returns:

DataFrame – Charge(s) stored in a DataFrame.

convert_df_to_array()[source]#

Convert charge dataframe to an array.

Charge in the detector volume is collected and assigned to the nearest pixel.

static convert_array_to_df(array, num_rows, num_cols, pixel_vertical_size, pixel_horizontal_size)[source]#

Convert charge array to a dataframe placing charge packets in pixels to the center and on top of pixels.

Parameters:
Returns:

DataFrame

add_charge_dataframe(new_charges)[source]#

Add new charge(s) or group of charge(s) to the charge dataframe.

Parameters:

new_charges (DataFrame) – Charges as a DataFrame

add_charge(*, particle_type, particles_per_cluster, init_energy, init_ver_position, init_hor_position, init_z_position, init_ver_velocity, init_hor_velocity, init_z_velocity)[source]#

Add new charge(s) or group of charge(s) inside the detector.

Parameters:
  • particle_type (str) – Type of particle. Valid values: ‘e’ for an electron or ‘h’ for a hole.

  • particles_per_cluster (array)

  • init_energy (array)

  • init_ver_position (array)

  • init_hor_position (array)

  • init_z_position (array)

  • init_ver_velocity (array)

  • init_hor_velocity (array)

  • init_z_velocity (array)

add_charge_array(array)[source]#

Add charge to the charge array. Add to charge dataframe if not empty instead.

Parameters:

array (ndarray)

property array#

Get charge in a numpy array.

to_xarray()[source]#

Convert into a DataArray object.

property frame#

Get charge in a pandas dataframe.

empty()[source]#

Empty all data stored in Charge class.

frame_empty()[source]#

Return True if frame is empty and False otherwise.

validate_type(value)[source]#

Validate a value.

Parameters:

value

validate_shape(value)[source]#

TBW.

get_frame_values(quantity, id_list=None)[source]#

Get quantity values of particles defined with id_list. By default it returns values of all particles.

Parameters:
  • quantity (str) – Name of the quantity: number, energy, position_ver, velocity_hor, etc.

  • id_list (sequence of int) – List of particle ids: [0, 12, 321]

Returns:

array

set_frame_values(quantity, new_value_list, id_list=None)[source]#

Update quantity values of particles defined with id_list. By default it updates all.

Parameters:
  • quantity (str) – Name of the quantity: number, energy, position_ver, velocity_hor, etc.

  • new_value_list (sequence of int) – List of values [1.12, 2.23, 3.65]

  • id_list (sequence of int) – List of particle ids: [0, 12, 321]

remove_from_frame(id_list=None)[source]#

Remove particles defined with id_list. By default it removes all particles from DataFrame.

Parameters:

id_list (sequence of int) – List of particle ids: [0, 12, 321]

class pyxel.data_structure.Pixel(geo)[source]#

Bases: ArrayBase

Pixel class defining and storing information of charge packets within pixel (unit: e⁻).

Accepted array types: np.int32, np.int64, np.uint32, np.uint64, np.float16, np.float32, np.float64.

TYPE_LIST = (dtype('float16'), dtype('float32'), dtype('float64'))#
NAME = 'Pixel'#
UNIT = 'e⁻'#
empty()[source]#

Empty the array by setting the array to zero array in detector shape.

update(data)[source]#

Update ‘array’ attribute.

This method updates ‘array’ attribute of this object with new data. If the data is None, then the object is empty.

Parameters:

data (array_like, Optional)

Examples

>>> from pyxel.data_structure import Photon
>>> obj = Photon(...)
>>> obj.update([[1, 2], [3, 4]])
>>> obj.array
array([[1, 2], [3, 4]])
>>> obj.update(None)  # Equivalent to obj.empty()
property array#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property dtype#

Return array data type.

property ndim#

Return number of dimensions of the array.

property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

plot(robust=True)#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png
property shape#

Return array shape.

to_xarray(dtype=None)#

Convert into a 2D DataArray object with dimensions ‘y’ and ‘x’.

Parameters:

dtype (data-type, optional) – Force a data-type for the array.

Returns:

DataArray – 2D DataArray objects.

Examples

>>> detector.photon.to_xarray(dtype=float)
<xarray.DataArray 'photon' (y: 100, x: 100)>
array([[15149., 15921., 15446., ..., 15446., 15446., 16634.],
       [15149., 15446., 15446., ..., 15921., 16396., 17821.],
       [14555., 14971., 15446., ..., 16099., 16337., 17168.],
       ...,
       [16394., 16334., 16334., ..., 16562., 16325., 16325.],
       [16334., 15978., 16215., ..., 16444., 16444., 16206.],
       [16097., 15978., 16215., ..., 16681., 16206., 16206.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
  * x        (x) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
Attributes:
    units:      Ph
class pyxel.data_structure.Phase(geo)[source]#

Bases: ArrayBase

Phase class.

Accepted array types: np.float16, np.float32 and np.float64.

TYPE_LIST = (dtype('float16'), dtype('float32'), dtype('float64'))#
NAME = 'Phase'#
UNIT = ''#
property array#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property dtype#

Return array data type.

empty()#

Empty the array by setting the array to None.

property ndim#

Return number of dimensions of the array.

property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

plot(robust=True)#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png
property shape#

Return array shape.

to_xarray(dtype=None)#

Convert into a 2D DataArray object with dimensions ‘y’ and ‘x’.

Parameters:

dtype (data-type, optional) – Force a data-type for the array.

Returns:

DataArray – 2D DataArray objects.

Examples

>>> detector.photon.to_xarray(dtype=float)
<xarray.DataArray 'photon' (y: 100, x: 100)>
array([[15149., 15921., 15446., ..., 15446., 15446., 16634.],
       [15149., 15446., 15446., ..., 15921., 16396., 17821.],
       [14555., 14971., 15446., ..., 16099., 16337., 17168.],
       ...,
       [16394., 16334., 16334., ..., 16562., 16325., 16325.],
       [16334., 15978., 16215., ..., 16444., 16444., 16206.],
       [16097., 15978., 16215., ..., 16681., 16206., 16206.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
  * x        (x) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
Attributes:
    units:      Ph
update(data)#

Update ‘array’ attribute.

This method updates ‘array’ attribute of this object with new data. If the data is None, then the object is empty.

Parameters:

data (array_like, Optional)

Examples

>>> from pyxel.data_structure import Photon
>>> obj = Photon(...)
>>> obj.update([[1, 2], [3, 4]])
>>> obj.array
array([[1, 2], [3, 4]])
>>> obj.update(None)  # Equivalent to obj.empty()
class pyxel.data_structure.Signal(geo)[source]#

Bases: ArrayBase

Signal class defining and storing information of detector signal (unit: Volt).

Accepted array types: np.float16, np.float32, np.float64.

TYPE_LIST = (dtype('float16'), dtype('float32'), dtype('float64'))#
NAME = 'Signal'#
UNIT = 'Volt'#
property array#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property dtype#

Return array data type.

empty()#

Empty the array by setting the array to None.

property ndim#

Return number of dimensions of the array.

property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

plot(robust=True)#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png
property shape#

Return array shape.

to_xarray(dtype=None)#

Convert into a 2D DataArray object with dimensions ‘y’ and ‘x’.

Parameters:

dtype (data-type, optional) – Force a data-type for the array.

Returns:

DataArray – 2D DataArray objects.

Examples

>>> detector.photon.to_xarray(dtype=float)
<xarray.DataArray 'photon' (y: 100, x: 100)>
array([[15149., 15921., 15446., ..., 15446., 15446., 16634.],
       [15149., 15446., 15446., ..., 15921., 16396., 17821.],
       [14555., 14971., 15446., ..., 16099., 16337., 17168.],
       ...,
       [16394., 16334., 16334., ..., 16562., 16325., 16325.],
       [16334., 15978., 16215., ..., 16444., 16444., 16206.],
       [16097., 15978., 16215., ..., 16681., 16206., 16206.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
  * x        (x) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
Attributes:
    units:      Ph
update(data)#

Update ‘array’ attribute.

This method updates ‘array’ attribute of this object with new data. If the data is None, then the object is empty.

Parameters:

data (array_like, Optional)

Examples

>>> from pyxel.data_structure import Photon
>>> obj = Photon(...)
>>> obj.update([[1, 2], [3, 4]])
>>> obj.array
array([[1, 2], [3, 4]])
>>> obj.update(None)  # Equivalent to obj.empty()
class pyxel.data_structure.Image(geo)[source]#

Bases: ArrayBase

Image class defining and storing information of detector image (unit: adu).

Accepted array types: np.uint16, np.uint32, np.uint64

TYPE_LIST = (dtype('uint8'), dtype('uint16'), dtype('uint32'), dtype('uint64'))#
NAME = 'Image'#
UNIT = 'adu'#
property array#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

property dtype#

Return array data type.

empty()#

Empty the array by setting the array to None.

property ndim#

Return number of dimensions of the array.

property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

plot(robust=True)#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png
property shape#

Return array shape.

to_xarray(dtype=None)#

Convert into a 2D DataArray object with dimensions ‘y’ and ‘x’.

Parameters:

dtype (data-type, optional) – Force a data-type for the array.

Returns:

DataArray – 2D DataArray objects.

Examples

>>> detector.photon.to_xarray(dtype=float)
<xarray.DataArray 'photon' (y: 100, x: 100)>
array([[15149., 15921., 15446., ..., 15446., 15446., 16634.],
       [15149., 15446., 15446., ..., 15921., 16396., 17821.],
       [14555., 14971., 15446., ..., 16099., 16337., 17168.],
       ...,
       [16394., 16334., 16334., ..., 16562., 16325., 16325.],
       [16334., 15978., 16215., ..., 16444., 16444., 16206.],
       [16097., 15978., 16215., ..., 16681., 16206., 16206.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
  * x        (x) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
Attributes:
    units:      Ph
update(data)#

Update ‘array’ attribute.

This method updates ‘array’ attribute of this object with new data. If the data is None, then the object is empty.

Parameters:

data (array_like, Optional)

Examples

>>> from pyxel.data_structure import Photon
>>> obj = Photon(...)
>>> obj.update([[1, 2], [3, 4]])
>>> obj.array
array([[1, 2], [3, 4]])
>>> obj.update(None)  # Equivalent to obj.empty()
class pyxel.data_structure.ArrayBase(shape)[source]#

Bases: object

Base Array class.

TYPE_LIST = ()#
NAME = ''#
UNIT = ''#
property shape#

Return array shape.

property ndim#

Return number of dimensions of the array.

property dtype#

Return array data type.

empty()[source]#

Empty the array by setting the array to None.

property array#

Two-dimensional numpy array storing the data.

Only accepts an array with the right type and shape.

Raises:

ValueError – Raised if ‘array’ is not initialized.

update(data)[source]#

Update ‘array’ attribute.

This method updates ‘array’ attribute of this object with new data. If the data is None, then the object is empty.

Parameters:

data (array_like, Optional)

Examples

>>> from pyxel.data_structure import Photon
>>> obj = Photon(...)
>>> obj.update([[1, 2], [3, 4]])
>>> obj.array
array([[1, 2], [3, 4]])
>>> obj.update(None)  # Equivalent to obj.empty()
property numbytes#

Recursively calculates object size in bytes using Pympler library.

Returns:

int – Size of the object in bytes.

to_xarray(dtype=None)[source]#

Convert into a 2D DataArray object with dimensions ‘y’ and ‘x’.

Parameters:

dtype (data-type, optional) – Force a data-type for the array.

Returns:

DataArray – 2D DataArray objects.

Examples

>>> detector.photon.to_xarray(dtype=float)
<xarray.DataArray 'photon' (y: 100, x: 100)>
array([[15149., 15921., 15446., ..., 15446., 15446., 16634.],
       [15149., 15446., 15446., ..., 15921., 16396., 17821.],
       [14555., 14971., 15446., ..., 16099., 16337., 17168.],
       ...,
       [16394., 16334., 16334., ..., 16562., 16325., 16325.],
       [16334., 15978., 16215., ..., 16444., 16444., 16206.],
       [16097., 15978., 16215., ..., 16681., 16206., 16206.]])
Coordinates:
  * y        (y) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
  * x        (x) int64 0 1 2 3 4 5 6 7 8 9 10 ... 90 91 92 93 94 95 96 97 98 99
Attributes:
    units:      Ph
plot(robust=True)[source]#

Plot the array using Matplotlib.

Parameters:

robust (bool, optional) – If True, the colormap is computed with 2nd and 98th percentile instead of the extreme values.

Examples

>>> detector.photon.plot()
../../_images/photon_plot.png