comet.snmr.kernel package

Submodules

comet.snmr.kernel.kernel_bib module

Part of comet/snmr/kernel

This file contents parts of the MRSmatlab Kernel function part

class comet.snmr.kernel.kernel_bib.Kernel(survey=None, fid=0, dimension=1, name=None)[source]

Bases: object

Basic class to solve the NMR kernel computation.

Parameters:
  • name (string [ None ]) – If kernel is loaded from file.
  • survey (survey class instance [ None ]) – Calls setSurvey to define underlaying survey class. Holds important attributes like pulse moments and the loops for tx and rx.
  • tx (integer [ 0 ]) – Transmitter index in corresponding survey.
  • rx (integer [ 0 ]) – Receiver index in corresponding survey.
  • fid (interger [ 0 ]) – Sounding index in corresponding survey.
  • dimension (integer [1]) – Defines the kernel integration.

Example

>>> from comet.snmr import kernel as k
>>> from comet.snmr import survey
>>> site = survey.Survey()
>>> kernel = k.kernel(site)
>>> kernel.calculate()
>>> kernel.save('savename')
>>> kernel.show()
BFieldCalculation(loop_mesh=None, dipole_mesh=None, interpolate=False, just_loop_fields=False, recalc_loop_fields=False, recalc_primary=False, num_cpu=12, **kwargs)[source]

Calculates the Bfield for the kernel function for tx and rx.

internal call of loop.calculate() including decision if cell based or node based Bfield is needed. All optional parameters are piped to the loop.calculate() call. Based on the desired dimension of the kernel a specialised mesh may be automatically generated for the calculation.

Part 1/3 of the kernel calculation. Called automatically if kernel.calculate is called.

calcInterpolationMatrix()[source]
calcMagnetization()[source]

Creates 3D mesh and calcualtes magnetization vector after excitation. Returns magnetization vector of shape (num_pulses, num_cells_3d, 3)

calculate(loop_mesh=None, dipole_mesh=None, interpolate=False, savename=None, forceNew=False, slices=True, slice_name=None, **kwargs)[source]

All three parts of the kernel calulation are called here.

All given kwargs are directed to BfieldCalculation(), see function info for details about possible keyword arguments.

>>> self.BFieldCalculation(**kwargs)
>>> self.ellipticalDecomposition()
>>> self.kernelIntegration()
>>> if savename is not None:
        self.save(savename)
Keyword Arguments:
 destinations – none for now, with exception of “num_cpu”, [12] which is directed to BfieldCalculation and/or sliceKernel
coincident
create1DInterpolationSlices()[source]
create1DKernelMesh(max_length=0.1, area=100.0, quality=32, zvec=None, size_factor=2.5, z_factor=2.5, export_xyplane=None, max_dipoles=2000, calc_3D_stats=True, xmin=None, xmax=None, ymin=None, ymax=None)[source]

In order to integrate the kernel to a 1D structure without interpolation errors, a special mesh consisting of triangular zylinders has to be defined.

Parameters:
  • max_length (float [ 0.1 ]) – Defines the smallest edge length for the discretisation of the loop . In order to get admirable kernel results a value of 0.1 meters should be the maximum.
  • area (float [ 100. ]) – Defines the maximum Area a triangle in the loop slice can have.
  • quality (float [ 32. ]) – Defines the smallest angle inside a triangle. Be careful with values above 35.
  • zvec (array_like [ None ]) – Usualy the zvec is defined automatically, this flag gives the user the optional possibility to give a zvec from outside the funktion.
  • size_factor (float [ 2.5 ]) – Extension of the kernel mesh (and therefore integration volume) in the x and y direction. Should be at least 2 times the loop diameter or shortest edge length. This value defines the multipier.
  • z_factor (float [ 2.5 ]) – Maximum depth of the Kernel. Should be at least 2 times the loop diameter or shortest edge length. This value defines the multipier.
  • export_xyplane (string [ None ]) – Filename for the resulting kernel mesh plane in 2D can be exported for debugging or simply to check the mesh (vtk).
  • max_dipoles (interger [ 2000 ]) – Fallback for high node density loops. This sets an overall maximum for the number of dipoles used for the loop discretization. However this only comes into account in rare cases.
