C++ API#

Fast code in skmap is written in C++ and exposed to Python with thin wrappers. This is a lower-level api. All Python-exposed functions handle setting up GDAL and then call a function on an array class.

Input-Output#

group IO / Storage

Functions for reading/writing raster or matrix data to/from disk.

This group contains all I/O related functions, including GDAL-backed reading/writing and Python-friendly wrapper functions.

Functions

void readData(Eigen::Ref<MatFloat> data, const uint_t n_threads, const std::vector<std::string> &file_locs, const std::vector<uint_t> perm_vec, const uint_t x_off, const uint_t y_off, const uint_t x_size, const uint_t y_size, const std::vector<int> bands_list, py::dict conf_GDAL, std::optional<float_t> value_to_mask, std::optional<float_t> value_to_set)#

see IoArray::readData

Currently only supports Float32 data

void writeData(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::vector<std::string> base_files, std::string base_folder, std::vector<std::string> file_names, std::vector<uint_t> data_indices, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, float_t no_data_value, std::string gdal_data_type_str, std::optional<std::string> bash_compression_command, std::optional<std::vector<std::string>> seaweed_path)#

Python-friendly wrapper for writing raster data using GDAL and Eigen.

This function wraps the IoArray::writeData method, handling GDAL setup, nodata type dispatching, and optional post-processing (compression or remote storage).

Other parameters are passed to IoArray::WriteData: IoArray::WriteData

Note

For full details on the writing behavior, memory handling, and parallelization, see IoArray::writeData: IoArray::writeData

Parameters:
  • data – Eigen::matrix reference containing the data to write.

  • n_threads – Number of threads for parallel I/O operations.

  • conf_GDAL – Python dictionary of GDAL configuration options.

void extractOverlay(Eigen::Ref<MatFloat> data, const uint_t n_threads, const std::vector<uint_t> pix_block_ids, const std::vector<uint_t> pix_inblock_idxs, const std::vector<uint_t> unique_blocks_ids_comb, const std::vector<uint_t> key_layer_ids_comb, Eigen::Ref<MatFloat> data_overlay)#

wrapper for IoArray::extractOverlay

void warpTile(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::string tile_path, std::string mosaic_path, std::string resample)#

wrapper for IoArray::warpTile

Deprecated:

Please use gdal VRTs

void getLatLonArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::string file_loc, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size)#
void readData(std::vector<std::string> file_locs, std::vector<uint_t> perm_vec, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, GDALDataType read_type, std::vector<int> bands_list, std::optional<float_t> value_to_mask, std::optional<float_t> value_to_set)#

Reads multiple raster datasets into the internal matrix in parallel.

This function reads specified regions from a set of raster files (file_locs) into the internal m_data matrix, using a permutation vector to control the row order. Reading is performed in parallel via the parRowPerm helper, which calls readDataCore for each file.

The function defines a lambda readTiff that wraps readDataCore, which reads a single raster into a given row buffer. parRowPerm executes this lambda in parallel over the permutation vector, efficiently filling the internal matrix m_data.

Note

Both value_to_mask and value_to_set are forwarded to readDataCore for per-pixel masking and replacement.

Parameters:
  • file_locs – Vector of file paths or URLs of raster datasets to read.

  • perm_vec – Vector of row indices specifying the permutation order for writing into m_data.

  • x_off – Horizontal offset of the reading window in the rasters.

  • y_off – Vertical offset of the reading window in the rasters.

  • x_size – Width of the window to read.

  • y_size – Height of the window to read.

  • read_type – GDALDataType specifying the type to read (e.g., GDT_Float32).

  • bands_list – List of band indices to read from each dataset.

  • value_to_mask – Optional value to treat as masked (e.g., nodata). If not provided, no masking is performed.

  • value_to_set – Optional value to replace masked values with.

Throws:

skmapAssertIfTrue – if the number of columns in m_data is smaller than the requested reading window (x_size * y_size).

template<typename T>
inline void writeData(std::vector<std::string> base_files, std::string base_folder, std::vector<std::string> file_names, std::vector<uint_t> data_indices, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, GDALDataType write_type, T no_data_value, std::optional<std::string> bash_compression_command, std::optional<std::vector<std::string>> seaweed_path)#

