magali.NonlinearMagneticDipoleBz#

class magali.NonlinearMagneticDipoleBz(initial_location, max_iter=100, tol=0.0001, alpha_init=1, alpha_scale=10.0)[source]#

Estimate the position and magnetic dipole moment vector from Bz measurements.

Uses the Bz component of the magnetic field to estimate both the position and the moment of a magnetic dipole through a nonlinear inversion based on the Levenberg-Marquardt algorithm. Returns the location and dipole moment vector that best fit the data in a least-squares sense. Requires an initial guess for the dipole location.

Parameters:
initial_locationtuple of float

Initial guess for the coordinates (x, y, z) of the dipole location, in µm.

max_iterint

Maximum number of iterations for both the outer and inner loops of the nonlinear inversion.

tolfloat

Convergence tolerance for the relative change in misfit between iterations.

alpha_initfloat

Initial damping parameter for the Levenberg-Marquardt algorithm.

alpha_scalefloat

Multiplicative factor used to increase or decrease the damping parameter during the optimization.

References

[Souza-Junior2024]

Attributes:
location_1d-array

Estimated location (x, y, z) of the dipole, in µm.

dipole_moment_1d-array

Estimated dipole moment vector (mx, my, mz), in A.m².

misfit_list of float

Norm of the residual vector at each outer iteration of the inversion.

Methods

fit(coordinates, data)

Fit the nonlinear magnetic dipole model to Bz data.

jacobian(coordinates, location, moment, jacobian)

Compute the Jacobian matrix for the linear point dipole model.

predict(coordinates)

Predict the Bz magnetic field from the estimated dipole parameters.

Method documentation#

NonlinearMagneticDipoleBz.fit(coordinates, data)[source]#

Fit the nonlinear magnetic dipole model to Bz data.

Performs nonlinear inversion using the Levenberg-Marquardt method to estimate both the dipole location and its magnetic moment. The method alternates between a nonlinear update of the dipole location (inner loop) and a linear least-squares estimate of the dipole moment (outer loop). The Jacobian matrix with respect to the location is computed numerically using JIT-accelerated code.

The Jacobian matrix used in the nonlinear step contains partial derivatives of the Bz field with respect to the dipole location: \(\frac{\partial B_z}{\partial x_0}\), \(\frac{\partial B_z}{\partial y_0}\), and \(\frac{\partial B_z}{\partial z_0}\). These are computed assuming a fixed moment vector.

At each inner iteration:

  • The forward field is computed using the trial location and fixed moment.

  • A trial update is accepted if it reduces the data misfit.

  • The damping parameter (alpha) is adapted based on success/failure.

At each outer iteration:

  • A new linear estimate of the dipole moment is computed for the current location.

  • Convergence is assessed based on relative reduction in residual norm.

Parameters:
coordinatestuple of array_like

Arrays with the x, y, z coordinates of the observation points. The arrays can have any shape as long as they all have the same shape.

dataarray_like

Observed Bz component of the magnetic field (in nT) at the observation points. Must have the same shape as the coordinate arrays.

Returns:
selfobject

This instance with updated location_ and dipole_moment_.

Notes

Internally uses:

  • jacobian_nonlinear_jit: JIT-compiled function that fills the Jacobian matrix \(\frac{\partial B_z}{\partial (x_0, y_0, z_0)}\) for a fixed moment.

  • dipole_bz: forward model for the Bz field of a dipole at given coordinates.

  • MagneticMomentBz: linear inversion for estimating moment given a fixed location.

NonlinearMagneticDipoleBz.jacobian(coordinates, location, moment, jacobian)[source]#

Compute the Jacobian matrix for the linear point dipole model.

The Jacobian is a matrix with derivatives of the forward modeling function (the magnetic field of a dipole) with regard to the parameters (the 3 dipole moment components) for each data point.

Parameters:
coordinatestuple = (x, y, z)

Arrays with the x, y, and z coordinates of the observations points. The arrays can have any shape as long as they all have the same shape.

locationarray_like

Dipole location (x, y, z), in µm.

momentarray_like

Dipole moment vector (mx, my, mz), in A.m².

jacobian(n, 3) ndarray

Preallocated array to store the resulting Jacobian matrix, where n is the number of observation points.

Returns:
jacobian(n, 3) ndarray

The Jacobian matrix of Bz with respect to the dipole moment, in nT/(A·m²). Each row corresponds to an observation point, and each column to the partial derivative with respect to mx, my, and mz, respectively.

NonlinearMagneticDipoleBz.predict(coordinates)[source]#

Predict the Bz magnetic field from the estimated dipole parameters.

Uses the estimated dipole location and moment vector to compute the Bz component of the magnetic field at the given observation points.

Parameters:
coordinatestuple of array_like

Arrays with the x, y, z coordinates of the observation points, in µm. The arrays must have the same shape.

Returns:
predictedarray

Array with the predicted Bz values (in nT) at the observation points. Has the same shape as the input coordinate arrays.

Raises:
AttributeError

If fit has not been called yet, and location_ or dipole_moment_ are not set.

See also

fit

Estimate the dipole location and moment from Bz measurements.