swirlspy.qpf package

Quantitative precipitation forecast methods.

swirlspy.qpf.constant(frames) DataArray

Wrapper of pysteps.motion.constant.constant.

Parameters

frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

Returns

motion – Three-dimensional array of shape (Dimensions, m, n) containing the dense x- and y-components of the motion field in units of pixels / timestep as given by frames. Dimensions - [u, v].

Return type

xarray.DataArray

swirlspy.qpf.darts(frames, **kwargs) DataArray

Wrapper of pysteps.motion.darts.DARTS.

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • N_x (int) – Number of DFT coefficients to use for the input images, x-axis (default=50).

  • N_y (int) – Number of DFT coefficients to use for the input images, y-axis (default=50).

  • N_t (int) – Number of DFT coefficients to use for the input images, time axis (default=4). N_t must be strictly smaller than T.

  • M_x (int) – Number of DFT coefficients to compute for the output advection field, x-axis (default=2).

  • M_y (int) – Number of DFT coefficients to compute for the output advection field, y-axis (default=2).

  • fft_method (str) – A string defining the FFT method to use, see utils.fft.get_method. Defaults to ‘numpy’.

  • output_type ({"spatial", "spectral"}) – The type of the output: “spatial”=apply the inverse FFT to obtain the spatial representation of the advection field, “spectral”=return the (truncated) DFT representation.

  • n_threads (int) – Number of threads to use for the FFT computation. Applicable if fft_method is ‘pyfftw’.

  • print_info (bool) – If True, print information messages.

  • lsq_method ({1, 2}) – The method to use for solving the linear equations in the least squares sense: 1=numpy.linalg.lstsq, 2=explicit computation of the Moore-Penrose pseudoinverse and SVD.

  • verbose (bool) – if set to True, it prints information about the program

Returns

motion – Three-dimensional array of shape (Dimensions, m, n) containing the dense x- and y-components of the motion field in units of pixels / timestep as given by frames. Dimensions - [u, v].

Return type

xarray.DataArray

swirlspy.qpf.dense_lucaskanade(frames, **kwargs) DataArray

Wrapper of pysteps.motion.lucaskanade.dense_lucaskanade.

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • buffer_mask (int, optional) – A mask buffer width in pixels. This extends the input mask (if any) to help avoiding the erroneous interpretation of velocities near the maximum range of the radars (0 by default).

  • max_corners_ST (int, optional) – The maxCorners parameter in the `Shi-Tomasi`_ corner detection method. It represents the maximum number of points to be tracked (corners), by default this is 500. If set to zero, all detected corners are used.

  • quality_level_ST (float, optional) – The qualityLevel parameter in the `Shi-Tomasi`_ corner detection method. It represents the minimal accepted quality for the points to be tracked (corners), by default this is set to 0.1. Higher quality thresholds can lead to no detection at all.

  • min_distance_ST (int, optional) – The minDistance parameter in the `Shi-Tomasi`_ corner detection method. It represents minimum possible Euclidean distance in pixels between corners, by default this is set to 3 pixels.

  • block_size_ST (int, optional) – The blockSize parameter in the `Shi-Tomasi`_ corner detection method. It represents the window size in pixels used for computing a derivative covariation matrix over each pixel neighborhood, by default this is set to 15 pixels.

  • winsize_LK (tuple of int, optional) – The winSize parameter in the `Lucas-Kanade`_ optical flow method. It represents the size of the search window that it is used at each pyramid level, by default this is set to (50, 50) pixels.

  • nr_levels_LK (int, optional) – The maxLevel parameter in the `Lucas-Kanade`_ optical flow method. It represents the 0-based maximal pyramid level number, by default this is set to 3.

  • nr_std_outlier (int, optional) – Maximum acceptable deviation from the mean/median in terms of number of standard deviations. Any anomaly larger than this value is flagged as outlier and excluded from the interpolation. By default this is set to 3.

  • multivariate_outlier (bool, optional) – If true (the default), the outlier detection is computed in terms of the Mahalanobis distance. If false, the outlier detection is simply computed in terms of velocity.

  • k_outlier (int, optional) – The number of nearest neighbours used to localize the outlier detection. If set equal to 0, it employs all the data points. The default is 30.

  • size_opening (int, optional) – The size of the structuring element kernel in pixels. This is used to perform a binary morphological opening on the input fields in order to filter isolated echoes due to clutter. By default this is set to 3. If set to zero, the fitlering is not perfomed.

  • decl_grid (int, optional) – The cell size in pixels of the declustering grid that is used to filter out outliers in a sparse motion field and get more representative data points before the interpolation. This simply computes new sparse vectors over a coarser grid by taking the median of all vectors within one cell. By default this is set to 20 pixels. If set to less than 2 pixels, the declustering is not perfomed.

  • min_nr_samples (int, optional) – The minimum number of samples necessary for computing the median vector within given declustering cell, otherwise all sparse vectors in that cell are discarded. By default this is set to 2.

  • rbfunction (string, optional) – The name of the radial basis function used for the interpolation of the sparse vectors. This is based on the Euclidian norm d. By default this is set to “inverse” and the available names are “nearest”, “inverse”, “gaussian”.

  • k (int, optional) – The number of nearest neighbours used for fast interpolation, by default this is set to 20. If set equal to zero, it employs all the neighbours.

  • epsilon (float, optional) – The adjustable constant used in the gaussian and inverse radial basis functions. by default this is computed as the median distance between the sparse vectors.

  • nchunks (int, optional) – Split the grid points in n chunks to limit the memory usage during the interpolation. By default this is set to 5, if set to 1 the interpolation is computed with the whole grid.

  • extra_vectors (numpy.ndarray, optional) – Additional sparse motion vectors as 2d array (columns: x,y,u,v; rows: nbr. of vectors) to be integrated with the sparse vectors from the Lucas-Kanade local tracking. x and y must be in pixel coordinates, with (0,0) being the upper-left corner of the field R. u and v must be in pixel units. By default this is set to None.

  • verbose (bool, optional) – If set to True, it prints information about the program (True by default).

