epygram.fields.D3Field — 3-D Field class

Contains the class that handle a Horizontal 2D field.


class epygram.fields.D3Field.D3CommonField(*args, **kwargs)[source]

Bases: epygram.base.Field

3-Dimensions common field class.

Note

This class is managed by footprint.

  • info: Not documented
  • priority: PriorityLevel::DEFAULT (rank=1)

Automatic parameters from the footprint:

  • comment (str) - rwd - Not documented, sorry.
    • Optional. Default is None.
  • fid (footprints.stdtypes.FPDict) - rwx - Not documented, sorry.
  • structure (str) - rxx - Type of Field geometry.
    • Values: set([‘3D’])
  • units (str) - rwd - Not documented, sorry.
    • Optional. Default is ‘’.

Constructor. See its footprint for arguments.

as_dicts(subzone=None)[source]

Export values as a list of dicts.

Parameters:subzone – defines the LAM subzone to be included, in LAM case, among: ‘C’, ‘CI’.
as_lists(order=u'C', subzone=None)[source]

Export values as a dict of lists (in fact numpy arrays).

Parameters:
  • order – whether to flatten arrays in ‘C’ (row-major) or ‘F’ (Fortran, column-major) order.
  • subzone – defines the LAM subzone to be included, in LAM case, among: ‘C’, ‘CI’.
as_points(subzone=None)[source]

Export values as a fieldset of points.

Parameters:subzone – defines the LAM subzone to be included, in LAM case, among: ‘C’, ‘CI’.
as_profiles(subzone=None)[source]

Export values as a fieldset of profiles.

Parameters:subzone – defines the LAM subzone to be included, in LAM case, among: ‘C’, ‘CI’.
dctspectrum(subzone=None, level_index=None, validity_index=None)[source]

Returns the DCT spectrum of the field, as a epygram.spectra.Spectrum instance.

Parameters:
  • subzone – LAM zone on which to compute the DCT, among (None, ‘C’, ‘CI’, ‘CIE’). Defaults to ‘C’.
  • level_index – the level index on which to compute the DCT
  • validity_index – the validity index on which to compute the DCT
decumulate(center=False)[source]

Decumulate cumulated fields (for a field with temporal dimension !).

Parameters:center
  • if False, values at t are mean values between t-1 and t, except at t=0 where it remains untouched.
  • if True, values at t are then (v[t-1, t] + v[t,t+1])/2., except for t=0 where it becomes v[0, 1] and t=last where it remains the same.
dump_to_nc(filename, variablename=None, fidkey=None)[source]

Dumps the field in a netCDF file.

Parameters:
  • filename – filename to dump in
  • variablename – variable name in netCDF
  • fidkey – forces variablename to self.fid[fidkey]; defaults to ‘netCDF’
extend(another_field_with_time_dimension)[source]

Extend the field with regard to time dimension with the field given as argument.

Be careful no check is done for consistency between the two fields geometry (except that dimensions match) nor their validities.

extract_subarray(first_i, last_i, first_j, last_j)[source]

Extract a rectangular sub-array from the field, given the i,j index limits of the sub-array, and return the extracted field.

extract_subdomain(geometry, interpolation=u'nearest', external_distance=None, exclude_extralevels=True, getdata=True)[source]

Extracts a subdomain from a field, given a new geometry.

Parameters:
  • geometry – defines the geometry on which extract data
  • interpolation
    defines the interpolation function used to compute
    the profile at requested lon/lat from the fields grid:
    • if ‘nearest’ (default), extracts profile at the horizontal nearest neighboring gridpoint;
    • if ‘linear’, computes profile with horizontal linear interpolation;
    • if ‘cubic’, computes profile with horizontal cubic interpolation.
  • external_distance – can be a dict containing the target point value and an external field on the same grid as self, to which the distance is computed within the 4 horizontally nearest points; e.g. {‘target_value’:4810, ‘external_field’:an_H2DField_with_same_geometry}. If so, the nearest point is selected with distance = abs(target_value - external_field.data)
  • exclude_extralevels – if True levels with no physical meaning are suppressed.
  • getdata – if False returns a field without data
extract_zoom(zoom, extra_10th=False)[source]

