cerr.radiomics package

Submodules

cerr.radiomics.first_order module

This module contains routine for calculation of 1st order statistical features.

cerr.radiomics.first_order.stats(planC, structNum, offsetForEnergy=0, binWidth=None, binNum=None)[source]

This routine calculates 1st order statistical features from the input image and segmentation based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#intensity-based-statistical-features.

Parameters:
  • planC (plan_container.planC or np.ndarray(dtype=flat)) – planC containing the image and segmentation or scan np.ndarray.

  • structNum (int or np.ndarray(dtype=int)) – index of structure number in planC or binary mask for segmentation that matches the scan dimensions

  • offsetForEnergy (float) – optional, value to add to scan for computing the Energy, TotalEnergy and RMS features

  • binWidth (float) – optional, bin width for discretizing the input scan. Required when binNum is None. The default value is None.

  • binNum (int) – optional, number of bins to discretize the input scan. Required when binWidth is None. The default value is None.

Returns:

dictionary of features

Return type:

dict

cerr.radiomics.gray_level_cooccurence module

This module contains routines for calculation of Gray Level Co-occurrence matrix and features.

cerr.radiomics.gray_level_cooccurence.calcCooccur(quantizedM, offsetsM, nL, cooccurType=1)[source]

This function calculates the cooccurrence matrix for the passed quantized image based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-co-occurrence-based-features

Parameters:
  • quantizedM (np.ndarray(dtype=int)) – quantized 3d matrix obtained, for example, from radiomics.preprocess.imquantize_cerr

  • offsetsM (np.ndarray(dtype=int)) – Offsets for directionality/neighbors, obtained from radiomics.ibsi1.getDirectionOffsets

  • nL (int) – Number of gray levels. nL must be less than 65535.

  • cooccurType (int) – flag, 1 or 2. 1: returns a single cooccurrence matrix by combining contributions from all offsets into one cooccurrence matrix. 2: returns cooccurM with each column containing cooccurrence matrix for the row of offsetsM.

Returns:

cooccurrence matrix of size (nL*nL) x 1 for cooccurType = 1,

or (nL*nL) x offsetsM.shape[0] for cooccurType = 2.

The output can be passed to cooccurToScalarFeatures to get GLCM texture features.

Return type:

np.ndarray

cerr.radiomics.gray_level_cooccurence.cooccurToScalarFeatures(cooccurM)[source]

This function calculates scalar texture features from cooccurrence matrix as per IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-co-occurrence-based-features

Parameters:

cooccurM (np.ndarray(dtype=int)) – cooccurrence matrix of size (nL*nL) x 1 for combined cooccurrrences from all the directions, or (nL*nL) x offsetsM.shape[0] for individual directions or patches.

Returns:

dictionary with scalar texture features as its

fields. Each field’s value is a vector containing the feature values for each column cooccurM.

Return type:

dict

cerr.radiomics.ibsi1 module

This module contains routines for IBSI1 compatible radiomics calculation

cerr.radiomics.ibsi1.getDirectionOffsets(direction)[source]
cerr.radiomics.ibsi1.calcRadiomicsForImgType(volToEval, maskBoundingBox3M, morphMask3M, gridS, paramS)[source]
cerr.radiomics.ibsi1.computeScalarFeatures(scanNum, structNum, settingsFile, planC)[source]
cerr.radiomics.ibsi1.getIBSINameMap()[source]
cerr.radiomics.ibsi1.createFieldNameFromParameters(imageType, settingS)[source]

Create unique fieldname for radiomics features returned by calcRadiomicsForImgType :param imageType: ‘Original’ or filtered image type (see processImage.m for valid options) :param settingS: Parameter dictionary for radiomics feature extraction :return: fieldName

cerr.radiomics.ibsi1.createFlatFeatureDict(featDict, imageType, avgType, directionality, mapToIBSI=False)[source]
cerr.radiomics.ibsi1.writeFeaturesToFile(featList, csvFileName, writeHeader=True)[source]

