wolfhece.wolf_array._mb_model

WolfArrayMBModel - Data-only multi-block array class (no GUI/OpenGL dependency).

For GUI support (OpenGL rendering, wx dialogs), use WolfArrayMB from _mb.py.

Extracted from _mb.py to follow the same Model/GUI separation pattern as WolfArrayModel / WolfArray.

Author: HECE - University of Liege, Pierre Archambeau Date: 2024

Copyright (c) 2024 University of Liege. All rights reserved.

This script and its content are protected by copyright law. Unauthorized copying or distribution of this file, via any medium, is strictly prohibited.

Module Contents

class wolfhece.wolf_array._mb_model.WolfArrayMBModel(fname=None, mold=None, masknull=True, crop=None, whichtype=WOLF_ARRAY_MB_SINGLE, preload=True, nullvalue=0.0, srcheader=None, idx: str = '')[source]

Bases: wolfhece.wolf_array._base.WolfArrayModel

Inheritance diagram of wolfhece.wolf_array._mb_model.WolfArrayMBModel

Data-only multi-block array class (no GUI/OpenGL dependency).

A multi-block array is composed of several independent rectangular blocks (each a WolfArrayModel) stored in the myblocks dictionary, keyed by block identifiers produced by getkeyblock().

The global header (header_wolf) describes the bounding box that encloses all blocks. Individual block headers are kept in head_blocks (inherited from header_wolf).

For GUI support (OpenGL rendering, wx dialogs, palette management), use WolfArrayMB from _mb.py.

Variables:
  • myblocks – Dictionary mapping block keys to WolfArrayModel instances. Keys are produced by getkeyblock().

  • _active_blocks – Number of currently active blocks (0 = all).

  • mngselection – Selection manager (SelectionDataMB).

myblocks: dict[str, wolfhece.wolf_array._base.WolfArrayModel][source]
mngselection[source]
_active_blocks = 0[source]
classmethod make_multiblocks(others: list, abs: bool = True, copy: bool = True) WolfArrayMBModel[source]

Merge several single-block arrays into a new multi-block array.

The global header is recomputed from the bounding box of all blocks.

Parameters:
  • others – List of WolfArrayModel instances to merge.

  • abs – If True, use global (absolute) world coordinates.

  • copy – If True, copy data from each source array. Set to False only if the source headers can be modified in place (saves memory).

Returns:

A new WolfArrayMBModel containing all blocks.

property zmin: float[source]

Minimum unmasked value across all blocks.

property zmax: float[source]

Maximum unmasked value across all blocks.

property zmin_global: float[source]

Minimum value across all blocks (including masked cells).

property zmax_global: float[source]

Maximum value across all blocks (including masked cells).

property nullvalue: float[source]

No-data / null value shared by all blocks.

add_block(arr: wolfhece.wolf_array._base.WolfArrayModel, force_idx: bool = False, copyarray=False)[source]

Add a block to this multi-block array.

After adding all blocks, call set_header_from_added_blocks() to recompute the global header.

Parameters:
  • arr – The block to add.

  • force_idx – If True, the block key is set automatically based on the current number of blocks. If False, the key must already be set on arr.

  • copyarray – If True, a deep copy of arr is stored instead of the original reference.

Raises:

AssertionError – If force_idx is False and the key on arr is missing, duplicated, or out of sequence.

count()[source]

Count the number of unmasked (non-null) cells across all blocks.

The result is stored in nbnotnull (total) and in each block’s nbnotnull attribute.

check_consistency(other)[source]

Check that other has the same block structure as self.

Two multi-block arrays are consistent if they have the same number of blocks and each pair of corresponding blocks has identical spatial headers.

Parameters:

other – Another WolfArrayMBModel to compare.

Returns:

True if the structures match.

Return type:

bool

allocate_ressources()[source]

Allocate memory for every block defined in head_blocks.

Creates WolfArrayModel instances for each block header if myblocks is empty. Does nothing if blocks are already populated.

set_header_from_added_blocks()[source]

Recompute the global header from the bounding box of all blocks.

Must be called after the last add_block() call. Updates origx, origy, dx, dy and the per-block sub-headers in head_blocks.

reset()[source]

