Issue with wrl.georef.georeference_dataset function, and CAPPI from netCDF !!

131 views
Skip to first unread message

Syed Hamid Ali

unread,
Jul 13, 2021, 6:34:15 AM7/13/21
to wradlib-users
Hello Kai,
I was trying to create CAPPI from netCDF, and after applying wrl.georef.georeference_dataset function; the newly created dimensions (x,y,z) are giving the same shape, which is ((643, 2500), (643, 2500), (643, 2500)), but I think z should differ from x and y as we are considering z up to 15 km only for the vertical profile. 
and the plot which I am getting is like (attached)
download (12).png

Kai Mühlbauer

unread,
Jul 13, 2021, 7:56:34 AM7/13/21
to wradli...@googlegroups.com
Hi Syed,

first, the x,y and z arrays have azimuth * range dimensions. So every
radar bin gets a distinct x, y and z coordinate. That is they have the
same dimensionality.

Try to plot these arrays:

swp.x.plot(x="x", y="y")
swp.y.plot(x="x", y="y")
swp.z.plot(x="x", y="y")

and you will immediately see what is going on. You can specify other
projections (AEQD is default) in the call to
wrl.georef.georeference_dataset. If you have any questions on this
please let me know.

Second, you are only evaluating the lowest sweep (vol[0]). For CAPPI you
would need to extract all sweeps and feed their coordinates/data to the
CAPPI function(s). Then your result looks something like the attached image.

Cheers,
Kai



Am 13.07.21 um 12:34 schrieb Syed Hamid Ali:
> <https://drive.google.com/file/d/1_nqDaUqzalYbSqAky-T86oA714aLLaWX/view?usp=sharing>
>
> --
> You received this message because you are subscribed to the Google
> Groups "wradlib-users" group.
> To unsubscribe from this group and stop receiving emails from it, send
> an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-user...@googlegroups.com>.
> To view this discussion on the web, visit
> https://groups.google.com/d/msgid/wradlib-users/a0fb6e7f-c90d-499d-82ac-31e0236c1518n%40googlegroups.com
> <https://groups.google.com/d/msgid/wradlib-users/a0fb6e7f-c90d-499d-82ac-31e0236c1518n%40googlegroups.com?utm_medium=email&utm_source=footer>.

MAXCAPPI.png

Syed Hamid Ali

unread,
Jul 13, 2021, 8:34:09 AM7/13/21
to wradli...@googlegroups.com
Dear Kai, would you please attach the corresponding code for the above picture of cappi, which you created. 
I am still struggling o extract all sweeps and feed their coordinates/data to the
CAPPI function(s).

To unsubscribe from this group and stop receiving emails from it, send an email to wradlib-user...@googlegroups.com.
To view this discussion on the web, visit https://groups.google.com/d/msgid/wradlib-users/005ea2d3-bf74-fa55-b193-d62c8315a35e%40uni-bonn.de.


--

Syed Hamid Ali

Master of Technology Student

Department of Atmospheric and Space Sciences,

Savitribai Phule Pune University &

Indian Institute of Tropical Meteorology,

Pune, India

Mobile : +91 9622222937

Syed Hamid Ali  

Kai Mühlbauer

unread,
Jul 13, 2021, 8:48:32 AM7/13/21
to wradli...@googlegroups.com
Sure!

# Open volume
vol = wrl.io.open_cfradial1_dataset("MDV-20180602-101927-PPIVol.nc")
# Iterate over sweeps and process
swp_list = []
for v in vol:
# georeference and stack dimensions
swp =
v.pipe(wrl.georef.georeference_dataset).stack(points=["azimuth", "range"])
swp_list.append(swp)
# concat sweeps to volume
vol0 = xr.concat(swp_list, dim="points")
# Create XYZ Coordinate DataArray
xyz = xr.concat([vol0.x, vol0.y, vol0.z], dim="xyz").transpose()
# Create Target 3D Grid
trgx = np.linspace(xyz[:, 0].min(),xyz[:, 0].max(), 100)
trgy = np.linspace(xyz[:, 1].min(),xyz[:, 1].max(), 100)
trgz = np.linspace(0, 20000, 200)
yy, hh, xx = np.meshgrid(trgy,trgz,trgx)
trgxyz = np.stack([xx.flatten(), yy.flatten(), hh.flatten()]).T
# Create Gridder/Interpolator
trgshape=xx.shape
gridder = wrl.vpr.CAPPI(polcoords=xyz,
gridcoords=trgxyz,
gridshape=trgshape,
maxrange=200000,
minelev=0,
maxelev=90,
ipclass=wrl.ipol.Nearest)


