Grid Classes

plateflex defines the following Grid classes:

These classes can be initialized with a grid of the corresponding data type, and contain methods for the following functionality:

  • Extracting contours at some level of the grid

  • Performing a wavelet transform using a Morlet wavelet

  • Obtaining the wavelet scalogram from the wavelet transform

  • Plotting the input grids, wavelet transform components, and scalograms

Grid

class plateflex.classes.Grid(grid, dx, dy)

An object containing a 2D array of data and Cartesian coordinates specifying the bounding box of the array. Contains methods to calculate the wavelet transform, wavelet scalogram and to plot those quantities at a specified wavenumber index.

Parameters
  • grid (ndarray) – 2D array of of topography/gravity data

  • dx (float) – Grid spacing in the x-direction (km)

  • dy (float) – Grid spacing in the y-direction (km)

Grid must be projected in km.

Default Attributes

datandarray

2D array of topography/gravity data (shape (nx,ny))

dxfloat

Grid spacing in the x-direction in km

dyfloat

Grid spacing in the y-direction in km

nxint

Number of grid cells in the x-direction

nyint

Number of grid cells in the y-direction

unitsstr

Units of data set

sg_unitsstr

Units of power-spectral density of data set

logsg_unitsstr

Units of power-spectral density of data set in log

titlestr

Descriptor for data set - used in title of plots

ns int

Number of wavenumber samples

knp.ndarray

1D array of wavenumbers

>>> import numpy as np
>>> from plateflex import Grid
>>> # Create zero-valued square grid
>>> nn = 200; dd = 10.
>>> xmin = ymin = 0.
>>> xmax = ymax = (nn-1)*dd
>>> data = np.zeros((nn, nn))
>>> grid = Grid(data, dx, dy)
>>> grid
<plateflex.grids.Grid object at 0x10613fe10>
make_contours(level=0.0)

This function returns the contours as a List of coordinate positions at one given level - run this more than once with different levels if desired

Parameters

level (float) – Level (z value) of grid to contour

Returns

contours: List of contours with coordinate positions

wlet_transform()

This method uses the module cpwt to calculate the wavelet transform of the grid. By default the method mirrors the grid edges, tapers them (10 points on each edge) and pads the array with zeros to the next power of 2. The wavelet transform is stored as an attribute of the object.

Additional Attributes

wl_transndarray

Wavelet transform of the grid (shape (nx,ny,na,ns))

Example

>>> import numpy as np
>>> from plateflex import Grid
>>> # Create zero-valued square array
>>> nn = 200; dd = 10.
>>> xmin = ymin = 0.
>>> xmax = ymax = (nn-1)*dd
>>> data = np.zeros((nn, nn))
>>> grid = Grid(data, dx, dy)
>>> grid.wlet_transform()
 #loops = 13:  1  2  3  4  5  6  7  8  9 10 11 12 13
>>> grid.wl_trans.shape
(100, 100, 11, 13)
wlet_scalogram()

This method uses the module cpwt to calculate the wavelet scalogram of the grid.

Note

If no wl_trans attribute is found, the method automatically calculates the wavelet transform first.

The wavelet scalogram is stored as an attribute of the object.

Additional Attributes

wl_sgndarray

Wavelet scalogram of the grid (shape (nx,ny,na,ns))

Examples

>>> import numpy as np
>>> from plateflex import Grid
>>> # Create random-valued square grid
>>> nn = 200; dd = 10.
>>> xmin = ymin = 0.
>>> xmax = ymax = (nn-1)*dd
>>> # set random seed
>>> np.random.seed(0)
>>> data = np.random.randn(nn, nn)
>>> grid1 = Grid(data, dx, dy)
>>> grid1.wlet_scalogram()
 #loops = 17:  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17
>>> grid1.wl_sg.shape
(200, 200, 17)
>>> # Perform wavelet transform first
>>> grid2 = Grid(data, dx, dy)
>>> grid2.wlet_transform()
 #loops = 17:  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17
>>> grid2.wlet_scalogram()
>>> np.allclose(grid1.wl_sg, grid2.wl_sg)
True
plot_transform(kindex=None, aindex=None, log=False, mask=None, save=None, clabel=None, contours=None, **kwargs)