Extract an unstructured field with the gridpoints contained in zoom.

Parameters:
  • zoom – a dict(lonmin=, lonmax=, latmin=, latmax=).
  • extra_10th – if True, add 1/10th of the X/Y extension of the zoom (regular_lonlat grid case only).
getvalue_ll(lon=None, lat=None, level=None, validity=None, interpolation=u'nearest', neighborinfo=False, one=True, external_distance=None)[source]

Returns the value of the field on point of coordinates (lon, lat, level). lon and lat may be list/numpy.array.

Parameters:
  • interpolation

    if:

    • ‘nearest’ (default), returns the value of the nearest neighboring gridpoint;
    • ‘linear’, computes and returns the field value with linear spline interpolation;
    • ‘cubic’, computes and returns the field value with cubic spline interpolation.
  • level – is the True level not the index of the level. Depending on the vertical coordinate, it could be expressed in Pa, m.
  • validity – is a FieldValidity or a FieldValidityList instance
  • neighborinfo – if set to True, returns a tuple (value, (lon, lat)), with (lon, lat) being the actual coordinates of the neighboring gridpoint (only for interpolation == ‘nearest’).
  • one – if False and len(lon) is 1, returns [value] instead of value.
  • external_distance – can be a dict containing the target point value and an external field on the same grid as self, to which the distance is computed within the 4 horizontally nearest points; e.g. {‘target_value’:4810, ‘external_field’:an_H2DField_with_same_geometry}. If so, the nearest point is selected with distance = abs(target_value - external_field.data)

Warning: for interpolation on Gauss geometries, requires the pyproj module.

global_shift_center(longitude_shift)[source]

Shifts the center of the geometry (and the data accordingly) by longitude_shift (in degrees).

Parameters:longitude_shift – has to be a multiple of the grid’s resolution in longitude.

For global RegLLGeometry grids only.

histogram(subzone=None, over=(None, None), title=None, get_n_bins_patches=False, together_with=[], mask_threshold=None, center_hist_on_0=False, minmax_in_title=True, figsize=None, **hist_kwargs)[source]

Build an histogram of the field data.

Parameters:
  • subzone – LAM zone on which to compute the DCT, among (None, ‘C’, ‘CI’, ‘CIE’).
  • over – any existing figure and/or ax to be used for the plot, given as a tuple (fig, ax), with None for missing objects. fig is the frame of the matplotlib figure, containing eventually several subplots (axes); ax is the matplotlib axes on which the drawing is done. When given (is not None), these objects must be coherent, i.e. ax being one of the fig axes.
  • title – title for the plot
  • get_n_bins_patches – if True, returns n, bins, patches (cf. matplotlib hist()) in addition
  • together_with – another field or list of fields to plot on the same histogram
  • mask_threshold – dict with min and/or max value(s) to mask outside.
  • center_hist_on_0 – to center the histogram on 0.
  • minmax_in_title – if True and range is not None, adds min and max values in title.
  • hist_kwargs – any keyword argument to be passed to matplotlib’s hist()
  • figsize – figure sizes in inches, e.g. (5, 8.5). Default figsize is config.plotsizes.
max(subzone=None)[source]

Returns the maximum value of data.

mean(subzone=None)[source]

Returns the mean value of data.

min(subzone=None)[source]

Returns the minimum value of data.

nonzero(subzone=None)[source]

Returns the number of non-zero values (whose absolute value > config.epsilon).

quadmean(subzone=None)[source]

Returns the quadratic mean of data.

resample(target_geometry, weighting=u'nearest', subzone=None, neighbour_info=None, radius_of_influence=None, neighbours=8, epsilon=0, reduce_data=True, nprocs=1, segments=None, fill_value=None, sigma=None, with_uncert=False)[source]

Resample data (interpolate) on a target geometry, using pyresample.

Parameters:
  • weighting

    among [‘nearest’, ‘gauss’ (in which case sigma can be specified,

    cf. sigma option below),

    lambda r:1/r**2 (or other function, r in m)]

  • subzone – restrain source field to LAMzone in LAM geometry case
  • neighbour_info
    • if True: skips the resampling, only compute
      neighbours and return this information as a tuple;
    • if given as a tuple, skips the computation of neighbours, taking those information given in argument (from a former call with neighbour_info=True)

