Fieldset object
- class Fieldset
Metview’s Fieldset object represents GRIB data. It is a container-like object with each entry representing a GRIB message.
Construction
Fieldsets can be directly constructed either as empty, with a path to a GRIB file or using read():
import metview as mv # empty Fieldset f1 = mv.Fieldset() # create from GRIB file f2 = mv.read("test.grib") # create from GRIB file f3 = mv.Fieldset(path="test.grib") # create from a set of GRIB files using globbing # (new in metview-python version 1.8.0) f4 = mv.Fieldset(path="test_*.grib") # create from multiple GRIB files f5 = mv.Fieldset(path=["a.grib", "b.grib", "/a/b/c.grib"])
Concatenation
Concatenation can be performed in these ways:
f4 = mv.Fieldset(fields=[f1, f2, f3]) # create from list of Fieldsets f4.append(f2) # append f2 onto the end of f4 f5 = mv.merge(f2, f3)The ‘list of files’ method of constructing a Fieldset shown in the previous section effectively creates a concatenation of those files into a Fieldset.
Indexing
Indexing and slicing works in the standard Python way. There is no such thing as a single field object in Metview, only a Fieldset with a single field, so all of the following return a Fieldset:
f[0] # first field f[1] # second field f[-1] # last field f[0:6] # the first 5 fields f[::2] # every second field my_fields = fs[np.array([1, 2, 0, 5])] # numpy array of indicesIt is also possible to assign fields into given locations in a Fieldset, for example:
grib = mv.read("t_for_xs.grib") grib[0] = grib[0] * 10 grib[4] = mv.sqrt(grib[3])Slicing is done with standard Python notation, e.g.
# select fields 4 to 7, step 2 my_slice = data[4:8:2]For more examples of indexing and slicing Fieldsets, see Slicing and dicing GRIB data.
Iteration
A Fieldset is iterable, with each iteration returning a single-field Fieldset, e.g.
fs = mv.Fieldset(path="test.grib") field_maxes = [f.maxvalue() for f in fs]len(fs) and fs.count() both return the number of fields in the Fieldset.
Inspection
The contents of a Fieldset can be easily inspected using the
ls()
anddescribe()
methods. See the Inspecting GRIB data notebook for some examples.
Filtering
A set of fields from a Fieldset can be extracted using the
read()
andselect()
functions. See the Filtering GRIB data notebook for the comparison of these two methods.For simple data extractions with
select()
a shorthand notation with the [] operator is also available. E.g. instead ofg = fs.select(shortName="t", level=500, typeOfLevel="isobaricInhPa")we can say:
g = fs["t500hPa"]See more examples here.
Note
select()
and its [] operator are only available from metview-python version 1.8.0.
Functions vs methods
Functions that work with Fieldsets can also be used as methods, provided that their first argument is a Fieldset. For example, the following two operations are shown in two equivalent ways:
a = mv.abs(fs) a = fs.abs() b = mv.bitmap(fs, 0) b = fs.bitmap(0)
Per-point methods
Unary functions and methods on Fieldsets act on each grid point of each field. For example, the
abs()
method will return a new Fieldset where all the grid values of all the fields have the absolute of their original value.Operations between Fieldsets act on corresponding grid points in the corresponding fields in each Fieldset. Both Fieldsets must have the same number of fields and the same number of points in their corresponding fields. For example, if we have one Fieldset containing analysis data for 99 vertical levels, and another with forecast data for the same 99 levels (stored in the same order) then we can easily compute the difference Fieldset like this:
diff = forecast_fs - analysis_fs # contains 99 fields of differencesSimilarly, operations work between Fieldsets and numbers, for example:
temperature_in_K = mv.read("temp.grib") temperature_in_C = temperature_in_K - 273.15The following list of operators are valid when acting between two Fieldsets and also when acting between a Fieldset and a number. Of these, the logical operators return a Fieldset containing values of 1 where they pass the test, or 0 where they fail.
Operator |
Decription |
---|---|
+ |
addition |
- |
subtraction |
* |
multiplication |
/ |
division |
** |
power |
& |
logical and (result is 1s and 0s) |
| |
logical or (result is 1s and 0s) |
~ |
logical not (result is 1s and 0s) |
> |
greater than (result is 1s and 0s) |
>= |
greater than or equal to (result is 1s and 0s) |
< |
less than (result is 1s and 0s) |
<= |
less than or equal to (result is 1s and 0s) |
== |
equal (result is 1s and 0s) |
!= |
not equal (result is 1s and 0s) |
Data extraction
A Fieldset can return an xarray by calling its to_dataset()
method.
A Fieldset can return a numpy array of values by calling its values()
method.
A Fieldset can return a Geopoints
object by calling the grib_to_geo()
function.
Methods and functions
List of methods
Computes the absolute value |
|
Computes the absolute vorticity from a relative vorticity |
|
Adds up the values per field in a |
|
Computes the arc cosine |
|
Computes the arc sine |
|
Computes the arc tangent |
|
Computes the arc tangent of 2 |
|
Computes the average per field in a |
|
Computes the zonal averages for each field in a |
|
Computes the meridional averages for each field in a |
|
Returns the base date(s) of a given |
|
Computes the bearings with respect to a reference in a |
|
Converts numbers to missing values in a |
|
Performs spatial convolution on a |
|
Computes the area-weighted correlation for each field in a |
|
Computes the cosine |
|
Generates a field with the cosine of the latitudes in a |
|
Returns the covariance of two |
|
Computes the area-weighted covariance for each field in a |
|
Returns information on missing values in a |
|
Converts an xndarray dataset to a |
|
De-accumulate values in a |
|
Prints summary of the contents of a |
|
Computes the wind direction |
|
Computes the distances in a |
|
Computes the integer part of a divison |
|
Computes the horizontal divergence of a vector |
|
Computes the exponential |
|
|
Fills missing values along the horizontal line |
Find locations of values in a |
|
Computes first West-East derivative of a |
|
Computes first South-North derivative of a |
|
Converts int GRIB to float GRIB |
|
Computes the frequencies of a |
|
Computes the geostrophic wind on pressure levels in a |
|
Finds values in field and returns the result as |
|
Computes horizontal gradient of a |
|
Reads GRIB headers using ecCodes keys |
|
Reads GRIB headers using ecCodes keys |
|
Reads GRIB headers using ecCodes keys |
|
Reads GRIB headers using ecCodes keys |
|
Reads GRIB headers using ecCodes keys |
|
Reads GRIB headers using ecCodes keys |
|
Returns list of indexes into the fields of the Fieldset |
|
Writes GRIB headers using ecCodes keys |
|
Writes GRIB headers using ecCodes keys |
|
Writes GRIB headers using ecCodes keys |
|
Writes GRIB headers using ecCodes keys |
|
Writes GRIB headers using ecCodes keys |
|
Sets GRIB packing bit width |
|
Computes the grid cell area in a |
|
Builds a |
|
Integer part |
|
Converts float GRIB to int GRIB |
|
Computes the surface integral of a |
|
Computes the average weighted by the gridcell area for each field in |
|
Interpolates |
|
Computes the horizontal Laplacian of |
|
Computes the natural logarithm |
|
Computes the base 10 logarithm |
|
Builds an output |
|
Dumps metadata for each message in a |
|
Maximum |
|
Maximum value of a |
|
Generates a |
|
Minimum |
|
Interpolates a model level |
|
Computes the integer remainder of a divison |
|
Computes the geopotential on model levels for a |
|
Interpolates a |
|
Returns the nearest grid point value from a |
|
Returns the nearest grid point value from a |
|
Converts missing values to numbers in a |
|
Interpolates pressure level fields to pressure levels |
|
Generates polygon masks for a |
|
Computes the pressure on model levels in a |
|
Computes the verical pressure derivative |
|
Generates masks based on a radius around a point for |
|
Computes he pointwise root mean square for |
|
Computes the area-weighted root mean square for |
|
Computes the second West-East derivative of a |
|
Computes the econd South-North derivative of a |
|
Filters |
|
Computes the sign |
|
Compute the shear deformation of a vector |
|
Computes the sine |
|
Generates a field with the cosine of the latitudes in a |
|
Performs spatial smoothing on a |
|
|
Performs spatial smoothing on a |
Computes the solar zenith angle for a |
|
Sorts a |
|
Computes the square root |
|
Returns the standard deviation of all the fields in a |
|
Computes the area-weighted standard deviation for each field in a |
|
Compute the stretch deformation of a vector |
|
Returns the indexes of the four surrounding grid points in a |
|
Computes the tangent |
|
Generates a field with the tangent of the latitudes in a |
|
Computes the pressure thickness on model levels in a |
|
Converts a |
|
Computes the pressure on model levels in a |
|
Computes the pressure thickness of model levels in a |
|
Performs a vertical integration for a |
|
Returns the valid date(s) of a given |
|
Returns the values from a data object |
|
Returns the variance of all the fields in a |
|
Computes the area-weighted variance for each field in a |
|
Performs a vertical integration for a |
|
Computes the relative vorticity of a vector |