create2DInterpolationSlices()[source]
create2DKernelMesh(area=15.0, quality=34, yvec=None, x_factor=5, z_factor=2, savename=None, export_xzplane=None, calc_3D_stats=True, order=0)[source]

Similary to the mesh in the 1D case a special mesh consisting of triangluar zylinders is generated. The Zylinders are pointing in the y direction to allow a perfect integration to the x-z plane.

Parameters:
  • area (float [15.]) – Affects the maximum area a triangle in the 2D slice is allowed to have. Higher Values lead to bigger cells.
  • quality (float [34]) – Defines the smallest angle inside a triangle. Be careful with values above 34.5. Higher values = more cells.
  • yvec (ndarray, list [None]) – Usualy the y vector is defined automatically, this flag gives the user the optional possibility to give a YVec from outside the function.
  • x_factor (float [2]) – Extension of the kernel mesh (and therefore integration volume) in the x direction. Should be at least 2 times the loop diameter or shortest edge length. This value defines the multipier.
  • z_factor (float [2]) – Extension of the kernel mesh (and therefore integration volume) in the z direction. Should be at least 2 times the loop diameter or shortest edge length. This value defines the multipier.
  • savename (string [None]) – If a savename is given, the resulting 2D Mesh is saved in the .bms format for later use.
  • export_xyplane (string [ None ]) – Filename for the resulting kernel mesh plane in 2D can be exported for debugging or simply to check the mesh (vtk).
createMagnetizationMesh()[source]

Creates full 3D mesh for display and calcualtion of magnetization vectors. Not needed for normal kernel calculation routine and big, therefore separate.

createSeperatedLoopMesh(name='SepLoopMesh', dipole=True, exportVTK=False, refinement_para=1.0, max_area_factor=1.0)[source]

Creates a mesh that contains the receiver and the transmitter loop.

createYVec(max_length=0.2, max_num=300, y_factor=2.0, calc_3D_stats=True)[source]

Creates the y vector discretization for the 2D kernel mesh.

The y vector represents the y values of the 3D Kernel mesh before the integration to 2D.

Parameters:
  • max_length (float [ 0.2 ]) – Maximum distance between two slices inbetween the source dipoles.
  • max_num (integer [ 300 ]) – Maximum number of slices. Overrides max_length if they conflict.
  • y_factor (float [ 2. ]) – Extension of the kernel mesh (and therefore integration volume) in the y direction. Should be at least 2 times the loop diameter or shortest edge length. This value defines the multipier.
createZVector(numz, minz, min_thk=0.5)[source]

Creates a sinus hyperbolicus shaped Z discretisation in numz steps between 0 and minz.

ellipticalDecomposition()[source]

Computes the counter and corotating parts of the given magnetic fields with respect to a given earth magnetic field.

Parameters:
  • Bfield (complex field [3, n] or string) – Optional. Possibility to insert a pre calculated field.
  • Inclination (float) – Inclination of the earth magnetic field at the loop site in rad [0… 2pi]
  • Declination (float) – Declination of the magnetic field at the loop site in rad [0… 2pi]
  • B (np.array of shape (3, n)) – Magnetic field of the loop
  • Second part of the kernel calculation.
  • - mainly from Weichman et al. (2000)
static ellipticalDecomposition_multi(Bfield, earth)[source]

Computes the counter and corotating parts of the given magnetic fields with respect to a given earth magnetic field.

Parameters:
  • Bfield (complex field [3, n] or string) – Optional. Possibility to insert a pre calculated field.
  • Inclination (float) – Inclination of the earth magnetic field at the loop site in rad [0… 2pi]
  • Declination (float) – Declination of the magnetic field at the loop site in rad [0… 2pi]
  • B (np.array of shape (3, n)) – Magnetic field of the loop
  • Second part of the kernel calculation.
  • Literature
  • ———-
  • - Weichman et al. (2000)
  • - Hertrich (2005, Appendix)
  • - Hertrich (2008, eq. 6 ff.)