This method plots the real and imaginary components of the wavelet transform of a Grid object at wavenumber and angle indices (int). Raises Exception for the cases where:

  • no wavenumber OR angle index is specified (kindex and aindex)

  • wavenumber index is lower than 0 or larger than self.ns

  • angle index is lower than 0 or larger than 11 (hardcoded)

Note

If no wl_trans attribute is found, the method automatically calculates the wavelet transform first.

plot_scalogram(kindex=None, log=True, mask=None, title='Wavelet scalogram', save=None, clabel=None, contours=None, **kwargs)

This method plots the wavelet scalogram of a Grid object at a wavenumber index (int). Raises Exception for the cases where:

  • no wavenumber index is specified (kindex)

  • wavenumber index is lower than 0 or larger than self.ns - 1

Note

If no wl_sg attribute is found, the method automatically calculates the wavelet scalogram (and maybe also the wavelet transform) first.

TopoGrid

class plateflex.classes.TopoGrid(grid, dx, dy)

Basic grid class of plateflex for Topography data that inherits from Grid

Additional Attributes

unitsstr

Units of Topography (‘\(m\)’)

sg_unitsstr

Units of wavelet PSD (scalogram) (‘\(m^2/|k|\)’)

logsg_unitsstr

Units of log of wavelet PSD (log(scalogram)) (‘\(log(m^2/|k|)\)’)

title: str

Descriptor for Topography data

>>> import numpy as np
>>> from plateflex import Grid, TopoGrid
>>> nn = 200; dd = 10.
>>> topogrid = TopoGrid(np.random.randn(nn, nn), dd, dd)
>>> isinstance(topogrid, Grid)
True

Note

Automatically converts grid values to ‘meters’ if standard deviation is lower than 20

GravGrid

class plateflex.classes.GravGrid(grid, dx, dy)

Basic grid class of plateflex for gravity data that inherits from Grid

Additional Attributes

unitsstr

Units of Gravity anomaly (‘\(mGal\)’)

sg_unitsstr

Units of wavelet PSD (scalogram) (‘\(mGal^2/|k|\)’)

logsg_unitsstr

Units of log of wavelet PSD (log(scalogram)) (‘\(log(mGal^2/|k|)\)’)

title: str

Descriptor for Gravity data

>>> import numpy as np
>>> from plateflex import Grid, GravGrid
>>> nn = 200; dd = 10.
>>> gravgrid = GravGrid(np.random.randn(nn, nn), dd, dd)
>>> isinstance(gravgrid, Grid)
True

BougGrid

class plateflex.classes.BougGrid(grid, dx, dy)

Basic grid class of plateflex for Bouguer gravity data that inherits from GravGrid

Additional Attributes

title: str

Descriptor for Bouguer gravity data

>>> import numpy as np
>>> from plateflex import Grid, BougGrid, GravGrid
>>> nn = 200; dd = 10.
>>> bouggrid = BougGrid(np.random.randn(nn, nn), dd, dd)
>>> isinstance(bouggrid, GravGrid)
True
>>> isinstance(bouggrid, Grid)
True

FairGrid

class plateflex.classes.FairGrid(grid, dx, dy)

Basic grid class of plateflex for Free-air gravity data that inherits from GravGrid

Additional Attributes

title: str

Descriptor for Free-air gravity data

>>> import numpy as np
>>> from plateflex import Grid, FairGrid, GravGrid
>>> nn = 200; dd = 10.
>>> fairgrid = FaiorGrid(np.random.randn(nn, nn), dd, dd)
>>> isinstance(fairgrid, GravGrid)
True
>>> isinstance(fairgrid, Grid)
True

RhocGrid

class plateflex.classes.RhocGrid(grid, dx, dy)

Basic grid class of plateflex for crustal density data that inherits from Grid

Additional Attributes

unitsstr

Units of Density (‘\(kg/m^3\)’)

title: str

Descriptor for Density data

Note

This class should only be used to specify the density of the crust at each cell location. Although the Grid methods are still available, they are not useful in this context.

