toupy.io package

Submodules

toupy.io.dataio module

class toupy.io.dataio.LoadData(**params)[source]

Bases: PathName, Variables

Load projections from an HDF5 file.

Parameters:
  • **params – Keyword arguments forwarded to PathName via a params HDF5 file lookup. Additional recognised keys:

  • params["amponly"] (bool, optional) – Return only the amplitude of complex projections. Default False.

  • params["phaseonly"] (bool, optional) – Return only the phase of complex projections. Default False.

  • params["loadroi"] (bool, optional) – Load only a sub-region of interest. Default False.

  • params["roi"] (list of int, optional) – [proj_start, proj_end, row_start, row_end, col_start, col_end] when loadroi=True.

classmethod load(*args, **params)[source]

Load data from h5 file

Parameters:
  • h5name (str) – File name from which data is loaded

  • **params – Dictionary of additonal parameters

  • params["autosave"] (bool) – Save the projections once load without asking

  • params["phaseonly"] (bool) – Load only phase projections. Used when the projections are complex-valued.

  • params["amponly"] (bool) – Load only amplitude projections. Used when the projections are complex-valued.

  • params["pixtol"] (float) – Tolerance for alignment, which is also used as a search step

  • params["alignx"] (bool) – True or False to activate align x using center of mass (default= False, which means align y only)

  • params["shiftmeth"] (str) – Shift images with fourier method (default). The options are linear -> Shift images with linear interpolation (default); fourier -> Fourier shift or spline -> Shift images with spline interpolation.

  • params["circle"] (bool) – Use a circular mask to eliminate corners of the tomogram

  • params["filtertype"] (str) – Filter to use for FBP

  • params["freqcutoff"] (float) – Frequency cutoff for tomography filter (between 0 and 1)

  • params["cliplow"] (float) – Minimum value in tomogram

  • params["cliphigh"] (float) – Maximum value in tomogram

  • params["correct_bad"] (bool) – If true, it will interpolate bad projections. The numbers of projections to be corrected is given by params[“bad_projs”].

  • params["bad_projs"] (list of ints) – List of projections to be interpolated. It starts at 0.

Returns:

  • stack_projs (array_like) – Stack of projections

  • theta (array_like) – Stack of thetas

  • shiftstack (array_like) – Shifts in vertical (1st dimension) and horizontal (2nd dimension)

  • datakwargs (dict) – Dictionary with metadata information

classmethod load_olddata(*args, **params)[source]

Load old data from h5 file. It should disappear soon.

Parameters:
  • h5name (str) – File name from which data is loaded

  • **params – Dictionary of additonal parameters

  • params["autosave"] (bool) – Save the projections once load without asking

  • params["phaseonly"] (bool) – Load only phase projections. Used when the projections are complex-valued.

  • params["amponly"] (bool) – Load only amplitude projections. Used when the projections are complex-valued.

  • params["pixtol"] (float) – Tolerance for alignment, which is also used as a search step

  • params["alignx"] (bool) – True or False to activate align x using center of mass (default= False, which means align y only)

  • params["shiftmeth"] (str) – Shift images with fourier method (default). The options are linear -> Shift images with linear interpolation (default); fourier -> Fourier shift or spline -> Shift images with spline interpolation.

  • params["circle"] (bool) – Use a circular mask to eliminate corners of the tomogram

  • params["filtertype"] (str) – Filter to use for FBP

  • params["freqcutoff"] (float) – Frequency cutoff for tomography filter (between 0 and 1)

  • params["cliplow"] (float) – Minimum value in tomogram

  • params["cliphigh"] (float) – Maximum value in tomogram

  • params["correct_bad"] (bool) – If true, it will interpolate bad projections. The numbers of projections to be corrected is given by params[“bad_projs”].

  • params["bad_projs"] (list of ints) – List of projections to be interpolated. It starts at 0.