Options from pyresample, neighbouring:

Parameters:
  • radius_of_influence – Cut off distance in meters, default value computed from geometry
  • neighbours – The number of neigbours to consider for each grid point
  • epsilon – Allowed uncertainty in meters. Increasing uncertainty reduces execution time
  • reduce_data – Perform initial coarse reduction of source dataset in order to reduce execution time
  • nprocs – Number of processor cores to be used
  • segments – Number of segments to use when resampling. If set to None an estimate will be calculated

Options from pyresample, resampling:

Parameters:
  • fill_value – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked
  • sigma – sigma (in m) to use for the gauss weighting w = exp(-dist^2/sigma^2) Default is twice the resolution.
  • with_uncert – Calculate uncertainty estimates (returned as side fields)
resample_on_regularll(borders, resolution_in_degrees, **kwargs)[source]

Resample data (interpolate) on a target geometry defined as a regular lonlat grid of given resolution. Wrapper to resample() method, building target geometry.

Parameters:
  • borders (dict) – dict(lonmin=, latmin=, lonmax=, latmax=)
  • resolution_in_degrees – resolution of the lonlat grid in degrees. Supposed to be a divisor of lonmax-lonmin and latmax-latmin.

Other kwargs to be passed to resample() method.

shave(minval=None, maxval=None)[source]

Shave data: cut values greater (resp. lower) than maxval to maxval (resp. minval).

Parameters:
  • maxval – upper threshold
  • minval – lower threshold
spectral

Returns True if the field is spectral.

stats(subzone=None)[source]

Computes some basic statistics on the field, as a dict containing: {‘min’, ‘max’, ‘mean’, ‘std’, ‘quadmean’, ‘nonzero’}. See each of these methods for details.

Parameters:subzone – optional, among (‘C’, ‘CI’), for LAM fields only, plots the data resp. on the C or C+I zone. Default is no subzone, i.e. the whole field.
std(subzone=None)[source]

Returns the standard deviation of data.

time_reduce(reduce_function=u'mean')[source]

Make a time reduction of field, i.e. apply the requested reduce_function on the time dimension of the field. Replace data (and validity) in place.

Parameters:reduce_function – must be a numpy.ndarray method, e.g. ‘mean’, ‘sum’, ‘min’, ‘max’, ...
time_smooth(length, window=u'center')[source]

Smooth data in time.

Parameters:
  • length – time span to average on
  • window

    replace data[t] by:

    • data[t-length/2:t+length/2].mean() if ‘center’;
    • data[t-length:t].mean() if ‘left’;
    • data[t:t+length/2].mean() if ‘right’.
what(out=<open file '<stdout>', mode 'w'>, validity=True, vertical_geometry=True, cumulativeduration=True, arpifs_var_names=False, fid=True)[source]

Writes in file a summary of the field.

Parameters:
  • out – the output open file-like object.
  • vertical_geometry – if True, writes the validity of the field.
  • vertical_geometry – if True, writes the vertical geometry of the field.
  • cumulativeduration – if False, not written.
  • arpifs_var_names – if True, prints the equivalent ‘arpifs’ variable names.
  • fid – if True, prints the fid.

class epygram.fields.D3Field.D3Field(*args, **kwargs)[source]

Bases: epygram.fields.D3Field.D3CommonField

3-Dimensions field class. A field is defined by its identifier ‘fid’, its data, its geometry (gridpoint and optionally spectral), and its validity.

The natural being of a field is gridpoint, so that: a field always has a gridpoint geometry, but it has a spectral geometry only in case it is spectral.

Note

This class is managed by footprint.

  • info: Not documented
  • priority: PriorityLevel::DEFAULT (rank=1)

