Skip to content

Api

DHyDAMO

The API references of DHyDAMO can be found below:

CrossSectionsIO

Source code in hydrolib/dhydamo/converters/hydamo2df.py
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
class CrossSectionsIO:
    def __init__(self, crosssections):
        self.crosssections = crosssections

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def from_datamodel(
        self, crsdefs: pd.DataFrame = None, crslocs: pd.DataFrame = None
    ) -> None:
        """ "
        From parsed data models of crsdefs and crslocs
        """

        if crslocs is not None:
            for _, crsloc in crslocs.iterrows():
                # add location
                self.crosssections.add_crosssection_location(
                    branchid=crsloc["branch_id"],
                    chainage=crsloc["branch_offset"],
                    shift=crsloc["shift"],
                    definition=crsloc["crosssectiondefinitionid"],
                )

        if crsdefs is not None:
            crsdefs = crsdefs.drop_duplicates(subset=["crosssectiondefinitionid"])
            for _, crsdef in crsdefs.iterrows():
                # Set roughness value on default if cross-section has non defined (e.g. culverts)
                if isinstance(crsdef["frictionid"], str):
                    roughtype = crsdef["frictionid"].split("_")[0]
                else:
                    roughtype = "Chezy"

                if isinstance(crsdef["frictionid"], str):
                    roughval = float(crsdef["frictionid"].split("_")[-1])
                else:
                    roughval = 45

                # add definition
                if crsdef["type"] == "circle":
                    self.crosssections.add_circle_definition(
                        diameter=crsdef["diameter"],
                        roughnesstype=roughtype,
                        roughnessvalue=roughval,
                        name=crsdef["crosssectiondefinitionid"],
                    )
                elif crsdef["type"] == "rectangle":
                    self.crosssections.add_rectangle_definition(
                        height=crsdef["height"],
                        width=crsdef["width"],
                        closed=crsdef["closed"],
                        roughnesstype=roughtype,
                        roughnessvalue=roughval,
                        name=crsdef["crosssectiondefinitionid"],
                    )
                elif crsdef["type"] == "trapezium":
                    self.crosssections.add_trapezium_definition(
                        slope=(crsdef["t_width"] - crsdef["width"])
                        / 2
                        / crsdef["height"],
                        maximumflowwidth=crsdef["t_width"],
                        bottomwidth=crsdef["width"],
                        closed=crsdef["closed"],
                        roughnesstype=roughtype,
                        roughnessvalue=roughval,
                        name=crsdef["crosssectiondefinitionid"],
                    )
                elif crsdef["type"] == "zw":
                    self.crosssections.add_zw_definition(
                        numLevels=crsdef["numlevels"],
                        levels=crsdef["levels"],
                        flowWidths=crsdef["flowwidths"],
                        totalWidths=crsdef["totalwidths"],
                        roughnesstype=roughtype,
                        roughnessvalue=roughval,
                        name=crsdef["crosssectiondefinitionid"],
                    )

                else:
                    raise NotImplementedError

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def profiles(
        self,
        branches: ExtendedGeoDataFrame, # = None,
        roughness_variant: RoughnessVariant,# =None,
        crosssections: Optional[ExtendedGeoDataFrame] = None,
        crosssection_roughness: Optional[ExtendedDataFrame] = None,
        profile_groups: Optional[ExtendedDataFrame] = None,
        profile_lines: Optional[ExtendedGeoDataFrame] = None,
        param_profile: Optional[ExtendedDataFrame] = None,
        param_profile_values: Optional[ExtendedDataFrame] = None,
    ) -> None:
        """
        Method to add cross section from hydamo files. Two files
        can be handed to the function, the cross section file (dwarsprofiel) and the
        parametrised file (normgeparametriseerd). The
        hierarchical order is 1. dwarsprofiel, 2. normgeparametriseerd.
        Each branch will be assigned a profile following this order. If parametrised
        and standard are not given, branches can be wit hout cross section. In that case
        a standard profile should be assigned
        """
        dp_branches = None
        dp_structures = None
        if profile_groups is not None:
            # check for profile_groups items with valid brugid or stuwid. They need to be droppped from profiles.
            groupidx = [
                idx
                for idx, group in profile_groups.iterrows()
                if ("brugid" in profile_groups.columns) & (not pd.isnull(group.brugid))
            ]

            groupidx = groupidx + [
                idx
                for idx, group in profile_groups.iterrows()
                if ("stuwid" in profile_groups.columns) & (not pd.isnull(group.stuwid))
            ]

            # index of the lines that are associated to these groups
            lineidx = [
                profile_lines[
                    profile_lines["profielgroepid"]
                    == profile_groups.loc[grindex, "globalid"]
                ].index.values[0]
                for grindex in groupidx
            ]
            # index of the profiles associated to these lines
            profidx = [
                crosssections[
                    crosssections["profiellijnid"]
                    == profile_lines.loc[lindex, "globalid"]
                ].index.values[0]
                for lindex in lineidx
            ]
            # make a copy and drop the profiles corresponding to a structure
            dp_branches = crosssections.copy(deep=True)
            dp_branches.drop(profidx, axis=0, inplace=True)

            dp_structures = crosssections.copy(deep=True)
            dp_structures = dp_structures.loc[profidx, :]
        else:
            dp_branches = crosssections.copy(deep=True)

        # Assign cross-sections to branches
        nnocross = len(self.crosssections.get_branches_without_crosssection())
        print(
            f"Before adding the number of branches without cross section is: {nnocross}."
        )

        if dp_branches is not None:
            # 1. Collect cross sections from 'dwarsprofielen'
            yz_profiles = self.crosssections.crosssection_to_yzprofiles(
                dp_branches,
                crosssection_roughness,
                branches,
                roughness_variant=roughness_variant,
            )

            for name, css in yz_profiles.items():
                # Add definition
                self.crosssections.add_yz_definition(
                    yz=css["yz"],
                    thalweg=css["thalweg"],
                    name=name,
                    roughnesstype=css["typeruwheid"],
                    roughnessvalue=css["ruwheid"],
                )
                # Add location
                self.crosssections.add_crosssection_location(
                    branchid=css["branchid"], chainage=css["chainage"], definition=name
                )

        # Check the number of branches with cross sections
        no_crosssection_id = self.crosssections.get_branches_without_crosssection()
        no_crosssection = [
            b for b in branches.itertuples() if b.code in no_crosssection_id
        ]

        nnocross = len(no_crosssection)
        print(
            f"After adding 'dwarsprofielen' the number of branches without cross section is: {nnocross}."
        )
        if nnocross == 0:
            print("No further branches without a profile.")
        elif param_profile is None:
            print("No parametrised crossections available for branches.")
        else:
            # Derive norm cross sections for norm parametrised
            param_profiles_converted = self.crosssections.parametrised_to_profiles(
                param_profile,
                param_profile_values,
                no_crosssection,
                roughness_variant=roughness_variant,
            )
            # Get branch information
            branchdata = self.crosssections.hydamo.branches.loc[
                list(param_profiles_converted.keys())
            ]
            branchdata["chainage_upper"] = 0.05 * branchdata.length
            branchdata["chainage_lower"] = 0.95 * branchdata.length


            # Add cross sections
            for branchid, css in param_profiles_converted.items():
                chainage_upper = branchdata.at[branchid, "chainage_upper"]
                chainage_lower = branchdata.at[branchid, "chainage_lower"]

                if css["type"] == "rectangle":
                    name = self.crosssections.add_rectangle_definition(
                        height=css["height"],
                        width=css["width"],
                        closed=css["closed"],
                        roughnesstype=css["typeruwheid"],
                        roughnessvalue=css["ruwheid"],
                    )

                if css["type"] == "trapezium":
                    name = self.crosssections.add_trapezium_definition(
                        slope=css["slope"],
                        maximumflowwidth=css["maximumflowwidth"],
                        bottomwidth=css["bottomwidth"],
                        closed=css["closed"],
                        roughnesstype=css["typeruwheid"],
                        roughnessvalue=css["ruwheid"],
                    )

                # Add location
                self.crosssections.add_crosssection_location(
                    branchid=branchid,
                    chainage=chainage_upper,
                    definition=name,
                    shift=css["bottomlevel_upper"],
                )

                self.crosssections.add_crosssection_location(
                    branchid=branchid,
                    chainage=chainage_lower,
                    definition=name,
                    shift=css["bottomlevel_lower"],
                )

        nnocross = len(self.crosssections.get_branches_without_crosssection())
        print(
            f"After adding 'normgeparametriseerd' the number of branches without cross section is: {nnocross}."
        )

        if dp_structures is not None:
            # 1. Collect cross sections from 'dwarsprofielen'
            yz_profiles = self.crosssections.crosssection_to_yzprofiles(
                dp_structures,
                crosssection_roughness,
                None,
                roughness_variant=roughness_variant,
            )

            for name, css in yz_profiles.items():
                # Add definition
                self.crosssections.add_yz_definition(
                    yz=css["yz"],
                    thalweg=css["thalweg"],
                    name=name,
                    roughnesstype=css["typeruwheid"],
                    roughnessvalue=css["ruwheid"],
                )

from_datamodel(crsdefs: pd.DataFrame = None, crslocs: pd.DataFrame = None) -> None

" From parsed data models of crsdefs and crslocs

Source code in hydrolib/dhydamo/converters/hydamo2df.py
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def from_datamodel(
    self, crsdefs: pd.DataFrame = None, crslocs: pd.DataFrame = None
) -> None:
    """ "
    From parsed data models of crsdefs and crslocs
    """

    if crslocs is not None:
        for _, crsloc in crslocs.iterrows():
            # add location
            self.crosssections.add_crosssection_location(
                branchid=crsloc["branch_id"],
                chainage=crsloc["branch_offset"],
                shift=crsloc["shift"],
                definition=crsloc["crosssectiondefinitionid"],
            )

    if crsdefs is not None:
        crsdefs = crsdefs.drop_duplicates(subset=["crosssectiondefinitionid"])
        for _, crsdef in crsdefs.iterrows():
            # Set roughness value on default if cross-section has non defined (e.g. culverts)
            if isinstance(crsdef["frictionid"], str):
                roughtype = crsdef["frictionid"].split("_")[0]
            else:
                roughtype = "Chezy"

            if isinstance(crsdef["frictionid"], str):
                roughval = float(crsdef["frictionid"].split("_")[-1])
            else:
                roughval = 45

            # add definition
            if crsdef["type"] == "circle":
                self.crosssections.add_circle_definition(
                    diameter=crsdef["diameter"],
                    roughnesstype=roughtype,
                    roughnessvalue=roughval,
                    name=crsdef["crosssectiondefinitionid"],
                )
            elif crsdef["type"] == "rectangle":
                self.crosssections.add_rectangle_definition(
                    height=crsdef["height"],
                    width=crsdef["width"],
                    closed=crsdef["closed"],
                    roughnesstype=roughtype,
                    roughnessvalue=roughval,
                    name=crsdef["crosssectiondefinitionid"],
                )
            elif crsdef["type"] == "trapezium":
                self.crosssections.add_trapezium_definition(
                    slope=(crsdef["t_width"] - crsdef["width"])
                    / 2
                    / crsdef["height"],
                    maximumflowwidth=crsdef["t_width"],
                    bottomwidth=crsdef["width"],
                    closed=crsdef["closed"],
                    roughnesstype=roughtype,
                    roughnessvalue=roughval,
                    name=crsdef["crosssectiondefinitionid"],
                )
            elif crsdef["type"] == "zw":
                self.crosssections.add_zw_definition(
                    numLevels=crsdef["numlevels"],
                    levels=crsdef["levels"],
                    flowWidths=crsdef["flowwidths"],
                    totalWidths=crsdef["totalwidths"],
                    roughnesstype=roughtype,
                    roughnessvalue=roughval,
                    name=crsdef["crosssectiondefinitionid"],
                )

            else:
                raise NotImplementedError

profiles(branches: ExtendedGeoDataFrame, roughness_variant: RoughnessVariant, crosssections: Optional[ExtendedGeoDataFrame] = None, crosssection_roughness: Optional[ExtendedDataFrame] = None, profile_groups: Optional[ExtendedDataFrame] = None, profile_lines: Optional[ExtendedGeoDataFrame] = None, param_profile: Optional[ExtendedDataFrame] = None, param_profile_values: Optional[ExtendedDataFrame] = None) -> None

