Matlab supporting functions

Main calculation functions

MagFldParent()

functions.MagFldParent(planet, r_km, theta, phi, InternalFieldModel, ExternalFieldModel, magPhase_deg, SPHOUT, Nmaxin, ATTEN_SHEET)

Evaluate the magnetic field of the desired planet at specified locations according to the specified internal and external field models.

Note

The standard reference frames used for each body, except as noted below, are IAU_BODY. The only exceptions to this rule are:

  • Uranus uses the ULS frame: System III coordinates with +z along the spin pole, opposite the direction of IAU_URANUS, and with longitudes defined arbitrarily based on the Voyager 2 flyby.

  • Neptune uses the NLS frame: The Neptune Longitude System derived from magnetic data. IAU_NEPTUNE is an outdated implementation of the IAU frame definition (as of the 2009 report). The 2015 report redefined the IAU frame to be a System II frame, rotating at a different rate than that inferred from magnetic measurements. The IAU_NEPTUNE frame rotates at the rate inferred from magnetic measurements, but is oriented differently. See Implemented frames and coordinates for more information.

All standard frames rotate with the body (planet or moon).

Parameters:
  • planet (char, 1xC) – Name of the planet for which to evaluate the field.

  • r_km (double, 1xN) – Radius of evaluation points in km from planet barycenter.

  • theta (double, 1xN) – Planetocentric colatitude of evaluation points in radians. Must be the same length as r_km.

  • phi (double, 1xN) – Planetocentric east longitude of evaluation points in radians. Must be the same length as r_km.

  • InternalFieldModel (char, 1xD) – Code name for planetary intrinsic field model to evaluate. Refer to GetModelOpts() for a listing of valid options. Also accepts the special option 'None', which skips evaluation of the intrinsic field model.

  • ExternalFieldModel (char, 1xE) – Code name for current sheet field model to evaluate. Refer to GetModelOpts() for a listing of valid options. Also accepts the special option 'None', which skips evaluation of the current sheet field model.

  • magPhase_deg (double, default=0) – Arbitrary offset in degrees by which to rotate the magnetospheric field evaluation.

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • Nmaxin (int, default=Inf) – Maximum degree to which to limit intrinsic field models. Only has an effect if the value passed is less than the lesser of the model maximum degree and the maximum implemented in spherical harmonic calculations (currently 10).

  • ATTEN_SHEET (bool, default=0) – Whether to attenuate current sheet models using the analytical approximations of Connerney et al. (1981), which are C1981, C2020, and Cassini 11, beyond \(50R_P\).

Returns:

  • Bvec_nT (double, 3xN) – Magnetic field vector in standard coordinates (typically IAU or spherical coordinates if SPHOUT is true) at evaluation points in nT. Output rows are x, y, z respectively for cartesian or r, theta, phi for spherical.

  • Mdip_nT (double, 1x3) – Dipole magnetic moment in coordinates matching Bvec_nT, as surface-equivalent nT, i.e. this vector times the planet volume yields the magnetic moment in SI units.

  • Odip_km (double, 1x3) – Dipole magnetic moment offset in km from planet barycenter in coordinates matching Bvec_nT.

MpauseFld()

functions.MpauseFld(nSW_pcc, vSW_kms, ets, xyz_km, Mdip_nT, Odip_km, parentName, S3coords, MPmodel, coeffPath, SPHOUT, Nmax)

Evaluate the magnetic field contribution from a magnetopause current model at specific locations and times.

Parameters:
  • nSW_pcc (double, 1xN) – Ion density in the solar wind in particles per cubic centimeter.

  • vSW_kms (double, 1xN) – Speed of incident solar wind ions in km/s.

  • ets (double, 1xN) – Ephemeris times of each measurement in TDB seconds relative to J2000.

  • xyz_km (double, 3xN) – Cartesian coordinates of measurement locations in S3coords frame in km.

  • Mdip_nT (double, 1x3) – Dipole magnetic moment in standard cartesian coordinates, as surface-equivalent nT.

  • Odip_km (double, 1x3) – Dipole magnetic moment offset in km from planet barycenter in standard cartesian coordinates.

  • parentName (char, 1xC) – Name of parent planet the evaluated model applies to.

  • S3coords (char, 1xD) – Standard coordinate frame used for evaluation of magnetic fields.

  • MPmodel (char, 1xE) – Code name for magnetopause current field model to apply. See GetModelOpts() for more details and available options.

  • coeffPath (char, 1xF, default='modelCoeffs') – Directory containing model coefficients files.

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • Nmax (int, default=10) – Maximum degree to which to limit spherical harmonic models. Only has an effect if the value passed is less than the lesser of the model maximum degree.

Returns:

  • mpBvecOut (double, 3xN) – Contribution to magnetic field at measurement positions from magnetopause currents.

  • OUTSIDE_MP (bool, 1xN) – Whether each point is outside the magnetopause, and so should have the net magnetospheric field zeroed out.

MagFldParentSingle()

functions.MagFldParentSingle(g, h, r_km, theta, phi, PlanetEqRadius, InternalFieldModel, ExternalFieldModel, magPhase_deg, Nmax, NmaxExt, SPHOUT, ATTEN_SHEET)

Evaluate the magnetic field of the desired planet at specified locations according to the specified internal and external field models.

See MagFldPlanet for additional notes.

Parameters:
  • g (double, (Nmax)x(Nmax+1)) – Internal field g coefficients in Schmidt normalization in G.

  • h (double, (Nmax)x(Nmax+1)) – Internal field h coefficients in Schmidt normalization in G.

  • r_km (double) – Radius of evaluation point in km from planet barycenter.

  • theta (double) – Planetocentric colatitude of evaluation points in radians.

  • phi (double) – Planetocentric east longitude of evaluation point in radians.

  • PlanetEqRadius (double) – Equatorial radius of planet in km associated with input g and h values.

  • InternalFieldModel (char, 1xD) – Code name for planetary intrinsic field model to evaluate. Refer to GetModelOpts() for a listing of valid options. Also accepts the special option 'None', which skips evaluation of the intrinsic field model.

  • ExternalFieldModel (char, 1xE) – Code name for current sheet field model to evaluate. Refer to GetModelOpts() for a listing of valid options. Also accepts the special option 'None', which skips evaluation of the current sheet field model.

  • magPhase_deg (double, default=0) – Arbitrary offset in degrees by which to rotate the magnetospheric field evaluation.

  • Nmax (int, default=Inf) – Maximum degree to which to limit intrinsic field models. Only has an effect if the value passed is less than the lesser of the model maximum degree and the maximum implemented in spherical harmonic calculations (currently 10).

  • NmaxExt (int, default=Inf) – Maximum degree to which to limit external field models. Only has an effect if the value passed is less than the lesser of the model maximum degree and the maximum implemented in external field calculations (currently 1).

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • ATTEN_SHEET (bool, default=0) – Whether to attenuate current sheet models using the analytical approximations of Connerney et al. (1981), which are C1981, C2020, and Cassini 11, beyond \(50R_P\).

Returns:

Bvec_nT – Magnetic field vector in standard coordinates (typically IAU or spherical coordinates if SPHOUT is true) at evaluation points in nT. Output rows are x, y, z respectively for cartesian or r, theta, phi for spherical.

Return type:

double, 3x1

LLSdecomposition()

functions.LLSdecomposition(moonName, parentName, ets, Bvec, magModelDescrip, SPHOUT, PLOT_DIAGNOSTIC, COMPARE_SEUFERT, COMPARE_PHIO, LIVE_PLOTS, figDir, figXtn)

Decomposes the input magnetic field vector time series into complex excitation moments using linear least-squares (LLS) decomposition.

The methodology applied to invert the input magnetic field vector time series Bvec (\(\mathbf{B}_i\)) for its complex excitation moments is described in detail in the PlanetMag paper, Styczinski and Cochrane (2024), DOI TBD. Each time in the input time series \(t_i\) corresponds to a row in the sampling matrix \(X_{ik}\), and each frequency in the list returned by GetExcitations() corresponds to a column in \(X_{ik}\). The entries in the matrix are sinusoidal, as we are decomposing the input time series \(\mathbf{B}_i\) into complex Fourier components. The sampling matrix is

\[\begin{split}X_{ik} = \begin{bmatrix} \cos(\omega_1 t_1) & \sin(\omega_1 t_1) & \cos(\omega_2 t_1) & \sin(\omega_2 t_1) & \dots & \cos(\omega_F t_1) & \sin(\omega_F t_1) & 1 \\ \cos(\omega_1 t_2) & \sin(\omega_1 t_2) & \cos(\omega_2 t_2) & \sin(\omega_2 t_2) & \dots & \cos(\omega_F t_2) & \sin(\omega_F t_2) & 1 \\ \vdots & \vdots & \vdots & \vdots & \ddots & \vdots & \vdots & \vdots \\ \cos(\omega_1 t_N) & \sin(\omega_1 t_N) & \cos(\omega_2 t_N) & \sin(\omega_2 t_N) & \dots & \cos(\omega_F t_N) & \sin(\omega_F t_N) & 1 \\ \end{bmatrix}\end{split}\]

