Skip to content

Net/grid file

The computational grid for D-Flow FM is stored in the net file. It is an unstructured grid, where 1D and 2D can be combined. The file format is NetCDF, using the CF-, UGRID- and Deltares-conventions.

The net file is represented by the classes below.

Model

Branch

generate_nodes(self, mesh1d_edge_length: float, structure_chainage: Optional[List[float]] = None, max_dist_to_struc: Optional[float] = None)

Generate the branch offsets and the nodes.

Parameters:

Name Type Description Default
mesh1d_edge_length float

The edge length of the 1d mesh.

required
structure_chainage Optional[List[float]]

A list with the structure chainages. If not specified, calculation will not take it into account. Defaults to None.

None
max_dist_to_struc Optional[float]

The maximum distance from a node to a structure. If not specified, calculation will not take it into account. Defaults to None.

None

Exceptions:

Type Description
ValueError

Raised when any of the structure offsets, if specified, is smaller than zero.

ValueError

Raised when any of the structure offsets, if specified, is greater than the branch length.

Source code in hydrolib/core/dflowfm/net/models.py
def generate_nodes(
    self,
    mesh1d_edge_length: float,
    structure_chainage: Optional[List[float]] = None,
    max_dist_to_struc: Optional[float] = None,
):
    """Generate the branch offsets and the nodes.

    Args:
        mesh1d_edge_length (float): The edge length of the 1d mesh.
        structure_chainage (Optional[List[float]], optional): A list with the structure chainages. If not specified, calculation will not take it into account. Defaults to None.
        max_dist_to_struc (Optional[float], optional): The maximum distance from a node to a structure. If not specified, calculation will not take it into account. Defaults to None.

    Raises:
        ValueError: Raised when any of the structure offsets, if specified, is smaller than zero.
        ValueError: Raised when any of the structure offsets, if specified, is greater than the branch length.
    """
    # Generate offsets
    self.branch_offsets = self._generate_offsets(
        mesh1d_edge_length, structure_chainage, max_dist_to_struc
    )
    # Calculate node positions
    self.node_xy = self.interpolate(self.branch_offsets)
    # Add mask (all False)
    self.mask = np.full(self.branch_offsets.shape, False)

interpolate(self, distance: npt.ArrayLike) -> np.ndarray

Interpolate coordinates along branch by length

Parameters:

Name Type Description Default
distance npt.ArrayLike

Length

required
Source code in hydrolib/core/dflowfm/net/models.py
def interpolate(self, distance: npt.ArrayLike) -> np.ndarray:
    """Interpolate coordinates along branch by length

    Args:
        distance (npt.ArrayLike): Length
    """
    intpcoords = np.stack(
        [
            np.interp(distance, self._distance, self._x_coordinates),
            np.interp(distance, self._distance, self._y_coordinates),
        ],
        axis=1,
    )

    return intpcoords

Link1d2d (BaseModel) pydantic-model

Link1d2d defines the 1D2D Links of a model network.

Attributes:

Name Type Description
meshkernel Optional[mk.MeshKernel]

The MeshKernel used to interact with this Link1d2d

link1d2d_id np.ndarray

The id of this Link1d2d

link1d2d_long_name np.ndarray

The long name of this Link1d2d

link1d2d_contact_type np.ndarray

The contact type of this Link1d2d

link1d2d np.ndarray

The underlying data object of this Link1d2d

clear(self) -> None

Remove all saved links from the links administration

Source code in hydrolib/core/dflowfm/net/models.py
def clear(self) -> None:
    """Remove all saved links from the links administration"""
    self.link1d2d_id = np.empty(0, object)
    self.link1d2d_long_name = np.empty(0, object)
    self.link1d2d_contact_type = np.empty(0, np.int32)
    self.link1d2d = np.empty((0, 2), np.int32)
    # The meshkernel object needs to be resetted
    self.meshkernel._deallocate_state()
    self.meshkernel._allocate_state(self.meshkernel.is_geographic)
    self.meshkernel.contacts_get()

is_empty(self) -> bool

Whether this Link1d2d is currently empty.

Returns:

Type Description
bool

True if the Link1d2d is currently empty; False otherwise.

Source code in hydrolib/core/dflowfm/net/models.py
def is_empty(self) -> bool:
    """Whether this Link1d2d is currently empty.

    Returns:
        bool: True if the Link1d2d is currently empty; False otherwise.
    """
    return self.link1d2d.size == 0