ZcGrid

class plateflex.classes.ZcGrid(grid, dx, dy)

Basic grid class of plateflex for crustal thickness data that inherits from Grid

Additional Attributes

unitsstr

Units of Crustal thickness (‘\(m\)’)

title: str

Descriptor for Crustal thickness data

Note

This class should only be used to specify the thickness of the crust at each cell location. Although the Grid methods are still available, they are not useful in this context.

Project Class

This module further contains the class Project, which itself is a container of Grid objects (at least one each of TopoGrid and GravGrid). Methods are available to:

  • Add Grid or Project objects to the project

  • Iterate over Grid objects

  • Initialize the project

  • Perform the wavelet admittance and coherence between topography (TopoGrid object) and gravity anomalies (GravGrid object)

  • Plot the wavelet admnittance and coherence spectra

  • Estimate model parameters at single grid cell

  • Estimate model parameters at every (or decimated) grid cell

  • Plot the statistics of the estimated parameters at single grid cell

  • Plot the observed and predicted admittance and coherence functions at single grid cell

  • Plot the final grids of model parameters

  • Save the results to disk as .csv files

Project

class plateflex.classes.Project(grids=None)

Container for Grid objects, with methods to calculate the wavelet admittance and coherence and estimate flexural model parameters as well as plot various results.

Parameters

grids (list of Grid, optional) – Initial list of PlateFlex Grid objects.

Default Attributes

gridsList

List of Grid objects

inversestr

Type of inversion to perform. By default the type is ‘L2’ for non-linear least-squares. Options are: ‘L2’ or ‘bayes’

maskArray

2D array of boolean values determined independently

initializedBool

Whether or not the project has been initialized and is ready for the wavelet analysis and estimation steps. By default this parameter is False, unless the method init() has been executed.

Note

Can hold a list of any length with any type of Grid objects - however the wavelet calculations will only proceed if the project holds one each of TopoGrid and GravGrid (BougGrid or FairGrid) objects.

Examples

Create zero-valued square grid

>>> import numpy as np
>>> from plateflex import Grid, TopoGrid, BougGrid, GravGrid, Project
>>> nn = 200; dd = 10.
>>> topogrid = TopoGrid(np.random.randn(nn, nn), dd, dd)
>>> bouggrid = BougGrid(np.random.randn(nn, nn), dd, dd)
>>> isinstance(bouggrid, GravGrid)
True
>>> isinstance(bouggrid, Grid)
True

Assign project with list of grids

>>> Project(grids=[topogrid, bouggrid])
<plateflex.classes.Project object at 0x10c62b320>

Add Grid to project

>>> project = Project()
>>> project += topogrid
>>> project += bouggrid
>>> project.grids
[<plateflex.classes.TopoGrid object at 0x1176b4240>, <plateflex.classes.BougGrid object at 0x1176b42e8>]

Append Grid to project

>>> project = Project()
>>> project.append(topogrid)
<plateflex.classes.Project object at 0x1176b9400>
>>> project.grids[0]
<plateflex.classes.TopoGrid object at 0x1176b4240>

Initialize project

>>> project.init()

Exception if project does not contain exactly one TopoGrid and one GravGrid

>>> project = Project(grids=[topogrid, topogrid])
>>> project.init()
Traceback (most recent call last):
  File "<stdin>", line 1, in <module>
  File "/Users/pascalaudet/Softwares/Python/plateflex-dev/PlateFlex/plateflex/classes.py", line 492, in wlet_admit_coh
    grids = self.grids + other.grids
Exception: There needs to be one GravGrid object in Project

Calculate wavelet admittance and coherence

>>> project = Project(grids=[topogrid, bouggrid])
>>> project.wlet_admit_coh()
 #loops = 17:  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17
 #loops = 17:  1  2  3  4  5  6  7  8  9 10 11 12 13 14 15 16 17
 Calculation jackknife error on admittance and coherence
append(grid)