cerr.radiomics.neighbor_gray_level_dependence module

This module contains routines for claculation of Gray Level Dependence texture features

cerr.radiomics.neighbor_gray_level_dependence.calcNGLDM(scan_array, patch_size, num_grayscale_levels, a)[source]

This function calculates the Neighborhood Gray Level Dependence Matrix for the passed quantized image based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#neighbouring-grey-level-dependence-based-features

Parameters:
  • scan_array (np.ndarray(dtype=int)) – quantized 3d matrix obtained, for example, from radiomics.preprocess.imquantize_cerr

  • patch_size (list) – list of length 3 defining patch radius for row, col, slc dimensions.

  • num_grayscale_levels (int) – Number of gray levels.

  • a – coarseness parameter

Returns:

NGLDM matrix os size (num_grayscale_levels, max_nbhood_size)

The output can be passed to ngldmToScalarFeatures to get NGLDM texture features.

Return type:

np.ndarray

cerr.radiomics.neighbor_gray_level_dependence.ngldmToScalarFeatures(s, numVoxels)[source]

This function calculates scalar texture features from Neighborhood Gray Level Dependence matrix as per IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#neighbouring-grey-level-dependence-based-features

Parameters:
  • (np.ndarray (s) – NGLDM matrix of size (num_grayscale_levels, max_nbhood_size)

  • numVoxels (int) – number of voxels in the region of interest used for generating s.

Returns:

dictionary with scalar NGLDM texture features as its

fields.

Return type:

dict

cerr.radiomics.neighbor_gray_tone module

This module contains routines for claculation of Gray Tone Difference texture features

cerr.radiomics.neighbor_gray_tone.calcNGTDM(scan_array, patch_size, numGrLevels)[source]

This function calculates the Neighborhood Gray Tone Difference Matrix for the passed quantized image based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#neighbourhood-grey-tone-difference-based-features

Parameters:
  • scan_array (np.ndarray(dtype=int)) – quantized 3d matrix obtained, for example, from radiomics.preprocess.imquantize_cerr

  • patch_size (list) – list of length 3 defining patch radius for row, col, slc dimensions.

  • num_grayscale_levels (int) – Number of gray levels.

  • a – coarseness parameter

Returns:

NGLDM matrix os size (num_grayscale_levels, max_nbhood_size)

The output can be passed to ngldmToScalarFeatures to get NGLDM texture features.

Return type:

np.ndarray

cerr.radiomics.neighbor_gray_tone.ngtdmToScalarFeatures(s, p, Nvc)[source]

cerr.radiomics.preprocess module

This module contains routines for processing image and segmentation mask for radiomics calculation.

cerr.radiomics.preprocess.imquantize(x, num_level=None, xmin=None, xmax=None, binwidth=None)[source]

Function to quantize an image based on the number of bins (num_level) or the bin width (bin_width). The min and max are computed from the input image image when they are not provided.

Parameters:
  • x (np.ndarray) – Input image matrix.

  • num_level (int) – [optional, default: Use fixed bin_width] Number of quantization levels.

  • xmin (int) – [optional, default: Use min intensity in x] Minimum value for quantization.

  • xmax (int) – [optional, default: Use max intensity in x] Maximum value for quantization.

  • bin_width (int) – [optional, default: Use fixed num_level] Bin width for quantization.

Returns:

Quantized image.

Return type:

np.ndarray(dtype=int)

cerr.radiomics.preprocess.getResampledGrid(resampResolutionV, xValsV, yValsV, zValsV, gridAlignMethod='center', *args)[source]

Function to create x,y,z image grid by resampling the input x,y,z grid at the required resolution.

Parameters:
  • resampResolutionV (np.ndarray) – Input image matrix.

  • xValsV (np.array) – x-coordinates of the image grid.

  • yValsV (np.array) – y-coordinates od the image grid.

  • zValsV (np.array) – z-coordinates od the image grid.

  • gridAlignMethod – [optional, default=’center’] Currently, the only method supported is “center”.

Returns:

(xResampleV, yResampleV, zResampleV) Resamples x,y,z grid coordinates.

Return type:

Tuple

cerr.radiomics.preprocess.imgResample3D(img3M, xValsV, yValsV, zValsV, xResampleV, yResampleV, zResampleV, method, extrapVal=None, inPlane=False)[source]
Parameters:
  • img3M (np.ndarray) – 3D scan array e.g. planC.scan[scanNum].getScanArray()

  • xValsV (np.array) – x-coordinates i.e. along columns of img3M. Must be monotonically increasing order as per CERR coordinate system.

  • yValsV (np.array) – y-coordinates i.e. along rows of img3M. Must be monotonically decreasing order as per CERR coordinate system.

  • zValsV (np.array) – z-coordinates i.e. along slices of img3M. Must be monotonically increasing order as per CERR coordinate system.

  • xResampleV (np.array) – new x-coordinates i.e. along columns of resampled image

  • yResampleV (np.array) – new y-coordinates i.e. along rows of resampled image

  • zResampleV (np.array) – new z-coordinates i.e. along slices of resampled image

  • method (string) – Resampling methods supported by SimpleITK, viz. ‘sitkNearestNeighbor’, ‘sitkLinear’, ‘sitkBSpline’, ‘sitkGaussian’, ‘sitkLabelGaussian’,’sitkHammingWindowedSinc’,’sitkCosineWindowedSinc’, ‘sitkWelchWindowedSinc’,’sitkLanczosWindowedSinc’, ‘sitkBlackmanWindowedSinc’

  • extrapVal (float) – Value of extrapolated pixels. When not specified, the value of an extrapolated pixel is assigned from the nearest neighbor.

  • inPlane (bool) – [optional, default=False] Specify whether to restrict the interpolation to in-plane. (e.g. bi-linear). The default is 3D interpolation.

Returns:

(np.ndarray) 3D array resampled at xResampleV, yResampleV, zResampleV

cerr.radiomics.preprocess.padScan(scan3M, mask3M, method, marginV, cropFlag=True)[source]

Function to pad the input image using specified method and padding size.

Parameters:
  • scan3M (np.ndarray) – 3D scan.

  • mask3M (np.ndarray) – 3D mask.

  • method (string) – Padding method. The following methods are supported: ‘expand’ - Image is cropped around the mask and expanded using the specified margin. ‘padzeros’ - Image is padded by zeros. ‘periodic’ - Image is padded with periodic expansion. ‘nearest’ - Image is padded by using values from nearest neighbors. ‘mirror’ - Image is padded by mirrorig the boundary region as per the margin.

  • marginV (np.array, 1D) – [nRows, nCols, nSlices] in voxels.

  • cropFlag – [optional, default:True] bool for flag to crop around mask bounding box when set to True.

Returns:

) outScan3M
outMask3M