Automatic parameters from the footprint:

  • comment (str) - rwd - Not documented, sorry.
    • Optional. Default is None.
  • fid (footprints.stdtypes.FPDict) - rwx - Not documented, sorry.
  • geometry (epygram.geometries.D3Geometry.D3Geometry) - rxx - Geometry defining the position of the field gridpoints.
  • processtype (str) - rxx - Generating process.
    • Optional. Default is None.
  • spectral_geometry (epygram.geometries.SpectralGeometry.SpectralGeometry) - rxx - For a spectral field, its spectral geometry handles spectral transforms and dimensions.
    • Optional. Default is None.
  • structure (str) - rxx - Type of Field geometry.
    • Values: set([‘3D’])
  • units (str) - rwd - Not documented, sorry.
    • Optional. Default is ‘’.
  • validity (epygram.base.FieldValidityList) - rwx - Validity of the field.
    • Optional. Default is FieldValidityList::<<as_list:: [FieldValidity::<epygram.base.FieldValidity object at 0x7ffa493aa990>]>>.

Constructor. See its footprint for arguments.

compute_xy_spderivatives()[source]

Compute the derivatives of field in spectral space, then come back in gridpoint space. Returns the two derivative fields.

The spectral transform and derivatives subroutines are actually included in the spectral geometry’s compute_xy_spderivatives() method.

data

Accessor to the field data.

getdata(subzone=None, d4=False)[source]

Returns the field data, with 3D shape if the field is not spectral, 2D if spectral.

Parameters:
  • subzone – optional, among (‘C’, ‘CI’), for LAM fields only, returns the data resp. on the C or C+I zone. Default is no subzone, i.e. the whole field.
  • d4
    • if True, returned values are shaped in a 4 dimensions array
    • if False, shape of returned values is determined with respect to geometry.

Shape of 4D data:

  • Rectangular grids:

    grid[t,k,0,0] is SW, grid[t,k,-1,-1] is NE

    grid[t,k,0,-1] is SE, grid[t,k,-1,0] is NW

    with k the level, t the temporal dimension

  • Gauss grids:

    grid[t,k,0,:Nj] is first (Northern) band of latitude, masked after Nj = number of longitudes for latitude j

    grid[t,k,-1,:Nj] is last (Southern) band of latitude (idem).

    with k the level, t the temporal dimension

getlevel(level=None, k=None)[source]

Returns a level of the field as a new field.

Parameters:
  • level – the requested level expressed in coordinate value (Pa, m...)
  • k – the index of the requested level
getvalidity(index_or_validity)[source]

Returns the field restrained to one of its temporal validity as a new field.

Parameters:index_or_validity – can be either a epygram.base.FieldValidity instance or the index of the requested validity in the field’s FieldValidityList.
getvalue_ij(i=None, j=None, k=None, t=None, one=True)[source]

Returns the value of the field on point of indices (i, j, k, t). Take care (i, j, k, t) is python-indexing, ranging from 0 to dimension-1.

Parameters:
  • i – the index in X direction
  • j – the index in Y direction
  • k – the index of the level (not a value in Pa or m...)
  • t – the index of the temporal dimension (not a validity object)

k and t can be scalar even if i and j are arrays.

If one is False, returns [value] instead of value.

gp2sp(spectral_geometry)[source]

Transforms the gridpoint field into spectral space, according to the spectral_geometry mandatorily passed as argument. Replaces data in place.

Parameters:spectral_geometry – instance of SpectralGeometry, actually containing spectral transform subroutine (in in its own gp2sp() method).
select_subzone(subzone)[source]

If a LAMzone defines the field, select only the subzone from it.

Parameters:subzone – among (‘C’, ‘CI’).

Warning: modifies the field and its geometry in place !

setdata(data)[source]

Sets field data, checking data to have the good shape according to geometry.

Parameters:data – dimensions should in any case be ordered in a subset of

(t,z,y,x), or (t,z,n) if spectral (2D spectral coefficients must be 1D with ad hoc ordering, n being the total number of spectral coefficients).

data may be 4D (3D if spectral) even if the field is not, as long as the above dimensions ordering is respected.

sp2gp()[source]

Transforms the spectral field into gridpoint, according to its spectral geometry. Replaces data in place.

The spectral transform subroutine is actually included in the spectral geometry’s sp2gp() method.


class epygram.fields.D3Field.D3VirtualField(*args, **kwargs)[source]

Bases: epygram.fields.D3Field.D3CommonField