Reset every block’s data to its null value.

read_data()[source]

Read binary block data from filename.

One WolfArrayModel is created per block header and its binary payload is read sequentially from the file.

write_array()[source]

Write the raw binary data of every block to filename.

Blocks are written sequentially in index order.

check_bounds_ij(i: int, j: int, which_block: int = 1)[source]

Check whether indices (i, j) fall inside the global bounds.

Converts indices to world coordinates via get_xy_from_ij() and delegates to check_bounds_xy().

Parameters:
  • i – Column index (0-based).

  • j – Row index (0-based).

  • which_block – 1-based block index used for the coordinate conversion.

Returns:

True if the resulting coordinates are inside the global bounding box.

Return type:

bool

check_bounds_xy(x: float, y: float)[source]

Check whether world coordinates (x, y) are inside the global bounding box.

Parameters:
  • x – X world coordinate.

  • y – Y world coordinate.

Returns:

True if inside bounds.

Return type:

bool

get_ij_from_xy(x: float, y: float, z: float = 0.0, scale: float = 1.0, aswolf: bool = False, abs: bool = True, which_block: int = 1)[source]

Return cell indices (i, j) for the given world coordinates.

Delegates to the single-block WolfArrayModel.get_ij_from_xy() of block which_block.

Parameters:
  • x – X world coordinate.

  • y – Y world coordinate.

  • z – Z coordinate (unused, kept for compatibility).

  • scale – Scale factor applied before conversion.

  • aswolf – If True, return 1-based (Fortran) indices.

  • abs – If True, account for translation.

  • which_block – 1-based block index.

Returns:

Tuple (i, j).

get_xy_from_ij(i: int, j: int, which_block: int = 1, aswolf: bool = False, abs: bool = True)[source]

Return world coordinates (x, y) for cell indices (i, j).

Delegates to the single-block WolfArrayModel.get_xy_from_ij() of block which_block.

Parameters:
  • i – Column index.

  • j – Row index.

  • which_block – 1-based block index.

  • aswolf – If True, indices are 1-based (Fortran).

  • abs – If True, account for translation.

Returns:

Tuple (x, y).

get_value(x: float, y: float, abs: bool = True, nullvalue: float = -99999, convert_to_float: bool = True)[source]

Read the value at world coordinate (x, y).

Iterates through all blocks and returns the first unmasked value found. If no block covers the coordinate, nullvalue is returned.

Parameters:
  • x – X world coordinate.

  • y – Y world coordinate.

  • abs – If True, account for translation.

  • nullvalue – Value returned when no block covers (x, y).

  • convert_to_float – If True, cast the result to Python float before returning.

Returns:

Cell value, or nullvalue if not found.

get_values_as_wolf(i: int, j: int, which_block: int = 1)[source]

Return the value at 1-based Fortran indices (i, j).

Parameters:
  • i – 1-based column index.

  • j – 1-based row index.

  • which_block – 1-based block index.

Returns:

Cell value or np.nan if out of bounds.

get_blockij_from_xy(x: float, y: float, abs: bool = True)[source]

Find which block covers the world coordinate (x, y).

Parameters:
  • x – X world coordinate.

  • y – Y world coordinate.

  • abs – If True, account for translation.

Returns:

Tuple (i, j, k) where i and j are the 0-based cell indices within block k (1-based block index). Returns (-1, -1, -1) if no block covers the coordinate.

_make_result_array()[source]

Create an empty result array of the same type.

Subclasses override to return their own type (e.g. WolfArrayMB).

Returns:

A new WolfArrayMBModel with the same header.

_make_block_op(op, curblock, other_block)[source]

Apply a binary operator between two blocks (or a block and a scalar).

Parameters:
  • op – Operator string: '+', '-', '*', '/', or '**'.

  • curblock – Left operand (WolfArrayModel).

  • other_block – Right operand (WolfArrayModel or scalar).

Returns:

A new WolfArrayModel with the result.

_apply_op(other, op: str)[source]

Generic element-wise arithmetic on multi-block arrays.

Handles both WolfArrayMBModel + WolfArrayMBModel and WolfArrayMBModel + scalar cases.

