Modules#
dfm_tools post-processing#
- dfm_tools.get_nc.get_Dataset_atdepths(data_xr: UgridDataset, depths, reference: str = 'z0')[source]#
Lazily depth-slice a dataset with layers. Performance can be increased by using a subset of variables or subsetting the dataset in any dimension. This can be done for instance with ds.isel(time=-1) or uds.ugrid.sel(x=slice(),y=slice()) to subset a ugrid dataset in space. The return dataset only contains the sliced variables.
Parameters#
- data_xrxu.UgridDataset
has to be Dataset (not a DataArray), otherwise mesh2d_flowelem_zw etc are not available (interface z values) in case of zsigma/sigma layers (or fullgrid), it is advisable to .sel()/.isel() the time dimension first, because that is less computationally heavy
- depthsTYPE
int/float or list/array of int/float. Depths w.r.t. reference level. If reference==’waterlevel’, depth>0 returns only nans. If reference==’bedlevel’, depth<0 returns only nans. Depths are sorted and only uniques are kept.
- referencestr, optional
compute depth w.r.t. z0/waterlevel/bed. The default is ‘z0’.
- zlayer_z0_selnearestbool, optional
Use xr.interp() to interpolate zlayer model to z-value. Only possible for reference=’z’ (not ‘waterlevel’ or ‘bedlevel’). Only used if “mesh2d_layer_z” is present (zlayer model) This is faster but results in values interpolated between zcc (z cell centers), so it is different than slicing.. The default is False.
Raises#
- Exception
DESCRIPTION.
Returns#
- xu.UgridDataset
Dataset with the depth-sliced variables.
- dfm_tools.get_nc.plot_ztdata(data_xr_sel, varname, ax=None, only_contour=False, **kwargs)[source]#
Parameters#
- data_xrTYPE
DESCRIPTION.
- varnameTYPE
DESCRIPTION.
- axmatplotlib.axes._subplots.AxesSubplot, optional
the figure axis. The default is None.
- only_contourbool, optional
Wheter to plot contour lines of the dataset. The default is False.
- kwargsTYPE
properties to give on to the pcolormesh function.
Raises#
- Exception
DESCRIPTION.
Returns#
- pcmatplotlib.collections.QuadMesh
DESCRIPTION.
- dfm_tools.get_nc.polyline_mapslice(uds: UgridDataset, line_array: array) UgridDataset [source]#
Slice trough mapdata, combine: intersect_edges_withsort, calculation of distances and conversion to ugrid dataset.
Parameters#
- udsxu.UgridDataset
DESCRIPTION.
- line_arraynp.array
DESCRIPTION.
- calcdist_fromlatlonbool, optional
DESCRIPTION. The default is None.
Raises#
- Exception
DESCRIPTION.
Returns#
- xr_crs_ugridxu.UgridDataset
DESCRIPTION.
- dfm_tools.get_nc.rasterize_ugrid(uds: UgridDataset, ds_like: Dataset = None, resolution: float = None)[source]#
Rasterizing ugrid dataset to regular dataset. ds_like has higher priority than resolution. If both are not passed, a raster is generated of at least 200x200 inspired by xugrid.plot.imshow and xugrid.ugrid.ugrid2d.rasterize/rasterize_like.
Parameters#
- udsxu.UgridDataset
DESCRIPTION.
- ds_likexr.Dataset, optional
xr.Dataset with ed x and y variables to interpolate uds to. The default is None.
- resolutionfloat, optional
Only used if ds_like is not supplied. The default is None.
Raises#
- Exception
DESCRIPTION.
Returns#
- dsTYPE
DESCRIPTION.
- dfm_tools.get_nc.reconstruct_zw_zcc(uds: UgridDataset)[source]#
reconstruct full grid output (time/face-varying z-values) for all layertypes, calls the respective reconstruction function
Parameters#
- udsxu.UgridDataset
DESCRIPTION.
Raises#
- KeyError
DESCRIPTION.
Returns#
- udsTYPE
DESCRIPTION.
- dfm_tools.get_nc_helpers.rename_fouvars(ds: (<class 'xarray.core.dataset.Dataset'>, <class 'xugrid.core.wrap.UgridDataset'>), drop_tidal_times: bool = True)[source]#
Rename all fourier variables in a dataset (like mesh2d_fourier033_amp) to a unique name containing gridname/quantity/analysistype/tstart/tstop
Parameters#
- ds(xr.Dataset,xu.UgridDataset)
DESCRIPTION.
Returns#
- dsTYPE
DESCRIPTION.
- dfm_tools.get_nc_helpers.rename_waqvars(ds: (<class 'xarray.core.dataset.Dataset'>, <class 'xugrid.core.wrap.UgridDataset'>))[source]#
Rename all water quality variables in a dataset (like mesh2d_water_quality_output_24) to their long_name attribute (like mesh2d_DOscore)
Parameters#
- ds(xr.Dataset,xu.UgridDataset)
DESCRIPTION.
Returns#
- dsTYPE
DESCRIPTION.
- dfm_tools.xugrid_helpers.add_network_cellinfo(uds: UgridDataset)[source]#
Reads a UgridDataset with a minimal topology as it occurs in dflowfm netfiles, this contains only node_x, node_y and edge_node_connectivity. xugrid reads this as Ugrid1d It converts it to a Ugrid2d grid which includes the face_node_connectivity. It also couples the variables like NetNode_z again, but drops the edge_dims to avoid conflicts.
- dfm_tools.xugrid_helpers.enrich_rst_with_map(ds_rst: Dataset)[source]#
enriches rst dataset with topology from mapfile in the same folder. in order to merge, the bnd dimension has to be dropped. Also, the domain variable is currently not merged, so ghostcells are not removed.
Parameters#
- ds_rstxr.Dataset
DESCRIPTION.
Returns#
- ds_rstTYPE
DESCRIPTION.
- dfm_tools.xugrid_helpers.open_dataset_curvilinear(file_nc, varn_lon='longitude', varn_lat='latitude', varn_vert_lon='vertices_longitude', varn_vert_lat='vertices_latitude', ij_dims=['i', 'j'], convert_360to180=False, **kwargs)[source]#
This is a first version of a function that creates a xugrid UgridDataset from a curvilinear dataset like CMCC. Curvilinear means in this case 2D lat/lon variables and i/j indexing. The CMCC dataset does contain vertices, which is essential for conversion to ugrid. It also works for WAQUA files that are converted with getdata
- dfm_tools.xugrid_helpers.open_partitioned_dataset(file_nc: str, decode_fillvals: bool = False, remove_edges: bool = False, remove_ghost: bool = True, **kwargs)[source]#
using xugrid to read and merge partitions, including support for delft3dfm mapformat1 by renaming old layerdim. Furthermore some optional extensions like removal of hanging edges and removal of ghost cells.
Parameters#
- file_ncstr
DESCRIPTION.
- decode_fillvalsbool, optional
DESCRIPTION. The default is False.
- remove_edgesbool, optional
Remove hanging edges from the mapfile, necessary to generate contour and contourf plots with xugrid. Can be enabled always, but takes some additional computation time. The default is False.
- remove_ghostbool, optional
Remove ghostcells from the partitions. This is also done by xugrid automatically upon merging, but then the domain numbers are not taken into account so the result will be different. The default is True.
- file_ncTYPE
DESCRIPTION.
- kwargsTYPE, optional
arguments that are passed to xr.open_dataset. The chunks argument is set if not provided chunks={‘time’:1} increases performance significantly upon reading, but causes memory overloads when performing sum/mean/etc actions over time dimension (in that case 100/200 is better). The default is {‘time’:1}.
Raises#
- Exception
DESCRIPTION.
Returns#
- ds_merged_xuTYPE
DESCRIPTION.
- dfm_tools.xugrid_helpers.uda_interfaces_to_centers(uda_int: UgridDataArray) UgridDataArray [source]#
- dfm_tools.xugrid_helpers.uda_to_faces(uda: UgridDataArray) UgridDataArray [source]#
Interpolates a ugrid variable (xu.DataArray) with a node or edge dimension to the faces by averaging the 3/4 nodes/edges around each face.
Parameters#
- uda_nodexu.UgridDataArray
DESCRIPTION.
Raises#
- KeyError
DESCRIPTION.
Returns#
- uda_facexu.UgridDataArray
DESCRIPTION.
- dfm_tools.coastlines.get_borders_gdb(res: str = 'h', bbox: tuple = (-180, -90, 180, 90), crs: str = None) GeoSeries [source]#
GSHHS coastlines: https://www.ngdc.noaa.gov/mgg/shorelines/data/gshhg/latest/readme.txt geopandas docs https://geopandas.org/en/stable/docs/reference/api/geopandas.read_file.html
Parameters#
- resstr, optional
f(ull), h(igh), i(ntermediate), l(ow), c(oarse) resolution. The default is ‘h’.
- bboxtuple, optional
(minx, miny, maxx, maxy), also includes shapes that are partly in the bbox. The default is (-180, -90, 180, 90).
- crsstr, optional
coordinate reference system
Returns#
- coastlines_gdbTYPE
DESCRIPTION.
- dfm_tools.coastlines.get_coastlines_gdb(res: str = 'h', bbox: tuple = (-180, -90, 180, 90), min_area: float = 0, crs: str = None, columns: list = ['area']) GeoSeries [source]#
GSHHS coastlines: https://www.ngdc.noaa.gov/mgg/shorelines/data/gshhg/latest/readme.txt geopandas docs https://geopandas.org/en/stable/docs/reference/api/geopandas.read_file.html
Parameters#
- resstr, optional
f(ull), h(igh), i(ntermediate), l(ow), c(oarse) resolution. The default is ‘h’.
- bboxtuple, optional
(minx, miny, maxx, maxy), also includes shapes that are partly in the bbox. The default is (-180, -90, 180, 90).
- min_areafloat, optional
in km2, min_area>0 speeds up process. The default is 0.
- crsstr, optional
coordinate reference system
- columnslist, optional
which shapefile columns to include, None gives all. The default is [‘area’].
Returns#
- coastlines_gdbTYPE
DESCRIPTION.
- dfm_tools.coastlines.plot_borders(ax=None, res: str = 'h', crs=None, **kwargs)[source]#
get borders with get_borders_gdb and bbox depending on axlims, plot on ax and set axlims back to original values
- dfm_tools.coastlines.plot_coastlines(ax=None, res: str = 'h', min_area: float = 0, crs=None, **kwargs)[source]#
get coastlines with get_coastlines_gdb and bbox depending on axlims, plot on ax and set axlims back to original values
- class dfm_tools.linebuilder.LineBuilder(ax=None)[source]#
Bases:
object
To interactively draw a line in a figure axis, for instance to use as cross-section line for dfmt.polyline_mapslice().
ctrl+leftmouseclick to add a point to the line
ctrl+rightmouseclick to remove the last point of the line
ctrl+doublemouseclick to finish and let the script continue
- property line_array#
numpy array of the x/y coordinates of the interactively clicked line.
dfm_tools pre-processing#
- dfm_tools.modelbuilder.cmems_nc_to_bc(ext_bnd, list_quantities, tstart, tstop, file_pli, dir_pattern, dir_output, conversion_dict=None, refdate_str=None)[source]#
- dfm_tools.modelbuilder.cmems_nc_to_ini(ext_old, dir_output, list_quantities, tstart, dir_pattern, conversion_dict=None)[source]#
This makes quite specific 3D initial input files based on CMEMS data. delft3dfm is quite picky, so it works with CMEMS files because they have a depth variable with standard_name=’depth’. If this is not present, or dimensions are ordered differently, there will be an error or incorrect model results.
- dfm_tools.modelbuilder.create_model_exec_files(file_mdu, nproc=1, dimrset_folder=None, path_style=None)[source]#
creates a dimr_config.xml and if desired a batfile to run the model
- dfm_tools.modelbuilder.make_paths_relative(mdu_file: str)[source]#
Making paths on windows and unix relative by removing the dir_model prefix from paths in the mdu and ext files This is a temporary workaround until implemented in Deltares/HYDROLIB-core#532
Parameters#
- mdu_filestr
path to mdu_file.
Returns#
None.
- dfm_tools.modelbuilder.preprocess_merge_meteofiles_era5(ext_old, varkey_list, dir_data, dir_output, time_slice)[source]#
- dfm_tools.download.download_CMEMS(varkey, longitude_min, longitude_max, latitude_min, latitude_max, date_min, date_max, freq='D', dataset_id=None, buffer=None, dir_output='.', file_prefix='', overwrite=False)[source]#
https://help.marine.copernicus.eu/en/articles/8283072-copernicus-marine-toolbox-api-subset
- dfm_tools.download.download_ERA5(varkey, longitude_min, longitude_max, latitude_min, latitude_max, date_min, date_max, dir_output='.', overwrite=False)[source]#
empty docstring
- dfm_tools.download.download_OPeNDAP(dataset_url, varkey, longitude_min, longitude_max, latitude_min, latitude_max, date_min, date_max, freq='D', dir_output='.', file_prefix='', overwrite=False)[source]#
Parameters#
- dataset_urlTYPE
DESCRIPTION.
- varkeyTYPE
DESCRIPTION.
- longitude_minTYPE
DESCRIPTION.
- longitude_maxTYPE
DESCRIPTION.
- latitude_minTYPE
DESCRIPTION.
- latitude_maxTYPE
DESCRIPTION.
- date_minTYPE
DESCRIPTION.
- date_maxTYPE
DESCRIPTION.
- freqTYPE, optional
DESCRIPTION. The default is ‘D’.
- dir_outputTYPE, optional
DESCRIPTION. The default is ‘.’.
- file_prefixTYPE, optional
DESCRIPTION. The default is ‘’.
- overwriteTYPE, optional
DESCRIPTION. The default is False.
Raises#
- KeyError
DESCRIPTION.
- OutOfRangeError
DESCRIPTION.
Returns#
None.
- dfm_tools.interpolate_grid2bnd.get_conversion_dict(ncvarname_updates={})[source]#
get the conversion_dict, optionally with updated ncvarnames conversion_dict.keys() are the dflowfm quantities and the ncvarname the corresponding netcdf variable name/key alternative ncvarnames can be supplied via ncvarname_updates, e.g. get_conversion_dict(ncvarname_updates={‘temperaturebnd’:’tos’})
interpolate_nc_to_bc() renames netcdf variable like this: data_xr = data_xr.rename({ncvarname:quantity})
for CMCC: conversion_dict = { # conversion is phyc in mol/m3 to newvar in g/m3 ‘tracerbndOXY’ : {‘ncvarname’: ‘o2’, ‘unit’: ‘g/m3’, ‘conversion’: 32.0 }, ‘tracerbndNO3’ : {‘ncvarname’: ‘no3’, ‘unit’: ‘g/m3’, ‘conversion’: 14.0 }, ‘tracerbndPO4’ : {‘ncvarname’: ‘po4’, ‘unit’: ‘g/m3’, ‘conversion’: 30.97 }, ‘tracerbndSi’ : {‘ncvarname’: ‘si’, ‘unit’: ‘g/m3’, ‘conversion’: 28.08}, ‘tracerbndPON1’ : {‘ncvarname’: ‘phyc’, ‘unit’: ‘g/m3’, ‘conversion’: 14.0}, # Caution: this empirical relation might not be applicable to your use case ‘tracerbndPOP1’ : {‘ncvarname’: ‘phyc’, ‘unit’: ‘g/m3’, ‘conversion’: 30.97}, # Caution: this empirical relation might not be applicable to your use case ‘tracerbndPOC1’ : {‘ncvarname’: ‘phyc’, ‘unit’: ‘g/m3’, ‘conversion’: 14.0 * (106 / 16)}, # Caution: this empirical relation might not be applicable to your use case ‘salinitybnd’ : {‘ncvarname’: ‘sos’}, #’1e-3’ ‘temperaturebnd’ : {‘ncvarname’: ‘tos’}, #’degC’ ‘ux’ : {‘ncvarname’: ‘uo’}, #’m/s’ ‘uy’ : {‘ncvarname’: ‘vo’}, #’m/s’ ‘waterlevelbnd’ : {‘ncvarname’: ‘zos’}, #’m’ #steric ‘tide’ : {‘ncvarname’: ‘’}, #’m’ #tide (dummy entry) }
The below clarification about the CMEMS conversion factors in this function is by Jos van Gils. The use of boundary conditions from CMEMS requires some conversion (scaling). In the first place, there is the unit conversion from mmol/m3 in CMEMS to g/m3 in DFM. This involves a factor M/1000, where M is the molar weight of the constituent in question: M = 12 for constituents expressed in gC/m3, 14 for gN/m3, 30.97 for gP/m3, 28.08 for gSi/m3 and 32 for dissolved oxygen.
Next there is the issue that CMEMS does not provide all state variables in the water quality model. Dead organic matter, particulate and dissolved, are missing. In practice, we use phytoplankton biomass from CMEMS (phyc) to estimate these missing variables. For particulate organic matter (POM), often a carbon-based ratio of 2.0 to phytoplankton biomass is used. PON and POP are derived from POC by classical Redfield ratios (C:N:P = 106:16:1) or revised Redfield ratios (C:N:P = 117:14:1). For Opal (biogenic silica from diatoms), the example here is based on a factor 0.5 expressing the share of diatoms in phytoplankton and a 0.13 Si:C stoichiometric ratio by Brzezinski (1985). For dissolved organic matter (DOM) there is limited experience. This conversion_dict uses the approach used for the “Bays Eutrophication Model” (Massachusetts Bay, https://www.mwra.com/harbor/enquad/pdf/2021-02.pdf). This relies on local field data.
It is noted that such scale factors are anyhow inaccurate and typically water-system-specific. Therefore, the values included in this conversion_dict should be considered indicative only and should by no means be seen as a default approach. Ideally, the scale factors to link missing states to available CMEMS variables are derived from local field data. A validation is strongly recommended, for example by a comparison of simulated and measured concentrations of total N and total P in the study area. Similar considerations would hold for other data sources than CMEMS.
- dfm_tools.interpolate_grid2bnd.interp_hisnc_to_plipoints(data_xr_his, file_pli, kdtree_k=3, load=True)[source]#
interpolate stations in hisfile to points of polyfile
- dfm_tools.interpolate_grid2bnd.interp_regularnc_to_plipointsDataset(data_xr_reg, gdf_points, load=True)[source]#
- dfm_tools.interpolate_grid2bnd.interp_uds_to_plipoints(uds: UgridDataset, gdf: GeoDataFrame) Dataset [source]#
To interpolate an unstructured dataset (like a _map.nc file) read with xugrid to plipoint locations
Parameters#
- udsxu.UgridDataset
dfm model output read using dfm_tools. Dims: mesh2d_nLayers, mesh2d_nInterfaces, time, mesh2d_nNodes, mesh2d_nFaces, mesh2d_nMax_face_nodes, mesh2d_nEdges.
- gdfgeopandas.GeoDataFrame
gdf with location, geometry (Point) and crs corresponding to model crs.
Raises#
- Exception
DESCRIPTION.
Returns#
- dsTYPE
xr.Dataset with dims: node, time, z.
- dfm_tools.interpolate_grid2bnd.interpolate_tide_to_bc(ext_new: ExtModel, tidemodel, file_pli, component_list=None, load=True)[source]#
- dfm_tools.interpolate_grid2bnd.interpolate_tide_to_plipoints(tidemodel, gdf_points, component_list=None, load=True)[source]#
empty docstring
- dfm_tools.interpolate_grid2bnd.plipointsDataset_to_ForcingModel(plipointsDataset)[source]#
empty docstring
- dfm_tools.hydrolib_helpers.DataFrame_to_PolyObject(poly_pd, name, content=None)[source]#
convert a pandas dataframe with x/y columns (and optional others like z/data/comment) to a hydrolib PolyObject
- dfm_tools.hydrolib_helpers.DataFrame_to_TimModel(tim_pd, refdate: (<class 'datetime.datetime'>, <class 'pandas._libs.tslibs.timestamps.Timestamp'>, <class 'str'>))[source]#
converts data from tim_pd to TimModel and puts all headers as comments. Ignores the index, assumes first column is time in minutes since a refdate.
- dfm_tools.hydrolib_helpers.Dataset_to_Astronomic(datablock_xr)[source]#
convert an xarray.Dataset (with amplitude and phase data_vars) to a hydrolib Astronomic object
- dfm_tools.hydrolib_helpers.Dataset_to_T3D(datablock_xr)[source]#
convert an xarray.DataArray (is one data_var) or an xarray.Dataset (with one or two data_vars) with time and depth dimension to a hydrolib T3D object
- dfm_tools.hydrolib_helpers.Dataset_to_TimeSeries(datablock_xr)[source]#
convert an xarray.DataArray or xarray.Dataset with time dimension to a hydrolib TimeSeries object
- dfm_tools.hydrolib_helpers.ForcingModel_to_plipointsDataset(forcingmodel: ForcingModel, npoints=None, convertnan=False) Dataset [source]#
- dfm_tools.hydrolib_helpers.PolyFile_to_geodataframe_linestrings(polyfile_object, crs=None)[source]#
empty docstring
- dfm_tools.hydrolib_helpers.PolyFile_to_geodataframe_points(polyfile_object: PolyFile, crs: str = None, add_pointnames: bool = True)[source]#
Parameters#
- polyfile_objecthcdfm.PolyFile
get this object with hcdfm.PolyFile(path_to_plifile).
- crsstr, optional
DESCRIPTION. The default is None’.
- add_pointnamesbool, optional
DESCRIPTION. The default is True.
Returns#
- gdfTYPE
gdf of all the pli files with columns: location, geometry (Points). Crs = 4326 (decimal degrees).
- dfm_tools.hydrolib_helpers.TimModel_to_DataFrame(data_tim: ~hydrolib.core.dflowfm.tim.models.TimModel, parse_column_labels: bool = True, refdate: (<class 'datetime.datetime'>, <class 'pandas._libs.tslibs.timestamps.Timestamp'>, <class 'str'>) = None)[source]#
Parameters#
- data_timhcdfm.TimModel
DESCRIPTION.
- parse_column_labelsbool, optional
Parse column labels from comments. This might fail since there is no standard way of prescribing the columns in the comments. The default is True.
- refdate(dt.datetime, pd.Timestamp, str, int, float), optional
DESCRIPTION. The default is None.
Returns#
- tim_pdTYPE
DESCRIPTION.
- dfm_tools.hydrolib_helpers.forcinglike_to_Dataset(forcingobj, convertnan=False)[source]#
convert a hydrolib forcing like object (like Timeseries, T3D, Harmonic, etc) to an xarray Dataset with one or more variables.
convertnan: convert depths with the same values over time as the deepest layer to nan (these were created with .bfill() or .ffill()).
- dfm_tools.hydrolib_helpers.geodataframe_to_PolyFile(poly_gdf, name='L')[source]#
convert a geopandas geodataframe with x/y columns (and optional others like z/data/comment) to a hydrolib PolyFile
- dfm_tools.hydrolib_helpers.pointlike_to_DataFrame(pointlike, drop_emptycols=True)[source]#
convert a hydrolib object with points (like PolyObject, XYZModel and possibly others) to a pandas DataFrame.
Parameters#
- pointlikeTYPE
Hydrolib-core object with point objects.
Returns#
- pointlike_pdTYPE
DESCRIPTION.
- drop_emptycolsbool, optional
Drop empty (all-nan) columns automatically, like the z-column in ldb file. The default is True.
- Example:
polyfile_object = PolyFile(file_pli) data_pol_pd_list = [pointlike_to_DataFrame(polyobj) for polyobj in polyfile_object.objects]
- dfm_tools.hydrolib_helpers.pointlike_to_geodataframe_points(polyline_object, crs: str = None, add_pointnames=True, only_xy=False)[source]#
empty docstring
- dfm_tools.meshkernel_helpers.generate_bndpli_cutland(mk: ~meshkernel.meshkernel.MeshKernel, res: str = 'f', min_area: float = 0, crs: (<class 'int'>, <class 'str'>) = None, buffer: float = 0)[source]#
Generate a boundary polyline from the meshkernel object and cut away the landward part. Be sure to do this on the base/refined grid, not on the grid where the landward cells were already cut.
Parameters#
- mkmeshkernel.MeshKernel
DESCRIPTION.
- resstr, optional
DESCRIPTION. The default is ‘f’.
- min_areafloat, optional
DESCRIPTION. The default is 0.
- crs(int,str), optional
DESCRIPTION. The default is None.
- bufferfloat, optional
DESCRIPTION. The default is 0.
Returns#
- bnd_gdfTYPE
DESCRIPTION.
- dfm_tools.meshkernel_helpers.interpolate_bndpli(bnd_gdf, res)[source]#
interpolate bnd_gdf to a new resolution
- dfm_tools.meshkernel_helpers.make_basegrid(lon_min, lon_max, lat_min, lat_max, dx, dy, angle=0, crs=None, is_geographic=None)[source]#
empty docstring
- dfm_tools.meshkernel_helpers.meshkernel_delete_withcoastlines(mk: ~meshkernel.meshkernel.MeshKernel, res: str = 'f', min_area: float = 0, crs: (<class 'int'>, <class 'str'>) = None)[source]#
Wrapper around meshkernel_delete_withgdf, which automatically gets the bbox from the meshkernel object and retrieves the coastlines_gdf.
Parameters#
- mkmeshkernel.MeshKernel
DESCRIPTION.
- resstr, optional
DESCRIPTION. The default is ‘f’.
- min_areafloat, optional
DESCRIPTION. The default is 0.
- crs(int,str), optional
DESCRIPTION. The default is None.
Returns#
None.
- dfm_tools.meshkernel_helpers.meshkernel_delete_withgdf(mk: MeshKernel, coastlines_gdf: GeoDataFrame)[source]#
Delete parts of mesh that are inside the polygons/Linestrings in a GeoDataFrame.
Parameters#
- mkmeshkernel.MeshKernel
DESCRIPTION.
- coastlines_gdfgpd.GeoDataFrame
DESCRIPTION.
Returns#
None.
- dfm_tools.meshkernel_helpers.meshkernel_delete_withshp(mk: ~meshkernel.meshkernel.MeshKernel, coastlines_shp: str, crs: (<class 'int'>, <class 'str'>) = None)[source]#
Delete parts of mesh that are inside the shapefile polygon.
Parameters#
- mkmeshkernel.MeshKernel
DESCRIPTION.
- coastlines_shpstr
Path to the shp file.
- crs(int,str), optional
DESCRIPTION. The default is None.
Returns#
None.
- dfm_tools.meshkernel_helpers.meshkernel_to_UgridDataset(mk: ~meshkernel.meshkernel.MeshKernel, crs: (<class 'int'>, <class 'str'>) = None) UgridDataset [source]#
Convert a meshkernel object to a UgridDataset, including a variable with the crs (used by dflowfm to distinguish spherical/cartesian networks). The UgridDataset enables bathymetry interpolation and writing to netfile.
Parameters#
- mkmeshkernel.MeshKernel
DESCRIPTION.
- crs(int,str), optional
DESCRIPTION. The default is None.
Returns#
- TYPE
DESCRIPTION.
- dfm_tools.meshkernel_helpers.refine_basegrid(mk, data_bathy_sel, min_edge_size)[source]#
empty docstring
- dfm_tools.xarray_helpers.merge_meteofiles(file_nc: str, preprocess=None, time_slice: slice = slice(None, None, None), add_global_overlap: bool = False, zerostart: bool = False, **kwargs) Dataset [source]#
for merging for instance meteo files x/y and lon/lat are renamed to longitude/latitude #TODO: is this desireable?
Parameters#
- file_ncstr
DESCRIPTION.
- preprocessTYPE, optional
DESCRIPTION. The default is None.
- time_sliceslice, optional
DESCRIPTION. The default is slice(None,None).
- add_global_overlapbool, optional
GTSM specific: extend data beyond -180 to 180 longitude. The default is False.
- zerostartbool, optional
GTSM specific: extend data with 0-value fields 1 and 2 days before time_slice.start. The default is False.
- kwargsdict, optional
arguments for xr.open_mfdataset() like chunks to prevent large chunks and resulting memory issues.
Returns#
- data_xrxr.Dataset
Merged meteo dataset.
- dfm_tools.xarray_helpers.preprocess_ERA5(ds)[source]#
Reduces the expver dimension in some of the ERA5 data (mtpr and other variables), which occurs in files with very recent data. The dimension contains the unvalidated data from the latest month in the second index in the expver dimension. The reduction is done with mean, but this is arbitrary, since there is only one valid value per timestep and the other one is nan.
- dfm_tools.xarray_helpers.preprocess_hisnc(ds)[source]#
Look for dim/coord combination and use this for Dataset.set_index(), to enable station/gs/crs/laterals label based indexing. If duplicate labels are found (like duplicate stations), these are dropped to avoid indexing issues.
Parameters#
- dsxarray.Dataset
DESCRIPTION.
- drop_duplicate_stationsTYPE, optional
DESCRIPTION. The default is True.
Returns#
- dsTYPE
DESCRIPTION.
dfm_tools insitu data retrieval#
- dfm_tools.observations.ssh_catalog_subset(source=None, lon_min=-180, lon_max=180, lat_min=-90, lat_max=90, time_min=None, time_max=None, **kwargs)[source]#
- dfm_tools.observations.ssh_catalog_tokmlfile(ssc_catalog_gpd, file_kml)[source]#
converts a ssh_catalog to a kml file for easy visualisation in google earth
- dfm_tools.observations.ssh_catalog_toxynfile(ssc_catalog_gpd, file_xyn)[source]#
converts a ssh_catalog to a _obs.xyn file for easy visualisation in FM GUI and interacter
dfm_tools unimportant stuff#
- dfm_tools.bathymetry.read_asc(file_asc: str) Dataset [source]#
Reading asc file into a xarray.Dataset
Parameters#
- file_ascstr
asc file with header ncols, nrows, xllcenter, yllcenter, cellsize, NODATA_value.
Returns#
- ds_ascxr.Dataset
xarray.Dataset with the data from the asc file as an array and the lat/lon coordinates as separate coordinate variables.
- dfm_tools.bathymetry.write_bathy_toasc(filename_asc, lon_sel_ext, lat_sel_ext, elev_sel_ext, asc_fmt='%9.2f', nodata_val=32767)[source]#
empty docstring
@author: Kieran Hunt kieranmrhunt/curved-quivers https://stackoverflow.com/questions/51843313/flow-visualisation-in-python-using-curved-path-following-vectors aligned with matplotlib.streamplot and improved by Jelmer Veenstra (30-03-2023), veenstrajelmer matplotlib.streamplot available at https://raw.githubusercontent.com/matplotlib/matplotlib/main/lib/matplotlib/streamplot.py
Streamline plotting for 2D vector fields.
- dfm_tools.modplot.velovect(axes, x, y, u, v, density=1, linewidth=None, color=None, cmap=None, norm=None, arrowsize=1, arrowstyle='-|>', transform=None, zorder=None, start_points=None, integration_direction='both', grains=15, broken_streamlines=True)[source]#
Draw streamlines of a vector flow.
Parameters#
- x, y1D/2D arrays
Evenly spaced strictly increasing arrays to make a grid. If 2D, all rows of x must be equal and all columns of y must be equal; i.e., they must be as if generated by
np.meshgrid(x_1d, y_1d)
.- u, v2D arrays
x and y-velocities. The number of rows and columns must match the length of y and x, respectively.
- densityfloat or (float, float)
Controls the closeness of streamlines. When
density = 1
, the domain is divided into a 30x30 grid. density linearly scales this grid. Each cell in the grid can have, at most, one traversing streamline. For different densities in each direction, use a tuple (density_x, density_y).- linewidthfloat or 2D array
The width of the streamlines. With a 2D array the line width can be varied across the grid. The array must have the same shape as u and v.
- colorcolor or 2D array
The streamline color. If given an array, its values are converted to colors using cmap and norm. The array must have the same shape as u and v.
- cmap, norm
Data normalization and colormapping parameters for color; only used if color is an array of floats. See ~.Axes.imshow for a detailed description.
- arrowsizefloat
Scaling factor for the arrow size.
- arrowstylestr
Arrow style specification. See ~matplotlib.patches.FancyArrowPatch.
- minlengthfloat
Minimum length of streamline in axes coordinates.
- start_points(N, 2) array
Coordinates of starting points for the streamlines in data coordinates (the same coordinates as the x and y arrays).
- zorderfloat
The zorder of the streamlines and arrows. Artists with lower zorder values are drawn first.
- maxlengthfloat
Maximum length of streamline in axes coordinates.
- integration_direction{‘forward’, ‘backward’, ‘both’}, default: ‘both’
Integrate the streamline in forward, backward or both directions.
- dataindexable object, optional
DATA_PARAMETER_PLACEHOLDER
- broken_streamlinesboolean, default: True
If False, forces streamlines to continue until they leave the plot domain. If True, they may be terminated if they come too close to another streamline.
Returns#
- StreamplotSet
Container object with attributes
lines
: .LineCollection of streamlinesarrows
: .PatchCollection containing .FancyArrowPatch objects representing the arrows half-way along streamlines.This container will probably change in the future to allow changes to the colormap, alpha, etc. for both lines and arrows, but these changes should be backward compatible.
- dfm_tools.energy_dissipation.compute_energy_dissipation(data_xr_map, file_ED_computed)[source]#
- Example:
data_xr_map = dfmt.open_partitioned_dataset(file_nc_map,chunks={‘time’:100}) #TODO: important to have time>1, otherwise time-mean floods memory (100-200 seems optimal for GTSM 1month, but memory still floods so 1 year would probably be impossible). data_xr_map = data_xr_map.sel(time=slice(‘2014-01-01’,’2014-02-01’)) dfmt.compute_energy_dissipation(data_xr_map,file_ED_computed)
TODO: clean up these comments
Stel de bodemwrijving is (tau_x,tau_y) dan is het verlies aan energie in W/m^2 E_loss = tau_x*u_x + tau_y*u_y
Voor Chezy geldt (zo uit mijn hoofd): tau_x = - (rho*g)/(C^2) u_x * sqrt(u_x^2+u_y^2) tau_y = - (rho*g)/(C^2) u_y * sqrt(u_x^2+u_y^2)
En gecombineerd E_loss = - (rho*g)/(C^2) * (u_x^2+u_y^2)^(3/2)
Approximation: using velocities in cell centers, this is averaged but easiest
E_loss(area) = sum_cells [ - (rho*g)/(C^2) * (sqrt(u_x^2+u_y^2))^3 ] * cell_area
Dit is dan voor een tijdstip. Ik zou middelen over een iets langere periode zodat het de variatie over het getij verdwijnt, want de waarde is bij springtij waarschijnlijk groter.