Example on retrieving and plotting windsΒΆ

This is a simple example for how to retrieve and plot winds from 2 radars using PyDDA.

Author: Robert C. Jackson

  • ../../_images/sphx_glr_plot_examples_001.png
  • ../../_images/sphx_glr_plot_examples_002.png
  • ../../_images/sphx_glr_plot_examples_003.png

Out:

Calculating weights for radars 0 and 1
Calculating weights for models...
Starting solver
rmsVR = 40.299436641264535
Total points:172054.0
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel | Max w
|  13.5607|  87.6626|   0.0000|   0.0000|   0.0000|   0.0000|  15.1602
Norm of gradient: 0.02024532293823782
Iterations before filter: 10
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel | Max w
|   8.5925|  90.8552|   0.0000|   0.0000|   0.0000|   0.0000|  16.3626
Norm of gradient: 0.021533456958078654
Iterations before filter: 20
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel | Max w
|   8.5925|  90.8552|   0.0000|   0.0000|   0.0000|   0.0000|  16.3626
Norm of gradient: 0.021533456958078654
Iterations before filter: 30
Applying low pass filter to wind field...
Iterations after filter: 1
Iterations after filter: 2
Done! Time = 91.7

import pyart
import pydda
from matplotlib import pyplot as plt


berr_grid = pyart.io.read_grid(pydda.tests.EXAMPLE_RADAR0)
cpol_grid = pyart.io.read_grid(pydda.tests.EXAMPLE_RADAR1)

sounding = pyart.io.read_arm_sonde(
    pydda.tests.SOUNDING_PATH)


# Load sounding data and insert as an intialization
u_init, v_init, w_init = pydda.initialization.make_wind_field_from_profile(
        cpol_grid, sounding[1], vel_field='VT')

# Start the wind retrieval. This example only uses the mass continuity
# and data weighting constraints.
Grids = pydda.retrieval.get_dd_wind_field([berr_grid, cpol_grid], u_init,
                                          v_init, w_init, Co=10.0, Cm=1500.0,
                                          Cz=0, vel_name='VT', refl_field='DT',
                                          frz=5000.0, filt_iterations=2,
                                          mask_outside_opt=True, upper_bc=1)
# Plot a horizontal cross section
plt.figure(figsize=(9, 9))
pydda.vis.plot_horiz_xsection_barbs(Grids, background_field='DT', level=6,
                                    w_vel_contours=[3, 6, 9, 12, 15],
                                    barb_spacing_x_km=5.0,
                                    barb_spacing_y_km=15.0)
plt.show()

# Plot a vertical X-Z cross section
plt.figure(figsize=(9, 9))
pydda.vis.plot_xz_xsection_barbs(Grids, background_field='DT', level=40,
                                 w_vel_contours=[3, 6, 9, 12, 15],
                                 barb_spacing_x_km=10.0,
                                 barb_spacing_z_km=2.0)
plt.show()

# Plot a vertical Y-Z cross section
plt.figure(figsize=(9, 9))
pydda.vis.plot_yz_xsection_barbs(Grids, background_field='DT', level=40,
                                 w_vel_contours=[3, 6, 9, 12, 15],
                                 barb_spacing_y_km=10.0,
                                 barb_spacing_z_km=2.0)
plt.show()

Total running time of the script: ( 1 minutes 32.681 seconds)

Gallery generated by Sphinx-Gallery