Returns:

  • stack_projs (array_like) – Stack of projections

  • theta (array_like) – Stack of thetas

  • shiftstack (array_like) – Shifts in vertical (1st dimension) and horizontal (2nd dimension)

  • datakwargs (dict) – Dictionary with metadata information

classmethod loadmasks(*args, **params)[source]

Load masks from previous h5 file

Parameters:

h5name (str) – File name from which data is loaded

Returns:

masks – Array with the masks

Return type:

array_like

classmethod loadshiftstack(*args, **params)[source]

Load the shift-stack from a previous HDF5 file.

Parameters:

h5name (str) – File name from which data is loaded.

Returns:

shiftstack – Shifts in the vertical (row 0) and horizontal (row 1) directions.

Return type:

ndarray, shape (2, nprojs)

classmethod loadtheta(*args, **params)[source]

Load theta angles from a previous HDF5 file.

Parameters:

h5name (str) – File name from which data is loaded.

Returns:

theta – Array of tomographic angles in degrees.

Return type:

ndarray

class toupy.io.dataio.LoadProjections(**params)[source]

Bases: PathName, Variables

Load the reconstructed projections from ptyr, cxi, or edf files.

Parameters:
  • **params – Keyword arguments forwarded to PathName. See that class for the full list. Additional recognised keys:

  • params["showrecons"] (bool, optional) – Show each projection interactively as it is loaded. Default False.

  • params["border_crop_x"] (int or None, optional) – Pixels to crop from each horizontal border. Default None.

  • params["border_crop_y"] (int or None, optional) – Pixels to crop from each vertical border. Default None.

  • params["checkextraprojs"] (bool, optional) – Remove projections acquired at or beyond 180 degrees. Default True.

  • params["missingprojs"] (bool, optional) – Interpolate missing projections. Default False.

  • params["missingnum"] (list of int, optional) – Indices of missing projections to interpolate. Default None.

  • params["cxientry"] (str or None, optional) – HDF5 path inside CXI file. Default None.

  • params["detector"] (str, optional) – Detector name used to find angles in the raw HDF5 header. Default "Frelon".

check_angles()[source]

Find projection angles from ICA metadata and plot them for verification.

Specific to the ID16A beamline at ESRF. Reads the theta motor position from the ICA HDF5 metadata file, removes duplicate angles, and prompts the user to confirm before returning.

Returns:

  • angles (list of float) – Sorted list of unique tomographic angles in degrees.

  • thetas (dict) – Mapping of scan-key strings to raw theta values.

check_angles_new()[source]

Find projection angles from raw HDF5 scan files and plot them for verification.

Reads theta directly from the raw acquisition HDF5 files rather than from the ICA metadata file. Specific to the ID16A beamline at ESRF.

Returns:

  • angles (list of float) – Sorted list of unique tomographic angles in degrees.

  • thetas (dict) – Mapping of scan-key strings to raw theta values.

static insert_missing(stack_objs, theta, missingnum)[source]

Insert missing projections by interpolation of their neighbours.

Parameters:
  • stack_objs (ndarray, shape (nprojs, nr, nc)) – Stack of projections.

  • theta (ndarray, shape (nprojs,)) – Corresponding angle array.

  • missingnum (list of int) – Zero-based indices of the missing projections to be inserted.

Returns:

  • stack_objs (ndarray) – Updated stack with interpolated projections inserted.

  • theta (ndarray) – Updated angle array with interpolated angles inserted.

classmethod load(**params)[source]

Load the reconstructed projections from phase-retrieved files.

