Example on retrieving and plotting winds on a distributed clusterΒΆ

This is a simple example for how to retrieve winds using the nested grid features of PyDDA.

Author: Robert C. Jackson

  • ../../_images/sphx_glr_plot_examples_nested_001.png
  • ../../_images/sphx_glr_plot_examples_nested_002.png
  • ../../_images/sphx_glr_plot_examples_nested_003.png

Out:

/home/travis/miniconda3/envs/testenv/lib/python3.6/site-packages/pyart/io/cfradial.py:384: UserWarning: WARNING: valid_min not used since it
cannot be safely cast to variable data type
  data = self.ncvar[:]
/home/travis/miniconda3/envs/testenv/lib/python3.6/site-packages/pyart/io/cfradial.py:384: UserWarning: WARNING: valid_max not used since it
cannot be safely cast to variable data type
  data = self.ncvar[:]
LocalCluster('tcp://127.0.0.1:33523', workers=2, threads=2, memory=7.84 GB)
<Client: 'tcp://127.0.0.1:33523' processes=2 threads=2, memory=7.84 GB>
/home/travis/build/openradar/PyDDA/pydda/retrieval/angles.py:24: RuntimeWarning: invalid value encountered in arccos
  elev = np.arccos((Re**2 + slantrsq - rh**2)/(2 * Re * slantr))
Calculating weights for radars 0 and 1
/home/travis/build/openradar/PyDDA/pydda/retrieval/wind_retrieve.py:544: RuntimeWarning: invalid value encountered in arccos
  theta_2 = np.arccos((x-rad2[1])/b)
Calculating weights for models...
Starting solver
rmsVR = 6.763780749907824
Total points:23059.0
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   7.3255|  31.7513|   0.0000|   0.0000|   0.0000|   0.0000|  11.4154
Norm of gradient: 0.0650530360051775
Iterations before filter: 10
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.7352|  14.2164|   0.0000|   0.0000|   0.0000|   0.0000|  25.2106
Norm of gradient: 0.024580790779265987
Iterations before filter: 20
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.6733|   9.3867|   0.0000|   0.0000|   0.0000|   0.0000|  48.6620
Norm of gradient: 0.026013982166250275
Iterations before filter: 30
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.2137|   7.3447|   0.0000|   0.0000|   0.0000|   0.0000|  56.3472
Norm of gradient: 0.02548859675514071
Iterations before filter: 40
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0640|   6.2785|   0.0000|   0.0000|   0.0000|   0.0000|  55.0991
Norm of gradient: 0.007875486945673972
Iterations before filter: 50
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0618|   5.6034|   0.0000|   0.0000|   0.0000|   0.0000|  47.7570
Norm of gradient: 0.007646483245098699
Iterations before filter: 60
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0476|   5.0452|   0.0000|   0.0000|   0.0000|   0.0000|  46.2112
Norm of gradient: 0.015678441265664804
Iterations before filter: 70
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0446|   4.7581|   0.0000|   0.0000|   0.0000|   0.0000|  41.1122
Norm of gradient: 0.013913955624800934
Iterations before filter: 80
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0601|   4.4503|   0.0000|   0.0000|   0.0000|   0.0000|  38.1186
Norm of gradient: 0.024220033810314595
Iterations before filter: 90
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.1290|   4.2479|   0.0000|   0.0000|   0.0000|   0.0000|  35.5431
Norm of gradient: 0.04007884802140612
Iterations before filter: 100
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0396|   4.2146|   0.0000|   0.0000|   0.0000|   0.0000|  35.0184
Norm of gradient: 0.013595536310328973
Iterations before filter: 110
| Jvel    | Jmass   | Jsmooth |   Jbg   | Jvort   | Jmodel  | Max w
|   0.0167|   4.2138|   0.0000|   0.0000|   0.0000|   0.0000|  35.0012
Norm of gradient: 0.006098645497136905
Iterations before filter: 120
Applying low pass filter to wind field...
Iterations after filter: 1
Iterations after filter: 2
Done! Time = 129.7
Waiting for nested grid to be retrieved...
/home/travis/build/openradar/PyDDA/pydda/vis/barb_plot.py:175: UserWarning: linewidths is ignored by contourf
  alpha=contour_alpha)
/home/travis/build/openradar/PyDDA/pydda/retrieval/wind_retrieve.py:544: RuntimeWarning: invalid value encountered in arccos
  theta_2 = np.arccos((x-rad2[1])/b)
/home/travis/build/openradar/PyDDA/pydda/vis/barb_plot.py:214: UserWarning: The following kwargs were not used by contour: 'color'
  levels=[bca_min, bca_max], color='k')
/home/travis/build/openradar/PyDDA/pydda/retrieval/wind_retrieve.py:544: RuntimeWarning: invalid value encountered in arccos
  theta_2 = np.arccos((x-rad2[1])/b)
/home/travis/build/openradar/PyDDA/pydda/vis/barb_plot.py:214: UserWarning: The following kwargs were not used by contour: 'color'
  levels=[bca_min, bca_max], color='k')
/home/travis/build/openradar/PyDDA/pydda/vis/barb_plot.py:637: UserWarning: linewidths is ignored by contourf
  alpha=contour_alpha)
/home/travis/build/openradar/PyDDA/pydda/vis/barb_plot.py:825: UserWarning: linewidths is ignored by contourf
  alpha=contour_alpha)

import pyart
import pydda
from matplotlib import pyplot as plt
from distributed import LocalCluster, Client

# Needed so that distributed doesn't run all of your code when the worker
# starts!
if __name__ == '__main__':

    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='corrected_velocity')

    # Start our dask distributed cluster. You can use any distributed cluster
    # for this...a LocalCluster is used here for the sake of being able to run
    # this example locally.
    cluster = LocalCluster(n_workers=2)
    print(cluster)
    client = Client(cluster)
    print(client)

    # Start the wind retrieval. This example only uses the mass continuity
    # and data weighting constraints.
    Grids = pydda.retrieval.get_dd_wind_field_nested(
        [berr_grid, cpol_grid], u_init,  v_init, w_init, client, Co=1.0,
        Cm=1500.0, Cz=0, 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='reflectivity',
                                        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='reflectivity',
                                     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='reflectivity',
                                     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: ( 6 minutes 52.250 seconds)

Gallery generated by Sphinx-Gallery