xWRF overview#

xWRF is a package designed to make the post-processing of WRF output data more pythonic. It’s aim is to smooth the rough edges around the unique, non CF-compliant WRF output data format and make the data accessible to utilities like dask and the wider Pangeo universe.

It is built as an Accessor on top of xarray, providing a very simple user interface.

Examining the data#

When opening up a normal WRF output file with the simple xarray netcdf backend, one can see that it does not provide a lot of useful information.

import xwrf

ds_old = xwrf.tutorial.open_dataset("wrfout")
ds_old
<xarray.Dataset> Size: 89MB
Dimensions:   (Time: 1, south_north: 340, west_east: 270, bottom_top: 39,
               west_east_stag: 271, south_north_stag: 341)
Coordinates:
    XLAT      (Time, south_north, west_east) float32 367kB ...
    XLONG     (Time, south_north, west_east) float32 367kB ...
    XTIME     (Time) datetime64[ns] 8B ...
    XLAT_U    (Time, south_north, west_east_stag) float32 369kB ...
    XLONG_U   (Time, south_north, west_east_stag) float32 369kB ...
    XLAT_V    (Time, south_north_stag, west_east) float32 368kB ...
    XLONG_V   (Time, south_north_stag, west_east) float32 368kB ...
Dimensions without coordinates: Time, south_north, west_east, bottom_top,
                                west_east_stag, south_north_stag
Data variables:
    ZNU       (Time, bottom_top) float32 156B ...
    Times     (Time) |S19 19B ...
    U         (Time, bottom_top, south_north, west_east_stag) float32 14MB ...
    V         (Time, bottom_top, south_north_stag, west_east) float32 14MB ...
    SINALPHA  (Time, south_north, west_east) float32 367kB ...
    COSALPHA  (Time, south_north, west_east) float32 367kB ...
    T         (Time, bottom_top, south_north, west_east) float32 14MB ...
    QVAPOR    (Time, bottom_top, south_north, west_east) float32 14MB ...
    P         (Time, bottom_top, south_north, west_east) float32 14MB ...
    PB        (Time, bottom_top, south_north, west_east) float32 14MB ...
    PSN       (Time, south_north, west_east) float32 367kB ...
Attributes: (12/149)
    TITLE:                            OUTPUT FROM WRF V4.1.3 MODEL
    START_DATE:                      2099-08-01_00:00:00
    SIMULATION_START_DATE:           2099-08-01_00:00:00
    WEST-EAST_GRID_DIMENSION:        271
    SOUTH-NORTH_GRID_DIMENSION:      341
    BOTTOM-TOP_GRID_DIMENSION:       40
    ...                              ...
    ISLAKE:                          21
    ISICE:                           15
    ISURBAN:                         13
    ISOILWATER:                      14
    HYBRID_OPT:                      0
    ETAC:                            0.0

While all variables are present, e.g. the information about the projection is still in the metadata and also for some fields, there are non-metpy compliant units attributes.

So let’s try to use the standard xWRF.postprocess() function in order to make this information useable.

ds_new = xwrf.tutorial.open_dataset("wrfout").xwrf.postprocess()
ds_new
<xarray.Dataset> Size: 104MB
Dimensions:                    (y: 340, x: 270, Time: 1, z: 39, x_stag: 271,
                                y_stag: 341)
Coordinates: (12/13)
    XLAT                       (y, x) float32 367kB 22.27 22.31 ... 57.54 57.58
    XLONG                      (y, x) float32 367kB -116.7 -116.6 ... -110.9
    XTIME                      (Time) datetime64[ns] 8B ...
    XLAT_U                     (y, x_stag) float32 369kB 22.25 22.29 ... 57.6
    XLONG_U                    (y, x_stag) float32 369kB -116.7 ... -110.8
    XLAT_V                     (y_stag, x) float32 368kB 22.23 22.28 ... 57.62
    ...                         ...
  * z                          (z) float32 156B 0.9969 0.9899 ... 0.002948
  * Time                       (Time) datetime64[ns] 8B 2099-10-01
  * y_stag                     (y_stag) float64 3kB -3.386e+05 ... 2.721e+06
  * y                          (y) float64 3kB -3.341e+05 ... 2.717e+06
  * x                          (x) float64 2kB -4.728e+06 ... -2.307e+06
  * x_stag                     (x_stag) float64 2kB -4.733e+06 ... -2.303e+06
