pysteps.blending.pca_ens_kalman_filter.forecast#
- pysteps.blending.pca_ens_kalman_filter.forecast(obs_precip, obs_timestamps, nwp_precip, nwp_timestamps, velocity, forecast_horizon, issuetime, n_ens_members, precip_mask_dilation=1, smooth_radar_mask_range=0, n_cascade_levels=6, precip_thr=-10.0, norain_thr=0.01, extrap_method='semilagrangian', decomp_method='fft', bandpass_filter_method='gaussian', noise_method='nonparametric', enkf_method='masked_enkf', enable_combination=True, noise_stddev_adj=None, ar_order=1, callback=None, return_output=True, seed=None, num_workers=1, fft_method='numpy', domain='spatial', extrap_kwargs=None, filter_kwargs=None, noise_kwargs=None, combination_kwargs=None, measure_time=False)#
Generate a combined nowcast ensemble by using the reduced-space ensemble Kalman filter method described in Nerini et al. 2019.
- Parameters:
obs_precip (np.ndarray) – Array of shape (ar_order+1,m,n) containing the observed input precipitation fields ordered by timestamp from oldest to newst. The time steps between the inputs are assumed to be regular.
obs_timestamps (np.ndarray) – Array of shape (ar_order+1) containing the corresponding time stamps of observed input precipitation fields as datetime objects.
nwp_precip (np.ndarray) – Array of shape (n_ens,n_times,m,n) containing the (NWP) ensemble model forecast.
nwp_timestamps (np.ndarray) – Array of shape (n_times) containing the corresponding time stamps of the (NWP) ensemble model forecast as datetime objects.
velocity (np.ndarray) – Array of shape (2,m,n) containing the x- and y-components of the advection field. The velocities are based on the observed input precipitation fields and are assumed to represent one time step between the inputs. All values are required to be finite.
forecast_horizon (int) – The length of the forecast horizon (the length of the forecast) in minutes.
issuetime (datetime object) – Issuetime of the combined forecast to compute.
n_ens_members (int) – The number of ensemble members to generate. This number should always be equal to or larger than the number of NWP ensemble members / number of NWP models.
precip_mask_dilation (int) – Range by which the precipitation mask within the forecast step should be extended per time step. Defaults to 1.
smooth_radar_mask_range (int, Default is 0.) – Method to smooth the transition between the radar-NWP-noise blend and the NWP-noise blend near the edge of the radar domain (radar mask), where the radar data is either not present anymore or is not reliable. If set to 0 (grid cells), this generates a normal forecast without smoothing. To create a smooth mask, this range should be a positive value, representing a buffer band of a number of pixels by which the mask is cropped and smoothed. The smooth radar mask removes the hard edges between NWP and radar in the final blended product. Typically, a value between 50 and 100 km can be used. 80 km generally gives good results.
n_cascade_levels (int, optional) – The number of cascade levels to use. Defaults to 6, see issue #385 on GitHub.
precip_thr (float, optional) – pecifies the threshold value for minimum observable precipitation intensity. Required if mask_method is not None or conditional is True. Defaults to -10.0.
norain_thr (float) – Specifies the threshold value for the fraction of rainy (see above) pixels in the radar rainfall field below which we consider there to be no rain. Depends on the amount of clutter typically present. Defaults to -15.0.
extrap_method (str, optional) – Name of the extrapolation method to use. See the documentation of
pysteps.extrapolation.interface. Defaults to ‘semilagrangian’.decomp_method ({'fft'}, optional) – Name of the cascade decomposition method to use. See the documentation of
pysteps.cascade.interface. Defaults to ‘fft’.bandpass_filter_method ({'gaussian', 'uniform'}, optional) – Name of the bandpass filter method to use with the cascade decomposition. See the documentation of
pysteps.cascade.interface. Defaults to ‘guassian’.noise_method ({'parametric','nonparametric','ssft','nested',None}, optional) – Name of the noise generator to use for perturbating the precipitation field. See the documentation of
pysteps.noise.interface. If set to None, no noise is generated. Defaults to ‘nonparametric’.enkf_method ({'masked_enkf}, optional) – Name of the ensemble Kalman filter method to use for the correction step. Currently, only ‘masked_enkf’ method is implemented that corresponds to the reduced-space ensemble Kalman filter technique described in Nerini et al. 2019. Defaults to ‘masked_enkf’.
enable_combination (bool, optional) – Flag to specify whether the correction step should be applied or a pure nowcasting ensemble should be computed. Defaults to True.
noise_stddev_adj ({'auto','fixed',None}, optional) – Optional adjustment for the standard deviations of the noise fields added to each cascade level. This is done to compensate incorrect std. dev. estimates of casace levels due to presence of no-rain areas. ‘auto’=use the method implemented in
pysteps.noise.utils. compute_noise_stddev_adjs(). ‘fixed’= use the formula given in [BPS06] (eq. 6), None=disable noise std. dev adjustment.ar_order (int, optional) – The order of the autoregressive model to use. Must be 1, since only this order is currently implemented.
callback (function, optional) – Optional function that is called after computation of each time step of the nowcast. The function takes one argument: a three-dimensional array of shape (n_ens_members,h,w), where h and w are the height and width of the input field precip, respectively. This can be used, for instance, writing the outputs into files.
return_output (bool, optional) – Set to False to disable returning the outputs as numpy arrays. This can save memory if the intermediate results are written to output files using the callback function.
num_workers (int, optional) – The number of workers to use for parallel computation. Applicable if dask is enabled or pyFFTW is used for computing the FFT. When num_workers>1, it is advisable to disable OpenMP by setting the environment variable OMP_NUM_THREADS to 1. This avoids slowdown caused by too many simultaneous threads.
fft_method (str, optional) – A string defining the FFT method to use (see FFT methods in
pysteps.utils.interface.get_method()). Defaults to ‘numpy’ for compatibility reasons. If pyFFTW is installed, the recommended method is ‘pyfftw’.domain ({"spatial", "spectral"}) – If “spatial”, all computations are done in the spatial domain (the classical STEPS model). If “spectral”, the AR(2) models and stochastic perturbations are applied directly in the spectral domain to reduce memory footprint and improve performance [PCH19].
extrap_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the extrapolation method. See the documentation of
pysteps.extrapolation.interface().filter_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the filter method. See the documentation of
pysteps.cascade.bandpass_filters.noise_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the initializer of the noise generator. See the documentation of
pysteps.noise. fftgenerators.combination_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the initializer of the ensemble Kalman filter method. See the documentation of
pysteps.blending.ens_kalman_filter_methods.measure_time (bool) – If set to True, measure, print and return the computation time.
- Returns:
out – If return_output is True, a four-dimensional array of shape (n_ens_members,num_timesteps,m,n) containing a time series of forecast precipitation fields for each ensemble member. Otherwise, a None value is returned. The time series starts from t0. The timestep is taken from the input precipitation fields precip.
- Return type:
np.ndarray
See also
pysteps.extrapolation.interface,pysteps.cascade.interface,pysteps.noise.interface,pysteps.noise.utils. compute_noise_stddev_adjs()References
[NFL+19]
Notes
1. The combination method currently supports only an AR(1) process for the forecast step.