export1D(savename, loop_layout=True, title='{0.survey.earth!r}')[source]
export2DKernel(fig=None, ax=None, savename=None, png_dpi=300, noYLabel=False, index=0, colorBar=True, size=13, pdf=None, fixed_cbar=False, **kwargs)[source]

Exports 2D Kernel for given pulse moment. Kwargs are redirected to show.

export2DKernel2PDF(name, fixed_cbar=False, **kwargs)[source]

Export 2D Kernel for all pulse moments as stiched pdf. Kwargs are redirected to export2DKernel.

exportMagnetization(name, vtk_export=False, pulse=0)[source]

Export a previously calculated magnetization vector as numpy vector and optionally vtk file.

exportVTK(savename, save=['abs'], save_in_log=False)[source]
fid

Reference to sounding (FID) class instance in survey.

getKernel(reduced=True)[source]
getSliceCoords()[source]

Returns input coordinates for custEM Slice interpolation of magnetic fields to the kernel slices.

getZVector(reduced=True)[source]
interpolateBFieldToKernel(recalc_prim_on_kernel=False, recalc_primary=False, num_cpu=32, calc_3D_stats=True)[source]

Takes the rx Bfield and interpolates it to the kernel mesh.

static kernelCalculation_multi(fid, earth, txalpha, txbeta, txzeta, txperpend, rxalpha=None, rxbeta=None, rxzeta=None, rxperpend=None, calc_theta=False)[source]
kernelIntegration(calc_theta=False)[source]

Computes the integration of the kernel with respect to the desired dimension.

Parameters:
  • decomposition ((alpha, beta, zeta)) – Bfield_part essentially consists of the output from the elliptical decomposition of the magnetic field.
  • measurement (class) – An instance of a measurement class has to be given in order to keep the number of input arguments manageable.
  • earthmagnitude (float) – Magnitude of the earth magnetic field [Tesla]. Aproximatly about 30000 to 65000 nT (1 nT = 1e-9 Tesla).
  • Third part of the kernel calculation.
larmor

Larmor frequency [Hz] from earth defined in survey.

load(savename, load_loopmesh=True, kernelmesh2d=None, load_kernelmesh=True, use_order_refinement=True)[source]

Load a previously saved kernel (.npz-format).

magnetization_magnitude
pulses

Reference to pulse moments from sounding (FID).

release_memory()[source]

Calling this function is releasing some attributes that are using a fairly big amount of memory.

Sets the following attributes back to None:

  • The interpolation matrix between the loop meshes and the kernel mesh

interpolationMatrix

  • local copies of the magnetic fields (fields in tx and rx are not

effected) txBfield, rxBfield

  • the 3D kernel mesh cell center and volumes

kernelMeshCellVolume, kernelMeshCellCenter

  • the elliptical decomposition of the tx and rx bfields

txalpha, txbeta, txzeta, txperpend, rxalpha, rxbeta, rxzeta, rxperpend

Note: a recalculation of the kernel will take about the same amount of time as the first call, as all cached variables are gone, however apart from a recalculation, the other purposes of the kernel class (export, figures, inversion(without recalculation)) are not effected.

Another note: If you want to use this method only for saving disk space in case you save the kernel class, then you might consider the light flag of the .save method instead.

rx

Reference to receiver class instance in survey.

rx_area

Area of the receiver loop.

rx_index
save(savename=None, save_interpolation_mat=False, save_loopmesh=False, light=True, kernelmesh_name=None)[source]

Save the basic information to restore the Kernel class later.

set1DKernelMesh(mesh, calc_3D_stats=True)[source]

Sets the 1D kernel mesh.

Parameters:
  • mesh (stirng or pygimli.Mesh) – Filename or mesh instance of a 2D mesh in the x-y plane.
  • Need
  • —-
  • z discretization – Can be setted via createZVector, setZVector or direct use of create1DKernelMesh. However the needed information to do that may not be available on the fly, therefore no default z vector is created.