Parameters:
  • **params – Container with parameters to load the files.

  • params["account"] (str) – User experiment number at ESRF.

  • params["samplename"] (str) – Sample name

  • params["pathfilename"] (str) – Path to the first projection file.

  • params["regime"] (str) – Imaging regime. The options are: nearfield, farfield, holoct.

  • params["showrecons"] (bool) – To show or not the projections once loaded

  • params["autosave"] (bool) – Save the projections once load without asking

  • params["phaseonly"] (bool) – Load only phase projections. Used when the projections are complex-valued.

  • params["amponly"] (bool) – Load only amplitude projections. Used when the projections are complex-valued.

  • params["border_crop_x"] (int, None) – Amount of pixels to crop at each border in x.

  • params["border_crop_y"] (int, None) – Amount of pixels to crop at each border in y.

  • params["checkextraprojs"] (bool) – Check for the projections acquired at and over 180 degrees.

  • params["missingprojs"] (bool) – Allow to interpolate for missing projections. The numbers of the projections need to be provided in params[“missingnum”].

  • params["missingnum"] (list of ints) – Numbers of the missing projections to be interpolated.

Returns:

  • stack_objs (array_like) – Array containing the projections

  • stack_angles (array_like) – Array containing the thetas

  • pxsize (list of floats) – List containing the pixel size in the vertical and horizontal directions. Typically, the resolution is isotropic and the two values are the same

  • paramsload (dict) – Parameters of the loading

classmethod loadedf(**params)[source]

Load the reconstructed projections from the edf files This is adapted for the phase-contrast imaging generating projections as edf files

Parameters:
  • **params – Container with parameters to load the files.

  • params["account"] (str) – User experiment number at ESRF.

  • params["samplename"] (str) – Sample name

  • params["pathfilename"] (str) – Path to the first projection file.

  • params["regime"] (str) – Imaging regime. The options are: nearfield, farfield, holoct.

  • params["showrecons"] (bool) – To show or not the projections once loaded

  • params["autosave"] (bool) – Save the projections once load without asking

Returns:

  • stack_objs (array_like) – Array containing the projections

  • stack_angles (array_like) – Array containing the thetas

  • pxsize (list of floats) – List containing the pixel size in the vertical and horizontal directions. Typically, the resolution is isotropic and the two values are the same

  • paramsload (dict) – Parameters of the loading

class toupy.io.dataio.LoadTomogram(**params)[source]

Bases: LoadData

Load a tomogram from an HDF5 file.

Parameters:

**params – Keyword arguments forwarded to LoadData.

classmethod load(*args, **params)[source]

Load tomographic data from h5 file

Parameters:

args[0] (str) – HDF5 file name from which data is loaded

Returns:

  • tomogram (array_like) – Stack of tomographic slices

  • theta (array_like) – Stack of thetas

  • shiftstack (array_like) – Shifts in vertical (1st dimension) and horizontal (2nd dimension)

  • datakwargs (dict) – Dictionary with metadata information

class toupy.io.dataio.PathName(**params)[source]

Bases: object

Class to manage file location and paths.

Parameters:
  • **params – Dictionary of parameters.

  • params["pathfilename"] (str) – Path to first file and filename.

  • params["account"] (str) – User account number.

  • params["samplename"] (str) – Sample name.

  • params["scanprefix"] (str) – Prefix of the scan folder names.

  • params["regime"] (str) – Imaging regime. Options are nearfield, farfield, holoct.

  • params["pycudaprojs"] (bool, optional) – Whether reconstructions were done with PyCUDA. Default True.

  • params["legacy"] (bool, optional) – Whether the raw-data filenames follow the legacy naming scheme. Default True.

  • params["thetameta"] (bool, optional) – Whether to read angles from the ESRF ICA metadata HDF5 file. Default True.

datafilewcard()[source]

Create file wildcard to search for files.

Returns:

file_wcard – Wildcard string for globbing data files.

Return type:

str

metadatafilewcard()[source]

Create file wildcard to search for metafiles.

Returns:

metafile_wcard – Wildcard string for globbing metadata files.

Return type:

str

Raises:

IOError – If the ICA HDF5 metadata file does not exist, or if the projection file extension is not .ptyr, .cxi, or .edf.

results_datapath(h5name)[source]

Create the full path for an HDF5 file inside the result folder.

Parameters:

h5name (str) – HDF5 file name (basename only).

Returns:

h5file – Absolute path to the HDF5 file.

Return type:

str

results_folder()[source]