Returns

motion – Three-dimensional array of shape (Dimensions, m, n) containing the dense x- and y-components of the motion field in units of pixels / timestep as given by frames. Dimensions - [u, v].

Return type

xarray.DataArray

swirlspy.qpf.extrapolation(frames, motion, nowcast_steps, **kwargs) DataArray

Wrapper of pysteps.nowcasts.extrapolation.forecast (semi-Lagrangian).

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • motion (xarray.DataArray) – x, y component of motion field. Dimensions - [components, y, x]. components - [u, v].

  • nowcast_steps (int) – Number of time steps to forecast.

  • outval (float, optional) – Optional argument for specifying the value for pixels advected from outside the domain. If outval is set to ‘min’, the value is taken as the minimum value of R. Default : np.nan

  • xy_coords (numpy.ndarray, optional) – Array with the coordinates of the grid dimension (2, m, n ). * xy_coords[0] : x coordinates * xy_coords[1] : y coordinates By default, the xy_coords are computed for each extrapolation.

  • allow_nonfinite_values (bool, optional) – If True, allow non-finite values in the precipitation and advection fields. This option is useful if the input fields contain a radar mask (i.e. pixels with no observations are set to nan).

  • D_prev (array-like) – Optional initial displacement vector field of shape (2,m,n) for the extrapolation. Default : None

  • n_iter (int) – Number of inner iterations in the semi-Lagrangian scheme. If n_iter > 0, the integration is done using the midpoint rule. Otherwise, the advection vectors are taken from the starting point of each interval. Default : 1

Returns

Three-dimensional array of shape (timesteps, m, n) containing a time series of nowcast precipitation fields. The time series starts from base time, where timestep is taken from the advection field velocity.

Return type

xarray.DataArray

swirlspy.qpf.persistence(frames, nowcast_steps) DataArray

Forecast by persistence, accepting an observation xarray and returning a forecast xarray. Persistence means that the forecast is exactly same as the observation. While it sounds trivial, this is usefulf for benchmarking the performance of another QPF algorithm.

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • nowcast_steps (int) – Number of time steps to forecast.

Returns

Three-dimensional array of shape (timesteps, m, n) containing a time series of nowcast precipitation fields. The time series starts from base time, where timestep is taken from the advection field velocity.

Return type

xarray.DataArray

swirlspy.qpf.rover(frames, **kwargs) DataArray

Running rover to generate u and v motion fields for QPF.