where \(N\) is the number of points in the time series and \(F\) is the number of frequencies for which to invert the excitation moments. The eigenvectors of \(X_{ik}\) are the columns of the weight matrix \(W_{kk'}\), such that

\[W = \left( X_{ik}^T X_{ik} \right)^{-1}\]

and the excitation moments are

\[B^e_{j,k} = (\mathbf{B}_i\cdot\mathbf{e}_j) X_{ik}W,\]

with each element \(B_{j,k}^e\) corresponding to the real (cos) or imaginary (sin) part of the excitation moment for frequency \(f_k\). The complex excitation moments are constructed by combining the real and imaginary parts for each \(f_k\):

\[\mathbf{B}^e_k = \sum_j \left( B^e_{j,k,\mathrm{Re}} + \mathrm{i}B^e_{j,k,\mathrm{Im}} \right) \mathbf{e}_j\]

and the reproduced time series is

\[\mathrm{Re}\{\mathbf{B}_\mathrm{exc}(t)\} = \sum_k \left[\mathbf{B}_{k,\mathrm{Re}}^e \cos(\omega_k t) + \mathbf{B}_{k,\mathrm{Im}}^e \sin(\omega_k t) \right].\]

From the input time series,

\[\mathrm{Re}\{\mathbf{B}_\mathrm{exc}(t_i)\} = \sum_{j,k} X_{ik}(B_{i,j}^e)^T \mathbf{e}_j.\]
Parameters:
  • moonName (char, 1xC) – Name of moon for which to evaluate excitation moments.

  • parentName (char, 1xD) – Name of parent planet for the moon for which to evaluate excitation moments.

  • ets (double, 1xN) – Ephemeris times associated with input time series magnetic field vectors in TDB seconds relative to J2000.

  • Bvec (double, 3xN) – Magnetic field vectors at each measurement time.

  • magModelDescrip (char, 1xE) – Text description of the magnetic field model that was evalauted for the input time series.

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • PLOT_DIAGNOSTIC (bool, default=1) – Whether to plot comparisons between the input and reproduced time series.

  • COMPARE_SEUFERT (bool, default=1) – Whether to plot coordinate directions aligned to the axes presented in Seufert et al. (2011) https://doi.org/10.1016/j.icarus.2011.03.017. Only has an effect when SPHOUT is true, because Seufert et al. used spherical coordinates for evaluating the excitation spectra of Jupiter’s moons.

  • COMPARE_PHIO (bool, default=1) – Whether to plot components in the same orientation as PhiO axes, as in past work, and label the axes with the comparisons (true), or just plot IAU frame axes without comparisons (false).

  • LIVE_PLOTS (bool, default=0) – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • figDir (char, 1xF, default='figures') – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • figXtn (char, 1xG, default='pdf') – Extension to use for figures, which determines the file type.

Returns:

BD

Contains fields:

  • f (double, 1xP) – Frequencies of estimated excitation moments in Hz.

  • fNames (string, 1xP`) – Descriptive names for each excitation.

  • BexcVec1i (double, 1xP) – Real part (i for “in-phase”) of excitation moments for vector component 1 (x or r) in nT.

  • BexcVec1q (double, 1xP) – Imaginary part (q for “quadrature”) of excitation moments for vector component 1 in nT.

  • BexcVec1o (double) – Vector component 1 of uniform background field in nT.

  • BexcVec2i (double, 1xP) – Real part (i for “in-phase”) of excitation moments for vector component 2 (y or theta) in nT.

  • BexcVec2q (double, 1xP) – Imaginary part (q for “quadrature”) of excitation moments for vector component 2 in nT.

  • BexcVec2o (double) – Vector component 2 of uniform background field in nT.

  • BexcVec3i (double, 1xP) – Real part (i for “in-phase”) of excitation moments for vector component 3 (z or phi) in nT.

  • BexcVec3q (double, 1xP) – Imaginary part (q for “quadrature”) of excitation moments for vector component 3 in nT.

  • BexcVec3o (double) – Vector component 3 of uniform background field in nT.

  • rmse (double) – Normalized root-mean-squared error calculated from comparing input and reproduced time series data.

  • Rsquared (double) – \(R^2\) value for the reproduced time series.

Return type:

struct

LegendreS()

functions.LegendreS(n, m, theta, UNNORM)

Evaluate Schmidt semi-normalized associated Legendre functions for a single n,m.

Returns the value of Schmidt semi-normalized associated Legendre functions \(S_n^m\) of degree n and order m at a list of colatitudes theta. The Condon–Shortley phase is omitted. The Schmidt normalization is

\[S_n^m(\theta) = \sqrt{(2 - \delta_{m0})\frac{(n-m)!}{(n+m)!}} P_n^m(\cos\theta)\]

where \(P_n^m\) are the unnormalized associated Legendre functions without the Condon–Shortley phase and \(\delta_{m,m'}\) is the Kronecker delta function, which is 0 unless \(m = m'\).

Parameters:
  • n (int) – Degree n and order m for which to evaluate Legendre functions.

  • m (int) – Degree n and order m for which to evaluate Legendre functions.

  • theta (double, 1xN) – Colatitudes in radians at which to evaluate Legendre functions.

  • UNNORM (bool, default=0) – Whether to return un-normalized associated Legendre functions \(P_n^m(\cos\theta)\). As with the default, semi-normalized \(S_n^m\), the Condon–Shortley phase is omitted.

Returns:

Snm – Schmidt semi-normalized associated Legendre function evaluation at given colatitudes for the specified n,m.

Return type:

double, 1xN

dLegendreS()

functions.dLegendreS(n, m, theta, UNNORM)

Evaluate the derivative of Schmidt semi-normalized associated Legendre functions for a single n,m.

Returns the value of the \(\theta\) derivative of Schmidt semi-normalized associated Legendre functions \(\frac{dS_n^m}{d\theta}\) of degree n and order m at a list of colatitudes theta. The Condon–Shortley phase is omitted. See LegendreS() for more information on the normalization.

Parameters:
  • n (int) – Degree n and order m for which to evaluate Legendre function derivatives.

  • m (int) – Degree n and order m for which to evaluate Legendre function derivatives.

  • theta (double, 1xN) – Colatitudes in radians at which to evaluate Legendre function derivatives.

  • UNNORM (bool, default=0) – Whether to return un-normalized associated Legendre function derivatives \(\frac{dP_n^m(\cos\theta)}{d\theta}\). As with the default, the Condon–Shortley phase is omitted.

Returns:

dSnm – Schmidt semi-normalized associated Legendre function derivative evaluation at given colatitudes for the specified n,m.

Return type:

double, 1xN

Coordinate transformation functions

Bcyl2Bxyz()

functions.coordinates.Bcyl2Bxyz(Brho, Bphi, Bzin, phi)

Convert cylindrical (rho, phi, z) vector to cartesian (x, y, z).

Convert vector components aligned to cylindrical coordinates into vector components aligned to cartesian axes. Source: Arfken, Weber, Harris, Mathematical Methods for Physicists, 7th ed, pg. 197 for the unit vectors.

Parameters:
  • Brho (double, 1xN) – x-aligned component of vectors to be converted.

  • Bphi (double, 1xN) – Azimuthal (\(\hat{\phi}\)) component of vectors to be converted.

  • Bzin (double, 1xN) – Axial (\(\hat{z}\)) component of vectors to be converted.

  • phi (double, 1xN) – East longitude for each location associated with each (Brho, Bphi, Bzin) vector in radians.

Returns:

Bx, By, Bz – x-, y-, and z-aligned components of converted vectors.

Return type:

double, 1xN

Bsph2Bxyz()

functions.coordinates.Bsph2Bxyz(Br, Bth, Bphi, theta, phi)

Convert spherical (r, theta, phi) vector to cartesian (x, y, z) vector.

Convert vector components aligned to spherical coordinates into vector components aligned to cartesian axes. Source: Arfken, Weber, Harris, Mathematical Methods for Physicists, 7th ed, pg. 199 for the unit vectors.

Parameters:
  • Br (double, 1xN) – Radial (\(\hat{r}\)) component of vectors to be converted.

  • Bth (double, 1xN) – Colatitudinal (\(\hat{\theta}\)) component of vectors to be converted.

  • Bphi (double, 1xN) – Azimuthal (\(\hat{\phi}\)) component of vectors to be converted.

  • theta (double, 1xN) – Colatitude for each location associated with each (Br, Bth, Bphi) vector in radians.

  • phi (double, 1xN) – East longitude for each location associated with each (Br, Bth, Bphi) vector in radians.

Returns:

Bx, By, Bz – x-, y-, and z-aligned components of converted vectors.

Return type:

double, 1xN

Bxyz2Bcyl()

functions.coordinates.Bxyz2Bcyl(Bx, By, Bz, phi)

Convert cartesian (x, y, z) vector to cylindrical (rho, phi, z) vector.

Convert vector components aligned to cartesian coordinates into vector components aligned to cylindrical axes. Source: Arfken, Weber, Harris, Mathematical Methods for Physicists, 7th ed, pg. 197 for the unit vectors.

Parameters:
  • Bx (double, 1xN) – x-aligned component of vectors to be converted.

  • By (double, 1xN) – y-aligned component of vectors to be converted.

  • Bz (double, 1xN) – z-aligned component of vectors to be converted.

  • phi (double, 1xN) – East longitude for each location associated with each (Bx, By, Bz) vector in radians.

Returns:

Brho, Bphi, Bzout – Equatorial projection of radial (\(\hat{\rho}\)), azimuthal (\(\hat{\phi}\)), and axial (\(\hat{z}\)) components of converted vectors.

Return type:

double, 1xN

Bxyz2Bsph()

functions.coordinates.Bxyz2Bsph(Bx, By, Bz, theta, phi)

Convert cartesian (x, y, z) vector to spherical (r, theta, phi) vector.

Convert vector components aligned to cartesian axes into vector components aligned to spherical coordinates. Source: Arfken, Weber, Harris, Mathematical Methods for Physicists, 7th ed, pg. 199 for the unit vectors.

Parameters:
  • Bx (double, 1xN) – x-aligned component of vectors to be converted.

  • By (double, 1xN) – y-aligned component of vectors to be converted.

  • Bz (double, 1xN) – z-aligned component of vectors to be converted.

  • theta (double, 1xN) – Colatitude for each location associated with each (Bx, By, Bz) vector in radians.

  • phi (double, 1xN) – East longitude for each location associated with each (Bx, By, Bz) vector in radians.

Returns:

Br, Bth, Bphi – Radial (\(\hat{r}\)), colatitudinal (\(\hat{\theta}\)), and azimuthal (\(\hat{\phi}\)) components of converted vectors.

Return type:

double, 1xN

cyl2xyz()

functions.coordinates.cyl2xyz(rho, phi, zin)

Convert position vectors from cylindrical coordinates to cartesian.

Parameters:
  • rho (double, 1xN) – rho coordinate to be converted.

  • phi (double, 1xN) – phi coordinate to be converted in radians.

  • zin (double, 1xN) – z coordinate to be converted.

Returns:

x, y, z – Cartesian coordinates of converted position vectors with the same units as rho and zin.

Return type:

double, 1xN

KS_BJSMtoBS3C()

functions.coordinates.KS_BJSMtoBS3C(BxJSM, ByJSM, BzJSM, ctimes, AS_CODED)

Rotate vectors aligned to Jupiter–Sun–magnetic coordinates to System III cartesian.

In Jupiter–Sun–magnetic (JSM) coordinates, \(\hat{x}\) points toward the Sun, the xz plane contains the magnetic dipole moment, and \(\hat{y}\) completes the right-handed set. This means \(\hat{z}\) is directed toward the \(x=0\) plane projection of the magnetic dipole moment. For Jupiter, \(\hat{y}\) nods roughly above and below the direction opposite the orbital velocity vector as the planet rotates.

Parameters:
  • BxJSM (double, 1xN) – x-aligned vector component to be converted in the JSM frame.

  • ByJSM (double, 1xN) – y-aligned vector component to be converted in the JSM frame.

  • BzJSM (double, 1xN) – z-aligned vector component to be converted in the JSM frame.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

BxS3, ByS3, BzS3 – x-, y-, and z-aligned components of converted vectors in System III frame.

Return type:

double, 1xN

KS_BpJSOtoBxyz()

functions.coordinates.KS_BpJSOtoBxyz(BxpJSO, BypJSO, BzpJSO, sphi)

Rotate cartesian vector aligned with pseudo-JSO coordinates to System III.

In psuedo-Jupiter–Sun–Orbital (pJSO) coordinates, \(\hat{z}\) is parallel to the jovian spin axis and so \(\hat{x}\) does not point toward the Sun as in JSO, but the Sun does lie in the xz plane. pJSO is equivalent to the Jupiter–De-Spun–Sun frame (JUNO_JSS) in the Juno frames kernel. See KS_S3CtoJSO() for a definition of the JSO frame.

Parameters:
  • BxpJSO (double, 1xN) – x-aligned vector component to be converted in the pseudo-JSO frame.

  • BypJSO (double, 1xN) – y-aligned vector component to be converted in the pseudo-JSO frame.

  • BzpJSO (double, 1xN) – z-aligned vector component to be converted in the pseudo-JSO frame.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in the System III frame in radians.

Returns:

Bx, By, Bz – x-, y-, and z-aligned components of converted vectors in System III frame.

Return type:

double, 1xN

KS_S3CtoJSM()

functions.coordinates.KS_S3CtoJSM(xS3C, yS3C, zS3c, ctimes, AS_CODED)

Rotate position vectors from System III cartesian to Jupiter–Sun–magnetic coordinates.

See KS_BJSMtoBS3C() for a definition of the JSM frame.

Parameters:
  • x (double, 1xN) – x coordinate to be converted in the System III frame.

  • y (double, 1xN) – y coordinate to be converted in the System III frame.

  • z (double, 1xN) – z coordinate to be converted in the System III frame.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

xJSM, yJSM, zJSM – Cartesian coordinates of converted vectors in JSM frame.

Return type:

double, 1xN

KS_S3CtoJSO()

functions.coordinates.KS_S3CtoJSO(xS3C, yS3C, zS3c, ctimes, AS_CODED)

Rotate position vectors from System III cartesian to Jupiter–Sun–Orbital coordinates.

In the Jupiter–Sun–Orbital (JSO) frame, \(\hat{x}\) points toward the Sun, \(\hat{y}\) points in the direction of the Sun’s velocity vector as observed from the planet center of mass (opposite the orbital velocity vector of Jupiter), and \(\hat{z}\) completes the right-handed set, roughly toward the ecliptic normal.

Parameters:
  • xS3C (double, 1xN) – x coordinate to be converted in the System III frame.

  • yS3C (double, 1xN) – y coordinate to be converted in the System III frame.

  • zS3C (double, 1xN) – z coordinate to be converted in the System III frame.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

xJSO, yJSO, zJSO – Cartesian coordinates of converted vectors in JSO frame.

Return type:

double, 1xN

KS_xyz2pJSO()

functions.coordinates.KS_xyz2pJSO(x, y, z, sphi)

Rotate position vectors from System III cartesian to pseudo-JSO coordinates.

See KS_BpJSOtoBxyz() for a definition of the pseudo-JSO frame.

Parameters:
  • x (double, 1xN) – x coordinate to be converted in the System III frame.

  • y (double, 1xN) – y coordinate to be converted in the System III frame.

  • z (double, 1xN) – z coordinate to be converted in the System III frame.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in the System III frame in radians.

Returns:

xpJSO, ypJSO, zpJSO – Cartesian coordinates of converted vectors in pseudo-JSO frame.

Return type:

double, 1xN

RotateVecSpice()

functions.coordinates.RotateVecSpice(Vx, Vy, Vz, ets, coordSource, coordDest)

Rotate vectors from one SPICE frame to another.

Parameters:
  • Vx (double, 1xN) – x-aligned component of vectors to be rotated from the source frame to the destination frame.

  • Vy (double, 1xN) – y-aligned component of vectors to be rotated from the source frame to the destination frame.

  • Vz (double, 1xN) – z-aligned component of vectors to be rotated from the source frame to the destination frame.

  • ets (double, 1xN) – Ephemeris times (ETs) in TDB seconds relative to J2000 for each vector to rotate.

  • coordSource (char, 1xC) – Source frame, with coordinate axes to which the input vectors are referenced.

  • coordDest (char, 1xC) – Destination frame, with coordinate axes to which the output vectors are referenced.

Returns:

VxMoon, VyMoon, VzMoon – Cartesian components of vectors in the destination frame for each input ET.

Return type:

double, 1xN

sph2xyz()

functions.coordinates.sph2xyz(r, theta, phi)

Convert position vectors from spherical coordinates to cartesian.

Parameters:
  • r (double, 1xN) – r coordinate to be converted.

  • theta (double, 1xN) – theta coordinate to be converted in radians.

  • phi (double, 1xN) – phi coordinate to be converted in radians.

Returns:

x, y, z – Cartesian coordinates of converted position vectors, with units matching r.

Return type:

double, 1xN

xyz2cyl()

functions.coordinates.xyz2cyl(x, y, zin)

Convert position vectors from cartesian coordinates to cylindrical.

Parameters:
  • x (double, 1xN) – x coordinate to be converted.

  • y (double, 1xN) – y coordinate to be converted.

  • zin (double, 1xN) – z coordinate to be converted.

Returns:

rho, phi, z – Cylindrical coordinates of converted position vectors. rho has the same units as the input coordinates. phi is in radians.

Return type:

double, 1xN

xyz2sph()

functions.coordinates.xyz2sph(x, y, z)

Convert position vectors from cartesian coordinates to spherical.

Parameters:
  • x (double, 1xN) – x coordinate to be converted.

  • y (double, 1xN) – y coordinate to be converted.

  • z (double, 1xN) – z coordinate to be converted.

Returns:

r, theta, phi – Spherical coordinates of converted position vectors. r has the same units as the input coordinates. theta and phi are in radians.

Return type:

double, 1xN

Khurana and Schwarzl (2005) Jupiter model functions

MagFldJupiterKS2005()

functions.KS2005functions.MagFldJupiterKS2005(r_km, theta, phi, ets, SPHOUT, coeffPath, nHeadLines, AS_CODED)

Calculate Jupiter magnetic field at evaluation points and times based on the Khurana and Schwarzl (2005) model.

Model implementation based on Fortran code from K. Khurana and H. Schwarzl after Khurana and Schwarzl (2005) https://doi.org/10.1029/2004JA010757. Much of the original software was refactored to create this module, although some functions (e.g. KS_VIP4noDipole(), KS_Upot()) were difficult to interpret and so were translated near-exactly without refactoring. These functions are identifiable from their use of non-descript variable naming, i.e. lots of single-letter variable names. In several places, more precise model coefficients or corrections to certain parameters have been implemented, but can be toggled off by setting AS_CODED to true, thereby matching the behavior of the original Fortran software.

See GetModelOpts() or refer to Khurana and Schwarzl (2005) for a description of the model. Note that although the embedded magnetopause current model matches the Fortran implementation, most of this model was not documented in the code or the publication, so some variable meanings are based on speculation.

Parameters:
  • r_km (double, 1xN) – Radius of evaluation points in km from planet barycenter.

  • theta (double, 1xN) – Planetocentric colatitude of evaluation points in radians. Must be the same length as r_km.

  • phi (double, 1xN) – Planetocentric east longitude of evaluation points in radians. Must be the same length as r_km.

  • ets (double, 1xN) – Ephemeris times associated with evaluation points in TDB seconds relative to J2000. Must be the same length as r_km.

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory where coefficients files are located.

  • nHeadLines (int, default=2) – Number of header lines in model coefficient files. This is typically 2: One for a description of the file contents and one for column headers.

  • AS_CODED (bool, default=0) –

    Whether to use model coefficients and implementation matching the original Fortran code (true) or with increased precision and corrected parameters (false). Specifically, unless AS_CODED is true:

    • Convert distances to units of RJ = 71323 km as in VIP4 paper (Connerney et al., 1998 https://doi.org/10.1029/97JA03726)

    • Use more precise VIP4 tables of g and h

    • Use corrected values for Jupiter rotation period and rate

Returns:

  • Bvec_nT (double, 3xN) – Magnetic field vector in System III coordinates at evaluation points in nT. Output rows are x, y, z respectively for cartesian or r, theta, phi for spherical.

  • Mdip_nT (double, 1x3) – Dipole magnetic moment in coordinates matching Bvec_nT, as surface-equivalent nT, i.e. this vector times the planet volume yields the magnetic moment in SI units.

  • Odip_km (double, 1x3) – Dipole magnetic moment offset in km from planet barycenter in coordinates matching Bvec_nT.

KS_BMPperp()

functions.KS2005functions.KS_BMPperp(rho, phi, x, Nmodes)

Get magnetic field components perpendicular to the magnetopause boundary.

The coordinates passed to this function are essentially cylindrical coordinates with the Jupiter–Sun–Orbital (JSO) x axis as the z axis of the coordinates and \(phi = 0\) measured from the JSO y axis. See KS_S3CtoJSO() for a definition of the JSO frame.

Parameters:
  • rho (double, 1xN) – Distance from JSO x axis in planetary radii.

  • phi (double, 1xN) – Azimuthal angle between JSO y axis and evaluation point in radians.

  • x (double, 1xN) – Distance from JSO yz plane in planetary radii.

  • Nmodes (int) – “Number of dipole modes” is how the parameter was labeled in the original Fortran code. Its usage in KS_BMPperp() suggests it’s some method of indexing Bessel functions and spherical harmonics in that function.

Returns:

Brds, Bpds, Bzds – Magnetic field contribution from shielded dipole aligned to cylindrical System III axes in nT.

Return type:

double, 1xN

KS_BtailShield()

functions.KS2005functions.KS_BtailShield(M, Mode, xJSO, yJSO, zJSO)

Get magnetotail contribution to shielded dipole field.

See KS_S3CtoJSO() for a definition of the Jupiter–Sun–Orbital (JSO) frame.

Parameters:
  • M (int) – Dimension (MxM) of magnetotail coefficients a and c, typically Mode + 1.

  • Mode (int) – Parameter selection for magnetotail model. Passing 7 or greater will include all coefficients.

  • xJSO (double, 1xN) – x coordinate of evaluation points in JSO frame.

  • yJSO (double, 1xN) – y coordinate of evaluation points in JSO frame.

  • zJSO (double, 1xN) – z coordinate of evaluation points in JSO frame.

Returns:

BxJSO, ByJSO, BzJSO – Magnetic field contribution from magnetotail in nT in JSO frame aligned with cartesian axes.

Return type:

double, 1xN

KS_CheckIfInsideMappedMP()

functions.KS2005functions.KS_CheckIfInsideMappedMP(xS3, yS3, zS3, zNS3, ctimes, AS_CODED)

Return logical array indicating whether each point is inside the mapped magnetopause.

Parameters:
  • xS3 (double, 1xN) – x coordinate of evaluation points in System III frame in planetary radii.

  • yS3 (double, 1xN) – y coordinate of evaluation points in System III frame in planetary radii.

  • zS3 (double, 1xN) – z coordinate of evaluation points in System III frame in planetary radii.

  • zNS3 (double, 1xN) – z distance of current sheet from equatorial plane in the System III frame in planetary radii.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

ptInsideMpause – Whether each evaluation point is inside (true) or outside (false) of the magnetopause.

Return type:

bool, 1xN

KS_CheckIfInsideMP()

functions.KS2005functions.KS_CheckIfInsideMP(xJSM, yJSM, zJSM)

Return logical array indicating whether each point is inside the unmapped magnetopause.

Parameters:
  • xJSM (double, 1xN) – x coordinate of evaluation point in Jupiter–Sun–magnetic (JSM) frame in planetary radii. See KS_BJSMtoBxyz for a definition of the JSM frame.

  • yJSM (double, 1xN) – y coordinate of evaluation point in JSM frame in planetary radii.

  • zJSM (double, 1xN) – z coordinate of evaluation point in JSM frame in planetary radii.

Returns:

ptInsideMpause – Whether each evaluation point is inside (true) or outside (false) of the magnetopause.

Return type:

bool, 1xN

KS_CsheetN()

functions.KS2005functions.KS_CsheetN(xS3, yS3, zS3, localTime_rads, stheta, ctimes, AS_CODED)

Determine current sheet normal vectors at the given evaluation points for the Khurana and Schwarzl (2005) model.

Parameters:
  • xS3 (double, 1xN) – x component of evaluation point position vector in the System III frame.

  • yS3 (double, 1xN) – y component of evaluation point position vector in the System III frame.

  • zS3 (double, 1xN) – z component of evaluation point position vector in the System III frame.

  • localTime_rads (double, 1xN) – Local time in radians, i.e. east longitude in the JSO frame, which is the azimuthal angle between the plane containing the Sun and the evaluation point.

  • stheta (double, 1xN) – (speculated) Angle between the planet–Sun direction and the evaluation point in radians.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

RNx, RNy, RNz – Orientation vector of unit normal for the current sheet in System III coordinates at each evaluation point.

Return type:

double, 1xN

KS_CsheetStruc()

functions.KS2005functions.KS_CsheetStruc(rho, phi, xJSO, yJSO, localTime_rads, stheta, AS_CODED)

Get distance of the current sheet from the System III equatorial plane.

Parameters:
  • rho (double, 1xN) – Distance from z axis of equatorial projection of evaluation point position vector in the System III frame in planetary radii.

  • phi (double, 1xN) – East longitude of evalaution point in radians.

  • xJSO (double, 1xN) – x-aligned vector component in JSO frame (see KS_S3CtoJSO()) in planetary radii.

  • yJSO (double, 1xN) – y-aligned vector component in JSO frame (see KS_S3CtoJSO()) in planetary radii.

  • localTime_rads (double, 1xN) – Local time in radians, i.e. east longitude in the JSO frame, which is the azimuthal angle between the plane containing the Sun and the evaluation point.

  • stheta (double, 1xN) – (speculated) Angle between the planet–Sun direction and the evaluation point in radians.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

zNS3 – z distance of current sheet from equatorial plane in the System III frame in planetary radii.

Return type:

double, 1xN

KS_ctime2et()

functions.KS2005functions.KS_ctime2et(ctimes)

Converts times from ctime to TDB seconds relative to J2000.

Parameters:

ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

Returns:

ets – Ephemeris times in TDB seconds relative to J2000.

Return type:

double, 1xN

KS_ctimer()

functions.KS2005functions.KS_ctimer(ets)

Converts ephemeris times from TDB seconds relative to J2000 to ctime.

MJS note: I don’t understand why the reference is against 1/1/1966. The JS3(1965) epoch is 1965-01-01 at midnight (see https://doi.org/10.1029/GL004i002p00065), but the ctime columns in the trajectory files provided by K. Khurana for validation are consistent with a 1966 reference time.

Parameters:

ets (double, 1xN) – Ephemeris times in TDB seconds relative to J2000.

Returns:

ctimes – Seconds past midnight Jan 1, 1966. See ctimer for more details.

Return type:

double, 1xN

KS_DipoleFld()

functions.KS2005functions.KS_DipoleFld(x, y, z, B0x, B0y, B0z, AS_CODED)

Calculate dipole field as a function of position and dipole components

Parameters:
  • x (double, 1xN) – x coordinate of evalaution points in the System III frame.

  • y (double, 1xN) – y coordinate of evalaution points in the System III frame.

  • z (double, 1xN) – z coordinate of evalaution points in the System III frame.

  • B0x (double, 1xN) – x-aligned component of dipole moment vector in surface-equivalent magnitude in the System III frame in nT.

  • B0y (double, 1xN) – y-aligned component of dipole moment vector in surface-equivalent magnitude in the System III frame in nT.

  • B0z (double, 1xN) – z-aligned component of dipole moment vector in surface-equivalent magnitude in the System III frame in nT.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

Bxd, Byd, Bzd – Magnetic field contribution from magnetic dipole moment aligned to cartesian System III axes in nT.

Return type:

double, 1xN

KS_DipoleShieldCylS3()

functions.KS2005functions.KS_DipoleShieldCylS3(parmod, rmap, pmap, zmapin, sphi, AS_CODED)

Get dipole shield field in cylindrical System III coordinates.

Parameters:
  • parmod (struct) –

    Must contain fields:

    • B0 (double) – Magnetic dipole moment in nT in surface-equivalent magnitude.

    • dipTilt_deg (double) – Tilt of dipole moment vector relative to JSM z axis in degrees. This angle varies as Jupiter rotates, but is treated as a constant in the KS2005 model.

    • Nmodes (int) – “Number of dipole modes” is how the parameter was labeled in the original Fortran code. Its usage in KS_BMPperp() suggests it’s some method of indexing Bessel functions and spherical harmonics in that function.

  • rmap (double, 1xN) – Distance from z axis in current sheet coordinates in planetary radii.

  • pmap (double, 1xN) – Azimuthal angle in current sheet coordinates in radians.

  • zmapin (double, 1xN) – Distance from current sheet normal plane in planetary radii.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in System III coordinates in radians.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

Brds, Bpds, Bzds – Magnetic field contribution from shielded dipole aligned to cylindrical System III axes in nT.

Return type:

double, 1xN

KS_DipoleShielded()

functions.KS2005functions.KS_DipoleShielded(parmod, xpJSO, ypJSO, zpJSO, AS_CODED)

Get shielded dipole field in pseudo-JSO coordinates.

See KS_BpJSOtoBxyz() for a definition of the pseudo-Jupiter–Sun–Orbital (pJSO) frame.

Parameters:
  • parmod (struct) –

    Must contain fields:

    • B0 (double) – Magnetic dipole moment in nT in surface-equivalent magnitude.

    • dipTilt_deg (double) – Tilt of dipole moment vector relative to JSM z axis in degrees. This angle varies as Jupiter rotates, but is treated as a constant in the KS2005 model.

    • Nmodes (int) – “Number of dipole modes” is how the parameter was labeled in the original Fortran code. Its usage in KS_BMPperp() suggests it’s some method of indexing Bessel functions and spherical harmonics in that function.

  • xpJSO (double, 1xN) – x coordinate of evaluation points in pseudo-JSO frame in planetary radii.

  • ypJSO (double, 1xN) – y coordinate of evaluation points in pseudo-JSO frame in planetary radii.

  • zpJSO (double, 1xN) – z coordinate of evaluation points in pseudo-JSO frame in planetary radii.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

BxpJSO, BypJSO, BzpJSO – Magnetic field contribution from shielded dipole in pseudo-JSO frame aligned with cartesian axes in nT.

Return type:

double, 1xN

KS_etimer()

functions.KS2005functions.KS_etimer(ctimes)

Converts times from ctime to etime, which accounts for leap seconds.

Parameters:

ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

Returns:

etime – ctimes adjusted for leap seconds.

Return type:

double, 1xN

KS_GetBIMF()

functions.KS2005functions.KS_GetBIMF(theta, phi, ctimes, AS_CODED)

Determine IMF field contribution in System III spherical coordinates.

Uses a parameterization from the original Fortran implmenetation of the Khurana and Schwarzl (2005) model to evaluate the magnetic field contribution from the interplanetary magnetic field (IMF).

Parameters:
  • theta (double, 1xN) – Colatitude of evaluation points in the System III frame in radians.

  • phi (double, 1xN) – East longitude of evaluation points in the System III frame in radians.

  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

BrS3IMF, BthS3IMF, BphiS3IMF – IMF contribution to magnetic field vectors at the evaluation point in System III spherical coordinates in nT.

Return type:

double, 1xN

KS_GetMappedSunAngle()

functions.KS2005functions.KS_GetMappedSunAngle(stheta, sphi, xpx, xpy, xpz, ypx, ypy, ypz, zpx, zpy, zpz)

Get angle to Sun in current sheet coordinates.

Parameters:
  • stheta (double, 1xN) – (speculated) Angle between the planet–Sun direction and the evaluation point in radians.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in System III coordinates in radians.

  • xpx (double, 1xN) – (speculated) Unit mapping vector x component for current sheet x axis in System III frame.

  • xpy (double, 1xN) – (speculated) Unit mapping vector y component for current sheet x axis in System III frame.

  • xpz (double, 1xN) – (speculated) Unit mapping vector z component for current sheet x axis in System III frame.

  • ypx (double, 1xN) – (speculated) Unit mapping vector x component for current sheet y axis in System III frame.

  • ypy (double, 1xN) – (speculated) Unit mapping vector y component for current sheet y axis in System III frame.

  • ypz (double, 1xN) – (speculated) Unit mapping vector z component for current sheet y axis in System III frame.

  • zpx (double, 1xN) – Current sheet unit normal x component at evaluation point as determined by KS_CsheetN().

  • zpy (double, 1xN) – Current sheet unit normal y component at evaluation point as determined by KS_CsheetN().

  • zpz (double, 1xN) – Current sheet unit normal z component at evaluation point as determined by KS_CsheetN().

Returns:

  • sthetaOut (double, 1xN) – (speculated) Angle between the planet–Sun direction and the evaluation point in radians.

  • sphiOut (double, 1xN) – East longitude of the plane containing the Sun in current sheet coordinates in radians.

KS_JSun()

functions.KS2005functions.KS_JSun(ctimes, AS_CODED)

Calculate the direction and phase angle from Jupiter to the Sun as a function of ctime.

Parameters:
  • ctimes (double, 1xN) – Seconds past midnight Jan 1, 1966. See ctimer for more details.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

  • stheta (double, 1xN) – (speculated) Angle between the planet–Sun direction and the evaluation point in radians.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in System III coordinates in radians.

  • sphase (double, 1xN) – (speculated) Angle between the planet–Sun direction and orbital velocity vector in radians.

KS_TailMagNoTilt()

functions.KS2005functions.KS_TailMagNoTilt(Mode, rho, phi, z, localTime_rads)

Calculate magnetotail field in current sheet coordinates.

Parameters:
  • Mode (int) – Parameter selection for magnetotail model. Passing 7 or greater will include all coefficients.

  • rho (double, 1xN) – Distance from z axis in current sheet coordinates in planetary radii.

  • phi (double, 1xN) – Azimuthal angle in current sheet coordinates in radians.

  • z (double, 1xN) – Distance from current sheet normal plane in planetary radii.

  • localTime_rads (double, 1xN) – Local time in radians, i.e. east longitude in the JSO frame, which is the azimuthal angle between the plane containing the Sun and the evaluation point.

Returns:

Brcs, Bpcs, Bzcs – Magnetic field contribution from magnetotail in current sheet cylindrical coordinates in nT.

Return type:

double, 1xN

KS_TailMagShieldCylS3()

functions.KS2005functions.KS_TailMagShieldCylS3(M, Mode, rmap, pmap, zmap, sphi)

Get magnetotail shield field in cylindrical System III coordinates.

Parameters:
  • M (int) – Dimension (MxM) of magnetotail coefficients a and c, typically Mode + 1.

  • Mode (int) – Parameter selection for magnetotail model. Passing 7 or greater will include all coefficients.

  • rmap (double, 1xN) – Distance from z axis in current sheet coordinates in planetary radii.

  • pmap (double, 1xN) – Azimuthal angle in current sheet coordinates in radians.

  • zmap (double, 1xN) – Distance from current sheet normal plane in planetary radii.

  • sphi (double, 1xN) – East longitude of the plane containing the Sun in System III coordinates in radians.

Returns:

Brcss, Bpcss, Bzcss – Magnetic field contribution from magnetotail in current sheet cylindrical coordinates in nT.

Return type:

double, 1xN

KS_Upot()

functions.KS2005functions.KS_Upot(x, y, z, a, c, p, r, M)

Calculate magnetic potential terms from model coefficients for magnetotail field contribution.

Parameters:
  • x (double, 1xN) – x component of evaluation point in the System III frame.

  • y (double, 1xN) – y component of evaluation point in the System III frame.

  • z (double, 1xN) – z component of evaluation point in the System III frame.

  • a (double, 8x64) – Magnetotail parameterization coefficients.

  • c (double, 8x64) – Magnetotail parameterization coefficients.

  • p (double, 8x16) – Magnetotail parameterization coefficients.

  • r (double, 8x16) – Magnetotail parameterization coefficients.

  • M (int) – Dimension (MxM) of magnetotail coefficients a and c, typically Mode + 1.

Returns:

outputName – Description.

Return type:

type, dims

KS_VIP4noDipole()

functions.KS2005functions.KS_VIP4noDipole(r_Rp, theta, phi, gVIP4, hVIP4, AS_CODED)

Calculate magnetic field vectors at evaluation points for the VIP4 model, sans the dipole moment.

Get magnetic field vectors in System III spherical coordinates for the non-dipole multipole moments of the VIP4 model. We omit the dipole moments because we have already accounted for them in evaluating the shielded dipole magnetosphere model.

Parameters:
  • r_Rp (double, 1xN) – Distance from Jupiter center of mass for evaluation points in terms of planet radius.

  • theta (double, 1xN) – Colatitude of evaluation points in System III coordinates in radians.

  • phi (double, 1xN) – East longitude of evaluation points in System III coordinates in radians.

  • gVIP4 (double, 4x5) – Schmidt semi-normalized g Gauss coefficient of VIP4 model in gauss.

  • hVIP4 (double, 4x5) – Schmidt semi-normalized h Gauss coefficient of VIP4 model in gauss.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

BrVIP4, BthVIP4, BphiVIP4 – Evalauted magnetic field vector components in System III spherical coordinates.

Return type:

double, 1xN

Additional supporting functions

ExcitationSpectrum()

functions.ExcitationSpectrum(moonName, nOsc, rate, Tinterest_h, coordType, outData, coeffPath, fPatternFT, fPatternTseries, magPhase, DO_EPO, SPHOUT, DO_MPAUSE)

Calculate a Fast Fourier Transform (FFT) spectrum of magnetic oscillations for a given moon.

Note

Although this function returns degree-1 complex excitation moments, these moments cannot be used for reproducing the input time series for any period except Tinterest_h. This is because the phase for each peak in the FFT can be tuned only based on the time step of the input time series—if a given oscillation period does not cover a whole number of time steps, the calculated phase of the complex peak amplitude in the FFT will drift.

For calculating excitation moments, instead use LLSdecomposition().

Parameters:
  • moonName (char, 1xC) – Name of target moon for which to calculate an excitation spectrum.

  • nOsc (int, default=5000) – Number of oscillations of Tinterest_h over which to span the time series to transform.

  • rate (int, default=100) – Number of time steps to use for each oscillation period.

  • Tinterest_h (double, optional) – Period of key magnetic oscillation to which to tune the FFT. The excitation moment for this period is the one most accurately reproduced from the FFT inversion. The default value is the synodic period, as evaluated from values returned by GetBodyParams().

  • outData (char, 1xF, default='out') – Directory to use for output of complex spectrum amplitudes.

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory containing model coefficients files.

  • fPatternFT (char, 1xG, default='FTdata') – Pattern for file names of FFT spectrum data saved to disk.

  • fPatternTseries (char, 1xH, default='TseriesData') – Pattern for file names of time series data saved to disk.

  • magPhase (double, default=0) – Arbitrary offset in radians by which to rotate the magnetospheric field evaluation.

  • DO_EPO (bool, default=0) – Whether to perform calculations for Europa in EphiO coordinates (true) or cartesian (false).

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • DO_MPAUSE (bool, default=0) – Whether to include a magnetopause current magnetic field model in spectrum evaluation.

Returns:

  • Tpeak_h (double, 1xP) – Oscillation periods in hours associated with peaks in the excitation spectrum greater than 1 nT.

  • B0vec (double, 1x3) – Background uniform magnetic field vector from evaluated time series.

  • B1vec (complex double, 3xP) – Degree-1 complex excitation moments in the selected coordinates cooresponding to each peak in Tpeak_h.

GetBodyParams()

functions.GetBodyParams(moonName)

Return physical properties of the target moon, its orbit, and the parent planet.

Parameters:
  • moonName (char, 1xC) –

    Name of the target moon. Currently implemented options are:

    • 'Moon'

    • 'Io'

    • 'Europa'

    • 'Ganymede'

    • 'Callisto'

    • 'Mimas'

    • 'Enceladus'

    • 'Tethys'

    • 'Dione'

    • 'Rhea'

    • 'Titan'

    • 'Iapetus'

    • 'Miranda'

    • 'Ariel'

    • 'Umbriel'

    • 'Titania'

    • 'Oberon'

    • 'Triton'

  • planets. (and their parent) –

Returns:

  • RparentEq_km (double) – Equatorial radius of the parent planet in km.

  • RparentPol_km (double) – Polar radius of the parent planet in km.

  • RmoonEq_km (double) – Equatorial radius of the target moon in km.

  • a_AU (double) – Semimajor axis of the parent planet’s solar orbit in AU.

  • omegaParent_radps (double) – Angular rotation rate of the parent planet in radians/s.

  • omegaMoon_radps (double) – Angular rotation rate of the target moon in radians/s.

  • Tparent_s (double) – Sidereal rotation period of the parent planet in s.

  • Tmoon_s (double) – Sidereal rotation period of the target moon in s.

  • nutPrecParent (double, 1xN’) – Nutation/precession rate coefficients as specified in the latest IAU report implemented in the loaded PCK file. Size depends on the body and which rates are defined.

  • nutPrecMoon (double, 1xM’) – Nutation/precession rate coefficients as specified in the latest IAU report implemented in the loaded PCK file. Size depends on the body and which rates are defined.

GetDespun()

functions.GetDespun(xyz_Rp, despin)

Rotate cartesian coordinates through some angle.

Intended to “de-spin” a rotating coordinate system for convenience in plotting, e.g. in the case of a trajectory around a rotating planet that is evaluated in planetocentric coordinates.

Parameters:
  • xyz_Rp (double, 3xN) – Cartesian position vectors in planetary radii.

  • despin (double, 1xN) – Angles by which to rotate each position vector, prograde (toward \(+\hat{\phi}\)). Thus, to despin the coordinates by \(\Delta\phi\), negate the angles passed.

Returns:

x, y, z – Rotated cartesian position vectors in planetary radii.

Return type:

double, 1xN

GetExcitations()

functions.GetExcitations(moonName, etMid_day)

Retrieve a list of excitation frequencies to invert.

Returned excitation frequencies consist of combinations of precise orbital periods, including those adjusted for nutation and precession. For most bodies, one or more empirically determined oscillation periods are also included, such as those corresponding to the true anomaly period. Only one of the parameters marked “adjusted” should be used for each body.

Parameters:
  • moonName (char, 1xC) – Name of moon for which to return excitation frequencies.

  • etMid_day (double, default=0) – Ephemeris time (ET) near the center of time series to be evaluated in days. Default is J2000.

Returns:

fList – Contains columns:

  • double: Frequencies of excitations in Hz.

  • string: Descriptive names for each excitation.

Return type:

table, P’x2

GetModelOpts()

functions.GetModelOpts(planet, opt, MPopt)

Convert from an index number to a compatible combination of planetary field models.

For each planet, only specific internal and external magnetic field models are mutually compatible, as some are derived using parallel assumptions or outputs from other models. For ease of implementation, this function provides an index-based lookup to select from among valid model combinations.

Parameters:
Returns:

  • MagModel (char, 1xD) – Planet intrinsic field model.

  • CsheetModel (char, 1xE) – Current sheet magnetic field model to apply, if applicable.

  • MPmodel (char, 1xF) – Magnetopause current magnetic field model to apply, if applicable.

  • magModelDescrip (char, 1xG) – Text description of magnetic field model combination for plot labels.

  • fEnd (char, 1xH) – File name descriptor to append for model-specific outputs, such as excitation moments.

GetMinMoonDist()

functions.GetMinMoonDist(sc, parentName, ets)

Get distance from a target to the nearest major moon of a planetary system.

Uses loaded SPICE kernels to determine the distance from the target to the major moons of the specified parent body at given ephemeris times (ETs). ETs are in terms of seconds relative to J2000. Position vectors are evaluated in IAU coordinates for the parent body since we only consider line-of-sight distance. Returns the distance to the nearest moon at each ET and the radial distance to the parent body center of mass, both in km.

Parameters:
  • sc (char, 1xC) – Name of target body (e.g. spacecraft) understood by SPICE when upper-cased.

  • parentName (char, 1xD) –

    Name of parent body from which to determine relative distance to the target and for the moons of this body. Must be one of the following:

    • Earth

    • Jupiter

    • Saturn

    • Uranus

    • Neptune

    and have SPICE kernels present in the loaded kernel pool sufficient to determine all relative locations.

  • ets (double, 1xN) – Ephemeris times in TDB seconds relative to J2000.

Returns:

  • rMinMoon_km (double, 1xN) – Radial distance from the target to the nearest major moon of the parent body in km.

  • rPar_km (double, 1xN) – Radial distance from the target to the parent body in km.

GetPosSpice()

functions.GetPosSpice(moonName, parentName, t_h, S3coords)

Retrieve locations of a target relative to a parent body at specified ephemeris times.

Uses loaded SPICE kernels to determine the position of one body relative to another in a specific coordinate system. If the coordinates are not passed, standard cartesian coordinates are selected. Ephemeris times are passed in terms of hours relative to J2000.

Parameters:
  • moonName (char, 1xC) – Name of target body, which must match a code name understood by SPICE and present in the loaded kernel pool when upper-cased. Can be any body loaded in the kernel pool for which sufficient data is present to determine location relative to the parent.

  • parentName (char, 1xD) – Name of parent body understood by SPICE when upper-cased.

  • t_h (double, 1xN) – Ephemeris times in TDB hours relative to J2000.

  • S3coords (char, 1xE, optional) – Standard (System III) coordinates to use for the parent body. Defaults to ULS for Uranus, NLS for Neptune, or IAU_PARENT for any other parent body. Any coordinate system understood by SPICE will return a result, but the first 3 return values from this function will only be converted correctly into spherical coordinates if cartesian coordinates are specified. Magnetic field model evaluations will only work correctly with the default coordinates.

Returns:

  • r_km (double, 1xN) – Radii from parent body center of mass in km.

  • theta_rad (double, 1xN) – Colatitude of positions relative to parent body spin pole in radians.

  • phi_rad (double, 1xN) – East longitude of positions relative to IAU-defined parent body prime meridian in radians.

  • xyz_km (double, 3xN) – Position vectors of target body for each time in the selected coordinates.

  • S3coords (char, 1xE) – Coordinate system of position vectors.

GetTargetMoonDist()

functions.GetTargetMoonDist(sc, moonName, parentName, ets)

Get distance between two targets at given ephemeris times.

Uses loaded SPICE kernels to determine the distance from the first target to the second in terms of planetary equatorial radii of the second target. Planetary radius is as defined in the last-loaded text PCK file. Ephemeris times are passed in terms of seconds relative to J2000. Position vectors are evaluated in IAU coordinates for the parent body of the second target because these coordinates are always defined in SPICE and we are only finding the line-of-sight distance, so the coordinates don’t matter.

Parameters:
  • sc (char, 1xC) – Name of target body (e.g. spacecraft) understood by SPICE when upper-cased.

  • moonName (char, 1xD) – Name of target body from which to determine relative distance. Must match a code name understood by SPICE and present in the loaded kernel pool when upper-cased. Can be any body loaded in the kernel pool for which sufficient data is present to determine sc location relative to this body AND for which a radius can be determined from loaded PCK files.

  • parentName (char, 1xE) – Name of parent body understood by SPICE when upper-cased to use for evaluation of position vectors.

  • ets (double, 1xN) – Ephemeris times in TDB seconds relative to J2000.

Returns:

r_RM – Distance between target bodies in equatorial radii of the second target body (moonName).

Return type:

double, 1xN

HodogramSingle()

functions.HodogramSingle(moonName, parentName, magModelDescrip, Bvec, SPHOUT, COMPARE_SEUFERT, COMPARE_PHIO)

Organizes data and labels for plotting a hodogram, showing the projection of the tip of the magnetic field vector in a relevant plane, usually the equatorial plane.

Must be plotted with a separate function.

Parameters:
  • moonName (char, 1xC) – Name of moon for which to evaluate excitation moments.

  • parentName (char, 1xD) – Name of parent planet for the moon for which to evaluate excitation moments.

  • magModelDescrip (char, 1xE) – Text description of the magnetic field model that was evalauted for the input time series.

  • Bvec (double, 3xN) – Magnetic field vectors at each measurement time.

  • SPHOUT (bool, default=0) – Whether to return vectors aligned to spherical coordinate axes (true) or cartesian (false).

  • COMPARE_SEUFERT (bool, default=1) – Whether to plot coordinate directions aligned to the axes presented in Seufert et al. (2011) https://doi.org/10.1016/j.icarus.2011.03.017. Only has an effect when SPHOUT is true, because Seufert et al. used spherical coordinates for evaluating the excitation spectra of Jupiter’s moons.

  • COMPARE_PHIO (bool, default=1) – Whether to plot components in the same orientation as PhiO axes, as in past work, and label the axes with the comparisons (true), or just plot IAU frame axes without comparisons (false).

Returns:

  • xx (double, 1xN) – x axis data.

  • yy (double, 1xN) – y axis data.

  • windowName (char, 1xC) – Name to use for the interactive figure window, which is shown when LIVE_PLOTS is true.

  • titleInfo (char, 1xD) – Description of plot contents to print at the top of the figure.

  • xInfo (char, 1xE) – Label for x axis of plot.

  • yInfo (char, 1xF) – Label for y axis of plot.

  • fName (char, 1xG) – File name to use for output figures, sans extension.

  • xlims (double, 1x2) – Minimum and maximum values to display on x axis, respectively.

  • ylims (double, 1x2) – Minimum and maximum values to display on y axis, respectively.

LoadSpice()

functions.LoadSpice(moonName, sc, spiceDir)

Load SPICE kernels needed to determine positions for an input spacecraft and target moon.

Loads into the SPICE kernel pool all those kernels needed to evaluate spacecraft and moon positions relative to the parent planet.

Parameters:
  • moonName (char, 1xC) –

    Name of target body. Currently implemented options are:

    • 'Moon'

    • 'Io'

    • 'Europa'

    • 'Ganymede'

    • 'Callisto'

    • 'Mimas'

    • 'Enceladus'

    • 'Dione'

    • 'Rhea'

    • 'Titan'

    • 'Miranda'

    • 'Ariel'

    • 'Umbriel'

    • 'Titania'

    • 'Oberon'

    • 'Triton'

    and their parent planets. The same binary kernels are loaded for any body within a particular planetary system, as these kernels contain information for all listed satellites and the planet.

  • sc (char, 1xD) –

    Spacecraft name for which trajectories will be loaded from relevant binary SPICE kernels. Currently implemented options are:

    • 'Swarm'

    • 'Galileo'

    • 'Cassini'

    • 'Juno'

    • 'Voyager 1'

    • 'Voyager 2'

  • spiceDir (char, 1xE, default='SPICE') – Directory where SPICE kernel files are located.

Returns:

parentName – Name of the parent body.

Return type:

char, 1xF

LshellTrace()

functions.LshellTrace(parentName, opt, MPopt, DO_MPAUSE, sc, t_h, Nmaxin, coeffPath)

Determine the L shell of the target body at specified ephemeris times.

Trace magnetic field line from the location of a target (e.g. a spacecraft) to the magnetic dipole equator to determine the radial distance in planetary radii to the planet’s center of mass at that location, also known as the L shell.

Note

This function is a work in progress and does not work correctly yet.

Parameters:
  • parentName (char, 1xC) – Parent body name recognized by SPICE. Sufficient kernels must be loaded in the pool to determine the location of the target relative to this body at each ephemeris time.

  • opt (int) – Magnetospheric field model option selection. 0 selects the default model—see GetModelOpts() for more details.

  • MPopt (int) – Magnetopause current magnetic field model option selection. 0 selects the default model—see GetModelOpts() for more details.

  • DO_MPAUSE (bool) – Whether to include magnetopause currents in the magnetic field model.

  • sc (char, 1xD) – Target body for which to evaluate L shells.

  • t_h (double, 1xN) – Ephemeris times to evaluate in TDB hours relative to J2000.

  • Nmaxin (int, default=Inf) – Spherical harmonic degree to which to limit magnetic field model evaluation in tracing field lines. Using the default (Inf) makes use of the highest degree both implemented in MagFldParent() and present in the model coefficients.

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory containing model coefficients files.

Returns:

L – L shell of each evaluation spacetime point for the target.

Return type:

double, 1xN

WriteTimeSeries()

functions.WriteTimeSeries(scName, opt, MPopt, moonName, orbNum, outFstr, FULLORBITS)

Evaluate a given magnetic field model and optional magnetopause model at spacecraft locations and save to disk.

Currently implemented for Galileo, Juno, and Cassini.

Datasets:

Parameters:
  • scName (string, default="Galileo") – Spacecraft name for which magnetic field data will be compared against implemented models. A directory must exist with this name in the MAG directory within the top-level PlanetMag directory. This directory will be searched for body-specific data files, and each of these files will be loaded.

  • opt (int, default=0) – Magnetic field model option ID number to evaluate, as defined in GetModelOpts().

  • MPopt (int, default=-1) – Magnetopause model option ID number to evaluate, as defined in GetModelOpts(). Passing -1 skips magnetopause modeling.

  • moonName (char, 1xC, default='Io') – Name of moon for which scName has flyby data. Flyby data is not used if FULLORBITS is true, but a name must be passed for kernel loading.

  • orbNum (int, default=-1) –

    Orbit number (or year for Cassini) to use for data comparison. Accepts the following special options:

    • -1: Use all orbits/years for which data is available.

    • -2: Use all orbits for which flyby data is available for moonName.

  • outFstr (char, 1xE, default='modelDescrip') – String to append to PDS file names to make output file name. Passing ‘modelDescrip’ adds the magModelDescrip file name code string.

  • FULLORBITS (bool, default=0) – Whether to use full-orbit MAG data files or flyby-only files. Only applies to Galileo data.

Model coefficient functions

modelCoeffs.printIGRFcoeffs.printIGRFcoeffs(fNameIn='IGRF13.shc', fNameOut='coeffsEarthIGRF13', datDir='modelCoeffs')

Convert IGRF-13 spherical harmonic coefficients described in IGRF13.shc to a PlanetMag-compatible format.

Parameters:
  • fNameIn (str, default='IGRF13.shc') – File name for table of spherical harmonics

  • fNameOut (str, default='coeffsEarthIGRF13') – File name base for output .csv files (one for g, one for h)

  • datDir (str, default='modelCoeffs') – Relative directory path where fNameIn is found and to which fNameOut will be printed

coeffsBode1994G()

modelCoeffs.coeffsBode1994G(alpha)

Calculate Jupiter magnetopause field coefficients as fit by Bode (1994).

For use with Engle (1992) method as a function of dipole “precession angle” alpha for Jupiter. See https://apps.dtic.mil/sti/pdfs/ADA284857.pdf.

Parameters:

alpha (double, 1xN) – Nod longitude for planetary dipole moment in degrees. This is the angle between the equatorial projection of the magnetic dipole moment vector and the direction of the Sun. alpha = 0 is when the the dipole nod longitude is at local noon.

Returns:

Gnm – Numerical coefficients that describe the fit by Bode (1994) for each dipole angle alpha in the time series to be modeled.

Return type:

double, 10x11xN

ConvertSwarmCDF()

modelCoeffs.ConvertSwarmCDF(fNameIn, fNameOut, fDir)

Convert magnetic field data files from the Swarm mission from CDF to ASCII text.

Converts CDF file containing Swarm data to .tab text file with positions and field vector components aligned to the IAU_EARTH frame. Spice kernels must already be loaded (use LoadSpice() with 'Earth', 'Swarm'). Note that all indices pertaining to timestamps, locations, and measurement data and the corresponding coordinate systems must be set bespoke to each CDF file. Converting this function to work with other missions’ data products will require a careful examination of the second output from cdfread (info) and the ordering of the contents in info.Variables.

Parameters:
  • fNameIn (char, 1xC) – Name of .cdf input file.

  • fNameOut (char, 1xD) – Name of .tab output file.

  • fDir (char, 1xE, default=fullfile('MAG', 'Swarm')) – Directory where input file is located and output file will be printed.

GetGaussCoeffs()

modelCoeffs.GetGaussCoeffs(planet, InternalFieldModel, coeffPath, nHeadLines)

Read spherical harmonic coefficients from formatted ASCII files on disk.

Spherical harmonic components for each model are returned in gauss, as are the maximum degree implemented in the reported model and in PlanetMag and the planet reference radius in km used in the spherical harmonic expansion, as reported in listed publication from which the coefficients have been collected.

Parameters:
  • planet (char, 1xC) – Planet for which to load model coefficients. Options are the planets listed under InternalFieldModel, as well as the special name 'PureHarmonic', which reads a register file from disk that has a single g or h represented for a specific n,m pair, for testing and validation purposes.

  • InternalFieldModel (char, 1xD) –

    Published spherical harmonic model for the desired planet. Currently implemented are:

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory where coefficients files are located.

  • nHeadLines (int, default=2) – Number of header lines in model coefficient files. This is typically 2: One for a description of the file contents and one for column headers.

Returns:

  • g, h (double, (Nmax)x(Nmax+1)) – Real spherical harmonic coefficients (Gauss coefficients) in Schmidt semi-normalized form for the internal magnetic field of the planet in terms of surface field strength in gauss. Coefficients are each referenced to a specific System III coordinate system: ITRF08 for Earth, ULS for Uranus, NLS for Neptune, and IAU for Jupiter and Saturn. Refer to the listed publications for more details.

  • G, H (double, (NmaxExt)x(NmaxExt+1)) – Real spherical harmonic coefficients (Gauss coefficients) in Schmidt semi-normalized form for the external magnetic field applied to the planet in gauss.

  • PlanetEqRadius (double) – Planet reference (equatorial) radius used in reporting spherical harmonic coefficients from each reference.

  • Nmax (int) – Maximum degree of spherical harmonic coefficients in the selected internal field model.

  • NmaxExt (int) – Maximum degree of spherical harmonic coefficients in the external field model associated with the selected internal field model, if applicable. Set to zero if the model does not include external field coefficients.

KS_coeffsBIMF()

modelCoeffs.KS_coeffsBIMF()

Coefficients describing the interplanetary magnetic field (IMF) in calculations for the Khurana and Schwarzl (2005) model.

Returns:

  • eps1, eps2 (double) – Scaling factors for IMF parameters.

  • ByConst, BzConst (double) – Magnetic field components to use for the IMF.

KS_coeffsCSdeform()

modelCoeffs.KS_coeffsCSdeform()

Import model parameters for current sheet deformation in the magnetotail in the Khurana and Schwarzl (2005) model.

Returns:

  • C (double, 1x6) – Current sheet deformation model parameters.

  • P, D0, D1, D (double) – Current sheet deformation model parameters.

  • f, beta (double, 6x6) – Current sheet deformation model parameters.

KS_coeffsCsheet()

modelCoeffs.KS_coeffsCsheet()

Coefficients characterizing the current sheet in the Khurana and Schwarzl (2005) model.

Returns:

CS – Current sheet model parameters.

Return type:

double, 1x6

KS_coeffsCsheetStruc()

modelCoeffs.KS_coeffsCsheetStruc(AS_CODED)

Import coefficients describing the structure of the current sheet in the Khurana and Schwarzl (2005) model.

Imports coefficients used to determine the current sheet structure used in calculating current sheet location relative to System III equatorial plane.

Parameters:

AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

  • period_hr (double) – Planetary rotation period in hours.

  • omegaJ_deghr (double) – Angular velocity of planetary rotation in degrees per hour.

  • X0 (double) – Distance from JSO \(x=0\) plane beyond which the current sheet is hinged due to solar wind forcing.

  • phip0 (double) – Current sheet parameterization coefficient.

  • obliq (double) – Obliquity of Jupiter relative to the solar spin axis in radians.

  • incl_solar (double) – Orbital inclination of Jupiter relative to the solar spin equator in radians.

  • C (double, 1x5) – Current sheet parameterization coefficient.

  • psi2, psi4 (double) – Current sheet parameterization coefficients.

KS_coeffsJSMdipole()

modelCoeffs.KS_coeffsJSMdipole()

Import parameters describing the orientation of the magnetic dipole moment in the JSM frame.

Includes a 202 degree rotation about z and a 9.6 degree rotation about y. According to the following document, this functions is therefore using “outdated” O4 model parameters: https://lasp.colorado.edu/home/mop/files/2015/02/CoOrd_systems12.pdf

Returns:

dipole – Rotation matrix used in converting between the System III frame and the Jupiter–Sun–magnetic frame.

Return type:

double, 3x3

KS_coeffsJSun()

modelCoeffs.KS_coeffsJSun(AS_CODED)

Import coefficients for Jupiter–Sun angle calculations.

Parameters:

AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

  • yrJup_s (double) – Jupiter’s orbital period around the Sun in seconds.

  • omegaJ_degps (double) – Sidereal rotation rate of Jupiter in degrees per second.

  • omegayr_radps (double) – Angular speed of orbital motion for Jupiter in radians per second.

  • etime1 (double) – Reference time from which orbital/rotation parameters are referenced.

  • obliq (double) – Obliquity of Jupiter relative to the solar spin axis in radians.

  • tan_ob (double) – tan(obliq).

  • aa, bb (double, 1x7) – Solar phase angle parameterization coefficients.

  • deltaPhi (double) – Sun direction offset in degrees.

KS_coeffsMpause()

modelCoeffs.KS_coeffsMpause()

Import parameter constants for KS_BMPperp().

These coefficients describe the structure of the magnetopause in the Khurana and Schwarzl (2005) model.

Returns:

  • A (double, 1x30) – Magnetopause parameterization coefficients.

  • B (double, 1x15) – Magnetopause parameterization coefficients.

KS_coeffsMtail()

modelCoeffs.KS_coeffsMtail()

Import parameter constants for KS_BtailShield().

These coefficients describe the structure of a shielded dipole’s magnetotail in the Khurana and Schwarzl (2005) model.

Returns:

  • a, c (double, 8x64) – Magnetotail parameterization coefficients.

  • p, r (double, 8x16) – Magnetotail parameterization coefficients.

KS_coeffsVIP4()

modelCoeffs.KS_coeffsVIP4(coeffPath, nHeadLines, AS_CODED)

Import spherical harmonic coefficients for the VIP4 model.

See GetGaussCoeffs() for a more detailed description of the VIP4 model. This function includes the option to use less-precise values for the coefficients as originally implemented in the Khurana and Schwarzl (2005) model.

Parameters:
  • coeffPath (char, 1xE, default='modelCoeffs') – Directory where coefficients files are located.

  • nHeadLines (int, default=2) – Number of header lines in model coefficient files. This is typically 2: One for a description of the file contents and one for column headers.

  • AS_CODED (bool, default=0) – Whether to match the original Fortran code (true) or with increased precision and corrected parameters (false). See MagFldJupiterKS2005() for more details.

Returns:

g, h – Real spherical harmonic coefficients (Gauss coefficients) in Schmidt semi-normalized form for the internal magnetic field of the planet in terms of surface field strength in gauss in the VIP4 model.

Return type:

double, 4x5

Magnetopause functions

GetMPsurfAB2005()

functions.magnetopauseModels.GetMPsurfAB2005(Rss_km, xyzPSM_km)

Determine paraboloid magnetopause surface shape based on Alexeev and Belenkaya (2005).

Determines shape of magnetopause screening surface current shape from the model detailed by Alexeev and Belenkaya (2005) https://doi.org/10.5194/angeo-23-809-2005. This is a simple paraboloid with an input sub-solar magnetopause standoff distance.

Parameters:
  • Rss_km (double, 1xN) – Sub-solar point magnetopause standoff distance in km.

  • xyzPSM_km (double, 3xN) – Cartesian coordinates of measurement locations in planet–solar–magnetospheric frame in km.

Returns:

xMP_km – Distance from body planet–solar–magnetospheric yz-plane to magnetopause surface for each (y,z) point, i.e. the x coordinate of the magnetopause in the PSM frame.

Return type:

double, 1xN

GetMPsurfS1997()

functions.magnetopauseModels.GetMPsurfS1997(Rss_km, xi, thDSZ)

Determine paraboloid magnetopause surface shape based on Shue et al (1997) model.

Determines shape of magnetopause screening surface current shape from the model detailed by Shue et al. (1997) https://doi.org/10.1029/97JA00196. This is a simple paraboloid with an input sub-solar magnetopause standoff distance.

Parameters:
  • Rss_km (double, 1xN) – Sub-solar point magnetopause standoff distance in km.

  • xi (double, 1x1 or 1xN) – Exponent to scale paraboloid shape.

  • thDSZ (double, 1xN) – Colatitude in dipole–solar–zenith coordinates in radians, i.e. the angle between the planet–Sun direction and the evaluation point location vector.

Returns:

rMP_km – Distance from body center of mass to magnetopause surface for each evaluation point.

Return type:

double, 1xN

GetMPsurfSM1996()

functions.magnetopauseModels.GetMPsurfSM1996(Rss_km, xyzPSMdipO_km)

Determine paraboloid magnetopause surface shape based on Schulz and McNab (1996)

Determines shape of magnetopause screening surface current shape from the model detailed by Schulz and McNabb (1996) https://doi.org/10.1029/95JA02987. This is a tanh paraboloid with an input sub-solar magnetopause standoff distance and dipole offset.

Parameters:
  • Rss_km (double, 1xN) – Sub-solar point magnetopause standoff distance in km.

  • xyzPSMdipO_km (double, 3xN) – Dipole offset vector in planet–solar–magnetospheric coordinates in km.

Returns:

rhoMP_km – Distance from planet–Sun line in km, i.e. \(\rho\) in dipole–solar–zenith cylindrical coordinates.

Return type:

double, 1xN

MPboxHarmonic()

functions.magnetopauseModels.MPboxHarmonic(xyzPSM_Rp, Psi, a, b, c, d, pqrs)

Calculate the magnetic field contributuon at specific points based on a box harmonix model.

Applies a box harnomic model based on Tsyganenko (2002) https://doi.org/10.1029/2001JA000219 as formulated by Arridge and Eggington (2021) https://doi.org/10.1016/j.icarus.2021.114562.

Parameters:
  • xyzPSM_Rp (double, 3xN) – Cartesian coordinates of measurement locations in planet–Sun–magnetic frame in units of planetary radii.

  • Psi (double, 1xN) – Dipole angle relative to perpendicular to pointing into the solar wind in radians, i.e. Psi = 0 is perpendicular to the solar wind direction, Psi = pi/2 is directed into the solar wind, and Psi = -pi/2 is directed along the solar wind velocity vector.

  • a (double, 3x3) – Model coefficients fit from Arridge and Eggington (2021).

  • b (double, 3x3) – Model coefficients fit from Arridge and Eggington (2021).

  • c (double, 4x4) – Model coefficients fit from Arridge and Eggington (2021).

  • d (double, 4x4) – Model coefficients fit from Arridge and Eggington (2021).

  • pqrs (double, 4x4) – Model coefficients fit from Arridge and Eggington (2021). p and r are 3x1, q and s are 4x1.

Returns:

Bx, By, Bz – Magnetic field contribution in nT from magnetopause surface currents in planet–Sun–magnetic coordinates.

Return type:

double, 1xN

MPsphericalHarmonic()

functions.magnetopauseModels.MPsphericalHarmonic(r_Rss, thMP, phiMP, Gnm, Nmax)

Evaluate external-source spherical harmonics for magnetopause surface currents, i.e. all hnm = 0.

Parameters:
  • r_Rss (double, 1xN) – Radial distance from planet center of mass in units of the sub-solar magnetopause standoff distance.

  • thMP (double, 1xN) – Colatitude in planet–Sun–magnetic coordinates in radians. The \(\hat{z}\) axis of this coordinate system is derived from \(\hat{x}\times\hat{y}\), where \(\hat{x}\) is along the planet–Sun direction and \(\hat{y}\) is along \(\mathbf{m}\times\hat{x}\), with \(\mathbf{m}\) the instantaneous magnetic dipole moment vector for the body.

  • phiMP (double, 1xN) – Azimuthal angle in planet–Sun–magnetic coordinates in radians, i.e. the angle between the planet–Sun direction and the position vector’s projection into the dipole equatorial plane.

  • Gnm (double, (Nmax)x(Nmax+1)) – Coefficients for external spherical harmonics in nT at body center.

  • Nmax (int) – Maximum degree to which to limit spherical harmonic models. Only has an effect if the value passed is less than the lesser of the model maximum degree.

Returns:

Br, Bth, Bphi – External-source magnetic field vector components aligned to spherical axes in the System III frame in nT.

Return type:

double, 1xN

Figure plotting functions

ApplyPlotDefaults()

functions.plotting.ApplyPlotDefaults(fig, interpreter, font, fontsize)

Applies default settings to plot labels.

Parameters:
  • fig (figure) – Figure object for the plot to adjust.

  • interpreter (char, 1xC) – The default interpreter to use for rendering text labels on plots.

  • font (char, 1xD) – The default font to use for text labels on plots.

  • fontsize (double) – The default font size in pt to use for text labels on plots.

PlotBandLsq()

functions.plotting.PlotBandLsq(ets, t_h, r_km, theta, phi, xyz_km, BrSC, BthSC, BphiSC, scName, parentName, S3coords, orbStr, opt, MPopt, SEQUENTIAL, coeffPath, figDir, figXtn, LIVE_PLOTS, jt_h, RELATIVE_t, RELATIVE_r)

Plots and calculates comparisons between modeled and measured magnetic fields.

Generates time series data of a specified combination of magnetic field models implemented in PlanetMag and compares against spacecraft measurements of the planetary magnetic field. Comparisons are plotted and least-squares differences are calculated and printed to the terminal. Intended as a final step in model validation.

Parameters:
  • ets (double, 1xN) – Ephemeris times of measurements to compare in TDB seconds relative to J2000.

  • t_h (double, 1xN) – Ephemeris times of measurements to compare in TDB hours relative to J2000 (ets/3600).

  • r_km (double, 1xN) – Radial distance of measurement locations from planet center of mass in km.

  • theta (double, 1xN) – Colatitude of measurement locations from planet center of mass in radians.

  • phi (double, 1xN) – East longitude of measurement locations from planet center of mass in radians.

  • xyz_km (double, 3xN) – Cartesian coordinates of measurement locations in S3coords frame in km.

  • BrSC (double, 1xN) – Radial component (\(B_r\)) of magnetic field measurements in S3coords frame in nT.

  • BthSC (double, 1xN) – Colatitudinal component (\(B_\theta\)) of magnetic field measurements in S3coords frame in nT.

  • BphiSC (double, 1xN) – Azimuthal component (\(B_\phi\)) of magnetic field measurements in S3coords frame in nT.

  • scName (string, 1xS') – Name(s) of spacecraft for measurement comparisons. Accepts a lone string or a list of strings.

  • parentName (char, 1xC) – Name of parent planet the evaluated model(s) and spacecraft measurements apply to.

  • S3coords (char, 1xD) – Standard coordinate frame used for evaluation of magnetic fields.

  • orbStr (char, 1xE) – Description of orbit(s) covered to place in legend labels.

  • opt (int) – Index of planetary magnetic field model. See GetModelOpts() for more details and available options.

  • MPopt (int) – Index of magnetopause current magnetic field model. See GetModelOpts() for more details and available options.

  • SEQUENTIAL (bool, default=0) – Whether to plot points by index or hours relative to a reference time.

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory containing model coefficients files.

  • figDir (char, 1xG, default='figures') – Directory to use for output figures.

  • figXtn (char, 1xH, default='pdf') – Extension to use for figures, which determines the file type.

  • LIVE_PLOTS (bool, default=0) – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • jt_h (double, 1xM, default=[]) – Ephemeris times of Juno measurements in TDB hours relative to J2000 to use in data comparisons.

  • RELATIVE_t (bool, default=0) – Whether to plot points relative to the start of the input time series. Only has an effect when SEQUENTIAL = 0.

  • RELATIVE_r (bool, default=0) – Whether to plot points relative to distance from the parent planet. Overrides SEQUENTIAL and RELATIVE_t.

PlotBandLsqMoon()

functions.plotting.PlotBandLsqMoon(ets, t_h, r_km, theta, phi, xyz_km, r_RM, BxSC, BySC, BzSC, scName, parentName, S3coords, moonName, era, fbStr, opt, MPopt, SEQUENTIAL, dataDir, figDir, figXtn, LIVE_PLOTS, coeffPath, jt_h)

Plots and calculates comparisons between modeled and measured magnetic fields in the vicinity of a target moon.

Generates time series data of a specified combination of magnetic field models implemented in PlanetMag and compares against spacecraft measurements of the planetary magnetic field. Comparisons are plotted and least-squares differences are calculated and printed to the terminal. Intended as a final step in model validation; unlike PlotBandLsq(), this function uses flyby data in the vicinity of the target moon and the recorded excitation moments to model the magnetic field applied by the parent planet.

Note

Because this function uses excitation moments saved from PlanetMag(), that function must be run for each of the input settings (moonName, opt, MPopt, and era) before this one.

Parameters:
  • ets (double, 1xN) – Ephemeris times of measurements to compare in TDB seconds relative to J2000.

  • t_h (double, 1xN) – Ephemeris times of measurements to compare in TDB hours relative to J2000 (ets/3600).

  • r_km (double, 1xN) – Radial distance of measurement locations from planet center of mass in km.

  • theta (double, 1xN) – Colatitude of measurement locations from planet center of mass in radians.

  • phi (double, 1xN) – East longitude of measurement locations from planet center of mass in radians.

  • xyz_km (double, 3xN) – Cartesian coordinates of measurement locations in S3coords frame in km.

  • r_RM (double, 1xN) – Radial distance of measurement locations from moon center of mass in moon equatorial radii.

  • BxSC (double, 1xN) – x component (\(B_x\)) of magnetic field measurements in S3coords frame in nT.

  • BySC (double, 1xN) – y component (\(B_y\)) of magnetic field measurements in S3coords frame in nT.

  • BzSC (double, 1xN) – z component (\(B_z\)) of magnetic field measurements in S3coords frame in nT.

  • scName (string, 1xS') – Name(s) of spacecraft for measurement comparisons. Accepts a lone string or a list of strings.

  • parentName (char, 1xC) – Name of parent planet the evaluated model(s) and spacecraft measurements apply to.

  • S3coords (char, 1xD) – Standard coordinate frame (for the parent planet) used for evaluation of magnetic fields.

  • moonName (char, 1xE) – Name of target moon the evaluated model(s) and spacecraft measurements apply to.

  • era (char, 1xF) – Name of spacecraft mission time frame to use for excitation moments. Must be

  • fbStr (char, 1xG) – Description of flyby(s) covered to place in legend labels.

  • opt (int) – Index of planetary magnetic field model. See GetModelOpts() for more details and available options.

  • MPopt (int) – Index of magnetopause current magnetic field model. See GetModelOpts() for more details and available options.

  • SEQUENTIAL (bool, default=0) – Whether to plot points by index or hours relative to a reference time (closest approach).

  • dataDir (char, 1xH, default='out') – Name of directory where excitation moments are printed to disk.

  • figDir (char, 1xI, default='figures') – Directory to use for output figures.

  • figXtn (char, 1xF, default='pdf') – Extension to use for figures, which determines the file type.

  • LIVE_PLOTS (bool, default=0) – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • coeffPath (char, 1xJ, default='modelCoeffs') – Directory containing model coefficients files.

  • jt_h (double, 1xM, default=[]) – Ephemeris times of Juno measurements in TDB hours relative to J2000 to use in data comparisons.

PlotBandLsqNeptune()

functions.plotting.PlotBandLsqNeptune(ets, t_h, r_km, theta, phi, xyz_km, BrSC, BthSC, BphiSC, scName, S3coords, orbStr, opt, MPopt, SEQUENTIAL, coeffPath, figDir, figXtn, LIVE_PLOTS, INC_PDS, INC_O8, INC_PDS_NM3, INC_NM3)

Plots and calculates comparisons between modeled and measured magnetic fields for Neptune.

Generates time series data of a specified combination of magnetic field models implemented in PlanetMag and compares against spacecraft measurements of the planetary magnetic field. Comparisons are plotted and least-squares differences are calculated and printed to the terminal. Intended as a final step in model validation; because Voyager 2 trajectories reported in PDS files do not perfectly match the trajectories reconstructed from any available SPICE kernels, and because the IAU coordinate system for Neptune uses a System II frame (rotating with a stable atmospheric feature) instead of a System III frame (rotating with the intrinsic magnetic field), this function includes comparisons between evaluation with various relevant coordinate systems.

Parameters:
  • ets (double, 1xN) – Ephemeris times of measurements to compare in TDB seconds relative to J2000.

  • t_h (double, 1xN) – Ephemeris times of measurements to compare in TDB hours relative to J2000 (ets/3600).

  • r_km (double, 1xN) – Radial distance of measurement locations from planet center of mass in km.

  • theta (double, 1xN) – Colatitude of measurement locations from planet center of mass in radians.

  • phi (double, 1xN) – East longitude of measurement locations from planet center of mass in radians.

  • xyz_km (double, 3xN) – Cartesian coordinates of measurement locations in S3coords frame in km.

  • BrSC (double, 1xN) – Radial component (\(B_r\)) of magnetic field measurements in S3coords frame in nT.

  • BthSC (double, 1xN) – Colatitudinal component (\(B_\theta\)) of magnetic field measurements in S3coords frame in nT.

  • BphiSC (double, 1xN) – Azimuthal component (\(B_\phi\)) of magnetic field measurements in S3coords frame in nT.

  • scName (string, 1xS') – Name(s) of spacecraft for measurement comparisons. Accepts a lone string or a list of strings.

  • S3coords (char, 1xD) – Standard coordinate frame used for evaluation of magnetic fields.

  • orbStr (char, 1xE) – Description of orbit(s) covered to place in legend labels.

  • opt (int) – Index of planetary magnetic field model. See GetModelOpts() for more details and available options.

  • MPopt (int) – Index of magnetopause current magnetic field model. See GetModelOpts() for more details and available options.

  • SEQUENTIAL (bool, default=0) – Whether to plot points by index or hours relative to a reference time.

  • coeffPath (char, 1xE, default='modelCoeffs') – Directory containing model coefficients files.

  • figDir (char, 1xG, default='figures') – Directory to use for output figures.

  • figXtn (char, 1xH, default='pdf') – Extension to use for figures, which determines the file type.

  • LIVE_PLOTS (bool, default=0) – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • INC_PDS (bool, default=1) – Whether to include the full O8 model, including unresolved coefficients, with PDS trajectory in NLS frame in coordinate/frame comparison.

  • INC_O8 (bool, default=1) – Whether to include the full O8 model, including unresolved coefficients, in the NLS frame as implemented in PlanetMag in coordinate/frame comparison.

  • INC_PDS_NM3 (bool, default=1) – Whether to include O8 model, limited to Nmax=3, with PDS trajectory in NLS frame in coordinate/frame comparison.

  • INC_NM3 (bool, default=1) – Whether to include O8 model, limited to Nmax=3, in NLS frame as implemented in PlanetMag in coordinate/frame comparison. This implementation uses the Neptune spin pole as referenced at the J2000 epoch and the System III rotation rate as reported in Ness et al. (1989) https://doi.org/10.1126/science.246.4936.1473.

PlotSpectrum()

functions.plotting.PlotSpectrum(moonName, LIVE_PLOTS, figNumber, dataDir, figDir, fPattern, figXtn)

Plot a pre-saved spectrum of magnetic oscillations as derived from an FFT.

Parameters:
  • moonName (char, 1xC) – Name of the moon for which to plot a spectrum. An FFT data file must be present.

  • LIVE_PLOTS (bool, default=0) – Whether to load interactive figure windows for plots (true) or print them to disk (false).

  • figNumber (int, default=0) – Figure number to use for assigning/reusing figure windows. If 0 is passed, the default (no assigned number, use increment) is used.

  • dataDir (char, 1xD, default='out') – Directory where data files are saved.

  • figDir (char, 1xE, default='figures') – Directory to use for output figures.

  • fPattern (char, 1xF, default='FTdata') – Pattern used for saved FFT data file names.

  • figXtn (char, 1xG, default='pdf') – Extension to use for figures, which determines the file type.

SetPlotDefaults()

functions.plotting.SetPlotDefaults()

Sets global Matlab plotting defaults, including some global font switches.

Returns:

  • interpreter (char, 1xC) – The default interpreter to use for rendering text labels on plots.

  • font (char, 1xD) – The default font to use for text labels on plots.

  • fontsize (double) – The default font size in pt to use for text labels on plots.

  • Sets globals

  • ————

  • nmTxt (char, 1xE (global)) – Font switch to use in LaTeX-formatted labels with “normal” text.

  • bnmTxt (char, 1xF (global)) – Font switch to use in LaTeX-formatted title labels (with bold text).

  • mathTxt (char, 1xG (global)) – Font switch to use in LaTeX-formatted labels with math text.

  • bmathTxt (char, 1xH (global)) – Font switch to use in LaTeX-formatted labels with bolded math text.

vlines()

functions.plotting.vlines(x, color)

Plot vertical lines at specific x locations with the named color.

Vertical lines are plotted on the currently focused axes object.

Parameters:
  • x (double, 1xN) – x coordinates at which to plot vertical lines in data units.

  • color (char, 1x7) – Object interpretable by Matlab as a color (e.g. HTML color string '#FF0000').