Create path for the result folder.

Returns:

results_path – Absolute path to the results folder. The folder is created if it does not already exist.

Return type:

str

Raises:

ValueError – If self.regime is not one of nearfield, farfield, or holoct.

search_projections()[source]

Search for projection files given the filenames.

Returns:

scan_wcard – Wildcard string for globbing projection files.

Return type:

str

Raises:

IOError – If the projection file extension is not .ptyr, .cxi, or .edf.

class toupy.io.dataio.SaveData(**params)[source]

Bases: PathName, Variables

Save projections to an HDF5 file.

Parameters:
  • **params – Keyword arguments forwarded to PathName.

  • params["autosave"] (bool, optional) – If True save without prompting the user. Default False.

  • params["cxientry"] (str or None, optional) – CXI entry path. Default None.

classmethod save(*args, **params)[source]

Save data to HDF5 File

Parameters:
  • *args – positional arguments

  • args[0] (str) – H5 file name

  • args[1] (array_like) – Array containing the stack of projections

  • args[2] (array_like) – Values of theta

  • args[3] (array_like) – Array containing the shifts for each projection in the stack. If not provided, it will be initialized with zeros

  • args[4] (array_like or None) – Array containing the projection masks

classmethod saveFSC(*args, **params)[source]

Save FSC data to HDF5 file

Parameters:
  • *args – positional arguments

  • args[0] (str) – H5 file name

  • args[1] (array_like) – Normalized frequencies

  • args[2] (array_like) – Value of the threshold for each frequency

  • args[3] (array_like) – The FSC curve

  • args[4] (array_like) – The first tomogram

  • args[5] (array_like) – The second tomogram

  • args[6] (array_like) – The array of theta values

  • args[7] (float) – Pixel size

savecheck()[source]

Decorator for save data

classmethod savemasks(*args, **params)[source]
class toupy.io.dataio.SaveTomogram(**params)[source]

Bases: SaveData

Save a tomogram volume to an HDF5 file.

Parameters:

**params – Keyword arguments forwarded to SaveData.

classmethod convert_to_tiff(*args, **params)[source]

Convert the HDF5 file with the tomogram to tiff

Parameters:
  • *args – positional arguments

  • args[0] (str) – H5 file name

  • args[1] (array_like) – Array containing the stack of slices (tomogram)

  • args[2] (array_like) – Values of theta

  • args[3] (array_like) – Array containing the shifts for each projection in the stack

classmethod save(*args, **params)[source]
Parameters:
  • *args – positional arguments

  • args[0] (str) – H5 file name

  • args[1] (array_like) – Array containing the stack of slices (tomogram)

  • args[2] (array_like) – Values of theta

  • args[3] (array_like) – Array containing the shifts for each projection in the stack

classmethod save_vol_to_h5(*args, **params)[source]
savecheck()[source]

Decorator for save data

tiff_folderpath(foldername)[source]

Create the path to the folder in which the TIFF files will be stored.

Parameters:

foldername (str) – Name of the sub-folder to be created inside the results folder.

Returns:

folderpath – Absolute path to the TIFF sub-folder.

Return type:

str

toupy.io.dataio.remove_extraprojs(stack_projs, theta, n_extra=0)[source]

Remove extra projections appended at the end of a tomographic scan (e.g. the 180°, 90° and 0° bookend angles).

Parameters:
  • stack_projs (array_like) – Stack of projections, first index = projection number.

  • theta (array_like) – Array of theta values.

  • n_extra (int, optional) – Number of trailing projections to remove. When 0 (default) the function just prints the last 5 angles and returns unchanged.

Returns:

  • stack_projs (array_like)

  • theta (array_like)

toupy.io.filesrw module

Files read and write

toupy.io.filesrw.convert16bitstiff(tiffimage, low_cutoff, high_cutoff)[source]

Convert a 16-bit TIFF image back to quantitative values.