Append a single Grid object to the current :class:`~plateflex.classes.Project object.

Parameters

grid (Grid) – object to append to project

Example

>>> import numpy as np
>>> from plateflex import Grid, Project
>>> nn = 200; dd = 10.
>>> grid = Grid(np.random.randn(nn, nn), dd, dd)
>>> project = Project()
>>> project.append(grid)
extend(grid_list)

Extend the current Project object with a list of Grid objects.

Parameters

trace_list – list of Grid objects or Project.

Example

>>> import numpy as np
>>> from plateflex import Grid, Project
>>> nn = 200; dd = 10.
>>> grid1 = Grid(np.random.randn(nn, nn), dd, dd)
>>> grid2 = Grid(np.random.randn(nn, nn), dd, dd)
>>> project = Project()
>>> project.extend(grids=[grid1, grid2])
init()

Method to initialize a project. This step is required before calculating the wavelet admittance and coherence. The method checks that the project contains one each of TopoGrid and GravGrid. It also ensures that all grids have the same shape and sampling intervals. If Grids of type RhocGrid and/or ZcGrid are present, the project attributes will be updated with data from those grids to be used in the estimation part. The method further sets the water depth at each cell location (grid point) from the TopoGrid object.

Additional Attributes

water_depthndarray

Grid of water depth from topography data (shape (nx,ny))

nxint

Number of grid cells in the x-direction

nyint

Number of grid cells in the y-direction

ns int

Number of wavenumber samples

knp.ndarray

1D array of wavenumbers

initializedbool

Set to True when method is called successfully

Optional Attributes

rhocndarray

Grid of crustal density data (shape (nx,ny))

zcndarray

Grid of crustal thickness data (shape (nx,ny))

Example

>>> project = Project[grids=[topogrid, fairgrid, zcgrid]]
>>> project.init()
wlet_admit_coh()

This method uses the module cpwt to calculate the wavelet admittance and coherence. The object needs to contain exactly two Grid objects, one of each of TopoGrid and GravGrid objects.

Note

If no wl_trans attribute is found for individual Grid objects, the method automatically calculates them first.

Stores the wavelet admittance, coherence and their errors as attributes

Additional Attributes

wl_admitndarray

Wavelet admittance (shape (nx,ny,ns))

ewl_admitndarray

Error of wavelet admittance (shape (nx,ny,ns))

wl_cohndarray

Wavelet coherence (shape (nx,ny,ns))

ewl_cohndarray

Error of wavelet coherence (shape (nx,ny,ns))

plot_admit_coh(kindex=None, mask=None, save=None, contours=None, **kwargs)

Method to plot grids of wavelet admittance and coherence at a given wavenumber index.

Parameters
  • kindex (int) – Index of wavenumber array

  • mask (ndarray, optional) – Array of booleans for masking data points

  • title (str, optional) – Title of plot

  • save (str, optional) – Name of file for to save figure

  • clabel (str, optional) – Label for colorbar

  • contours (List) – List of contour lines obtained separately

  • kwargs (Keyword arguments) – Keyword arguments allowing more control on plots

estimate_cell(cell=(0, 0), alph=False, atype='joint', returned=False)

Method to estimate the parameters of the flexural model at a single cell location of the input grids. The type of estimation performed is set by the project attribute inverse. See Project` and estimate for details. The model parameters rhoc and zc are extracted from the corresponding RhocGrid and ZcGrid objects if they were initialized in the project.

Parameters
  • cell (tuple) – Indices of cell location within grid

  • alph (bool, optional) – Whether or not to estimate parameter alpha

  • atype (str, optional) – Whether to use the admittance (‘admit’), coherence (‘coh’) or both (‘joint’)

  • returned (bool, optional) – Whether or not to return the estimates

Additional Attributes

traceMultiTrace

Posterior samples from the MCMC chains

summaryDataFrame

Summary statistics from Posterior distributions

map_estimatedict, optional

Container for Maximum a Posteriori (MAP) estimates

celltuple

Indices of cell location within grid

Results are stored as attributes of Project object.

estimate_grid(nn=10, alph=False, atype='joint', parallel=False)

Method to estimate the parameters of the flexural model at all grid point locations. It is also possible to decimate the number of grid cells at which to estimate parameters.

