
PyCWR:天气雷达数据处理工具库¶
- Release
V0.3.4
- Date
Mar 01, 2023
PyCWR (Python based China Weather Radar Toolkit)主要面向中国天气雷达基数据的科研和应用, 广泛支持中国业务天气雷达基数据(WSR98D、CINRAD/SA/SB/CB、CINRAD/CC/CCJ和CINRAD/SC/CD等), 方便用于科研和业务上数据处理和分析,PyCWR现阶段主要由NJIAS团队进行维护和更新。
使用文档¶
开始使用
PyCWR简介
安装方法
图形化界面显示
常见问题
PyCWR简介¶
中国天气雷达基数据的格式由于历史遗留问题多种多样,文件复杂,雷达数据的处理和应用往往需要重要的先验知识和技能,各种算法难以统一实施。 基于Python的中国天气雷达(PyCWR)库旨在解决这些问题并使雷达数据变得易于使用。
Important
它支持几乎所有中国业务天气雷达基数据格式,包含WSR98D、CINRAD/SA/SB/CB、CINRAD/CC/CCJ、CINRAD/SC/CD和相控阵等雷达基数据格式;
定义用于表示雷达数据的核心Python对象,而该对象是使用xarray的Dataset构建的,Dataset是类似dict的具有对齐尺寸的带标签数组的容器;
提供了经过良好测试的算法(数据质量控制,水凝物分类,定量降水估计等),并使用Cython加速核心算法;
建立开放的平台,供研究人员从天气雷达观测中开发复杂的分析和应用程序,它还提供了图形用户界面工具,以方便进行水文学和气象学分析。
天气雷达PPI扫描显示:
天气雷达RHI扫描显示:
天气雷达CAPPI插值显示:
天气雷达组合反射率因子显示:
天气雷达VCS垂直剖面图:
水凝物分类算法:
更多功能更新中…
安装方法¶
Anaconda¶
目前 PyCWR 仅支持Python3,推荐使用Python3.8及以上版本的Anaconda。
由于国内墙下载速度的限制,推荐在 清华镜像站 下载并安装Anaconda,提高下载速度。
Anaconda安装完毕后请按照 conda源修改教程 修改conda源为清华源,按照以下教程修改pip源为豆瓣源以提高库的安装速度,节省您的时间。
Linux/Mac平台
Linux/Mac用户将它命名为pip.conf
[global] timeout = 60 index-url = http://pypi.douban.com/simple trusted-host = pypi.douban.com
然后将该文件放在
$HOME/.pip/pip.conf
位置Windows平台
Windows用户将它命名为pip.ini
[global] timeout = 60 index-url = http://pypi.douban.com/simple trusted-host = pypi.douban.com
然后将该文件放在
%HOME%\pip\pip.ini
位置
直接安装 PyCWR 在base环境¶
由于Cartopy采用pip安装容易出错,但 PyCWR 画图部分需要依赖 Cartopy, 因此需要先采用conda来安装 Cartopy,再使用pip安装 PyCWR。
首先使用conda安装Cartopy
conda install cartopy -c conda-forge --yes
然后使用pip安装PyCWR
pip install pycwr==0.3.6
隔离环境安装 PyCWR (推荐)¶
为防止 PyCWR 环境与base环境冲突,推荐使用隔离环境的方式安装。
首先使用conda隔离环境并安装Cartopy
conda create -n cwr cartopy -c conda-forge --yes
再切换到隔离的cwr环境
conda activate cwr
最后在cwr环境中使用pip安装PyCWR
pip install pycwr==0.3.6
图形化界面显示¶
用Python执行此脚本可以启动图形化界面程序:
1import warnings
2warnings.filterwarnings('ignore')
3import os, sys
4sys.path.append(os.path.dirname(os.path.dirname(os.path.abspath(__file__))))
5from pycwr.GraphicalInterface.RadarInterface import MainWindow
6from PyQt5 import QtWidgets
7app = QtWidgets.QApplication(sys.argv)
8ui = MainWindow()
9ui.show()
10sys.exit(app.exec_())
图形化界面打开后如下图:

在右侧可以选择要显示的 产品 和 仰角 ,在 设置->显示设置->叠加地图 可以选择显示的图像是否叠加地图。
常见问题¶
用户指南
数据读取
数据结构
绘图
导出数据
选取数据
水凝物分类
衰减订正
雷达数据插值
组合反射率产品
CAPPI产品
风场反演
数据读取¶
双偏振雷达普及之前,我国新一代天气雷达数据主要有CINRAD/SA/SB/CB、CINRAD/CC/CCJ、CINRAD/SC/CD等格式, 随着双偏振雷达数据的普及,我国逐渐将数据统一为WSR98D格式。
文档中采用的数据可通过百度网盘下载,链接: https://pan.baidu.com/s/1MhW9jE7r3LYHO7enUOfRsw 提取码: 8ba9
自动读取¶
可以使用read_auto函数实现雷达数据类型的识别,并转化为pycwr内置的雷达数据存储格式PRD(Polarimetry Radar Data)类型:
1from pycwr.io import read_auto
2PRD = read_auto(r"./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2")
3print(PRD.scan_info)
4print(PRD.fields)
read_auto函数在 pycwr.io
下,可以传入站点的经纬度高度参数,如不传入该参数,则会根据站号或者站点信息识别站点
的经纬度数据。
read_auto函数 read_auto(filename, station_lon=None, station_lat=None, station_alt=None)
共有四个参数:
filename : 雷达基数据路径,基数据可以是bz2和gz格式的压缩文件(不需要解压)。
station_lon : 雷达站点的经度信息,范围是-180 ~ 180。(units: degree east)
station_lat : 雷达站点的纬度信息,范围是-90 ~ 900。(units: degree north)
station_alt : 雷达站点的高度信息。(units: meters)
指定类型读取¶
如果已经知道雷达数据格式的类型,可以直接指定格式进行读取。
数据类型是WSR98D:
1from pycwr.io import WSR98DFile 2filename = "./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2" 3PRD = WSR98DFile.WSR98D2NRadar(WSR98DFile.WSR98DBaseData(filename, station_lon=None,\ 4 station_lat=None, station_alt=None)).ToPRD() 5print(PRD.scan_info) 6print(PRD.fields)
数据类型是SA/SB:
1from pycwr.io import SABFile 2filename = "./data/Z9396_BASE_SB_20180724_055400.bin.bz2" 3PRD = SABFile.SAB2NRadar(SABFile.SABBaseData(filename, station_lon=None,\ 4 station_lat=None, station_alt=None)).ToPRD() 5print(PRD.scan_info) 6print(PRD.fields)
数据类型是CC/CCJ:
1from pycwr.io import CCFile 2filename = "./data/2016070817.48V.gz" 3PRD = CCFile.CC2NRadar(CCFile.CCBaseData(filename, station_lon=None,\ 4 station_lat=None, station_alt=None)).ToPRD() 5print(PRD.scan_info) 6print(PRD.fields)
数据类型是SC/CD:
1from pycwr.io import SCFile 2filename = "./data/Z_RADR_I_Z9240_20190703101340_O_DOR_SC_CAP.bin.bz2" 3PRD = SCFile.SC2NRadar(SCFile.SCBaseData(filename, station_lon=None,\ 4 station_lat=None, station_alt=None)).ToPRD() 5print(PRD.scan_info) 6print(PRD.fields)
数据类型是Phase Array:
1from pycwr.io import read_PA 2filename = "./data/Z_RADR_I_ZGZ01_20200820220246_O_DOR_DXK_CAR.bin.bz2" 3PRD = read_PA(filename) 4print(PRD.scan_info) 5print(PRD.fields)
配合 Py-ART 库对雷达数据进行处理¶
1from pycwr.io import read_auto
2PRD = read_auto(r"./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2")
3print(PRD.scan_info)
4print(PRD.fields)
5PyartRadar = PRD.ToPyartRadar()
数据结构¶
PyCWR的主要模块及作用:

为了增加雷达数据的易用性,读取后的雷达数据会转换为偏振雷达数据(PRD, Polarimetric Radar Data)类,该类由两部分构成, 一部分是扫描信息(scan_info),它主要包含了雷达站点的经纬度、高度、扫描类型、雷达频率、不模糊速度等信息, 另一部分是观测数据(fields),它是一个列表容器,存储不同仰角的雷达扫描数据,列表的每个元素都是xarray的DataArray类, DataArray类提供了一个围绕numpy的ndarray类的包装器,该包装器使用标注的维度和坐标来支持元数据操作, 极大的提高了数据的可读性和可理解性。其结构如下:

绘图¶
PPI绘图不叠加地图:
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import Graph
4
5filename = r"./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2"
6PRD = read_auto(filename)
7fig, ax = plt.subplots()
8graph = Graph(PRD)
9graph.plot_ppi(ax, 0, "dBZ", cmap="CN_ref") ## 0代表第一层, dBZ代表反射率产品
10graph.add_rings(ax, [0, 50, 100, 150, 200, 250, 300])
11ax.set_title("PPI Plot", fontsize=16)
12ax.set_xlabel("Distance From Radar In East (km)", fontsize=14)
13ax.set_ylabel("Distance From Radar In North (km)", fontsize=14)
14plt.show()

PPI绘图叠加地图:
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import GraphMap
4import cartopy.crs as ccrs
5filename = r"./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2"
6PRD = read_auto(filename)
7
8ax = plt.axes(projection=ccrs.PlateCarree())
9graph = GraphMap(PRD, ccrs.PlateCarree())
10graph.plot_ppi_map(ax, 0, "dBZ", cmap="CN_ref") ## 0代表第一层, dBZ代表反射率产品,cmap
11ax.set_title("PPI Plot with Map", fontsize=16)
12plt.tight_layout()
13plt.show()

雷达RHI绘图:
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import Graph
4
5filename = r"./data/NUIST.20170323.142921.AR2"
6PRD = read_auto(filename)
7
8fig, ax = plt.subplots()
9graph = Graph(PRD)
10graph.plot_rhi(ax, 0, field_name="dBZ", cmap="CN_ref", clabel="Radar Reflectivity")
11ax.set_ylim([0, 10]) #设置rhi的高度范围 (units:km)
12ax.set_xlabel("distance from radar (km)", fontsize=14)
13ax.set_ylabel("Height (km)", fontsize=14)
14plt.tight_layout()
15plt.show()

天气雷达剖面图:
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import Graph
4
5filename = r"./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2"
6PRD = read_auto(filename)
7
8fig, ax = plt.subplots()
9graph = Graph(PRD)
10graph.plot_vcs(ax, (0,0), (150, 0), "dBZ", cmap="copy_pyart_NWSRef") #起点,终点 (units: km)
11ax.set_ylim([0, 15])
12ax.set_xlim([0, 80])
13ax.set_ylabel("Height (km)", fontsize=14)
14ax.set_xlabel("Distance From Section Start (Uints:km)", fontsize=14)
15ax.set_title("VCS Plot", fontsize=16)
16plt.tight_layout()
17plt.show()