Parameters
  • frames (xarray.DataArray) – Contains two continuous frames data. Dimensions - [time, y, x].

  • start_level (int) – Starting level used as base in multigrid algorithm. Higher=coarser. Defaults to 1.

  • max_level (int) – Maximum level to be reached in multigrid algorithm. Higher=coarser. Defaults to 7.

  • rho (float) – Gaussian smoothing parameter. Defaults to 9.0.

  • alpha (float) – Regularisation parameter in the energy functional. Defaults to 2000.0.

  • sigma (float) – Gaussian smoothing parameter. Defaults to 1.5.

  • atan_dBZ_min (float) – Minimum arc tangent value. Defaults to -1.482.

  • atan_dBZ_max (float) – Maximum arc tangent value. Defaults to 1.412.

  • ashift (float) – Point of inflection. Defaults to 33 dBZ.

  • afact (float) – Sharpness of inflection. Defaults to 4.0 dBZ.

Returns

motion – Three-dimensional array of shape (Dimensions, m, n) containing the dense x- and y-components of the motion field in units of pixels / timestep as given by frames. Dimensions - [u, v].

Return type

xarray.DataArray

Notes

Default values of start_level, max_level, rho, alpha and sigma are those of ROVER-A.

swirlspy.qpf.sla(frames: DataArray, motion: DataArray, nowcast_steps: int) DataArray

Calculates forecast reflectivity using Semi-Lagrangian Advection from horizontal and vertical motion fields.

Parameters
  • frames (xarray.DataArray) – Contains two continuous frames data. Dimensions - [time, y, x].

  • motion (xarray.DataArray) – x, y component of motion field. Dimensions - [components, y, x]. components - [u, v].

Returns

forecast – A xarray.DataArray populated with forecast reflectivity and metadata.

Return type

xarray.DataArray

swirlspy.qpf.sprog(frames, motion, nowcast_steps, **kwargs) DataArray

Wrapper of pysteps.nowcasts.sprog.forecast (Spectral Prognosis (S-PROG)).

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • motion (xarray.DataArray) – x, y component of motion field. Dimensions - [components, y, x]. components - [u, v].

  • nowcast_steps (int) – Number of time steps to forecast.

  • n_cascade_levels (int, optional) – The number of cascade levels to use.

  • R_thr (float, optional) – Specifies the threshold value for minimum observable precipitation intensity. Required if mask_method is not None or conditional is True.

  • decomp_method ({'fft'}, optional) – Name of the cascade decomposition method to use. See the documentation of pysteps.cascade.interface.

  • 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.

  • ar_order (int, optional) – The order of the autoregressive model to use. Must be >= 1.

  • conditional (bool, optional) – If set to True, compute the statistics of the precipitation field conditionally by excluding pixels where the values are below the threshold R_thr.

  • probmatching_method ({'cdf','mean',None}, optional) – Method for matching the conditional statistics of the forecast field (areas with precipitation intensity above the threshold R_thr) with those of the most recently observed one. ‘cdf’=map the forecast CDF to the observed one, ‘mean’=adjust only the mean value, None=no matching applied.

  • 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 utils.fft.get_method). Defaults to ‘numpy’ for compatibility reasons. If pyFFTW is installed, the recommended method is ‘pyfftw’.

  • extrap_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the extrapolation method. See the documentation of pysteps.extrapolation.

  • filter_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the filter method. See the documentation of pysteps.cascade.bandpass_filters.py.

Returns

Three-dimensional array of shape (timesteps, m, n) containing a time series of nowcast precipitation fields. The time series starts from base time, where timestep is taken from the advection field velocity.

Return type

xarray.DataArray

swirlspy.qpf.sseps(frames, motion, nowcast_steps, **kwargs) DataArray

