Maps

That’s not a map...?!

A SimPEG Map operates on a vector and transforms it to another space. We will use an example commonly applied in electromagnetics (EM) of the log-conductivity model.

\[m = \log(\sigma)\]

Here we require a mapping to get from \(m\) to \(\sigma\), we will call this map \(\mathcal{M}\).

\[\sigma = \mathcal{M}(m) = \exp(m)\]

In SimPEG, we use a (SimPEG.Maps.ExpMap) to describe how to map back to conductivity. This is a relatively trivial example (we are just taking the exponential!) but by defining maps we can start to combine and manipulate exactly what we think about as our model, \(m\). In code, this looks like

1
2
3
4
5
6
7
8
9
M = Mesh.TensorMesh([100]) # Create a mesh
expMap = Maps.ExpMap(M)    # Create a mapping
m = np.zeros(M.nC)         # Create a model vector
m[M.vectorCCx>0.5] = 1.0   # Set half of it to 1.0
sig = expMap * m           # Apply the mapping using *
print m
# [ 0.    0.    0.    1.     1.     1. ]
print sig
# [ 1.    1.    1.  2.718  2.718  2.718]

Combining Maps

We will use an example where we want a 1D layered earth as our model, but we want to map this to a 2D discretization to do our forward modeling. We will also assume that we are working in log conductivity still, so after the transformation we want to map to conductivity space. To do this we will introduce the vertical 1D map (SimPEG.Maps.SurjectVertical1D), which does the first part of what we just described. The second part will be done by the SimPEG.Maps.ExpMap described above.

1
2
3
4
5
6
M = Mesh.TensorMesh([7,5])
v1dMap = Maps.SurjectVertical1D(M)
expMap = Maps.ExpMap(M)
myMap = expMap * v1dMap
m = np.r_[0.2,1,0.1,2,2.9] # only 5 model parameters!
sig = myMap * m
../../_images/sphx_glr_plot_combo_001.png

If you noticed, it was pretty easy to combine maps. What is even cooler is that the derivatives also are made for you (if everything goes right). Just to be sure that the derivative is correct, you should always run the test on the mapping that you create.

Taking Derivatives

Now that we have wrapped up the mapping, we can ensure that it is easy to take derivatives (or at least have access to them!). In the SimPEG.Maps.ExpMap there are no dependencies between model parameters, so it will be a diagonal matrix:

\[\left(\frac{\partial \mathcal{M}(m)}{\partial m}\right)_{ii} = \frac{\partial \exp(m_i)}{\partial m} = \exp(m_i)\]

Or equivalently:

\[\frac{\partial \mathcal{M}(m)}{\partial m} = \text{diag}( \exp(m) )\]

The mapping API makes this really easy to test that you have got the derivative correct. When these are used in the inverse problem, this is extremely important!!

import numpy as np
from SimPEG import Mesh, Maps
import matplotlib.pyplot as plt
M = Mesh.TensorMesh([100])
expMap = Maps.ExpMap(M)
m = np.zeros(M.nC)
m[M.vectorCCx>0.5] = 1.0
expMap.test(m, plotIt=True)

(Source code, png, hires.png, pdf)

../../_images/api_Maps-1.png

Common Maps

Exponential Map

Electrical conductivity varies over many orders of magnitude, so it is a common technique when solving the inverse problem to parameterize and optimize in terms of log conductivity. This makes sense not only because it ensures all conductivities will be positive, but because this is fundamentally the space where conductivity lives (i.e. it varies logarithmically).

class SimPEG.Maps.ExpMap(mesh=None, nP=None, **kwargs)

Electrical conductivity varies over many orders of magnitude, so it is a common technique when solving the inverse problem to parameterize and optimize in terms of log conductivity. This makes sense not only because it ensures all conductivities will be positive, but because this is fundamentally the space where conductivity lives (i.e. it varies logarithmically).

Changes the model into the physical property.

A common example of this is to invert for electrical conductivity in log space. In this case, your model will be log(sigma) and to get back to sigma, you can take the exponential:

\[ \begin{align}\begin{aligned}m = \log{\sigma}\\\exp{m} = \exp{\log{\sigma}} = \sigma\end{aligned}\end{align} \]
deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model

The transform changes the model into the physical property. The transformDeriv provides the derivative of the transform.

If the model transform is:

\[ \begin{align}\begin{aligned}m = \log{\sigma}\\\exp{m} = \exp{\log{\sigma}} = \sigma\end{aligned}\end{align} \]

Then the derivative is:

\[\frac{\partial \exp{m}}{\partial m} = \text{sdiag}(\exp{m})\]
inverse(D)
Parameters:D (numpy.array) – physical property
Return type:numpy.array
Returns:model

The transformInverse changes the physical property into the model.

\[m = \log{\sigma}\]

Vertical 1D Map

class SimPEG.Maps.SurjectVertical1D(mesh, **kwargs)

SurjectVertical1DMap

Given a 1D vector through the last dimension of the mesh, this will extend to the full model space.

deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

Map 2D Cross-Section to 3D Model

class SimPEG.Maps.Surject2Dto3D(mesh, **kwargs)

Map2Dto3D

Given a 2D vector, this will extend to the full 3D model space.

deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

normal = u'Y'

Mesh to Mesh Map

class SimPEG.Maps.Mesh2Mesh(meshes, **kwargs)
Takes a model on one mesh are translates it to another mesh.

Required Properties:

  • indActive (Array): active indices on target mesh, a list or numpy array of <type ‘bool’> with shape (*)
P
deriv(m, v=None)
indActive

indActive (Array) – active indices on target mesh, a list or numpy array of <type ‘bool’> with shape (*)

nP

Number of parameters in the model.

shape

Number of parameters in the model.

../../_images/sphx_glr_plot_mesh2mesh_001.png

Under the Hood

Combo Map

The ComboMap holds the information for multiplying and combining maps. It also uses the chain rule to create the derivative. Remember, any time that you make your own combination of mappings be sure to test that the derivative is correct.

class SimPEG.Maps.ComboMap(maps, **kwargs)

Combination of various maps.

The ComboMap holds the information for multiplying and combining maps. It also uses the chain rule to create the derivative. Remember, any time that you make your own combination of mappings be sure to test that the derivative is correct.

deriv(m, v=None)
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

shape

The API

The IdentityMap is the base class for all mappings, and it does absolutely nothing.

class SimPEG.Maps.ActiveCells(mesh, indActive, valInactive, nC=None)

Bases: SimPEG.Maps.InjectActiveCells

ActiveCells is depreciated. Use InjectActiveCells instead

class SimPEG.Maps.BaseParametric(mesh, **kwargs)

Bases: SimPEG.Maps.IdentityMap

