The fatiando package has been deprecated. Please check out the new tools in the Fatiando a Terra website: www.fatiando.org

Forward modeling with 3D polygonal prisms (fatiando.gravmag.polyprism)

Calculate the potential fields of the 3D prism with polygonal crossection using the formula of Plouff (1976).

Gravity

First and second derivatives of the gravitational potential:

Magnetic

There are functions to calculate the total-field anomaly and the 3 components of magnetic induction:

Auxiliary Functions

Calculates the second derivatives of the function

\[\phi(x,y,z) = \int\int\int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

with respect to the variables \(x\), \(y\), and \(z\). In this equation,

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}\]

and \(\nu\), \(\eta\), \(\zeta\) are the Cartesian coordinates of an element inside the volume of a 3D prism with polygonal crossection. These second derivatives are used to calculate the total field anomaly and the gravity gradient tensor components produced by a 3D prism with polygonal crossection.

References

Plouff, D. , 1976, Gravity and magnetic fields of polygonal prisms and applications to magnetic terrain corrections, Geophysics, 41(4), 727-741.


fatiando.gravmag.polyprism.bx(xp, yp, zp, prisms)[source]

Calculates the x component of the magnetic induction produced by 3D prisms with polygonal crosssection.

Note

Input units are SI. Output is in nT

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates where the anomaly will be calculated

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the total field anomaly. Prisms without the physical property 'magnetization' will be ignored. The 'magnetization' must be a vector.

Returns:

  • bx: array

    The x component of the magnetic induction

fatiando.gravmag.polyprism.by(xp, yp, zp, prisms)[source]

Calculates the y component of the magnetic induction produced by 3D prisms with polygonal crosssection.

Note

Input units are SI. Output is in nT

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates where the anomaly will be calculated

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the total field anomaly. Prisms without the physical property 'magnetization' will be ignored. The 'magnetization' must be a vector.

Returns:

  • by: array

    The y component of the magnetic induction

fatiando.gravmag.polyprism.bz(xp, yp, zp, prisms)[source]

Calculates the z component of the magnetic induction produced by 3D prisms with polygonal crosssection.

Note

Input units are SI. Output is in nT

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates where the anomaly will be calculated

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the total field anomaly. Prisms without the physical property 'magnetization' will be ignored. The 'magnetization' must be a vector.

Returns:

  • bz: array

    The z component of the magnetic induction

fatiando.gravmag.polyprism.gxx(xp, yp, zp, prisms)[source]

Calculates the \(g_{xx}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gxy(xp, yp, zp, prisms)[source]

Calculates the \(g_{xy}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gxz(xp, yp, zp, prisms)[source]

Calculates the \(g_{xz}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gyy(xp, yp, zp, prisms)[source]

Calculates the \(g_{yy}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gyz(xp, yp, zp, prisms)[source]

Calculates the \(g_{yz}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gz(xp, yp, zp, prisms)[source]

Calculates the \(g_{z}\) gravity acceleration component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in mGal!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.gzz(xp, yp, zp, prisms)[source]

Calculates the \(g_{zz}\) gravity gradient tensor component.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input values in SI units and output in Eotvos!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the field. Prisms must have the physical property 'density' will be ignored.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelxx(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial x^2},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelxy(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial x \partial y},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelxz(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial x \partial z},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelyy(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial y^2},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelyz(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial y \partial z},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.kernelzz(xp, yp, zp, prism)[source]

Calculates the function

\[\frac{\partial^2 \phi(x,y,z)}{\partial z^2},\]

where

\[\phi(x,y,z) = \int \int \int \frac{1}{r} \mathrm{d}\nu \mathrm{d}\eta \mathrm{d}\zeta\]

and

\[r = \sqrt{(x - \nu)^2 + (y - \eta)^2 + (z - \zeta)^2}.\]

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

All input and output values in SI!

Parameters:

  • xp, yp, zp
    : arrays

    The x, y, and z coordinates of the computation points.

  • prisms
    : object of fatiando.mesher.PolygonalPrism

    The model used to calculate the function.

Returns:

  • res
    : array

    The effect calculated on the computation points.

fatiando.gravmag.polyprism.tf(xp, yp, zp, prisms, inc, dec, pmag=None)[source]

Calculate the total-field anomaly of polygonal prisms.

Note

The coordinate system of the input parameters is to be x -> North, y -> East and z -> Down.

Note

Input units are SI. Output is in nT

Parameters:

  • xp, yp, zp
    : arrays

    Arrays with the x, y, and z coordinates of the computation points.

  • prisms
    : list of fatiando.mesher.PolygonalPrism

    The model used to calculate the total field anomaly. Prisms without the physical property 'magnetization' will be ignored.

  • inc
    : float

    The inclination of the regional field (in degrees)

  • dec
    : float

    The declination of the regional field (in degrees)

  • pmag
    : [mx, my, mz] or None

    A magnetization vector. If not None, will use this value instead of the 'magnetization' property of the prisms. Use this, e.g., for sensitivity matrix building.

Returns:

  • res
    : array

    The field calculated on xp, yp, zp