Wrapper of pysteps.nowcasts.sseps.forecast (Short-space ensemble prediction system (SSEPS)).

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • motion (xarray.DataArray) – x, y component of motion field. Dimensions - [components, y, x]. components - [u, v].

  • nowcast_steps (int) – Number of time steps to forecast.

  • win_size (int or two-element sequence of ints) – Size-length of the localization window.

  • overlap (float [0,1[) – A float between 0 and 1 prescribing the level of overlap between successive windows. If set to 0, no overlap is used.

  • war_thr (float) – Threshold for the minimum fraction of rain in a given window.

  • n_ens_members (int) – The number of ensemble members to generate.

  • n_cascade_levels (int) – The number of cascade levels to use.

  • extrap_method ({'semilagrangian'}) – Name of the extrapolation method to use. See the documentation of pysteps.extrapolation.interface.

  • decomp_method ({'fft'}) – Name of the cascade decomposition method to use. See the documentation of pysteps.cascade.interface.

  • bandpass_filter_method ({'gaussian', 'uniform'}) – Name of the bandpass filter method to use with the cascade decomposition.

  • noise_method ({'parametric','nonparametric','ssft','nested',None}) – 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.

  • ar_order (int) – The order of the autoregressive model to use. Must be >= 1.

  • vel_pert_method ({'bps',None}) – Name of the noise generator to use for perturbing the advection field. See the documentation of pysteps.noise.interface. If set to None, the advection field is not perturbed.

  • mask_method ({'incremental', None}) – The method to use for masking no precipitation areas in the forecast field. The masked pixels are set to the minimum value of the observations. ‘incremental’ = iteratively buffer the mask with a certain rate (currently it is 1 km/min), None=no masking.

  • probmatching_method ({'cdf', None}) – Method for matching the statistics of the forecast field with those of the most recently observed one. ‘cdf’=map the forecast CDF to the observed one, None=no matching applied. Using ‘mean’ requires that mask_method is not None.

  • callback (function) – 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 R, respectively. This can be used, for instance, writing the outputs into files.

  • return_output (bool) – 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.

  • seed (int) – Optional seed number for the random generators.

  • num_workers (int) – 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) – A string defining the FFT method to use (see utils.fft.get_method). Defaults to ‘numpy’ for compatibility reasons. If pyFFTW is installed, the recommended method is ‘pyfftw’.

  • extrap_kwargs (dict) – Optional dictionary containing keyword arguments for the extrapolation method. See the documentation of pysteps.extrapolation.

  • filter_kwargs (dict) – Optional dictionary containing keyword arguments for the filter method. See the documentation of pysteps.cascade.bandpass_filters.py.

  • noise_kwargs (dict) – Optional dictionary containing keyword arguments for the initializer of the noise generator. See the documentation of pysteps.noise.fftgenerators.

  • vel_pert_kwargs (dict) – Optional dictionary containing keyword arguments “p_pert_par” and “p_pert_perp” for the initializer of the velocity perturbator. See the documentation of pysteps.noise.motion.

  • mask_kwargs (dict) – Optional dictionary containing mask keyword arguments ‘mask_f’ and ‘mask_rim’, the factor defining the the mask increment and the rim size, respectively. The mask increment is defined as mask_f*timestep/kmperpixel.

  • measure_time (bool) – If set to True, measure, print and return the computation time.

Returns

Three-dimensional array of shape (timesteps, m, n) containing a time series of nowcast precipitation fields. The time series starts from base time, where timestep is taken from the advection field velocity.

Return type

xarray.DataArray

swirlspy.qpf.steps(frames, motion, nowcast_steps, statistical_operations=None, **kwargs) DataArray

Wrapper of pysteps.nowcasts.steps.forecast (Short-Term Ensemble Prediction System (STEPS)).

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • motion (xarray.DataArray) – x, y component of motion field. Dimensions - [components, y, x]. components - [u, v].

  • nowcast_steps (int) – Number of time steps to forecast.

  • statistical_operations ({'mean', 'median', 'max', 'min'}, optional) – The statistics method applied to ensemble members.

  • n_ens_members (int, optional) – The number of ensemble members to generate.

  • n_cascade_levels (int, optional) – The number of cascade levels to use.

  • R_thr (float, optional) – Specifies the threshold value for minimum observable precipitation intensity. Required if mask_method is not None or conditional is True.

  • kmperpixel (float, optional) – Spatial resolution of the input data (kilometers/pixel). Required if vel_pert_method is not None or mask_method is ‘incremental’.

  • timestep (float, optional) – Time step of the motion vectors (minutes). Required if vel_pert_method is not None or mask_method is ‘incremental’.

  • extrap_method (str, optional) – Name of the extrapolation method to use. See the documentation of pysteps.extrapolation.interface.

  • decomp_method ({'fft'}, optional) – Name of the cascade decomposition method to use. See the documentation of pysteps.cascade.interface.

  • 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.

  • 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.

  • 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 :cite:`BPS2006` (eq. 6), None=disable noise std. dev adjustment.

  • ar_order (int, optional) – The order of the autoregressive model to use. Must be >= 1.

  • vel_pert_method ({'bps',None}, optional) – Name of the noise generator to use for perturbing the advection field. See the documentation of pysteps.noise.interface. If set to None, the advection field is not perturbed.

  • conditional (bool, optional) – If set to True, compute the statistics of the precipitation field conditionally by excluding pixels where the values are below the threshold R_thr.

  • mask_method ({'obs','sprog','incremental',None}, optional) – The method to use for masking no precipitation areas in the forecast field. The masked pixels are set to the minimum value of the observations. ‘obs’ = apply R_thr to the most recently observed precipitation intensity field, ‘sprog’ = use the smoothed forecast field from S-PROG, where the AR(p) model has been applied, ‘incremental’ = iteratively buffer the mask with a certain rate (currently it is 1 km/min), None=no masking.

  • probmatching_method ({'cdf','mean',None}, optional) – Method for matching the statistics of the forecast field with those of the most recently observed one. ‘cdf’=map the forecast CDF to the observed one, ‘mean’=adjust only the conditional mean value of the forecast field in precipitation areas, None=no matching applied. Using ‘mean’ requires that mask_method is not None.

  • 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 R, 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.

  • seed (int, optional) – Optional seed number for the random generators.

  • 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 utils.fft.get_method). Defaults to ‘numpy’ for compatibility reasons. If pyFFTW is installed, the recommended method is ‘pyfftw’.

  • extrap_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the extrapolation method. See the documentation of pysteps.extrapolation.

  • filter_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the filter method. See the documentation of pysteps.cascade.bandpass_filters.py.

  • noise_kwargs (dict, optional) – Optional dictionary containing keyword arguments for the initializer of the noise generator. See the documentation of pysteps.noise.fftgenerators.

  • vel_pert_kwargs (dict, optional) – Optional dictionary containing keyword arguments ‘p_par’ and ‘p_perp’ for the initializer of the velocity perturbator. The choice of the optimal parameters depends on the domain and the used optical flow method. See pysteps.noise.motion for additional documentation.

  • mask_kwargs (dict) – Optional dictionary containing mask keyword arguments ‘mask_f’ and ‘mask_rim’, the factor defining the the mask increment and the rim size, respectively. The mask increment is defined as mask_f*timestep/kmperpixel.

  • measure_time (bool) – If set to True, measure, print and return the computation time.