Method to add cross section from hydamo files. Two files can be handed to the function, the cross section file (dwarsprofiel) and the parametrised file (normgeparametriseerd). The hierarchical order is 1. dwarsprofiel, 2. normgeparametriseerd. Each branch will be assigned a profile following this order. If parametrised and standard are not given, branches can be wit hout cross section. In that case a standard profile should be assigned

Source code in hydrolib/dhydamo/converters/hydamo2df.py
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def profiles(
    self,
    branches: ExtendedGeoDataFrame, # = None,
    roughness_variant: RoughnessVariant,# =None,
    crosssections: Optional[ExtendedGeoDataFrame] = None,
    crosssection_roughness: Optional[ExtendedDataFrame] = None,
    profile_groups: Optional[ExtendedDataFrame] = None,
    profile_lines: Optional[ExtendedGeoDataFrame] = None,
    param_profile: Optional[ExtendedDataFrame] = None,
    param_profile_values: Optional[ExtendedDataFrame] = None,
) -> None:
    """
    Method to add cross section from hydamo files. Two files
    can be handed to the function, the cross section file (dwarsprofiel) and the
    parametrised file (normgeparametriseerd). The
    hierarchical order is 1. dwarsprofiel, 2. normgeparametriseerd.
    Each branch will be assigned a profile following this order. If parametrised
    and standard are not given, branches can be wit hout cross section. In that case
    a standard profile should be assigned
    """
    dp_branches = None
    dp_structures = None
    if profile_groups is not None:
        # check for profile_groups items with valid brugid or stuwid. They need to be droppped from profiles.
        groupidx = [
            idx
            for idx, group in profile_groups.iterrows()
            if ("brugid" in profile_groups.columns) & (not pd.isnull(group.brugid))
        ]

        groupidx = groupidx + [
            idx
            for idx, group in profile_groups.iterrows()
            if ("stuwid" in profile_groups.columns) & (not pd.isnull(group.stuwid))
        ]

        # index of the lines that are associated to these groups
        lineidx = [
            profile_lines[
                profile_lines["profielgroepid"]
                == profile_groups.loc[grindex, "globalid"]
            ].index.values[0]
            for grindex in groupidx
        ]
        # index of the profiles associated to these lines
        profidx = [
            crosssections[
                crosssections["profiellijnid"]
                == profile_lines.loc[lindex, "globalid"]
            ].index.values[0]
            for lindex in lineidx
        ]
        # make a copy and drop the profiles corresponding to a structure
        dp_branches = crosssections.copy(deep=True)
        dp_branches.drop(profidx, axis=0, inplace=True)

        dp_structures = crosssections.copy(deep=True)
        dp_structures = dp_structures.loc[profidx, :]
    else:
        dp_branches = crosssections.copy(deep=True)

    # Assign cross-sections to branches
    nnocross = len(self.crosssections.get_branches_without_crosssection())
    print(
        f"Before adding the number of branches without cross section is: {nnocross}."
    )

    if dp_branches is not None:
        # 1. Collect cross sections from 'dwarsprofielen'
        yz_profiles = self.crosssections.crosssection_to_yzprofiles(
            dp_branches,
            crosssection_roughness,
            branches,
            roughness_variant=roughness_variant,
        )

        for name, css in yz_profiles.items():
            # Add definition
            self.crosssections.add_yz_definition(
                yz=css["yz"],
                thalweg=css["thalweg"],
                name=name,
                roughnesstype=css["typeruwheid"],
                roughnessvalue=css["ruwheid"],
            )
            # Add location
            self.crosssections.add_crosssection_location(
                branchid=css["branchid"], chainage=css["chainage"], definition=name
            )

    # Check the number of branches with cross sections
    no_crosssection_id = self.crosssections.get_branches_without_crosssection()
    no_crosssection = [
        b for b in branches.itertuples() if b.code in no_crosssection_id
    ]

    nnocross = len(no_crosssection)
    print(
        f"After adding 'dwarsprofielen' the number of branches without cross section is: {nnocross}."
    )
    if nnocross == 0:
        print("No further branches without a profile.")
    elif param_profile is None:
        print("No parametrised crossections available for branches.")
    else:
        # Derive norm cross sections for norm parametrised
        param_profiles_converted = self.crosssections.parametrised_to_profiles(
            param_profile,
            param_profile_values,
            no_crosssection,
            roughness_variant=roughness_variant,
        )
        # Get branch information
        branchdata = self.crosssections.hydamo.branches.loc[
            list(param_profiles_converted.keys())
        ]
        branchdata["chainage_upper"] = 0.05 * branchdata.length
        branchdata["chainage_lower"] = 0.95 * branchdata.length


        # Add cross sections
        for branchid, css in param_profiles_converted.items():
            chainage_upper = branchdata.at[branchid, "chainage_upper"]
            chainage_lower = branchdata.at[branchid, "chainage_lower"]

            if css["type"] == "rectangle":
                name = self.crosssections.add_rectangle_definition(
                    height=css["height"],
                    width=css["width"],
                    closed=css["closed"],
                    roughnesstype=css["typeruwheid"],
                    roughnessvalue=css["ruwheid"],
                )

            if css["type"] == "trapezium":
                name = self.crosssections.add_trapezium_definition(
                    slope=css["slope"],
                    maximumflowwidth=css["maximumflowwidth"],
                    bottomwidth=css["bottomwidth"],
                    closed=css["closed"],
                    roughnesstype=css["typeruwheid"],
                    roughnessvalue=css["ruwheid"],
                )

            # Add location
            self.crosssections.add_crosssection_location(
                branchid=branchid,
                chainage=chainage_upper,
                definition=name,
                shift=css["bottomlevel_upper"],
            )

            self.crosssections.add_crosssection_location(
                branchid=branchid,
                chainage=chainage_lower,
                definition=name,
                shift=css["bottomlevel_lower"],
            )

    nnocross = len(self.crosssections.get_branches_without_crosssection())
    print(
        f"After adding 'normgeparametriseerd' the number of branches without cross section is: {nnocross}."
    )

    if dp_structures is not None:
        # 1. Collect cross sections from 'dwarsprofielen'
        yz_profiles = self.crosssections.crosssection_to_yzprofiles(
            dp_structures,
            crosssection_roughness,
            None,
            roughness_variant=roughness_variant,
        )

        for name, css in yz_profiles.items():
            # Add definition
            self.crosssections.add_yz_definition(
                yz=css["yz"],
                thalweg=css["thalweg"],
                name=name,
                roughnesstype=css["typeruwheid"],
                roughnessvalue=css["ruwheid"],
            )

ExternalForcingsIO

