Documentation
Spatial pre-processing functions
- wavebreaking.processing.spatial.calculate_momentum_flux(u, v, *args, **kwargs)[source]
Calculate the momentum flux derived from the product of the deviations of both wind components from the zonal mean. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value arguments.
Parameters
- uxarray.DataArray
zonal (x) component of the wind
- vxarray.DataArray
meridional (y) component of the wind
Returns
- momentum flux: xarray.DataArray
Data containing the momentum flux
- wavebreaking.processing.spatial.calculate_smoothed_field(data, passes, weights=array([[0, 1, 0], [1, 2, 1], [0, 1, 0]]), mode='wrap', *args, **kwargs)[source]
Calculate smoothed field based on a two-dimensional weight kernel and multiple smoothing passes. Default weight kernel is a 3x3 5-point smoothing with double-weighted centre. The arguments “weight” and “mode” must be accepted by scipy.ndimage.convolve. Values at the latitude border are always set to NaN. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value arguments.
Parameters
- dataxarray.DataArray
data to smooth
- passesint or float
number of smoothing passes of the 5-point smoothing
- weigthsarray_like, optional
array of weight, two-dimensional (see scipy.ndimage.convolve function)
- modestring, optional
defines how the array is extended at boundaries (see scipy.ndimage.convolve function)
Returns
- smoothed data: xarray.DataArray
Data containing the smoothed field
Indices
- wavebreaking.indices.contour_index.calculate_contours(data, contour_levels, periodic_add=120, original_coordinates=True, *args, **kwargs)[source]
Calculate contour lines for a set of contour levels. The calculations are based on the measure.find_contours module. If periodic_add is provided, the data array is expanded in the longitudinal direction and undulations at the date border are captured correctly. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument.
Parameters
- dataxarray.DataArray
data for the contour calculation
- contour_levelsarray_like
levels for contour calculation
- periodic_add: int or float, optional
number of longitudes in degrees to expand the dataset to correctly capture undulations at the date border if the input field is not periodic, use periodic_add = 0
- original_coordinates: bool, optional
if False, the array indices of the contour lines are returned instead of the original coordinates
Returns
- contours: geopandas.GeoDataFrame
- GeoDataFrame containing different characteristics of the contours:
“date”: date of the contour line
“level”: level of the contour line
“closed”: True if contour line is closed
“exp_lon”: expansion in degrees of the contours in the longitudinal direction
“mean_lat”: mean latitude of the contours
“geometry”: LineString object with the contour coordinates in the format (x,y)
- wavebreaking.indices.streamer_index.calculate_streamers(data, contour_level, contours=None, geo_dis=800, cont_dis=1500, intensity=None, periodic_add=120, *args, **kwargs)[source]
Identify streamer structures based on the Streamer Index developed by Wernli and Sprenger (2007). The default parameters for the streamer identification are based on the study by Wernli and Sprenger (2007). Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument. Before the index calculation, the contour lines are calculated if not provided.
Parameters
- dataxarray.DataArray
data for the contour and streamer calculation
- contour_levelsarray_like
levels for contour calculation
- contoursgeopandas.GeoDataFrame, optional
contours calculated with wavebreaking.calculate_contours(…, original_coordinates=False)
- geo_disint or float, optional
Maximal geographic distance between two contour points that describe a streamer
- cont_disint or float, optional
Minimal distance along the contour line between two points that describe a streamer
- intensityxarray.DataArray, optional
data for the intensity calculation (hint: use wb_spatial.calculate_momentum_flux)
- periodic_add: int or float, optional
number of longitudes in degrees to expand the dataset to correctly capture undulations at the date border if the input field is not periodic, use periodic_add = 0
Returns
- streamers: geopandas.GeoDataFrame:
- GeoDataFrame containing different characteristics of the streamers
“date”: date of the streamers
“level”: level of the contour line
“com”: center of mass in the format (x,y)
“mean_var”: mean of the variable used for the streamer calculations
“event_area”: area of a streamer
“intensity”: sum of the intensity (momentum flux)
“geometry”: (Multi)Polygon with coordinates in the format (x,y)
- wavebreaking.indices.overturning_index.calculate_overturnings(data, contour_levels, contours=None, range_group=5, min_exp=5, intensity=None, periodic_add=120, *args, **kwargs)[source]
Identify overturning structures based on the Overturning Index developed by Barnes and Hartmann (2012). The default parameters for the overturning identification are based on the study by Barnes and Hartmann (2012). The overturning region is represented by a rectangle enclosing all overturning contour points. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument. Before the index calculation, the contour lines are calculated if not provided.
Parameters
- dataxarray.DataArray
data for the contour and overturning calculation
- contour_levelsarray_like
levels for contour calculation
- contoursgeopandas.GeoDataFrame, optional
contours calculated with wavebreaking.calculate_contours(…, original_coordinates=False)
- range_groupint or float, optional
Maximal degrees in the longitudinal direction in which two overturning are grouped
- min_expint or float, optional
Minimal longitudinal expansion of an overturning event
- intensityxarray.DataArray, optional
data for the intensity calculation (hint: use wb_spatial.calculate_momentum_flux)
- periodic_add: int or float, optional
number of longitudes in degrees to expand the dataset to correctly capture undulations at the date border if the input field is not periodic, use periodic_add = 0
Returns
- overturnings: geopandas.GeoDataFrame
- GeoDataFrame containing different characteristics of the overturning events:
“date”: date of the overturning
“level”: level of the contour line
“com”: center of mass in the format (x,y)
“mean_var”: mean of the variable used for the overturning calculations
“event_area”: area of a overturning event;
“intensity”: sum of the intensity (momentum flux)
“orientation”: orientation of the most west- and eastward point
“geometry”: (Multi)Polygon with coordinates in the format (x,y)
- wavebreaking.indices.cutoff_index.calculate_cutoffs(data, contour_level, contours=None, min_exp=5, intensity=None, periodic_add=120, *args, **kwargs)[source]
Identify cutoff structures. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument. Before the index calculation, the contour lines are calculated if not provided.
Parameters
- dataxarray.DataArray
data for the contour and cutoff calculation
- contour_levelsarray_like
levels for contour calculation
- contoursgeopandas.GeoDataFrame, optional
contours calculated with wavebreaking.calculate_contours(…, original_coordinates=False)
- min_expint or float, optional
Minimal longitudinal expansion of a cutoff event
- intensityxarray.DataArray, optional
data for the intensity calculation (hint: use wb_spatial.calculate_momentum_flux)
- periodic_add: int or float, optional
number of longitudes in degrees to expand the dataset to correctly capture undulations at the date border if the input field is not periodic, use periodic_add = 0
Returns
- cutoffs: geopandas.GeoDataFrame
- GeoDataFrame containing different characteristics of the cutoffs events:
“date”: date of the cutoffs
“level”: level of the contour line
“com”: center of mass in the format (x,y)
“mean_var”: mean of the variable used for the contour calculation
“event_area”: area of a cutoff event
“intensity”: sum of the intensity (momentum flux)
“geometry”: (Multi)Polygon with the coordinates in the format (x,y)
Events post-processing functions
- wavebreaking.processing.events.to_xarray(data, events, flag='ones', name='flag', *args, **kwargs)[source]
Create xarray.DataArray from events stored in a geopandas.GeoDataFrame. Grid cells where an event is present are flagged with the value 1. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument.
Parameters
- dataxarray.DataArray
data used for the index calculation
- eventsgeopandas.GeoDataFrame
GeoDataFrame with the date and geometry for each event
- flagstring, optional
column name of the events geopandas.GeoDataFrame flag is set where an event is present default value is “ones”
- namestring, optional
name of the xarray variable that is created
Returns
- flag: xarray.DataArray
Data with events flagged with the value 1
- wavebreaking.processing.events.track_events(events, time_range=None, method='by_overlap', buffer=0, overlap=0, distance=1000)[source]
Temporal tracking of events. Events receive the same label if they spatially overlap at step t and t + time_range.
Parameters
- eventsgeopandas.GeoDataFrame
GeoDataFrame with the date and coordinates of each identified event
- time_range: int or float, optional
Time range for temporally tracking the events. The units of time_range is hours if the type of the time dimension is np.datetime64. If not specified, the smallest time difference larger than zero is used.
- method{“by_overlap”, “by_distance”}, optional
- Method for temporally tracking the events:
- “by_overlap”: Events receive the same label if they spatially
overlap at step t and t + time_range.
- “by_distance”: Events receive the same label if their centre of mass
is closer than “distance”
- bufferfloat, optional
buffer around event polygon in degrees for the ‘by_overlap’ method
- overlapfloat, optional
minimum percentage of overlapping for the ‘by_overlap’ method
- distanceint or float, optional
maximum distance in km between two events for the ‘by_distance’ method
Returns
- events: geopandas.GeoDataFrame
GeoDataFrame with label column showing the temporal coherence
Plotting functions
- wavebreaking.processing.plots.plot_clim(flag_data, seasons=None, proj=None, size=None, smooth_passes=5, periodic=True, labels=True, levels=None, cmap=None, title='', *args, **kwargs)[source]
Creates a simple climatological plot showing the occurrence frequency of the detected events. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument.
Parameters
- flag_dataxarray.DataArray
data containing the locations of the events flagged with the value 1
- seasonslist or array, optional
months of the seasons for occurrence frequency calculation (e.g. [12, 1, 2])
- projcartopy.crs, optional
cartopy projection object
- sizetuple of integers, optional
size of the figure in the format (width, height)
- smooth_passesint or float, optional
number of smoothing passes of the 5-point smoothing of the occurrence frequencies
- periodicbool, optional
If True, the first longitude is added at the end to close the gap in a polar projection
- labelsbool, optional
If False, no labels are added to the plot
- levelslist or array, optional
Colorbar levels. If not provided, default levels are used.
- cmapstring, optional
Name of a valid cmap. If not provided, a default cmap is used.
- titlestring, optional
Title of the plot
Returns
- plotmatplotlib.pyplot
Climatological plot of the occurrence frequencies.
- wavebreaking.processing.plots.plot_step(flag_data, step, data=None, contour_levels=None, proj=None, size=None, periodic=True, labels=True, levels=None, cmap='RdBu_r', color_events='gold', title='', *args, **kwargs)[source]
Creates a plot showing the events at one time step. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument.
Parameters
- flag_dataxarray.DataArray
Data containing the locations of the events flagged with the value 1
- stepint or string
index or name of a time step in the xarray.Dataset
- dataxarray.DataArray
Data that has been used to calculate the contours and the indices
- contour_levelarray_like
contour levels that are shown in the plot
- projcartopy.crs, optional
cartopy projection object
- sizetuple of integers, optional
size of the figure in the format (width, height)
- periodicbool, optional
If True, the first longitude is added at the end to close the gap in a polar projection
- labelsbool, optional
If False, no labels are added to the plot
- levelslist or array, optional
Colorbar levels. If not provided, default levels are used.
- cmapstring, optional
Name of a valid cmap. If not provided, a default cmap is used.
- color_eventsstring, optional
Color of the events
- titlestring, optional
Title of the plot
Returns
- plotmatplotlib.pyplot
Plot of one time step.
- wavebreaking.processing.plots.plot_tracks(data, events, proj=None, size=None, min_path=0, plot_events=False, labels=True, title='', *args, **kwargs)[source]
Creates a plot showing the tracks of the temporally coherent events. Dimension names (“time_name”, “lon_name”, “lat_name”), size (“ntime”, “nlon”, “nlat”) and resolution (“dlon”, “dlat”) can be passed as key=value argument.
Parameters
- dataxarray.DataArray
Data that has been used to calculate the contours and the indices
- events: pd.DataFrame
DataFrame with the date, coordinates and label of the identified events
- projcartopy.crs, optional
cartopy projection object
- sizetuple of integers, optional
size of the figure in the format (width, height)
- min_path: int, optional
Minimal number of time steps an event has to be tracked
- plot_events: bool, optional
If True, the events are also plotted by a shaded area
- labels: bool, optional
If False, no labels are added to the plot
- title: string, optional
Title of the plot
Returns
- plotmatplotlib.pyplot
Plot of the tracks