outLimitsV

np.ndarray(dtype=int): , maskBoundingBox3M, morphmask3M, gridS, paramS, diagS

Return type:

np.ndarray(dtype=

cerr.radiomics.preprocess.unpadScan(padScan3M, marginV)[source]

This function return the image array after stripping off specified padding margin

Parameters:
  • padScan3M (np.ndarray) – image to be unpadded

  • marginV (list) – padding margin for row, col, slc

Returns:

unpadded image

Return type:

np.ndarray

cerr.radiomics.preprocess.preProcessForRadiomics(scanNum, structNum, paramS, planC)[source]

This function applies pre-processing to scanNum as per settings specified in paramS dictionary.

Parameters:
  • scanNum (int)

  • structNum (int)

  • paramS (dict)

  • planC (plan_container.planC)

Returns:

unpadded image np.ndarray(dtype=int): maskBoundingBox3M np.ndarray(dtype=int): tuple: (x,y,z) grid corresponding to the processed image paramS (dictionary): Parameters for pre-processing. diagS (dictionary): Diagnostic features from the region of interest

Return type:

np.ndarray(dtype=float)

cerr.radiomics.run_length module

This module contains routines for calculation of Run Length texture features

cerr.radiomics.run_length.calcRLM(quantizedM, offsetsM, nL, rlmType=1)[source]

This function calculates the run-length matrix for the passed quantized image based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-run-length-based-features

Parameters:
  • quantizedM (np.ndarray(dtype=int)) – quantized 3d matrix obtained, for example, from radiomics.preprocess.imquantize_cerr

  • offsetsM (np.ndarray(dtype=int)) – Offsets for directionality/neighbors, obtained from radiomics.ibsi1.getDirectionOffsets

  • nL (int) – Number of gray levels. nL must be less than 65535.

  • rlmType (int) – flag, 1 or 2. 1: returns a single run length matrix by combining contributions from all offsets into one cooccurrence matrix. 2: returns a list of run length matrices, per row row of offsetsM.

Returns:

run-length matrix of size (nL x L) for rlmType = 1,

list of size equal to the number of directions for rlmType = 2.

The output can be passed to rlmToScalarFeatures to get RLM texture features.

Return type:

np.ndarray

cerr.radiomics.run_length.rlmToScalarFeatures(rlmM, numVoxels)[source]

This function calculates scalar texture features from run length matrix as per IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-run-length-based-features

Parameters:

rlmM (np.ndarray(dtype=int)) – run-length matrix of size (nL x L) for rlmType = 1, list of size equal to the number of directions for rlmType = 2.

Returns:

dictionary with scalar texture features as its

fields. Each field’s value is a vector containing the feature values for each list element of rlmM.

Return type:

dict

cerr.radiomics.shape module

shape module.

The shape module contains routines for calculation of shape features.

cerr.radiomics.shape.trimeshSurfaceArea(v, f)[source]

Routine to calculate surface area from vertices and faces of triangular mesh

Parameters:
  • v (numpy.array) – (numPoints x 3) vertices of triangular mesh

  • f (numpy.array) – (numFaces x 3) faces of triangular mesh

Returns:

Surface area

Return type:

float

cerr.radiomics.shape.vectorNorm3d(v)[source]
cerr.radiomics.shape.eig(a)[source]
cerr.radiomics.shape.sepsq(a, b)[source]
cerr.radiomics.shape.calcMaxDistBetweenPts(ptsM, distType)[source]

This routine calculates the maximum distance between the input points

Parameters:
  • ptsM (numpy.array) – (nunPoints x 3) coordinates of points

  • distType (str or Callable) – Type of distance. E.g. ‘euclidean’ as supported by scipy.spatial.distance.cdist

Returns:

Maximum distance between the input points

Return type:

float

cerr.radiomics.shape.calcShapeFeatures(mask3M, xValsV, yValsV, zValsV)[source]

Routine to calculate shape features for the inout mask and grid

Parameters:
  • mask3M (numpy.nparray) – Binary mask where 1s represent the segmentation

  • xValsV (numpy.nparray) – x-values i.e. coordinates of columns of input mask

  • yValsV (numpy.nparray) – y-values i.e. coordinates of rows of input mask

  • zValsV (numpy.nparray) – z-values i.e. coordinates of slices of input mask

Returns:

Dictionary containing shape features

Return type:

dict

cerr.radiomics.size_zone module

This module contains routines for claculation of Size Zone texture features

cerr.radiomics.size_zone.calcSZM(quantized3M, nL, szmType)[source]

This function calculates the Size Zone Matrix for the passed quantized image based on IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-size-zone-based-features

Parameters:
  • quantized3M (np.ndarray(dtype=int)) – quantized 3d matrix obtained, for example, from radiomics.preprocess.imquantize_cerr

  • nL (int) – Number of gray levels.

  • szmType – flag, 1 or 2. 1: 3D zones 2: 2D zones

Returns:

size-zone matrix of size (nL x L)

The output can be passed to szmToScalarFeatures to get SZM texture features.

Return type:

np.ndarray

cerr.radiomics.size_zone.szmToScalarFeatures(szmM, numVoxels)[source]

This function calculates scalar texture features from Size Zone Matrix as per IBSI definitions https://ibsi.readthedocs.io/en/latest/03_Image_features.html#grey-level-size-zone-based-features

Parameters:
  • szmM (np.ndarray(dtype=int)) – size-zone matrix of size (nL x L)

  • numVoxels (int) – number of voxels in the region of interest for szmM calculation

Returns:

dictionary with scalar texture features as its

fields. Each field’s value is a vector containing the feature values for each list element of rlmM.

Return type:

dict

cerr.radiomics.texture_filters module

This module contains definitions of image texture filters and a wrapper function to apply any of them. Supported filters include: “mean”, “sobel”, “LoG”, “gabor”, “gabor3d”, “laws”, “lawsEnergy” “rotationInvariantLaws”, “rotationInvariantLawsEnergy”

cerr.radiomics.texture_filters.meanFilter(scan3M, kernelSize, absFlag=False)[source]

Function to compute patchwise mean on input image using specified kernel size.

Parameters:
  • scan3M (np.ndarray) – Input image matrix.

  • kernelSize (np.array) – Size of filter kernel [nRow, nCol, nSlc].

  • absFlag (bool) – [optional, default:False] Flag to use absolute intensities.

Returns:

Mean filter response.

Return type:

np.ndarray(dtype=float)

cerr.radiomics.texture_filters.sobelFilter(scan3M)[source]

Function to compute gradient magnitude and direction using Sobel filter.

Parameters:

scan3M (np.ndarray) – 3D input scan matrix.

Returns:

gradient magnitude np.ndarray(dtype=float): gradient direction

Return type:

np.ndarray(dtype=float)

cerr.radiomics.texture_filters.LoGFilter(scan3M, sigmaV, cutoffV, voxelSizeV)[source]

Function to apply IBSI standard Laplacian of Gaussian (LoG) filter

Parameters:
  • scan3M (np.ndarray) – 3D scan matrix.

  • sigmaV (np.array, 1-D) – Gaussian smoothing widths, [sigmaRows,sigmaCols,sigmaSlc] in mm.

  • cutoffV (np.array 1-D) – Filter cutoffs [cutOffRow, cutOffCol cutOffSlc] in mm. Note: Filter size = 2.*cutoffV+1

  • voxelSizeV (np.array, 1-D) – Scan voxel dimensions [dx, dy, dz] in mm.

Returns:

IBSI-compatible Laplacian of Gaussian filter response.

Return type:

np.ndarray(dtype=float)

cerr.radiomics.texture_filters.gaborFilter(scan3M, sigma, wavelength, gamma, thetaV, aggS=None, radius=None, paddingV=None)[source]

Function to apply IBSI standard 2D Gabor filter

Parameters:
  • scan3M (np.ndarray) – 3D scan matrix.

  • sigma (float) – Std. deviation Gaussian envelope in no. voxels.

  • lambda (float) – wavelength in no. voxels.

  • gamma (float) – Spatial aspect ratio

  • thetaV (np.array, 1D) – Orientations in degrees

  • aggS (dict) – [optional, default=None] Parameters for averaging responses across orientations.

  • radius (np.array(dtype=int)) – [optional, default=None] Kernel radius in voxels [nRows nCols].

  • paddingV (np.array, 1D) – [optional, default=None] Amount of padding applied to scan in voxels [nRows nCols].

Returns:

IBSI-compatible 2D Gabor filter response.

Return type:

np.ndarray(dtype=float)

cerr.radiomics.texture_filters.gaborFilter3d(scan3M, sigma, wavelength, gamma, thetaV, aggS, radius=None, paddingV=None)[source]

Function to return Gabor filter responses aggregated across the 3 orthogonal planes (IBSI-compatible)

Parameters:
  • scan3M (np.ndarray) – 3D scan array.

  • sigma (int) – Std. dev. of Gaussian envelope in no. voxels.

  • lambda (int) – Wavelength in no. voxels.

  • gamma (float) – Spatial aspect ratio.

  • thetaV (list(dtype=float)) – Orientations in degrees.

  • aggS (dict) – Parameters for aggregation of responses across orientations and/or planes.

  • radius (np.array, 1D) – [optional, default=None] Kernel radii in voxels [nRows, nCols].

  • paddingV (np.aray, 1D) – [optional, default=None] Amount of padding applied to scan in voxels [nRows nCols].

Returns:

Gabor filter responses. hGaborPlane (dict): Gabor filter kernels for each plane.

Return type:

gaborOut (dict)

cerr.radiomics.texture_filters.get3dLawsMask(x, y, z)[source]

Function to get 3D Laws’ filter kernel (IBSI-compatible)

Parameters:
  • x (np.array, 1D) – Supported Laws filter coefficients applied along rows.

  • y (np.array, 1D) – Supported Laws filter coefficients applied along cols.

  • z (np.array, 1D) – Supported Laws filter coefficients applied along slices.

Returns:

3D Laws’ kernel.

Return type:

conved3M (np.ndarray)

cerr.radiomics.texture_filters.getLawsMasks(direction='all', filterType='all', normFlag=False)[source]

Function to get Laws filter kernels.

Parameters:
  • direction (string) – [optional, default=’all’] specifying ‘2d’, ‘3d’ or ‘All’

  • filterType (string) – [optional, default=’all’] specifying ‘3’, ‘5’, ‘all’, or a combination of any 2 (if 2d) or 3 (if 3d) of E3, L3, S3, E5, L5, S5.

  • normFlag (bool) – [optional, default=False] Flag to normalize filter coefficients when True. Normalization ensures average pixel in filtered image is as bright as the average pixel in the original image.

Returns:

Laws kernels for specified filter types and directions.

Return type:

lawsMasks (dict)

cerr.radiomics.texture_filters.lawsFilter(scan3M, direction, filterDim, normFlag)[source]

Function to compute Laws’ filter response (IBSI-compatible)

Parameters:
  • scan3M (np.ndarray) – 3D scan.

  • direction (string) – specifying ‘2d’, ‘3d’ or ‘All’

  • filterDim (string) – specifying ‘3’, ‘5’, ‘all’, or a combination of any 2 (if 2d) or 3 (if 3d) of E3, L3, S3, E5, L5, S5.

  • normFlag (bool) – Flag to normalize filter coefficients if True. Normalization ensures average pixel in filtered image is as bright as the average pixel in the original image)

