API reference

Climatology

class pyfesom.climatology.climatology(path, climname='woa05')

Bases: object

Class that contains information from ocean 1 degree annual climatology. Presently options are WOA2005 and PHC3.0

Parameters:

path : str

Path to the directory with climatology files

climname : str

Name of the climatology (‘woa05’ or ‘phc’)

Returns:

object with climatology fields

Attributes

x (2d array) longitudes
y (2d array) latitudes
T (3d array) temperatures
S (3d array) salinity
z (1d array) depths
Tyz (2d array) zonal mean of temperature
Syz (3d array) zonal mean of salinity

Load mesh

pyfesom.load_mesh_data.cut_region(mesh, box=[13, 30, 53, 66], depth=0)

Cut region from the mesh.

Parameters:

mesh : object

FESOM mesh object

box : list

Coordinates of the box in [-180 180 -90 90] format. Default set to [13, 30, 53, 66], Baltic Sea.

depth : float

depth

Returns:

elem_no_nan : array

elements that belong to the region defined by box.

no_nan_triangles : array

boolian array of size elem2d with True for elements that belong to the region defines by box.

pyfesom.load_mesh_data.fesom2depth(depth, data3, mesh, verbose=True)

Return 2d array of the 2d mesh shape with data from the model level closest to the desired depth. There is no interpolation to the depth.

Parameters:

depth : int

desired depth

data3 : array

complete 3d data (vector) for one timestep

mesh : fesom_mesh object

mesh representation

verbose : bool

flag to turn off information about which level will be used.

Returns:

data2 : array

2d array (actually vector) with data from the desired level.

class pyfesom.load_mesh_data.fesom_mesh(path, abg=[50, 15, -90], get3d=True)

Bases: object

Creates instance of the FESOM mesh.

This class creates instance that contain information about FESOM mesh. At present the class works with ASCII representation of the FESOM grid, but should be extended to be able to read also netCDF version (probably UGRID convention).

Minimum requirement is to provide the path to the directory, where following files should be located (not nessesarely all of them will be used):

  • nod2d.out
  • nod3d.out
  • elem2d.out
  • aux3d.out
Parameters:

path : str

Path to the directory with mesh files

abg : list

alpha, beta and gamma Euler angles. Default [50, 15, -90]

get3d : bool

do we load complete 3d mesh or only 2d nodes.

Returns:

mesh : object

fesom_mesh object

Attributes

path (str) Path to the directory with mesh files
x2 (array) x position (lon) of the surface node
y2 (array) y position (lat) of the surface node
n2d (int) number of 2d nodes
e2d (int) number of 2d elements (triangles)
n3d (int) number of 3d nodes
nlev (int) number of vertical levels
zlevs (array) array of vertical level depths
voltri (array) array with 2d volume of triangles
alpha (float) Euler angle alpha
beta (float) Euler angle beta
gamma (float) Euler angle gamma
read2d()

Reads only surface part of the mesh. Useful if your mesh is large and you want to visualize only surface.

read3d()

Reads 3d part of the mesh.

pyfesom.load_mesh_data.get_data(data, mesh, depth=0, verbose=True)

Show data from the model level that is closest to the desired depth.

Parameters:

data : array

complete 3d data for one timestep

mesh : fesom_mesh object

mesh representation

depth : int

desired depth

verbose : bool

flag to turn off information about which level will be used.

Returns:

level_data : array

2d array (actually vector) with data from model level that is closest to the desired depth.

elem_no_nan : array

array with triangles (defined as triplets of node indexes) with not NaN elements.

pyfesom.load_mesh_data.get_layer_mean(data, depth, mesh, timeslice=None)

Return mean over the model depth that is closest to specified depth.

pyfesom.load_mesh_data.ind_for_depth(depth, mesh)

Return indexes that belong to certain depth.

Parameters:

depth : float

desired depth. Note thre will be no interpolation, the model level that is closest to desired depth will be selected.

mesh : object

FESOM mesh object

Returns:

ind_depth : 1d array

vector with the size equal to the size of the surface nodes with index values where we have data values and missing values where we don’t have data values.

ind_noempty : 1d array

vector with indexes of the ind_depth that have data values.

ind_empty : 1d array