Write portions of the internal matrix m_data to multiple GeoTIFF files.

This function writes selected rows of the m_data matrix into GeoTIFF files based on a set of base files and offsets. It supports in-place row casting, NoData masking, optional bash compression commands, and optional uploading to a remote storage (e.g., SeaweedFS).

For each row specified in data_indices, the function:

  1. Opens the corresponding base file and copies geotransform, projection, and spatial reference.

  2. Masks NaN values with no_data_value.

  3. Writes the row to a temporary raster if compression is requested, otherwise directly to the final file.

  4. Applies optional bash compression.

  5. Uploads the file to seaweed_path if provided and removes local copies as needed.

Template Parameters:

T – The data type to write to disk (e.g., float, double).

Parameters:
  • base_files – Vector of file paths to base raster files to copy metadata from.

  • base_folder – Folder where output files will be written.

  • file_names – Names of the output files (without folder path).

  • data_indices – Indices of rows in m_data corresponding to each output file.

  • x_off – X-offset within the raster to start writing data.

  • y_off – Y-offset within the raster to start writing data.

  • x_size – Width of the region to write.

  • y_size – Height of the region to write.

  • write_type – GDAL data type to write (e.g., GDT_Float32).

  • no_data_value – Value to use for missing or NaN cells.

  • bash_compression_command – Optional: shell command to compress the output files (e.g., gdal_translate or gzip).

  • seaweed_path – Optional: vector of remote storage paths to upload each output file.

void extractOverlay(std::vector<uint_t> pix_block_ids, std::vector<uint_t> pix_inblock_idxs, std::vector<uint_t> unique_blocks_ids_comb, std::vector<uint_t> key_layer_ids_comb, Eigen::Ref<MatFloat> data_overlay)#

Extract overlay values from raster blocks for a set of points.

For each pixel (given by pix_block_ids and pix_inblock_idxs), this function finds the matching block and copies the corresponding raster values from m_data into data_overlay.

The operation is parallelized across pixels.

Parameters:
  • pix_block_ids – Vector of block IDs corresponding to each pixel.

  • pix_inblock_idxs – Vector of linear indices of each pixel within its block.

  • unique_blocks_ids_comb – Vector of all unique block IDs for the current chunk.

  • key_layer_ids_comb – Vector of layer indices corresponding to each block-layer combination.

  • data_overlay – Output array (layers x pixels) where extracted values will be stored.

void warpTile(std::string ref_tile_path, std::string mosaic_path, std::string resample)#

Warp a single-band mosaic to match a reference tile.

Deprecated:

This function is deprecated. Use GDAL VRTs instead, which provide on-the-fly mosaicking and warping, support multiple bands, and avoid loading full rasters into memory.

This function opens the reference tile and mosaic, computes the target dimensions, sets up GDALWarpOptions, and performs an in-memory warp to match the reference. The resulting raster is read into m_data.

Note

Only works for single-band float32 rasters.

Parameters:
  • ref_tile_path – File path to reference raster tile.

  • mosaic_path – File path to mosaic raster.

  • resample – Resampling method (“nearest”, “bilinear”, etc.)

Data mangling#

group Data mangling

Functions for changing the shape of data for efficient parallel processing

For canonical usage examples, see unit tests

Functions

void copyVecInMatrixRow(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<VecFloat> in_vec, uint_t row_idx)#

see TransArray::copyVecInMatrixRow

void selArrayRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> row_select)#

see TransArray::selArrayRows

void selArrayCols(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_select)#

see TransArray::selArrayCols

void expandArrayRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> row_select)#

see TransArray::expandArrayRows

void expandArrayCols(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_select)#

see TransArray::expandArrayCols

void reorderArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> indices_matrix)#

see TransArray::reorderArray

void transposeArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data)#

simple transpose, see TransArray::transposeArray for details

void transposeReorderArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> permutation_matrix)#

see TransArray::transposeReorderArray

void inverseReorderArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> indices_matrix)#

see TransArray::inverseReorderArray

void copyVecInMatrixRow(Eigen::Ref<VecFloat> in_vec, uint_t row_idx)#

Copy a vector into the matrix.

If you created a vector in Python, this is the way to add it to the matrix

Parameters:
  • in_vec – Vector to copy

  • row_idx – Matrix row to copy into

