epygram.fields.H2DVectorField — Horizontal 2-D Vector Field class

Contains the class for a Horizontal 2D field.

Plus a function to create a Vector field from 2 scalar fields.


epygram.fields.H2DVectorField.make_vector_field(fX, fY)[source]

Creates a new epygram.H2DVectorField from two epygram.H2DField fX, fY representing resp. the X and Y components of the vector in the field geometry.


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

Bases: epygram.base.Field

Horizontal 2-Dimensions Vector field class.

This is a wrapper to a list of H2DField(s), representing the components of a vector projected on its geometry (the grid axes).

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.
  • components (footprints.stdtypes.FPList) - rxx - List of Fields that each compose a component of the vector.
    • Optional. Default is FPList::<<as_list:: []>>.
  • fid (footprints.stdtypes.FPDict) - rwx - Not documented, sorry.
  • processtype (str) - rxx - Generating process.
    • Optional. Default is None.
  • structure (str) - rxx - Type of Field geometry.
    • Values: set([‘H2D’])
  • 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 0x7ffa4939b210>]>>.
  • vector (bool) - rxx - Intrinsic vectorial nature of the field.
    • Values: set([True])

Constructor. See its footprint for arguments.

attach_components(*components)[source]

Attach components of the vector to the VectorField. components must be a series of H2DField.

compute_direction()[source]

Returns a epygram.H2DField whose data is the direction of the Vector field, in degrees.

compute_vordiv(divide_by_m=False)[source]

Compute vorticity and divergence fields from the vector field.

Parameters:divide_by_m – if True, apply f = f/m beforehand, where m is the map factor.
data

Accessor to the field data.

deldata()[source]

Empties the data.

extract_subarray(*args, **kwargs)[source]

Cf. D3Field.extract_subarray()

extract_subdomain(*args, **kwargs)[source]

Cf. D3Field.extract_subdomain()

extract_zoom(*args, **kwargs)[source]

Cf. D3Field.extract_zoom()

getdata(subzone=None, **kwargs)[source]

Returns the field data, with 2D shape if the field is not spectral, 1D if spectral, as a tuple with data for each component.

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.

Shape of 2D data: (x (0) being the X component, y (1) the Y one)

  • Rectangular grids:

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

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

  • Gauss grids:

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

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

getvalue_ij(*args, **kwargs)[source]

Returns the value of the different components of the field from indexes.

getvalue_ll(*args, **kwargs)[source]

Returns the value of the different components of the field from coordinates.

global_shift_center(longitude_shift)[source]

Shifts the center of the geometry (and the data accordingly) by longitude_shift (in degrees). longitude_shift has to be a multiple of the grid’s resolution in longitude.

For global RegLLGeometry grids only.

gp2sp(spectral_geometry=None)[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).
map_factorize(reverse=False)[source]

Multiply the field by its map factor.

Parameters:reverse – if True, divide.
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).

plotanimation(title=u'__auto__', repeat=False, interval=1000, **kwargs)[source]

Plot the field with animation with regards to time dimension. Returns a matplotlib.animation.FuncAnimation.

In addition to those specified below, all plotfield() method arguments can be provided.

Parameters:
  • title – title for the plot. ‘__auto__’ (default) will print the current validity of the time frame.
  • repeat – to repeat animation
  • interval – number of milliseconds between two validities
plotfield(over=(None, None), subzone=None, title=None, gisquality=u'i', specificproj=None, zoom=None, use_basemap=None, drawcoastlines=True, drawcountries=True, drawrivers=False, departments=False, boundariescolor=u'0.25', parallels=u'auto', meridians=u'auto', subsampling=1, symbol=u'barbs', symbol_options={u'color': u'k'}, plot_module=True, plot_module_options=None, bluemarble=0.0, background=False, quiverkey=None, quiver_options=None, components_are_projected_on=u'grid', map_factor_correction=True, mask_threshold=None, figsize=None, rcparams=None)[source]

Makes a simple plot of the field, with a number of options.