Parameters:
  • other – Another WolfArrayMBModel (same structure) or a float scalar.

  • op – Operator string ('+', '-', etc.).

Returns:

A new multi-block array with the result, or None on inconsistency.

mask_data(value)[source]

Mask cells whose data equals value in every block.

For integer array types the value is cast to int first. Handles NaN comparisons correctly.

Parameters:

value – The value to mask. None is a no-op.

mask_reset()[source]

Reset the mask of every block to False (all cells unmasked).

mask_lower(value)[source]

Mask cells whose value is strictly less than value.

Parameters:

value – Threshold value.

mask_lowerequal(value)[source]

Mask cells whose value is less than or equal to value.

Parameters:

value – Threshold value.

mask_greater(value)[source]

Mask cells whose value is strictly greater than value.

Parameters:

value – Threshold value.

mask_greaterequal(value)[source]

Mask cells whose value is greater than or equal to value.

Parameters:

value – Threshold value.

mask_union(source: WolfArrayMBModel)[source]

Compute the mask union (logical AND) with another multi-block array.

The intersection of masked regions is applied block by block.

Parameters:

source – Another WolfArrayMBModel with the same block structure.

copy_mask(source: WolfArrayMBModel, forcenullvalue: bool = False)[source]

Copy the mask from source to self, block by block.

Parameters:
  • source – Source WolfArrayMBModel whose mask is copied.

  • forcenullvalue – If True, also set masked cells to the null value.

interpolate_on_polygon(working_vector: wolfhece.PyVertexvectors.vector, method: Literal['nearest', 'linear', 'cubic'] = 'linear', keep: Literal['all', 'below', 'above'] = 'all')[source]

Interpolate elevation values under a polygon for every block.

For each block, cells inside the polygon (or the current selection if one exists) are filled with values interpolated from the 3-D vertices of working_vector using scipy.interpolate.griddata().

Parameters:
  • working_vector – Polygon whose 3-D vertices provide the interpolation source.

  • method – Interpolation method ('nearest', 'linear', or 'cubic').

  • keep – Which interpolated values to keep: 'all', 'below' (only values lower than existing), or 'above'.

interpolate_on_polygons(working_zone: wolfhece.PyVertexvectors.zone, method: Literal['nearest', 'linear', 'cubic'] = 'linear', keep: Literal['all', 'below', 'above'] = 'all')[source]

Interpolate elevation values under every polygon in a zone.

Calls interpolate_on_polygon() for each vector in working_zone.

Parameters:
  • working_zone – Zone containing one or more polygon vectors.

  • method – Interpolation method.

  • keep – Which interpolated values to keep.

interpolate_on_polyline(working_vector: wolfhece.PyVertexvectors.vector, usemask=True)[source]

Interpolate elevation values under a polyline for every block.

Cells under the polyline are filled with values interpolated along the 3-D vertices of working_vector using Shapely’s line.interpolate().

Parameters:
  • working_vector – Polyline whose 3-D vertices provide the interpolation source.

  • usemask – If True, only consider unmasked cells.

interpolate_on_polylines(working_zone: wolfhece.PyVertexvectors.zone, usemask=True)[source]

Interpolate elevation values under every polyline in a zone.

Calls interpolate_on_polyline() for each vector in working_zone.

Parameters:
  • working_zone – Zone containing one or more polyline vectors.

  • usemask – If True, only consider unmasked cells.

filter_zone(set_null: bool = False)[source]

Filter zones, keeping only those that contain selected cells.

Applied independently to every block.

Parameters:

set_null – If True, set filtered-out cells to the null value instead of just masking them.

labelling()[source]

Label connected components in every block using SciPy.

Each block is labelled independently via scipy.ndimage.label().

as_WolfArray(abs: bool = True, forced_header: wolfhece.wolf_array._header_wolf.header_wolf = None) wolfhece.wolf_array._base.WolfArrayModel[source]

Convert to a single-block WolfArrayModel.

When blocks have different resolutions, the finer blocks are preserved and coarser blocks are rebinned to match the finest resolution.

Parameters:
  • abs – If True, use absolute coordinates (translation included).

  • forced_header – If provided, override the target header instead of computing it from blocks.

Returns:

A new single-block WolfArrayModel.