void selArrayRows(Eigen::Ref<MatFloat> out_data, std::vector<uint_t> row_select)#

Select row_select rows to put in out_data in-order.

_images/sel_expand_rows.svg

out_data.rows() MUST be >=row_select.len() and out_data.cols() MUST match the width.

Parameters:
  • out_data – output matrix

  • row_select – rows to select from initial matrix

void selArrayCols(Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_select)#

Select col_select cols to put in out_data in-order.

out_data.cols() MUST be >=row_select.len() and out_data.rows() MUST match the height.

Parameters:
  • out_data – output matrix

  • col_select – columns to select from initial matrix

void expandArrayRows(Eigen::Ref<MatFloat> out_data, std::vector<uint_t> row_select)#

Inserts rows of self into out_data at row_select indices.

_images/sel_expand_rows.svg

length of row_select MUST equal m_data.rows()

Parameters:
  • out_data – Output array

  • row_select – indices at which to put m_data

void expandArrayCols(Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_select)#

Inserts columns of self into out_data at col_select indices.

length of col_select MUST equal m_data.cols()

Parameters:
  • out_data – Output array

  • col_select – indices at which to put m_data

void reorderArray(Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> indices_matrix)#

Expands the indices of indices_matrix to get a wider matrix.

_images/reorder_transpose_inverse.svg

This function is almost always slower than transposeReorderArray

Parameters:
  • out_data – The wider output matrix

  • indices_matrix – a rectangular matrix of indices to expand

void transposeArray(Eigen::Ref<MatFloat> out_data)#

Transposes the array.

_images/reorder_transpose_inverse.svg

Transposes the array

Parameters:

out_data – Array to be filled with transposed values

void transposeReorderArray(Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> permutation_matrix)#

reorders and transposes in one step

_images/reorder_transpose_inverse.svg

Warning

Currently not working as intended

Parameters:
  • out_data – output data matrix

  • permutation_matrix – permutation matrix

void inverseReorderArray(Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> indices_matrix)#

Inverse of TransArray::reorderArray

_images/reorder_transpose_inverse.svg

Selects blocks from the array based on out_data width and indices_matrix

Parameters:
  • out_data – Narrower output matrix

  • indices_matrix – matrix of indices to be appended

Data manipulation#

group Data manipulation

parallel operations on data arrays

Functions

void fillArray(Eigen::Ref<MatFloat> data, const uint_t n_threads, float_t val)#

see TransArray::fillArray

void maskNan(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, float_t new_value_in_data)#

see TransArray::maskNan

void maskNanRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, Eigen::Ref<VecFloat> new_value_vec)#

see TransArray::maskNanRows

void maskData(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, Eigen::Ref<MatFloat> mask, float_t value_of_mask_to_mask, float_t new_value_in_data)#

see TransArray::maskData

void maskDataRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, Eigen::Ref<MatFloat> mask, float_t value_of_mask_to_mask, float_t new_value_in_data)#

see TransArray::maskDataRows

void swapRowsValues(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, float_t value_to_mask, float_t new_value)#

see TransArray::swapRowsValues