Data variables:
    Times                      (Time) |S19 19B b'2099-10-01_00:00:00'
    U                          (Time, z, y, x_stag) float32 14MB ...
    V                          (Time, z, y_stag, x) float32 14MB ...
    SINALPHA                   (Time, y, x) float32 367kB 0.5507 ... 0.489
    COSALPHA                   (Time, y, x) float32 367kB 0.8347 ... 0.8723
    QVAPOR                     (Time, z, y, x) float32 14MB ...
    PSN                        (Time, y, x) float32 367kB ...
    air_potential_temperature  (Time, z, y, x) float32 14MB 302.2 ... 501.3
    air_pressure               (Time, z, y, x) float32 14MB 1.009e+05 ... 5.2...
    wind_east                  (Time, z, y, x) float32 14MB -0.8541 ... 12.5
    wind_north                 (Time, z, y, x) float32 14MB -4.625 ... -0.8308
    wrf_projection             object 8B +proj=lcc +x_0=0 +y_0=0 +a=6370000 +...
Attributes: (12/149)
    TITLE:                            OUTPUT FROM WRF V4.1.3 MODEL
    START_DATE:                      2099-08-01_00:00:00
    SIMULATION_START_DATE:           2099-08-01_00:00:00
    WEST-EAST_GRID_DIMENSION:        271
    SOUTH-NORTH_GRID_DIMENSION:      341
    BOTTOM-TOP_GRID_DIMENSION:       40
    ...                              ...
    ISLAKE:                          21
    ISICE:                           15
    ISURBAN:                         13
    ISOILWATER:                      14
    HYBRID_OPT:                      0
    ETAC:                            0.0

As you see, xWRF added some coordinate data, reassigned some dimensions and generally increased the amount of information available in the dataset.

Projection treatment#

xWRF has determined the correct projection from the netCDF metadata and has added a wrf_projection variable to the dataset storing the pyproj projection object.

ds_new['wrf_projection'].item()
<Projected CRS: +proj=lcc +x_0=0 +y_0=0 +a=6370000 +b=6370000 +lat ...>
Name: unknown
Axis Info [cartesian]:
- E[east]: Easting (metre)
- N[north]: Northing (metre)
Area of Use:
- undefined
Coordinate Operation:
- name: unknown
- method: Lambert Conic Conformal (2SP)
Datum: unknown
- Ellipsoid: unknown
- Prime Meridian: Greenwich

Using this projection xWRF has also calculated the regular model grid, which the WRF simulation is performed on and which can be used for e.g. bilinear interpolation.

ds_new[['x', 'y']]
<xarray.Dataset> Size: 739kB
Dimensions:  (x: 270, y: 340)
Coordinates:
  * x        (x) float64 2kB -4.728e+06 -4.719e+06 ... -2.316e+06 -2.307e+06
  * y        (y) float64 3kB -3.341e+05 -3.251e+05 ... 2.708e+06 2.717e+06
    XLAT     (y, x) float32 367kB 22.27 22.31 22.35 22.4 ... 57.5 57.54 57.58
    XLONG    (y, x) float32 367kB -116.7 -116.6 -116.6 ... -111.2 -111.0 -110.9
Data variables:
    *empty*
Attributes: (12/149)
    TITLE:                            OUTPUT FROM WRF V4.1.3 MODEL
    START_DATE:                      2099-08-01_00:00:00
    SIMULATION_START_DATE:           2099-08-01_00:00:00
    WEST-EAST_GRID_DIMENSION:        271
    SOUTH-NORTH_GRID_DIMENSION:      341
    BOTTOM-TOP_GRID_DIMENSION:       40
    ...                              ...
    ISLAKE:                          21
    ISICE:                           15
    ISURBAN:                         13
    ISOILWATER:                      14
    HYBRID_OPT:                      0
    ETAC:                            0.0

Attribute changes#

xWRF adds additional attributes to variables in order to make them CF- and COMODO-compliant. It also amends unit attributes to work with metpy units, enabling a seamless integration with the Pangeo software stack.

Here, for example the x-wind component gets the correct CF standard_name and a COMODO grid_mapping attribute indicating the respective projection.

ds_old['U'].attrs, ds_new['U'].attrs
({'FieldType': 104,
  'MemoryOrder': 'XYZ',
  'description': 'x-wind component',
  'units': 'm s-1',
  'stagger': 'X'},
 {'FieldType': 104,
  'MemoryOrder': 'XYZ',
  'description': 'x-wind component',
  'units': 'm s-1',
  'stagger': 'X',
  'standard_name': 'x_wind',
  'grid_mapping': 'wrf_projection'})

Also, the units attribute of the PSN variable was cleaned up to conform to metpy unit conventions.

Note

As of now, unit translations are implemented on a manual basis, so please raise an issue with us if you encounter any problems in this regard. In the future, this will be implemented in a more structured manner.