导出数据¶
目前可以通过pyart将数据导出为cfradial/UF格式
1from pycwr.io import read_auto
2import pyart
3import warnings
4warnings.filterwarnings("ignore")
5
6file = "./data/Z_RADR_I_Z9898_20190828181529_O_DOR_SAD_CAP_FMT.bin.bz2"
7PRD = read_auto(file)
8pyart_radar = PRD.ToPyartRadar()
9pyart.io.write_cfradial("./cfradial.nc", pyart_radar)
选取数据¶
fields属性存放了雷达体扫数据:
>>> from pycwr.io import read_auto
>>> PRD = read_auto("./data/NUIST.20150627.002438.AR2.bz2")
通过[层数]方法索引对应产品, 层数从0开始:
>>> print(PRD.fields[0])
<xarray.Dataset>
Dimensions: (range: 1933, time: 668)
Coordinates:
azimuth (time) float64 9.149 9.69 10.22 10.77 ... 8.094 8.611 9.168 9.693
elevation (time) float64 0.5933 0.5933 0.5933 ... 0.5054 0.5054 0.5054
x (time, range) float64 11.92 23.85 35.77 ... 2.439e+04 2.44e+04
y (time, range) float64 74.04 148.1 222.1 ... 1.428e+05 1.429e+05
z (time, range) float64 174.8 175.6 176.3 ... 2.688e+03 2.689e+03
lat (time, range) float64 32.21 32.21 32.21 ... 33.49 33.49 33.49
lon (time, range) float64 118.7 118.7 118.7 ... 119.0 119.0 119.0
* range (range) float64 75.0 150.0 225.0 ... 1.448e+05 1.449e+05 1.45e+05
* time (time) datetime64[ns] 2015-06-27T00:24:38.640231 ... 2015-06-2...
Data variables:
dBZ (time, range) float64 -13.6 17.34 24.1 ... 18.62 18.24 19.32
dBT (time, range) float64 -13.6 17.34 24.1 ... 18.62 18.24 19.32
V (time, range) float64 0.0 -0.08 -1.12 -1.55 ... -1.76 -3.22 -3.41
W (time, range) float64 0.01 0.01 0.01 1.12 ... 1.63 1.43 1.5 1.68
ZDR (time, range) float64 -0.81 13.73 0.41 -0.39 ... 1.1 0.67 0.43
KDP (time, range) float64 7.03 6.13 5.44 4.78 ... 0.5 0.47 0.48 0.51
PhiDP (time, range) float64 100.0 33.3 1.02 2.86 ... 103.2 104.5 102.6
CC (time, range) float64 1.0 0.57 0.938 0.994 ... 0.923 0.938 0.96
SQI (time, range) float64 1.0 1.0 1.0 0.966 ... 0.913 0.905 0.898
SNRH (time, range) float64 67.77 81.21 81.95 ... 14.46 14.08 15.16
通过[层数][关键词]方法索引对应产品,关键词详见表1:
>>> print(PRD.fields[0]["dBZ"])
<xarray.DataArray 'dBZ' (time: 668, range: 1933)>
array([[-13.60000038, 17.34000015, 24.10000038, ..., 16.87999916,
16.76000023, 17.52000046],
[-13.60000038, 16.65999985, 24.09000015, ..., 18.37000084,
16.42000008, 17.86000061],
[-13.60000038, 18.03000069, 23.57999992, ..., 16.98999977,
20.04000092, 20.85000038],
...,
[-13.60000038, 12.61999989, 24. , ..., 15.63000011,
14.63000011, 13.81000042],
[-13.60000038, 15.75 , 25.35000038, ..., 12.85999966,
17.29999924, 14.15999985],
[-13.60000038, 14.97999954, 24.55999947, ..., 18.62000084,
18.23999977, 19.31999969]])
Coordinates:
azimuth (time) float64 9.149 9.69 10.22 10.77 ... 8.094 8.611 9.168 9.693
elevation (time) float64 0.5933 0.5933 0.5933 ... 0.5054 0.5054 0.5054
x (time, range) float64 11.92 23.85 35.77 ... 2.439e+04 2.44e+04
y (time, range) float64 74.04 148.1 222.1 ... 1.428e+05 1.429e+05
z (time, range) float64 174.8 175.6 176.3 ... 2.688e+03 2.689e+03
lat (time, range) float64 32.21 32.21 32.21 ... 33.49 33.49 33.49
lon (time, range) float64 118.7 118.7 118.7 ... 119.0 119.0 119.0
* range (range) float64 75.0 150.0 225.0 ... 1.448e+05 1.449e+05 1.45e+05
* time (time) datetime64[ns] 2015-06-27T00:24:38.640231 ... 2015-06-2...
Attributes:
units: dBZ
standard_name: equivalent_reflectivity_factor
long_name: Reflectivity
valid_max: 80.0
valid_min: -30.0
coordinates: elevation azimuth range
scan_info存放了雷达位置信息以及扫描方式:
>>> print(PRD.scan_info)
<xarray.Dataset>
Dimensions: (sweep: 15)
Coordinates:
* sweep (sweep) int64 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14
Data variables:
latitude float64 32.21
longitude float64 118.7
altitude int64 87
scan_type <U3 'ppi'
frequency float64 5.592
start_time datetime64[ns] 2015-06-27T00:24:38.640231
end_time datetime64[ns] 2015-06-27T00:31:47.207645
nyquist_velocity (sweep) float32 13.4 13.4 13.4 13.4 ... 26.81 26.81 26.81
unambiguous_range (sweep) int32 145000 145000 145000 ... 70000 70000 10000
rays_per_sweep (sweep) int64 668 668 668 668 668 ... 669 668 668 668 668
fixed_angle (sweep) float32 0.5 1.5 2.4 3.4 ... 14.0 16.7 19.5 90.0
beam_width (sweep) float64 0.5389 0.5389 0.5389 ... 0.5389 0.5389
表1,不同变量对应索引的关键词:
key |
variable |
---|---|
dBT |
total_power |
dBZ |
reflectivity |
V |
velocity |
W |
spectrum_width |
SQI |
normalized_coherent_power |
CPA |
clutter_phase_alignment |
ZDR |
differential_reflectivity |
LDR |
linear_depolarization_ratio |
CC |
cross_correlation_ratio |
PhiDP |
differential_phase |
KDP |
specific_differential_phase |
CP |
clutter_probability |
Flag |
flag_of_rpv_data |
HCL |
hydro_class |
CF |
clutter_flag |
Zc |
corrected_reflectivity |
Vc |
corrected_velocity |
Wc |
spectrum_width_corrected |
SNRH |
horizontal_signal_noise_ratio |
SNRV |
vertical_signal_noise_ratio |
水凝物分类¶
水凝物分类示例程序:
1# -*- coding: utf-8 -*-
2from pycwr.io import read_auto
3from pycwr.retrieve.HID import fhc_HCL
4import matplotlib.pyplot as plt
5import numpy as np
6from pycwr.draw.RadarPlot import plot_xy, add_rings
7import pandas as pd
8
9file = r"./data/NUIST.20150627.002438.AR2.bz2"
10file_t = r"./data/20150627.csv"
11temp = pd.read_csv(file_t, index_col=0, header=None, names=['temp'])
12
13NRadar = read_auto(file)
14num = 3
15dBZ = np.where(NRadar.fields[num].CC>0.9, NRadar.fields[num].dBZ, np.nan)
16KDP = np.where(NRadar.fields[num].CC>0.9, NRadar.fields[num].KDP, np.nan)
17ZDR = np.where(NRadar.fields[num].CC>0.9, NRadar.fields[num].ZDR, np.nan)
18CC = np.where(NRadar.fields[num].CC>0.9, NRadar.fields[num].CC, np.nan)
19temp_2d = np.interp(NRadar.fields[num].z.values/1000., temp.index, temp['temp'])
20dBZ[:,0] = np.nan
21ticks = np.arange(1, 11, 1)
22ticklabels = ['Drizzle', 'Rain', 'Ice Crystals', 'Aggregates', 'Wet Snow', 'Vertical Ice',
23 'LD Graupel', 'HD Graupel', 'Hail', 'Big Drops']
24
25hcl = fhc_HCL(dBZ=dBZ, KDP=KDP, ZDR=ZDR, CC = CC, T=temp_2d)
26fig, ax = plt.subplots()
27plot_xy(ax, NRadar.fields[num].x, NRadar.fields[num].y, hcl,
28 cmap="CN_hcl", bounds=np.arange(0.5,10.6,1),
29 cbar_ticks=ticks, cbar_ticklabels=ticklabels)
30add_rings(ax=ax, rings=[0, 50, 100, 150])
31ax.set_xlim([-150, 150])
32ax.set_ylim([-150, 150])
33ax.set_xlabel("Distance From Radar In East (km)", fontsize=12)
34ax.set_ylabel("Distance From Radar In North (km)", fontsize=12)
35ax.set_title("Hydrometeor classification, El : 3.4", fontsize=14)
36plt.savefig(r"./201506270024_HC.png", dpi=600)
37plt.show()
水凝物分类效果:

