{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Map a single PPI sweep to a Cartesian grid using 2D weighting\n\nMap the reflectivity field of a single PPI sweep from Antenna (polar) coordinates\nto a Cartesian grid, while using a 2D weighting.\nThis solution is valid only for this case (single PPI sweep), yet it is a useful one\n(an arguably global solution since it overlooks the z-dimension grid limits).\nThe exclusion of the z-dimension from the RoI and weighting calculations results in\nminor errors, especially considering the high co-variance of neighboring radar\nvolumes and the typically small Cartesian grid separation. Errors are effectively\nnegligible at low elevation angles common to PPI sweeps.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print(__doc__)\n\n# =====================\n# Author: Israel Silber\n# =====================\n\nimport cartopy.crs as ccrs\nimport matplotlib.pyplot as plt\nfrom open_radar_data import DATASETS\n\nimport pyart\n\n# file, and fields\n# =======================\nfile_name = DATASETS.fetch(\"example_plot_ppi_single_sweep.nc\")\nprocessed_field = \"reflectivity_at_cor\"\n\n# load file\n# =======================\nradar_obj = pyart.io.read(file_name)\nprint(f\"Total sweeps = {radar_obj.nsweeps}\")\n\n# Plot polar coordinate (raw) data\n# =======================\nprint(\n \"===\\n\\nnow displaying Ze data for 2nd and 4th sweeps in polar coordinates\\n\"\n \"(both approaching the 40 km range denoted by the purple ring)\\n\\n===\"\n)\nfig = plt.figure(figsize=(12, 6), tight_layout=True)\nfig.suptitle(\"polar\")\nfor sweep, ax_ind in zip([1, 3], range(2)):\n ax = fig.add_subplot(1, 2, ax_ind + 1, projection=ccrs.Mercator())\n sweep_obj = radar_obj.extract_sweeps((sweep,))\n display = pyart.graph.RadarDisplay(sweep_obj)\n display.plot(\n processed_field,\n sweep=0,\n ax=ax,\n )\n display.plot_range_rings([10, 20, 30])\n display.plot_range_rings([40], col=\"purple\")\nplt.show()\n\nprint(\n \"===\\n\\nnow displaying gridded Ze data for 2nd and 4th (final) sweeps; note the \"\n \"truncated max range in the case of the 4th sweep\\n\\n===\"\n)\nfig2 = plt.figure(figsize=(12, 6), tight_layout=True)\nfig2.suptitle(\"Cartesian gridded\")\nfor sweep, ax_ind in zip([1, 3], range(2)):\n sweep_obj = radar_obj.extract_sweeps((sweep,))\n grid = pyart.map.grid_from_radars(\n (sweep_obj,),\n grid_shape=(1, 1601, 1601),\n grid_limits=((0, 10000.0), [-40000, 40000], [-40000, 40000]),\n fields=[processed_field],\n )\n ax = fig2.add_subplot(1, 2, ax_ind + 1)\n ax.imshow(\n grid.fields[processed_field][\"data\"][0],\n origin=\"lower\",\n extent=(-40, 40, -40, 40),\n )\n ax.set_title(f\"sweep #{sweep + 1}\")\nplt.show()\n\nprint(\n \"===\\n\\nUsing 2D weighting by \"\n \"setting h_factor and dist_factor z component to 0.0, the max range looks OK now\\n\\n===\"\n)\nfig2 = plt.figure(figsize=(12, 6), tight_layout=True)\nfig2.suptitle(\"Cartesian gridded\")\nfor sweep, ax_ind in zip([1, 3], range(2)):\n sweep_obj = radar_obj.extract_sweeps((sweep,))\n grid = pyart.map.grid_from_radars(\n (sweep_obj,),\n grid_shape=(1, 1601, 1601),\n grid_limits=((0, 10000.0), [-40000, 40000], [-40000, 40000]),\n fields=[processed_field],\n h_factor=(0.0, 1.0, 1.0),\n dist_factor=(0.0, 1.0, 1.0),\n )\n ax = fig2.add_subplot(1, 2, ax_ind + 1)\n ax.imshow(\n grid.fields[processed_field][\"data\"][0],\n origin=\"lower\",\n extent=(-40, 40, -40, 40),\n )\n ax.set_title(f\"sweep #{sweep + 1}\")\nplt.show()" ] } ], "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.13.0" } }, "nbformat": 4, "nbformat_minor": 0 }