3-Dimensions Virtual field class.

Data is taken from other fields, either: - a given fieldset - a resource in which are stored fields defined by resource_fids.

Note

This class is managed by footprint.

  • info: Not documented
  • priority: PriorityLevel::DEFAULT (rank=1)

Automatic parameters from the footprint:

  • comment (str) - rwd - Not documented, sorry.
    • Optional. Default is None.
  • fid (footprints.stdtypes.FPDict) - rwx - Not documented, sorry.
  • fieldset (epygram.base.FieldSet) - rxx - Set of real fields that can compose the Virtual Field.
    • Optional. Default is FieldSet::<<as_list:: []>>.
  • resource (epygram.base.Resource) - rxx - Resource in which is stored the fields defined by resource_fids.
    • Optional. Default is None.
  • resource_fids (footprints.stdtypes.FPList) - rxx - Definition of the fields in resource that compose the virtual field.
    • Optional. Default is FPList::<<as_list:: []>>.
  • structure (str) - rxx - Type of Field geometry.
    • Values: set([‘3D’])
  • units (str) - rwd - Not documented, sorry.
    • Optional. Default is ‘’.
data

Returns the field data, with 3D shape if the field is not spectral, 2D if spectral.

Parameters:
  • subzone – optional, among (‘C’, ‘CI’), for LAM fields only, returns the data resp. on the C or C+I zone. Default is no subzone, i.e. the whole field.
  • d4
    • if True, returned values are shaped in a 4 dimensions array
    • if False, shape of returned values is determined with respect to geometry

Shape of 3D data:

  • Rectangular grids:

    grid[k,0,0] is SW, grid[k,-1,-1] is NE

    grid[k,0,-1] is SE, grid[k,-1,0] is NW

    with k the level

  • Gauss grids:

    grid[k,0,:Nj] is first (Northern) band of latitude, masked after Nj = number of longitudes for latitude j

    grid[k,-1,:Nj] is last (Southern) band of latitude (idem).

    with k the level

deldata()[source]

deldata() not implemented on virtual fields.

getdata(subzone=None, d4=False)[source]

Returns the field data, with 3D shape if the field is not spectral, 2D if spectral.

Parameters:
  • subzone – optional, among (‘C’, ‘CI’), for LAM fields only, returns the data resp. on the C or C+I zone. Default is no subzone, i.e. the whole field.
  • d4
    • if True, returned values are shaped in a 4 dimensions array
    • if False, shape of returned values is determined with respect to geometry

Shape of 3D data:

  • Rectangular grids:

    grid[k,0,0] is SW, grid[k,-1,-1] is NE

    grid[k,0,-1] is SE, grid[k,-1,0] is NW

    with k the level

  • Gauss grids:

    grid[k,0,:Nj] is first (Northern) band of latitude, masked after Nj = number of longitudes for latitude j

    grid[k,-1,:Nj] is last (Southern) band of latitude (idem).

    with k the level

getlevel(level=None, k=None)[source]

Returns a level of the field as a new field.

Parameters:
  • level – the requested level expressed in coordinate value (Pa, m...)
  • k – the index of the requested level
getvalue_ij(i=None, j=None, k=None, t=None, one=True)[source]

Returns the value of the field on point of indices (i, j, k, t). Take care (i, j, k, t) is python-indexing, ranging from 0 to dimension-1.

Parameters:
  • i – the index in X direction
  • j – the index in Y direction
  • k – the index of the level (not a value in Pa or m...)
  • t – the index of the temporal dimension (not a validity object)

k and t can be scalar even if i and j are arrays.

If one is False, returns [value] instead of value.

gp2sp(spectral_geometry)[source]

Transforms the gridpoint field into spectral space, according to the spectral_geometry mandatorily passed as argument. Replaces data in place.

Parameters:spectral_geometry – instance of SpectralGeometry, actually containing spectral transform subroutine (in in its own gp2sp() method).
setdata(data)[source]

setdata() not implemented on virtual fields.

sp2gp()[source]

Transforms the spectral field into gridpoint, according to its spectral geometry. Replaces data in place.

The spectral transform subroutine is actually included in the spectral geometry’s sp2gp() method.