Source code in hydrolib/dhydamo/converters/hydamo2df.py
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
class ExternalForcingsIO:
    def __init__(self, external_forcings):
        self.external_forcings = external_forcings

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def boundaries(
        self, boundary_conditions: ExtendedGeoDataFrame, mesh1d: Network = None
    ) -> None:
        """
        Generate boundary conditions from hydamo 'randvoorwaarden' file. The file format does not allow for timeseries.

        To add a time series, use 'add_boundary' from the workflow.

        Parameters
        ----------
        boundary_conditions: gpd.GeoDataFrame
            geodataframe with the locations and properties of the boundary conditions

        Returns
        -------
        dictionary
            Dictionary with attributes of boundary conditions, usable for dflowfm
        """
        # Read from Hydamo
        bcdct = {}

        for bndcnd in boundary_conditions.itertuples():
            if "waterstand" in bndcnd.typerandvoorwaarde:
                quantity = "waterlevelbnd"
            elif "debiet" in bndcnd.typerandvoorwaarde:
                quantity = "dischargebnd"

            # Add boundary condition
            bcdct[bndcnd.code] = {
                "code": bndcnd.code,
                "quantity": quantity,
                "value": bndcnd.waterstand
                if not np.isnan(bndcnd.waterstand)
                else bndcnd.debiet,
                "time": None,
                "geometry": bndcnd.geometry,
            }

        # Add all items
        for key, item in bcdct.items():
            self.external_forcings.add_boundary_condition(
                key, item["geometry"], item["quantity"], item["value"], mesh1d=mesh1d
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def laterals(
        self,
        locations: ExtendedGeoDataFrame,
        overflows: Optional[ExtendedGeoDataFrame] = None,
        greenhouse_laterals: Optional[ExtendedGeoDataFrame] = None,
        lateral_discharges: Optional[Union[pd.DataFrame, pd.Series]]=None,
        rr_boundaries: Optional[dict] = None,
    ) -> None:
        """
        Process laterals

        Parameters
        ----------
        locations: gpd.GeoDataFrame
            GeoDataFrame with at least 'geometry' (Point) and the column 'code'
        lateral_discharges: pd.DataFrame
            DataFrame with lateral discharges. The index should be a time object (datetime or similar).
        rr_boundaries: pd.DataFrame
            DataFrame with RR-catchments that are coupled
        """

        if rr_boundaries is None:
            rr_boundaries = []

        # in case of 3d points, remove the 3rd dimension
        locations["geometry2"] = [
            Point([point.geometry.x, point.geometry.y])
            for _, point in locations.iterrows()
        ]
        locations.drop("geometry", inplace=True, axis=1)
        locations.rename(columns={"geometry2": "geometry"}, inplace=True)

        latdct = {}
        if overflows is not None:
            locations = pd.concat([locations, overflows], ignore_index=True)
        if greenhouse_laterals is not None:
            locations = pd.concat([locations, greenhouse_laterals], ignore_index=True)

        # Get time series and add to dictionary
        # for nidx, lateral in zip(nearest_idx, locations.itertuples()):
        for lateral in locations.itertuples():
            # Check if a time is provided for the lateral
            if lateral.code in rr_boundaries:
                # Add to dictionary
                latdct[lateral.code] = {
                    "branchid": lateral.branch_id,
                    "chainage": str(lateral.branch_offset),
                    "discharge": "realtime",
                }

            else:
                if lateral_discharges is None:
                    logger.warning(
                        f"No lateral_discharges provided. {lateral.code} expects them. Skipping."
                    )
                else:
                    if isinstance(lateral_discharges, pd.Series):
                        series = lateral_discharges.loc[lateral.code]

                        # Add to dictionary
                        latdct[lateral.code] = {
                            "branchid": lateral.branch_id,
                            "chainage": str(lateral.branch_offset),
                            "discharge": series,
                        }

                    else:
                        if lateral.code not in lateral_discharges.columns:
                            logger.warning(
                                f"No data found for {lateral.code}. Skipping."
                            )
                            continue

                        # Get timeseries
                        series = lateral_discharges.loc[:, lateral.code]

                        # Add to dictionary
                        latdct[lateral.code] = {
                            "branchid": lateral.branch_id,
                            "chainage": str(lateral.branch_offset),
                            "discharge": series,
                        }

        # Add all items
        for key, item in latdct.items():
            self.external_forcings.add_lateral(
                id=key,
                branchid=item["branchid"],
                chainage=item["chainage"],
                discharge=item["discharge"],
            )

boundaries(boundary_conditions: ExtendedGeoDataFrame, mesh1d: Network = None) -> None

Generate boundary conditions from hydamo 'randvoorwaarden' file. The file format does not allow for timeseries.

To add a time series, use 'add_boundary' from the workflow.

Parameters

boundary_conditions: gpd.GeoDataFrame geodataframe with the locations and properties of the boundary conditions

Returns

dictionary Dictionary with attributes of boundary conditions, usable for dflowfm

Source code in hydrolib/dhydamo/converters/hydamo2df.py
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def boundaries(
    self, boundary_conditions: ExtendedGeoDataFrame, mesh1d: Network = None
) -> None:
    """
    Generate boundary conditions from hydamo 'randvoorwaarden' file. The file format does not allow for timeseries.

    To add a time series, use 'add_boundary' from the workflow.

    Parameters
    ----------
    boundary_conditions: gpd.GeoDataFrame
        geodataframe with the locations and properties of the boundary conditions

    Returns
    -------
    dictionary
        Dictionary with attributes of boundary conditions, usable for dflowfm
    """
    # Read from Hydamo
    bcdct = {}

    for bndcnd in boundary_conditions.itertuples():
        if "waterstand" in bndcnd.typerandvoorwaarde:
            quantity = "waterlevelbnd"
        elif "debiet" in bndcnd.typerandvoorwaarde:
            quantity = "dischargebnd"

        # Add boundary condition
        bcdct[bndcnd.code] = {
            "code": bndcnd.code,
            "quantity": quantity,
            "value": bndcnd.waterstand
            if not np.isnan(bndcnd.waterstand)
            else bndcnd.debiet,
            "time": None,
            "geometry": bndcnd.geometry,
        }

    # Add all items
    for key, item in bcdct.items():
        self.external_forcings.add_boundary_condition(
            key, item["geometry"], item["quantity"], item["value"], mesh1d=mesh1d
        )

laterals(locations: ExtendedGeoDataFrame, overflows: Optional[ExtendedGeoDataFrame] = None, greenhouse_laterals: Optional[ExtendedGeoDataFrame] = None, lateral_discharges: Optional[Union[pd.DataFrame, pd.Series]] = None, rr_boundaries: Optional[dict] = None) -> None

Process laterals

Parameters

locations: gpd.GeoDataFrame GeoDataFrame with at least 'geometry' (Point) and the column 'code' lateral_discharges: pd.DataFrame DataFrame with lateral discharges. The index should be a time object (datetime or similar). rr_boundaries: pd.DataFrame DataFrame with RR-catchments that are coupled

Source code in hydrolib/dhydamo/converters/hydamo2df.py
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def laterals(
    self,
    locations: ExtendedGeoDataFrame,
    overflows: Optional[ExtendedGeoDataFrame] = None,
    greenhouse_laterals: Optional[ExtendedGeoDataFrame] = None,
    lateral_discharges: Optional[Union[pd.DataFrame, pd.Series]]=None,
    rr_boundaries: Optional[dict] = None,
) -> None:
    """
    Process laterals

    Parameters
    ----------
    locations: gpd.GeoDataFrame
        GeoDataFrame with at least 'geometry' (Point) and the column 'code'
    lateral_discharges: pd.DataFrame
        DataFrame with lateral discharges. The index should be a time object (datetime or similar).
    rr_boundaries: pd.DataFrame
        DataFrame with RR-catchments that are coupled
    """

    if rr_boundaries is None:
        rr_boundaries = []

    # in case of 3d points, remove the 3rd dimension
    locations["geometry2"] = [
        Point([point.geometry.x, point.geometry.y])
        for _, point in locations.iterrows()
    ]
    locations.drop("geometry", inplace=True, axis=1)
    locations.rename(columns={"geometry2": "geometry"}, inplace=True)

    latdct = {}
    if overflows is not None:
        locations = pd.concat([locations, overflows], ignore_index=True)
    if greenhouse_laterals is not None:
        locations = pd.concat([locations, greenhouse_laterals], ignore_index=True)

    # Get time series and add to dictionary
    # for nidx, lateral in zip(nearest_idx, locations.itertuples()):
    for lateral in locations.itertuples():
        # Check if a time is provided for the lateral
        if lateral.code in rr_boundaries:
            # Add to dictionary
            latdct[lateral.code] = {
                "branchid": lateral.branch_id,
                "chainage": str(lateral.branch_offset),
                "discharge": "realtime",
            }

        else:
            if lateral_discharges is None:
                logger.warning(
                    f"No lateral_discharges provided. {lateral.code} expects them. Skipping."
                )
            else:
                if isinstance(lateral_discharges, pd.Series):
                    series = lateral_discharges.loc[lateral.code]

                    # Add to dictionary
                    latdct[lateral.code] = {
                        "branchid": lateral.branch_id,
                        "chainage": str(lateral.branch_offset),
                        "discharge": series,
                    }

                else:
                    if lateral.code not in lateral_discharges.columns:
                        logger.warning(
                            f"No data found for {lateral.code}. Skipping."
                        )
                        continue

                    # Get timeseries
                    series = lateral_discharges.loc[:, lateral.code]

                    # Add to dictionary
                    latdct[lateral.code] = {
                        "branchid": lateral.branch_id,
                        "chainage": str(lateral.branch_offset),
                        "discharge": series,
                    }

    # Add all items
    for key, item in latdct.items():
        self.external_forcings.add_lateral(
            id=key,
            branchid=item["branchid"],
            chainage=item["chainage"],
            discharge=item["discharge"],
        )

StorageNodesIO

Source code in hydrolib/dhydamo/converters/hydamo2df.py
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
class StorageNodesIO:
    def __init__(self, storagenodes):
        self.storagenodes = storagenodes

    def storagenodes_from_datamodel(self, storagenodes):
        """ "From parsed data model of storage nodes"""
        for storagenode_idx, storagenode in storagenodes.iterrows():
            self.add_storagenode(
                id=storagenode.id,
                name=storagenode.name if "name" in storagenode.index else np.nan,
                usestreetstorage=storagenode.usestreetstorage,
                nodetype="unspecified",
                nodeid=storagenode.nodeid,
                usetable="false",
                bedlevel=storagenode.bedlevel,
                area=storagenode.area,
                streetlevel=storagenode.streetlevel,
                streetstoragearea=storagenode.streetstoragearea,
                storagetype=storagenode.storagetype,
            )

    def storagenodes_from_input(self, storagenodes=None, nodeid=None, xy=None, storagedata=None, usestreetstorage= True, network=None):
        """ "From parsed data model of storage nodes"""

        for storagenode_idx, storagenode in storagenodes.iterrows():
            data = storagedata[storagedata.code == storagenode_idx]
            name = storagenode.name
            if pd.isna(name):
                name = storagenode.code
            self.storagenodes.add_storagenode(
                id=storagenode_idx,
                name=name,
                usestreetstorage=usestreetstorage,
                nodetype="unspecified",
                nodeid=nodeid,
                xy=xy,
                branchid=storagenode.branch_id,
                chainage=storagenode.branch_offset,
                usetable="true",
                storagearea=' '.join(data.area.astype(str)),
                levels=' '.join(data.level.astype(str)),
                network=network
            )

storagenodes_from_datamodel(storagenodes)

"From parsed data model of storage nodes

Source code in hydrolib/dhydamo/converters/hydamo2df.py
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
def storagenodes_from_datamodel(self, storagenodes):
    """ "From parsed data model of storage nodes"""
    for storagenode_idx, storagenode in storagenodes.iterrows():
        self.add_storagenode(
            id=storagenode.id,
            name=storagenode.name if "name" in storagenode.index else np.nan,
            usestreetstorage=storagenode.usestreetstorage,
            nodetype="unspecified",
            nodeid=storagenode.nodeid,
            usetable="false",
            bedlevel=storagenode.bedlevel,
            area=storagenode.area,
            streetlevel=storagenode.streetlevel,
            streetstoragearea=storagenode.streetstoragearea,
            storagetype=storagenode.storagetype,
        )

storagenodes_from_input(storagenodes=None, nodeid=None, xy=None, storagedata=None, usestreetstorage=True, network=None)

"From parsed data model of storage nodes

Source code in hydrolib/dhydamo/converters/hydamo2df.py
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
def storagenodes_from_input(self, storagenodes=None, nodeid=None, xy=None, storagedata=None, usestreetstorage= True, network=None):
    """ "From parsed data model of storage nodes"""

    for storagenode_idx, storagenode in storagenodes.iterrows():
        data = storagedata[storagedata.code == storagenode_idx]
        name = storagenode.name
        if pd.isna(name):
            name = storagenode.code
        self.storagenodes.add_storagenode(
            id=storagenode_idx,
            name=name,
            usestreetstorage=usestreetstorage,
            nodetype="unspecified",
            nodeid=nodeid,
            xy=xy,
            branchid=storagenode.branch_id,
            chainage=storagenode.branch_offset,
            usetable="true",
            storagearea=' '.join(data.area.astype(str)),
            levels=' '.join(data.level.astype(str)),
            network=network
        )

StructuresIO

Source code in hydrolib/dhydamo/converters/hydamo2df.py
 427
 428
 429
 430
 431
 432
 433
 434
 435
 436
 437
 438
 439
 440
 441
 442
 443
 444
 445
 446
 447
 448
 449
 450
 451
 452
 453
 454
 455
 456
 457
 458
 459
 460
 461
 462
 463
 464
 465
 466
 467
 468
 469
 470
 471
 472
 473
 474
 475
 476
 477
 478
 479
 480
 481
 482
 483
 484
 485
 486
 487
 488
 489
 490
 491
 492
 493
 494
 495
 496
 497
 498
 499
 500
 501
 502
 503
 504
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 546
 547
 548
 549
 550
 551
 552
 553
 554
 555
 556
 557
 558
 559
 560
 561
 562
 563
 564
 565
 566
 567
 568
 569
 570
 571
 572
 573
 574
 575
 576
 577
 578
 579
 580
 581
 582
 583
 584
 585
 586
 587
 588
 589
 590
 591
 592
 593
 594
 595
 596
 597
 598
 599
 600
 601
 602
 603
 604
 605
 606
 607
 608
 609
 610
 611
 612
 613
 614
 615
 616
 617
 618
 619
 620
 621
 622
 623
 624
 625
 626
 627
 628
 629
 630
 631
 632
 633
 634
 635
 636
 637
 638
 639
 640
 641
 642
 643
 644
 645
 646
 647
 648
 649
 650
 651
 652
 653
 654
 655
 656
 657
 658
 659
 660
 661
 662
 663
 664
 665
 666
 667
 668
 669
 670
 671
 672
 673
 674
 675
 676
 677
 678
 679
 680
 681
 682
 683
 684
 685
 686
 687
 688
 689
 690
 691
 692
 693
 694
 695
 696
 697
 698
 699
 700
 701
 702
 703
 704
 705
 706
 707
 708
 709
 710
 711
 712
 713
 714
 715
 716
 717
 718
 719
 720
 721
 722
 723
 724
 725
 726
 727
 728
 729
 730
 731
 732
 733
 734
 735
 736
 737
 738
 739
 740
 741
 742
 743
 744
 745
 746
 747
 748
 749
 750
 751
 752
 753
 754
 755
 756
 757
 758
 759
 760
 761
 762
 763
 764
 765
 766
 767
 768
 769
 770
 771
 772
 773
 774
 775
 776
 777
 778
 779
 780
 781
 782
 783
 784
 785
 786
 787
 788
 789
 790
 791
 792
 793
 794
 795
 796
 797
 798
 799
 800
 801
 802
 803
 804
 805
 806
 807
 808
 809
 810
 811
 812
 813
 814
 815
 816
 817
 818
 819
 820
 821
 822
 823
 824
 825
 826
 827
 828
 829
 830
 831
 832
 833
 834
 835
 836
 837
 838
 839
 840
 841
 842
 843
 844
 845
 846
 847
 848
 849
 850
 851
 852
 853
 854
 855
 856
 857
 858
 859
 860
 861
 862
 863
 864
 865
 866
 867
 868
 869
 870
 871
 872
 873
 874
 875
 876
 877
 878
 879
 880
 881
 882
 883
 884
 885
 886
 887
 888
 889
 890
 891
 892
 893
 894
 895
 896
 897
 898
 899
 900
 901
 902
 903
 904
 905
 906
 907
 908
 909
 910
 911
 912
 913
 914
 915
 916
 917
 918
 919
 920
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
class StructuresIO:
    def __init__(self, structures):
        self.structures = structures

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def generalstructures_from_datamodel(self, generalstructures: pd.DataFrame) -> None:
        """From parsed data model of orifices

        Args:
            generalstructures (pd.DataFrame): dataframe containing the data
        """

        for generalstructure_idx, generalstructure in generalstructures.iterrows():
            self.structures.add_generalstructure(
                id=generalstructure.id,
                name=generalstructure.name
                if "name" in generalstructure.index
                else np.nan,
                branchid=generalstructure.branch_id,
                chainage=generalstructure.branch_offset,
                allowedflowdir="both",
                upstream1width=generalstructure.upstream1width
                if "upstream1width" in generalstructure.index
                else np.nan,
                upstream1level=generalstructure.upstream1level
                if "upstream1level" in generalstructure.index
                else np.nan,
                upstream2width=generalstructure.upstream2width
                if "upstream2width" in generalstructure.index
                else np.nan,
                upstream2level=generalstructure.upstream2level
                if "upstream2level" in generalstructure.index
                else np.nan,
                crestwidth=generalstructure.crestwidth
                if "crestwidth" in generalstructure.index
                else np.nan,
                crestlevel=generalstructure.crestlevel
                if "crestlevel" in generalstructure.index
                else np.nan,
                crestlength=generalstructure.crestlength
                if "crestlength" in generalstructure.index
                else np.nan,
                downstream1width=generalstructure.downstream1width
                if "downstream1width" in generalstructure.index
                else np.nan,
                downstream1level=generalstructure.downstream1level
                if "downstream1level" in generalstructure.index
                else np.nan,
                downstream2width=generalstructure.downstream2width
                if "downstream2width" in generalstructure.index
                else np.nan,
                downstream2level=generalstructure.downstream2level
                if "downstream2level" in generalstructure.index
                else np.nan,
                gateloweredgelevel=generalstructure.gateloweredgelevel
                if "gateloweredgelevel" in generalstructure.index
                else np.nan,
                posfreegateflowcoeff=generalstructure.posfreegateflowcoeff
                if "posfreegateflowcoeff" in generalstructure.index
                else np.nan,
                posdrowngateflowcoeff=generalstructure.posdrowngateflowcoeff
                if "posdrowngateflowcoeff" in generalstructure.index
                else np.nan,
                posfreeweirflowcoeff=generalstructure.posfreeweirflowcoeff
                if "posfreeweirflowcoeff" in generalstructure.index
                else np.nan,
                posdrownweirflowcoeff=generalstructure.posdrownweirflowcoeff
                if "posdrownweirflowcoeff" in generalstructure.index
                else np.nan,
                poscontrcoeffreegate=generalstructure.poscontrcoeffreegate
                if "poscontrcoeffreegate" in generalstructure.index
                else np.nan,
                negfreegateflowcoeff=generalstructure.negfreegateflowcoeff
                if "negfreegateflowcoeff" in generalstructure.index
                else np.nan,
                negdrowngateflowcoeff=generalstructure.negdrowngateflowcoeff
                if "negdrowngateflowcoeff" in generalstructure.index
                else np.nan,
                negfreeweirflowcoeff=generalstructure.negfreeweirflowcoeff
                if "negfreeweirflowcoeff" in generalstructure.index
                else np.nan,
                negdrownweirflowcoeff=generalstructure.negdrownweirflowcoeff
                if "negdrownweirflowcoeff" in generalstructure.index
                else np.nan,
                negcontrcoeffreegate=generalstructure.negcontrcoeffreegate
                if "negcontrcoeffreegate" in generalstructure.index
                else np.nan,
                extraresistance=generalstructure.extraresistance
                if "extraresistance" in generalstructure.index
                else np.nan,
                gateheight=generalstructure.gateheight
                if "gateheight" in generalstructure.index
                else np.nan,
                gateopeningwidth=generalstructure.gateopeningwidth
                if "gateopeningwidth" in generalstructure.index
                else np.nan,
                gateopeninghorizontaldirection=generalstructure.gateopeninghorizontaldirection
                if "gateopeninghorizontaldirection" in generalstructure.index
                else np.nan,
                usevelocityheight=generalstructure.usevelocityheight
                if "usevelocityheight" in generalstructure.index
                else np.nan,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def weirs(
        self,
        weirs: ExtendedGeoDataFrame = None,
        profile_groups = None,
        profile_lines = None,
        profiles: Optional[ExtendedGeoDataFrame] = None,
        opening: ExtendedDataFrame = None,
        management_device: ExtendedDataFrame = None,
        usevelocityheight: Optional[str] = "true",
    ) -> None:
        """
        Method to convert HyDAMO weirs to DFlowFM structures: regular weirs, orifices and universal weirs.
        If a weir ID corresponds to a 'stuwid' in the profile_groups object, a universal weir is created.
        If the management_device corresponding to a weir has the field "overlaatonderlaat" set to "onderlaat', an orifice is schematized. In all other cases, a regular (rectangular) weir is created.

        Parameters corrspond to the HyDAMO DAMO2.2 objects. DFlowFM keyword 'usevelocityheight' can be specificied as a string, default is 'true'.
        """
        # bypass HyDAMO and add stuwid directly to managment for use in RTC
        if not self.structures.hydamo.management.empty:
            self.structures.hydamo.management['stuwid'] = None

        index = np.zeros((len(weirs.code)))
        if profile_groups is not None and hasattr(profile_groups, "stuwid"):
            index[np.isin(weirs.globalid, np.asarray(profile_groups.stuwid))] = 1

        rweirs = weirs[index == 0]
        for weir in rweirs.itertuples():
            weir_opening = opening[opening.stuwid == weir.globalid]

            # check if a separate name field is present
            if "naam" in weirs:
                name = weir.naam
                if not name:
                    name = weir.code
            else:
                name = weir.code

            if weir_opening.shape[0] > 1:
                print(f'Weir {weir.code} contains {weir_opening.shape[0]} openings. Creating a compound structure with a fictional weir for each one.')
                cmp_list = []
                for num_op, (_, op_row) in enumerate(weir_opening.iterrows()):
                    weir_mandev = management_device[
                        management_device.kunstwerkopeningid
                        == op_row.globalid
                    ]
                    weir_id = f'{weir.code}_{num_op+1}'
                    if (not self.structures.hydamo.management.empty) & (hasattr(self.structures.hydamo.management, 'regelmiddelid')):
                        if weir_mandev.globalid.isin(self.structures.hydamo.management.regelmiddelid).item():
                            idx = self.structures.hydamo.management[self.structures.hydamo.management.regelmiddelid == weir_mandev.globalid.squeeze()].index.values[0]
                            self.structures.hydamo.management.loc[idx, 'stuwid'] =weir_id
                    if weir_mandev.overlaatonderlaat.squeeze().lower() == 'overlaat':
                        cmp_list.append(weir_id)
                        self.structures.add_rweir(id=weir_id,
                                                    name=name,
                                                    branchid=weir.branch_id,
                                                    chainage=weir.branch_offset,
                                                    crestlevel=op_row.laagstedoorstroomhoogte,
                                                    crestwidth=op_row.laagstedoorstroombreedte,
                                                    corrcoeff=weir.afvoercoefficient,
                                                    allowedflowdir="both",
                                                    usevelocityheight=usevelocityheight,
                                                 )
                    elif weir_mandev.overlaatonderlaat.squeeze().lower() == 'onderlaat':
                        cmp_list.append(weir_id)
                        if "maximaaldebiet" not in weir_mandev or pd.isnull(weir_mandev.maximaaldebiet.values[0]):
                            limitflow = "false"
                            maxq = 0.0
                        else:
                            limitflow = "true"
                            maxq = float(weir_mandev.maximaaldebiet.values[0])
                        self.structures.add_orifice(
                            id=weir_id,
                            name=name,
                            branchid=weir.branch_id,
                            chainage=weir.branch_offset,
                            crestlevel=float(weir_opening.laagstedoorstroomhoogte.values[0]),
                            crestwidth=float(weir_opening.laagstedoorstroombreedte.values[0]),
                            corrcoeff=weir.afvoercoefficient,
                            allowedflowdir="both",
                            usevelocityheight=usevelocityheight,
                            gateloweredgelevel=float(weir_opening.laagstedoorstroomhoogte.values[0])
                            + float(weir_mandev.hoogteopening.values[0]),
                            uselimitflowpos=limitflow,
                            limitflowpos=maxq,
                            uselimitflowneg=limitflow,
                            limitflowneg=maxq,
                        )
                    else:
                        print(f'Skipping {weir.code} - from "overlaatonderlaat" {weir_mandev.overlaatonderlaat} the type of structure could not be determined.')
                self.structures.add_compound(id=f'cmp_{weir.code}', structureids =cmp_list)
                # self.structures.rweirs_df.drop(weir.code)
            else:
                if weir_opening.empty:
                    print(f'Skipping {weir.code} because there is no associated opening.')
                    continue

                weir_id = weir.code
                weir_mandev = management_device[
                        management_device.kunstwerkopeningid
                        == weir_opening.globalid.values[0]
                    ]

                if weir_mandev.empty:
                    print(f'Skipping {weir.code} because there is no associated management device.')
                    continue

                if (not self.structures.hydamo.management.empty) & hasattr(self.structures.hydamo.management, 'regelmiddelid'):
                    if weir_mandev.globalid.isin(self.structures.hydamo.management.regelmiddelid).item():
                        idx = self.structures.hydamo.management[self.structures.hydamo.management.regelmiddelid == weir_mandev.globalid.squeeze()].index.values[0]
                        self.structures.hydamo.management.loc[idx, 'stuwid'] = weir_id


                if isinstance(weir_mandev.overlaatonderlaat, pd.Series):
                    overlaatonderlaat = weir_mandev.overlaatonderlaat.squeeze()
                else:
                    overlaatonderlaat = weir_mandev.overlaatonderlaat

                if (
                    overlaatonderlaat.lower()
                    == "overlaat"
                ):
                    self.structures.add_rweir(
                        id=weir_id,
                        name=name,
                        branchid=weir.branch_id,
                        chainage=weir.branch_offset,
                        crestlevel=weir_opening.laagstedoorstroomhoogte.values[0],
                        crestwidth=weir_opening.laagstedoorstroombreedte.values[0],
                        corrcoeff=weir.afvoercoefficient,
                        allowedflowdir="both",
                        usevelocityheight=usevelocityheight,
                    )

                elif (
                    overlaatonderlaat.lower()
                    == "onderlaat"
                ):
                    if "maximaaldebiet" not in weir_mandev or pd.isnull(weir_mandev.maximaaldebiet.values[0]):
                        limitflow = "false"
                        maxq = 0.0
                    else:
                        limitflow = "true"
                        maxq = float(weir_mandev.maximaaldebiet.values[0])
                    self.structures.add_orifice(
                        id=weir_id,
                        name=name,
                        branchid=weir.branch_id,
                        chainage=weir.branch_offset,
                        crestlevel=float(weir_opening.laagstedoorstroomhoogte.values[0]),
                        crestwidth=float(weir_opening.laagstedoorstroombreedte.values[0]),
                        corrcoeff=weir.afvoercoefficient,
                        allowedflowdir="both",
                        usevelocityheight=usevelocityheight,
                        gateloweredgelevel=float(weir_opening.laagstedoorstroomhoogte.values[0])
                        + float(weir_mandev.hoogteopening.values[0]),
                        uselimitflowpos=limitflow,
                        limitflowpos=maxq,
                        uselimitflowneg=limitflow,
                        limitflowneg=maxq,
                    )
                else:
                    print(f'Skipping {weir.code} - from "overlaatonderlaat" {weir_mandev.overlaatonderlaat} the type of structure could not be determined.')

        uweirs = weirs[index == 1]
        for uweir in uweirs.itertuples():
            # check if a separate name field is present
            if "naam" in uweirs:
                name = uweir.naam
            else:
                name = uweir.code

            prof = np.empty(0)
            if (profiles is not None) & ("stuwid" in profile_groups):
                group = profile_groups[profile_groups["stuwid"] == uweir.globalid]
                line = profile_lines[
                    profile_lines["profielgroepid"] == group["globalid"].values[0]
                ]
                prof = profiles[profiles["profiellijnid"] == line["globalid"].values[0]]
                if not prof.empty:
                    counts = len(prof.geometry.iloc[0].coords[:])
                    xyz = np.vstack(prof.geometry.iloc[0].coords[:])
                    length = np.r_[
                        0,
                        np.cumsum(np.hypot(np.diff(xyz[:, 0]), np.diff(xyz[:, 1]))),
                    ]
                    yzvalues = np.c_[length, xyz[:, -1] - np.min(xyz[:, -1])]

            if not hasattr(uweir, 'laagstedoorstroomhoogte') or pd.isnull(uweir.laagstedoorstroomhoogte):
                kruinhoogte = np.min(xyz[:,-1])
            else:
                kruinhoogte = uweir.laagstedoorstroomhoogte

            if len(prof) == 0:
                # return an error it is still not found
                raise ValueError(f"{uweir.code} is not found in any cross-section.")
            self.structures.add_uweir(
                id=uweir.code,
                name=name,
                branchid=uweir.branch_id,
                chainage=uweir.branch_offset,
                crestlevel=kruinhoogte,
                dischargecoeff=uweir.afvoercoefficient,
                allowedflowdir="both",
                numlevels=counts,
                yvalues=" ".join([f"{yz[0]:7.3f}" for yz in yzvalues]),
                zvalues=" ".join([f"{yz[1]:7.3f}" for yz in yzvalues]),
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def weirs_from_datamodel(self, weirs: pd.DataFrame) -> None:
        """ "From parsed data model of weirs"""
        for weir_idx, weir in weirs.iterrows():
            self.structures.add_weir(
                id=weir.id,
                name=weir.name if "name" in weir.index else np.nan,
                branchid=weir.branch_id,
                chainage=weir.branch_offset,
                crestlevel=weir.crestlevel,
                crestwidth=weir.crestwidth,
                corrcoeff=weir.corrcoeff,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def orifices_from_datamodel(self, orifices: pd.DataFrame) -> None:
        """ "From parsed data model of orifices"""
        for orifice_idx, orifice in orifices.iterrows():
            self.structures.add_orifice(
                id=orifice.id,
                name=orifice.name if "name" in orifice.index else np.nan,
                branchid=orifice.branch_id,
                chainage=orifice.branch_offset,
                allowedflowdir="both",
                crestlevel=orifice.crestlevel,
                crestwidth=orifice.crestwidth,
                gateloweredgelevel=orifice.gateloweredgelevel,
                corrcoeff=orifice.corrcoef,
                uselimitflowpos=orifice.uselimitflowpos,
                limitflowpos=orifice.limitflowpos,
                uselimitflowneg=orifice.uselimitflowneg,
                limitflowneg=orifice.limitflowneg,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def uweirs_from_datamodel(self, uweirs: pd.DataFrame) -> None:
        """ "From parsed data model of universal weirs"""
        for uweir_idx, uweir in uweirs.iterrows():
            self.structures.add_uweir(
                id=uweir.id,
                name=uweir.name if "name" in uweir.index else np.nan,
                branchid=uweir.branch_id,
                chainage=uweir.branch_offset,
                crestlevel=uweir.crestlevel,
                yvalues=uweir.yvalues,
                zvalues=uweir.zvalues,
                allowedflowdir="both",
                dischargecoeff=uweir.dischargecoeff,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def bridges(
        self,
        bridges: ExtendedGeoDataFrame,
        profile_groups: ExtendedDataFrame = None,
        profile_lines: ExtendedGeoDataFrame = None,
        profiles: ExtendedGeoDataFrame = None,
    ) -> None:
        """
        Method to convert HyDAMO bridges to DFlowFM bridges. Every bridge needs an associated YZ-profile through profile_groups ('brugid'), profile_lines and profiles.

        Parameters corrspond to the HyDAMO DAMO2.2 objects.
        """
        for bridge in bridges.itertuples():
            # first search in yz-profiles
            group = profile_groups[profile_groups["brugid"] == bridge.globalid]
            line = profile_lines[
                profile_lines["profielgroepid"] == group["globalid"].values[0]
            ]
            prof = profiles[profiles["profiellijnid"] == line["globalid"].values[0]]

            if len(prof) == 0:
                raise ValueError(f"{bridge.code} is not found in any cross-section.")

            if "naam" in bridges:
                name = bridge.naam
            else:
                name = bridge.code

            profile_id = prof.code.values[0]
            self.structures.add_bridge(
                id=bridge.code,
                name=name,
                branchid=bridge.branch_id,
                chainage=bridge.branch_offset,
                csdefid=profile_id,
                shift=0.0,
                allowedflowdir="both",
                inletlosscoeff=bridge.intreeverlies,
                outletlosscoeff=bridge.uittreeverlies,
                length=bridge.lengte,
                frictiontype=bridge.typeruwheid,
                friction=bridge.ruwheid,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def bridges_from_datamodel(self, bridges: pd.DataFrame) -> None:
        """ "From parsed data model of bridges"""
        for bridge_idx, bridge in bridges.iterrows():
            self.structures.add_bridge(
                id=bridge.code,
                name=bridge.name if "name" in bridge.index else np.nan,
                branchid=bridge.branch_id,
                chainage=bridge.branch_offset,
                csdefid=bridge.csdefid,
                shift=0.0,
                allowedflowdir="both",
                inletlosscoeff=bridge.intreeverlies,
                outletlosscoeff=bridge.uittreeverlies,
                length=bridge.lengte,
                frictiontype=bridge.typeruwheid,
                friction=bridge.ruwheid,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def culverts(
        self,
        culverts: ExtendedGeoDataFrame,
        management_device: Optional[ExtendedDataFrame] = None,
    ) -> None:
        """
        Method to convert HyDAMO culverts to DFlowFM culverts. Devices like a valve and a slide can be schematized from the management_device object.
        According to HyDAMO DAMO2.2 a closing_device ('afsluitmiddel') could also be used but this is not supported.

        Parameters corrspond to the HyDAMO DAMO2.2 objects.
        """
        if management_device is not None:
            if 'soortafsluitmiddel' not in management_device.columns:
               management_device['soortafsluitmiddel'] = management_device['soortregelmiddel']

        for culvert in culverts.itertuples():
            # Generate cross section definition name
            if culvert.vormkoker.lower() == "rond" or culvert.vormkoker.lower() == "ellipsvormig":
                crosssection = {"shape": "circle", "diameter": culvert.hoogteopening}
            elif (
                culvert.vormkoker.lower() == "rechthoekig"
                or culvert.vormkoker.lower() == "onbekend"
                or culvert.vormkoker.lower() == "eivormig"
                or culvert.vormkoker.lower() == "muilprofiel"
                or culvert.vormkoker.lower() == "heulprofiel"
            ):
                crosssection = {
                    "shape": "rectangle",
                    "height": culvert.hoogteopening,
                    "width": culvert.breedteopening,
                    "closed": 1,
                }
            else:
                crosssection = {"shape": "circle", "diameter": 0.40}
                print(
                    f"Culvert {culvert.code} has an unknown shape: {culvert.vormkoker}. Applying a default profile (round - 40cm)"
                )

            # check whether an afsluitmiddel is present and take action dependent on its settings
            if management_device is not None:
                mandev = management_device[
                    management_device.duikersifonhevelid == culvert.globalid
                ]
                if 'soortafsluitmiddel' not in mandev:
                    mandev.loc[mandev.index,'soortafsluitmiddel'] = mandev['soortregelmiddel']
            else:
                mandev = pd.DataFrame()

            if mandev.empty:
                allowedflowdir = "both"
                valveonoff = 0
                numlosscoeff = None
                valveopeningheight = 0
                relopening = None
                losscoeff = None
            else:
                for _, i in mandev.iterrows():
                    if i["soortafsluitmiddel"] == "terugslagklep":
                        allowedflowdir = "positive"
                        valveonoff = 0
                        numlosscoeff = None
                        valveopeningheight = 0
                        relopening = None
                        losscoeff = None
                    elif i["soortafsluitmiddel"] == "schuif":
                        allowedflowdir = "positive"
                        valveonoff = 1
                        valveopeningheight = float(i["hoogteopening"])
                        numlosscoeff = 1
                        relopening = [float(i["hoogteopening"]) / culvert.hoogteopening]
                        losscoeff = [float(i["afvoercoefficient"])]
                    else:
                        raise NotImplementedError(
                            f'Type of management device for culvert {culvert.code} is not implemented; only "schuif" and "terugslagklep" are allowed.'
                        )

            # check if a separate name field is present
            if "naam" in culverts:
                name = culvert.naam
            else:
                name = culvert.code

            self.structures.add_culvert(
                id=culvert.code,
                name=name,
                branchid=culvert.branch_id,
                chainage=culvert.branch_offset,
                leftlevel=culvert.hoogtebinnenonderkantbov,
                rightlevel=culvert.hoogtebinnenonderkantbene,
                length=culvert.lengte,
                inletlosscoeff=culvert.intreeverlies,
                outletlosscoeff=culvert.uittreeverlies,
                crosssection=crosssection,
                allowedflowdir=allowedflowdir,
                valveonoff=valveonoff,
                numlosscoeff=numlosscoeff,
                valveopeningheight=valveopeningheight,
                relopening=relopening,
                losscoeff=losscoeff,
                bedfrictiontype=culvert.typeruwheid,
                bedfriction=culvert.ruwheid,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def culverts_from_datamodel(self, culverts: pd.DataFrame) -> None:
        """
        From parsed model of culverts
        """

        # Add to dict
        for culvert_idx, culvert in culverts.iterrows():
            self.structures.add_culvert(
                id=culvert.id,
                name=culvert.name if "name" in culvert.index else np.nan,
                branchid=culvert.branch_id,
                chainage=culvert.branch_offset,
                leftlevel=culvert.leftlevel,
                rightlevel=culvert.rightlevel,
                crosssection=culvert.crosssectiondefinitionid,
                length=culvert.geometry.length
                if "geometry" in culvert.index
                else culvert.length,
                inletlosscoeff=culvert.inletlosscoeff,
                outletlosscoeff=culvert.outletlosscoeff,
                allowedflowdir="both",
                valveonoff=0,
                numlosscoeff=0,
                valveopeningheight=np.nan,
                relopening=np.nan,
                losscoeff=np.nan,
                frictiontype=culvert.frictiontype,
                frictionvalue=culvert.frictionvalue,
            )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def pumps(
        self,
        pumpstations: ExtendedGeoDataFrame,
        pumps: ExtendedDataFrame = None,
        management: ExtendedDataFrame = None,
    ) -> None:
        """
        Method to convert HyDAMO pumps to DFlowFM pumps. Three objects are required: pumpstations, pumps and management ('sturing').

        Parameters corrspond to the HyDAMO DAMO2.2 objects.
        """

        # Add sturing to pumps
        for pumpstation in pumpstations.itertuples():

            # find pumps for gemaal
            pumps_subset = pumps[pumps.gemaalid == pumpstation.globalid]
            if pumps_subset.empty:
                print(f'Skipping {pumpstation.code} because there is no associated pump.')
                continue

            if "naam" in pumpstation:
                name = pumpstation.name
            else:
                name = pumpstation.code
            if pumps_subset.shape[0] > 1:
                # more than one pump
                cmp_list = []
                print(f'Pumpstation {pumpstation.code} contains {pumps_subset.shape[0]} openings. Creating a compound structure with a fictional pump for each one.')
                for ipump, (_,pump) in enumerate(pumps_subset.iterrows()):

                    pump_control = management[management.pompid== pump.globalid]
                    if pump_control.empty:
                        print(f'No management found for {pump.code}')
                        continue

                    startlevelsuctionside = [pump_control["bovengrens"]]
                    stoplevelsuctionside = [pump_control["ondergrens"]]

                    pumpid = f'{pumpstation.code}_{ipump+1}'
                    cmp_list.append(pumpid)

                    self.structures.add_pump(
                        id=pumpid,
                        name=name,
                        branchid=pumpstation.branch_id,
                        chainage=pumpstation.branch_offset,
                        orientation="positive",
                        numstages=1,
                        controlside="suctionside",
                        capacity=pump.maximalecapaciteit/60.,
                        startlevelsuctionside=startlevelsuctionside,
                        stoplevelsuctionside=stoplevelsuctionside,
                        startleveldeliveryside=startlevelsuctionside,
                        stopleveldeliveryside=stoplevelsuctionside,
                    )
                self.structures.add_compound(id=f'cmp_{pumpstation.code}', structureids =cmp_list)

            else:
                #  only one pump
                pump_control = management[management.pompid== pumps_subset.globalid.values[0]]
                if pump_control.empty:
                    print(f'Skipping {pumpstation.code} because there is no associated management.')

                startlevelsuctionside = [pump_control["bovengrens"]]
                stoplevelsuctionside = [pump_control["ondergrens"]]


                # the pumpstation has only one pump
                self.structures.add_pump(
                    id=pumpstation.code,
                    name=name,
                    branchid=pumpstation.branch_id,
                    chainage=pumpstation.branch_offset,
                    orientation="positive",
                    numstages=1,
                    controlside="suctionside",
                    capacity=pumps_subset.maximalecapaciteit.values[0]/60.,
                    startlevelsuctionside=startlevelsuctionside,
                    stoplevelsuctionside=stoplevelsuctionside,
                    startleveldeliveryside=startlevelsuctionside,
                    stopleveldeliveryside=stoplevelsuctionside,
                )

    @validate_arguments(config=dict(arbitrary_types_allowed=True))
    def pumps_from_datamodel(self, pumps: pd.DataFrame) -> None:
        """From parsed data model of pumps"""

        for pump_idx, pump in pumps.iterrows():
            self.structures.add_pump(
                id=pump.id,
                name=pump.name if "name" in pump.index else np.nan,
                branchid=pump.branch_id,
                chainage=pump.branch_offset,
                orientation="positive",
                numstages=1,
                controlside=pump.controlside,
                capacity=pump.maximumcapacity,
                startlevelsuctionside=pump.startlevelsuctionside,
                stoplevelsuctionside=pump.stoplevelsuctionside,
                startleveldeliveryside=pump.startleveldeliveryside,
                stopleveldeliveryside=pump.stopleveldeliveryside,
            )

    @staticmethod
    def move_structure(struc, struc_dict, branch, offset):
        """
        Function the move a structure if needed for a compound structure.

        Parameters
        ----------
        struc : string
            current sub-structure id
        struc_dict : dict
            dict with all structures of a certain type
        branch : string
            branch id of the first structure in the compound
        offset : float
            chainage of the first structure in the compound

        Returns
        -------
        Dict with shifted coordinates.

        """
        branch2 = struc_dict[struc]["branchid"]
        if branch2 != branch:
            logger.warning(
                f"Structures of not on the same branche. Moving structure {struc} to branch {branch}."
            )
        struc_dict[struc]["branchid"] = branch
        struc_dict[struc]["chainage"] = offset
        return struc_dict

    def compound_structures(self, idlist, structurelist):
        # probably the coordinates should all be set to those of the first structure (still to do)
        for c_i, c_id in enumerate(idlist):
            # check the substructure coordinates. If they do not coincide, move subsequent structures to the coordinates of the first
            for s_i, struc in enumerate(structurelist[c_i]):
                if s_i == 0:
                    # find out what type the first structure it is and get its coordinates
                    if not self.structures.pumps_df.empty:
                        if struc in list(self.structures.pumps_df.id):
                            branch = self.structures.pumps_df[
                                self.structures.pumps_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.pumps_df[
                                self.structures.pumps_df.id == struc
                            ].chainage.values[0]
                    if not self.structures.rweirs_df.empty:
                        if struc in list(self.structures.rweirs_df.id):
                            branch = self.structures.rweirs_df[
                                self.structures.rweirs_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.rweirs_df[
                                self.structures.rweirs_df.id == struc
                            ].chainage.values[0]
                    if not self.structures.uweirs_df.empty:
                        if struc in list(self.structures.uweirs_df.id):
                            branch = self.structures.uweirs_df[
                                self.structures.uweirs_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.uweirs_df[
                                self.structures.uweirs_df.id == struc
                            ].chainage.values[0]
                    if not self.structures.culverts_df.empty:
                        if struc in list(self.structures.culverts_df.id):
                            branch = self.structures.culverts_df[
                                self.structures.culverts_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.culverts_df[
                                self.structures.culverts_df.id == struc
                            ].chainage.values[0]
                    if not self.structures.bridges_df.empty:
                        if struc in list(self.structures.bridges_df.id):
                            branch = self.structures.bridges_df[
                                self.structures.bridges_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.bridges_df[
                                self.structures.bridges_df.id == struc
                            ].chainage.values[0]
                    if not self.structures.orifices_df.empty:
                        if struc in list(self.structures.orifices_df.id):
                            branch = self.structures.orifices_df[
                                self.structures.orifices_df.id == struc
                            ].branchid.values[0]
                            offset = self.structures.orifices_df[
                                self.structures.orifices_df.id == struc
                            ].chainage.values[0]
                else:
                    # move a subsequent structure to the location of the first
                    if not self.structures.pumps_df.empty:
                        if struc in list(self.structures.pumps_df.id):
                            self.structures.pumps_df.loc[
                                self.structures.pumps_df[
                                    self.structures.pumps_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.pumps_df.loc[
                                self.structures.pumps_df[
                                    self.structures.pumps_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset
                    if not self.structures.rweirs_df.empty:
                        if struc in list(self.structures.rweirs_df.id):
                            self.structures.rweirs_df.loc[
                                self.structures.rweirs_df[
                                    self.structures.rweirs_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.rweirs_df.loc[
                                self.structures.rweirs_df[
                                    self.structures.rweirs_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset
                    if not self.structures.uweirs_df.empty:
                        if struc in list(self.structures.uweirs_df.id):
                            self.structures.uweirs_df.loc[
                                self.structures.uweirs_df[
                                    self.structures.uweirs_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.uweirs_df.loc[
                                self.structures.uweirs_df[
                                    self.structures.uweirs_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset
                    if not self.structures.culverts_df.empty:
                        if struc in list(self.structures.culverts_df.id):
                            self.structures.culverts_df.loc[
                                self.structures.culverts_df[
                                    self.structures.culverts_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.culverts_df.loc[
                                self.structures.culverts_df[
                                    self.structures.culverts_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset
                    if not self.structures.bridges_df.empty:
                        if struc in list(self.structures.bridges_df.id):
                            self.structures.bridges_df.loc[
                                self.structures.bridges_df[
                                    self.structures.bridges_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.bridges_df.loc[
                                self.structures.bridges_df[
                                    self.structures.bridges_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset
                    if not self.structures.orifices_df.empty:
                        if struc in list(self.structures.orifices_df.id):
                            self.structures.orifices_df.loc[
                                self.structures.orifices_df[
                                    self.structures.orifices_df.id == struc
                                ].index,
                                "branchid",
                            ] = branch
                            self.structures.orifices_df.loc[
                                self.structures.orifices_df[
                                    self.structures.orifices_df.id == struc
                                ].index,
                                "chainage",
                            ] = offset

            self.structures.add_compound(id=c_id, structureids=structurelist[c_i])

bridges(bridges: ExtendedGeoDataFrame, profile_groups: ExtendedDataFrame = None, profile_lines: ExtendedGeoDataFrame = None, profiles: ExtendedGeoDataFrame = None) -> None

Method to convert HyDAMO bridges to DFlowFM bridges. Every bridge needs an associated YZ-profile through profile_groups ('brugid'), profile_lines and profiles.

Parameters corrspond to the HyDAMO DAMO2.2 objects.

Source code in hydrolib/dhydamo/converters/hydamo2df.py
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def bridges(
    self,
    bridges: ExtendedGeoDataFrame,
    profile_groups: ExtendedDataFrame = None,
    profile_lines: ExtendedGeoDataFrame = None,
    profiles: ExtendedGeoDataFrame = None,
) -> None:
    """
    Method to convert HyDAMO bridges to DFlowFM bridges. Every bridge needs an associated YZ-profile through profile_groups ('brugid'), profile_lines and profiles.

    Parameters corrspond to the HyDAMO DAMO2.2 objects.
    """
    for bridge in bridges.itertuples():
        # first search in yz-profiles
        group = profile_groups[profile_groups["brugid"] == bridge.globalid]
        line = profile_lines[
            profile_lines["profielgroepid"] == group["globalid"].values[0]
        ]
        prof = profiles[profiles["profiellijnid"] == line["globalid"].values[0]]

        if len(prof) == 0:
            raise ValueError(f"{bridge.code} is not found in any cross-section.")

        if "naam" in bridges:
            name = bridge.naam
        else:
            name = bridge.code

        profile_id = prof.code.values[0]
        self.structures.add_bridge(
            id=bridge.code,
            name=name,
            branchid=bridge.branch_id,
            chainage=bridge.branch_offset,
            csdefid=profile_id,
            shift=0.0,
            allowedflowdir="both",
            inletlosscoeff=bridge.intreeverlies,
            outletlosscoeff=bridge.uittreeverlies,
            length=bridge.lengte,
            frictiontype=bridge.typeruwheid,
            friction=bridge.ruwheid,
        )

bridges_from_datamodel(bridges: pd.DataFrame) -> None

"From parsed data model of bridges

Source code in hydrolib/dhydamo/converters/hydamo2df.py
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def bridges_from_datamodel(self, bridges: pd.DataFrame) -> None:
    """ "From parsed data model of bridges"""
    for bridge_idx, bridge in bridges.iterrows():
        self.structures.add_bridge(
            id=bridge.code,
            name=bridge.name if "name" in bridge.index else np.nan,
            branchid=bridge.branch_id,
            chainage=bridge.branch_offset,
            csdefid=bridge.csdefid,
            shift=0.0,
            allowedflowdir="both",
            inletlosscoeff=bridge.intreeverlies,
            outletlosscoeff=bridge.uittreeverlies,
            length=bridge.lengte,
            frictiontype=bridge.typeruwheid,
            friction=bridge.ruwheid,
        )

culverts(culverts: ExtendedGeoDataFrame, management_device: Optional[ExtendedDataFrame] = None) -> None

Method to convert HyDAMO culverts to DFlowFM culverts. Devices like a valve and a slide can be schematized from the management_device object. According to HyDAMO DAMO2.2 a closing_device ('afsluitmiddel') could also be used but this is not supported.

Parameters corrspond to the HyDAMO DAMO2.2 objects.

Source code in hydrolib/dhydamo/converters/hydamo2df.py
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def culverts(
    self,
    culverts: ExtendedGeoDataFrame,
    management_device: Optional[ExtendedDataFrame] = None,
) -> None:
    """
    Method to convert HyDAMO culverts to DFlowFM culverts. Devices like a valve and a slide can be schematized from the management_device object.
    According to HyDAMO DAMO2.2 a closing_device ('afsluitmiddel') could also be used but this is not supported.

    Parameters corrspond to the HyDAMO DAMO2.2 objects.
    """
    if management_device is not None:
        if 'soortafsluitmiddel' not in management_device.columns:
           management_device['soortafsluitmiddel'] = management_device['soortregelmiddel']

    for culvert in culverts.itertuples():
        # Generate cross section definition name
        if culvert.vormkoker.lower() == "rond" or culvert.vormkoker.lower() == "ellipsvormig":
            crosssection = {"shape": "circle", "diameter": culvert.hoogteopening}
        elif (
            culvert.vormkoker.lower() == "rechthoekig"
            or culvert.vormkoker.lower() == "onbekend"
            or culvert.vormkoker.lower() == "eivormig"
            or culvert.vormkoker.lower() == "muilprofiel"
            or culvert.vormkoker.lower() == "heulprofiel"
        ):
            crosssection = {
                "shape": "rectangle",
                "height": culvert.hoogteopening,
                "width": culvert.breedteopening,
                "closed": 1,
            }
        else:
            crosssection = {"shape": "circle", "diameter": 0.40}
            print(
                f"Culvert {culvert.code} has an unknown shape: {culvert.vormkoker}. Applying a default profile (round - 40cm)"
            )

        # check whether an afsluitmiddel is present and take action dependent on its settings
        if management_device is not None:
            mandev = management_device[
                management_device.duikersifonhevelid == culvert.globalid
            ]
            if 'soortafsluitmiddel' not in mandev:
                mandev.loc[mandev.index,'soortafsluitmiddel'] = mandev['soortregelmiddel']
        else:
            mandev = pd.DataFrame()

        if mandev.empty:
            allowedflowdir = "both"
            valveonoff = 0
            numlosscoeff = None
            valveopeningheight = 0
            relopening = None
            losscoeff = None
        else:
            for _, i in mandev.iterrows():
                if i["soortafsluitmiddel"] == "terugslagklep":
                    allowedflowdir = "positive"
                    valveonoff = 0
                    numlosscoeff = None
                    valveopeningheight = 0
                    relopening = None
                    losscoeff = None
                elif i["soortafsluitmiddel"] == "schuif":
                    allowedflowdir = "positive"
                    valveonoff = 1
                    valveopeningheight = float(i["hoogteopening"])
                    numlosscoeff = 1
                    relopening = [float(i["hoogteopening"]) / culvert.hoogteopening]
                    losscoeff = [float(i["afvoercoefficient"])]
                else:
                    raise NotImplementedError(
                        f'Type of management device for culvert {culvert.code} is not implemented; only "schuif" and "terugslagklep" are allowed.'
                    )

        # check if a separate name field is present
        if "naam" in culverts:
            name = culvert.naam
        else:
            name = culvert.code

        self.structures.add_culvert(
            id=culvert.code,
            name=name,
            branchid=culvert.branch_id,
            chainage=culvert.branch_offset,
            leftlevel=culvert.hoogtebinnenonderkantbov,
            rightlevel=culvert.hoogtebinnenonderkantbene,
            length=culvert.lengte,
            inletlosscoeff=culvert.intreeverlies,
            outletlosscoeff=culvert.uittreeverlies,
            crosssection=crosssection,
            allowedflowdir=allowedflowdir,
            valveonoff=valveonoff,
            numlosscoeff=numlosscoeff,
            valveopeningheight=valveopeningheight,
            relopening=relopening,
            losscoeff=losscoeff,
            bedfrictiontype=culvert.typeruwheid,
            bedfriction=culvert.ruwheid,
        )

culverts_from_datamodel(culverts: pd.DataFrame) -> None

From parsed model of culverts

Source code in hydrolib/dhydamo/converters/hydamo2df.py
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def culverts_from_datamodel(self, culverts: pd.DataFrame) -> None:
    """
    From parsed model of culverts
    """

    # Add to dict
    for culvert_idx, culvert in culverts.iterrows():
        self.structures.add_culvert(
            id=culvert.id,
            name=culvert.name if "name" in culvert.index else np.nan,
            branchid=culvert.branch_id,
            chainage=culvert.branch_offset,
            leftlevel=culvert.leftlevel,
            rightlevel=culvert.rightlevel,
            crosssection=culvert.crosssectiondefinitionid,
            length=culvert.geometry.length
            if "geometry" in culvert.index
            else culvert.length,
            inletlosscoeff=culvert.inletlosscoeff,
            outletlosscoeff=culvert.outletlosscoeff,
            allowedflowdir="both",
            valveonoff=0,
            numlosscoeff=0,
            valveopeningheight=np.nan,
            relopening=np.nan,
            losscoeff=np.nan,
            frictiontype=culvert.frictiontype,
            frictionvalue=culvert.frictionvalue,
        )

generalstructures_from_datamodel(generalstructures: pd.DataFrame) -> None

From parsed data model of orifices

Parameters:

Name Type Description Default
generalstructures DataFrame

dataframe containing the data

required
Source code in hydrolib/dhydamo/converters/hydamo2df.py
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def generalstructures_from_datamodel(self, generalstructures: pd.DataFrame) -> None:
    """From parsed data model of orifices

    Args:
        generalstructures (pd.DataFrame): dataframe containing the data
    """

    for generalstructure_idx, generalstructure in generalstructures.iterrows():
        self.structures.add_generalstructure(
            id=generalstructure.id,
            name=generalstructure.name
            if "name" in generalstructure.index
            else np.nan,
            branchid=generalstructure.branch_id,
            chainage=generalstructure.branch_offset,
            allowedflowdir="both",
            upstream1width=generalstructure.upstream1width
            if "upstream1width" in generalstructure.index
            else np.nan,
            upstream1level=generalstructure.upstream1level
            if "upstream1level" in generalstructure.index
            else np.nan,
            upstream2width=generalstructure.upstream2width
            if "upstream2width" in generalstructure.index
            else np.nan,
            upstream2level=generalstructure.upstream2level
            if "upstream2level" in generalstructure.index
            else np.nan,
            crestwidth=generalstructure.crestwidth
            if "crestwidth" in generalstructure.index
            else np.nan,
            crestlevel=generalstructure.crestlevel
            if "crestlevel" in generalstructure.index
            else np.nan,
            crestlength=generalstructure.crestlength
            if "crestlength" in generalstructure.index
            else np.nan,
            downstream1width=generalstructure.downstream1width
            if "downstream1width" in generalstructure.index
            else np.nan,
            downstream1level=generalstructure.downstream1level
            if "downstream1level" in generalstructure.index
            else np.nan,
            downstream2width=generalstructure.downstream2width
            if "downstream2width" in generalstructure.index
            else np.nan,
            downstream2level=generalstructure.downstream2level
            if "downstream2level" in generalstructure.index
            else np.nan,
            gateloweredgelevel=generalstructure.gateloweredgelevel
            if "gateloweredgelevel" in generalstructure.index
            else np.nan,
            posfreegateflowcoeff=generalstructure.posfreegateflowcoeff
            if "posfreegateflowcoeff" in generalstructure.index
            else np.nan,
            posdrowngateflowcoeff=generalstructure.posdrowngateflowcoeff
            if "posdrowngateflowcoeff" in generalstructure.index
            else np.nan,
            posfreeweirflowcoeff=generalstructure.posfreeweirflowcoeff
            if "posfreeweirflowcoeff" in generalstructure.index
            else np.nan,
            posdrownweirflowcoeff=generalstructure.posdrownweirflowcoeff
            if "posdrownweirflowcoeff" in generalstructure.index
            else np.nan,
            poscontrcoeffreegate=generalstructure.poscontrcoeffreegate
            if "poscontrcoeffreegate" in generalstructure.index
            else np.nan,
            negfreegateflowcoeff=generalstructure.negfreegateflowcoeff
            if "negfreegateflowcoeff" in generalstructure.index
            else np.nan,
            negdrowngateflowcoeff=generalstructure.negdrowngateflowcoeff
            if "negdrowngateflowcoeff" in generalstructure.index
            else np.nan,
            negfreeweirflowcoeff=generalstructure.negfreeweirflowcoeff
            if "negfreeweirflowcoeff" in generalstructure.index
            else np.nan,
            negdrownweirflowcoeff=generalstructure.negdrownweirflowcoeff
            if "negdrownweirflowcoeff" in generalstructure.index
            else np.nan,
            negcontrcoeffreegate=generalstructure.negcontrcoeffreegate
            if "negcontrcoeffreegate" in generalstructure.index
            else np.nan,
            extraresistance=generalstructure.extraresistance
            if "extraresistance" in generalstructure.index
            else np.nan,
            gateheight=generalstructure.gateheight
            if "gateheight" in generalstructure.index
            else np.nan,
            gateopeningwidth=generalstructure.gateopeningwidth
            if "gateopeningwidth" in generalstructure.index
            else np.nan,
            gateopeninghorizontaldirection=generalstructure.gateopeninghorizontaldirection
            if "gateopeninghorizontaldirection" in generalstructure.index
            else np.nan,
            usevelocityheight=generalstructure.usevelocityheight
            if "usevelocityheight" in generalstructure.index
            else np.nan,
        )

move_structure(struc, struc_dict, branch, offset) staticmethod

Function the move a structure if needed for a compound structure.

Parameters

struc : string current sub-structure id struc_dict : dict dict with all structures of a certain type branch : string branch id of the first structure in the compound offset : float chainage of the first structure in the compound

Returns

Dict with shifted coordinates.

Source code in hydrolib/dhydamo/converters/hydamo2df.py
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
@staticmethod
def move_structure(struc, struc_dict, branch, offset):
    """
    Function the move a structure if needed for a compound structure.

    Parameters
    ----------
    struc : string
        current sub-structure id
    struc_dict : dict
        dict with all structures of a certain type
    branch : string
        branch id of the first structure in the compound
    offset : float
        chainage of the first structure in the compound

    Returns
    -------
    Dict with shifted coordinates.

    """
    branch2 = struc_dict[struc]["branchid"]
    if branch2 != branch:
        logger.warning(
            f"Structures of not on the same branche. Moving structure {struc} to branch {branch}."
        )
    struc_dict[struc]["branchid"] = branch
    struc_dict[struc]["chainage"] = offset
    return struc_dict

orifices_from_datamodel(orifices: pd.DataFrame) -> None

"From parsed data model of orifices

Source code in hydrolib/dhydamo/converters/hydamo2df.py
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def orifices_from_datamodel(self, orifices: pd.DataFrame) -> None:
    """ "From parsed data model of orifices"""
    for orifice_idx, orifice in orifices.iterrows():
        self.structures.add_orifice(
            id=orifice.id,
            name=orifice.name if "name" in orifice.index else np.nan,
            branchid=orifice.branch_id,
            chainage=orifice.branch_offset,
            allowedflowdir="both",
            crestlevel=orifice.crestlevel,
            crestwidth=orifice.crestwidth,
            gateloweredgelevel=orifice.gateloweredgelevel,
            corrcoeff=orifice.corrcoef,
            uselimitflowpos=orifice.uselimitflowpos,
            limitflowpos=orifice.limitflowpos,
            uselimitflowneg=orifice.uselimitflowneg,
            limitflowneg=orifice.limitflowneg,
        )

pumps(pumpstations: ExtendedGeoDataFrame, pumps: ExtendedDataFrame = None, management: ExtendedDataFrame = None) -> None

Method to convert HyDAMO pumps to DFlowFM pumps. Three objects are required: pumpstations, pumps and management ('sturing').

Parameters corrspond to the HyDAMO DAMO2.2 objects.

Source code in hydrolib/dhydamo/converters/hydamo2df.py
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def pumps(
    self,
    pumpstations: ExtendedGeoDataFrame,
    pumps: ExtendedDataFrame = None,
    management: ExtendedDataFrame = None,
) -> None:
    """
    Method to convert HyDAMO pumps to DFlowFM pumps. Three objects are required: pumpstations, pumps and management ('sturing').

    Parameters corrspond to the HyDAMO DAMO2.2 objects.
    """

    # Add sturing to pumps
    for pumpstation in pumpstations.itertuples():

        # find pumps for gemaal
        pumps_subset = pumps[pumps.gemaalid == pumpstation.globalid]
        if pumps_subset.empty:
            print(f'Skipping {pumpstation.code} because there is no associated pump.')
            continue

        if "naam" in pumpstation:
            name = pumpstation.name
        else:
            name = pumpstation.code
        if pumps_subset.shape[0] > 1:
            # more than one pump
            cmp_list = []
            print(f'Pumpstation {pumpstation.code} contains {pumps_subset.shape[0]} openings. Creating a compound structure with a fictional pump for each one.')
            for ipump, (_,pump) in enumerate(pumps_subset.iterrows()):

                pump_control = management[management.pompid== pump.globalid]
                if pump_control.empty:
                    print(f'No management found for {pump.code}')
                    continue

                startlevelsuctionside = [pump_control["bovengrens"]]
                stoplevelsuctionside = [pump_control["ondergrens"]]

                pumpid = f'{pumpstation.code}_{ipump+1}'
                cmp_list.append(pumpid)

                self.structures.add_pump(
                    id=pumpid,
                    name=name,
                    branchid=pumpstation.branch_id,
                    chainage=pumpstation.branch_offset,
                    orientation="positive",
                    numstages=1,
                    controlside="suctionside",
                    capacity=pump.maximalecapaciteit/60.,
                    startlevelsuctionside=startlevelsuctionside,
                    stoplevelsuctionside=stoplevelsuctionside,
                    startleveldeliveryside=startlevelsuctionside,
                    stopleveldeliveryside=stoplevelsuctionside,
                )
            self.structures.add_compound(id=f'cmp_{pumpstation.code}', structureids =cmp_list)

        else:
            #  only one pump
            pump_control = management[management.pompid== pumps_subset.globalid.values[0]]
            if pump_control.empty:
                print(f'Skipping {pumpstation.code} because there is no associated management.')

            startlevelsuctionside = [pump_control["bovengrens"]]
            stoplevelsuctionside = [pump_control["ondergrens"]]


            # the pumpstation has only one pump
            self.structures.add_pump(
                id=pumpstation.code,
                name=name,
                branchid=pumpstation.branch_id,
                chainage=pumpstation.branch_offset,
                orientation="positive",
                numstages=1,
                controlside="suctionside",
                capacity=pumps_subset.maximalecapaciteit.values[0]/60.,
                startlevelsuctionside=startlevelsuctionside,
                stoplevelsuctionside=stoplevelsuctionside,
                startleveldeliveryside=startlevelsuctionside,
                stopleveldeliveryside=stoplevelsuctionside,
            )

pumps_from_datamodel(pumps: pd.DataFrame) -> None

From parsed data model of pumps

Source code in hydrolib/dhydamo/converters/hydamo2df.py
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def pumps_from_datamodel(self, pumps: pd.DataFrame) -> None:
    """From parsed data model of pumps"""

    for pump_idx, pump in pumps.iterrows():
        self.structures.add_pump(
            id=pump.id,
            name=pump.name if "name" in pump.index else np.nan,
            branchid=pump.branch_id,
            chainage=pump.branch_offset,
            orientation="positive",
            numstages=1,
            controlside=pump.controlside,
            capacity=pump.maximumcapacity,
            startlevelsuctionside=pump.startlevelsuctionside,
            stoplevelsuctionside=pump.stoplevelsuctionside,
            startleveldeliveryside=pump.startleveldeliveryside,
            stopleveldeliveryside=pump.stopleveldeliveryside,
        )

uweirs_from_datamodel(uweirs: pd.DataFrame) -> None

"From parsed data model of universal weirs

Source code in hydrolib/dhydamo/converters/hydamo2df.py
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def uweirs_from_datamodel(self, uweirs: pd.DataFrame) -> None:
    """ "From parsed data model of universal weirs"""
    for uweir_idx, uweir in uweirs.iterrows():
        self.structures.add_uweir(
            id=uweir.id,
            name=uweir.name if "name" in uweir.index else np.nan,
            branchid=uweir.branch_id,
            chainage=uweir.branch_offset,
            crestlevel=uweir.crestlevel,
            yvalues=uweir.yvalues,
            zvalues=uweir.zvalues,
            allowedflowdir="both",
            dischargecoeff=uweir.dischargecoeff,
        )

weirs(weirs: ExtendedGeoDataFrame = None, profile_groups=None, profile_lines=None, profiles: Optional[ExtendedGeoDataFrame] = None, opening: ExtendedDataFrame = None, management_device: ExtendedDataFrame = None, usevelocityheight: Optional[str] = 'true') -> None

Method to convert HyDAMO weirs to DFlowFM structures: regular weirs, orifices and universal weirs. If a weir ID corresponds to a 'stuwid' in the profile_groups object, a universal weir is created. If the management_device corresponding to a weir has the field "overlaatonderlaat" set to "onderlaat', an orifice is schematized. In all other cases, a regular (rectangular) weir is created.

Parameters corrspond to the HyDAMO DAMO2.2 objects. DFlowFM keyword 'usevelocityheight' can be specificied as a string, default is 'true'.

Source code in hydrolib/dhydamo/converters/hydamo2df.py
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def weirs(
    self,
    weirs: ExtendedGeoDataFrame = None,
    profile_groups = None,
    profile_lines = None,
    profiles: Optional[ExtendedGeoDataFrame] = None,
    opening: ExtendedDataFrame = None,
    management_device: ExtendedDataFrame = None,
    usevelocityheight: Optional[str] = "true",
) -> None:
    """
    Method to convert HyDAMO weirs to DFlowFM structures: regular weirs, orifices and universal weirs.
    If a weir ID corresponds to a 'stuwid' in the profile_groups object, a universal weir is created.
    If the management_device corresponding to a weir has the field "overlaatonderlaat" set to "onderlaat', an orifice is schematized. In all other cases, a regular (rectangular) weir is created.

    Parameters corrspond to the HyDAMO DAMO2.2 objects. DFlowFM keyword 'usevelocityheight' can be specificied as a string, default is 'true'.
    """
    # bypass HyDAMO and add stuwid directly to managment for use in RTC
    if not self.structures.hydamo.management.empty:
        self.structures.hydamo.management['stuwid'] = None

    index = np.zeros((len(weirs.code)))
    if profile_groups is not None and hasattr(profile_groups, "stuwid"):
        index[np.isin(weirs.globalid, np.asarray(profile_groups.stuwid))] = 1

    rweirs = weirs[index == 0]
    for weir in rweirs.itertuples():
        weir_opening = opening[opening.stuwid == weir.globalid]

        # check if a separate name field is present
        if "naam" in weirs:
            name = weir.naam
            if not name:
                name = weir.code
        else:
            name = weir.code

        if weir_opening.shape[0] > 1:
            print(f'Weir {weir.code} contains {weir_opening.shape[0]} openings. Creating a compound structure with a fictional weir for each one.')
            cmp_list = []
            for num_op, (_, op_row) in enumerate(weir_opening.iterrows()):
                weir_mandev = management_device[
                    management_device.kunstwerkopeningid
                    == op_row.globalid
                ]
                weir_id = f'{weir.code}_{num_op+1}'
                if (not self.structures.hydamo.management.empty) & (hasattr(self.structures.hydamo.management, 'regelmiddelid')):
                    if weir_mandev.globalid.isin(self.structures.hydamo.management.regelmiddelid).item():
                        idx = self.structures.hydamo.management[self.structures.hydamo.management.regelmiddelid == weir_mandev.globalid.squeeze()].index.values[0]
                        self.structures.hydamo.management.loc[idx, 'stuwid'] =weir_id
                if weir_mandev.overlaatonderlaat.squeeze().lower() == 'overlaat':
                    cmp_list.append(weir_id)
                    self.structures.add_rweir(id=weir_id,
                                                name=name,
                                                branchid=weir.branch_id,
                                                chainage=weir.branch_offset,
                                                crestlevel=op_row.laagstedoorstroomhoogte,
                                                crestwidth=op_row.laagstedoorstroombreedte,
                                                corrcoeff=weir.afvoercoefficient,
                                                allowedflowdir="both",
                                                usevelocityheight=usevelocityheight,
                                             )
                elif weir_mandev.overlaatonderlaat.squeeze().lower() == 'onderlaat':
                    cmp_list.append(weir_id)
                    if "maximaaldebiet" not in weir_mandev or pd.isnull(weir_mandev.maximaaldebiet.values[0]):
                        limitflow = "false"
                        maxq = 0.0
                    else:
                        limitflow = "true"
                        maxq = float(weir_mandev.maximaaldebiet.values[0])
                    self.structures.add_orifice(
                        id=weir_id,
                        name=name,
                        branchid=weir.branch_id,
                        chainage=weir.branch_offset,
                        crestlevel=float(weir_opening.laagstedoorstroomhoogte.values[0]),
                        crestwidth=float(weir_opening.laagstedoorstroombreedte.values[0]),
                        corrcoeff=weir.afvoercoefficient,
                        allowedflowdir="both",
                        usevelocityheight=usevelocityheight,
                        gateloweredgelevel=float(weir_opening.laagstedoorstroomhoogte.values[0])
                        + float(weir_mandev.hoogteopening.values[0]),
                        uselimitflowpos=limitflow,
                        limitflowpos=maxq,
                        uselimitflowneg=limitflow,
                        limitflowneg=maxq,
                    )
                else:
                    print(f'Skipping {weir.code} - from "overlaatonderlaat" {weir_mandev.overlaatonderlaat} the type of structure could not be determined.')
            self.structures.add_compound(id=f'cmp_{weir.code}', structureids =cmp_list)
            # self.structures.rweirs_df.drop(weir.code)
        else:
            if weir_opening.empty:
                print(f'Skipping {weir.code} because there is no associated opening.')
                continue

            weir_id = weir.code
            weir_mandev = management_device[
                    management_device.kunstwerkopeningid
                    == weir_opening.globalid.values[0]
                ]

            if weir_mandev.empty:
                print(f'Skipping {weir.code} because there is no associated management device.')
                continue

            if (not self.structures.hydamo.management.empty) & hasattr(self.structures.hydamo.management, 'regelmiddelid'):
                if weir_mandev.globalid.isin(self.structures.hydamo.management.regelmiddelid).item():
                    idx = self.structures.hydamo.management[self.structures.hydamo.management.regelmiddelid == weir_mandev.globalid.squeeze()].index.values[0]
                    self.structures.hydamo.management.loc[idx, 'stuwid'] = weir_id


            if isinstance(weir_mandev.overlaatonderlaat, pd.Series):
                overlaatonderlaat = weir_mandev.overlaatonderlaat.squeeze()
            else:
                overlaatonderlaat = weir_mandev.overlaatonderlaat

            if (
                overlaatonderlaat.lower()
                == "overlaat"
            ):
                self.structures.add_rweir(
                    id=weir_id,
                    name=name,
                    branchid=weir.branch_id,
                    chainage=weir.branch_offset,
                    crestlevel=weir_opening.laagstedoorstroomhoogte.values[0],
                    crestwidth=weir_opening.laagstedoorstroombreedte.values[0],
                    corrcoeff=weir.afvoercoefficient,
                    allowedflowdir="both",
                    usevelocityheight=usevelocityheight,
                )

            elif (
                overlaatonderlaat.lower()
                == "onderlaat"
            ):
                if "maximaaldebiet" not in weir_mandev or pd.isnull(weir_mandev.maximaaldebiet.values[0]):
                    limitflow = "false"
                    maxq = 0.0
                else:
                    limitflow = "true"
                    maxq = float(weir_mandev.maximaaldebiet.values[0])
                self.structures.add_orifice(
                    id=weir_id,
                    name=name,
                    branchid=weir.branch_id,
                    chainage=weir.branch_offset,
                    crestlevel=float(weir_opening.laagstedoorstroomhoogte.values[0]),
                    crestwidth=float(weir_opening.laagstedoorstroombreedte.values[0]),
                    corrcoeff=weir.afvoercoefficient,
                    allowedflowdir="both",
                    usevelocityheight=usevelocityheight,
                    gateloweredgelevel=float(weir_opening.laagstedoorstroomhoogte.values[0])
                    + float(weir_mandev.hoogteopening.values[0]),
                    uselimitflowpos=limitflow,
                    limitflowpos=maxq,
                    uselimitflowneg=limitflow,
                    limitflowneg=maxq,
                )
            else:
                print(f'Skipping {weir.code} - from "overlaatonderlaat" {weir_mandev.overlaatonderlaat} the type of structure could not be determined.')

    uweirs = weirs[index == 1]
    for uweir in uweirs.itertuples():
        # check if a separate name field is present
        if "naam" in uweirs:
            name = uweir.naam
        else:
            name = uweir.code

        prof = np.empty(0)
        if (profiles is not None) & ("stuwid" in profile_groups):
            group = profile_groups[profile_groups["stuwid"] == uweir.globalid]
            line = profile_lines[
                profile_lines["profielgroepid"] == group["globalid"].values[0]
            ]
            prof = profiles[profiles["profiellijnid"] == line["globalid"].values[0]]
            if not prof.empty:
                counts = len(prof.geometry.iloc[0].coords[:])
                xyz = np.vstack(prof.geometry.iloc[0].coords[:])
                length = np.r_[
                    0,
                    np.cumsum(np.hypot(np.diff(xyz[:, 0]), np.diff(xyz[:, 1]))),
                ]
                yzvalues = np.c_[length, xyz[:, -1] - np.min(xyz[:, -1])]

        if not hasattr(uweir, 'laagstedoorstroomhoogte') or pd.isnull(uweir.laagstedoorstroomhoogte):
            kruinhoogte = np.min(xyz[:,-1])
        else:
            kruinhoogte = uweir.laagstedoorstroomhoogte

        if len(prof) == 0:
            # return an error it is still not found
            raise ValueError(f"{uweir.code} is not found in any cross-section.")
        self.structures.add_uweir(
            id=uweir.code,
            name=name,
            branchid=uweir.branch_id,
            chainage=uweir.branch_offset,
            crestlevel=kruinhoogte,
            dischargecoeff=uweir.afvoercoefficient,
            allowedflowdir="both",
            numlevels=counts,
            yvalues=" ".join([f"{yz[0]:7.3f}" for yz in yzvalues]),
            zvalues=" ".join([f"{yz[1]:7.3f}" for yz in yzvalues]),
        )

weirs_from_datamodel(weirs: pd.DataFrame) -> None

"From parsed data model of weirs

Source code in hydrolib/dhydamo/converters/hydamo2df.py
740
741
742
743
744
745
746
747
748
749
750
751
752
@validate_arguments(config=dict(arbitrary_types_allowed=True))
def weirs_from_datamodel(self, weirs: pd.DataFrame) -> None:
    """ "From parsed data model of weirs"""
    for weir_idx, weir in weirs.iterrows():
        self.structures.add_weir(
            id=weir.id,
            name=weir.name if "name" in weir.index else np.nan,
            branchid=weir.branch_id,
            chainage=weir.branch_offset,
            crestlevel=weir.crestlevel,
            crestwidth=weir.crestwidth,
            corrcoeff=weir.corrcoeff,
        )