Returns

Three-dimensional array of shape (timesteps, m, n) containing a time series of nowcast precipitation fields. The time series starts from base time, where timestep is taken from the advection field velocity.

Return type

xarray.DataArray

swirlspy.qpf.vet(frames, **kwargs) DataArray

Wrapper of pysteps.motion.vet.vet.

Parameters
  • frames (xarray.DataArray) – Contains continuous frames data. Dimensions - [time, y, x].

  • sectors (list or array, optional) – The number of sectors for each dimension used in the scaling procedure. If dimension is 1, the same sectors will be used both image dimensions (x and y). If is 2D, the each row determines the sectors of the each dimension.

  • smooth_gain (float, optional) – Smooth gain factor

  • first_guess (numpy.ndarray, optional) –

    The shape of the first guess should have the same shape as the initial sectors shapes used in the scaling procedure. If first_guess is not present zeros are used as first guess. E.g.:

    If the first sector shape in the scaling procedure is (ni,nj), then the first_guess should have (2, ni, nj ) shape.

  • intermediate_steps (bool, optional) – If True, also return a list with the first guesses obtained during the scaling procedure. False, by default.

  • verbose (bool, optional) – Verbosity enabled if True (default).

  • indexing (str, optional) – Input indexing order.’ij’ and ‘xy’ indicates that the dimensions of the input are (time, longitude, latitude), while ‘yx’ indicates (time, latitude, longitude). The displacement field dimensions are ordered accordingly in a way that the first dimension indicates the displacement along x (0) or y (1). That is, UV displacements are always returned.

  • padding (int) – Padding width in grid points. A border is added to the input array to reduce the effects of the minimization at the border.

  • options (dict, optional) – A dictionary of solver options. See `scipy minimization`_ function for more details.

Returns

motion – Three-dimensional array of shape (Dimensions, m, n) containing the dense x- and y-components of the motion field in units of pixels / timestep as given by frames. Dimensions - [u, v].

Return type

xarray.DataArray

Subpackages

Submodules

swirlspy.qpf.persistence module

swirlspy.qpf.persistence.persistence(xr2d, timesteps, timesteps_size)

Forecast by persistence, accepting an observation xarray and returning a forecast xarray. Persistence means that the forecast is exactly same as the observation. While it sounds trivial, this is usefulf for benchmarking the performance of another QPF algorithm.

