Grid Classes
This platecurie module contains the following Grid classes:
These classes can be initiatlized with a grid of corresponding data type,
and contain methods inherited from Grid
for the following functionality:
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
MagGrid
- class platecurie.classes.MagGrid(grid, dx, dy)
Basic grid class of
platecuriefor Magnetic anomaly data that inherits fromGridContains method to plot the grid data with default title and units using function
plot_real_grid().Additional Attributes
unitsstrUnits of Magnetic data (’\(nTesla\)’)
sg_unitsstrUnits of wavelet PSD (scalogram) (’\(nTesla^2/|k|\)’)
logsg_unitsstrUnits of log of wavelet PSD (log(scalogram)) (’\(log(nTesla^2/|k|)\)’)
title: strDescriptor for Magnetic anomaly data
>>> import numpy as np >>> from plateflex import Grid >>> from platecurie import MagGrid >>> maggrid = MagGrid(np.zeros((100, 100)), 10, 10) >>> isinstance(maggrid, Grid) True
ZtGrid
- class platecurie.classes.ZtGrid(grid, dx, dy)
Basic grid class of
platecuriefor the top of the magnetic layer that inherits fromGridContains method to plot the grid data with default title and units using function
plot_real_grid().Additional Attributes
unitsstrUnits of depth (‘m’)
title: strDescriptor for Depth to top of magnetic layer
Note
This class should only be used to specify the depth to the top of the magnetic layer at each cell location obtained from independent data. Although the
Gridmethods 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
(with at least one MagGrid). Methods are
available to:
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 fitted power-spectral density of the magnetic anomaly at single grid cell
Plot the final grids of model parameters
Save the results to disk as .csv files
Project
- class platecurie.classes.Project(grids=None)
Container for
MagGridand/orZtGridobjects, with methods to estimate model parameters of the magnetic layer as well as plot various results.- Parameters:
grids (list of
Grid, optional) – Initial list of grid objects.
Default Attributes
gridsListList of
MagGridobjectsmaskArray2D array of boolean values determined independently
initializedBoolWhether or not the project has been initialized and is ready for the estimation steps. By default this parameter is
False, unless the methodinit()has been executed.
Note
A Project can hold a list of any length with any type of
Gridor those defined inplatecurie- however the estimation will only proceed if the project holds exactly oneMagGridobject.Examples
- append(grid)
Append a single grid object to the current :class:`~platecurie.classes.Project object.
- Parameters:
grid (
Grid) – object to append to project
Example
>>> import numpy as np >>> from platecurie import MagGrid, Project >>> nn = 200; dd = 10. >>> grid = MagGrid(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 platecurie import MagGrid, ZtGrid, Project >>> nn = 200; dd = 10. >>> maggrid = MagGrid(np.random.randn(nn, nn), dd, dd) >>> ztgrid = ZtGrid(np.random.randn(nn, nn), dd, dd) >>> project = Project() >>> project.extend(grids=[maggrid, ztgrid])
- init(inverse='L2')
Method to initialize a project. This step is required before estimating the model parameters. The method checks that the project contains exactly one
MagGrid. It also ensures that all grids have the same shape and sampling intervals. If a grid of typeZtGrid(and possiblySigZtGrid) is present, the project attributes will be updated with data from the grid to be used in the estimation part.- Type:
str, optional
- Param:
Type of inversion to perform. By default the type is ‘L2’ for non-linear least-squares. Options are: ‘L2’ or ‘bayes’
Additional Attributes
nxintNumber of grid cells in the x-direction
nyintNumber of grid cells in the y-direction
nsintNumber of wavenumber samples
knp.ndarray1D array of wavenumbers
initializedboolSet to
Truewhen method is called successfully
Optional Attributes
ztndarrayGrid of depth to top of magnetic anomaly (km) (shape (nx,ny))
sztndarrayGrid of uncertainty in depth to top of magnetic anomaly (km) (shape (nx,ny))
Example
>>> from platecurie import MagGrid, ZtGrid, Project >>> nn = 200; dd = 10. >>> maggrid = MagGrid(np.random.randn(nn, nn), dd, dd) >>> ztgrid = ZtGrid(np.random.randn(nn, nn), dd, dd) >>> project = Project[grids=[maggrid, ztgrid]] >>> project.init()
- estimate_cell(cell=(0, 0), fix_beta=None, returned=False)
Method to estimate the parameters of the magnetized layer at a single cell location of the input grid.
- Parameters:
cell (tuple) – Indices of cell location within grid
fix_beta (float, optional) – If not None, fix
betaparameter - otherwise estimate itreturned (bool, optional) – Whether to use to return the estimates
Additional Attributes
traceMultiTracePosterior samples from the MCMC chains
summaryDataFrameSummary statistics from Posterior distributions
map_estimatedictContainer for Maximum a Posteriori (MAP) estimates
celltupleIndices of cell location within grid
Results are stored as attributes of
Projectobject.
- estimate_grid(nn=10, fix_beta=None, parallel=False)
Method to estimate the parameters of the magnetized layer 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)).fix_beta (float, optional) – If not None, fix
betaparameter - otherwise estimate it
Additional Attributes
mean_A_gridndarrayGrid with mean A estimates (shape
(nx, ny))MAP_A_gridndarrayGrid with MAP A estimates (shape
(nx, ny))std_A_gridndarrayGrid with std A estimates (shape
(nx, ny))mean_dz_gridndarrayGrid with mean dz estimates (shape
(nx, ny))MAP_dz_gridndarrayGrid with MAP dz estimates (shape
(nx, ny))std_dz_gridndarrayGrid with std dz estimates (shape
(nx, ny))
Optional Additional Attributes
mean_zt_gridndarrayGrid with mean zt estimates (shape
(nx, ny))MAP_zt_gridndarrayGrid with MAP zt estimates (shape
(nx, ny))std_zt_gridndarrayGrid with std zt estimates (shape
(nx, ny))mean_beta_gridndarrayGrid with mean beta estimates (shape
(nx, ny))MAP_beta_gridndarrayGrid with MAP beta estimates (shape
(nx, ny))std_beta_gridndarrayGrid with std beta estimates (shape
(nx, ny))
- plot_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 PSD function using one of
MAPormeanestimates. Calls the functionplot_fitted()with attributes as arguments.- Parameters:
est (str, optional) – Type of inference estimate to use for predicting th PSD
title (str, optional) – Title of plot
save (str, optional) – Name of file for to save figure
- plot_results(mean_A=False, MAP_A=False, std_A=False, mean_zt=False, MAP_zt=False, std_zt=False, mean_dz=False, MAP_dz=False, std_dz=False, mean_beta=False, MAP_beta=False, std_beta=False, chi2=False, mask=False, contours=None, filter=False, sigma=1, save=None, **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_A/zt/dz/beta (bool) – Type of plot to produce. All variables default to False (no plot generated)
- save_results(mean_A=False, MAP_A=False, std_A=False, mean_dz=False, MAP_dz=False, std_dz=False, mean_beta=False, MAP_beta=False, std_beta=False, mean_zt=False, MAP_zt=False, std_zt=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_A/dz/beta/zt (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
Trueprefix (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.