ds_old['PSN'].attrs['units'], ds_new['PSN'].attrs['units']
('umol co2/m2/s', 'umol m-2 s-1')
import metpy

try:
    ds_old['PSN'].metpy.quantify()
except metpy.units.UndefinedUnitError as e:
    print(e)
ds_new['PSN'].metpy.quantify()
'co' is not defined in the unit registry
<xarray.DataArray 'PSN' (Time: 1, y: 340, x: 270)> Size: 367kB
<Quantity([[[0.        0.        0.        ... 0.        0.        0.       ]
  [0.        0.        0.        ... 0.        0.        0.       ]
  [0.        0.        0.        ... 0.        0.        0.       ]
  ...
  [0.        0.        0.        ... 2.0479841 2.014872  2.0163262]
  [0.        0.        0.        ... 2.122164  2.3794565 2.0393817]
  [0.        0.        0.        ... 2.0216792 1.9838098 1.9889725]]], 'micromole / meter ** 2 / second')>
Coordinates:
    XLAT     (y, x) float32 367kB 22.27 22.31 22.35 22.4 ... 57.5 57.54 57.58
    XLONG    (y, x) float32 367kB -116.7 -116.6 -116.6 ... -111.2 -111.0 -110.9
    XTIME    (Time) datetime64[ns] 8B ...
  * Time     (Time) datetime64[ns] 8B 2099-10-01
  * y        (y) float64 3kB -3.341e+05 -3.251e+05 ... 2.708e+06 2.717e+06
  * x        (x) float64 2kB -4.728e+06 -4.719e+06 ... -2.316e+06 -2.307e+06
Attributes:
    FieldType:     104
    MemoryOrder:   XY 
    description:   total photosynthesis
    stagger:       
    grid_mapping:  wrf_projection

Diagnostic variables#

Because some WRF output fields are quite raw, essential diagnostic variables like air_pressure or air_potential_temperature are missing. Also wind fields are natively in grid-space and not oriented in west-east, north-south direction, so the wind_east and wind_north fields are added. These useful variables get added by xWRF by default. Users can choose to keep the fields after the computation of diagnostics is done by using .xwrf.postprocess(drop_diagnostic_variable_components=False), however grid-relative winds are always retained.

ds_new['air_pressure']
<xarray.DataArray 'air_pressure' (Time: 1, z: 39, y: 340, x: 270)> Size: 14MB
array([[[[100893.734 , 100889.61  , 100888.09  , ...,  99661.97  ,
           99665.875 ,  99651.086 ],
         [100893.55  , 100891.305 , 100887.91  , ...,  99460.77  ,
           99437.125 ,  99507.88  ],
         [100896.195 , 100891.875 , 100890.65  , ...,  99259.664 ,
           99260.1   ,  99446.16  ],
         ...,
         [101696.73  , 101671.125 , 101662.305 , ...,  95169.67  ,
           95047.93  ,  94998.36  ],
         [101684.664 , 101665.67  , 101655.04  , ...,  95912.805 ,
           95763.6   ,  95652.57  ],
         [101672.53  , 101665.94  , 101659.2   , ...,  96404.76  ,
           96312.66  ,  96160.69  ]],

        [[100203.97  , 100202.63  , 100201.086 , ...,  98989.78  ,
           98996.55  ,  98979.    ],
         [100206.67  , 100205.6   , 100202.766 , ...,  98787.086 ,
           98764.5   ,  98832.36  ],
         [100209.414 , 100205.04  , 100203.88  , ...,  98584.98  ,
           98589.43  ,  98770.89  ],
...
         [  5886.9263,   5886.67  ,   5886.6123, ...,   5827.5776,
            5825.9956,   5821.3135],
         [  5886.796 ,   5886.5864,   5886.468 , ...,   5833.852 ,
            5831.8096,   5818.641 ],
         [  5886.6597,   5886.6133,   5886.5464, ...,   5813.8115,
            5815.558 ,   5817.2427]],

        [[  5281.6   ,   5281.608 ,   5281.594 , ...,   5275.668 ,
            5275.5195,   5275.705 ],
         [  5281.6064,   5281.692 ,   5281.7773, ...,   5277.55  ,
            5277.2354,   5277.1577],
         [  5281.6157,   5281.8145,   5281.881 , ...,   5276.886 ,
            5276.828 ,   5277.4326],
         ...,
         [  5284.9956,   5284.882 ,   5284.917 , ...,   5266.0586,
            5265.268 ,   5260.8257],
         [  5284.9634,   5284.8643,   5284.846 , ...,   5267.7246,
            5266.5693,   5254.022 ],
         [  5284.9243,   5284.8877,   5284.871 , ...,   5244.3447,
            5246.701 ,   5249.4062]]]], dtype=float32)