set2DKernelMesh(inmesh, yvec=None, order=0, integration_mat=None, calc_3D_stats=True)[source]

kwargs to createYVec if YVec is None

setDebug(debug: bool)[source]
setModel(*args, **kwargs)[source]

Pipes args and kwargs to self.tx.setModel. Same for rx.

setPulsesDirectly(pulses)[source]

Set pulse moment vector manually if not supported by survey + fid. (This is called when loading a kernel from the harddisk, mainly for plotting reasons). For all calculation purposes a survey and fid class is recommended.

setRx(rx, **kwargs)[source]

Sets initialized loop or pipe arg and kwargs to loadLoop.

setSurvey(survey, fid=0)[source]

Sets survey class containing necessary information for the kernel.

Parameters:
  • survey (comet.snmr.survey.Survey or None) – Sets given survey class instance or create empty class instance.
  • fid (integer [ 0 ]) – Index of corresponding sounding in the survey.
setTx(tx, **kwargs)[source]

Sets initialized loop or pipe arg and kwargs to loadLoop.

setYVector(vector)[source]
setZVector(vector, min_thk=0.5)[source]

Defines the attribute zvec.

Sets the given vector as z discretization. Attention: the value for min_thk defines the minimum thickness of the discretization used in the end. For all thicknesses in vector smaller than min_thk, the Kernel is integrated to match the min_thk. For calulation of the kernel function the original given vector is used.

Parameters:
  • vector (array_like) – Z discretization in m to be used for the kernel calculation. If a new vector is to be created, please also take a look at the method createZVector.
  • min_thk (float) – Minimum thickness te kernel and zvec is integrated if returned. This leads to higher accuracy in the vicinity of the loop.
shape
show(toplot=['real', 'imag', 'amp', 'phase', '0D'], indices=None, savename=None, normed=False, suptitle=None, ax=None, pulse_in_log=False, kernel_absolute_values=False, cbar_percentage=0.99, fixed_cbar=False, lut=33, show_marked_edges=False, **kwargs)[source]

Visualise the Kernel with respect to the desired dimension.

Automatically defined within the kernel class via the parameter kernel.dimension = [0…3]. Plotting of a kernel in the desired dimension is only possible if the kernel is also calculated with respect to that dimension. It’s not possible to calculate the kernel with kernel.dimension = 1 and then plot the kernel with kernel.dimension = 2.

0D :
Simple Graph plotting kernel-values over pulsemoments
1D :
Graph with 1D integrated kernels over the depth of the model
2D :
Slice of the x-z-plane with triangle mesh containing the 2D
3D :
Export of the kernel in vtk format for visualising.

none so far

Plots the 1D integrated Kernel with a given z discretisation over the measured pulse sequences.

toplot: list [ [‘real’, ‘imag’, ‘amp’, ‘phase’, ‘1D’] ]
There are different possibilities to plot the kernel. This parameter defines which part of the kernel is shown. Possible options are: ‘real’, ‘imag’, ‘amp’, ‘phase’, ‘0D’ (integrated over z). All strings in the toplot variable will be plotted in the same order given in the list.
cMap: string [‘viridis’]
Defines the colormap used to display the kernel. In order to get a good contrast between the max and min as well as being useful in comparison with MRSMatlab, ‘viridis’ is the default colormap. Any colormap reachable by the plt.get_cmap(…) method can be chosen.
normed: bool [True]
A on the dimension based normalisation of the plot permits a better assessment of the kernel distribution.
ax: plotting ax or list of axes [None]
Plot on a predefined ax and gives back the ax. A onedimensionla list of axes is also accepted, if the number of items in ‘toplot’ is the same as the available axes.
lut: None or int [None]
Number of colors for the colorbar. If lut is not None it must be an integer giving the number of entries desired in the lookup table, and name must be a standard mpl colormap name.
indices: list
By default one 2D plot is created for each pulsemoment. In order to limit the number of plots the optional paramter indices can be given as a list of indices referring to the pulse moments to be shown.
cMap: string [‘viridis’]
See Parameter 1D.
normed: bool [True]
A on the dimension based normalisation of the plot permits a better assessment of the kernel distribution.
show_marked_edges: boolean [ False ]
Whether or not marked edges gets drawn.
possible kwargs for matplotlib:
cMin, cMax for range of the colorbar. All other kwargs are reaching matplotlib functions.
default label 2D:
‘integrated kernel (2D) [nV/$m^2$] pulsemoment: {:.3f} As’ .format(self.pulses[i])

