{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# 1.0 Basic Variable Extraction" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from __future__ import (absolute_import, division, print_function, unicode_literals)\n", "\n", "from wrf import getvar\n", "from netCDF4 import Dataset as nc\n", "#ncfile = nc(\"/Users/ladwig/Documents/wrf_files/wrfout_d01_2016-02-25_18_00_00\")\n", "ncfile = nc(\"/Users/ladwig/Documents/wrf_files/wrfout_d01_2016-10-07_00_00_00\")\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "p = getvar(ncfile, \"P\")\n", "print(p)\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.0.1 DataArray attributes: 'dims', 'coords', 'attrs'" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print(\"dims: \", p.dims)\n", "print(\"coords: \", p.coords) \n", "print(\"attrs: \", p.attrs)\n", "del p" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.0.2 Removing implicit 'squeeze' behavior to preserve single sized dimensions" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "p_nosqueeze = getvar(ncfile, \"P\", timeidx=0, squeeze=False)\n", "print (p_nosqueeze)\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.0.3 Single element metadata" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print (p_nosqueeze[0,0,100,200])\n", "del p_nosqueeze" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 1.0.4 Disabling/Enabling xarray" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf import disable_xarray, enable_xarray\n", "\n", "# Disable xarray completely\n", "disable_xarray()\n", "p_no_meta = getvar(ncfile, \"P\")\n", "print(type(p_no_meta))\n", "print (p_no_meta)\n", "del p_no_meta\n", "enable_xarray()\n", "\n", "# Disable on extraction\n", "p_no_meta = getvar(ncfile, \"P\", meta=False)\n", "print(\"\\n\")\n", "print(type(p_no_meta))\n", "print(p_no_meta)\n", "del p_no_meta\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 2.0 Sequences of Input Files " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 2.0.1 Combining via the 'cat' method" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf import ALL_TIMES\n", "import numpy as np\n", "wrflist = [ncfile, ncfile, ncfile]\n", "p_cat = getvar(wrflist, \"P\", timeidx=ALL_TIMES, method=\"cat\")\n", "\n", "print(p_cat)\n", "del p_cat" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 2.0.2 Combining via the 'join' method" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "p_join = getvar(wrflist, \"P\", timeidx=ALL_TIMES, method=\"join\")\n", "print(p_join)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note how the Time dimension was replaced with the file dimension, due to the 'squeezing' of the Time dimension.\n", "\n", "\n", "To maintain the Time dimension, set squeeze to False." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf import ALL_TIMES\n", "p_join = getvar(wrflist, \"P\", timeidx=ALL_TIMES, method=\"join\", squeeze=False)\n", "print(p_join)\n", "del p_join" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 2.0.3 Dictionary Sequences" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "wrf_dict = {\"label1\" : [ncfile, ncfile],\n", " \"label2\" : [ncfile, ncfile]}\n", "p_dict = getvar(wrf_dict, \"P\", timeidx=ALL_TIMES)\n", "print(p_dict)\n", "del p_dict" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 2.0.4 Generator Sequences" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "def gen_seq():\n", " wrfseq = [ncfile, ncfile, ncfile]\n", " for wrf in wrfseq:\n", " yield wrf\n", " \n", "p_gen = getvar(gen_seq(), \"P\", method=\"join\")\n", "print(p_gen)\n", "del p_gen" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 2.0.5 Custom Iterable Classes" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "class FileGen(object):\n", " def __init__(self, ncfile, count=3):\n", " self._total = count\n", " self._i = 0\n", " self.ncfile = [ncfile]*count\n", " \n", " def __iter__(self):\n", " return self\n", " \n", " def next(self):\n", " if self._i >= self._total:\n", " raise StopIteration\n", " else:\n", " val = self.ncfile[self._i]\n", " self._i += 1\n", " return val\n", " \n", " # Python 3\n", " def __next__(self):\n", " return self.next()\n", "\n", "obj_gen = FileGen(ncfile, 3)\n", "\n", "p_obj_gen = getvar(gen_seq(), \"P\", method=\"join\", squeeze=False)\n", "print(p_obj_gen)\n", "\n", "del p_obj_gen\n", " " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "# 3.0 WRF Variable Computational Routines" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "wrf_vars = [\"avo\", \"eth\", \"cape_2d\", \"cape_3d\", \"ctt\", \"dbz\", \"mdbz\", \n", " \"geopt\", \"helicity\", \"lat\", \"lon\", \"omg\", \"p\", \"pressure\", \n", " \"pvo\", \"pw\", \"rh2\", \"rh\", \"slp\", \"ter\", \"td2\", \"td\", \"tc\",\n", " \"theta\", \"tk\", \"tv\", \"twb\", \"updraft_helicity\", \"ua\", \"va\", \n", " \"wa\", \"uvmet10\", \"uvmet\", \"z\", \"ctt\"]\n", "wrf_vars = [\"slp\"]\n", "\n", "vard = {varname: getvar(ncfile, varname, method=\"cat\", squeeze=True) for varname in wrf_vars}\n", "for varname in wrf_vars:\n", " print(vard[varname])\n", " print (\"\\n\")\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "(Note all of the NaNs in the above routines which produce missing values (e.g. cape_2d). xarray always converts all masked_array missing values to NaN in order to work with pandas. To get back the original missing values in a numpy masked_array, you need to use the 'to_np' method from wrf.)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf import to_np\n", "masked_ndarray = to_np(vard[\"slp\"])\n", "print(type(masked_ndarray))\n", "del masked_ndarray" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "keys = [x for x in vard.keys()]\n", "for key in keys:\n", " del vard[key]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 3.1 Interpolation Routines" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 3.1.1 Horizontal Level Interpolation" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# 500 MB Heights\n", "from wrf import getvar, interplevel\n", "\n", "z = getvar(ncfile, \"z\")\n", "p = getvar(ncfile, \"pressure\")\n", "ht_500mb = interplevel(z, p, 500)\n", "\n", "print(ht_500mb)\n", "del ht_500mb, z, p" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 3.1.2 Vertical Cross Section Interpolation" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Pressure using pivot and angle\n", "from wrf import getvar, vertcross, CoordPair\n", "\n", "z = getvar(ncfile, \"z\")\n", "p = getvar(ncfile, \"pressure\")\n", "pivot_point = CoordPair((z.shape[-1]-1) // 2, (z.shape[-2] - 1) // 2) \n", "angle = 90.0\n", "\n", "p_vert = vertcross(p, z, pivot_point=pivot_point, angle=angle, latlon=True)\n", "print(p_vert)\n", "print (\"\\n\")\n", "del p_vert\n", "\n", "# Pressure using start_point and end_point\n", "start_point = CoordPair(0, (z.shape[-2]-1) // 2)\n", "end_point = CoordPair(-1, (z.shape[-2]-1) // 2)\n", "\n", "p_vert = vertcross(p, z, start_point=start_point, end_point=end_point, latlon=True)\n", "print(p_vert)\n", "del p_vert, p, z" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Pressure using pivot and angle\n", "from wrf import getvar, vertcross, CoordPair, xy_to_ll\n", "\n", "z = getvar(ncfile, \"z\")\n", "p = getvar(ncfile, \"pressure\")\n", "lats = getvar(ncfile, \"lat\")\n", "lons = getvar(ncfile, \"lon\")\n", "\n", "#print ((lats.shape[-2]-1) / 2)\n", "#print ((lats.shape[-1]-1) / 2)\n", "\n", "#print (to_np(lats[529, 899]))\n", "#print (to_np(lons[529, 899]))\n", "\n", "#print (to_np(lats[529, 0]))\n", "#print (to_np(lons[529, 0]))\n", "\n", "#print (to_np(lats[529, -1]))\n", "#print (to_np(lons[529, -1]))\n", "\n", "pivot_point = CoordPair(lat=38.5, lon=-97.5) \n", "angle = 90.0\n", "\n", "p_vert = vertcross(p, z, wrfin=ncfile, pivot_point=pivot_point, angle=angle, latlon=True)\n", "print (p_vert)\n", "print (\"\\n\")\n", "\n", "start_lat = lats[(lats.shape[-2]-1)//2, 0]\n", "end_lat = lats[(lats.shape[-2]-1)//2, -1]\n", "start_lon = lons[(lats.shape[-2]-1)//2, 0]\n", "end_lon = lons[(lats.shape[-2]-1)//2, -1]\n", "\n", "print (start_lat)\n", "print (end_lat)\n", "print (start_lon)\n", "print (end_lon)\n", "\n", "# Pressure using start_point and end_point\n", "start_point = CoordPair(lat=start_lat, lon=start_lon)\n", "end_point = CoordPair(lat=end_lat, lon=end_lon)\n", "\n", "p_vert = vertcross(p, z, wrfin=ncfile, start_point=start_point, end_point=end_point, latlon=True)\n", "print(p_vert)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# Pressure using pivot and angle\n", "from wrf import getvar, vertcross, CoordPair, xy_to_ll\n", "\n", "z = getvar(ncfile, \"z\")\n", "p = getvar(ncfile, \"pressure\")\n", "lats = getvar(ncfile, \"lat\")\n", "lons = getvar(ncfile, \"lon\")\n", "\n", "#print ((lats.shape[-2]-1) / 2)\n", "#print ((lats.shape[-1]-1) / 2)\n", "\n", "#print (to_np(lats[529, 899]))\n", "#print (to_np(lons[529, 899]))\n", "\n", "#print (to_np(lats[529, 0]))\n", "#print (to_np(lons[529, 0]))\n", "\n", "#print (to_np(lats[529, -1]))\n", "#print (to_np(lons[529, -1]))\n", "\n", "pivot_point = CoordPair(lat=38.5, lon=-97.5) \n", "angle = 90.0\n", "\n", "p_vert = vertcross(p, z, wrfin=ncfile, pivot_point=pivot_point, angle=angle, latlon=True)\n", "print (p_vert)\n", "print (\"\\n\")\n", "\n", "start_lat = lats[(lats.shape[-2]-1)//2, 0]\n", "end_lat = lats[(lats.shape[-2]-1)//2, -1]\n", "start_lon = lons[(lats.shape[-2]-1)//2, 0]\n", "end_lon = lons[(lats.shape[-2]-1)//2, -1]\n", "\n", "print (start_lat)\n", "print (end_lat)\n", "print (start_lon)\n", "print (end_lon)\n", "\n", "# Pressure using start_point and end_point\n", "start_point = CoordPair(lat=start_lat, lon=start_lon)\n", "end_point = CoordPair(lat=end_lat, lon=end_lon)\n", "\n", "levels = [1000., 2000., 3000.]\n", "\n", "p_vert = vertcross(p, z, wrfin=ncfile, levels=levels, start_point=start_point, end_point=end_point, latlon=True)\n", "print(p_vert)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 3.1.3 Interpolate 2D Variable to a Line" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# T2 using pivot and angle\n", "from wrf import interpline, getvar, CoordPair\n", "\n", "t2 = getvar(ncfile, \"T2\")\n", "pivot_point = CoordPair((t2.shape[-1]-1)//2, (t2.shape[-2]-1)//2) \n", "angle = 0.0\n", "\n", "t2_line = interpline(t2, pivot_point=pivot_point, angle=angle, latlon=True)\n", "print(t2_line, \"\\n\")\n", "\n", "del t2_line\n", "\n", "# T2 using start_point and end_point\n", "start_point = CoordPair((t2.shape[-1]-1)//2, 0)\n", "end_point = CoordPair((t2.shape[-1]-1)//2, -1)\n", "\n", "t2_line = interpline(t2, start_point=start_point, end_point=end_point, latlon=True)\n", "print(t2_line, \"\\n\")\n", "\n", "del t2_line\n", "\n", "t2 = getvar(ncfile, \"T2\")\n", "lats = getvar(ncfile, \"lat\")\n", "lons = getvar(ncfile, \"lon\")\n", "\n", "start_lat = lats[0, (lats.shape[-1]-1)//2]\n", "end_lat = lats[-1, (lats.shape[-1]-1)//2]\n", "start_lon = lons[0, (lons.shape[-1]-1)//2]\n", "end_lon = lons[-1, (lons.shape[-1]-1)//2]\n", "\n", "start_point = CoordPair(lat=start_lat, lon=start_lon)\n", "end_point = CoordPair(lat=end_lat, lon=end_lon)\n", "\n", "t2_line = interpline(t2, wrfin=ncfile, start_point=start_point, end_point=end_point, latlon=True)\n", "print (t2_line)\n", "\n", "del t2_line, t2\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### 3.1.4 Vertical Coordinate Interpolation" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf import vinterp, getvar\n", "\n", "# Interpolate tk to theta levels \n", "tk = getvar(ncfile, \"tk\") \n", "interp_levels = [200, 300, 500, 1000]\n", "\n", "interp_field = vinterp(ncfile, \n", " field=tk, \n", " vert_coord=\"theta\", \n", " interp_levels=interp_levels, \n", " extrapolate=True, \n", " field_type=\"tk\", \n", " log_p=True)\n", "\n", "print(interp_field)\n", "del interp_field\n", "\n", "# Interpolate tk to theta-e levels \n", " \n", "interp_levels = [200, 300, 500, 1000]\n", "\n", "interp_field = vinterp(ncfile, \n", " field=tk, \n", " vert_coord=\"eth\", \n", " interp_levels=interp_levels, \n", " extrapolate=True, \n", " field_type=\"tk\", \n", " log_p=True)\n", "\n", "print(interp_field)\n", "del interp_field\n", "\n", "# Interpolate tk to geopotential height (MSL) levels \n", " \n", "interp_levels = [30, 60, 90]\n", "\n", "interp_field = vinterp(ncfile, \n", " field=tk, \n", " vert_coord=\"ght_msl\", \n", " interp_levels=interp_levels, \n", " extrapolate=True, \n", " field_type=\"tk\", \n", " log_p=True)\n", "\n", "print(interp_field)\n", "del interp_field\n", "\n", "# Interpolate tk to geopotential height (MSL) levels \n", " \n", "interp_levels = [30, 60, 90]\n", "\n", "interp_field = vinterp(ncfile, \n", " field=tk, \n", " vert_coord=\"ght_agl\", \n", " interp_levels=interp_levels, \n", " extrapolate=True, \n", " field_type=\"tk\", \n", " log_p=True)\n", "\n", "print(interp_field)\n", "del interp_field\n", "\n", "# Interpolate tk to pressure levels\n", "interp_levels = [850, 500]\n", " \n", "interp_field = vinterp(ncfile, \n", " field=tk, \n", " vert_coord=\"pressure\", \n", " interp_levels=interp_levels, \n", " extrapolate=True, \n", " field_type=\"tk\", \n", " log_p=True)\n", "\n", "print(interp_field)\n", "del interp_field, tk\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## 3.2 Lat/Lon to X/Y Routines" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf.latlon import xy_to_ll, ll_to_xy \n", "\n", "a = xy_to_ll(ncfile, 400, 200)\n", "a1 = ll_to_xy(ncfile, a[0], a[1])\n", "\n", "#print(a)\n", "#print(\"\\n\")\n", "#print(a1)\n", "#print(\"\\n\")\n", "\n", "a = xy_to_ll(ncfile, [400,105], [200,205])\n", "a1 = ll_to_xy(ncfile, a[0,:], a[1,:])\n", "\n", "b = ll_to_xy(ncfile, 45.5, -110.8, as_int=True)\n", "\n", "# Note: Lists/Dictionaries of files will add a new dimension ('domain') only if the domain is moving\n", "c = xy_to_ll([ncfile, ncfile, ncfile], [400,105], [200,205])\n", "d = xy_to_ll({\"label1\" : [ncfile, ncfile],\n", " \"label2\" : [ncfile, ncfile]}, \n", " [400,105], [200,205])\n", "\n", "print(a)\n", "print(\"\\n\")\n", "print(a1)\n", "print(\"\\n\")\n", "print(b)\n", "print(\"\\n\")\n", "print(c)\n", "print(\"\\n\")\n", "print(d)" ] }, { "cell_type": "markdown", "metadata": { "collapsed": true }, "source": [ "# 4.0 Plotting with Cartopy" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# SLP\n", "import matplotlib.pyplot as plt\n", "from matplotlib.cm import get_cmap\n", "import cartopy.crs as crs\n", "from cartopy.feature import NaturalEarthFeature\n", "\n", "from wrf import to_np, getvar, smooth2d, get_cartopy, cartopy_xlim, cartopy_ylim, latlon_coords\n", "\n", "slp = getvar(ncfile, \"slp\")\n", "smooth_slp = smooth2d(slp, 3)\n", "lats, lons = latlon_coords(slp)\n", "\n", "cart_proj = get_cartopy(slp)\n", "\n", "fig = plt.figure(figsize=(10,10))\n", "ax = plt.axes(projection=cart_proj)\n", "\n", "states = NaturalEarthFeature(category='cultural', scale='50m', facecolor='none',\n", " name='admin_1_states_provinces_shp')\n", "ax.add_feature(states, linewidth=.5)\n", "ax.coastlines('50m', linewidth=0.8)\n", "\n", "# Can only get this to work if I manually transform the lat/lon points to projected space.\n", "xform_coords = cart_proj.transform_points(crs.PlateCarree(), to_np(lons), to_np(lats))\n", "x = xform_coords[:,:,0]\n", "y = xform_coords[:,:,1]\n", "\n", "plt.contour(x, y, to_np(smooth_slp), 10, colors=\"black\")\n", "plt.contourf(x, y, to_np(smooth_slp), 10)\n", "plt.colorbar(ax=ax, shrink=.47)\n", "\n", "ax.set_xlim(cartopy_xlim(slp))\n", "ax.set_ylim(cartopy_ylim(slp))\n", "ax.gridlines()\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# SLP\n", "from __future__ import (absolute_import, division, print_function, unicode_literals)\n", " \n", "from netCDF4 import Dataset \n", "import matplotlib.pyplot as plt\n", "from matplotlib.cm import get_cmap\n", "import cartopy.crs as crs\n", "from cartopy.feature import NaturalEarthFeature\n", "\n", "from wrf import to_np, getvar, smooth2d, get_cartopy, cartopy_xlim, cartopy_ylim, latlon_coords\n", "\n", "ncfile = Dataset(\"/Users/ladwig/Documents/wrf_files/wrfout_d01_2016-10-07_00_00_00\")\n", "\n", "# Get the sea level pressure\n", "slp = getvar(ncfile, \"slp\")\n", "\n", "# Smooth the sea level pressure since it tends to be noisey near the mountains\n", "smooth_slp = smooth2d(slp, 3)\n", "\n", "# Get the numpy array from the XLAT and XLONG coordinates\n", "lats, lons = latlon_coords(slp, as_np=True)\n", "\n", "# The cartopy() method returns a cartopy.crs projection object\n", "cart_proj = get_cartopy(slp)\n", "\n", "# Create a figure that's 10x10\n", "fig = plt.figure(figsize=(10,10))\n", "# Get the GeoAxes set to the projection used by WRF\n", "ax = plt.axes(projection=cart_proj)\n", "\n", "# Download and add the states and coastlines\n", "states = NaturalEarthFeature(category='cultural', scale='50m', facecolor='none',\n", " name='admin_1_states_provinces_shp')\n", "ax.add_feature(states, linewidth=.5)\n", "ax.coastlines('50m', linewidth=0.8)\n", "\n", "# Make the contour outlines and filled contours for the smoothed sea level pressure.\n", "# The transform keyword indicates that the lats and lons arrays are lat/lon coordinates and tells \n", "# cartopy to transform the points in to grid space.\n", "plt.contour(lons, lats, to_np(smooth_slp), 10, colors=\"black\", transform=crs.PlateCarree())\n", "plt.contourf(lons, lats, to_np(smooth_slp), 10, transform=crs.PlateCarree())\n", "\n", "# Add a color bar\n", "plt.colorbar(ax=ax, shrink=.47)\n", "\n", "# Set the map limits\n", "ax.set_xlim(cartopy_xlim(slp))\n", "ax.set_ylim(cartopy_ylim(slp))\n", "\n", "# Add the gridlines\n", "ax.gridlines()" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "# 500 MB Heights and Winds\n", "\n", "import matplotlib.pyplot as plt\n", "from matplotlib.cm import get_cmap\n", "import cartopy.crs as crs\n", "from cartopy.feature import NaturalEarthFeature\n", "\n", "from wrf import getvar, interplevel, to_np, get_cartopy, cartopy_xlim, cartopy_ylim, latlon_coords\n", "\n", "\n", "p = getvar(ncfile, \"pressure\")\n", "z = getvar(ncfile, \"z\", units=\"dm\")\n", "ua = getvar(ncfile, \"ua\", units=\"kts\")\n", "va = getvar(ncfile, \"va\", units=\"kts\")\n", "\n", "ht_500 = interplevel(z, p, 500)\n", "u_500 = interplevel(ua, p, 500)\n", "v_500 = interplevel(va, p, 500)\n", "\n", "lats, lons = latlon_coords(ht_500)\n", "\n", "cart_proj = get_cartopy(slp)\n", "\n", "fig = plt.figure(figsize=(20,20))\n", "ax = plt.axes([0.1,0.1,0.8,0.8], projection=cart_proj)\n", "\n", "states = NaturalEarthFeature(category='cultural', scale='50m', facecolor='none',\n", " name='admin_1_states_provinces_shp')\n", "ax.add_feature(states, linewidth=0.5)\n", "ax.coastlines('50m', linewidth=0.8)\n", "\n", "# Can only get this to work if I manually transform the lat/lon points to projected space.\n", "xform_coords = cart_proj.transform_points(crs.PlateCarree(), to_np(lons), to_np(lats))\n", "x = xform_coords[:,:,0]\n", "y = xform_coords[:,:,1]\n", "\n", "plt.contour(x, y, to_np(ht_500), 20, cmap=get_cmap(\"plasma\"))\n", "plt.barbs(x[::50,::50], y[::50,::50], to_np(u_500[::50, ::50]), to_np(v_500[::50, ::50]))\n", "plt.colorbar(ax=ax, shrink=.7)\n", "\n", "ax.set_xlim(cartopy_xlim(slp))\n", "ax.set_ylim(cartopy_ylim(slp))\n", "ax.gridlines()\n", "\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false, "scrolled": false }, "outputs": [], "source": [ "# Cross-section of pressure using xarray's builtin plotting\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "from matplotlib.cm import get_cmap\n", "\n", "from wrf import getvar, vertcross, to_np, CoordPair\n", "\n", "p = getvar(ncfile, \"pressure\")\n", "z = getvar(ncfile, \"z\", units=\"dm\")\n", "\n", "pivot_point = CoordPair(z.shape[-1] / 2, z.shape[-2] / 2) \n", "angle = 90.0\n", "\n", "p_vert = vertcross(p, z, pivot_point=pivot_point, angle=angle, levels=[1000,850,500])\n", "\n", "fig = plt.figure(figsize=(20,8))\n", "ax = plt.axes([0.1,0.1,0.8,0.8])\n", "\n", "p_vert.plot.contour(ax=ax, levels=[0 + 50*n for n in range(20)], cmap=get_cmap(\"viridis\"))\n", "\n", "\n" ] }, { "cell_type": "markdown", "metadata": { "collapsed": false }, "source": [ "# Multi-time Moving Domain Files" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import os\n", "from wrf import getvar, ALL_TIMES\n", "from netCDF4 import Dataset as nc\n", "\n", "dir = \"/Users/ladwig/Documents/wrf_files/wrf_vortex_multi\"\n", "ncfilenames = [os.path.join(dir, x) for x in os.listdir(dir) if x.find(\"_d02_\") > 0]\n", "ncfiles = [nc(x) for x in ncfilenames]\n", "\n", "\n", "#print (ncfiles[0].variables[\"XLONG\"][0,0,-1], ncfiles[0].variables[\"XLONG\"][-1,0,-1])\n", "#print (ncfiles[1].variables[\"XLONG\"][0,0,-1], ncfiles[1].variables[\"XLONG\"][-1,0,-1])\n", "#print (ncfiles[-1].variables[\"XLONG\"][0,0,-1], ncfiles[-1].variables[\"XLONG\"][-1,0,-1])\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "p = getvar(ncfiles, \"ctt\", timeidx=ALL_TIMES)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print (p)\n", "#print (p.attrs[\"projection\"].shape)\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print (p.attrs[\"projection\"])\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "ncfiles[2].variables[\"XTIME\"][:]\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "p = getvar(ncfiles, \"P\", timeidx=None, method=\"cat\", meta=True, squeeze=True)\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print (p)" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print (type(p.coords[\"Time\"]))" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import datetime\n", "import pandas\n", "print (type(p.coords[\"Time\"].values.astype(datetime.datetime)))\n", "print (repr(datetime.datetime.utcfromtimestamp(p.coords[\"Time\"][0].values.astype(int) * 1E-9)))\n", "print (pandas.to_datetime(p.coords[\"Time\"].values))" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "wrf_vars = [\"avo\", \"eth\", \"cape_2d\", \"cape_3d\", \"ctt\", \"dbz\", \"mdbz\", \n", " \"geopt\", \"helicity\", \"lat\", \"lon\", \"omg\", \"p\", \"pressure\", \n", " \"pvo\", \"pw\", \"rh2\", \"rh\", \"slp\", \"ter\", \"td2\", \"td\", \"tc\",\n", " \"theta\", \"tk\", \"tv\", \"twb\", \"updraft_helicity\", \"ua\", \"va\", \n", " \"wa\", \"uvmet10\", \"uvmet\", \"z\", \"ctt\"]\n", "#wrf_vars = [\"cape_2d\"]\n", "\n", "vard = {}\n", "for varname in wrf_vars:\n", " print (varname)\n", " vard[varname] = getvar(ncfiles, varname, timeidx=None, method=\"cat\", squeeze=False)\n", " \n", "#vard = {varname: getvar(ncfiles, varname, method=\"join\", squeeze=False) for varname in wrf_vars}\n", "for varname in wrf_vars:\n", " print(vard[varname])\n", " " ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import os\n", "from wrf import getvar\n", "from netCDF4 import Dataset as nc\n", "\n", "dir = \"/Users/ladwig/Documents/wrf_files/wrf_vortex_multi\"\n", "ncfilenames = [os.path.join(dir, x) for x in os.listdir(dir) if x.find(\"_d02_\") > 0]\n", "ncfiles = [nc(x) for x in ncfilenames]\n", "\n", "# Pressure using pivot and angle\n", "from wrf import getvar, vertcross, CoordPair\n", "\n", "timeidx = 0\n", "z = getvar(ncfiles, \"z\", timeidx, method=\"join\")\n", "p = getvar(ncfiles, \"pressure\", timeidx, method=\"join\")\n", "pivot_point = CoordPair(z.shape[-1] / 2, z.shape[-2] / 2) \n", "angle = 40.0\n", "\n", "p_vert = vertcross(p, z, pivot_point=pivot_point, angle=angle)\n", "print(p_vert)\n", "print (\"\\n\")\n", "del p_vert\n", "\n", "# Pressure using start_point and end_point\n", "start_point = CoordPair(0, z.shape[-2]/2)\n", "end_point = CoordPair(-1, z.shape[-2]/2)\n", "\n", "p_vert = vertcross(p, z, start_point=start_point, end_point=end_point)\n", "print(p_vert)\n", "del p_vert, p, z" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import os\n", "from wrf import getvar\n", "from netCDF4 import Dataset as nc\n", "\n", "dir = \"/Users/ladwig/Documents/wrf_files/wrf_vortex_multi\"\n", "ncfilenames = [os.path.join(dir, x) for x in os.listdir(dir) if x.find(\"_d02_\") > 0]\n", "ncfiles = [nc(x) for x in ncfilenames]\n", "\n", "timeidx = None\n", "\n", "# T2 using pivot and angle\n", "from wrf import interpline, getvar, to_np, CoordPair\n", "\n", "t2 = getvar(ncfiles, \"T2\", timeidx)\n", "pivot_point = CoordPair(t2.shape[-2] / 2, t2.shape[-1] / 2) \n", "angle = 90.0\n", "\n", "t2_line = interpline(t2, pivot_point=pivot_point, angle=angle, latlon=True)\n", "print(t2_line)\n", "print(\"\\n\")\n", "\n", "\n", "del t2_line\n", "\n", "# T2 using start_point and end_point\n", "start_point = CoordPair(t2.shape[-2]/2, 0)\n", "end_point = CoordPair(t2.shape[-2]/2, -1)\n", "\n", "t2_line = interpline(t2, start_point=start_point, end_point=end_point, latlon=True)\n", "print(t2_line)\n", "print(\"\\n\")\n", "\n", "del t2_line, t2" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "from wrf import getvar\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [ "getvar?" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from wrf.latlon import xy_to_ll, ll_to_xy \n", "\n", "a = xy_to_ll(ncfiles, 400, 200)\n", "a = xy_to_ll(ncfiles, [400,105], [200,205])\n", "b = ll_to_xy(ncfiles, 45.5, -110.8, as_int=True)\n", "\n", "# Note: Lists/Dictionaries of files will add a new dimension ('domain') only if the domain is moving\n", "c = xy_to_ll(ncfiles, [400,105], [200,205])\n", "d = xy_to_ll({\"label1\" : ncfiles,\n", " \"label2\" : ncfiles}, \n", " [400,105], [200,205])\n", "\n", "print(a)\n", "print(\"\\n\")\n", "print(b)\n", "print(\"\\n\")\n", "print(c)\n", "print(\"\\n\")\n", "print(d)\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": true }, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Python 2", "language": "python", "name": "python2" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 2 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython2", "version": "2.7.12" } }, "nbformat": 4, "nbformat_minor": 0 }