Coordinates:
    XLAT     (y, x) float32 367kB 22.27 22.31 22.35 22.4 ... 57.5 57.54 57.58
    XLONG    (y, x) float32 367kB -116.7 -116.6 -116.6 ... -111.2 -111.0 -110.9
    XTIME    (Time) datetime64[ns] 8B ...
  * z        (z) float32 156B 0.9969 0.9899 0.981 ... 0.0161 0.009174 0.002948
  * Time     (Time) datetime64[ns] 8B 2099-10-01
  * y        (y) float64 3kB -3.341e+05 -3.251e+05 ... 2.708e+06 2.717e+06
  * x        (x) float64 2kB -4.728e+06 -4.719e+06 ... -2.316e+06 -2.307e+06
Attributes:
    units:          Pa
    standard_name:  air_pressure
    grid_mapping:   wrf_projection
ds_new['wind_east']
<xarray.DataArray 'wind_east' (Time: 1, z: 39, y: 340, x: 270)> Size: 14MB
array([[[[-0.8540844 , -0.738251  , -0.545264  , ...,  1.4950651 ,
           1.7660816 ,  1.9395261 ],
         [-0.8818562 , -0.74259746, -0.5333221 , ...,  1.348903  ,
           1.5153437 ,  1.7764862 ],
         [-0.9169271 , -0.7991271 , -0.6141951 , ...,  1.3527424 ,
           1.460252  ,  1.6192899 ],
         ...,
         [ 9.255407  ,  9.28659   ,  9.1721    , ...,  4.580781  ,
           4.4629016 ,  3.488285  ],
         [ 9.394391  ,  9.397248  ,  9.32838   , ...,  5.3224483 ,
           4.4839735 ,  3.1068704 ],
         [ 9.494457  ,  9.411124  ,  9.363975  , ...,  5.2005925 ,
           4.8733425 ,  3.7735646 ]],

        [[-0.87797356, -0.7614403 , -0.59394455, ...,  1.882653  ,
           2.3743167 ,  2.5777702 ],
         [-0.90974176, -0.77123356, -0.5847976 , ...,  1.7839468 ,
           2.272645  ,  2.497768  ],
         [-0.9420047 , -0.8337804 , -0.67352176, ...,  1.8757367 ,
           2.2196953 ,  2.3774064 ],
...
         [14.407449  , 14.719514  , 14.945457  , ..., 11.257124  ,
          11.583581  , 14.111986  ],
         [14.562296  , 14.829544  , 15.028374  , ..., 11.101351  ,
          11.73059   , 14.226654  ],
         [14.707831  , 14.862129  , 15.002537  , ..., 12.704432  ,
          12.850094  , 14.030186  ]],

        [[-1.1588178 , -1.0782948 , -1.0017282 , ...,  1.8423607 ,
           1.7980661 ,  1.6273091 ],
         [-1.1129029 , -1.0681975 , -0.9338307 , ...,  3.0397072 ,
           2.6984587 ,  1.9549323 ],
         [-1.0991675 , -1.181081  , -1.0825561 , ...,  3.7225795 ,
           3.3059304 ,  2.319018  ],
         ...,
         [12.773672  , 13.064896  , 13.278723  , ...,  9.394521  ,
           9.826231  , 12.684486  ],
         [12.933495  , 13.184     , 13.374461  , ...,  9.11593   ,
           9.9368305 , 12.795868  ],
         [13.088163  , 13.237101  , 13.3752    , ..., 10.920462  ,
          11.129587  , 12.496285  ]]]], dtype=float32)
Coordinates:
    XLAT     (y, x) float32 367kB 22.27 22.31 22.35 22.4 ... 57.5 57.54 57.58
    XLONG    (y, x) float32 367kB -116.7 -116.6 -116.6 ... -111.2 -111.0 -110.9
    XTIME    (Time) datetime64[ns] 8B ...
  * z        (z) float32 156B 0.9969 0.9899 0.981 ... 0.0161 0.009174 0.002948
  * Time     (Time) datetime64[ns] 8B 2099-10-01
  * y        (y) float64 3kB -3.341e+05 -3.251e+05 ... 2.708e+06 2.717e+06
  * x        (x) float64 2kB -4.728e+06 -4.719e+06 ... -2.316e+06 -2.307e+06
Attributes:
    description:    earth-relative x-wind component
    standard_name:  eastward_wind
    units:          m s-1
    grid_mapping:   wrf_projection