read_file(self, file_path: Path) -> None

Read the Link1d2d data from the specified netCDF file at file_path into this

Parameters:

Name Type Description Default
file_path Path

Path to the netCDF file.

required
Source code in hydrolib/core/dflowfm/net/models.py
def read_file(self, file_path: Path) -> None:
    """Read the Link1d2d data from the specified netCDF file at file_path into this

    Args:
        file_path (Path): Path to the netCDF file.
    """

    reader = UgridReader(file_path)
    reader.read_link1d2d(self)

Mesh1d (BaseModel) pydantic-model

get_node_mask(self, branchids: List[str] = None)

Get node mask, give a mask with True for each node that is in the given branchid list

Source code in hydrolib/core/dflowfm/net/models.py
def get_node_mask(self, branchids: List[str] = None):
    """Get node mask, give a mask with True for each node that is in the given branchid list"""

    mask = np.full(self.mesh1d_node_id.shape, False, dtype=bool)
    if branchids is None:
        mask[:] = True
        return mask

    # Get number (index) of given branches
    idx = np.where(np.isin(self.network1d_branch_id, branchids))[0]
    if idx.size == 0:
        raise KeyError("No branches corresponding to the given keys were found.")

    mask[np.isin(self.mesh1d_node_branch_id, idx)] = True

    return mask

Mesh2d (BaseModel) pydantic-model

Mesh2d defines a single two dimensional grid.

Attributes:

Name Type Description
meshkernel mk.MeshKernel

The meshkernel used to manimpulate this Mesh2d.

mesh2d_node_x np.ndarray

The node positions on the x-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_node_y np.ndarray

The node positions on the y-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_node_z np.ndarray

The node positions on the z-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_edge_x np.ndarray

The edge positions on the x-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_edge_y np.ndarray

The edge positions on the y-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_edge_z np.ndarray

The edge positions on the z-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_edge_nodes np.ndarray

The mapping of edges to node indices. Defaults to np.empty((0, 2), dtype=np.int32).

mesh2d_face_x np.ndarray

The face positions on the x-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_face_y np.ndarray

The face positions on the y-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_face_z np.ndarray

The face positions on the z-axis. Defaults to np.empty(0, dtype=np.double).

mesh2d_face_nodes np.ndarray

The mapping of faces to node indices. Defaults to np.empty((0, 0), dtype=np.int32)

clip(self, geometrylist: mk.GeometryList, deletemeshoption: int = 1, inside = False) -> None

Clip the 2D mesh by a polygon. Both outside the exterior and inside the interiors is clipped

Parameters:

Name Type Description Default
geometrylist GeometryList

Polygon stored as GeometryList

required
deletemeshoption int

[description]. Defaults to 1.

1
Source code in hydrolib/core/dflowfm/net/models.py
def clip(
    self,
    geometrylist: mk.GeometryList,
    deletemeshoption: int = 1,
    inside=False,
) -> None:
    """Clip the 2D mesh by a polygon. Both outside the exterior and inside the interiors is clipped

    Args:
        geometrylist (GeometryList): Polygon stored as GeometryList
        deletemeshoption (int, optional): [description]. Defaults to 1.
    """

    # Add current mesh to Mesh2d instance
    self._set_mesh2d()

    deletemeshoption = mk.DeleteMeshOption(deletemeshoption)

    # For clipping outside
    if not inside:
        # Check if a multipolygon was provided when clipping outside
        if geometrylist.geometry_separator in geometrylist.x_coordinates:
            raise NotImplementedError(
                "Deleting outside more than a single exterior (MultiPolygon) is not implemented."
            )

        # Get exterior and interiors
        parts = split_by(geometrylist, geometrylist.inner_outer_separator)

        exteriors = [parts[0]]
        interiors = parts[1:]

    # Inside
    else:
        # Check if any polygon contains holes, when clipping inside
        if geometrylist.inner_outer_separator in geometrylist.x_coordinates:
            raise NotImplementedError(
                "Deleting inside a (Multi)Polygon with holes is not implemented."
            )

        # Get exterior and interiors
        parts = split_by(geometrylist, geometrylist.geometry_separator)

        exteriors = parts[:]
        interiors = []

    # Check if parts are closed
    for part in exteriors + interiors:
        if (part.x_coordinates[0], part.y_coordinates[0]) != (
            part.x_coordinates[-1],
            part.y_coordinates[-1],
        ):
            raise ValueError(
                "First and last coordinate of each GeometryList part should match."
            )

    # Delete everything outside the (Multi)Polygon
    for exterior in exteriors:
        self.meshkernel.mesh2d_delete(
            geometry_list=exterior,
            delete_option=deletemeshoption,
            invert_deletion=not inside,
        )

    # Delete all holes.
    for interior in interiors:
        self.meshkernel.mesh2d_delete(
            geometry_list=interior,
            delete_option=deletemeshoption,
            invert_deletion=inside,
        )

    # Process
    self._process(self.meshkernel.mesh2d_get())