Parameters:
  • tiffimage (array_like) – Image read from a 16-bit TIFF file.

  • low_cutoff (float) – Low cutoff of the gray level used during the original normalisation.

  • high_cutoff (float) – High cutoff of the gray level used during the original normalisation.

Returns:

tiffimage – Array containing the image restored to quantitative values.

Return type:

ndarray, float

toupy.io.filesrw.convert8bitstiff(tiffimage, low_cutoff, high_cutoff)[source]

Convert an 8-bit TIFF image back to quantitative values.

Parameters:
  • tiffimage (array_like) – Image read from an 8-bit TIFF file.

  • low_cutoff (float) – Low cutoff of the gray level used during the original normalisation.

  • high_cutoff (float) – High cutoff of the gray level used during the original normalisation.

Returns:

tiffimage – Array containing the image restored to quantitative values.

Return type:

ndarray, float

toupy.io.filesrw.convertimageto16bits(input_image, low_cutoff, high_cutoff)[source]

Convert image gray-level to 16 bits with normalization

Parameters:
  • input_image (array_like) – Input image to be converted.

  • low_cutoff (float) – Low cutoff of the gray level.

  • high_cutoff (float) – High cutoff of the gray level.

Returns:

tiffimage – Array containing the image at 16 bits.

Return type:

array_like

toupy.io.filesrw.convertimageto8bits(input_image, low_cutoff, high_cutoff)[source]

Convert image gray-level to 8 bits with normalization.

Parameters:
  • input_image (array_like) – Input image to be converted.

  • low_cutoff (float) – Low cutoff of the gray level.

  • high_cutoff (float) – High cutoff of the gray level.

Returns:

tiffimage – Array containing the image at 8 bits.

Return type:

array_like

toupy.io.filesrw.create_paramsh5(**params)[source]

Create parameter file in HDF5 format

Parameters:

params (dict) – Dictionary containing the parameters to be saved

toupy.io.filesrw.crop_array(input_array, delcropx, delcropy)[source]

Crop borders from 2D arrays

Parameters:
  • input_array (array_like) – Input array to be cropped

  • delcropx (int) – Number of pixels to be crop from borders in x and y directions

  • delcropy (int) – Number of pixels to be crop from borders in x and y directions

Returns:

cropped_array – Cropped array

Return type:

array_like

toupy.io.filesrw.load_paramsh5(**params)[source]

Load parameters from the sample’s HDF5 parameter file.

Parameters:

**params – Must contain params["samplename"] (str) to derive the filename <samplename>_params.h5.

Returns:

out_params – Parameters read from the HDF5 file, updated with any keys in params that were not present in the file.

Return type:

dict

toupy.io.filesrw.memmap_volfile(filename)[source]

Memory map the tomogram from .vol file

Parameters:

filename (str) – filename to be read

Returns:

  • tomogram (numpy.memmap, float32, shape (z_size, x_size, y_size)) – Memory-mapped 3D array containing the tomogram.

  • voxelsize (float) – Voxel size in metres.

  • arrayshape (tuple of int) – The array shape (x_size, y_size, z_size).

Examples

>>> volpath = 'volfilename.vol'
>>> tomogram, voxelsize, arrayshape = memmap_volfile(volpath)

Notes

The volume info file (.vol.info) containing the metadata of the volume must reside in the same folder as the .vol file.

toupy.io.filesrw.read_cxi(pathfilename, correct_orientation=True)[source]

Read reconstruction files .cxi from PyNX

Parameters:
  • pathfilename (str) – Path to file

  • correct_orientation (bool) – True for correcting the image orientation and False to keep as it is. The default value is True.

Returns:

  • data1 (array_like, complex) – Object image

  • probe1 (array_like, complex) – Probe images

  • pixelsize (list of floats) – List with pixelsizes in vertical and horizontal directions

  • energy (float) – Energy of the incident photons

Examples

>>> imgpath = 'filename.cxi'
>>> objdata, probedata, pixel, energy = read_cxi(imgpath)
toupy.io.filesrw.read_edf(fname)[source]

