{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Dealias doppler velocities for an RHI file\n\nIn this example doppler velocities from range height indicator (RHI) scans are dealiased,\nusing custom parameters with the region-based dealiasing algorithm in Py-ART.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print(__doc__)\n\n# Author: Max Grover (mgrover@anl.gov)\n\nimport matplotlib.pyplot as plt\nfrom open_radar_data import DATASETS\n\nimport pyart\n\n# Read data from the open-radar-data package\nfile = DATASETS.fetch(\"cfrad.20211011_223602.712_to_20211011_223612.091_DOW8_RHI.nc\")\nradar = pyart.io.read(file)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Apply the default parameters for dealiasing**\n\nNow let's apply the default region-based technique\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Retrieve the nyquist velocity from the first sweep\nnyq = radar.instrument_parameters[\"nyquist_velocity\"][\"data\"][0]\n\n# Apply the dealiasing algorithm with default parameters\nvelocity_dealiased = pyart.correct.dealias_region_based(\n radar,\n vel_field=\"VEL\",\n nyquist_vel=nyq,\n)\n\n# Add the field to the radar object\nradar.add_field(\"corrected_velocity\", velocity_dealiased, replace_existing=True)\n\n# Visualize the output\nfig = plt.figure(figsize=(12, 4))\nax1 = fig.add_subplot(121)\ndisplay = pyart.graph.RadarDisplay(radar)\ndisplay.plot(\n \"VEL\",\n ax=ax1,\n vmin=-30,\n vmax=30,\n cmap=\"balance\",\n title=\"Uncorrected Radial Velocity\",\n)\ncbar = display.cbs[0]\n# Modify the colorbar label and size\ncbar.set_label(label=\"Raw Radial Velocity ($V_{R}$) [m/s]\", fontsize=12)\nplt.ylim(0, 10)\nplt.xlim(20, 60)\n\nax2 = fig.add_subplot(122)\ndisplay = pyart.graph.RadarDisplay(radar)\ndisplay.plot(\n \"corrected_velocity\",\n ax=ax2,\n vmin=-30,\n vmax=30,\n cmap=\"balance\",\n title=\"Corrected Radial Velocity\",\n)\ncbar = display.cbs[0]\n# Modify the colorbar label and size\ncbar.set_label(label=\"Corrected Radial Velocity ($V_{R}$) [m/s]\", fontsize=12)\nplt.ylim(0, 10)\nplt.xlim(20, 60)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Refining the technique**\n\nIf we use the default configuration, notice how extreme\nthe values are for the output. The algorithm keys on the\nsmall variations in the clear air returns, and\nresults in much larger values than would be expected.\nIt is considered best practice to apply a base-level of\nquality control here to improve the results. We can use velocity\ntexture as a base level, where the noisier data will be removed.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Calculate the velocity texture and add it to the radar object\nvel_texture = pyart.retrieve.calculate_velocity_texture(\n radar,\n vel_field=\"VEL\",\n nyq=nyq,\n)\nradar.add_field(\"velocity_texture\", vel_texture, replace_existing=True)\n\n# Set up the gatefilter to be based on the velocity texture.\ngatefilter = pyart.filters.GateFilter(radar)\ngatefilter.exclude_above(\"velocity_texture\", 6)\n\n# Dealias with the gatefilter and add the corrected field to the radar object\nvelocity_dealiased = pyart.correct.dealias_region_based(\n radar,\n vel_field=\"VEL\",\n nyquist_vel=nyq,\n gatefilter=gatefilter,\n)\nradar.add_field(\"corrected_velocity\", velocity_dealiased, replace_existing=True)\n\n# Visualize the output\nfig = plt.figure(figsize=(12, 4))\nax1 = fig.add_subplot(121)\ndisplay = pyart.graph.RadarDisplay(radar)\ndisplay.plot(\n \"VEL\",\n ax=ax1,\n vmin=-30,\n vmax=30,\n cmap=\"balance\",\n title=\"Uncorrected Radial Velocity\",\n)\ncbar = display.cbs[0]\n\n# Modify the colorbar label and size\ncbar.set_label(label=\"Raw Radial Velocity ($V_{R}$) [m/s]\", fontsize=12)\nplt.ylim(0, 10)\nplt.xlim(20, 60)\n\nax2 = fig.add_subplot(122)\ndisplay = pyart.graph.RadarDisplay(radar)\ndisplay.plot(\n \"corrected_velocity\",\n ax=ax2,\n vmin=-30,\n vmax=30,\n cmap=\"balance\",\n title=\"Corrected Radial Velocity\",\n)\ncbar = display.cbs[0]\n\n# Modify the colorbar label and size\ncbar.set_label(label=\"Corrected Radial Velocity ($V_{R}$) [m/s]\", fontsize=12)\nplt.ylim(0, 10)\nplt.xlim(20, 60)" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.14.3" } }, "nbformat": 4, "nbformat_minor": 0 }