Contour

Docstrings for the Contour class

Objects

Contour()
Contour.get_contours()
Contour.plot_contour()
Contour.get_contour_segment()
Contour.process_contour()
Contour.gen_z_levels()
ContourF()
ContourF.calc_cross_contour_flow()
ContourF._update_cross_flow_vars()
ContourF._update_cross_flow_latlon()
ContourF._pressure_gradient_fpoint2()
ContourF.calc_geostrophic_flow()
ContourT()
ContourT.construct_pressure()
ContourT.calc_along_contour_flow()
ContourT.calc_along_contour_flow_2d()
ContourT._update_flow_vars()
ContourT._update_along_flow_latlon()

Contour classes

Contour()

class Contour():
None
Contour.get_contours()
@staticmethod
def Contour.get_contours(gridded, contour_depth):

A method to obtain the continuous isbobath contours within a supplied
gridded domain as a set of y indices and x indices for the model grid.

Parameters
----------
gridded : Coast
                The gridded object containing the dataset with the 'bathymetry' variable
contour_depth : int
                Depth of desired contours

Returns
-------
List of 2d ndarrays
                Each item of the list contains a different continuous isobath contour as a
                2d ndarray of indicies, i.e. for each list item:
                contour[:,0] contains the y indices for the contour on the model grid
                contour[:,1] contains the x indices for the contour on the model grid

Contour.plot_contour()
@staticmethod
def Contour.plot_contour(gridded, contour):

Quick plot method to plot a contour over a pcolormesh of the
model bathymetry

Parameters
----------
gridded : Coast
                The gridded object containing the dataset with the 'bathymetry' variable
contour : 2d array
                contour[:,0] contains the y indices for the contour on the model grid
                contour[:,1] contains the x indices for the contour on the model grid
                i.e. contour = np.vstack((y_indices,x_indices)).T

Returns
-------
None

Contour.get_contour_segment()
@staticmethod
def Contour.get_contour_segment(gridded, contour, start_coords, end_coords):

Method that will take a contour from the list of contours generated by
coast.Contour.get_contours() and trim it to start at supplied (lat,lon)
coordinates and end at supplied (lat, lon) coordinates.

Parameters
----------
gridded : Coast
                The gridded object containing the dataset with the 'bathymetry' variable
contour : numpy.ndarray
                contour[:,0] contains the y indices for the contour on the model grid
                contour[:,1] contains the x indices for the contour on the model grid
start_coords : numpy.ndarray
                1d array containing [latitude,longitude] of the start point of the contour
end_coords : numpy.ndarray
                1d array containing [latitude,longitude] of the end point of the contour

Returns
-------
y_ind : numpy.ndarray
                y indices of the contour on the model grid
x_ind : numpy.ndarray
                x indices of the contour on the model grid
contour : numpy.ndarray
                For the convenience of plotting using coast.Contour.plot_contour()
                contour[:,0] = y_ind
                contour[:,1] = x_ind

Contour.process_contour()

def Contour.process_contour(self, dataset, y_ind, x_ind):

Redefine contour so that each point on the contour defined by
y_ind and x_ind is seperated from its neighbours by a single index change
in y or x, but not both.
example: convert y_ind = [10,11], x_ind = [1,2] to y_ind = [10,10], x_ind = [1,2]
or y_ind = [10,11], x_ind = [1,1]

Parameters
----------
dataset : xarray.Dataset
                xarray Dataset from supplied gridded object
y_ind : numpy.ndarray
                1d array of y indices defining the contour on the model grid
x_ind : numpy.ndarray
                1d array of x indices defining the contour on the model grid

Returns
-------
y_ind : numpy.ndarray
                processed y indices of the contour on the model grid
x_ind : numpy.ndarray
                processed x indices of the contour on the model grid

Contour.gen_z_levels()
@staticmethod
def Contour.gen_z_levels(max_depth):

Generates a pre-defined 1d vertical depth coordinates,
i.e. horizontal z-level vertical coordinates up to a supplied
maximum depth, 'max_depth'

ContourF()

class ContourF(Contour):
Class defining a Contour type on the f-grid, which is a 3d dataset of points between a point A and
a point B defining an isobath contour. The dataset has a time, depth and contour dimension.
The contour dimension defines the points along the contour.
The supplied model f-grid Data is subsetted in its entirety along these dimensions
within Contour_f.data_contour of type xarray.Dataset


Parameters
----------
gridded_f : Coast
    f-grid gridded object containing the model dataset.
y_ind : numpy.ndarray
    1d array of y indices defining the contour on the model grid