indActive = None
slope = None
slopeFact = 100.0
x
y
z
class SimPEG.Maps.ChiMap(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Chi Map

Convert Magnetic Susceptibility to Magnetic Permeability.

\[\mu(m) = \mu_0 (1 + \chi(m))\]
deriv(m, v=None)
inverse(m)
class SimPEG.Maps.CircleMap(mesh, logSigma=True)

Bases: SimPEG.Maps.ParametricCircleMap

CircleMap is depreciated. Use ParametricCircleMap instead

class SimPEG.Maps.ComboMap(maps, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Combination of various maps.

The ComboMap holds the information for multiplying and combining maps. It also uses the chain rule to create the derivative. Remember, any time that you make your own combination of mappings be sure to test that the derivative is correct.

deriv(m, v=None)
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

shape
class SimPEG.Maps.ComplexMap(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

default nP is nC in the mesh times 2 [real, imag]

deriv(m, v=None)
nP
shape
class SimPEG.Maps.ExpMap(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Electrical conductivity varies over many orders of magnitude, so it is a common technique when solving the inverse problem to parameterize and optimize in terms of log conductivity. This makes sense not only because it ensures all conductivities will be positive, but because this is fundamentally the space where conductivity lives (i.e. it varies logarithmically).

Changes the model into the physical property.

A common example of this is to invert for electrical conductivity in log space. In this case, your model will be log(sigma) and to get back to sigma, you can take the exponential:

\[ \begin{align}\begin{aligned}m = \log{\sigma}\\\exp{m} = \exp{\log{\sigma}} = \sigma\end{aligned}\end{align} \]
deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model

The transform changes the model into the physical property. The transformDeriv provides the derivative of the transform.

If the model transform is:

\[ \begin{align}\begin{aligned}m = \log{\sigma}\\\exp{m} = \exp{\log{\sigma}} = \sigma\end{aligned}\end{align} \]

Then the derivative is:

\[\frac{\partial \exp{m}}{\partial m} = \text{sdiag}(\exp{m})\]
inverse(D)
Parameters:D (numpy.array) – physical property
Return type:numpy.array
Returns:model

The transformInverse changes the physical property into the model.

\[m = \log{\sigma}\]
class SimPEG.Maps.FullMap(mesh, **kwargs)

Bases: SimPEG.Maps.SurjectFull

FullMap is depreciated. Use SurjectVertical1DMap instead

class SimPEG.Maps.IdentityMap(mesh=None, nP=None, **kwargs)

Bases: properties.base.base.HasProperties

SimPEG Map

deriv(m, v=None)

The derivative of the transformation.

Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
inverse(D)

Changes the physical property into the model.

Note

The transformInverse may not be easy to create in general.

Parameters:D (numpy.array) – physical property
Return type:numpy.array
Returns:model
nP

rtype – int :return: number of parameters that the mapping accepts

shape

The default shape is (mesh.nC, nP) if the mesh is defined. If this is a meshless mapping (i.e. nP is defined independently) the shape will be the the shape (nP,nP).

Return type:tuple
Returns:shape of the operator as a tuple (int,int)
test(m=None, num=4, **kwargs)

Test the derivative of the mapping.

Parameters:
Return type:

bool

Returns:

passed the test?

testVec(m=None, **kwargs)

Test the derivative of the mapping times a vector.

Parameters:
Return type:

bool

Returns:

passed the test?

class SimPEG.Maps.InjectActiveCells(mesh, indActive, valInactive, nC=None)

Bases: SimPEG.Maps.IdentityMap

Active model parameters.

deriv(m, v=None)
indActive = None
inverse(D)
nP

Number of parameters in the model.

shape
valInactive = None
class SimPEG.Maps.LogMap(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Changes the model into the physical property.

If (p) is the physical property and (m) is the model, then

\[p = \log(m)\]

and

\[m = \exp(p)\]

NOTE: If you have a model which is log conductivity (ie. (m = log(sigma))), you should be using an ExpMap

deriv(m, v=None)
inverse(m)
class SimPEG.Maps.Map2Dto3D(mesh, **kwargs)

Bases: SimPEG.Maps.Surject2Dto3D

Map2Dto3D is depreciated. Use Surject2Dto3D instead

class SimPEG.Maps.Mesh2Mesh(meshes, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Takes a model on one mesh are translates it to another mesh.

Required Properties:

  • indActive (Array): active indices on target mesh, a list or numpy array of <type ‘bool’> with shape (*)
P
deriv(m, v=None)
indActive

indActive (Array) – active indices on target mesh, a list or numpy array of <type ‘bool’> with shape (*)

nP

Number of parameters in the model.

shape

Number of parameters in the model.

class SimPEG.Maps.MuRelative(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Invert for relative permeability

\[\mu(m) = \mu_0 * \mathbf{m}\]
deriv(m, v=None)
inverse(m)
class SimPEG.Maps.ParametricBlock(mesh, **kwargs)

Bases: SimPEG.Maps.BaseParametric

Parametric Block in a Homogeneous Space

For 2D:

m = [
    val_background,
    val_block,
    block_x0,
    block_dx,
    block_y0,
    block_dy
]

For 3D:

m = [
    val_background,
    val_block,
    block_x0,
    block_dx,
    block_y0,
    block_dy
    block_z0,
    block_dz
]

Required

Parameters:mesh (discretize.BaseMesh.BaseMesh) – SimPEG Mesh, 2D or 3D

Optional

Parameters:
  • slopeFact (float) – arctan slope factor - divided by the minimum h spacing to give the slope of the arctan functions
  • slope (float) – slope of the arctan function
  • indActive (numpy.ndarray) – bool vector with active indices
deriv(m)
mDict(m)
nP
shape
xleft(mDict)
xright(mDict)
yleft(mDict)
yright(mDict)
zleft(mDict)
zright(mDict)
class SimPEG.Maps.ParametricBlockInLayer(mesh, **kwargs)

Bases: SimPEG.Maps.ParametricLayer

Parametric Block in a Layered Space

For 2D:

m = [val_background,
     val_layer,
     val_block,
     layer_center,
     layer_thickness,
     block_x0,
     block_dx
]

For 3D:

m = [val_background,
     val_layer,
     val_block,
     layer_center,
     layer_thickness,
     block_x0,
     block_y0,
     block_dx,
     block_dy
]

Required

Parameters:mesh (discretize.BaseMesh.BaseMesh) – SimPEG Mesh, 2D or 3D

Optional

Parameters:
  • slopeFact (float) – arctan slope factor - divided by the minimum h spacing to give the slope of the arctan functions
  • slope (float) – slope of the arctan function
  • indActive (numpy.ndarray) – bool vector with
deriv(m)
mDict(m)
nP
shape
xleft(mDict)
xright(mDict)
yleft(mDict)
yright(mDict)
class SimPEG.Maps.ParametricCasingAndLayer(mesh, **kwargs)

Bases: SimPEG.Maps.ParametricLayer

Parametric layered space with casing.

m = [val_background,
     val_layer,
     val_casing,
     val_insideCasing,
     layer_center,
     layer_thickness,
     casing_radius,
     casing_thickness,
     casing_bottom,
     casing_top
]
casing_a(mDict)
casing_b(mDict)
deriv(m)
layer_cont(mDict)
mDict(m)
nP
shape
class SimPEG.Maps.ParametricCircleMap(mesh, logSigma=True)

Bases: SimPEG.Maps.IdentityMap

Parameterize the model space using a circle in a wholespace.

\[\sigma(m) = \sigma_1 + (\sigma_2 - \sigma_1)\left( \arctan\left(100*\sqrt{(\vec{x}-x_0)^2 + (\vec{y}-y_0)}-r \right) \pi^{-1} + 0.5\right)\]

Define the model as:

\[m = [\sigma_1, \sigma_2, x_0, y_0, r]\]
deriv(m, v=None)
nP
slope = 0.1
class SimPEG.Maps.ParametricLayer(mesh, **kwargs)

Bases: SimPEG.Maps.BaseParametric

Parametric Layer Space

m = [
    val_background,
    val_layer,
    layer_center,
    layer_thickness
]

Required

Parameters:mesh (discretize.BaseMesh.BaseMesh) – SimPEG Mesh, 2D or 3D

Optional

Parameters:
  • slopeFact (float) – arctan slope factor - divided by the minimum h spacing to give the slope of the arctan functions
  • slope (float) – slope of the arctan function
  • indActive (numpy.ndarray) – bool vector with
deriv(m)
layer_cont(mDict)
mDict(m)
nP
shape
class SimPEG.Maps.ParametricPolyMap(mesh, order, logSigma=True, normal=u'X', actInd=None)

Bases: SimPEG.Maps.IdentityMap

PolyMap

Parameterize the model space using a polynomials in a wholespace.

\[y = \mathbf{V} c\]

Define the model as:

\[m = [\sigma_1, \sigma_2, c]\]

Can take in an actInd vector to account for topography.

deriv(m, v=None)
nP
shape
slope = 10000.0
class SimPEG.Maps.ParametricSplineMap(mesh, pts, ptsv=None, order=3, logSigma=True, normal=u'X')

Bases: SimPEG.Maps.IdentityMap

SplineMap

Parameterize the boundary of two geological units using a spline interpolation

\[g = f(x)-y\]

Define the model as:

\[m = [\sigma_1, \sigma_2, y]\]
deriv(m, v=None)
nP
slope = 10000.0
class SimPEG.Maps.PolyMap(mesh, order, logSigma=True, normal=u'X', actInd=None)

Bases: SimPEG.Maps.ParametricPolyMap

PolyMap is depreciated. Use ParametricSplineMap instead

class SimPEG.Maps.Projection(nP, index, **kwargs)

Bases: SimPEG.Maps.IdentityMap

A map to rearrange / select parameters

Parameters:
  • nP (int) – number of model parameters
  • index (numpy.array) – indices to select
deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
shape

Shape of the matrix operation (number of indices x nP)

class SimPEG.Maps.ReciprocalMap(mesh=None, nP=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Reciprocal mapping. For example, electrical resistivity and conductivity.

\[\rho = \frac{1}{\sigma}\]
deriv(m, v=None)
inverse(D)
class SimPEG.Maps.SelfConsistentEffectiveMedium(mesh=None, nP=None, sigstart=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap, properties.base.base.HasProperties

Two phase self-consistent effective medium theory mapping for ellipsoidal inclusions. The inversion model is the concentration (volume fraction) of the phase 2 material.

The inversion model is \(\varphi\). We solve for \(\sigma\) given \(\sigma_0\), \(\sigma_1\) and \(\varphi\) . Each of the following are implicit expressions of the effective conductivity. They are solved using a fixed point iteration.

Spherical Inclusions

If the shape of the inclusions are spheres, we use

\[\sum_{j=1}^N (\sigma^* - \sigma_j)R^{j} = 0\]

where \(j=[1,N]\) is the each material phase, and N is the number of phases. Currently, the implementation is only set up for 2 phase materials, so we solve

\[(1-\varphi)(\sigma - \sigma_0)R^{(0)} + \varphi(\sigma - \sigma_1)R^{(1)} = 0.\]

Where \(R^{(j)}\) is given by

\[R^{(j)} = \left[1 + \frac{1}{3}\frac{\sigma_j - \sigma}{\sigma} \right]^{-1}.\]

Ellipsoids

If the inclusions are aligned ellipsoids, we solve

\[\sum_{j=1}^N \varphi_j (\Sigma^* - \sigma_j\mathbf{I}) \mathbf{R}^{j, *} = 0\]

where

\[\mathbf{R}^{(j, *)} = \left[ \mathbf{I} + \mathbf{A}_j {\Sigma^{*}}^{-1}(\sigma_j \mathbf{I} - \Sigma^*) \right]^{-1}\]

and the depolarization tensor \(\mathbf{A}_j\) is given by

\[\begin{split}\mathbf{A}^* = \left[\begin{array}{ccc} Q & 0 & 0 \\ 0 & Q & 0 \\ 0 & 0 & 1-2Q \end{array}\right]\end{split}\]

for a spheroid aligned along the z-axis. For an oblate spheroid (\(\alpha < 1\), pancake-like)

\[Q = \frac{1}{2}\left( 1 + \frac{1}{\alpha^2 - 1} \left[ 1 - \frac{1}{\chi}\tan^{-1}(\chi) \right] \right)\]

where

\[\chi = \sqrt{\frac{1}{\alpha^2} - 1}\]

For reference, see Torquato (2002), Random Heterogeneous Materials

Required Properties:

  • alpha0 (Float): aspect ratio of the phase-0 ellipsoids, a float, Default: 1.0
  • alpha1 (Float): aspect ratio of the phase-1 ellipsoids, a float, Default: 1.0
  • maxIter (Integer): maximum number of iterations for the fixed point iteration calculation, an integer, Default: 50
  • orientation0 (Vector3): orientation of the phase-0 inclusions, a 3D Vector of <type ‘float’> with shape (3), Default: Z
  • orientation1 (Vector3): orientation of the phase-1 inclusions, a 3D Vector of <type ‘float’> with shape (3), Default: Z
  • random (Boolean): are the inclusions randomly oriented (True) or preferentially aligned (False)?, a boolean, Default: True
  • rel_tol (Float): relative tolerance for convergence for the fixed-point iteration, a float, Default: 0.001
  • sigma0 (Float): physical property value for phase-0 material, a float in range [0.0, inf]
  • sigma1 (Float): physical property value for phase-1 material, a float in range [0.0, inf]
alpha0

alpha0 (Float) – aspect ratio of the phase-0 ellipsoids, a float, Default – 1.0

alpha1

alpha1 (Float) – aspect ratio of the phase-1 ellipsoids, a float, Default – 1.0

deriv(m)

Derivative of the effective conductivity with respect to the volume fraction of phase 2 material

getA(alpha, orientation)

Depolarization tensor

getQ(alpha)

Geometric factor in the depolarization tensor

getR(sj, se, alpha, orientation=None)

Electric field concentration tensor

getdR(sj, se, alpha, orientation=None)

Derivative of the electric field concentration tensor with respect to the concentration of the second phase material.

hashin_shtrikman_bounds(phi1)

Hashin Shtrikman bounds

See Torquato, 2002

hashin_shtrikman_bounds_anisotropic(phi1)

Hashin Shtrikman bounds for anisotropic media

See Torquato, 2002

inverse(sige)

Compute the concentration given the effective conductivity

maxIter

maxIter (Integer) – maximum number of iterations for the fixed point iteration calculation, an integer, Default – 50

orientation0

orientation0 (Vector3) – orientation of the phase-0 inclusions, a 3D Vector of <type ‘float’> with shape (3), Default – Z

orientation1

orientation1 (Vector3) – orientation of the phase-1 inclusions, a 3D Vector of <type ‘float’> with shape (3), Default – Z

random

random (Boolean) – are the inclusions randomly oriented (True) or preferentially aligned (False)?, a boolean, Default – True

rel_tol

rel_tol (Float) – relative tolerance for convergence for the fixed-point iteration, a float, Default – 0.001

sigma0

sigma0 (Float) – physical property value for phase-0 material, a float in range [0.0, inf]

sigma1

sigma1 (Float) – physical property value for phase-1 material, a float in range [0.0, inf]

sigstart

first guess for sigma

tol

absolute tolerance for the convergence of the fixed point iteration calc

wiener_bounds(phi1)

Define Wenner Conductivity Bounds

See Torquato, 2002

class SimPEG.Maps.SplineMap(mesh, pts, ptsv=None, order=3, logSigma=True, normal=u'X')

Bases: SimPEG.Maps.ParametricSplineMap

SplineMap is depreciated. Use ParametricSplineMap instead

class SimPEG.Maps.Surject2Dto3D(mesh, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Map2Dto3D

Given a 2D vector, this will extend to the full 3D model space.

deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

normal = u'Y'
class SimPEG.Maps.SurjectFull(mesh, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Given a scalar, the SurjectFull maps the value to the full model space.

deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:numpy.array
Returns:derivative of transformed model
nP
class SimPEG.Maps.SurjectVertical1D(mesh, **kwargs)

Bases: SimPEG.Maps.IdentityMap

SurjectVertical1DMap

Given a 1D vector through the last dimension of the mesh, this will extend to the full model space.

deriv(m, v=None)
Parameters:m (numpy.array) – model
Return type:scipy.sparse.csr_matrix
Returns:derivative of transformed model
nP

Number of model properties.

The number of cells in the last dimension of the mesh.

class SimPEG.Maps.Vertical1DMap(mesh, **kwargs)

Bases: SimPEG.Maps.SurjectVertical1D

Vertical1DMap is depreciated. Use SurjectVertical1D instead

class SimPEG.Maps.Weighting(mesh=None, nP=None, weights=None, **kwargs)

Bases: SimPEG.Maps.IdentityMap

Model weight parameters.

P
deriv(m, v=None)
inverse(D)
shape
class SimPEG.Maps.Wires(*args)

Bases: object

nP