Welcome to wrf-python’s documentation!

Contents:

Extraction Routine

wrf.getvar(wrfnc, varname, timeidx=0, method=u'cat', squeeze=True, cache=None, meta=True, **kargs)

Returns basic diagnostics from the WRF ARW model output.

Below is a table of available diagnostics.

Variable Name Description Units Additional Keyword Arguments
avo Absolute Vorticity 10-5 s-1  
eth/theta_e Equivalent Potential Temperature K  
cape_2d 2D cape (mcape/mcin/lcl/lfc) J/kg / J/kg / m / m missing: Fill value for output only (float)
cape_3d 3D cape and cin J/kg missing: Fill value for output only (float)
ctt Cloud Top Temperature C  
dbz Reflectivity dBz do_variant: Set to True to enable variant calculation. Default is False. do_liqskin : Set to True to enable liquid skin calculation. Default is False.
mdbz Maximum Reflectivity dBz do_variant: Set to True to enable variant calculation. Default is False. do_liqskin: Set to True to enable liquid skin calculation. Default is False.
geopt/geopotential Full Model Geopotential m2 s-2  
helicity Storm Relative Helicity m-2/s-2 top: The top level for the calculation in meters (float). Default is 3000.0.
lat Latitude decimal degrees  
lon Longitude decimal degrees  
omg/omega Omega Pa/s  
p/pres Full Model Pressure Pa  
pressure Full Model Pressure hPa  
pvo Potential Vorticity PVU  
pw Precipitable Water kg m-2  
rh2 2m Relative Humidity %  
slp Sea Level Pressure hPa  
ter Model Terrain Height m  
td2 2m Dew Point Temperature C  
td Dew Point Temperature C  
tc Temperature C  
th/theta Potential Temperature K  
tk Temperature K  
times Times in the File or Sequence    
tv Virtual Temperature K  
twb Wet Bulb Temperature K  
updraft_helicity Updraft Helicity m-2/s-2 bottom: The bottom level for the calculation in meters (float). Default is 2000.0. top: The top level for the calculation in meters (float). Default is 5000.0.
ua U-component of Wind on Mass Points m/s  
va V-component of Wind on Mass Points m/s  
wa W-component of Wind on Mass Points m/s  
uvmet10 10 m U and V Components of Wind Rotated to Earth Coordinates m/s  
uvmet U and V Components of Wind Rotated to Earth Coordinates m/s  
z/height Full Model Height m msl: Set to False to return AGL values. Otherwise, MSL. Default is True.
Parameters:
  • wrfnc (netCD4F.Dataset, Nio.NioFile, or a sequence) – Input WRF ARW NetCDF data as a netCDF4.Dataset, Nio.NioFile or an iterable sequence of the aforementioned types.
  • varname (str) – The variable name.
  • timeidx (int or wrf.ALL_TIMES, optional) – The desired time index. This value can be a positive integer, negative integer, or wrf.ALL_TIMES (an alias for None) to return all times in the file or sequence. The default is 0.
  • method ({'cat', 'join'}, optional) – The aggregation method to use for sequences, either ‘cat’ or ‘join’. ‘cat’ combines the data along the Time index. ‘join’ is creates a new index for the file. The default is ‘cat’.
  • squeeze ({True, False}, optional) – Set to False to prevent dimensions with a size of 1 from being removed from the shape of the output. Default is True.
  • cache (dict, optional) – A dictionary of (varname, ndarray) which can be used to supply pre-extracted NetCDF variables to the computational routines. This can be used to prevent the repeated variable extraction from large sequences of data files. Default is None.
  • meta ({True, False}, optional) – Set to False to disable metadata and return numpy.ndarray instead of xarray.DataArray. Default is True.
Returns:

Returns the specified diagnostics output. If xarray is enabled and the meta parameter is True, then the result will be an xarray.DataArray object. Otherwise, the result will be a numpy.ndarray object with no metadata.

Return type:

xarray.DataArray or numpy.ndarray

Interpolation Routines

wrf.interplevel(field3d, z, desiredlev, missing=9.969209968386869e+36, meta=True)

Interpolates a three-dimensional field specified pressure or height level.

Parameters:
  • field3d (xarray.DataArray or numpy.ndarray) – A three-dimensional field.
  • z (xarray.DataArray or numpy.ndarray) – A three-dimensional array for the vertical coordinate, typically pressure or height.
  • desiredlev (float) – The desired vertical level. Must be in the same units as the z parameter.
  • missing (float) – The fill value to use for the output. Default is wrf.Constants.DEFAULT_FILL.
  • meta ({True, False}) – Set to False to disable metadata and return numpy.ndarray instead of xarray.DataArray. Default is True.
Returns:

Returns the interpolated variable. If xarray is enabled and the meta parameter is True, then the result will be an xarray.DataArray object. Otherwise, the result will be a numpy.ndarray object with no metadata.

Return type:

xarray.DataArray or numpy.ndarray

wrf.vertcross(field3d, z, missing=9.969209968386869e+36, pivot_point=None, angle=None, start_point=None, end_point=None, include_latlon=False, cache=None, meta=True)

Return the vertical cross section for a 3D field, interpolated to a verical plane defined by a horizontal line.

The horizontal line is defined by either including the pivot_point and angle parameters, or the start_point and end_point parameters.