void hadamardProduct(Eigen::Ref<MatFloat> out, const uint_t n_threads, Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#

see TransArray::hadamardProduct

void offsetAndScale(Eigen::Ref<MatFloat> data, const uint_t n_threads, float_t offset, float_t scaling)#

see TransArray::offsetAndScale

void scaleAndOffset(Eigen::Ref<MatFloat> data, const uint_t n_threads, float_t offset, float_t scaling)#

see TransArray::scaleAndOffset

void offsetsAndScales(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> row_select, Eigen::Ref<VecFloat> offsets, Eigen::Ref<VecFloat> scalings)#

see TransArray::offsetsAndScales

void castFloat32ToFloat64(Eigen::Ref<Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> data, const uint_t n_threads, Eigen::Ref<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> out_data)#

Cast a float array to a double array.

Will error on a size mismatch

Parameters:
  • data – Input data

  • n_threads – number of threads to use

  • out_data – output array

void castFloat64ToFloat32(Eigen::Ref<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> data, const uint_t n_threads, Eigen::Ref<Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> out_data)#

Cast a double array to a float array.

Will error on a size mismatch

Parameters:
  • data – Input data

  • n_threads – number of threads to use

  • out_data – output array

void fillArray(float_t val)#

Fill array with val

Parameters:

val – The fill value

void maskNan(std::vector<uint_t> row_select, float_t new_value_in_data)#

Fill NANs in row_select with new_value_in_data

If the fill value should be different per row, see TransArray::maskNanRows

Parameters:
  • row_select – vec of row indices to mask

  • new_value_in_data – mask value

void maskNanRows(std::vector<uint_t> row_select, Eigen::Ref<VecFloat> new_value_vec)#

Fill NANs in row_select with corresponding value in new_value_vec

The fill value will correspond to the index in row_select

Parameters:
  • row_select – vec of row indices to mask

  • new_value_vec – mask values per index

void maskData(std::vector<uint_t> row_select, Eigen::Ref<MatFloat> masks, float_t value_of_mask_to_mask, float_t new_value_in_data)#

fills data with new_value_in_data where mask==value_of_mask_to_mask

The mask row corresponds to the index of row_select. For using the same mask for all rows, see TransArray::maskDataRows.

Warning

maskDataRows and maskData may swap names

Parameters:
  • row_select – vec of rows to apply masks to

  • masks – masks, one for each row

  • value_of_mask_to_mask – if the mask cell equals this value, overwrite data

  • new_value_in_data – the value to overwrite data with

void maskDataRows(std::vector<uint_t> row_select, Eigen::Ref<MatFloat> mask, float_t value_of_mask_to_mask, float_t new_value_in_data)#

fills data with new_value_in_data where mask==value_of_mask_to_mask

This function applies the same mask to all rows. To apply a different mask per row, see TransArray::maskData

Warning

maskDataRows and maskData may swap names

Parameters:
  • row_select – vec of rows to apply masks to

  • mask – mask for all rows

  • value_of_mask_to_mask – if the mask cell equals this value, overwrite data

  • new_value_in_data – the value to overwrite data with

void swapRowsValues(std::vector<uint_t> row_select, float_t value_to_mask, float_t new_value)#

swap value_to_mask with new_value in row_select rows

Warning

This may be renamed to swapValues

Parameters:
  • row_select – vec of rows to swap values on

  • value_to_mask – value to filter out (noData value?)

  • new_value – fill value

void hadamardProduct(Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#

element-wise product of two matrices

will fill this array with the result. Therefore, in1, in2 and this array’s dimensions MUST match exactly.

Parameters:

in1, in2 – The two matrices to multiply element-wise

void offsetAndScale(float_t offset, float_t scaling)#

Applies offset and scale to all rows.

v_new = (v_old+offset)*scaling

Parameters:
  • offset – Offset to apply

  • scaling – Scale factor

void scaleAndOffset(float_t offset, float_t scaling)#

Applies scale and offset to all rows.

v_new = v_old*scaling+offset

Parameters:
  • offset – Offset to add

  • scaling – Scale factor

void offsetsAndScales(std::vector<uint_t> row_select, Eigen::Ref<VecFloat> offsets, Eigen::Ref<VecFloat> scalings)#

Applies multiple offsets and scales.

v_new = (v_old+offset)*scaling

Warning

UNSTABLE Currently index of offsets and scalings is the same as in m_data. This may be updated so that it is the same as row_select.

Parameters:
  • row_select – permutation matrix with rows to select

  • offsets – Offsets to apply.

  • scalings – Scale factors

Data Processing#

group Data processing

data processing and statistics

Functions

void nanMean(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<VecFloat> out_data)#

see TransArray::nanMean

void blocksAverage(Eigen::Ref<MatFloat> out, const uint_t n_threads, Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2, uint_t n_pix, uint_t y)#

see TransArray::blocksAverage

void blocksAverageVecs(Eigen::Ref<MatFloat> out, const uint_t n_threads, Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2, uint_t n_pix, uint_t y, uint_t row_offset)#
void elementwiseAverage(Eigen::Ref<MatFloat> out, const uint_t n_threads, Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#
void extractIndicators(Eigen::Ref<MatFloat> data_in, const uint_t n_threads, Eigen::Ref<MatFloat> data_out, uint_t col_in_select, std::vector<uint_t> col_out_select, std::vector<uint_t> classes)#
void fitPercentage(Eigen::Ref<MatFloat> out, const uint_t n_threads, Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#
void texturesBwTransform(Eigen::Ref<MatFloat> texture_1, const uint_t n_threads, Eigen::Ref<MatFloat> texture_2, float_t k, float_t a, Eigen::Ref<MatFloat> sand, Eigen::Ref<MatFloat> silt, Eigen::Ref<MatFloat> clay)#
void linearRegression(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<VecFloat> x, Eigen::Ref<VecFloat> beta_0, Eigen::Ref<VecFloat> beta_1)#
void computeMannKendallPValues(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<VecFloat> out_data)#
void averageAggregate(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, uint_t agg_factor)#
void maskDifference(Eigen::Ref<MatFloat> data, const uint_t n_threads, float_t diff_th, uint_t count_th, Eigen::Ref<MatFloat> ref_data, Eigen::Ref<MatFloat> mask_out)#
void computeNormalizedDifference(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> positive_indices, std::vector<uint_t> negative_indices, std::vector<uint_t> result_indices, float_t positive_scaling, float_t negative_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeNirv(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> nir_indices, std::vector<uint_t> red_indices, std::vector<uint_t> result_indices, float_t nir_scaling, float_t red_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeBsi(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> swir1_indices, std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> blue_indices, std::vector<uint_t> result_indices, float_t swir1_scaling, float_t red_scaling, float_t nir_scaling, float_t blue_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeEvi(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> blue_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t blue_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeFapar(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeSavi(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeGeometricTemperature(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> latitude, Eigen::Ref<MatFloat> elevation, float_t elevation_scaling, float_t a, float_t b, float_t result_scaling, std::vector<uint_t> result_indices, std::vector<float_t> days_of_year)#
void computePercentiles(Eigen::Ref<MatFloat> data, const uint_t n_threads, std::vector<uint_t> col_in_select, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_out_select, std::vector<float_t> percentiles)#

see TransArray::computePercentiles

void fitProbabilities(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, float_t input_scaling, uint_t target_scaling, Eigen::Ref<MatFloat> best_classes_data, uint_t n_best_classes)#
void applyTsirf(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, uint_t out_index_offset, float_t w_0, Eigen::Ref<VecFloat> w_p, Eigen::Ref<VecFloat> w_f, bool keep_original_values, const std::string &version, const std::string &backend)#
void nanMeanAggregatePattern(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> &agg_pattern)#
void convolveRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, float_t w_0, Eigen::Ref<VecFloat> w_p, Eigen::Ref<VecFloat> w_f)#
void slidingWindowClassMode(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, size_t window_size)#
void nanMean(Eigen::Ref<VecFloat> out_data)#

Computes the per-row mean of non-Nan values.

Processes chunks in parallel

Parameters:

out_data – Vec of rows length

void elementwiseAverage(Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#

Fills this array with the element-wise average of in1 and in2

Parameters:
  • in1 – The first matrix

  • in2 – The second matrix

void computePercentiles(std::vector<uint_t> col_in_select, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_out_select, std::vector<float_t> percentiles)#

Computes per-pixel percentiles across layers.

This works on transposed data where layers are columns and rows are pixels

_images/compute_percentiles.svg

Parameters:
  • col_in_select – The columns to use for computing percentiles

  • out_data – Output array, should have n_pix rows

  • col_out_select – columns in Output array to fill

  • percentiles – List of percentiles [0-100]

IoArray Class#

This class is used by TiledDataLoader, TiledDataExporter and SpaceOverlay

class IoArray : public skmap::ParArray#

Public Functions

void readDataCore(Eigen::Ref<MatFloat::RowXpr> row, std::string file_loc, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, GDALDataType read_type, int band, std::optional<float_t> value_to_mask, std::optional<float_t> value_to_set)#

Reads a portion of a raster dataset into a row buffer.

This function opens a raster dataset using GDAL, reads a rectangular window defined by offsets and sizes into an Eigen row expression, and optionally replaces a specific masked value with another value.

Opens the GDAL dataset in read-only mode, reads the specified rectangular portion into the provided Eigen row buffer, and closes the dataset. If both value_to_mask and value_to_set are provided and are different, the function replaces all occurrences of value_to_mask in the buffer with value_to_set.

Note

If only value_to_set is provided, the function automatically determines the dataset’s nodata value to use as value_to_mask.

Parameters:
  • row – Reference to an Eigen row expression where the data will be stored.

  • file_loc – File path or URL of the raster dataset to read.

  • x_off – Horizontal offset of the window to read.

  • y_off – Vertical offset of the window to read.

  • x_size – Width of the window to read.

  • y_size – Height of the window to read.

  • read_type – GDAL data type to read (e.g., GDT_Float32).

  • bands_list – List of band indices to read from the dataset.

  • value_to_mask – Optional value in the dataset to treat as “masked” (e.g., nodata).

  • value_to_set – Optional value to replace the masked values with.

Throws:

skmapAssertIfTrue – if the dataset cannot be opened or read.

void setupGdal(dict_t dict)#

Initialize GDAL with custom configuration options and error handling.

Sets GDAL runtime options from the input dictionary, registers all GDAL drivers, and configures error logging to suppress console output if possible.

Each key-value pair in dict is applied via CPLSetConfigOption(). GDALAllRegister() is called to ensure all drivers are available. Error logging is redirected to “/dev/null” if possible; otherwise, a quiet error handler is used to suppress warnings and errors.

Note

This function modifies global GDAL state and should be called before performing any GDAL I/O operations. All Python-exposed functions call it already.

Parameters:

dict – Dictionary of GDAL configuration options (key-value pairs) where keys and values are strings recognized by GDAL.

TransArray Class#

class TransArray : public skmap::ParArray#

Public Functions

TransArray(Eigen::Ref<MatFloat> data, const uint_t n_threads)#
void linearRegression(Eigen::Ref<VecFloat> x, Eigen::Ref<VecFloat> beta_0, Eigen::Ref<VecFloat> beta_1)#
void averageAggregate(Eigen::Ref<MatFloat> out_data, uint_t agg_factor)#
void maskDifference(float_t diff_th, uint_t count_th, Eigen::Ref<MatFloat> ref_data, Eigen::Ref<MatFloat> mask_out)#
void fitPercentage(Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2)#
void computeMannKendallPValues(Eigen::Ref<VecFloat> out_data)#
void blocksAverage(Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2, uint_t n_pix, uint_t y)#
void blocksAverageVecs(Eigen::Ref<MatFloat> in1, Eigen::Ref<MatFloat> in2, uint_t n_pix, uint_t y, uint_t row_offset)#
void texturesBwTransform(Eigen::Ref<MatFloat> texture_2, float_t k, float_t a, Eigen::Ref<MatFloat> sand, Eigen::Ref<MatFloat> silt, Eigen::Ref<MatFloat> clay)#
void slidingWindowClassMode(Eigen::Ref<MatFloat> out_data, size_t window_size)#
void computeNormalizedDifference(std::vector<uint_t> positive_indices, std::vector<uint_t> negative_indices, std::vector<uint_t> result_indices, float_t positive_scaling, float_t negative_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeBsi(std::vector<uint_t> swir1_indices, std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> blue_indices, std::vector<uint_t> result_indices, float_t swir1_scaling, float_t red_scaling, float_t nir_scaling, float_t blue_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeEvi(std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> blue_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t blue_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeNirv(std::vector<uint_t> nir_indices, std::vector<uint_t> red_indices, std::vector<uint_t> result_indices, float_t nir_scaling, float_t red_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeFapar(std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeSavi(std::vector<uint_t> red_indices, std::vector<uint_t> nir_indices, std::vector<uint_t> result_indices, float_t red_scaling, float_t nir_scaling, float_t result_scaling, float_t result_offset, std::vector<float_t> clip_value)#
void computeGeometricTemperature(Eigen::Ref<MatFloat> latitude, Eigen::Ref<MatFloat> elevation, float_t elevation_scaling, float_t a, float_t b, float_t result_scaling, std::vector<uint_t> result_indices, std::vector<float_t> days_of_year)#
void fitProbabilities(Eigen::Ref<MatFloat> out_data, float_t input_scaling, uint_t target_scaling, Eigen::Ref<MatFloat> best_classes_data, uint_t n_best_classes)#
void applyTsirf(Eigen::Ref<MatFloat> out_data, uint_t out_index_offset, float_t w_0, Eigen::Ref<VecFloat> w_p, Eigen::Ref<VecFloat> w_f, bool keep_original_values, const std::string &version, const std::string &backend)#
void convolveRows(Eigen::Ref<MatFloat> out_data, float_t w_0, Eigen::Ref<VecFloat> w_p, Eigen::Ref<VecFloat> w_f)#
void extractIndicators(Eigen::Ref<MatFloat> data_out, uint_t col_in_select, std::vector<uint_t> col_out_select, std::vector<uint_t> classes)#
void nanMeanAggregatePattern(Eigen::Ref<MatFloat> out_data, std::vector<std::vector<uint_t>> &agg_pattern)#

ParArray Class#

class ParArray#

Subclassed by skmap::IoArray, skmap::TransArray

Public Functions

ParArray(Eigen::Ref<MatFloat> data, const uint_t n_threads)#
void printData()#
template<typename F>
inline void parForRange(F f, uint_t n_max)#
template<typename F>
inline void parRowPerm(F f_in, std::vector<uint_t> perm_vec)#

process f_in parallelly on rows, based on the order in perm_vec

perm_vec is allowed to be any length less than rows().

Warning

If an index appears twice, this will cause a data race

Parameters:
  • f_in – The function to apply to each row

  • perm_vec – vector with row indices

template<typename F>
inline void parChunk(F f_in)#

Process row chunks in parallel by f_in

Will divide the matrix into approximately equal blocks, where the remainder is shared over the first number of chunks

Parameters:

f_in – The function to apply to each chunk

Bindings#

Functions

GDALDataType GetGDALDataTypeFromString(const std::string &gdal_data_type_str)#
std::optional<NodataVariant> getNodataVariant(const std::string &gdal_data_type_str, float_t no_data_value)#
dict_t convPyDict(py::dict in_dict)#
map_t convPyMap(py::dict in_map)#
void readDataBlocks(Eigen::Ref<MatFloat> data, const uint_t n_threads, const std::vector<std::string> &file_locs, const std::vector<uint_t> perm_vec, const std::vector<uint_t> x_off_vec, const std::vector<uint_t> y_off_vec, const std::vector<uint_t> x_size_vec, const std::vector<uint_t> y_size_vec, const std::vector<int> bands_list, py::dict conf_GDAL, std::optional<std::vector<float_t>> value_to_mask_vec, std::optional<float_t> value_to_set)#
void readDataCore(Eigen::Ref<MatFloat> data, const uint_t n_threads, const std::string file_loc, const uint_t x_off, const uint_t y_off, const uint_t x_size, const uint_t y_size, int band, py::dict conf_GDAL, std::optional<float_t> value_to_mask, std::optional<float_t> value_to_set)#

core implementation, see IoArray::readDataCore

void writeInt16Data(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::vector<std::string> base_files, std::string base_folder, std::vector<std::string> file_names, std::vector<uint_t> data_indices, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, int16_t no_data_value, std::optional<std::string> bash_compression_command, std::optional<std::vector<std::string>> seaweed_path)#

Deprecated:

IoArray::writeData

void writeUInt16Data(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::vector<std::string> base_files, std::string base_folder, std::vector<std::string> file_names, std::vector<uint_t> data_indices, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, uint16_t no_data_value, std::optional<std::string> bash_compression_command, std::optional<std::vector<std::string>> seaweed_path)#

Deprecated:

use IoArray::writeData

void writeByteData(Eigen::Ref<MatFloat> data, const uint_t n_threads, py::dict conf_GDAL, std::vector<std::string> base_files, std::string base_folder, std::vector<std::string> file_names, std::vector<uint_t> data_indices, uint_t x_off, uint_t y_off, uint_t x_size, uint_t y_size, byte_t no_data_value, std::optional<std::string> bash_compression_command, std::optional<std::vector<std::string>> seaweed_path)#

Deprecated:

use IoArray::writeData

void extractArrayRows(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> row_select)#

Deprecated:

use selArrayRows in stead

void extractArrayCols(Eigen::Ref<MatFloat> data, const uint_t n_threads, Eigen::Ref<MatFloat> out_data, std::vector<uint_t> col_select)#

Deprecated:

use selArrayCols in stead

void checkSimdInstructionSetsInUse()#
PYBIND11_MODULE(skmap_bindings, m)#