A self-sufficient plot of the kernel without any integration would result in a set of 3D Cubes and is not implemented for now.

Instead the kernel will be saved in vtk format which can be easily handled.

savename: string
A String defining the relative path to the vtk-file the kernel will be saved in. If not given the default savename will be flagged with the string ‘_default_’ and contain some information about the kernel.

Example

2D:

>>> ax, cbar = kernel.show(indices=[16], cMin=-1,
>>>                        cMax=2, size=20, pad=0.7)
>>> ax.set_ylim(-50, 0)
show2DMesh()[source]
showLoopLayout(ax=None, **kwargs)[source]
sliceKernel1D(num_cpu=None, loop_mesh=None, new_bfield=False, interpolate_bfield=True, slice_name=None)[source]
sliceKernel2D(savename=None, forceNew=False, loopSaveName=None, num_cpu=None, new_bfield=False, loop_mesh=None, slice_name=None, **kwargs)[source]

2D Kernel in a memory saving parallel computation approach.

tx

Reference to transmitter class instance in survey.

tx_area

Area of the transmitter loop.

tx_index
updatable
zvec

z discretisation

comet.snmr.kernel.kernel_bib.calcInterpolationMatrix_para(source_mesh, target_pos, num_cpu=8)[source]
comet.snmr.kernel.kernel_bib.calculateKernelFromSlices(survey_name, invmesh, cfgname, max_length=0.05, max_num=400, path_name='kernel', num_cpu=48, h_order=1, json_name=None, force_new_paths=True, kernel_name='kernel/kern_{}', slice_export_name='kernel_slice')[source]
survey: string
Filepath of the survey containing the FIDs.
invmesh: string
Filepath for the inversion mesh the kernel is calculated on.
cfgname: string
Filepath of the secondary config containingin formation for custEM. The same file should have been used for the field calculation.
max_length: float [ 0.05 ]
Minimum distance between two slices in y direction.
max_num: integer [ 400 ]
Number of slices for y discretization.
path_name: string [ ‘kernel’ ]
Final name for the slices will be {mdir}/paths/{path_name}_{number}_path.xml “mdir” is defined in the secondary config.
num_cpu: integer [ 48 ]
Number of cores used.
h_order: integer [ 1 ]
Order of h-refinement when setting the invmesh in the kernelclass.
json_name: string [ None ]
Optional name for the json file. Alternatively a tempory file is created.
force_new_path: boolean [ True ]
Deletes old pathfiles (slices) and forces the generation of new ones.
kernel_name: string [ ‘kernel/kern_{}’ ]
Filepath of used for export of the kernel functions. Need to contain a “{}” which is filled with the index of the corresponding FID in the survey class.
slice_export_name: string [ ‘kernel_slice’ ]
Filepath of the interpolated magnetic fields for the individual kernel slices. Full slice path contains: “{r_dir}/{approach}/{mesh_name}/{mod_name}_interpolated/ tx_{tx_number}_{slice_export_name}_imesh_{slice_number}.npy”
kernel_name (added a “_{}” if not in original string)
comet.snmr.kernel.kernel_bib.checkForKernel(name, mkdir=False)[source]
comet.snmr.kernel.kernel_bib.create1DInterpolationSlices(kern)[source]
comet.snmr.kernel.kernel_bib.create2DInterpolationSlices(kern)[source]
comet.snmr.kernel.kernel_bib.integrateKernelH2(mat, array)[source]
comet.snmr.kernel.kernel_bib.simpleZVec(numz, minz, reduced=False)[source]

Module contents

Module comet/snmr/kernel