Vertical Profile of Reflectivity (VPR)¶
Precipitation is 3-dimensional in space. The vertical distribution of precipitation (and thus reflectivity) is typically non-uniform. As the height of the radar beam increases with the distance from the radar location (beam elevation, earth curvature), one sweep samples from different heights. The effects of the non-uniform VPR and the different sampling heights need to be accounted for if we are interested in the precipitation near the ground or in defined heights. This module is intended to provide a set of tools to account for these effects.
The first step will normally be to reference the polar volume data in a
3-dimensional Cartesian coordinate system. The three dimensional Cartesian
coordinates of the original polar volume data can be computed using
wradlib.vpr.volcoords_from_polar
.
Then, we can create regular 3-D grids in order to analyse the vertical profile
of reflectivity or rainfall intensity. For some applications you might want
to create so-called Constant Altitude Plan Position Indicators (CAPPI)
in order to make radar observations at different distances from the radar more
comparable. Basically, a CAPPI is simply one slice out of a 3-D volume grid.
Analoguous, we will refer to the elements in a three dimensional Cartesian grid
as voxels. In wradlib, you can create
CAPPIS (CAPPI
) and Pseudo CAPPIs
(PseudoCAPPI
) for different altitudes at once.
Here’s an example how a set of CAPPIs can be created from synthetic polar volume data:
import wradlib
import numpy as np
# define elevation and azimuth angles, ranges, radar site coordinates,
# projection
elevs = np.array([0.5,1.5,2.4,3.4,4.3,5.3,6.2,7.5,8.7,10,12,14,16.7,19.5])
azims = np.arange(0., 360., 1.)
ranges = np.arange(0., 120000., 1000.)
sitecoords = (14.924218,120.255547,500.)
proj = osr.SpatialReference()
proj.ImportFromEPSG(32651)
# create Cartesian coordinates corresponding the location of the
# polar volume bins
polxyz = wradlib.vpr.volcoords_from_polar(sitecoords, elevs,
azims, ranges, proj) # noqa
poldata = wradlib.vpr.synthetic_polar_volume(polxyz)
# this is the shape of our polar volume
polshape = (len(elevs),len(azims),len(ranges))
# now we define the coordinates for the 3-D grid (the CAPPI layers)
x = np.linspace(polxyz[:,0].min(), polxyz[:,0].max(), 120)
y = np.linspace(polxyz[:,1].min(), polxyz[:,1].max(), 120)
z = np.arange(500.,10500.,500.)
xyz = wradlib.util.gridaspoints(x, y, z)
gridshape = (len(x), len(y), len(z))
# create an instance of the CAPPI class and
# use it to create a series of CAPPIs
gridder = wradlib.vpr.CAPPI(polxyz, xyz, maxrange=ranges.max(),
gridshape=gridshape, ipclass=wradlib.ipol.Idw)
gridded = np.ma.masked_invalid( gridder(poldata) ).reshape(gridshape)
# plot results
levels = np.linspace(0,100,25)
wradlib.vis.plot_max_plan_and_vert(x, y, z, gridded, levels=levels,
cmap=pl.cm.viridis)
volcoords_from_polar |
Create Cartesian coordinates for regular polar volumes |
make_3d_grid |
Generate Cartesian coordinates for a regular 3-D grid based on radar specs. |
CartesianVolume |
Create 3-D regular volume grid in Cartesian coordinates from polar data with multiple elevation angles |
CAPPI |
Create a Constant Altitude Plan Position Indicator (CAPPI) |
PseudoCAPPI |
Create a Pseudo-CAPPI Constant Altitude Plan Position Indicator (CAPPI) |