衰减订正¶
雷达数据插值¶
组合反射率产品¶
以雷达为中心,生成笛卡尔坐标的组合反射率网格产品
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import plot_xy
4import numpy as np
5
6filename = r"../../data/NUIST.20150627.002438.AR2.bz2"
7PRD = read_auto(filename)
8
9x1d = np.arange(-150000, 150001, 1000) ##x方向1km等间距, -150km~150km范围
10y1d = np.arange(-150000, 150001, 1000) ##y方向1km等间距, -150km~150km范围
11PRD.add_product_CR_xy(XRange=x1d, YRange=y1d)
12print(PRD.product)
13grid_x, grid_y = np.meshgrid(x1d, y1d, indexing="ij")
14fig, ax = plt.subplots()
15plot_xy(ax, grid_x, grid_y, PRD.product.CR) ##画图显示
16ax.set_xlabel("Distance From Radar In East (km)", fontsize=14)
17ax.set_ylabel("Distance From Radar In North (km)", fontsize=14)
18plt.tight_layout()
19plt.show()

利用经纬度坐标信息生成组合反射率网格产品
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import plot_lonlat_map
4import cartopy.crs as ccrs
5import numpy as np
6
7filename = r"./data/NUIST.20150627.002438.AR2.bz2"
8PRD = read_auto(filename)
9
10lon1d = np.arange(117, 120.001, 0.01) ##lon方向0.01等间距,117-120范围
11lat1d = np.arange(31, 34.001, 0.01) ##lat方向0.01等间距, 31-34度范围
12PRD.add_product_CR_lonlat(XLon=lon1d, YLat=lat1d)
13# XLon:np.ndarray, 1d, units:degrees
14# YLat:np.ndarray, 1d, units:degrees
15# level_height:常量,要计算的高度 units:meters
16grid_lon, grid_lat = np.meshgrid(lon1d, lat1d, indexing="ij")
17ax = plt.axes(projection=ccrs.PlateCarree())
18plot_lonlat_map(ax, grid_lon, grid_lat, PRD.product.CR_geo, transform=ccrs.PlateCarree())
19ax.set_extent([117, 120, 31, 34], crs = ccrs.PlateCarree()) #设置显示范围
20plt.tight_layout()
21plt.show()