x_ind : numpy.ndarray
    1d array of x indices defining the contour on the model grid
depth : int
    Depth of contour isobath
ContourF.calc_cross_contour_flow()

def ContourF.calc_cross_contour_flow(self, gridded_u, gridded_v):

Method that will calculate the flow across the contour and store this data
within Contour_f.data_cross_flow, which is an xarray.Dataset. Specifically
Contour_f.normal_velocities are the velocities across the contour
(time, depth, position along contour) in m/s
Contour_f.depth_integrated_normal_transport are the depth integrated
volume transports across the contour (time, position along contour) in Sv

If the time dependent cell thicknesses (e3) on the u and v grids are
present in the gridded_u and gridded_v datasets they will be used, if they
are not then the initial cell thicknesses (e3_0) will be used.

Parameters
----------
gridded_u : Coast
                The gridded object containing the model data on the u-grid.
gridded_v : Coast
                The gridded object containing the model data on the v-grid.

Returns
-------
None.

ContourF._update_cross_flow_vars()

def ContourF._update_cross_flow_vars(self, var, u_var, v_var, dr_n, dr_s, dr_e, dr_w, pos):

This method will pull variable data at specific points along the contour
from the u and v grid datasets and put them into the self.data_cross_flow dataset

ContourF._update_cross_flow_latlon()

def ContourF._update_cross_flow_latlon(self, ds_u, ds_v, dr_n, dr_s, dr_e, dr_w):

This method will pull the latitude and longitude data at specific points along the
contour from the u and v grid datasets and put them into the self.data_cross_flow dataset

ContourF._pressure_gradient_fpoint2()
@staticmethod
def ContourF._pressure_gradient_fpoint2(ds_t, ds_t_j1, ds_t_i1, ds_t_j1i1, r_ind, velocity_component):

Calculates the hydrostatic and surface pressure gradients at a set of f-points
along the contour, i.e. at a set of specific values of r_dim (but for all time and depth).
The caller must supply four datasets that contain the variables which define
the hydrostatic and surface pressure at all vertical z_levels and all time
on the t-points around the contour i.e. for a set of f-points on the contour
defined each defined at (j+1/2, i+1/2), we want t-points at (j,i), (j+1,i), (j,i+1), (j+1,i+1),
corresponding to ds_t, ds_t_j1, ds_t_i1, ds_t_j1i1, respectively.
ds_t, ds_t_j1, ds_t_i1, ds_t_j1i1 will have dimensions in time and depth.

The velocity_component defines whether u or v is normal to the contour
for the segments of the contour. A segment of contour is
defined as being r_dim to r_dim+1 where r_dim is the along contour dimension.


Returns
-------
hpg_f : DataArray with dimensions in time and depth and along contour
                hydrostatic pressure gradient at a set of f-points along the contour
                for all time and depth
spg_f : DataArray with dimensions in time and depth and along contour
                surface pressure gradient at a set of f-points along the contour

ContourF.calc_geostrophic_flow()

def ContourF.calc_geostrophic_flow(self, gridded_t, ref_density=None, config_u=config/example_nemo_grid_u.json, config_v=config/example_nemo_grid_v.json):