Read EDF files of tomographic datasets

Parameters:

fname (str) – Path to file

Returns:

  • projs (array_like) – Array of projections

  • pixelsize (list of floats) – List with pixelsizes in vertical and horizontal directions

  • energy (float) – Energy of the incident photons

  • nvue (int) – Number of projections

toupy.io.filesrw.read_ptyr(pathfilename, correct_orientation=True)[source]

Read reconstruction files .ptyr from Ptypy

Parameters:
  • pathfilename (str) – Path to file

  • correct_orientation (bool, optional) – True for correcting the image orientation and False to keep as it is. The default value is True.

Returns:

  • data1 (array_like, complex) – Object image

  • probe1 (array_like, complex) – Probe images

  • pixelsize (list of floats) – List with pixelsizes in vertical and horizontal directions

  • energy (float) – Energy of the incident photons

Examples

>>> imgpath = 'filename.ptyr'
>>> objdata, probedata, pixel, energy = read_ptyr(imgpath)
toupy.io.filesrw.read_recon(filename, correct_orientation=False)[source]

Wrapper for choosing the function to read recon file

Parameters:
  • pathfilename (str) – Path to file

  • correct_orientation (bool, optional) – True for correcting the image orientation and False to keep as it is. The default value is False.

Returns:

  • data1 (array_like, complex) – Object image

  • probe1 (array_like, complex) – Probe images

  • pixelsize (list of floats) – List with pixelsizes in vertical and horizontal directions

  • energy (float) – Energy of the incident photons

Examples

>>> imgpath = 'filename.ptyr'
>>> objdata, probedata, pixel, energy = read_recon(imgpath)
toupy.io.filesrw.read_theta_raw(pathfilename, detector='Frelon')[source]

Read the tomographic angle from a raw acquisition HDF5 file at ID16A.

Parameters:
  • pathfilename (str) – Path to the raw HDF5 file.

  • detector (str, optional) – Detector name used to construct the HDF5 path to the motor header. Default is "Frelon".

Returns:

theta – Tomographic angle in degrees (value of the somega motor).

Return type:

float

Examples

>>> imgpath = 'filename.h5'
>>> theta = read_theta_raw(imgpath)
toupy.io.filesrw.read_theta_recon(reconfile)[source]

Auxiliary function to read theta from recon files

Parameters:

reconfile (str) – Path to recon file

Returns:

theta – Tomographic angle

Return type:

float

Examples

>>> imgpath = 'filename.ptyr'
>>> theta = read_theta_recon(imgpath)
toupy.io.filesrw.read_tiff(imgpath)[source]

Read tiff files using skimage.io.imread

Parameters:

imgpath (str) – Path to tiff file with extension

Returns:

imgout – Array containing the image

Return type:

array_like

Examples

>>> imgpath = 'image.tiff'
>>> ar = read_tiff(imgpath)
>>> ar.dtype
dtype('uint16')
>>> np.max(ar)
65535
toupy.io.filesrw.read_tiff_info(tiff_info_file)[source]

Read info file from tiff slices of the reconstructed tomographic volume

Parameters:

tiff_info_file (str) – Info filename

Returns:

  • low_cutoff (float) – Low cutoff of the gray level

  • high_cutoff (float) – High cutoff of the gray level

  • pixelsize (float) – Pixelsize in nanometers

Notes

The info file here is the file that is saved when the volume is exported to TIFF files. It is not the info file saved by the volume reconstruction when saving in .vol format.

toupy.io.filesrw.read_volfile(filename)[source]

Read tomogram from .vol file

Parameters:

filename (str) – filename to be read

Returns:

  • tomogram (ndarray, float32, shape (z_size, x_size, y_size)) – 3D array containing the tomogram.

  • voxelsize (float) – Voxel size in metres.

  • arrayshape (tuple of int) – The array shape (x_size, y_size, z_size).

Examples

>>> volpath = 'volfilename.vol'
>>> tomogram, voxelsize, arrayshape = read_volfile(volpath)