Returns:

Filter responses for specified filter types and directions.

Return type:

lawsOut (dict)

cerr.radiomics.texture_filters.energyFilter(tex3M, mask3M, texPadFlag, texPadSizeV, texPadMethod, energyKernelSizeV, energyPadSizeV, energyPadMethod)[source]

Function to compute energy (local mean of absolute intensities)

Parameters:
  • tex3M (np.ndarray(dtype=float)) – Input filter response map

  • mask3M (np.ndarray(dtype=bool)) – Processed mask returned by radiomics.preprocess.

  • texPadFlag (bool) – Flag to indicate if padding was applied to compute Laws filter response.

  • texPadSizeV (np.array(dtype=int)) – Amount of padding applied to compute tex3M

  • texPadMethod (string) – Padding method applied prior to compute tex3M

  • energyKernelSizeV (np.array(dtype=int)) – Patch size used to calculate local energy [numRows numCols num_slc] in voxels.

  • energyPadMethod (string) – Padding method applied to calculate local energy.

  • energyPadSizeV – np.array(dtype=int) for amount padding applied to calculate local energy [numRows numCols num_slc] in voxels.

Returns:

Energy filter response.

Return type:

texEnergyPad3M (np.ndarray(dtype=float))

cerr.radiomics.texture_filters.lawsEnergyFilter(scan3M, mask3M, direction, filterDim, normFlag, lawsPadFlag, lawsPadSizeV, lawsPadMethod, energyKernelSizeV, energyPadSizeV, energyPadMethod)[source]