This method will calculate the geostrophic velocity and volume transport
(due to the geostrophic current) across the contour.
Four variables are added to the Contour.data_cross_flow dataset:
1. normal_velocity_hpg                (t_dim, depth_z_levels, r_dim)
This is the velocity due to the hydrostatic pressure gradient
2. normal_velocity_spg                (t_dim, r_dim)
This is the velocity due to the surface pressure gradient
3. transport_across_AB_hpg (t_dim, r_dim)
This is the volume transport due to the hydrostatic pressure gradient
4. transport_across_AB_spg (t_dim, r_dim
This is the volume transport due to the surface pressure gradient

This implementation works by regridding vertically onto horizontal z_levels in order
to perform the horizontal gradients. Currently s_level depths are
assumed fixed at their initial depths, i.e. at time zero.

Requirements: The gridded t-grid dataset, gridded_t, must contain the sea surface height,
Practical Salinity and the Potential Temperature variables. The depth_0
field must also be supplied. The GSW package is used to calculate
The Absolute Pressure, Absolute Salinity and Conservate Temperature.

Parameters
----------
gridded_t : Coast
                This is the gridded model data on the t-grid for the entire domain.
ref_density : TYPE, optional
                reference density value. If not supplied a mean in time, depth and
                along the contour will be used as the mean reference value.
config_u : file
                configuration file for u-grid object
config_v : file
                configuration file for v-grid object

Returns
-------
None.

ContourT()

class ContourT(Contour):
Class defining a Contour type on the t-grid, which is a 3d dataset of points between a point A and
a point B defining an isobath contour. The dataset has a time, depth and contour dimension.
The contour dimension defines the points along the contour.
The supplied model t-grid Data is subsetted in its entirety along these dimensions and
calculations can be performed on this dataset.

Parameters
----------
gridded_t : Coast
    t-grid gridded object containing the model dataset.
y_ind : numpy.ndarray
    1d array of y indices defining the contour on the model grid
x_ind : numpy.ndarray
    1d array of x indices defining the contour on the model grid
depth : int
    Depth of contour isobath
ContourT.construct_pressure()

def ContourT.construct_pressure(self, ref_density=None, z_levels=None, extrapolate=False):

                This method is for calculating the hydrostatic and surface pressure fields
                on horizontal levels in the vertical (z-levels). The motivation
                is to enable the calculation of horizontal gradients; however,
                the variables can quite easily be interpolated onto the original
                vertical grid.

                Requirements: The object's t-grid dataset must contain the sea surface height,
                Practical Salinity and the Potential Temperature variables.
                The GSW package is used to calculate the Absolute Pressure,
                Absolute Salinity and Conservate Temperature.

                Three new variables (density, hydrostatic pressure, surface pressure)
                are created and added to the Contour_t.data_contour dataset:
                                density_zlevels                (t_dim, depth_z_levels, r_dim)
                                pressure_h_zlevels                (t_dim, depth_z_levels, r_dim)
                                pressure_s                                                (t_dim, r_dim)

                Note that density is constructed using the EOS10
                equation of state.

Parameters
----------
ref_density: float
                reference density value, if None, then the Contour mean across time,
                depth and along contour will be used.
z_levels : (optional) numpy array
                1d array that defines the depths to interpolate the density and pressure
                on to.
extrapolate : boolean, default False
                If true the variables are extrapolated to the deepest z_level, if false,
                values below the bathymetry are set to NaN
Returns
-------
None.

ContourT.calc_along_contour_flow()

def ContourT.calc_along_contour_flow(self, gridded_u, gridded_v):

Function that will calculate the flow along the contour and store this data
within Contour_t.data_along_flow, which is an xarray.Dataset. Specifically
Contour_t.data_along_flow.velocities are the velocities along the contour with dimensions
(t_dim, z_dim, r_dim), where r_dim is the dimension along the contour.
Contour_t.data_along_flow.transport are the velocities along the contour multiplied by the
thickness of the cell (velocity * e3) with dimensions
(t_dim, z_dim, r_dim).

If the time dependent cell thicknesses (e3) on the u and v grids are
present in the gridded_u and gridded_v datasets they will be used, if they
are not then the initial cell thicknesses (e3_0) will be used.

Parameters
----------
gridded_u : Coast
                The nemo object containing the model data on the u-grid.
gridded_v : Coast
                The nemo object containing the model data on the v-grid.

Returns
-------
None.

ContourT.calc_along_contour_flow_2d()

def ContourT.calc_along_contour_flow_2d(self, gridded_u, gridded_v):

Function that will calculate the 2d flow (no vertical dimension
along the contour and store this data within Contour_t.data_along_flow,
which is an xarray.Dataset. Contour_t.data_along_flow.velocities are
the velocities along the contour with dimensions (t_dim, r_dim),
where r_dim is the dimension along the contour. e3 and
e3_0 are interpreted to be the water column thicknesses and
are included in the dataset as Contour_t.data_along_flow.e3 and
Contour_t.data_along_flow.e3_0

Parameters
----------
gridded_u : Coast
                The nemo object containing the model data on the u-grid.
gridded_v : Coast
                The nemo object containing the model data on the v-grid.

Returns
-------
None.

ContourT._update_flow_vars()

def ContourT._update_flow_vars(self, var, u_var, v_var, dr_n, dr_s, dr_e, dr_w, pos):

This method will pull variable data at specific points along the contour
from the u and v grid datasets and put them into the self.data_along_flow dataset

ContourT._update_along_flow_latlon()

def ContourT._update_along_flow_latlon(self, ds_u, ds_v, dr_n, dr_s, dr_e, dr_w):

This method will pull latitude and longitude data at specific points along the
contour from the u and v grid datasets and put them into the self.data_along_flow dataset


Last modified November 23, 2022: Updated docstrings from coast repo. (9daee18)