Parameters:
  • over – to plot the vectors over an existing figure (e.g. colorshades). 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.
  • subzone

    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.

  • gisquality – among (‘c’, ‘l’, ‘i’, ‘h’, ‘f’) – by increasing quality. Defines the quality for GIS elements (coastlines, countries boundaries...). Default is ‘i’. Cf. ‘basemap’ doc for more details.
  • specificproj

    enables to make basemap on the specified projection, among: ‘kav7’, ‘cyl’, ‘ortho’, (‘nsper’, {...}) (cf. Basemap doc).

    In ‘nsper’ case, the {} may contain:

    • ‘sat_height’ = satellite height in km;
    • ‘lon’ = longitude of nadir in degrees;
    • ‘lat’ = latitude of nadir in degrees.

    Overwritten by zoom.

  • zoom

    specifies the lon/lat borders of the map, implying hereby a ‘cyl’ projection. Must be a dict(lonmin=, lonmax=, latmin=, latmax=).

    Overwrites specificproj.

  • use_basemap – a basemap.Basemap object used to handle the projection of the map. If given, the map projection options (specificproj, zoom, gisquality ...) are ignored, keeping the properties of the use_basemap object. (because making Basemap is the most time-consuming step).
  • drawrivers – to add rivers on map.
  • departments – if True, adds the french departments on map (instead of countries).
  • boundariescolor – color of lines for boundaries (countries, departments, coastlines)
  • drawcoastlines – to add coast lines on map.
  • drawcountries – to add countries on map.
  • title – title for the plot. Default is field identifier.
  • meridians

    enable to fine-tune the choice of lines to plot, with either:

    • ‘auto’: automatic scaling to the basemap extents
    • ‘default’: range(0,360,10)
    • a list of values
    • a grid step, e.g. 5 to plot each 5 degree.
    • None: no one is plot
    • meridian == ‘greenwich’ // ‘datechange’ // ‘greenwich+datechange’ combination (,) will plot only these.
  • parallels

    enable to fine-tune the choice of lines to plot, with either:

    • ‘auto’: automatic scaling to the basemap extents
    • ‘default’: range(-90,90,10)
    • a list of values
    • a grid step, e.g. 5 to plot each 5 degree.
    • None: no one is plot
    • ‘equator’ // ‘polarcircles’ // ‘tropics’ or any combination (,) will plot only these.
  • subsampling – to subsample the number of gridpoints to plot. Ex: subsampling = 10 will only plot one gridpoint upon 10.
  • symbol – among (‘barbs’, ‘arrows’, ‘stream’)
  • symbol_options – a dict of options to be passed to barbs or quiver method.
  • plot_module – to plot module as colorshades behind vectors.
  • plot_module_options – options (dict) to be passed to module.plotfield().
  • bluemarble – if > 0.0 (and <=1.0), displays NASA’s “blue marble” as background. The numerical value sets its transparency.
  • background – if True, set a background color to continents and oceans.
  • quiverkey – to activate quiverkey; must contain arguments to be passed to pyplot.quiverkey(), as a dict.
  • components_are_projected_on – inform the plot on which axes the vector components are projected on (‘grid’ or ‘lonlat’).
  • map_factor_correction – if True, applies a correction of magnitude to vector due to map factor.
  • mask_threshold – dict with min and/or max value(s) to mask outside.
  • figsize – figure sizes in inches, e.g. (5, 8.5). Default figsize is config.plotsizes.
  • rcparams – list of (*args, **kwargs) to be passed to pyplot.rc() defaults to [((‘font’,), dict(family=’serif’)),]

This method uses (hence requires) ‘matplotlib’ and ‘basemap’ libraries.

quadmean(subzone=None)[source]

Returns the quadratic mean of data.

reproject_wind_on_lonlat(map_factor_correction=True, reverse=False)[source]

Reprojects a wind vector (u, v) from the grid axes onto real sphere, i.e. with components on true zonal/meridian axes.

Parameters:
  • map_factor_correction – if True, apply a correction of magnitude due to map factor.
  • reverse – if True, apply the reverse reprojection.
resample(*args, **kwargs)[source]

Cf. D3Field.resample()

resample_on_regularll(*args, **kwargs)[source]

Cf. D3Field.resample_on_regularll()

setdata(data)[source]

Sets data to its components.

Parameters:data – [data_i for i components]
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.

spectral

Returns True if the field is spectral.

std(subzone=None)[source]

Returns the standard deviation of data.

to_module()[source]

Returns a epygram.H2DField whose data is the module of the Vector field.

what(out=<open file '<stdout>', mode 'w'>, vertical_geometry=True, cumulativeduration=True)[source]

Writes in file a summary of the field.

Parameters:
  • out – the output open file-like object (duck-typing: out.write() only is needed).
  • vertical_geometry – if True, writes the vertical geometry of the field.