Parameters:
  • field3d (xarray.DataArray or numpy.ndarray) – A three-dimensional field.
  • z (xarray.DataArray or numpy.ndarray) – A three-dimensional array for the vertical coordinate, typically pressure or height
  • pivot_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the x,y location through which the plane will pass. Must also specify angle.
  • angle (float, optional) – Only valid for cross sections where a plane will be plotted through a given point on the model domain. 0.0 represents a S-N cross section and 90.0 a W-E cross section.
  • start_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the start x,y location through which the plane will pass.
  • end_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the end x,y location through which the plane will pass.
  • include_latlon ({True, False}, optional) – Set to True to also interpolate the two-dimensional latitude and longitude coordinates along the same horizontal line and include this information in the metadata (if enabled). This can be helpful for plotting. Default is False.
  • cache (dict, optional) – A dictionary of (varname, numpy.ndarray) pairs which can be used to supply pre-extracted NetCDF variables to the computational routines. This can be used to prevent the repeated variable extraction from large sequences of data files. Default is None.
  • meta ({True, False}, optional) – Set to False to disable metadata and return numpy.ndarray instead of xarray.DataArray. Default is True.
Returns:

Returns the interpolated variable. If xarray is enabled and the meta parameter is True, then the result will be an xarray.DataArray object. Otherwise, the result will be a numpy.ndarray object with no metadata.

Return type:

xarray.DataArray or numpy.ndarray

wrf.interpline(field2d, pivot_point=None, angle=None, start_point=None, end_point=None, include_latlon=False, cache=None, meta=True)

Return the two-dimensional field interpolated along a line.

Parameters:
  • field2d (xarray.DataArray or numpy.ndarray) – A two-dimensional field.
  • pivot_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the x,y location through which the plane will pass. Must also specify angle.
  • angle (float, optional) – Only valid for cross sections where a plane will be plotted through a given point on the model domain. 0.0 represents a S-N cross section and 90.0 a W-E cross section.
  • start_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the start x,y location through which the plane will pass.
  • end_point (tuple or list, optional) – A tuple or list with two entries, in the form of [x, y] (or [west_east, south_north]), which indicates the end x,y location through which the plane will pass.
  • include_latlon ({True, False}, optional) – Set to True to also interpolate the two-dimensional latitude and longitude coordinates along the same horizontal line and include this information in the metadata (if enabled). This can be helpful for plotting. Default is False.
  • cache (dict, optional) – A dictionary of (varname, numpy.ndarray) pairs which can be used to supply pre-extracted NetCDF variables to the computational routines. This can be used to prevent the repeated variable extraction from large sequences of data files. Default is None.
  • meta ({True, False}, optional) – Set to False to disable metadata and return numpy.ndarray instead of xarray.DataArray. Default is True.
Returns:

Returns the interpolated variable. If xarray is enabled and the meta parameter is True, then the result will be an xarray.DataArray object. Otherwise, the result will be a numpy.ndarray object with no metadata.

Return type:

xarray.DataArray or numpy.ndarray

wrf.vinterp(wrfnc, field, vert_coord, interp_levels, extrapolate=False, field_type=None, log_p=False, timeidx=0, method=u'cat', squeeze=True, cache=None, meta=True)

Return the field vertically interpolated to the given the type of surface and a set of new levels.

Parameters:
  • wrfnc (netCD4F.Dataset, Nio.NioFile, or a sequence) – Input WRF ARW NetCDF data as a netCDF4.Dataset, Nio.NioFile or an iterable sequence of the aforementioned types.
  • field2d (xarray.DataArray or numpy.ndarray) – A two-dimensional field.
  • vert_coord ({'pressure', 'pres', 'p', 'ght_msl', 'ght_agl', 'theta', 'th', 'theta-e', 'thetae', 'eth'}) –

    A string indicating the vertical coordinate type to interpolate to.

    Valid strings are:
    • ‘pressure’, ‘pres’, ‘p’: pressure [hPa]
    • ‘ght_msl’: grid point height msl [km]
    • ‘ght_agl’: grid point height agl [km]
    • ‘theta’, ‘th’: potential temperature [K]
    • ‘theta-e’, ‘thetae’, ‘eth’: equivalent potential temperature [K]
  • interp_levels (sequence) – A 1D sequence of vertical levels to interpolate to.
  • extrapolate ({True, False}, optional) – Set to True to extrapolate values below ground. Default is False.
  • field_type ({'none', 'pressure', 'pres', 'p', 'z', 'tc', 'tk', 'theta', 'th', 'theta-e', 'thetae', 'eth', 'ght'}, optional) – The type of field. Default is None.
  • log_p ({True, False}) – Use the log of the pressure for interpolation instead of just pressure. Default is False.
  • timeidx (int, optional) – The time index to use when extracting auxiallary variables used in the interpolation. This value must be set to match the same value used when the field variable was extracted. Default is 0.
  • method ({'cat', 'join'}, optional) – The aggregation method to use for sequences, either ‘cat’ or ‘join’. ‘cat’ combines the data along the Time index. ‘join’ is creates a new index for the file. This must be set to the same method used when extracting the field variable. The default is ‘cat’.
  • squeeze ({True, False}, optional) – Set to False to prevent dimensions with a size of 1 from being removed from the shape of the output. Default is True.
  • cache (dict, optional) – A dictionary of (varname, ndarray) which can be used to supply pre-extracted NetCDF variables to the computational routines. This can be used to prevent the repeated variable extraction from large sequences of data files. Default is None.
  • meta ({True, False}, optional) – Set to False to disable metadata and return numpy.ndarray instead of xarray.DataArray. Default is True.
Returns:

Returns the interpolated variable. If xarray is enabled and the meta parameter is True, then the result will be an xarray.DataArray object. Otherwise, the result will be a numpy.ndarray object with no metadata.

Return type:

xarray.DataArray or numpy.ndarray

Indices and tables