Notes

The volume info file (.vol.info) containing the metadata of the volume must reside in the same folder as the .vol file.

toupy.io.filesrw.write_edf(fname, data_array, hd=None)[source]

Write EDF files

Parameters:
  • fname (str) – File name

  • data_array (array_like) – Data to be saved as edf

  • hd (dict) – Dictionary with header information

toupy.io.filesrw.write_paramsh5(h5filename, **params)[source]

Writes params to HDF5 file

Parameters:
  • h5filename (str) – Filename of the params file

  • params (dict) – Dictionary containing the parameters to be saved

toupy.io.filesrw.write_tiff(input_array, pathfilename, plugin='tifffile')[source]

Write tiff files using skimag.io.imsave

Parameters:
  • input_array (array_like) – Input array to be saved

  • pathfilename (str) – Path and filename to save the file

toupy.io.filesrw.write_tiffmetadata(filename, low_cutoff, high_cutoff, factor, **params)[source]

Creates a txt file with the information about the Tiff normalization

Parameters:
  • filename (str) – Filename to save the file.

  • low_cutoff (float) – Low cutoff value for the tiff normalization.

  • high_cutoff (float) – High cutoff value for the tiff normalization.

  • factor (float) – Multiplicative factor in case it is needed.

  • params (dict) – Dictionary of additional parameters.

  • params["voxelsize"] (float) – Voxel size.

  • params["filtertype"] (str) – Filter used in the tomographic reconstruction.

  • params["freqcutoff"] (float) – Frequency cutoff used in the tomographic reconstruction.

  • params["bits"] (int) – The tiff type. Options: 8 for 8 bits or 16 for 16 bits.

toupy.io.h5chunk_shape_3D module

toupy.io.h5chunk_shape_3D.binlist(n, width=0)[source]

Return list of bits that represent a non-negative integer.

Parameters:
  • n (int) – Non-negative integer to convert.

  • width (int, optional) – Minimum number of bits in the returned zero-filled list. Default is 0.

Returns:

List of 0/1 values representing n in binary, zero-padded on the left to at least width elements.

Return type:

list of int

toupy.io.h5chunk_shape_3D.chunk_shape_3D(varShape, valSize=4, chunkSize=4096)[source]

Return a ‘good shape’ for a 3D variable, assuming balanced 1D/(n-1)D access [1]

Parameters:
  • varShape (sequence of ints) – length 3 list of variable dimension sizes

  • chunkSize (int, optional) – maximum chunksize desired, in bytes (default 4096)

  • valSize (int, optional) – size of each data value, in bytes (default 4)

Returns:

Returns integer chunk lengths of a chunk shape that provides balanced access of 1D subsets and 2D subsets of a netCDF or HDF5 variable var with shape (T, X, Y), where the 1D subsets are of the form var[:,x,y] and the 2D slices are of the form var[t,:,:], typically 1D time series and 2D spatial slices.

Return type:

tuple

Notes

‘Good shape’ for chunks means that the number of chunks accessed to read either kind of 1D or 2D subset is approximately equal, and the size of each chunk (uncompressed) is no more than chunkSize, which is often a disk block size. Code fetched from [2] and [3].

References

toupy.io.h5chunk_shape_3D.numVals(shape)[source]

Return number of values in chunk of specified shape, given by a list of dimension lengths.

Parameters:

shape (sequence of ints) – List of variable dimension sizes.

Returns:

Product of all dimension sizes, or 1 if shape is empty.

Return type:

int

toupy.io.h5chunk_shape_3D.perturbShape(shape, onbits)[source]

Return shape perturbed by adding 1 to elements corresponding to 1 bits in onbits.

Parameters:
  • shape (sequence of ints) – List of variable dimension sizes.

  • onbits (int) – Non-negative integer less than 2**len(shape).

Returns:

New shape with 1 added to each dimension whose corresponding bit in onbits is set.

Return type:

list of int