CAPPI产品¶
以雷达为中心,生成笛卡尔坐标的CAPPI网格产品
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import plot_xy
4import numpy as np
5
6filename = r"./data/NUIST.20150627.002438.AR2.bz2"
7PRD = read_auto(filename)
8
9x1d = np.arange(-150000, 150001, 1000) ##x方向1km等间距, -150km~150km范围
10y1d = np.arange(-150000, 150001, 1000) ##y方向1km等间距, -150km~150km范围
11PRD.add_product_CAPPI_xy(XRange=x1d, YRange=y1d, level_height=3000) ##level height units:meters
12print(PRD.product)
13grid_x, grid_y = np.meshgrid(x1d, y1d, indexing="ij")
14fig, ax = plt.subplots()
15plot_xy(ax, grid_x, grid_y, PRD.product.CAPPI_3000) ##画图显示
16ax.set_xlabel("Distance From Radar In East (km)", fontsize=14)
17ax.set_ylabel("Distance From Radar In North (km)", fontsize=14)
18plt.tight_layout()
19plt.show()

利用经纬度坐标信息生成CAPPI网格产品
1from pycwr.io import read_auto
2import matplotlib.pyplot as plt
3from pycwr.draw.RadarPlot import plot_lonlat_map
4import cartopy.crs as ccrs
5import numpy as np
6
7filename = r"./data/NUIST.20150627.002438.AR2.bz2"
8PRD = read_auto(filename)
9
10lon1d = np.arange(117, 120.001, 0.01) ##lon方向0.01等间距,117-120范围
11lat1d = np.arange(31, 34.001, 0.01) ##lat方向0.01等间距, 31-34度范围
12PRD.add_product_CAPPI_lonlat(XLon=lon1d, YLat=lat1d, level_height=3000) ##插值1500m高度的
13# XLon:np.ndarray, 1d, units:degrees
14# YLat:np.ndarray, 1d, units:degrees
15# level_height:常量,要计算的高度 units:meters
16grid_lon, grid_lat = np.meshgrid(lon1d, lat1d, indexing="ij")
17ax = plt.axes(projection=ccrs.PlateCarree())
18plot_lonlat_map(ax, grid_lon, grid_lat, PRD.product.CAPPI_geo_3000, transform=ccrs.PlateCarree())
19ax.set_extent([117, 120, 31, 34], crs = ccrs.PlateCarree()) #设置范围
20plt.tight_layout()
21plt.show()

风场反演¶
开发者信息¶
- 作者
- 地址
南京市建邺区雨顺路6号江苏省气象预警中心
- 邮箱
- 贡献
特别感谢李扬、张昕、贾鹏程、吕星超和樊丝慧等人对项目的支持。