API Reference
The simca
package is based on the work [Rou22]
classes
- class simca.CassiSystem.CassiSystem(system_config=None, system_config_path=None)[source]
Bases:
object
Class that contains the cassi system main attributes and methods
- apply_psf()[source]
Apply the PSF to the last measurement
- Returns:
last measurement cube convolved with by PSF (shape= R x C x W). Each slice of the 3D filtered scene is convolved with the PSF
- Return type:
numpy.ndarray
- create_coordinates_grid(nb_of_pixels_along_x, nb_of_pixels_along_y, delta_x, delta_y)[source]
Create a coordinates grid for a given number of samples along X and Y axis and a given pixel size
- Parameters:
nb_of_pixels_along_x (int) – number of samples along X axis
nb_of_pixels_along_y (int) – number of samples along Y axis
delta_x (float) – pixel size along X axis
delta_y (float) – pixel size along Y axis
- Returns:
X coordinates grid (numpy.ndarray) and Y coordinates grid (numpy.ndarray)
- Return type:
tuple
- generate_2D_pattern(config_pattern)[source]
Generate the coded aperture 2D pattern based on the “pattern” configuration file
- Parameters:
config_pattern (dict) – coded-aperture pattern configuration
- Returns:
coded-aperture 2D pattern based on the configuration file (shape = H x L)
- Return type:
numpy.ndarray
- generate_filtering_cube()[source]
Generate filtering cube : each slice of the cube is a propagated pattern interpolated on the detector grid
- Returns:
filtering cube generated according to the optical system & the pattern configuration (R x C x W)
- Return type:
numpy.ndarray
- generate_multiple_filtering_cubes(number_of_patterns)[source]
Generate multiple filtering cubes, each cube corresponds to a pattern, and for each pattern, each slice is a propagated coded apertureinterpolated on the detector grid
- Parameters:
number_of_patterns (int) – number of patterns to generate
- Returns:
filtering cubes generated according to the current optical system and the pattern configuration
- Return type:
list
- generate_multiple_patterns(config_pattern, number_of_patterns)[source]
Generate a list of coded aperture patterns based on the “pattern” configuration file
- Parameters:
config_pattern (dict) – pattern configuration
number_of_patterns (int) – number of patterns to generate
- Returns:
coded aperture patterns (numpy.ndarray) generated according to the configuration file
- Return type:
list
- image_acquisition(use_psf=False, chunck_size=50)[source]
Run the acquisition/measurement process depending on the cassi system type
- Parameters:
chunck_size (int) – default block size for the interpolation
- Returns:
compressed measurement (R x C)
- Return type:
numpy.ndarray
- interpolate_dataset_along_wavelengths(new_wavelengths_sampling, chunk_size)[source]
Interpolate the dataset cube along the wavelength axis to match the system sampling
- Parameters:
new_wavelengths_sampling (numpy.ndarray) – new wavelengths on which to interpolate the dataset (shape = W)
chunk_size (int) – chunk size for the multiprocessing
- Returns:
interpolated dataset cube along the wavelength axis (shape = R_dts x C_dts x W)
- Return type:
numpy.ndarray
- load_dataset(directory, dataset_name)[source]
Loading the dataset and related attributes
- Parameters:
directory (str) – name of the directory containing the dataset
dataset_name (str) – dataset name
- Returns:
a list containing the dataset (shape= R_dts x C_dts x W_dts), the corresponding wavelengths (shape= W_dts), the labeled dataset, the label names and the ignored labels
- Return type:
list
- multiple_image_acquisitions(use_psf=False, nb_of_filtering_cubes=1, chunck_size=50)[source]
Run the acquisition process depending on the cassi system type
- Parameters:
chunck_size (int) – default block size for the dataset
- Returns:
list of compressed measurements (list of numpy.ndarray of size R x C)
- Return type:
list
- propagate_coded_aperture_grid(X_input_grid=None, Y_input_grid=None)[source]
Propagate the coded_aperture pattern through one CASSI system
- Parameters:
X_input_grid (numpy.ndarray) – x coordinates grid
Y_input_grid (numpy.ndarray) – y coordinates grid
- Returns:
propagated coded aperture x coordinates grid in the detector plane (3D numpy.ndarray), propagated coded aperture y coordinates grid in the detector plane (3D numpy.ndarray), 1D array of propagated coded aperture x coordinates (numpy.ndarray), 1D array of system wavelengths (numpy.ndarray)
- Return type:
tuple
- save_acquisition(config_pattern, config_acquisition)[source]
Save the all data related to an acquisition
- Parameters:
config_pattern (dict) – configuration dictionary related to pattern generation
config_acquisition (dict) – configuration dictionary related to acquisition parameters
- set_up_system(system_config_path=None, system_config=None)[source]
Loading system config & initializing the grids coordinates for the coded aperture and the detector
- Parameters:
system_config_path (str) – path to the configs file
system_config (dict) – system configuration
- update_config(system_config_path=None, system_config=None)[source]
Update the system configuration file and re-initialize the grids for the coded aperture and the detector
- Parameters:
system_config_path (str) – path to the configs file
system_config (dict) – system configuration
- Returns:
updated system configuration
- Return type:
dict
- class simca.OpticalModel.OpticalModel(system_config)[source]
Bases:
object
Class that contains the optical model caracteristics and propagation models
- calculate_alpha_c()[source]
Calculate the relative angle of incidence between the lenses and the dispersive element
- Returns:
angle of incidence
- Return type:
float
- calculate_central_dispersion()[source]
Calculate the dispersion related to the central pixel of the coded aperture
- Returns:
spectral dispersion of the central pixel of the coded aperture
- Return type:
numpy.float
- calculate_minimum_deviation(n, A)[source]
minimum deviation angle of a prism of index n and apex angle A
- Parameters:
n (float or numpy.ndarray) – index of the prism – no units
A (float) – apex angle of the prism – in radians
- Returns:
minimum deviation angle – in radians
- Return type:
float or numpy.ndarray
- check_if_sampling_is_sufficiant()[source]
Check if the sampling is sufficiant to avoid aliasing.
- Returns:
number of sample points per pixel
- Return type:
float
- generate_2D_gaussian(radius, sample_size_x, sample_size_y, nb_of_samples)[source]
Generate a 2D Gaussian of a given radius
- Parameters:
radius (float) – radius of the Gaussian
sample_size_x (float) – size of each sample along the X axis
sample_size_y (float) – size of each sample along the Y axis
nb_of_samples (int) – number of samples along each axis
- Returns:
2D Gaussian shape array
- Return type:
numpy.ndarray
- generate_psf(type, radius)[source]
Generate a PSF
- Parameters:
type (str) – type of PSF to generate
radius (float) – radius of the PSF
- Returns:
PSF generated (shape = R x C)
- Return type:
numpy.ndarray
- get_incident_angle_min_dev(A, D_m)[source]
Calculate the angle of incidence corresponding to the minimum deviation angle
- Parameters:
A (float) – apex angle of the prism – in radians
D_m (float) – minimum deviation angle – in radians
- Returns:
angle of incidence corresponding to minimum of deviation – in radians
- Return type:
float
- model_Grating_angle_to_angle(k_in, lba, m, G)[source]
Model of the grating
- Parameters:
k_in (numpy.ndarray) – wave vector of the incident ray (shape = 3 x N)
lba (numpy.ndarray) – wavelengths (shape = N) – in nm
m (float) – diffraction order of the grating – no units
G (float) – lines density of the grating – in lines/mm
- Returns:
wave vector of the outgoing ray (shape = 3 x N)
- Return type:
numpy.ndarray
- model_Lens_angle_to_position(k_in, F)[source]
Model of the lens : angle to position
- Parameters:
k_in (numpy.ndarray) – wave vector of the incident ray (shape = 3 x N)
F (float) – focal length of the lens – in um
- Returns:
position in the image plane (X,Y) – in um
- Return type:
tuple
- model_Lens_pos_to_angle(x_obj, y_obj, F)[source]
Model of the lens : position to angle
- Parameters:
x_obj (numpy.ndarray) – position X in the image plane (shape = N) – in um
y_obj (numpy.ndarray) – position Y in the image plane (shape = N) – in um
F (float) – focal length of the lens – in um
- Returns:
wave vector of the outgoing ray (shape = 3 x N)
- Return type:
numpy.ndarray
- model_Prism_angle_to_angle(k0, n, A)[source]
Ray tracing through the prism
- Parameters:
k0 (numpy.ndarray) – wave vector of the incident ray (shape = 3 x N)
n (numpy.ndarray) – refractive index of the prism (shape = N)
A (float) – angle of the prism – in radians
- Returns:
wave vector of the outgoing ray (shape = 3 x N)
- Return type:
numpy.ndarray
- propagate_through_arm(X_vec_in, Y_vec_in, n, lba)[source]
Propagate the light through one system arm : (lens + dispersive element + lens)
- Parameters:
X_vec_in (numpy.ndarray) – X coordinates of the coded aperture pixels (1D array)
Y_vec_in (numpy.ndarray) – Y coordinates of the coded aperture pixels (1D array)
n (numpy.ndarray) – refractive indexes of the system (at the corresponding wavelength)
lba (numpy.ndarray) – wavelengths
- Returns:
flatten arrays corresponding to the propagated X and Y coordinates
- Return type:
tuple
- propagation_with_distorsions(X_input_grid, Y_input_grid)[source]
Propagate the coded aperture coded_aperture through one CASSI system
- Parameters:
X_input_grid (numpy.ndarray) – x coordinates grid
Y_input_grid (numpy.ndarray) – y coordinates grid
- Returns:
X coordinates of the propagated coded aperture grids, Y coordinates of the propagated coded aperture grids
- Return type:
tuple
- propagation_with_no_distorsions(X_input_grid, Y_input_grid)[source]
Vanilla Propagation model used in most cassi acquisitions simulation.
- Parameters:
X_input_grid (numpy.ndarray) – X coordinates of the grid to be propagated (2D)
Y_input_grid (numpy.ndarray) – Y coordinates of the grid to be propagated (2D)
- Returns:
X coordinates grids of the propagated coded apertures, Y coordinates grids of the propagated coded apertures
- Return type:
tuple
- sellmeier(lambda_, glass_type='BK7')[source]
Evaluating the refractive index value of a prism for a given lambda based on Sellmeier equation
- Parameters:
lambda (numpy.ndarray of float) – wavelength in nm
- Returns:
index value corresponding to the input wavelength
- Return type:
numpy.ndarray of float
- set_optical_config(config)[source]
Set the optical model configuration
- Parameters:
config (dict) – configuration file
- set_wavelengths(wavelength_min, wavelength_max, nb_of_spectral_samples)[source]
Set the wavelengths range of the optical system
- Parameters:
wavelength_min (float) – minimum wavelength of the system
wavelength_max (float) – maximum wavelength of the system
nb_of_spectral_samples (int) – number of spectral samples of the system
Returns:
- simplified_grating_in_out(alpha, lba, m, G)[source]
Model 1D of the grating in the dispersion direction
- Parameters:
alpha (numpy.ndarray or float) – angle of the incident ray (shape = N) – in radians
lba (numpy.ndarray or float) – wavelengths (shape = N) – in nm
m (float) – diffraction order of the grating – no units
G (float) – lines density of the grating – in lines/mm
- Returns:
angle of the outgoing ray (shape = N) – in radians
- Return type:
numpy.ndarray
functions
- simca.functions_scenes.explore_spectrums(img, complete_gt, class_names, ignored_labels=None)[source]
Plot sampled spectrums with mean + std for each class.
- Parameters:
img – 3D hyperspectral image
complete_gt – 2D array of labels
class_names – list of class names
ignored_labels (optional) – list of labels to ignore
- Returns:
dict of mean spectrum by class
- Return type:
mean_spectrums
- simca.functions_scenes.get_dataset(dataset_name, folder='./datasets/')[source]
Gets the dataset specified by name and return the related components. :param dataset_name: the name of the dataset :type dataset_name: str :param folder: folder where the datasets are stored, defaults to “./datasets/” :type folder: str
- Returns:
3D hyperspectral image (WxHxB) numpy.ndarray: 2D array of labels (integers) list: list of class names ignored_labels: list of int classes to ignore
- Return type:
numpy.ndarray
- simca.functions_scenes.interpolate_data_along_wavelength(data, current_sampling, new_sampling, chunk_size=50)[source]
Interpolate the input 3D data along a new sampling in the third axis.
- Parameters:
data (numpy.ndarray) – 3D data to interpolate
current_sampling (numpy.ndarray) – current sampling for the 3rd axis
new_sampling (numpy.ndarray) – new sampling for the 3rd axis
chunk_size (int) – size of the chunks to use for the interpolation
- simca.functions_patterns_generation.FindLargestVoid(BinaryPattern, StandardDeviation)[source]
- This function returns the indices of the largest void in the given binary
pattern as defined by Ulichney.
- param BinaryPattern A boolean array (should be two-dimensional although the
implementation works in arbitrary dimensions).
- param StandardDeviation The standard deviation used for the Gaussian filter
in pixels. This can be a single float for an isotropic Gaussian or a tuple with one float per dimension for an anisotropic Gaussian.
- eturn A flat index i such that BinaryPattern.flat[i] corresponds to the
largest void. By definition this is a majority pixel.
sa GetVoidAndClusterBlueNoise
- simca.functions_patterns_generation.FindTightestCluster(BinaryPattern, StandardDeviation)[source]
- Like FindLargestVoid() but finds the tightest cluster which is a minority
pixel by definition.
sa GetVoidAndClusterBlueNoise
- simca.functions_patterns_generation.GetVoidAndClusterBlueNoise(OutputShape, StandardDeviation=1.5, InitialSeedFraction=0.1)[source]
- Generates a blue noise dither array of the given shape using the method
proposed by Ulichney [1993] in “The void-and-cluster method for dither array generation” published in Proc. SPIE 1913.
- param OutputShape The shape of the output array. This function works in
arbitrary dimension, i.e. OutputShape can have arbitrary length. Though it is only tested for the 2D case where you should pass a tuple (Height,Width).
- param StandardDeviation The standard deviation in pixels used for the
Gaussian filter defining largest voids and tightest clusters. Larger values lead to more low-frequency content but better isotropy. Small values lead to more ordered patterns with less low-frequency content. Ulichney proposes to use a value of 1.5. If you want an anisotropic Gaussian, you can pass a tuple of length len(OutputShape) with one standard deviation per dimension.
- param InitialSeedFraction The only non-deterministic step in the algorithm
marks a small number of pixels in the grid randomly. This parameter defines the fraction of such points. It has to be positive but less than 0.5. Very small values lead to ordered patterns, beyond that there is little change.
- eturn An integer array of shape OutputShape containing each integer from 0
to np.prod(OutputShape)-1 exactly once.
- simca.functions_patterns_generation.generate_blue_noise_type_1_pattern(shape)[source]
Generate blue noise (high frequency pseudo-random) type pattern
- Parameters:
shape (tuple of int) – shape of the pattern
- Returns:
binary blue noise type pattern
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.generate_blue_noise_type_2_pattern(shape, std=1.5, initial_seed_fraction=0.1)[source]
Generate blue noise pattern according to the void-and-cluster method proposed by Ulichney [1993] in “The void-and-cluster method for dither array generation” published in Proc. SPIE 1913.
- Parameters:
shape (tuple) – size of the pattern
std (float) – standard deviation in pixels used for the Gaussian filter
initial_seed_fraction (float) – Initial fraction of marked pixels in the grid. Has to be less than 0.5. Very small values lead to ordered patterns
- Returns:
float blue noise pattern
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.generate_ln_orthogonal_pattern(size, W, N)[source]
Generate a Length-N orthogonal pattern according to https://hal.laas.fr/hal-02993037
- Parameters:
size (tuple) – size of the pattern
W (int) – number of wavelengths in the scene
N (int) – number of acquisitions
- Returns:
length-N orthogonal pattern of shape = size[0] x (size[1]+W-1) x N):
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.generate_orthogonal_pattern(size, W, N)[source]
Generate an orthogonal pattern according to https://hal.laas.fr/hal-02993037
- Parameters:
size (list of int) – size of the pattern
W (int) – number of wavelengths in the scene
N (int) – number of acquisitions
- Returns:
orthogonal pattern : shape = size[0] x (size[1]+W-1) x N):
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.generate_random_pattern(shape, ROM)[source]
Generate a random pattern with a given rate of open/close mirrors
- Parameters:
shape (tuple of int) – shape of the pattern
ROM (float) – ratio of open mirrors
- Returns:
random pattern
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.generate_slit_pattern(shape, slit_position, slit_width)[source]
Generate a slit pattern that starts at the center of the image and goes to the right as slit position increases.
- Parameters:
shape (tuple) – shape of the pattern
slit_position (int) – position of the slit in relation to the central column
slit_width (int) – width of the slit in pixels
- Returns:
slit pattern
- Return type:
numpy.ndarray
- simca.functions_patterns_generation.load_custom_list_of_patterns(shape, patterns_path)[source]
Load custom list of patterns from h5 file. If the patterns are not the same size as the coded aperture, crop from the center of the loaded patterns.
- Parameters:
shape (tuple) – size of the pattern
patterns_path (str) – path to the h5 file containing the patterns
- Returns:
list of patterns
- Return type:
list
- simca.functions_patterns_generation.load_custom_pattern(shape, pattern_path)[source]
Load custom pattern from h5 file. If the pattern is not the same size as the coded aperture, crop from the center of the loaded pattern.
- Parameters:
shape (tuple) – size of the pattern
pattern_path (str) – path to the h5 file containing the pattern
- Returns:
float blue noise pattern
- Return type:
numpy.ndarray
- simca.functions_acquisition.crop_center(array, nb_of_pixels_along_x, nb_of_pixels_along_y)[source]
Crop the given array to the given size, centered on the array
- Parameters:
array (numpy.ndarray) – 2D array to be cropped
nb_of_pixels_along_x (int) – number of samples to keep along the X axis
nb_of_pixels_along_y (int) – number of samples to keep along the Y axis
- Returns:
cropped array
- Return type:
numpy.ndarray
- simca.functions_acquisition.generate_dd_measurement(scene, filtering_cube, chunk_size)[source]
Generate DD-CASSI type system measurement from a scene and a filtering cube. ref : “Single-shot compressive spectral imaging with a dual-disperser architecture”, M.Gehm et al., Optics Express, 2007
- Parameters:
scene (numpy.ndarray) – observed scene (shape = R x C x W)
filtering_cube (numpy.ndarray) – filtering cube of the instrument for a given pattern (shape = R x C x W)
chunk_size (int) – size of the spatial chunks in which the Hadamard product is performed
- Returns:
filtered scene (shape = R x C x W)
- Return type:
numpy.ndarray
- simca.functions_acquisition.generate_sd_measurement_cube(filtered_scene, X_input, Y_input, X_target, Y_target, grid_type, interp_method)[source]
Generate SD measurement cube from the coded aperture and the scene. For Single-Disperser CASSI systems, the scene is filtered then propagated in the detector plane.
- Parameters:
filtered_scene (numpy.ndarray) – filtered scene (shape = R x C x W)
- Returns:
SD measurement cube (shape = R x C x W)
- Return type:
numpy.ndarray
- simca.functions_acquisition.interpolate_data_on_grid_positions(data, X_init, Y_init, X_target, Y_target, grid_type='unstructured', interp_method='linear')[source]
Interpolate data on a single 2D grid defined by X_target and Y_target
- Parameters:
data (numpy.ndarray) – data to interpolate (3D or 2D)
X_init (numpy.ndarray) – X coordinates of the initial grid (3D)
Y_init (numpy.ndarray) – Y coordinates of the initial grid (3D)
X_target (numpy.ndarray) – X coordinates of the target grid (2D)
Y_target (numpy.ndarray) – Y coordinates of the target grid (2D)
grid_type (str) – type of the target grid (default = “unstructured”, other option = “regular”)
interp_method (str) – interpolation method (default = “linear”)
- Returns:
3D data interpolated on the target grid
- Return type:
numpy.ndarray
- simca.functions_acquisition.match_dataset_labels_to_instrument(dataset_labels, filtering_cube)[source]
Match the size of the dataset labels to the size of the filtering cube. Either by padding or by cropping
- Parameters:
dataset_labels (numpy.ndarray) – dataset labels (shape = R_dts x C_dts)
filtering_cube (numpy.ndarray) – filtering cube of the instrument
- Returns:
scene labels (shape = R x C)
- Return type:
numpy.ndarray
- simca.functions_acquisition.match_dataset_to_instrument(dataset, filtering_cube)[source]
Match the size of the dataset to the size of the filtering cube. Either by padding or by cropping
- Parameters:
dataset (numpy.ndarray) – dataset
filtering_cube (numpy.ndarray) – filtering cube of the instrument
- Returns:
observed scene (shape = R x C x W)
- Return type:
numpy.ndarray
- simca.functions_acquisition.worker_regulargrid(args)[source]
Process to parallellize the structured griddata interpolation between the propagated grid (mask and the detector grid Note : For now it is identical to the unstructured method but it could be faster …
- Parameters:
args (tuple) – containing the following elements: X_init_2D, Y_init_2D, data_2D, X_target_2D, Y_target_2D
- Returns:
2D array of the data interpolated on the target grid
- Return type:
numpy.ndarray
- simca.functions_acquisition.worker_unstructured(args)[source]
Process to parallellize the unstructured griddata interpolation between the propagated grid (mask and the detector grid
- Parameters:
args (tuple) – containing the following elements: X_init_2D, Y_init_2D, data_2D, X_target_2D, Y_target_2D
- Returns:
2D array of the data interpolated on the target grid
- Return type:
numpy.ndarray
- simca.functions_general_purpose.initialize_acquisitions_directory(config)[source]
Initialize the directory where the results of the acquisition will be stored
- Parameters:
config (dict) – a configuration dictionary containing storing information
- Returns:
path to the directory where the results will be stored
- Return type:
str
- simca.functions_general_purpose.load_yaml_config(file_path)[source]
Load a YAML configuration file as a dictionary
- Parameters:
file_path (str) – path to the YAML configuration file
- Returns:
configuration dictionary
- Return type:
dict
- simca.functions_general_purpose.rotation_x(theta)[source]
Rotate 3D matrix around the X axis
- Parameters:
theta (float) – Input angle (in rad)
- Returns:
2D rotation matrix
- Return type:
numpy.ndarray
- simca.functions_general_purpose.rotation_y(theta)[source]
Rotate 3D matrix around the Y axis
- Parameters:
theta (float) – Input angle (in rad)
- Returns:
2D rotation matrix
- Return type:
numpy.ndarray
- simca.functions_general_purpose.rotation_z(theta)[source]
Rotate 3D matrix around the Z axis
- Parameters:
theta (float) – Input angle (in rad)
- Returns:
2D rotation matrix
- Return type:
numpy.ndarray
- simca.functions_general_purpose.save_config_file(config_file_name, config_file, result_directory)[source]
Save a configuration file in a YAML file
- Parameters:
config_file_name (str) – name of the file
config_file (dict) – configuration file to save
result_directory (str) – path to the directory where the results will be stored