pyewt.ewt2dvoronoi module

pyewt.ewt2dvoronoi.ewt2d_Voronoi_Filterbank(sizeImg, vorocells, tau, extH, extW)[source]

Builds the filter bank based on the detected Voronoi cells.

Parameters:
  • sizeImg (tuple) – Size of the image (height, width).

  • vorocells (list of ndarray) – List of Voronoi cell masks (binary arrays).

  • tau (float) – Transition width.

  • extH (bool or int) – 1 (or True) if there is vertical extension, 0 otherwise.

  • extW (bool or int) – 1 (or True) if there is horizontal extension, 0 otherwise.

Returns:

  • mfb (list of ndarray) – List of filters (same size as sizeImg).

  • Author (Jerome Gilles)

  • Institution (San Diego State University)

  • Version (1.0 (07/16/2025))

pyewt.ewt2dvoronoi.ewt2d_Voronoi_LP_function(x, tau)[source]

Calculate the Littlewood-Paley value based on the provided distance.

Parameters:
  • x (float) – Distance from the edge.

  • tau (float) – Transition width.

Returns:

  • y (float) – Value of the filter.

  • Author (Jerome Gilles)

  • Institution (San Diego State University)

  • Version (1.0 (07/16/2025))

pyewt.ewt2dvoronoi.ewt2d_Voronoi_Partition(centroids, sizeImg)[source]

Return the Voronoi partition corresponding to the provided centroids. The cells are properly assigned in order to guarantee a symmetric partition.

Parameters:
  • centroids (ndarray) – Array of shape (n_centroids, 2) containing the coordinates of each centroid.

  • sizeImg (tuple) – Size of the image (height, width).

Returns:

  • labelImage (ndarray) – Image containing the Voronoi partition tagged from 1 to the number of cells.

  • voronoi_cells (list of ndarray) – List of binary images containing the mask of each Voronoi cell.

  • Author (Jerome Gilles)

  • Institution (San Diego State University)

  • Version (1.0 (07/15/2025))

pyewt.ewt2dvoronoi.ewt2d_merge_symmetric(vor_cells, centers)[source]

Merge symmetric Voronoi cells based on their centers. :param vor_cells: List of Voronoi cell masks (binary arrays). :type vor_cells: list of ndarray :param centers: Array of shape (n_centers, 2) containing the coordinates of each center. :type centers: ndarray

Returns:

symvorcells – List of merged Voronoi cell masks.

Return type:

list of ndarray

pyewt.ewt2dvoronoi.ewt2d_voronoi(f, params)[source]

Compute the 2D EWT based on a Voronoi partitioning of the Fourier domain.

Parameters:
  • f (2D ndarray) – Input image.

  • params (dict) – Transform parameters.

Returns:

  • ewtc (list of 2D ndarray) – Collection of outputs of each EW filter.

  • mfb (list of 2D ndarray) – The built filter bank.

  • maxima (ndarray) – Coordinates of each meaningful detected maxima.

  • vorpartition (2D ndarray) – Detected Voronoi partition.

  • plane (ndarray) – Scale-space domain.

Notes

If params[“log”] is set to True, then the detection is performed on the logarithm of the spectrum of f. params[“typeDetect”] must be set to the wanted classification

method. The available options are: - “otsu” : uses Otsu’s technique - “halfnormal” : uses a half-normal law to model the problem - “empiricallaw” : uses the data itself to build a model of the problem - “kmeans” : uses kmeans to classify

params[“t”]: is the initial scale for the Gaussian kernel (0.8 is a good default value). params[“kn”]: is the kernel size for the Gaussian kernel (6 is a good default value). params[“niter”]: is the number of iterations through the scales (4 is a good default value). params[“edge”]: is the size (in pixels) of the strip to ignore at the edge of the image (0 is no strip). params[“includeCenter”]: if 1, the center of the image is included in the scale space maxima (0 is not included). params[“tau”]: is the half width of the transition area for the Voronoi partition (0.1 is a good default value). params[“complex”]: if 1, the Voronoi partition is complex, otherwise it is real (0 is real).

Though usually not providing better results, this function allows regularization of the spectrum of f using the selected method set in params[“reg”]. The available methods are: - “none” : does nothing, returns f - “gaussian” : convolve f with a Gaussian of length and standard

deviation given by params[“lengthfilter”] and params[“sigmafilter”], respectively

  • “average”: convolve f with a constant filter of length given

    by params[“lengthfilter”]

Author: Jerome Gilles Institution: San Diego State University Version: 1.0 (07/21/2025)

pyewt.ewt2dvoronoi.iewt2d_voronoi(ewtc, mfb)[source]

Performs the inverse Empirical Voronoi Wavelet Transform, returning a reconstruction of the image.

Parameters:
  • ewtc (list of 2D ndarray) – Empirical wavelet coefficients.

  • mfb (list of 2D ndarray) – Corresponding empirical wavelet filters.

Returns:

  • rec (2D ndarray) – Reconstructed image.

  • Author (Jerome Gilles)

  • Institution (San Diego State University)

  • Version (1.0 (07/16/2025))

pyewt.ewt2dvoronoi.plot_voronoi_filterbank(mfb)[source]

Plot 2D empirical Voronoi filters

Parameters:

  • mfb: list

    list containing each Voronoi filter

Author: Jerome Gilles Institution: San Diego State University Version: 1.0 (07/16/2025)

pyewt.ewt2dvoronoi.show_ewt2d_voronoi_boundaries(f, vor, color=None, logspec=0)[source]

Plots the edges of the Voronoi partition onto the magnitude of the Fourier spectrum of the input image.

Parameters:
  • f – ndarray Input image.

  • vor – ndarray Voronoi partition image (same size as f).

  • color – list or tuple RGB color for the partition edges, values in [0,1]. Default is red.

  • logspec – int If 1, plot the logarithm of the spectrum. Default is 0.

Author: Jerome Gilles Institution: San Diego State University Version: 1.0 (07/16/2025)