Parameters
  • nn (int) – Decimator. If grid shape is (nx, ny), resulting grids will have shape of (int(nx/nn), int(ny/nn)).

  • alph (bool, optional) – Whether or not to estimate parameter alpha

  • atype (str, optional) – Whether to use the admittance (‘admit’), coherence (‘coh’) or both (‘joint’)

Additional Attributes

mean_Te_gridndarray

Grid with mean Te estimates (shape (nx, ny))

MAP_Te_gridndarray

Grid with MAP Te estimates (shape (nx, ny))

std_Te_gridndarray

Grid with std Te estimates (shape (nx, ny))

mean_F_gridndarray

Grid with mean F estimates (shape (nx, ny))

MAP_F_gridndarray

Grid with MAP F estimates (shape (nx, ny))

std_F_gridndarray

Grid with std F estimates (shape (nx, ny))

Optional Attributes

mean_a_gridndarray

Grid with mean alpha estimates (shape (nx, ny)). Only present if alph=True.

MAP_a_gridndarray

Grid with MAP alpha estimates (shape (nx, ny)). Only present if alph=True.

std_a_gridndarray

Grid with std alpha estimates (shape (nx, ny)). Only present if alph=True.

chi2_gridndarray

Grid with reduced chi-squared estimates (shape (nx, ny)). Only present if project.inverse='L2'

new_maskArray

New grid of masked (boolean) values, corresponding to decimated mask.

plot_bayes_stats(title=None, save=None)

Method to plot the marginal and joint distributions of samples drawn from the posterior distribution as well as the extracted statistics. Calls the function plot_stats() with attributes as arguments.

Parameters
  • title (str, optional) – Title of plot

  • save (str, optional) – Name of file for to save figure

plot_functions(est='mean', title=None, save=None)

Method to plot observed and fitted admittance and coherence functions using one of mean or MAP estimates. The MAP is only available if the project attribute has been set to project='bayes'. Calls the function plot_functions() with attributes as arguments.

Parameters
  • est (str, optional) – Type of inference estimate to use for predicting admittance and coherence

  • title (str, optional) – Title of plot

  • save (str, optional) – Name of file for to save figure

plot_results(mean_Te=False, MAP_Te=False, std_Te=False, mean_F=False, MAP_F=False, std_F=False, mean_a=False, MAP_a=False, std_a=False, chi2=False, mask=False, contours=None, save=None, filter=True, sigma=1, **kwargs)

Method to plot grids of estimated parameters with fixed labels and titles. To have more control over the plot rendering, use the function plot_real_grid() with the relevant quantities and plotting options.

Parameters
  • mean/MAP/std_Te/F/a (bool) – Type of plot to produce. All variables default to False (no plot generated)

  • mask (bool) – Whether or not to plot the mask

  • contours (List) – List of contours with coordinate positions

  • filter (bool) – Whether or not to filter the resulting grid using a Gaussian filter

  • sigma (int) – Standard deviation of filter (in terms of adjacent cells), if set to True

  • kwargs (Keyword arguments) – Keyword arguments allowing more control on plots

Note

It is advisable to plot each grid separately using one call per grid, in order to have more control on plot parameters (e.g., colormap, min and max values of colorbar, etc.).

save_results(mean_Te=False, MAP_Te=False, std_Te=False, mean_F=False, MAP_F=False, std_F=False, mean_a=False, MAP_a=False, std_a=False, chi2=False, mask=False, contours=None, filter=True, sigma=1, prefix='PlateFlex_')

Method to save grids of estimated parameters with default file names.

Parameters
  • mean/MAP/std_Te/F/a (bool) – Type of grid to save. All variables default to False (no grid saved)

  • mask (bool) – Whether or not to save the mask

  • contours (List) – List of contours with coordinate positions

  • filter (bool) – Whether or not to filter the resulting grid using a Gaussian filter

  • sigma (int) – Standard deviation of filter (in terms of adjacent cells), if set to True

  • prefix (str) – Prefix of file name.

Note

It is advisable to save each grid separately using one call per grid, in order to have more control on file names.