# Interpolate Data into 3D Grid
crtd_ref = vol0.UH.where(((vol0.RHOHV < 1.0) & ((vol0.PHIDP > 150) |
(vol0.PHIDP < -100))))
vol_zh = np.ma.masked_invalid(gridder(crtd_ref.values).reshape(trgshape))

trgx = trgxyz[:, 0].reshape(trgshape)[0, 0, :]
trgy = trgxyz[:, 1].reshape(trgshape)[0, :, 0]
trgz = trgxyz[:, 2].reshape(trgshape)[:, 0, 0]
yy, hh, xx = np.meshgrid(trgy,trgz,trgx)

# Plot MAXCAPPI
wrl.vis.plot_max_plan_and_vert(trgx, trgy, trgz, vol_zh, unit="dBZH",
cmap="jet", levels=range(-10, 60))

For real CAPPI you would just need to plot the single wanted height layer.

Cheers,
Kai

Am 13.07.21 um 14:33 schrieb Syed Hamid Ali:
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>
> > <mailto:wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>>.
> <https://groups.google.com/d/msgid/wradlib-users/a0fb6e7f-c90d-499d-82ac-31e0236c1518n%40googlegroups.com?utm_medium=email&utm_source=footer
> <https://groups.google.com/d/msgid/wradlib-users/a0fb6e7f-c90d-499d-82ac-31e0236c1518n%40googlegroups.com?utm_medium=email&utm_source=footer>>.
>
> --
> You received this message because you are subscribed to the Google
> Groups "wradlib-users" group.
> To unsubscribe from this group and stop receiving emails from it,
> send an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>.
> To view this discussion on the web, visit
> https://groups.google.com/d/msgid/wradlib-users/005ea2d3-bf74-fa55-b193-d62c8315a35e%40uni-bonn.de
> <https://groups.google.com/d/msgid/wradlib-users/005ea2d3-bf74-fa55-b193-d62c8315a35e%40uni-bonn.de>.
>
>
>
> --
>
> Syed Hamid Ali
>
> Master of Technology Student
>
> Department of Atmospheric and Space Sciences,
> <http://www.unipune.ac.in/dept/science/atmospheric_and_space_sciences/default.htm>
>
> Savitribai Phule Pune University &
>
> Indian Institute of Tropical Meteorology,
>
> Pune, India
>
> Mobile : +91 9622222937
>
> <https://github.com/syedhamidali>Syed Hamid Ali
> <https://www.linkedin.com/in/rizvihamid/><https://www.researchgate.net/profile/Syed-Ali-148><https://twitter.com/hamidrixvi>
>
> --
> You received this message because you are subscribed to the Google
> Groups "wradlib-users" group.
> To unsubscribe from this group and stop receiving emails from it, send
> an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-user...@googlegroups.com>.
> To view this discussion on the web, visit
> https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com
> <https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com?utm_medium=email&utm_source=footer>.

Syed Hamid Ali

unread,
Jul 13, 2021, 12:46:49 PM7/13/21
to wradli...@googlegroups.com
Thanks, Kai,
Saved a lot of time. 
I was wondering if I could add text in the figure like this one, I tried to do that by fig.text() command, but somehow it is not working with wrl.vis.plot_max_plan_and_vert()
This is what i first tried manually, without interpolating the panel data., that's why it is looking like that. 
image.png

To unsubscribe from this group and stop receiving emails from it, send an email to wradlib-user...@googlegroups.com.
To view this discussion on the web, visit https://groups.google.com/d/msgid/wradlib-users/14f382e9-aab9-c17e-eec8-4c620e5b4b33%40uni-bonn.de.


--

Syed Hamid Ali

Savitribai Phule Pune University &

Indian Institute of Tropical Meteorology,

Pune, India

Mobile : +91 9622222937

Syed Hamid Ali  

Kai Mühlbauer

unread,
Jul 14, 2021, 9:55:10 AM7/14/21
to wradli...@googlegroups.com
Am 13.07.21 um 18:46 schrieb Syed Hamid Ali:
> Thanks, Kai,
> Saved a lot of time.
> I was wondering if I could add text in the figure like this one, I tried
> to do that by fig.text() command, but somehow it is not working
> with wrl.vis.plot_max_plan_and_vert()

No unfortunately not. This function is there almost from the beginning
of wradlib 10 years ago, when interactive use of pylab was very common.
Today where the object-oriented approach is heavily advertised, we
should try to revisit all these plotting functions and refactor them.

