.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/qpe_hk.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note Click :ref:`here ` to download the full example code .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_qpe_hk.py: QPE (Hong Kong) ==================================================== This example demonstrates how to perform QPE, using raingauge and radar data from Hong Kong. .. GENERATED FROM PYTHON SOURCE LINES 10-13 Definitions -------------------------------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 13-42 .. code-block:: default import os import time import tarfile import numpy as np import pandas as pd import copy import xarray import scipy import pyproj import matplotlib from matplotlib.colors import BoundaryNorm import matplotlib.pyplot as plt from PIL import Image from pyresample import utils from cartopy.io import shapereader from swirlspy.obs.rain import Rain from swirlspy.rad.irisref import read_irisref from swirlspy.core.resample import grid_resample from swirlspy.qpe.rfmap import rg_interpolate, comp_qpe, show_raingauge from swirlspy.qpe.utils import timestamps_ending, locate_file, dbz2rr, rr2rf, \ temporal_interp, acc plt.switch_backend('agg') THIS_DIR = os.getcwd() os.chdir(THIS_DIR) .. GENERATED FROM PYTHON SOURCE LINES 43-50 Initialising ----------------------------------------------------------------- This section demonstrates extracting raingauge and radar data. Step 1: Defining an end-time for accumulating rainfall. .. GENERATED FROM PYTHON SOURCE LINES 50-54 .. code-block:: default acctime = pd.Timestamp('20190420150000').floor('min') acctime_str = acctime.strftime('%Y%m%d%H%M') .. GENERATED FROM PYTHON SOURCE LINES 55-57 Step 2: Setting up directories for raingauge and radar files. .. GENERATED FROM PYTHON SOURCE LINES 57-61 .. code-block:: default rad_dir = THIS_DIR + '/../tests/samples/iris/' rg_dir = THIS_DIR + '/../tests/samples/rfmap/' .. GENERATED FROM PYTHON SOURCE LINES 62-64 Step 3: Generating timestamps and pattern for both radar and raingauge files. .. GENERATED FROM PYTHON SOURCE LINES 64-90 .. code-block:: default # Timestamps of raingauges rg_timestrings = timestamps_ending( acctime + pd.Timedelta(minutes=5), duration=pd.Timedelta(hours=1), interval=pd.Timedelta(minutes=5) ) # Raingauge pattern rg_pattern = ['rf5m_20'+ts for ts in rg_timestrings] # Finding time nearest radar file # to accumulation end time minute = acctime.minute nearest_6min = acctime.minute // 6 * 6 nearest_rad_timestamp = pd.Timestamp( acctime_str[:-2]+f'{nearest_6min:02}' ) rad_timestrings = timestamps_ending( nearest_rad_timestamp, duration=pd.Timedelta(hours=1), interval=pd.Timedelta(minutes=6) ) .. GENERATED FROM PYTHON SOURCE LINES 91-93 Step 4: Extracting raingauge and radar files from their respective directories. .. GENERATED FROM PYTHON SOURCE LINES 93-102 .. code-block:: default located_rg_files = [] for pat in rg_pattern: located_rg_files.append(locate_file(rg_dir, pat)) located_radar_files = [] for ts in rad_timestrings: located_radar_files.append(locate_file(rad_dir, ts)) .. GENERATED FROM PYTHON SOURCE LINES 103-105 Step 5: Read data from raingauge files into a Rain object. Coordinates are in WGS84, following that in the files. .. GENERATED FROM PYTHON SOURCE LINES 105-114 .. code-block:: default rg_object_wgs84 = Rain( located_rg_files, 'WGS84', duration=pd.Timedelta(minutes=5), NAN=[3276.7, 32767] ) .. GENERATED FROM PYTHON SOURCE LINES 115-118 Step 6: Read radar files into xarray.DataArrays. The data in the files are already in Cartesian Coordinates, in the Centered Azimuthal Projection. .. GENERATED FROM PYTHON SOURCE LINES 118-128 .. code-block:: default reflec_list = [] for file in located_radar_files: reflec = read_irisref( file ) reflec_list.append(reflec) .. rst-class:: sphx-glr-script-out .. code-block:: none /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/wradlib/io/iris.py:4095: FutureWarning: IRIS Cartesian Product is currently returned with ``origin='upper'``. From wradlib version 2.0 the images will be returned with ``origin='lower'``. To silence this warning set kwarg ``origin='upper'`` or ``origin='lower'``. filename, loaddata=loaddata, rawdata=rawdata, debug=debug, **kwargs /tmp/build/docs/swirlspy/swirlspy/rad/_utils.py:48: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() .. GENERATED FROM PYTHON SOURCE LINES 129-134 Data Reprojection --------------------------------------------------------------------- This section demonstrates the reprojection of extracted raingauge and radar data to a user-defined grid. .. GENERATED FROM PYTHON SOURCE LINES 136-137 Step 1: Define the target grid as a pyresample AreaDefinition. .. GENERATED FROM PYTHON SOURCE LINES 137-159 .. code-block:: default # Defining target grid area_id = "hk1980_250km" description = ("A 250 m resolution rectangular grid " "centred at HKO and extending to 250 km " "in each direction in HK1980 easting/northing coordinates") proj_id = 'hk1980' projection = ('+proj=tmerc +lat_0=22.31213333333334 ' '+lon_0=114.1785555555556 +k=1 +x_0=836694.05 ' '+y_0=819069.8 +ellps=intl +towgs84=-162.619,-276.959,' '-161.764,0.067753,-2.24365,-1.15883,-1.09425 +units=m ' '+no_defs') x_size = 1000 y_size = 1000 area_extent = (792000, 796000, 880000, 862000) area_def = utils.get_area_def( area_id, description, proj_id, projection, x_size, y_size, area_extent ) .. rst-class:: sphx-glr-script-out .. code-block:: none /tmp/build/docs/swirlspy/swirlspy/examples/qpe_hk.py:156: UserWarning: 'get_area_def' has moved, import it with 'from pyresample import get_area_def' area_id, description, proj_id, projection, x_size, y_size, area_extent /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() .. GENERATED FROM PYTHON SOURCE LINES 160-162 Step 2: Convert coordinates of raingauge object to desired projection. In this example, the desired projection is HK1980. .. GENERATED FROM PYTHON SOURCE LINES 162-214 .. code-block:: default # Write a function to make conversion neater def convert_rain( rg_object_wgs84, area_def ): """ Convert coordinates in object from WGS84 to HK1980. Parameters ----------- rg_object_wgs84: Rain Rain object with WGS84 coordinates. area_def: pyresample.geometry.AreaDefinition AreaDefinition of the target grid. Returns ------- rg_object: Rain Rain object in HK1980. """ # Defining Proj objects geod = pyproj.Proj(init="epsg:4326") outproj = pyproj.Proj(area_def.proj_str) gauge_x = [] gauge_y = [] data = [] for tup in rg_object_wgs84.data: lat = tup[0] lon = tup[1] easting, northing = pyproj.transform( geod, outproj, lon, lat ) gauge_x.append(easting) gauge_y.append(northing) data.append((northing, easting, tup[2], tup[3])) rg_object = copy.copy(rg_object_wgs84) rg_object.gauge_x = gauge_x rg_object.gauge_y = gauge_y rg_object.data = data rg_object.proj = 'HK80' return rg_object .. GENERATED FROM PYTHON SOURCE LINES 215-219 .. code-block:: default # Calling the function to convert projection system of coordinates rg_object = convert_rain(rg_object_wgs84, area_def) .. rst-class:: sphx-glr-script-out .. code-block:: none /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:53: FutureWarning: '+init=:' syntax is deprecated. ':' is the preferred initialization method. When making the change, be mindful of axis order changes: https://pyproj4.github.io/pyproj/stable/gotchas.html#axis-order-changes-in-proj-6 return _prepare_from_string(" ".join(pjargs)) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:294: FutureWarning: '+init=:' syntax is deprecated. ':' is the preferred initialization method. When making the change, be mindful of axis order changes: https://pyproj4.github.io/pyproj/stable/gotchas.html#axis-order-changes-in-proj-6 projstring = _prepare_from_string(" ".join((projstring, projkwargs))) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() .. GENERATED FROM PYTHON SOURCE LINES 220-223 Step 3: Regrid radar reflectivity from Centered Azimuthal Projection to HK1980. .. GENERATED FROM PYTHON SOURCE LINES 223-235 .. code-block:: default reproj_reflec_list = [] for reflec in reflec_list: reproj_reflec = grid_resample( reflec, reflec.attrs['area_def'], area_def, coord_label=['easting', 'northing'] ) reproj_reflec_list.append(reproj_reflec) reflectivity = xarray.concat(reproj_reflec_list, dim='time') .. rst-class:: sphx-glr-script-out .. code-block:: none /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/image.py:62: FutureWarning: Usage of ImageContainer is deprecated, please use NumpyResamplerBilinear class instead FutureWarning) /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() .. GENERATED FROM PYTHON SOURCE LINES 236-244 Accumulating and interpolating rainfall ----------------------------------------------------------------- Interpolate rainfall recorded by raingauges into the user-defined grid and accumulate radar rainfall over an hour after making the necessary adjustments. .. GENERATED FROM PYTHON SOURCE LINES 246-248 Step 1: Interpolate Rain object to user-defined grid. .. GENERATED FROM PYTHON SOURCE LINES 248-255 .. code-block:: default interpolated_rg = rg_interpolate( rg_object, area_def, 'rbf', coord_label=['easting', 'northing'] ) .. GENERATED FROM PYTHON SOURCE LINES 256-259 Step 2: Convert to radar reflectivity to rainrates, interpolate rainrates to times of raingauges, and convert to rainfalls accumulated every 5 minutes. .. GENERATED FROM PYTHON SOURCE LINES 259-279 .. code-block:: default rainrates = dbz2rr(reflectivity, a=58.53, b=1.56) # Convert time coordinates of rainrates from start time # to endtime rainrates_time_endtime = rainrates.copy() rainrates_time_endtime.coords['time'] = \ [ pd.Timestamp(t) + pd.Timedelta(minutes=6) for t in rainrates.coords['time'].values ] rainrates_5min = temporal_interp( rainrates_time_endtime, rg_object.start_time + pd.Timedelta(minutes=5), rg_object.end_time, intvl=pd.Timedelta(minutes=5), interp_type='quadratic' ) rainfalls = rr2rf(rainrates_5min, scan_duration=5) .. GENERATED FROM PYTHON SOURCE LINES 280-281 Step 3: Accumulate rainfall over an hour. .. GENERATED FROM PYTHON SOURCE LINES 281-288 .. code-block:: default acc_rf = acc( rainfalls, rg_object.end_time, acc_period=pd.Timedelta(minutes=60) ) .. GENERATED FROM PYTHON SOURCE LINES 289-294 Compositing rainfall ----------------------------------------------------------------- Perform compositing on radar and raingauge derived rainfall to obtain a composite QPE. .. GENERATED FROM PYTHON SOURCE LINES 294-302 .. code-block:: default comprf = comp_qpe( area_def, rg_object=rg_object, rg_interp=interpolated_rg, rad_rf=acc_rf ) .. GENERATED FROM PYTHON SOURCE LINES 303-308 Plotting --------------------------------------------------------------- Plot composited radar and raingauge rainfall. .. GENERATED FROM PYTHON SOURCE LINES 308-415 .. code-block:: default # Plotting function for neatness. def plot_map( da, rg_object, acctime, area_def, based='raingauge and radar', savepath='', ): """ A custom function for plotting a map. Parameters -------------- da: xarray.DataArray Contains data to be plotted. rg_object: Rain Contains raingauge data. acctime: pd.Timestamp Contains the endtime of the accumulation period. area_def: pyresample.geometry.AreaDefinition AreaDefinition of the grid. based: str Type of data plotted in the map. savepath: str Path to save the image to. """ # Defining the colour scheme levels = [ 0, 0.5, 2, 5, 10, 20, 30, 40, 50, 70, 100, 150, 200, 300, 400, 500, 600, 700 ] cmap = matplotlib.colors.ListedColormap([ '#ffffff', '#9bf7f7', '#00ffff', '#00d5cc', '#00bd3d', '#2fd646', '#9de843', '#ffdd41', '#ffac33', '#ff621e', '#d23211', '#9d0063', '#e300ae', '#ff00ce', '#ff57da', '#ff8de6', '#ffe4fd' ]) norm = BoundaryNorm(levels, ncolors=cmap.N, clip=True) # Plotting axes plt.figure(figsize=(28, 21)) crs = area_def.to_cartopy_crs() ax = plt.axes(projection=crs) # Plot data quadmesh = da.plot( cmap=cmap, norm=norm, extend='neither', cbar_kwargs={'ticks': levels} ) # Adjusting size of colorbar cb = quadmesh.colorbar cb.ax.set_ylabel( da.attrs['long_name']+'['+da.attrs['units']+']', fontsize=28 ) cb.ax.tick_params(labelsize=24) # Setting labels ax.xaxis.set_visible(True) ax.yaxis.set_visible(True) for tick in ax.xaxis.get_major_ticks(): tick.label.set_fontsize(24) for tick in ax.yaxis.get_major_ticks(): tick.label.set_fontsize(24) ax.xaxis.label.set_size(28) ax.yaxis.label.set_size(28) # Coastlines shpfile = THIS_DIR + \ '/../tests/samples/gadm36_HKG_shp/gadm36_HKG_0_hk1980.shp' shp = shapereader.Reader(shpfile) for record, geometry in zip(shp.records(), shp.geometries()): ax.add_geometries([geometry], crs, facecolor='None', edgecolor='black') # Show raingauge show_raingauge( rg_object, ax, show_value=False, color='red', markersize=16, fontsize=16, adjust_labels=False ) # Show title plt.title( (f"Last Hour Rainfall at {acctime.strftime('%H:%MH %d-%b-%Y')}\n" f"(based on {based} data)"), fontsize=32 ) # from cartopy import feature # hires = feature.GSHHSFeature( # levels=[1], # scale='h', # edgecolor='red' # ) # ax.add_feature(hires) plt.savefig(savepath, dpi=300) .. GENERATED FROM PYTHON SOURCE LINES 416-441 .. code-block:: default # Plotting maps # Raingauge only plot_map( interpolated_rg, rg_object, acctime, area_def, based='raingauge', savepath=THIS_DIR+f'/../tests/outputs/raingauge_{acctime_str}.png' ) # Radar only plot_map( acc_rf, rg_object, acctime, area_def, based='radar', savepath=THIS_DIR+f'/../tests/outputs/radar_{acctime_str}.png' ) # Composite raingauge and radar plot_map( comprf, rg_object, acctime, area_def, savepath=THIS_DIR+f'/../tests/outputs/comp_{acctime_str}.png' ) .. rst-class:: sphx-glr-horizontal * .. image-sg:: /auto_examples/images/sphx_glr_qpe_hk_001.png :alt: Last Hour Rainfall at 15:00H 20-Apr-2019 (based on raingauge data) :srcset: /auto_examples/images/sphx_glr_qpe_hk_001.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/images/sphx_glr_qpe_hk_002.png :alt: Last Hour Rainfall at 15:00H 20-Apr-2019 (based on radar data) :srcset: /auto_examples/images/sphx_glr_qpe_hk_002.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/images/sphx_glr_qpe_hk_003.png :alt: Last Hour Rainfall at 15:00H 20-Apr-2019 (based on raingauge and radar data) :srcset: /auto_examples/images/sphx_glr_qpe_hk_003.png :class: sphx-glr-multi-img .. rst-class:: sphx-glr-script-out .. code-block:: none /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyresample/geometry.py:1617: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_params = self.crs.to_proj4() /opt/conda/envs/swirlspy/lib/python3.6/site-packages/pyproj/crs/crs.py:543: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems proj_string = self.to_proj4() .. rst-class:: sphx-glr-timing **Total running time of the script:** ( 1 minutes 10.751 seconds) .. _sphx_glr_download_auto_examples_qpe_hk.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: qpe_hk.py ` .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: qpe_hk.ipynb ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_