create_rectilinear(self, extent: tuple, dx: float, dy: float) -> None

Create a rectilinear mesh within a polygon. A rectangular grid is generated within the polygon bounds

Parameters:

Name Type Description Default
extent tuple

Bounding box of mesh (left, bottom, right, top)

required
dx float

Horizontal distance

required
dy float

Vertical distance

required

Exceptions:

Type Description
NotImplementedError

MultiPolygons

Source code in hydrolib/core/dflowfm/net/models.py
def create_rectilinear(self, extent: tuple, dx: float, dy: float) -> None:
    """Create a rectilinear mesh within a polygon. A rectangular grid is generated within the polygon bounds

    Args:
        extent (tuple): Bounding box of mesh (left, bottom, right, top)
        dx (float): Horizontal distance
        dy (float): Vertical distance

    Raises:
        NotImplementedError: MultiPolygons
    """

    xmin, ymin, xmax, ymax = extent

    # Generate mesh
    mesh2d_input = mk.Mesh2dFactory.create(
        rows=int((ymax - ymin) / dy),
        columns=int((xmax - xmin) / dx),
        origin_x=xmin,
        origin_y=ymin,
        spacing_x=dx,
        spacing_y=dy,
    )

    # Process
    self._process(mesh2d_input)

create_triangular(self, geometry_list: mk.GeometryList) -> None

Create triangular grid within GeometryList object

Parameters:

Name Type Description Default
geometry_list mk.GeometryList

GeometryList represeting a polygon within which the mesh is generated.

required
Source code in hydrolib/core/dflowfm/net/models.py
def create_triangular(self, geometry_list: mk.GeometryList) -> None:
    """Create triangular grid within GeometryList object

    Args:
        geometry_list (mk.GeometryList): GeometryList represeting a polygon within which the mesh is generated.
    """
    # Call meshkernel
    self.meshkernel.mesh2d_make_mesh_from_polygon(geometry_list)

    # Process new mesh
    self._process(self.get_mesh2d())

get_mesh2d(self) -> mk.Mesh2d

Get the mesh2d as represented in the MeshKernel

Returns:

Type Description
(mk.Mesh2d)

The mesh2d as represented in the MeshKernel

Source code in hydrolib/core/dflowfm/net/models.py
def get_mesh2d(self) -> mk.Mesh2d:
    """Get the mesh2d as represented in the MeshKernel

    Returns:
        (mk.Mesh2d): The mesh2d as represented in the MeshKernel
    """
    return self.meshkernel.mesh2d_get()

is_empty(self) -> bool

Determine whether this Mesh2d is empty.

Returns:

Type Description
(bool)

Whether this Mesh2d is empty.

Source code in hydrolib/core/dflowfm/net/models.py
def is_empty(self) -> bool:
    """Determine whether this Mesh2d is empty.

    Returns:
        (bool): Whether this Mesh2d is empty.
    """
    return self.mesh2d_node_x.size == 0

read_file(self, file_path: Path) -> None

Read the Mesh2d from the file at file_path.

Parameters:

Name Type Description Default
file_path Path

Path to the file to be read.

required
Source code in hydrolib/core/dflowfm/net/models.py
def read_file(self, file_path: Path) -> None:
    """Read the Mesh2d from the file at file_path.

    Args:
        file_path (Path): Path to the file to be read.
    """
    reader = UgridReader(file_path)
    reader.read_mesh2d(self)

refine(self, polygon: mk.GeometryList, level: int)

Refine the mesh within a polygon, by a number of steps (level)

Parameters:

Name Type Description Default
polygon GeometryList

Polygon in which to refine

required
level int

Number of refinement steps