vector with indexes of the ind_depth that do not have data values.

pyfesom.load_mesh_data.load_mesh(path, abg=[50, 15, -90], get3d=True, usepickle=True, usejoblib=False)

Loads FESOM mesh

Parameters:

path : str

Path to the directory with mesh files

abg : list

alpha, beta and gamma Euler angles. Default [50, 15, -90]

get3d : bool

do we load complete 3d mesh or only 2d nodes.

Returns:

mesh : object

fesom_mesh object

pyfesom.load_mesh_data.read_fesom_2d(str_id, months, years, mesh, result_path, runid, ext, how='mean')

Note

Deprecated in pyfesom 0.1 read_fesom_2d will be removed in future, it is replaced by get_data.

pyfesom.load_mesh_data.read_fesom_3d(str_id, months, years, mesh, result_path, runid, ext, how='mean')

Note

Deprecated in pyfesom 0.1 read_fesom_3d will be removed in future, it is replaced by get_data.

pyfesom.load_mesh_data.read_fesom_mesh(path, alpha, beta, gamma, read_diag=True)

Note

Deprecated in pyfesom 0.1 read_fesom_mesh will be removed in future, it is replaced by load_mesh.

Regriding

pyfesom.regriding.clim2regular(climatology, param, olons, olats, levels=None, verbose=True, radius_of_influence=100000)

Interpolation of data on the regular grid to climatology for the set of levels.

Parameters:

climatology: climatology object

FESOM climatology object

param : str

name of the parameter to interpolate. Only ‘T’ for temperature and ‘S’ for salinity are currently supported.

olons : array

1d or 2d array of longitudes to interpolate climatology on.

olats : array

1d or 2d array of longitudes to interpolate climatology on.

levels : list like

list of levels to interpolate to.

Returns:

xx : 2d array

longitudes

yy : 2d array

latitudes

out_data : 2d array

array with climatology data interpolated to desired levels

pyfesom.regriding.create_indexes_and_distances(mesh, lons, lats, k=1, n_jobs=2)

Creates KDTree object and query it for indexes of points in FESOM mesh that are close to the points of the target grid. Also return distances of the original points to target points.

Parameters:

mesh : fesom_mesh object

pyfesom mesh representation

lons/lats : array

2d arrays with target grid values.

k : int

k-th nearest neighbors to return.

n_jobs : int, optional

Number of jobs to schedule for parallel processing. If -1 is given all processors are used. Default: 1.

Returns:

distances : array of floats

The distances to the nearest neighbors.

inds : ndarray of ints

The locations of the neighbors in data.

pyfesom.regriding.fesom2clim(data, mesh, climatology, levels=None, verbose=True, how='nn', k_neighbors=10, radius_of_influence=100000)

Interpolation of fesom data to grid of the climatology for set of levels.

Parameters:

data : array

1d array of FESOM 3d data for one time step

mesh : mesh object

FESOM mesh object

climatology: climatology object

FESOM climatology object

levels : list like

list of levels to interpolate. At present you can use only standard levels of WOA. If levels are not specified, all standard WOA levels will be used.

how : str

Interpolation method. Options are ‘nn’ (nearest neighbor) and ‘idist’ (inverce distance)

k : int

k-th nearest neighbors to use. Only used when how==idist

radius_of_influence : int

Cut off distance in meters.

Returns:

xx : 2d array

longitudes

yy : 2d array

latitudes

out_data : 2d array

array with data interpolated to climatology level

pyfesom.regriding.fesom2regular(data, mesh, lons, lats, distances=None, inds=None, how='nn', k=10, radius_of_influence=100000, n_jobs=2)

Interpolates data from FESOM mesh to target (usually regular) mesh.

Parameters:

data : array

1d array that represents FESOM data at one level.

mesh : fesom_mesh object

pyfesom mesh representation

lons/lats : array

2d arrays with target grid values.

distances : array of floats, optional

The distances to the nearest neighbors.

inds : ndarray of ints, optional

The locations of the neighbors in data.

how : str

Interpolation method. Options are ‘nn’ (nearest neighbor) and ‘idist’ (inverce distance)

k : int

k-th nearest neighbors to use. Only used when how==idist

radius_of_influence : int

Cut off distance in meters.