Function to compute local mean of absolute values of laws filter response

Args:

scan3M (np.ndarray): 3D scan. mask3M (np.ndarray(dtype=bool)): 3D mask. direction (string): Specifying ‘2d’, ‘3d’ or ‘All’. filterDim (string: Specifying ‘3’, ‘5’, ‘all’, or a combination

of any 2 (if 2d) or 3 (if 3d) of E3, L3, S3, E5, L5, S5.

normFlag (bool): Flag to normalize kernel coefficients if set to True.

Normalization ensures average pixel in filtered image is as bright as the average pixel in the original image.

lawsPadFlag (bool): Flag to indicate if padding was applied to compute

Laws filter response.

lawsPadSizeV (np.array(dtype=int)): Amount of padding applied to compute Laws

filter response.

lawsPadMethod (string): Padding method applied to compute Laws filter response. energyKernelSizeV (np.array(dtype=int)): Patch size used to calculate local energy

[numRows numCols num_slc] in voxels.

energyPadMethod (np.array(dtype=int)): Padding method applied to

calculate local energy.

energyPadSizeV (np.array(dtype=int)): Amount padding applied to

calculate local energy [numRows numCols num_slc] in voxels.

Returns:

Filter responses for specified directions and filter types.

Return type:

outS (dict)

cerr.radiomics.texture_filters.rotationInvariantLawsFilter(scan3M, direction, filterDim, normFlag, rotS)[source]

Function to return rotation-invariant Laws filter response.

Parameters:
  • scan3M (np.ndarray) – 3D scan.

  • direction (string) – Specifying ‘2d’, ‘3d’ or ‘All’.

  • filterDim (string) – Specifying ‘3’, ‘5’, ‘all’, or a combination of any 2 (if 2d) or 3 (if 3d) of E3, L3, S3, E5, L5, S5.

  • normFlag (bool) – Flag to normalize kernel coefficients if set to True. Normalization ensures average pixel in filtered image is as bright as the average pixel in the original image

  • rotS (dict) – Parameters for aggregating filter response across different orientations

Returns:

Laws filter response aggregated

across orientations as specified.

Return type:

out3M (np.ndarray(dtype=float))

cerr.radiomics.texture_filters.rotationInvariantLawsEnergyFilter(scan3M, mask3M, direction, filterDim, normFlag, lawsPadFlag, lawsPadSizeV, lawsPadMethod, energyKernelSizeV, energyPadSizeV, energyPadMethod, rotS)[source]

Function to return rotation-invariant Laws energy response.

Parameters:
  • scan3M (np.ndarray) – 3D scan.

  • direction (string) – Specifying ‘2d’, ‘3d’ or ‘All’.

  • filterDim (string) – Specifying ‘3’, ‘5’, ‘all’, or a combination of any 2 (if 2d) or 3 (if 3d) of E3, L3, S3, E5, L5, S5.

  • normFlag (bool) – Flag to normalize kernel coefficients if set to True. Normalization ensures average pixel in filtered image is as bright as the average pixel in the original image

  • lawsPadFlag (bool) – Flag to indicate if padding was applied to compute Laws filter response.

  • lawsPadSizeV (np.array(dtype=int)) – Amount of padding applied to compute Laws filter response.

  • lawsPadMethod (string) – Padding method applied to compute Laws filter response.

  • energyKernelSizeV (np.array(dtype=int)) – Patch size used to calculate local energy [numRows numCols num_slc] in voxels.

  • energyPadMethod (np.array(dtype=int)) – Padding method applied to calculate local energy.

  • energyPadSizeV (np.array(dtype=int)) – Amount padding applied to calculate local energy [numRows numCols num_slc] in voxels.

  • rotS (dict) – Parameters for aggregating filter response across different orientations

Returns:

Energy of Laws filter response

aggregated across orientations as specified.

Return type:

lawsEnergyAggPad3M (np.ndarray(dtype=float))

cerr.radiomics.texture_filters.getMatchFields(dictIn, *args)[source]

Return vals in list of dictionaries matching input key

cerr.radiomics.texture_filters.aggregateRotatedResponses(rotTexture, aggregationMethod)[source]

Function to aggregate textures from all orientations.

Parameters:
  • rotTexture – list of filter response maps to be aggregated, computed at different orientations.

  • aggregationMethod – string for aggregation method. May be ‘avg’,’max’, or ‘std’.

Returns:

Aggregated response map.

cerr.radiomics.texture_filters.rot3d90(arr3M, axis, angle)[source]

Function to rotate a 3D array 90 degrees around a specified axis.

Parameters:
  • arr3M (np.ndarray) – 3D scan.

  • axis (int) – Axis around which to rotate the array.

  • angle (float) – Angle of rotation in degrees.

Returns:

Rotated scan.

Return type:

rotArr3M (np.ndarray)

cerr.radiomics.texture_filters.rotate3dSequence(vol3M, index, sign)[source]

Function to apply pre-defined sequence of right angle rotations to 2D/3D arrays

Parameters:
  • vol3M (np.ndarray) – 3D scan.

  • index (int) – Multiplicative factor of 90 degrees.

  • sign (int) – Indicate direction of rotation (+1 or -1).

Returns:

Rotated scan.

Return type:

rotArr3M (np.ndarray)

cerr.radiomics.texture_filters.flipSequenceForWavelets()[source]
cerr.radiomics.texture_filters.rotationInvariantFilt(scan3M, mask3M, filter, *params)[source]

Function to apply filter over a range of orientations to approximate rotation invariant filter.

Parameters:
  • scan3M (np.ndarray) – 3D scan.

  • mask3M (np.ndarray(dtype=bool)) – 3D mask.

  • filter (dict) – Filter names and associated parameters.

  • params (dict)

Returns:

Rotation -invariant filter responses.

Return type:

aggS (dict)

cerr.radiomics.texture_utils module

cerr.radiomics.texture_utils.loadSettingsFromFile(settingsFile, scanNum=None, planC=None)[source]

Function to load filter parameters from user-input JSON.

Parameters:
  • settingsFile (string) – Path to JSON file.

  • scanNum (int) – [optional, default=None] Scan no. from which to extract additional parameters like voxel size.

  • planC (plan_container.planC) – [optional, default=None] pyCERR’s plan container object.

Returns:

Radiomics parameters parsed from JSON file. filterTypes (list): Texture filters specified in JSON file.

Return type:

paramS (dict)

cerr.radiomics.texture_utils.processImage(filterType, scan3M, mask3M, paramS)[source]

Function to process scan using selected filter and parameters :param filterType: Name of supported filter. :type filterType: string :param scan3M: 3D scan. :type scan3M: np.ndarray :param mask3M: 3D binary mask. :type mask3M: np.ndarray(dtype=bool) :param paramS: Parameters (read from JSON). :type paramS: dict

Returns:

Containing response maps for each filter type.

Return type:

outS(dict)

cerr.radiomics.texture_utils.generateTextureMapFromPlanC(planC, scanNum, strNum, configFilePath)[source]

Function to filter scan and import result to planC.

Parameters:
  • planC – pyCERR’s plan container object.

  • scanNum – int for index of scan to be filtered.

  • strNum – int for index of ROI.

  • configFilePath – string for path to JSON config file with filter parameters.

Returns:

pyCERR plan_container object with texture map as pseudo-scan.

Return type:

planC (plan_container.planC)

Module contents