required
Source code in hydrolib/core/dflowfm/net/models.py
def refine(self, polygon: mk.GeometryList, level: int):
    """Refine the mesh within a polygon, by a number of steps (level)

    Args:
        polygon (GeometryList): Polygon in which to refine
        level (int): Number of refinement steps
    """
    # Add current mesh to Mesh2d instance
    mesh2d_input = mk.Mesh2d(
        node_x=self.mesh2d_node_x,
        node_y=self.mesh2d_node_y,
        edge_nodes=self.mesh2d_edge_nodes.ravel(),
    )
    self.meshkernel.mesh2d_set(mesh2d_input)

    # Check if parts are closed
    # if not (polygon.x_coordinates[0], polygon.y_coordinates[0]) == (
    #     polygon.x_coordinates[-1],
    #     polygon.y_coordinates[-1],
    # ):
    #     raise ValueError("First and last coordinate of each GeometryList part should match.")

    parameters = mk.MeshRefinementParameters(
        refine_intersected=True,
        use_mass_center_when_refining=False,
        min_face_size=10.0,  # Does nothing?
        refinement_type=1,  # No effect?
        connect_hanging_nodes=True,
        account_for_samples_outside_face=False,
        max_refinement_iterations=level,
    )
    self.meshkernel.mesh2d_refine_based_on_polygon(polygon, parameters)

    # Process
    self._process(self.meshkernel.mesh2d_get())

Network

from_file(file_path: Path) -> Network classmethod

Read network from file. This classmethod checks what mesh components (mesh1d & network1d, mesh2d, link1d2d) are present, and loads them one by one.

Parameters:

Name Type Description Default
file_path Path

path to netcdf file with network data

required

Returns:

Type Description
Network

The instance of the class itself that is returned

Source code in hydrolib/core/dflowfm/net/models.py
@classmethod
def from_file(cls, file_path: Path) -> Network:
    """Read network from file. This classmethod checks what mesh components (mesh1d & network1d, mesh2d, link1d2d) are
    present, and loads them one by one.

    Args:
        file_path (Path): path to netcdf file with network data

    Returns:
        Network: The instance of the class itself that is returned
    """

    network = cls()
    ds = nc.Dataset(file_path)  # type: ignore[import]

    reader = UgridReader(file_path)

    reader.read_mesh1d_network1d(network._mesh1d)
    reader.read_mesh2d(network._mesh2d)
    reader.read_link1d2d(network._link1d2d)

    ds.close()

    return network

mesh2d_create_triangular_within_polygon(self, polygon: mk.GeometryList) -> None

Create triangular grid within GeometryList object. Calls _mesh2d.create_triangular directly, but is easier accessible for users.

Parameters:

Name Type Description Default
polygon mk.GeometryList

GeometryList representing a polygon within which the mesh is generated.

required
Source code in hydrolib/core/dflowfm/net/models.py
def mesh2d_create_triangular_within_polygon(self, polygon: mk.GeometryList) -> None:
    """Create triangular grid within GeometryList object. Calls _mesh2d.create_triangular
    directly, but is easier accessible for users.

    Args:
        polygon (mk.GeometryList): GeometryList representing a polygon within which the mesh is generated.
    """
    self._mesh2d.create_triangular(geometry_list=polygon)

to_file(self, file: Path) -> None

Write network to file

Parameters:

Name Type Description Default
file Path

File where _net.nc is written to.

required
Source code in hydrolib/core/dflowfm/net/models.py
def to_file(self, file: Path) -> None:
    """Write network to file

    Args:
        file (Path): File where _net.nc is written to.
    """

    writer = UgridWriter()
    writer.write(self, file)

NetworkModel (ParsableFileModel) pydantic-model

Network model representation.

split_by(gl: mk.GeometryList, by: float) -> list

Function to split mk.GeometryList by seperator.

Parameters:

Name Type Description Default
gl mk.GeometryList

The geometry list to split.

required
by float

The value by which to split the gl.

required

Returns:

Type Description
list

The split lists.

Source code in hydrolib/core/dflowfm/net/models.py
def split_by(gl: mk.GeometryList, by: float) -> list:
    """Function to split mk.GeometryList by seperator.

    Args:
        gl (mk.GeometryList): The geometry list to split.
        by (float): The value by which to split the gl.

    Returns:
        list: The split lists.
    """
    x, y = gl.x_coordinates.copy(), gl.y_coordinates.copy()
    idx = np.where(x == by)[0]

    xparts = np.split(x, idx)
    yparts = np.split(y, idx)

    lists = [
        mk.GeometryList(xp[min(i, 1) :], yp[min(i, 1) :])
        for i, (xp, yp) in enumerate(zip(xparts, yparts))
    ]

    return lists
Back to top