You can copy the two needed functions from wradlib codebase and add the
text there.

If you like you can open an issue on the github issue tracker so we do
not forget.

Cheers,
Kai

> This is what i first tried manually, without interpolating the panel
> data., that's why it is looking like that.
> <mailto:kai.mue...@uni-bonn.de
> >     <mailto:wradlib-users%2Bunsu...@googlegroups.com
> <mailto:wradlib-users%252Buns...@googlegroups.com>>
> >      > <mailto:wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>
> >     <mailto:wradlib-users%2Bunsu...@googlegroups.com
> <mailto:wradlib-users%252Buns...@googlegroups.com>>>.
> >     <mailto:wradlib-users%2Bunsu...@googlegroups.com
> <mailto:wradlib-users%252Buns...@googlegroups.com>>.
> <https://www.linkedin.com/in/rizvihamid/>><https://www.researchgate.net/profile/Syed-Ali-148
> <https://www.researchgate.net/profile/Syed-Ali-148>><https://twitter.com/hamidrixvi
> <https://twitter.com/hamidrixvi>>
> >
> > --
> > You received this message because you are subscribed to the Google
> > Groups "wradlib-users" group.
> > To unsubscribe from this group and stop receiving emails from it,
> send
> > an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>
> > <mailto:wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>>.
> > To view this discussion on the web, visit
> >
> https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com
> <https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com>
>
> >
> <https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com?utm_medium=email&utm_source=footer
> <https://groups.google.com/d/msgid/wradlib-users/CADqG3xDv1MOK9NGDcgCJWN4MAsuRzKRFHG5eQ2c%3DGwmM8Mkr9g%40mail.gmail.com?utm_medium=email&utm_source=footer>>.
>
> --
> You received this message because you are subscribed to the Google
> Groups "wradlib-users" group.
> To unsubscribe from this group and stop receiving emails from it,
> send an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-users%2Bunsu...@googlegroups.com>.
> To view this discussion on the web, visit
> https://groups.google.com/d/msgid/wradlib-users/14f382e9-aab9-c17e-eec8-4c620e5b4b33%40uni-bonn.de
> <https://groups.google.com/d/msgid/wradlib-users/14f382e9-aab9-c17e-eec8-4c620e5b4b33%40uni-bonn.de>.
>
>
>
> --
>
> Syed Hamid Ali
>
> Master of Technology Student
>
> Department of Atmospheric and Space Sciences,
> <http://www.unipune.ac.in/dept/science/atmospheric_and_space_sciences/default.htm>
>
> Savitribai Phule Pune University &
>
> Indian Institute of Tropical Meteorology,
>
> Pune, India
>
> Mobile : +91 9622222937
>
> <https://github.com/syedhamidali>Syed Hamid Ali
> <https://www.linkedin.com/in/rizvihamid/><https://www.researchgate.net/profile/Syed-Ali-148><https://twitter.com/hamidrixvi>
>
> --
> You received this message because you are subscribed to the Google
> Groups "wradlib-users" group.
> To unsubscribe from this group and stop receiving emails from it, send
> an email to wradlib-user...@googlegroups.com
> <mailto:wradlib-user...@googlegroups.com>.
> To view this discussion on the web, visit
> https://groups.google.com/d/msgid/wradlib-users/CADqG3xBze118c8TaOsJvxn6AtEYchSnj3QNH%3D4Mhzr%2BPO6zfzQ%40mail.gmail.com
> <https://groups.google.com/d/msgid/wradlib-users/CADqG3xBze118c8TaOsJvxn6AtEYchSnj3QNH%3D4Mhzr%2BPO6zfzQ%40mail.gmail.com?utm_medium=email&utm_source=footer>.

Syed Hamid Ali

unread,
Jul 14, 2021, 12:17:14 PM7/14/21
to wradli...@googlegroups.com
Sure, I will modify the codebase, and yes I'll also open an issue on wradlib github page. 

To unsubscribe from this group and stop receiving emails from it, send an email to wradlib-user...@googlegroups.com.
To view this discussion on the web, visit https://groups.google.com/d/msgid/wradlib-users/b516b3f8-0fb0-81b7-bd64-cd286e3a022e%40uni-bonn.de.


--

Syed Hamid Ali

Savitribai Phule Pune University &

Indian Institute of Tropical Meteorology,

Pune, India

Mobile : +91 9622222937

Syed Hamid Ali  

Reply all
Reply to author
Forward
0 new messages