Parameters
  • xr2d (xarray.DataArray (2-dimensional, in 'x' 'y')) – The xarray containing data to be forecast

  • timesteps (int) – Number of timesteps to be forecast

  • timesteps_size (int) – Time interval between each steps, in minutes.

Returns

persistence_xrnd – Data values in ‘xr2d’ being copied along the ‘end_time’ axis, coordinates of ‘end_time’ axis have a unique value for each timestep.

Return type

xarray.DataArray

swirlspy.qpf.rover module

swirlspy.qpf.rover.rover(qpf_xarray, start_level=1, max_level=7, rho=9.0, alpha=2000.0, sigma=1.5, gray_scale=True, atan_dBZ_min=-1.482, atan_dBZ_max=1.412, ashift=33.0, afact=4.0, zero_value=None)

Running rover to generate u and v motion fields for QPF. fall back wrapper

Parameters
  • qpf_xarray (xarray.DataArray) – Contains reflectivity data. First element of array corresponds to the geographic upper left corner.

  • start_level (int) – Starting level used as base in multigrid algorithm. Higher=coarser. Defaults to 1.

  • max_level (int) – Maximum level to be reached in multigrid algorithm. Higher=coarser. Defaults to 7.

  • rho (float) – Gaussian smoothing parameter. Defaults to 9.0.

  • alpha (float) – Regularisation parameter in the energy functional. Defaults to 2000.0.

  • sigma (float) – Gaussian smoothing parameter. Defaults to 1.5.

  • atan_dBZ_min (float) – Minimum arc tangent value. Defaults to -1.482.

  • atan_dBZ_max (float) – Maximum arc tangent value. Defaults to 1.412.

  • ashift (float) – Point of inflection. Defaults to 33 dBZ.

  • afact (float) – Sharpness of inflection. Defaults to 4.0 dBZ.

  • zero_value (int) – Value that assigned to be no data. Default: 9999 (If zero_value is None, numpy.nan will be used)

Returns

  • motion_u (xarray.DataArray) – x component of motion field.

  • motion_v (xarray.DataArray) – y component of motion field.

Notes

Default values of start_level, max_level, rho, alpha and sigma are those of ROVER-A.

swirlspy.qpf.sla module

swirlspy.qpf.sla.sla(qpf_xarray, motion_u, motion_v, steps=6, zero_value=9999)

Calculates forecast reflectivity using Semi-Lagrangian Advection from horizontal and vertical motion fields. fall back wrapper

Parameters
  • qpf_xarray (xarray.DataArray) – Contains reflectivity data. First element of array corresponds to the geographic upper left corner.

  • motion_u (xarray.DataArray) – x component of motion field.

  • motion_v (xarray.DataArray) – y component of motion field.

  • steps (int, optional) – Number of time steps to forecast. Default: 6

  • zero_value (int) – Value that assigned to be no data. Default: 9999 (If zero_value is None, numpy.nan will be used)

Returns

forecast – A xarray.DataArray populated with forecast reflectivity and metadata.

Return type

xarray.DataArray

Notes

The program will automatically use GPU to run extrapolation if your computer has a CUDA-capable GPU and has installed a NVIDIA driver.

Please follow CUDA Installation Instruction to install NVIDIA driver. If you use Ubuntu 18.04, execute the below commands to install the driver:

sudo apt-get install linux-headers-$(uname -r)
wget -P /tmp/sla_gpu https://developer.download.nvidia.com/compute/cuda/repos/ubuntu1804/x86_64/cuda-repo-ubuntu1804_10.2.89-1_amd64.deb
sudo dpkg -i /tmp/sla_gpu/cuda-repo-ubuntu1804_10.2.89-1_amd64.deb
wget -P /tmp/sla_gpu https://developer.download.nvidia.com/compute/cuda/repos/ubuntu1804/x86_64/cuda-ubuntu1804.pin
sudo mv /tmp/sla_gpu/cuda-ubuntu1804.pin /etc/apt/preferences.d/cuda-repository-pin-600
sudo apt-key adv --fetch-keys https://developer.download.nvidia.com/compute/cuda/repos/ubuntu1804/x86_64/7fa2af80.pub
sudo add-apt-repository "deb http://developer.download.nvidia.com/compute/cuda/repos/ubuntu1804/x86_64/ /"
sudo apt-get install cuda-drivers -y

After installing the driver, regardless of which linux distribution you are using, you need to reboot your computer for the driver to function properly.