API Reference

pySVA: A package for Salinity Variance Analysis in ocean modeling.

This package provides tools for analyzing salinity variance in hydrodynamic models, particularly for D-Flow FM simulations.

Core module for pySVA containing the main constructorSVA class.

This module defines the constructorSVA class which handles data loading, grid operations, and provides access to various hydrodynamic properties for salinity variance analysis.

class pySVA.sva_core.constructorSVA(input_file, data_description, **kwargs)[source]

Bases: object

Core class for Salinity Variance Analysis (SVA) on unstructured grids.

This class loads and processes hydrodynamic data from D-Flow FM simulations and provides access to grid topology, geometric properties, velocity fields, tracers, and derived quantities required for variance budget analysis.

The class supports multiple input formats (netCDF files, xarray.Dataset, xugrid.UgridDataset, or xarray.DataArray) and automatically extracts grid structure and coordinate systems following UGRID conventions.

Parameters:
  • input_file (str or xr.Dataset or xu.UgridDataset or xr.DataArray) – Input data source. Can be a file path or pre-loaded dataset.

  • data_description (dict) –

    Mapping between logical variable names and dataset variables.

    Required keys: - velx, vely, velz - viscosity - tracer - depth

    Optional keys: - horizontal_diffusivity - interfaces - bed_level - volume - flow_area

ds

Underlying hydrodynamic dataset.

Type:

xr.Dataset or xu.UgridDataset

grid

UGRID mesh object.

Type:

object

face_coords, edge_coords, node_coords

Cartesian coordinates of mesh elements.

Type:

xr.DataArray

dimn_*

Standardized dimension names for grid topology.

Type:

str

Notes

  • All geometric and topological quantities follow UGRID conventions.

  • Coordinates are expressed in projected Cartesian space (typically meters).

  • Many properties are computed lazily using functools.cached_property().

advection(integration='depth', depth_averaged=False)[source]

Compute the advection term in the tracer variance budget, computed as:

\[\text{advection} = \nabla_h \cdot (\mathbf{u} (c')^2)\]

where: - \(c' = c - \bar{c}\) is the tracer perturbation relative to the depth-averaged tracer; - \(\nabla_h\) indicates the horizontal divergence; - \(\mathbf{u}=(u,v)\) represents the full velocity vector.

Integration over depth or volume is optional.

Parameters:
  • integration ({"depth", "volume", "none"}, optional) – Integration method applied to the computed advection field. Default is “depth”.

  • depth_averaged (bool, optional) – If True, return depth-averaged advection.

Returns:

Advection term in the tracer variance budget. Units are \(c^2\, s^{-1}\).

Return type:

xarray.DataArray

property bed_level[source]

Compute bed level at each grid face.

The bed level is derived from the minimum interface elevation, optionally averaged over time.

Returns:

Bed elevation with dimensions (n_faces, ...).

Return type:

xr.DataArray

Raises:

ValueError – If interface data is not available.

Notes

  • Computed as:

\[z_b = \min(\text{interfaces})\]
  • If time-dependent, the result is averaged over time.

property cell_area[source]

Compute horizontal area of each grid cell.

The horizontal cell area is derived from edge geometry using edge lengths and outward unit normal vectors, consistent with finite-volume formulations on unstructured grids.

Returns:

Horizontal cell area defined on faces with dimension (n_faces,) or including additional dimensions if time-dependent.

Return type:

xr.DataArray

Notes

  • The area is computed by integrating contributions from all

edges surrounding a face.

  • In discrete form, this corresponds to a polygonal area

reconstruction using edge vectors and normals:

\[A = \frac{1}{2} \sum_e \mathbf{r}_e \cdot \mathbf{n}_e\]

where \(\mathbf{r}_e\) is the edge vector and \(\mathbf{n}_e\) is the outward unit normal.

  • The formulation is consistent with finite-volume flux

integration and ensures geometric conservation.

  • Assumes a planar (Cartesian) coordinate system.

property cell_thickness[source]

Compute vertical thickness of grid cells.

The cell thickness is defined as the vertical distance between consecutive interfaces in the water column.

Returns:

Cell thickness with dimensions including the vertical layer coordinate (typically n_layers or equivalent) and the horizontal grid dimension (e.g., faces or edges).

Return type:

xr.DataArray

Notes

  • The thickness is computed as:

\[\Delta z_k = z_{k+1} - z_k\]

where \(z_k\) and \(z_{k+1}\) are consecutive interface elevations.

  • The result represents layer thickness in physical space

(typically meters).

  • The sign convention follows the vertical coordinate system of

the input dataset (usually positive upward for interfaces, resulting in positive layer thickness).

compute_gradient_on_face(uda, scalar=False, add_to_dataset=False, mode_depth=None)[source]

Compute the horizontal gradient of a variable on grid faces.

The gradient is evaluated using Gauss’ divergence theorem by summing flux contributions across edges surrounding each face.

Parameters:
  • uda (xu.UgridDataArray) – Variable for which the gradient is computed.

  • scalar (bool, optional) – If True, treat the variable as scalar and use face areas instead of flow areas. Default is False.

  • add_to_dataset (bool, optional) – If True, store the resulting gradient in the dataset.

  • mode_depth ({"sum", "mean", None}, optional) – Controls how vertical layers are treated when computing the gradient.

Returns:

Gradient vector with dimensions (faces, cartesian_dim).

Return type:

xarray.DataArray

property depth_integrated_tracer_variance[source]

Compute the depth-integrated tracer variance.

The integration is performed using the trapezoidal rule.

Returns:

Depth-integrated tracer variance.

Return type:

xarray.DataArray

property distance_vectors[source]

Compute vectors from face centroids to edge midpoints.

These vectors are used in gradient and geometric calculations involving the unstructured mesh.

Returns:

Distance vectors with dimensions (faces, max_face_edges, cartesian_dim).

Return type:

xarray.DataArray

property edge_face_weights[source]

Compute interpolation weights from edges to faces.

For each edge, weights are computed for all connected faces using distance-based weighting.

Returns:

Weights with dimensions (n_edges, nMax_edge_faces).

Return type:

xr.DataArray

Notes

  • Based on distances between edge midpoints and face centroids.

  • Weights are normalized per edge.

  • Returned as a Dask-backed array for scalability.

property edge_faces[source]

Return the edge-to-face connectivity matrix.

This property exposes the UGRID edge-face connectivity stored in the dataset as an xarray.DataArray.

Returns:

Connectivity array with dimensions (edges, max_edge_faces) mapping each edge to its adjacent faces.

Return type:

xarray.DataArray

Notes

The array contains integer indices referencing faces in the grid. Missing neighbors are indicated using the dataset’s fill value.

property edge_length[source]

Compute geometric length of each edge.

The edge length is calculated as the Euclidean distance between the two nodes defining each edge in Cartesian space.

Returns:

Edge length with dimension (n_edges,).

Return type:

xr.DataArray

Notes

  • Computed using the Pythagorean distance:

\[L = \sqrt{(x_2 - x_1)^2 + (y_2 - y_1)^2}\]
property edge_node_coords[source]

Retrieve the coordinates of nodes associated with each edge.

Maps the node coordinates to the corresponding edge-node connectivity, returning a DataArray with NaNs where no node is connected (fill value).

Returns:

Node coordinates per edge with dimensions (n_edges, nMax_edge_nodes, nCartesian_coords).

Return type:

xr.DataArray

Notes

  • Derived from edge–node connectivity.

  • Missing nodes (fill values) are replaced with NaN.

  • Used for geometric calculations such as edge length and normals.

property edge_nodes[source]

Return the edge-to-node connectivity matrix.

Constructs an xarray.DataArray describing which nodes define each edge of the unstructured grid.

Returns:

Edge-node connectivity array with dimensions (edges, max_edge_nodes).

Return type:

xarray.DataArray

Notes

The returned array contains CF-UGRID compliant metadata including:

  • cf_role = "edge_node_connectivity"

  • start_index = 0

Missing entries are filled using the dataset’s configured fill value.

property face_area[source]

Compute face-associated cross-sectional area on edges.

The face area is defined as the product of interpolated cell thickness and edge length.

Returns:

Face area with dimensions matching edges (typically n_edges and vertical layers if present).

Return type:

xr.DataArray

Notes

  • Cell thickness is interpolated from faces to edges.

  • Edge length is computed geometrically from node coordinates.

  • The result is stored as {gridname}_face_area.

property face_edge_weights[source]

Compute interpolation weights from faces to edges.

For each face, weights are computed for all connected edges using distance-based weighting.

Returns:

Weights with dimensions (n_faces, nMax_face_edges).

Return type:

xr.DataArray

Notes

  • Based on distances between face centroids and edge midpoints.

  • Weights are normalized such that they sum to 1 per face.

  • Typically used for interpolation from face-centered variables

to edge-centered quantities.

property face_edges[source]

Return the face-to-edge connectivity matrix.

Constructs an xarray.DataArray describing which edges bound each face of the unstructured grid.

Returns:

Connectivity array with dimensions (faces, max_face_edges).

Return type:

xarray.DataArray

Notes

The returned object follows CF-UGRID conventions with cf_role = "face_edge_connectivity".

get_all_coordinates()[source]

Retrieve Cartesian coordinates of all mesh elements.

Returns:

Tuple containing:

  • face_coords : (n_faces, nCartesian_coords)

  • edge_coords : (n_edges, nCartesian_coords)

  • node_coords : (n_nodes, nCartesian_coords)

Return type:

tuple of xr.DataArray

Notes

  • Coordinates are derived from UGRID metadata.

  • Returned arrays include CF-compliant attributes.

  • Cartesian dimension is defined as {gridname}_nCartesian_coords.

horizontal_dissipation(integration='depth', depth_averaged=False)[source]

Compute the horizontal dissipation term in the tracer variance budget.

Horizontal dissipation is calculated using the horizontal diffusivity and squared horizontal gradients of the tracer:

\[\text{horizontal dissipation} = 2 D_h \left[\left(\frac{\partial c}{\partial x}\right)^2 + \left(\frac{\partial c}{\partial y}\right)^2\right]\]

where: - \(c\) is the tracer; - \(D_h\) is the horizontal diffusivity.

Parameters:
  • integration ({"depth", "volume", "none"}, optional) – Integration method applied to the dissipation field. Default is “depth”.

  • depth_averaged (bool, optional) – If True, return depth-averaged horizontal dissipation.

Returns:

Horizontal dissipation term in the tracer variance budget. Units are \(c^2\, s^{-1}\). Returns None if calculation cannot be performed (e.g., horizontal diffusivity undefined).

Return type:

xarray.DataArray or None

property kzz[source]

Compute the vertical turbulent eddy diffusivity.

The vertical diffusivity \(K_{zz}\) is calculated as the sum of molecular diffusivity, background diffusivity, and a tracer-dependent laminar contribution.

\[K_{zz} = \frac{\nu}{Pr} + D_{0} + k_l\]

where \(\nu\) is the kinematic viscosity, \(Pr\) is the Prandtl/Schmidt number, \(D_{0}\) is a constant background diffusivity, and \(k_l\) is a tracer-dependent laminar component.

Parameters:
  • dicoww (float, optional) – Background vertical diffusivity (\(m^2 s^{-1}\)). Default is 5e-5.

  • prandtl_schmidt (float, optional) – Prandtl/Schmidt number used to convert viscosity to diffusivity. Default is 0.7.

Returns:

Vertical turbulent eddy diffusivity defined on grid edges.

Return type:

xarray.DataArray

Notes

The computed diffusivity is added to the dataset as {gridname}_dicwws and contains CF-style metadata describing the grid location and units.

straining(integration='depth', depth_averaged=False)[source]

Compute the straining term from the tracer variance budget.

The straining term is defined as

\[\text{straining} = -2 \, c' \, \mathbf{u}' \cdot \nabla_h \bar{c}\]

where

  • \(c' = c - \bar{c}\) is the tracer perturbation relative to the depth-averaged tracer;

  • \(\mathbf{u}' = \mathbf{u} - \bar{\mathbf{u}}\) is the velocity perturbation relative to the depth-averaged velocity;

  • \(\bar{c}\) indicates a depth-averaged quantity;

  • \(\nabla_h\) denotes the horizontal gradient operator.

Parameters:
  • integration ({"depth", "volume", "none"}, optional) – If “depth”, integrate the straining term over the vertical layers. If “volume”, integrate over the total volume (not yet implemented). If “none”, no integration is performed.

  • depth_averaged (bool, optional) – If True, return the term already averaged over depth layers. Otherwise, return the full 3D field.

Returns:

Straining term from the tracer variance budget, with units consistent with \(c^2\,s^{-1}\) (tracer squared per second).

Return type:

xarray.DataArray

tendency(integration='depth', depth_averaged=False)[source]

Compute the temporal tendency term of the tracer variance budget.

The tendency term is calculated as the time derivative of the tracer variance:

\[\frac{\partial (c')^2}{\partial t}\]

where \(c' = c - \bar{c}\) is the tracer perturbation relative to the depth-averaged tracer.

Optionally, this term can be integrated over depth or volume. A depth-averaged output can also be requested.

Parameters:
  • integration ({"depth", "volume", "none"}, optional) – Integration method for the tracer variance prior to differentiation. - “depth”: integrate over the vertical coordinate - “volume”: integrate over cell volumes - “none”: no integration applied Default is “depth”.

  • depth_averaged (bool, optional) – If True, divide the integrated tendency by the local water depth to return a depth-averaged quantity.

Returns:

Tendency term in the tracer variance budget. Units are \(c^2\, s^{-1}\), where \(c\) is the tracer concentration.

Return type:

xarray.DataArray

property tracer_variance[source]

Compute tracer variance.

The tracer variance is defined as the squared deviation of the tracer from its depth-averaged value.

Returns:

Tracer variance with the same dimensions as the input tracer field (typically including vertical layers and horizontal grid dimensions).

Return type:

xr.DataArray

Notes

  • The variance is computed as:

\[c'^2 = (c - \overline{c})^2\]

where \(c\) is the tracer concentration and \(\overline{c}\) is the depth-averaged tracer.

  • The depth-average is taken over the vertical coordinate

(layer or interface dimension).

uda_to_edges(uda)[source]

Interpolate a face-centered variable to grid edges.

The interpolation uses inverse-distance weights derived from the adjacent faces connected to each edge.

Parameters:

uda (xu.UgridDataArray) – Face-centered variable defined on the grid.

Returns:

Edge-centered variable obtained through weighted interpolation.

Return type:

xu.UgridDataArray

Raises:

ValueError – If the provided variable does not contain the face dimension.

Notes

Boundary edges with only one valid neighboring face inherit the value of that face.

uda_to_faces(uda)[source]

Interpolate an edge-centered variable to grid faces.

The interpolation averages values from the edges that bound each face using predefined edge weights.

Parameters:

uda (xu.UgridDataArray) – Edge-centered variable.

Returns:

Face-centered variable.

Return type:

xu.UgridDataArray

Raises:

ValueError – If the provided variable does not contain the edge dimension.

property unvs[source]

Compute unit normal vectors for all grid edges.

The normal vectors are defined perpendicular to each edge and normalized to unit length. If the dataset uses geographic coordinates (WGS84), coordinates are first projected to the Dutch RD stereographic coordinate system before vector calculation.

Returns:

Unit normal vectors with dimensions (edges, cartesian_dim).

Return type:

xarray.DataArray

Notes

The vectors are oriented outward relative to adjacent cells. The resulting array is stored in the dataset under the variable {gridname}_unvs.

vertical_dissipation(integration='depth', depth_averaged=False)[source]

Compute the vertical dissipation term in the tracer variance budget.

Vertical dissipation is calculated using the vertical eddy diffusivity and squared vertical gradients of the tracer:

\[\text{vertical dissipation} = 2 K_{zz} \left(\frac{\partial c}{\partial z}\right)^2\]

where: - \(c\) is the tracer,; - \(K_{zz}\) is the vertical diffusivity.

Parameters:
  • integration ({"depth", "volume", "none"}, optional) – Integration method applied to the dissipation field. Default is “depth”.

  • depth_averaged (bool, optional) – If True, return depth-averaged vertical dissipation.

Returns:

Vertical dissipation term in the tracer variance budget. Units are \(c^2\, s^{-1}\).

Return type:

xarray.DataArray

property volume[source]

Compute the volume of each grid cell.

Returns:

Cell volume calculated as cell_area times the cell_thickness.

Return type:

xarray.DataArray

property water_depth[source]

Compute water depth at each grid face.

Water depth is defined as the difference between the maximum interface elevation (water surface) and the bed level.

Returns:

Water depth with dimensions (n_faces, ...).

Return type:

xr.DataArray

Notes

  • Positive downward convention is used.

  • Computed as:

\[H = \eta - z_b\]

where \(\eta\) is the water surface elevation and \(z_b\) is the bed level.

Helper functions for pySVA.

This module contains utility functions for grid operations, coordinate calculations, and data transformations used in salinity variance analysis.

pySVA.sva_helpers.build_edge_face_weights(constructorSVA, **kwargs)[source]

Compute inverse-distance interpolation weights from faces to edges.

This function constructs weights for interpolating face-centered variables onto edges using inverse-distance weighting based on geometric distances between edge centers and neighboring face centroids.

Parameters:
  • constructorSVA (object) – Constructor object containing the UGRID dataset.

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    coordstuple of xr.DataArray, optional

    Tuple (face_coords, edge_coords, node_coords) containing precomputed coordinate arrays. If not provided, coordinates are derived from the dataset.

Returns:

Weights for edge-face interpolation with dimensions (n_edges, nMax_edge_faces).

Return type:

xr.DataArray

Notes

  • Uses edge-face connectivity to associate edges with neighboring faces.

  • Weights are computed using inverse-distance weighting.

  • Fill values in connectivity arrays are masked before computation.

pySVA.sva_helpers.build_edge_node_connectivity(constructorSVA)[source]

Construct edge–node connectivity array from a UGRID dataset.

This function extracts the edge–node connectivity from the underlying grid and returns it as a properly dimensioned xarray.DataArray with CF-compliant metadata.

Parameters:

constructorSVA (object) – Object containing the UGRID dataset (.ds attribute).

Returns:

Edge–node connectivity with dimensions (n_edges, nMax_edge_nodes).

Return type:

xr.DataArray

Raises:

IOError – If the dataset is not a xu.UgridDataset.

Notes

  • The connectivity defines which nodes belong to each edge.

  • Missing entries are filled with the grid _FillValue.

  • Attributes follow CF/UGRID conventions:

    • cf_role = "edge_node_connectivity"

    • start_index = 0

  • Coordinates associated with edges (e.g., edge centroids) are attached for convenience.

pySVA.sva_helpers.build_face_edge_connectivity(constructorSVA)[source]

Construct face–edge connectivity array from a UGRID dataset.

This function extracts the face–edge connectivity from the underlying grid and returns it as a properly dimensioned xarray.DataArray with CF-compliant metadata.

Parameters:

constructorSVA (object) – Object containing the UGRID dataset (.ds attribute).

Returns:

Face–edge connectivity with dimensions (n_faces, nMax_face_edges).

Return type:

xr.DataArray

Raises:

IOError – If the dataset is not a xu.UgridDataset.

Notes

  • The connectivity defines which edges bound each face.

  • The second dimension corresponds to the maximum number of edges per face (padded with _FillValue where necessary).

  • Attributes follow CF/UGRID conventions:

    • cf_role = "face_edge_connectivity"

    • start_index = 0

  • Face coordinates are attached as auxiliary coordinates.

pySVA.sva_helpers.build_inverse_distance_weights(a, b)[source]

Compute distance-based interpolation weights between coordinate sets.

This function constructs weights for mapping data between two sets of coordinates using pairwise distances. Typically used for interpolation between grid elements (e.g., faces → edges or nodes → edges).

Parameters:
  • a (xr.DataArray) – Reference coordinates with dimensions (N, nCartesian_coords).

  • b (xr.DataArray) – Target coordinates associated with each reference point, with dimensions (N, M, nCartesian_coords).

Returns:

Weights with dimensions (N, M), normalized such that the weights sum to 1 along the second dimension for each reference point.

Return type:

xr.DataArray

Notes

  • Distances are computed as:

    \[d_{ij} = \| \mathbf{a}_i - \mathbf{b}_{ij} \|\]
  • Weights are currently defined as:

    \[w_{ij} = \frac{d_{ij}}{\sum_j d_{ij}}\]

    Note that this corresponds to distance weighting, not true inverse distance weighting. For inverse-distance weighting, use:

    \[w_{ij} = \frac{1 / d_{ij}}{\sum_j (1 / d_{ij})}\]
  • Missing values (NaNs) are ignored in the normalization.

  • Implemented using xarray.apply_ufunc() for vectorized operation.

pySVA.sva_helpers.calculate_distance_haversine(lat, lon, lat_or, lon_or)[source]

Compute geodesic distance between coordinates using the WGS84 ellipsoid.

Parameters:
  • lat (float or pandas.Series) – Latitude(s) of target location(s) in degrees.

  • lon (float or pandas.Series) – Longitude(s) of target location(s) in degrees.

  • lat_or (float) – Reference latitude in degrees.

  • lon_or (float) – Reference longitude in degrees.

Returns:

Geodesic distance(s) in meters.

Return type:

float or ndarray

Notes

  • Uses pyproj.Geod.inv for accurate ellipsoidal distance computation.

  • Supports both scalar inputs and pandas Series.

pySVA.sva_helpers.calculate_distance_pythagoras(x1, y1, x2, y2)[source]

Compute Euclidean distance in a Cartesian coordinate system.

Parameters:
  • x1 (float or array-like) – Coordinates of the starting point.

  • y1 (float or array-like) – Coordinates of the starting point.

  • x2 (float or array-like) – Coordinates of the ending point.

  • y2 (float or array-like) – Coordinates of the ending point.

Returns:

Euclidean distance(s) in meters.

Return type:

float or ndarray

Notes

  • Assumes planar (projected) coordinates.

  • Equivalent to the L2 norm of the coordinate difference.

pySVA.sva_helpers.calculate_distance_vectors(constructorSVA, **kwargs)[source]

Compute distance vectors from face centroids to edge centers.

Parameters:
  • constructorSVA (object) – Constructor object containing the UGRID dataset.

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    face_edgesxr.DataArray, optional

    Precomputed face-edge connectivity.

Returns:

Distance vectors with dimensions (n_faces, nMax_face_nodes, nCartesian_coords).

Return type:

xr.DataArray

Notes

  • Distance vectors are defined as: edge coordinate minus face centroid coordinate.

  • Used for orientation of normal vectors and flux computations.

pySVA.sva_helpers.calculate_unit_normal_vectors(constructorSVA, **kwargs)[source]

Compute unit normal vectors for all edges in the unstructured grid.

The normal vectors are constructed perpendicular to each edge using node coordinates and normalized to unit length. The resulting vectors are oriented consistently with the edge definition.

Parameters:
  • constructorSVA (object) – Object containing the UGRID dataset (.ds attribute).

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    edge_node_coordsxr.DataArray, optional

    Precomputed coordinates of nodes belonging to each edge with dimensions (n_edges, 2, nCartesian_coords). If not provided, they are reconstructed from node connectivity.

Returns:

Unit normal vectors with dimensions (n_edges, nCartesian_coords). Returns None if the input dataset is not a UGRID dataset.

Return type:

xr.DataArray or None

Notes

  • The normal vector is defined as:

    \[\mathbf{n} = \frac{(-\Delta y, \Delta x)}{\|(-\Delta y, \Delta x)\|}\]

    where \(\Delta x\) and \(\Delta y\) are edge direction components.

  • The result is stored in the dataset as {gridname}_unvs.

  • If already present, the existing variable is reused.

  • Missing connectivity entries (fill values) are masked before computation.

pySVA.sva_helpers.compute_divergence_on_face(constructorSVA, uda, **kwargs)[source]

Compute the divergence of a vector field on mesh faces.

The divergence is calculated using a finite-volume formulation by summing fluxes across edges and normalizing by cell volume.

Parameters:
  • constructorSVA (object) – Constructor object containing the UGRID dataset.

  • uda (xr.DataArray or xu.UgridDataArray) – Vector field defined on edges. Must be aligned with edge-normal directions.

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    unvsxr.DataArray, optional

    Precomputed unit normal vectors.

Returns:

Divergence field on faces.

Return type:

xr.DataArray

Raises:

ValueError – If the input variable is not defined on edges.

Notes

  • Accounts for edge orientation using connectivity-based sign correction.

  • Result is stored in the dataset as {varname}_divergence.

Deprecated since version 1.0.0: This function will be removed in version 1.2.0. Use :meth: constructorSVA.compute_gradient_on_face` instead.

pySVA.sva_helpers.compute_gradient_on_face(constructorSVA, uda, **kwargs)[source]

Compute the gradient of a scalar field on mesh faces.

The gradient is computed using a finite-volume formulation by integrating fluxes across edges and normalizing by cell volume.

Parameters:
  • constructorSVA (object) – Constructor object containing the UGRID dataset.

  • uda (xr.DataArray or xu.UgridDataArray) – Input scalar field defined on edges.

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    unvsxr.DataArray, optional

    Precomputed unit normal vectors on edges. If not provided, they are calculated internally.

Returns:

Gradient field on faces with an additional Cartesian dimension.

Return type:

xr.DataArray

Raises:
  • ValueError – If the input variable is not defined on edges.

  • NameError – If required variables (flow area or volume) are missing.

Notes

  • Uses face-edge connectivity and unit normal vectors.

  • Ensures outward-pointing normals via dot-product sign correction.

  • Result is stored in the dataset as {varname}_gradient.

Deprecated since version 1.0.0: This function will be removed in version 1.2.0. Use :meth: constructorSVA.compute_gradient_on_face` instead.

pySVA.sva_helpers.compute_kzz(uds, dicoww=5e-05, prandtl_schmidt=0.7, tracer='mesh2d_sa1')[source]

Compute vertical turbulent diffusivity (kzz).

Parameters:
  • uds (xu.UgridDataset) – Input dataset containing viscosity field.

  • dicoww (float, optional) – Background diffusivity [m² s⁻¹].

  • prandtl_schmidt (float, optional) – Turbulent Prandtl/Schmidt number.

  • tracer (str, optional) – Tracer name used to determine molecular diffusivity.

Returns:

Dataset with added variable mesh2d_dicwwu.

Return type:

xu.UgridDataset

Notes

  • Combines turbulent and molecular diffusivity contributions.

  • Molecular diffusivity depends on tracer type (salinity or temperature).

  • Settings based on the standard settings in D-FLOW FM.

pySVA.sva_helpers.deprecated(reason, version, removal)[source]
pySVA.sva_helpers.depth_int2volume_int(uda: xarray.DataArray, cell_area: xarray.DataArray, dimn_faces: str) xarray.DataArray[source]

Convert depth-integrated values to volume-integrated values.

Parameters:
  • uda (xr.DataArray) – Depth-integrated variable.

  • cell_area (xr.DataArray) – Horizontal cell area.

  • dimn_faces (str) – Face dimension name.

Returns:

Volume-integrated variable.

Return type:

xr.DataArray

Notes

  • Performs horizontal integration after depth integration.

pySVA.sva_helpers.differentiate_over_3d_coord(uda: xarray.DataArray, coord_var: str, axis: int = -1) xarray.DataArray[source]

Differentiate a DataArray along a 3D coordinate.

Parameters:
  • uda (xr.DataArray) – Input variable.

  • coord_var (str) – Name of the coordinate variable.

  • axis (int, optional) – Axis corresponding to the coordinate dimension.

Returns:

Differentiated variable.

Return type:

xr.DataArray

Notes

  • Uses finite differences with Dask-compatible operations.

  • Pads result to preserve original array shape.

pySVA.sva_helpers.get_all_coordinates(uds)[source]

Extract face, edge, and node coordinates from a UGRID dataset.

This function constructs standardized coordinate arrays for faces, edges, and nodes, each expressed in Cartesian coordinates and returned as xarray.DataArray objects with consistent metadata.

Parameters:

uds (xu.UgridDataset or xu.UgridDataArray) – Input dataset or data array containing an unstructured grid following UGRID conventions.

Returns:

Tuple containing:

  • face_coordsxr.DataArray

    Face centroid coordinates with dimensions (n_faces, nCartesian_coords).

  • edge_coordsxr.DataArray

    Edge midpoint coordinates with dimensions (n_edges, nCartesian_coords).

  • node_coordsxr.DataArray

    Node coordinates with dimensions (n_nodes, nCartesian_coords).

Return type:

tuple of xr.DataArray

Raises:

IOError – If the input is not a xu.UgridDataset or xu.UgridDataArray.

Notes

  • Coordinate names are inferred from the UGRID attribute mesh2d.node_coordinates.

  • Coordinates are returned in projected Cartesian space (typically meters).

  • Metadata follows CF conventions, including:

    • standard_name = "projection_x_coordinate, projection_y_coordinate"

    • units = "m"

  • For xu.UgridDataset:

    • Face coordinates are reconstructed manually from dataset variables.

    • Edge and node coordinates are taken from the grid object.

  • For xu.UgridDataArray:

    • All coordinates are obtained directly from the grid object.

  • The Cartesian dimension is defined as:

    {gridname}_nCartesian_coords

  • These coordinate arrays are commonly used for:

    • geometric calculations (distances, normals)

    • interpolation between grid elements

    • finite-volume operators

pySVA.sva_helpers.integrate_trapz(y: xarray.DataArray, x: xarray.DataArray, dim: str | None = None) xarray.DataArray[source]

Integrate a DataArray using the trapezoidal rule.

Parameters:
  • y (xr.DataArray) – Values to integrate.

  • x (xr.DataArray) – Coordinate values along the integration dimension.

  • dim (str, optional) – Dimension along which to integrate.

Returns:

Integrated result.

Return type:

xr.DataArray

Notes

  • Assumes x is monotonically increasing.

  • Drops coordinate variables to avoid alignment conflicts.

pySVA.sva_helpers.reconstruct_vector_form(constructorSVA, vectors_list, **kwargs)[source]

Reconstruct a vector field from its Cartesian components.

This function concatenates scalar component fields (e.g., u and v) into a single vector-valued DataArray along the Cartesian dimension. The resulting vector is stored in the underlying dataset and the original component variables are removed.

Parameters:
  • constructorSVA (object) – Constructor object containing the underlying UGRID dataset (.ds).

  • vectors_list (list of str or xr.DataArray or xu.UgridDataArray) –

    List of vector components. Can be:

    • List of variable names (str) in the dataset

    • List of xarray.DataArray objects

    • List of xu.UgridDataArray objects

    All elements must be of the same type.

  • **kwargs (dict, optional) –

    Additional keyword arguments.

    vector_namestr, optional

    Name of the reconstructed vector variable. Defaults to {gridname}_uc.

Returns:

Reconstructed vector field with an added Cartesian dimension.

Return type:

xr.DataArray or xu.UgridDataArray

Raises:

UserWarning – If the type of vectors_list is not recognized.

Notes

  • The resulting vector is stored in constructorSVA.ds.

  • Original component variables are removed from the dataset.

  • The Cartesian dimension is defined as {gridname}_nCartesian_coords.

pySVA.sva_helpers.reconstruct_vector_form_magnitude(constructorSVA, varname, **kwargs)[source]

Reconstruct a vector field from edge-normal magnitude values.

This function converts a scalar magnitude defined along edge-normal directions (e.g., velocity components stored as normal fluxes) into full Cartesian vector components by multiplying with unit normal vectors.

Parameters:
  • constructorSVA (object) – Object containing the UGRID dataset (.ds attribute).

  • varname (str) – Name of the magnitude variable defined on edges (e.g., velocity in normal direction).

  • **kwargs (dict, optional) –

    Optional keyword arguments.

    face_edgesxr.DataArray, optional

    Face-edge connectivity array.

    edge_facesxr.DataArray, optional

    Edge-face connectivity array.

    unvsxr.DataArray, optional

    Precomputed unit normal vectors.

    varname_unvsstr, optional

    Name of the unit normal vector variable in the dataset.

Returns:

Reconstructed vector field with dimensions including nCartesian_coords.

Return type:

xr.DataArray or xu.UgridDataArray

Raises:

IOError – If the input dataset is not a xu.UgridDataset.

Notes

  • The magnitude variable is assumed to represent values in the direction of the edge-normal vector.

  • The vector is reconstructed as:

    \[\mathbf{u} = u_n \, \mathbf{n}\]

    where \(u_n\) is the scalar magnitude and \(\mathbf{n}\) is the unit normal vector.

  • Edge orientation is corrected using edge-face connectivity to ensure consistent direction across faces.

  • The result is stored in the dataset as {varname}c.

  • Designed primarily for velocity variables such as u0 and u1.

pySVA.sva_helpers.uda_to_edges(constructorSVA, uda)[source]

Interpolate a face-centered variable to edges.

Performs inverse-distance weighted interpolation using edge-face connectivity.

Parameters:
  • constructorSVA (object) – Constructor object containing the UGRID dataset.

  • uda (xu.UgridDataArray) – Input variable defined on faces.

Returns:

Interpolated variable on edges.

Return type:

xu.UgridDataArray

Notes

  • Uses weights from build_edge_face_weights().

  • Missing neighbors are handled via NaN masking.

  • Result is assigned location attribute 'edge'.

Deprecated since version 1.0.0: This function will be removed in version 1.2.0. Use :meth: constructorSVA.uda_to_edges()` instead.