n_jobs : int, optional

Number of jobs to schedule for parallel processing. If -1 is given all processors are used. Default: 1.

Returns:

data_interpolated : 2d array

array with data interpolated to the target grid.

pyfesom.regriding.lon_lat_to_cartesian(lon, lat, R=6371000)

Calculates lon, lat coordinates of a point on a sphere with radius R. Taken from http://earthpy.org/interpolation_between_grids_with_ckdtree.html

Parameters:

lon : 1d array

longitudes

lat : 1d array

latitudes

R : float

radius of the sphere

Returns:

x,y,z : 1d arrays

cartesian coordinates

pyfesom.regriding.regular2clim(data, ilons, ilats, izlevs, climatology, levels=None, verbose=True)

Interpolation of data on the regular grid to climatology for the set of levels.

Parameters:

data : array

3d array of data on regular grid for one time step

ilons : array

1d or 2d array of longitudes

ilats : array

1d or 2d array of latitudes

izlevs : array

depths of the source data

climatology: climatology object

FESOM climatology object

levels : list like

list of levels to interpolate. At present you can use only standard levels of WOA. If levels are not specified, all standard WOA levels will be used.

Returns:

xx : 2d array

longitudes

yy : 2d array

latitudes

out_data : 2d array

array with data interpolated to climatology level

pyfesom.regriding.regular2regular(data, ilons, ilats, olons, olats, distances=None, inds=None, how='nn', k=10, radius_of_influence=100000, n_jobs=2)

Interpolates from regular to regular grid. It’s a wraper around fesom2regular that creates an object that mimic fesom mesh class and contain only coordinates and flatten the data.

Parameters:

data : array

1d or 2d array that represents gridded data at one level.

ilons : arrea

mesh : fesom_mesh object

pyfesom mesh representation

lons/lats : array

2d arrays with target grid values.

distances : array of floats, optional

The distances to the nearest neighbors.

inds : ndarray of ints, optional

The locations of the neighbors in data.

how : str

Interpolation method. Options are ‘nn’ (nearest neighbor) and ‘idist’ (inverce distance)

k : int

k-th nearest neighbors to use. Only used when how==idist

radius_of_influence : int

Cut off distance in meters.

n_jobs : int, optional

Number of jobs to schedule for parallel processing. If -1 is given all processors are used. Default: 1.

Returns:

data_interpolated : 2d array

array with data interpolated to the target grid.

Utilites

pyfesom.ut.scalar_g2r(al, be, ga, lon, lat)

Converts geographical coordinates to rotated coordinates.

Parameters:

al : float

alpha Euler angle

be : float

beta Euler angle

ga : float

gamma Euler angle

lon : array

1d array of longitudes in geographical coordinates

lat : array

1d array of latitudes in geographical coordinates

Returns:

rlon : array

1d array of longitudes in rotated coordinates

rlat : array

1d araay of latitudes in rotated coordinates

pyfesom.ut.scalar_r2g(al, be, ga, rlon, rlat)

Converts rotated coordinates to geographical coordinates.

Parameters:

al : float

alpha Euler angle

be : float

beta Euler angle

ga : float

gamma Euler angle

rlon : array

1d array of longitudes in rotated coordinates

rlat : array

1d araay of latitudes in rotated coordinates

Returns:

lon : array

1d array of longitudes in geographical coordinates

lat : array

1d array of latitudes in geographical coordinates

pyfesom.ut.vec_rotate_r2g(al, be, ga, lon, lat, urot, vrot, flag)

Rotate vectors from rotated coordinates to geographical coordinates.

Parameters:

al : float

alpha Euler angle

be : float

beta Euler angle

ga : float

gamma Euler angle

lon : array

1d array of longitudes in rotated or geographical coordinates (see flag parameter)

lat : array

1d array of latitudes in rotated or geographical coordinates (see flag parameter)

urot : array

1d array of u component of the vector in rotated coordinates

vrot : array

1d array of v component of the vector in rotated coordinates

flag : 1 or 0

flag=1 - lon,lat are in geographical coordinate flag=0 - lon,lat are in rotated coordinate

Returns:

u : array

1d array of u component of the vector in geographical coordinates

v : array

1d array of v component of the vector in geographical coordinates

Module contents