API: Functions
The following is a full listing of the public functions within the package.
MagNav.TL_mat2vec — MethodTL_mat2vec(TL_coef_p, TL_coef_i, TL_coef_e, terms; Bt_scale = 50000f0)Extract the vector form of Tolles-Lawson coefficients from the matrix form.
Arguments:
TL_coef_p: length-3vector of permanent field coefficientsTL_coef_i:3x3symmetric matrix of induced field coefficients, denormalizedTL_coef_e:3x3matrix of eddy current coefficients, denormalizedterms: Tolles-Lawson terms used {:permanent,:induced,:eddy}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
TL_coef: Tolles-Lawson coefficients
MagNav.TL_vec2mat — MethodTL_vec2mat(TL_coef::Vector, terms; Bt_scale = 50000f0)Extract the matrix form of Tolles-Lawson coefficients from the vector form.
Arguments:
TL_coef: Tolles-Lawson coefficients (must include:permanent&:induced)terms: Tolles-Lawson terms used {:permanent,:induced,:eddy}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
TL_coef_p: length-3vector of permanent field coefficientsTL_coef_i:3x3symmetric matrix of induced field coefficients, denormalizedTL_coef_e:3x3matrix of eddy current coefficients, denormalized
MagNav.bpf_data! — Methodbpf_data!(x::AbstractVecOrMat; bpf=get_bpf())Bandpass (or low-pass or high-pass) filter vector or columns of matrix.
Arguments:
x: data vector or matrixbpf: (optional) filter object
Returns:
nothing:xis mutated with filtered data
MagNav.bpf_data — Methodbpf_data(x::AbstractMatrix; bpf=get_bpf())Bandpass (or low-pass or high-pass) filter columns of matrix.
Arguments:
x: data matrix (e.g., Tolles-LawsonAmatrix)bpf: (optional) filter object
Returns:
x_f: data matrix, filtered
MagNav.bpf_data — Methodbpf_data(x::AbstractVector; bpf=get_bpf())Bandpass (or low-pass or high-pass) filter vector.
Arguments:
x: data vector (e.g., magnetometer measurements)bpf: (optional) filter object
Returns:
x_f: data vector, filtered
MagNav.chunk_data — Methodchunk_data(x, y, l_window::Int)Break data into non-overlapping sequences (vectors).
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorl_window: temporal window length
Returns:
x_seqs: sequence (vector) ofNfxl_windowdata matricesy_seqs: sequence (vector) of length-l_windowtarget vectors
MagNav.comp_m2bc_test — Methodcomp_m2bc_test(comp_params::NNCompParams, lines,
df_line::DataFrame, df_flight::DataFrame, df_map::DataFrame;
silent::Bool = false)Evaluate performance of neural network-based aeromagnetic compensation, model 2b or 2c with additional outputs for explainability.
Arguments:
comp_params:NNCompParamsneural network-based aeromagnetic compensation parameters structlines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
silent: (optional) if true, no print outs
Returns:
y_nn: length-Nneural network compensation portiony_TL: length-NTolles-Lawson compensation portiony: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_m3_test — Methodcomp_m3_test(comp_params::NNCompParams, lines,
df_line::DataFrame, df_flight::DataFrame, df_map::DataFrame;
temp_params::TempParams = TempParams(),
silent::Bool = false)Evaluate performance of neural network-based aeromagnetic compensation, model 3 with additional outputs for explainability.
Arguments:
comp_params:NNCompParamsneural network-based aeromagnetic compensation parameters structlines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
temp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
TL_perm:3xNmatrix of TL permanent vector fieldTL_induced:3xNmatrix of TL induced vector fieldTL_eddy:3xNmatrix of TL eddy current vector fieldTL_aircraft:3xNmatrix of TL aircraft vector fieldB_unit:3xNmatrix of normalized vector magnetometer measurementsB_vec:3xNmatrix of vector magnetometer measurementsy_nn:3xNmatrix of vector neural network correction (for scalar models, in direction ofBt)vec_aircraft:3xNmatrix of predicted aircraft vector fieldy: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_test — Functioncomp_test(comp_params::CompParams, xyz::XYZ, ind,
mapS::Union{MapS,MapSd,MapS3D} = mapS_null;
temp_params::TempParams = TempParams(),
silent::Bool = false)Evaluate performance of an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
xyz:XYZflight data structind: selected data indicesmapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct, only used fory_type = :b, :ctemp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
y: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Ncompensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_test — Methodcomp_test(comp_params::CompParams, lines,
df_line::DataFrame, df_flight::DataFrame, df_map::DataFrame;
temp_params::TempParams = TempParams(),
silent::Bool = false)Evaluate performance of an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
lines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
temp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
y: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_train — Functioncomp_train(comp_params::CompParams, xyz::XYZ, ind,
mapS::Union{MapS,MapSd,MapS3D} = mapS_null;
temp_params::TempParams = TempParams(),
xyz_test::XYZ = xyz,
ind_test = BitVector(),
silent::Bool = false)Train an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
xyz:XYZflight data structind: selected data indicesmapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct, only used fory_type = :b, :ctemp_params: (optional)TempParamstemporary temporal parameters structxyz_test: (optional)XYZheld-out test data structind_test: (optional) indices for test data structsilent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters structy: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Ncompensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_train — Functioncomp_train(comp_params::CompParams, xyz_vec::Vector, ind_vec::Vector,
mapS::Union{MapS,MapSd,MapS3D} = mapS_null;
temp_params::TempParams = TempParams(),
xyz_test::XYZ = xyz_vec[1],
ind_test = BitVector(),
silent::Bool = false)Train an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
xyz_vec: vector ofXYZflight data structsind_vec: vector of selected data indicesmapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct, only used fory_type = :b, :ctemp_params: (optional)TempParamstemporary temporal parameters structxyz_test: (optional)XYZheld-out test data structind_test: (optional) indices for test data structsilent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters structy: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Ncompensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_train — Methodcomp_train(comp_params::CompParams, lines,
df_line::DataFrame, df_flight::DataFrame, df_map::DataFrame;
temp_params::TempParams = TempParams(),
silent::Bool = false)Train an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
lines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
temp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters structy: length-Ntarget vectory_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_train_test — Functioncomp_train_test(comp_params::CompParams,
xyz_train::XYZ, xyz_test::XYZ, ind_train, ind_test,
mapS_train::Union{MapS,MapSd,MapS3D} = mapS_null,
mapS_test::Union{MapS,MapSd,MapS3D} = mapS_null;
temp_params::TempParams = TempParams(),
silent::Bool = false)Train & evaluate performance of an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
xyz_train:XYZflight data struct for trainingxyz_test:XYZflight data struct for testingind_train: selected data indices for trainingind_test: selected data indices for testingmapS_train: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct for training, only used fory_type = :b, :cmapS_test: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct for testing, only used fory_type = :b, :ctemp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters structy_train: length-N_traintraining target vectory_train_hat: length-N_traintraining prediction vectorerr_train: length-N_traintraining compensation errory_test: length-N_testtesting target vectory_test_hat: length-N_testtesting prediction vectorerr_test: length-N_testtesting compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.comp_train_test — Methodcomp_train_test(comp_params::CompParams, lines_train, lines_test,
df_line::DataFrame, df_flight::DataFrame, df_map::DataFrame;
temp_params::TempParams = TempParams(),
silent::Bool = false)Train & evaluate performance of an aeromagnetic compensation model.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
lines_train: selected line number(s) for traininglines_test: selected line number(s) for testingdf_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
temp_params: (optional)TempParamstemporary temporal parameters structsilent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters structy_train: length-N_traintraining target vectory_train_hat: length-N_traintraining prediction vectorerr_train: length-N_trainmean-corrected (per line) training compensation errory_test: length-N_testtesting target vectory_test_hat: length-N_testtesting prediction vectorerr_test: length-N_testmean-corrected (per line) testing compensation errorfeatures: length-Nffeature vector (including components of TLA, etc.)
MagNav.compare_fields — Methodcompare_fields(s1, s2; silent::Bool = false)Compare data for each data field in 2 structs of the same type.
Arguments:
s1: struct 1s2: struct 2silent: (optional) if true, no summary print out
Returns:
N_dif: ifsilent = false, number of different fields
MagNav.corrupt_mag — Methodcorrupt_mag(mag_c, Bx, By, Bz;
dt = 0.1,
cor_sigma = 1.0,
cor_tau = 600.0,
cor_var = 1.0^2,
cor_drift = 0.001,
cor_perm_mag = 5.0,
cor_ind_mag = 5.0,
cor_eddy_mag = 0.5)Corrupt compensated (clean) magnetometer measurements with random, FOGM, drift, and Tolles-Lawson noise to create uncompensated (corrupted) scalar magnetometer measurements. FOGM is a First-order Gauss-Markov stochastic process.
Arguments:
mag_c: compensated (clean) scalar magnetometer measurements [nT]Bx,By,Bz: vector magnetometer measurements [nT]dt: (optional) measurement time step [s]cor_sigma: (optional) corruption FOGM catch-all bias [nT]cor_tau: (optional) corruption FOGM catch-all time constant [s]cor_var: (optional) corruption measurement (white) noise variance [nT^2]cor_drift: (optional) corruption measurement linear drift [nT/s]cor_perm_mag: (optional) corruption permanent field TL coef std devcor_ind_mag: (optional) corruption induced field TL coef std devcor_eddy_mag: (optional) corruption eddy current TL coef std dev
Returns:
mag_uc: uncompensated (corrupted) scalar magnetometer measurements [nT]TL_coef: Tolles-Lawson coefficients (partially) used to createmag_uccor_fogm: corruption FOGM portion [nT]
MagNav.corrupt_mag — Methodcorrupt_mag(mag_c, flux;
dt = 0.1,
cor_sigma = 1.0,
cor_tau = 600.0,
cor_var = 1.0^2,
cor_drift = 0.001,
cor_perm_mag = 5.0,
cor_ind_mag = 5.0,
cor_eddy_mag = 0.5)Corrupt compensated (clean) magnetometer measurements with random, FOGM, drift, and Tolles-Lawson noise to create uncompensated (corrupted) scalar magnetometer measurements. FOGM is a First-order Gauss-Markov stochastic process.
Arguments:
mag_c: compensated (clean) scalar magnetometer measurements [nT]flux:MagVvector magnetometer measurement structdt: (optional) measurement time step [s]cor_sigma: (optional) corruption FOGM catch-all bias [nT]cor_tau: (optional) corruption FOGM catch-all time constant [s]cor_var: (optional) corruption measurement (white) noise variance [nT^2]cor_drift: (optional) corruption measurement linear drift [nT/s]cor_perm_mag: (optional) corruption permanent field TL coef std devcor_ind_mag: (optional) corruption induced field TL coef std devcor_eddy_mag: (optional) corruption eddy current TL coef std dev
Returns:
mag_uc: uncompensated (corrupted) scalar magnetometer measurements [nT]TL_coef: Tolles-Lawson coefficients (partially) used to createmag_uccor_fogm: corruption FOGM portion [nT]
MagNav.create_P0 — Functioncreate_P0(lat1 = deg2rad(45);
init_pos_sigma = 3.0,
init_alt_sigma = 0.001,
init_vel_sigma = 0.01,
init_att_sigma = deg2rad(0.01),
ha_sigma = 0.001,
a_hat_sigma = 0.01,
acc_sigma = 0.000245,
gyro_sigma = 0.00000000727,
fogm_sigma = 3.0,
vec_sigma = 1000.0,
vec_states::Bool = false,
fogm_state::Bool = true,
P0_TL = [])Create initial covariance matrix P0.
Arguments:
lat1: initial approximate latitude [rad]init_pos_sigma: (optional) initial position uncertainty [m]init_alt_sigma: (optional) initial altitude uncertainty [m]init_vel_sigma: (optional) initial velocity uncertainty [m/s]init_att_sigma: (optional) initial attitude uncertainty [rad]ha_sigma: (optional) barometer aiding altitude bias [m]a_hat_sigma: (optional) barometer aiding vertical accel bias [m/s^2]acc_sigma: (optional) accelerometer bias [m/s^2]gyro_sigma: (optional) gyroscope bias [rad/s]fogm_sigma: (optional) FOGM catch-all bias [nT]vec_sigma: (optional) vector magnetometer noise std devvec_states: (optional) if true, include vector magnetometer statesfogm_state: (optional) if true, include FOGM catch-all bias stateP0_TL: (optional) initial Tolles-Lawson covariance matrix
Returns:
P0: initial covariance matrix
MagNav.create_Qd — Functioncreate_Qd(dt = 0.1;
VRW_sigma = 0.000238,
ARW_sigma = 0.000000581,
baro_sigma = 1.0,
acc_sigma = 0.000245,
gyro_sigma = 0.00000000727,
fogm_sigma = 3.0,
vec_sigma = 1000.0,
TL_sigma = [],
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
vec_states::Bool = false,
fogm_state::Bool = true)Create the discrete time process/system noise matrix Qd.
Arguments:
dt: measurement time step [s]VRW_sigma: (optional) velocity random walk [m/s^2 /sqrt(Hz)]ARW_sigma: (optional) angular random walk [rad/s /sqrt(Hz)]baro_sigma: (optional) barometer bias [m]acc_sigma: (optional) accelerometer bias [m/s^2]gyro_sigma: (optional) gyroscope bias [rad/s]fogm_sigma: (optional) FOGM catch-all bias [nT]vec_sigma: (optional) vector magnetometer noise std devTL_sigma: (optional) Tolles-Lawson coefficients estimate std devbaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]vec_states: (optional) if true, include vector magnetometer statesfogm_state: (optional) if true, include FOGM catch-all bias state
Returns:
Qd: discrete time process/system noise matrix
MagNav.create_TL_A — Functioncreate_TL_A(flux::MagV, ind = trues(length(flux.x));
Bt = sqrt.(flux.x.^2+flux.y.^2+flux.z.^2)[ind],
terms = [:permanent,:induced,:eddy],
Bt_scale = 50000,
return_B = false)Create Tolles-Lawson A matrix using vector magnetometer measurements. Optionally returns the magnitude & derivatives of total field.
Arguments:
flux:MagVvector magnetometer measurement structind: (optional) selected data indicesBt: (optional) magnitude of vector magnetometer measurements or scalar magnetometer measurements for modified Tolles-Lawson [nT]terms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]return_B: (optional) if true, also returnBt&B_dot
Returns:
A: Tolles-LawsonAmatrixBt: ifreturn_B = true, magnitude of total field measurements [nT]B_dot: ifreturn_B = true, finite differences of total field vector [nT]
MagNav.create_TL_A — Methodcreate_TL_A(Bx, By, Bz;
Bt = sqrt.(Bx.^2+By.^2+Bz.^2),
terms = [:permanent,:induced,:eddy],
Bt_scale = 50000,
return_B = false)Create Tolles-Lawson A matrix using vector magnetometer measurements. Optionally returns the magnitude & derivatives of total field.
Arguments:
Bx,By,Bz: vector magnetometer measurements [nT]Bt: (optional) magnitude of vector magnetometer measurements or scalar magnetometer measurements for modified Tolles-Lawson [nT]terms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]return_B: (optional) if true, also returnBt&B_dot
Returns:
A: Tolles-LawsonAmatrixBt: ifreturn_B = true, magnitude of total field measurements [nT]B_dot: ifreturn_B = true, finite differences of total field vector [nT]
MagNav.create_TL_coef — Functioncreate_TL_coef(flux::MagV, B, ind = trues(length(flux.x));
Bt = sqrt.(flux.x.^2+flux.y.^2+flux.z.^2)[ind],
λ = 0,
terms = [:permanent,:induced,:eddy],
pass1 = 0.1,
pass2 = 0.9,
fs = 10.0,
pole::Int = 4,
trim::Int = 20,
Bt_scale = 50000,
return_var = false)Create Tolles-Lawson coefficients using vector & scalar magnetometer measurements with a bandpass, low-pass, or high-pass filter.
Arguments:
flux:MagVvector magnetometer measurement structB: scalar magnetometer measurements [nT]ind: (optional) selected data indicesBt: (optional) magnitude of vector magnetometer measurements or scalar magnetometer measurements for modified Tolles-Lawson [nT]λ: (optional) ridge parameterterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}pass1: (optional) first passband frequency [Hz]pass2: (optional) second passband frequency [Hz]fs: (optional) sampling frequency [Hz]pole: (optional) number of poles for Butterworth filtertrim: (optional) number of elements to trim after filteringBt_scale: (optional) scaling factor for induced & eddy current terms [nT]return_var: (optional) if true, also returnB_var
Returns:
coef: Tolles-Lawson coefficientsB_var: ifreturn_var = true, fit error variance
MagNav.create_TL_coef — Methodcreate_TL_coef(Bx, By, Bz, B;
Bt = sqrt.(Bx.^2+By.^2+Bz.^2),
λ = 0,
terms = [:permanent,:induced,:eddy],
pass1 = 0.1,
pass2 = 0.9,
fs = 10.0,
pole::Int = 4,
trim::Int = 20,
Bt_scale = 50000,
return_var = false)Create Tolles-Lawson coefficients using vector & scalar magnetometer measurements with a bandpass, low-pass, or high-pass filter.
Arguments:
Bx,By,Bz: vector magnetometer measurements [nT]B: scalar magnetometer measurements [nT]Bt: (optional) magnitude of vector magnetometer measurements or scalar magnetometer measurements for modified Tolles-Lawson [nT]λ: (optional) ridge parameterterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}pass1: (optional) first passband frequency [Hz]pass2: (optional) second passband frequency [Hz]fs: (optional) sampling frequency [Hz]pole: (optional) number of poles for Butterworth filtertrim: (optional) number of elements to trim after filteringBt_scale: (optional) scaling factor for induced & eddy current terms [nT]return_var: (optional) if true, also returnB_var
Returns:
coef: Tolles-Lawson coefficientsB_var: ifreturn_var = true, fit error variance
MagNav.create_XYZ0 — Functioncreate_XYZ0(mapS::Union{MapS,MapSd,MapS3D} = get_map(namad);
alt = 1000,
dt = 0.1,
t = 300,
v = 68,
ll1::Tuple = (),
ll2::Tuple = (),
N_waves::Int = 1,
attempts::Int = 10,
info::String = "Simulated data",
flight = 1,
line = 1,
year = 2023,
doy = 154,
mapV::MapV = get_map(emm720),
cor_sigma = 1.0,
cor_tau = 600.0,
cor_var = 1.0^2,
cor_drift = 0.001,
cor_perm_mag = 5.0,
cor_ind_mag = 5.0,
cor_eddy_mag = 0.5,
init_pos_sigma = 3.0,
init_alt_sigma = 0.001,
init_vel_sigma = 0.01,
init_att_sigma = deg2rad(0.01),
VRW_sigma = 0.000238,
ARW_sigma = 0.000000581,
baro_sigma = 1.0,
ha_sigma = 0.001,
a_hat_sigma = 0.01,
acc_sigma = 0.000245,
gyro_sigma = 0.00000000727,
fogm_sigma = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
save_h5::Bool = false,
xyz_h5::String = "xyz_data.h5",
silent::Bool = false)Create basic flight data. Assumes constant altitude (2D flight). No required arguments, though many are available to create custom data.
Trajectory Arguments:
mapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map structalt: (optional) altitude [m]dt: (optional) measurement time step [s]t: (optional) total flight time, ignored ifll2is set [s]v: (optional) approximate aircraft velocity [m/s]ll1: (optional) initial (lat,lon) point [deg]ll2: (optional) final (lat,lon) point [deg]N_waves: (optional) number of sine waves along pathattempts: (optional) maximum attempts at creating flight path onmapSinfo: (optional) flight data informationflight: (optional) flight numberline: (optional) line number, i.e., segment withinflightyear: (optional) yeardoy: (optional) day of yearmapV: (optional)MapVvector magnetic anomaly map structsave_h5: (optional) if true, savexyztoxyz_h5xyz_h5: (optional) path/name of flight data HDF5 file to save (.h5extension optional)
Compensated Measurement Corruption Arguments:
cor_var: (optional) corruption measurement (white) noise variance [nT^2]fogm_sigma: (optional) FOGM catch-all bias [nT]fogm_tau: (optional) FOGM catch-all time constant [s]
Uncompensated Measurement Corruption Arguments:
cor_sigma: (optional) corruption FOGM catch-all bias [nT]cor_tau: (optional) corruption FOGM catch-all time constant [s]cor_var: (optional) corruption measurement (white) noise variance [nT^2]cor_drift: (optional) corruption measurement linear drift [nT/s]cor_perm_mag: (optional) corruption permanent field TL coef std devcor_ind_mag: (optional) corruption induced field TL coef std devcor_eddy_mag: (optional) corruption eddy current TL coef std dev
INS Arguments:
init_pos_sigma: (optional) initial position uncertainty [m]init_alt_sigma: (optional) initial altitude uncertainty [m]init_vel_sigma: (optional) initial velocity uncertainty [m/s]init_att_sigma: (optional) initial attitude uncertainty [rad]VRW_sigma: (optional) velocity random walk [m/s^2 /sqrt(Hz)]ARW_sigma: (optional) angular random walk [rad/s /sqrt(Hz)]baro_sigma: (optional) barometer bias [m]ha_sigma: (optional) barometer aiding altitude bias [m]a_hat_sigma: (optional) barometer aiding vertical accel bias [m/s^2]acc_sigma: (optional) accelerometer bias [m/s^2]gyro_sigma: (optional) gyroscope bias [rad/s]baro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]
General Arguments:
silent: (optional) if true, no print outs
Returns:
xyz:XYZ0flight data struct
MagNav.create_flux — Functioncreate_flux(lat, lon, mapV::MapV = get_map(emm720);
Cnb = repeat(I(3),1,1,length(lat)),
alt = 1000,
dt = 0.1,
meas_var = 1.0^2,
fogm_sigma = 1.0,
fogm_tau = 600.0,
silent::Bool = false)Create compensated (clean) vector magnetometer measurements from a vector magnetic anomaly map.
Arguments:
lat: latitude [rad]lon: longitude [rad]mapV: (optional)MapVvector magnetic anomaly map structCnb: (optional) direction cosine matrix (body to navigation) [-]alt: (optional) altitude [m]dt: (optional) measurement time step [s]meas_var: (optional) measurement (white) noise variance [nT^2]fogm_sigma: (optional) FOGM catch-all bias [nT]fogm_tau: (optional) FOGM catch-all time constant [s]silent: (optional) if true, no print outs
Returns:
flux:MagVvector magnetometer measurement struct
MagNav.create_flux — Functioncreate_flux(path::Path, mapV::MapV = get_map(emm720);
meas_var = 1.0^2,
fogm_sigma = 1.0,
fogm_tau = 600.0,
silent::Bool = false)Create compensated (clean) vector magnetometer measurements from a vector magnetic anomaly map.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structmapV: (optional)MapVvector magnetic anomaly map structmeas_var: (optional) measurement (white) noise variance [nT^2]fogm_sigma: (optional) FOGM catch-all bias [nT]fogm_tau: (optional) FOGM catch-all time constant [s]silent: (optional) if true, no print outs
Returns:
flux:MagVvector magnetometer measurement struct
MagNav.create_informed_xyz — Methodcreate_informed_xyz(xyz::XYZ, ind, mapS::Union{MapS,MapSd,MapS3D},
use_mag::Symbol, use_vec::Symbol, TL_coef::Vector;
terms::Vector{Symbol} = [:permanent,:induced,:eddy],
disp_min = 100,
disp_max = 500,
Bt_disp = 50,
Bt_scale = 50000)Create knowledge-informed data from existing data. Given map information, an XYZ structure with magnetometer readings, and a fitted Tolles-Lawson model, this function creates a consistent, displaced XYZ structure representing what the use_mag and mag_1_c would have collected had the entire flight been laterally shifted by (disp_min,disp_max), assuming that the linear aircraft model is reasonably accurate. It makes use of a Taylor expansion of the map information to update the expected changes due to the alternative map location and due to the imputed aircraft field. The aircraft "noise" is then carried over into use_mag in the new XYZ data.
Arguments:
xyz:XYZflight data structind: selected data indicesmapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structuse_mag: uncompensated scalar magnetometer to use forytarget vector {:mag_1_uc, etc.}use_vec: vector magnetometer (fluxgate) to use for Tolles-LawsonAmatrix {:flux_a, etc.}TL_coef: Tolles-Lawson coefficientsterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy}disp_min: (optional) minimum trajectory displacement [m]disp_max: (optional) maximum trajectory displacement [m]Bt_disp: (optional) target total magnetic field magnitude displacement offset [nT]Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
xyz_disp:XYZflight data struct with displaced trajectory & modified magnetometer readings
MagNav.create_ins — Methodcreate_ins(traj::Traj;
init_pos_sigma = 3.0,
init_alt_sigma = 0.001,
init_vel_sigma = 0.01,
init_att_sigma = deg2rad(0.01),
VRW_sigma = 0.000238,
ARW_sigma = 0.000000581,
baro_sigma = 1.0,
ha_sigma = 0.001,
a_hat_sigma = 0.01,
acc_sigma = 0.000245,
gyro_sigma = 0.00000000727,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
save_h5::Bool = false,
ins_h5::String = "ins_data.h5")Creates an INS trajectory about a true trajectory. Propagates a 17-state Pinson error model to create INS errors.
Arguments:
traj:Trajtrajectory structinit_pos_sigma: (optional) initial position uncertainty [m]init_alt_sigma: (optional) initial altitude uncertainty [m]init_vel_sigma: (optional) initial velocity uncertainty [m/s]init_att_sigma: (optional) initial attitude uncertainty [rad]VRW_sigma: (optional) velocity random walk [m/s^2 /sqrt(Hz)]ARW_sigma: (optional) angular random walk [rad/s /sqrt(Hz)]baro_sigma: (optional) barometer bias [m]ha_sigma: (optional) barometer aiding altitude bias [m]a_hat_sigma: (optional) barometer aiding vertical accel bias [m/s^2]acc_sigma: (optional) accelerometer bias [m/s^2]gyro_sigma: (optional) gyroscope bias [rad/s]baro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]save_h5: (optional) if true, saveinstoins_h5ins_h5: (optional) path/name of INS data HDF5 file to save (.h5extension optional)
Returns:
ins:INSinertial navigation system struct
MagNav.create_mag_c — Functioncreate_mag_c(lat, lon, mapS::Union{MapS,MapSd,MapS3D} = get_map(namad);
alt = 1000,
dt = 0.1,
meas_var = 1.0^2,
fogm_sigma = 1.0,
fogm_tau = 600.0,
silent::Bool = false)Create compensated (clean) scalar magnetometer measurements from a scalar magnetic anomaly map.
Arguments:
lat: latitude [rad]lon: longitude [rad]mapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map structalt: (optional) altitude [m]dt: (optional) measurement time step [s]meas_var: (optional) measurement (white) noise variance [nT^2]fogm_sigma: (optional) FOGM catch-all bias [nT]fogm_tau: (optional) FOGM catch-all time constant [s]silent: (optional) if true, no print outs
Returns:
mag_c: compensated (clean) scalar magnetometer measurements [nT]
MagNav.create_mag_c — Functioncreate_mag_c(path::Path, mapS::Union{MapS,MapSd,MapS3D} = get_map(namad);
meas_var = 1.0^2,
fogm_sigma = 1.0,
fogm_tau = 600.0,
silent::Bool = false)Create compensated (clean) scalar magnetometer measurements from a scalar magnetic anomaly map.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structmapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map structmeas_var: (optional) measurement (white) noise variance [nT^2]fogm_sigma: (optional) FOGM catch-all bias [nT]fogm_tau: (optional) FOGM catch-all time constant [s]silent: (optional) if true, no print outs
Returns:
mag_c: compensated (clean) scalar magnetometer measurements [nT]
MagNav.create_model — Functioncreate_model(dt = 0.1, lat1 = deg2rad(45);
init_pos_sigma = 3.0,
init_alt_sigma = 0.001,
init_vel_sigma = 0.01,
init_att_sigma = deg2rad(0.01),
meas_var = 3.0^2,
VRW_sigma = 0.000238,
ARW_sigma = 0.000000581,
baro_sigma = 1.0,
ha_sigma = 0.001,
a_hat_sigma = 0.01,
acc_sigma = 0.000245,
gyro_sigma = 0.00000000727,
fogm_sigma = 3.0,
vec_sigma = 1000.0,
TL_sigma = [],
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
vec_states::Bool = false,
fogm_state::Bool = true,
P0_TL = [])Create a magnetic navigation filter model for use in an EKF or a MPF.
Arguments:
dt: measurement time step [s]lat1: initial approximate latitude [rad]init_pos_sigma: (optional) initial position uncertainty [m]init_alt_sigma: (optional) initial altitude uncertainty [m]init_vel_sigma: (optional) initial velocity uncertainty [m/s]init_att_sigma: (optional) initial attitude uncertainty [rad]meas_var: (optional) measurement (white) noise variance [nT^2]VRW_sigma: (optional) velocity random walk [m/s^2 /sqrt(Hz)]ARW_sigma: (optional) angular random walk [rad/s /sqrt(Hz)]baro_sigma: (optional) barometer bias [m]ha_sigma: (optional) barometer aiding altitude bias [m]a_hat_sigma: (optional) barometer aiding vertical accel bias [m/s^2]acc_sigma: (optional) accelerometer bias [m/s^2]gyro_sigma: (optional) gyroscope bias [rad/s]fogm_sigma: (optional) FOGM catch-all bias [nT]vec_sigma: (optional) vector magnetometer noise std devTL_sigma: (optional) Tolles-Lawson coefficients estimate std devbaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]vec_states: (optional) if true, include vector magnetometer statesfogm_state: (optional) if true, include FOGM catch-all bias stateP0_TL: (optional) initial Tolles-Lawson covariance matrix
Returns:
P0: initial covariance matrixQd: discrete time process/system noise matrixR: measurement (white) noise variance
MagNav.create_traj — Functioncreate_traj(mapS::Union{MapS,MapSd,MapS3D} = get_map(namad);
alt = 1000,
dt = 0.1,
t = 300,
v = 68,
ll1::Tuple = (),
ll2::Tuple = (),
N_waves::Int = 1,
attempts::Int = 10,
save_h5::Bool = false,
traj_h5::String = "traj_data.h5")Create Traj trajectory struct with a straight or sinusoidal flight path at constant altitude (2D flight).
Arguments:
mapS: (optional)MapS,MapSd, orMapS3Dscalar magnetic anomaly map structalt: (optional) altitude [m]dt: (optional) measurement time step [s]t: (optional) total flight time, ignored ifll2is set [s]v: (optional) approximate aircraft velocity [m/s]ll1: (optional) initial (lat,lon) point [deg]ll2: (optional) final (lat,lon) point [deg]N_waves: (optional) number of sine waves along pathattempts: (optional) maximum attempts at creating flight path onmapSsave_h5: (optional) if true, savetrajtotraj_h5traj_h5: (optional) path/name of trajectory data HDF5 file to save (.h5extension optional)
Returns:
traj:Trajtrajectory struct
MagNav.crlb — Methodcrlb(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, dt, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Cramér–Rao lower bound (CRLB) computed with classic Kalman Filter. Equations evaluated about true trajectory.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
P: non-linear covariance matrix
MagNav.crlb — Methodcrlb(traj::Traj, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Cramér–Rao lower bound (CRLB) computed with classic Kalman Filter. Equations evaluated about true trajectory.
Arguments:
traj:Trajtrajectory structitp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
P: non-linear covariance matrix
MagNav.dcm2euler — Functiondcm2euler(dcm, order::Symbol = :body2nav)Converts a DCM (direction cosine matrix) to yaw, pitch, and roll Euler angles. Yaw is synonymous with azimuth and heading here. There are 2 use cases:
- With
order = :body2nav, the provided DCM is assumed to rotate from the
body frame in the standard -roll, -pitch, -yaw sequence to the navigation frame. For example, if v1 is a 3x1 vector in the body frame [nose, right wing, down], then that vector rotated into the navigation frame [north, east, down] would be v2 = dcm * v1.
- With
order = :nav2body, the provided DCM is assumed to rotate from the
navigation frame in the standard yaw, pitch, roll sequence to the body frame. For example, if v1 is a 3x1 vector in the navigation frame [north, east, down], then that vector rotated into the body frame [nose, right wing, down] would be v2 = dcm * v1.
Reference: Titterton & Weston, Strapdown Inertial Navigation Technology, 2004, Section 3.6 (pg. 36-41 & 537).
Arguments:
dcm:3x3xNdirection cosine matrix [-]order: (optional) rotation order {:body2nav,:nav2body}
Returns:
roll: length-Nroll angle [rad], right-handed rotation about x-axispitch: length-Npitch angle [rad], right-handed rotation about y-axisyaw: length-Nyaw angle [rad], right-handed rotation about z-axis
MagNav.de2dlon — Methodde2dlon(de, lat)Convert east-west position (easting) difference to longitude difference.
Arguments:
de: east-west position (easting) difference [m]lat: nominal latitude [rad]
Returns:
dlon: longitude difference [rad]
MagNav.denorm_sets — Methoddenorm_sets(train_bias, train_scale, train, test)Denormalize (or destandardize) features (columns) of training & testing data.
Arguments:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalizedtest:N_testxNftesting data, normalized
Returns:
train:N_trainxNftraining data, denormalizedtest:N_testxNftesting data, denormalized
MagNav.denorm_sets — Methoddenorm_sets(train_bias, train_scale, train, val, test)Denormalize (or destandardize) features (columns) of training, validation, & testing data.
Arguments:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalizedval:N_valxNfvalidation data, normalizedtest:N_testxNftesting data, normalized
Returns:
train:N_trainxNftraining data, denormalizedval:N_valxNfvalidation data, denormalizedtest:N_testxNftesting data, denormalized
MagNav.denorm_sets — Methoddenorm_sets(train_bias, train_scale, train)Denormalize (or destandardize) features (columns) of training data.
Arguments:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalized
Returns:
train:N_trainxNftraining data, denormalized
MagNav.detrend — Functiondetrend(y, x = [eachindex(y);]; λ = 0, mean_only::Bool = false)Detrend signal (remove mean and optionally slope).
Arguments:
y: length-Nobserved data vectorx: (optional)NxNfinput data matrix (Nfis number of features)λ: (optional) ridge parametermean_only: (optional) if true, only remove mean (not slope)
Returns:
y: length-Nobserved data vector, detrended
MagNav.dlat2dn — Methoddlat2dn(dlat, lat)Convert latitude difference to north-south position (northing) difference.
Arguments:
dlat: latitude difference [rad]lat: nominal latitude [rad]
Returns:
dn: north-south position (northing) difference [m]
MagNav.dlon2de — Methoddlon2de(dlon, lat)Convert longitude difference to east-west position (easting) difference.
Arguments:
dlon: longitude difference [rad]lat: nominal latitude [rad]
Returns:
de: east-west position (easting) difference [m]
MagNav.dn2dlat — Methoddn2dlat(dn, lat)Convert north-south position (northing) difference to latitude difference.
Arguments:
dn: north-south position (northing) difference [m]lat: nominal latitude [rad]
Returns:
dlat: latitude difference [rad]
MagNav.downward_L — Methoddownward_L(map_map::Matrix, dx, dy, dz, α::Vector;
map_mask::BitMatrix = map_params(map_map)[2],
expand::Bool = true)Downward continuation using a sequence of regularization parameters to create a characteristic L-curve. The optimal regularization parameter is at a local minimum on the L-curve, which is a local maximum of curvature. The global maximum of curvature may or may not be the optimal regularization parameter.
Arguments:
map_map:nyxnx2D gridded map datadx: x-direction map step size [m]dy: y-direction map step size [m]dz: z-direction upward/downward continuation distance [m]α: (geometric) sequence of regularization parametersmap_mask: (optional)nyxnxmask for valid (not filled-in) map dataexpand: (optional) if true, expand map temporarily to reduce edge effects
Returns:
norms: L-infinity norm of difference between sequential D.C. solutions
MagNav.downward_L — Methoddownward_L(mapS::Union{MapS,MapSd,MapS3D}, alt, α::Vector;
expand::Bool = true)Downward continuation using a sequence of regularization parameters to create a characteristic L-curve. The optimal regularization parameter is at a local minimum on the L-curve, which is a local maximum of curvature. The global maximum of curvature may or may not be the optimal regularization parameter.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structalt: target downward continuation altitude [m]α: (geometric) sequence of regularization parametersexpand: (optional) if true, expand map temporarily to reduce edge effects
Returns:
norms: L-infinity norm of difference between sequential D.C. solutions
MagNav.ekf — Methodekf(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas, dt, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
der_mapS = nothing,
map_alt = 0)Extended Kalman filter (EKF) for airborne magnetic anomaly navigation.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementder_mapS: (optional) scalar map vertical derivative map interpolation function (f(lat,lon)or (f(lat,lon,alt))map_alt: (optional) map altitude [m]
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf — Methodekf(ins::INS, meas, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
der_mapS = map_itp(zeros(2,2),[-pi,pi],[-pi/2,pi/2]),
map_alt = 0)Extended Kalman filter (EKF) for airborne magnetic anomaly navigation.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementder_mapS: (optional) scalar map vertical derivative map interpolation function (f(lat,lon)or (f(lat,lon,alt))map_alt: (optional) map altitude [m]
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf_online — Methodekf_online(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas,
Bx, By, Bz, dt, itp_mapS, x0_TL, P0, Qd, R;
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
terms = [:permanent,:induced,:eddy,:bias],
Bt_scale = 50000)Extended Kalman filter (EKF) with online learning of Tolles-Lawson coefficients.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]Bx,By,Bz: vector magnetometer measurements [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x0_TL: initial Tolles-Lawson coefficient statesP0: initial covariance matrixQd: discrete time process/system noise matrixR: measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf_online — Methodekf_online(ins::INS, meas, flux::MagV, itp_mapS, x0_TL, P0, Qd, R;
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
terms = [:permanent,:induced,:eddy,:bias],
Bt_scale = 50000)Extended Kalman filter (EKF) with online learning of Tolles-Lawson coefficients.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]flux:MagVvector magnetometer measurement structitp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x0_TL: initial Tolles-Lawson coefficient statesP0: initial covariance matrixQd: discrete time process/system noise matrixR: measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}Bt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf_online_nn — Methodekf_online_nn(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas,
dt, itp_mapS, x_nn, m, y_norms, P0, Qd, R;
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
terms = [:permanent,:induced,:eddy])Extended Kalman filter (EKF) with online learning of neural network weights.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)m: neural network model, does not work with skip connectionsy_norms: length-2tuple ofynormalizations,(y_bias,y_scale)P0: initial covariance matrixQd: discrete time process/system noise matrixR: measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf_online_nn — Methodekf_online_nn(ins::INS, meas, itp_mapS, x_nn, m, y_norms, P0, Qd, R;
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Extended Kalman filter (EKF) with online learning of neural network weights.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)m: neural network model, does not work with skip connectionsy_norms: length-2tuple ofynormalizations,(y_bias,y_scale)P0: initial covariance matrixQd: discrete time process/system noise matrixR: measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
filt_res:FILTresfilter results struct
MagNav.ekf_online_nn_setup — Methodekf_online_nn_setup(x, y, m, y_norms; N_sigma::Int = 1000)Setup for extended Kalman filter (EKF) with online learning of neural network weights.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorm: neural network model, does not work with skip connectionsy_norms: tuple ofynormalizations, i.e.,(y_bias,y_scale)N_sigma: (optional) number of neural network weights sets to use to createnn_sigma
Returns:
P0_nn: initial neural network weights covariance matrixnn_sigma: initial neural network weights estimate std dev
MagNav.ekf_online_setup — Functionekf_online_setup(flux::MagV, meas, ind = trues(length(meas));
Bt = sqrt.(flux.x.^2+flux.y.^2+flux.z.^2)[ind],
λ = 0.025,
terms = [:permanent,:induced,:eddy,:bias],
pass1 = 0.1,
pass2 = 0.9,
fs = 10.0,
pole::Int = 4,
trim::Int = 20,
N_sigma::Int = 100,
Bt_scale = 50000)Setup for extended Kalman filter (EKF) with online learning of Tolles-Lawson coefficients.
Arguments:
flux:MagVvector magnetometer measurement structmeas: scalar magnetometer measurement [nT]ind: selected data indicesBt: (optional) magnitude of vector magnetometer measurements or scalar magnetometer measurements for modified Tolles-Lawson [nT]λ: (optional) ridge parameterterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}pass1: (optional) first passband frequency [Hz]pass2: (optional) second passband frequency [Hz]fs: (optional) sampling frequency [Hz]pole: (optional) number of poles for Butterworth filtertrim: (optional) number of elements to trim after filteringN_sigma: (optional) number of Tolles-Lawson coefficient sets to use to createTL_sigmaBt_scale: (optional) scaling factor for induced & eddy current terms [nT]
Returns:
x0_TL: initial Tolles-Lawson coefficient statesP0_TL: initial Tolles-Lawson covariance matrixTL_sigma: Tolles-Lawson coefficients estimate std dev
MagNav.elasticnet_fit — Functionelasticnet_fit(x, y, α::Real = 0.99, no_norm = falses(size(x,2));
λ::Real = -1,
data_norms::Tuple = (zeros(1,1),zeros(1,1),[0.0],[0.0]),
l_segs::Vector = [length(y)],
silent::Bool = false)Fit an elastic net (ridge regression and/or Lasso) model to data.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorα: (optional) ridge regression (α=0) vs Lasso (α=1) balancing parameter {0:1}no_norm: (optional) length-NfBoolean indices of features to not be normalizedλ: (optional) elastic net parameter,-1to ignore & determine with cross-validationdata_norms: (optional) length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
model: length-2tuple of elastic net-based model, (length-Nfcoefficients, bias)data_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.err_segs — Methoderr_segs(y_hat, y, l_segs; silent::Bool = true)Remove mean error from multiple individual flight lines within larger dataset.
Arguments:
y_hat: length-Nprediction vectory: length-Ntarget vectorl_segs: length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
err: length-Nmean-corrected (per line) error
MagNav.euler2dcm — Functioneuler2dcm(roll, pitch, yaw, order::Symbol = :body2nav)Converts a (Euler) roll-pitch-yaw (X-Y-Z) right-handed body to navigation frame rotation (or the opposite rotation), to a DCM (direction cosine matrix). Yaw is synonymous with azimuth and heading here. If frame 1 is rotated to frame 2, then the returned DCM, when pre-multiplied, rotates a vector in frame 1 into frame 2. There are 2 use cases:
- With
order = :body2nav, the body frame is rotated in the standard
-roll, -pitch, -yaw sequence to the navigation frame. For example, if v1 is a 3x1 vector in the body frame [nose, right wing, down], then that vector rotated into the navigation frame [north, east, down] would be v2 = dcm * v1.
- With
order = :nav2body, the navigation frame is rotated in the standard
yaw, pitch, roll sequence to the body frame. For example, if v1 is a 3x1 vector in the navigation frame [north, east, down], then that vector rotated into the body frame [nose, right wing, down] would be v2 = dcm * v1.
Reference: Titterton & Weston, Strapdown Inertial Navigation Technology, 2004, Section 3.6 (pg. 36-41 & 537).
Arguments:
roll: length-Nroll angle [rad], right-handed rotation about x-axispitch: length-Npitch angle [rad], right-handed rotation about y-axisyaw: length-Nyaw angle [rad], right-handed rotation about z-axisorder: (optional) rotation order {:body2nav,:nav2body}
Returns:
dcm:3x3xNdirection cosine matrix [-]
MagNav.eval_crlb — Methodeval_crlb(traj::Traj, crlb_P::Array)Extract Cramér–Rao lower bound (CRLB) results.
Arguments:
traj:Trajtrajectory structcrlb_P: Cramér–Rao lower bound non-linear covariance matrix
Returns:
crlb_out:CRLBoutCramér–Rao lower bound extracted output struct
MagNav.eval_filt — Methodeval_filt(traj::Traj, ins::INS, filt_res::FILTres)Extract filter results.
Arguments:
traj:Trajtrajectory structins:INSinertial navigation system structfilt_res:FILTresfilter results struct
Returns:
filt_out:FILToutfilter extracted output struct
MagNav.eval_gsa — Functioneval_gsa(m::Chain, x, n::Int = min(10000,size(x,1)))Global sensitivity analysis (GSA) with the Morris Method.
Reference: https://book.sciml.ai/notes/17-GlobalSensitivityAnalysis/
Reference: https://docs.sciml.ai/GlobalSensitivity/stable/methods/morris/
Arguments:
m: neural network modelx:NxNfdata matrix (Nfis number of features)n: (optional) number of samples (instances) to use for explanation
Returns:
means: means of elementary effects
MagNav.eval_ins — Methodeval_ins(traj::Traj, ins::INS)Extract INS results.
Arguments:
traj:Trajtrajectory structins:INSinertial navigation system struct
Returns:
ins_out:INSoutinertial navigation system extracted output struct
MagNav.eval_results — Methodeval_results(traj::Traj, ins::INS, filt_res::FILTres, crlb_P::Array)Extract CRLB, INS, & filter results.
Arguments:
traj:Trajtrajectory structins:INSinertial navigation system structfilt_res:FILTresfilter results structcrlb_P: Cramér–Rao lower bound non-linear covariance matrix
Returns:
crlb_out:CRLBoutCramér–Rao lower bound extracted output structins_out:INSoutinertial navigation system extracted output structfilt_out:FILToutfilter extracted output struct
MagNav.eval_shapley — Functioneval_shapley(m::Chain, x, features::Vector{Symbol},
N::Int = min(10000,size(x,1)),
num_mc::Int = 10)Compute stochastic Shapley effects for global feature importance.
Reference: https://nredell.github.io/ShapML.jl/dev/#Examples-1
Arguments:
m: neural network modelx:NxNfdata matrix (Nfis number of features)features: length-Nffeature vector (including components of TLA, etc.)N: (optional) number of samples (instances) to use for explanationnum_mc: (optional) number of Monte Carlo simulations
Returns:
df_shap: DataFrame of Shapley effectsbaseline_shap: intercept of Shapley effects
MagNav.fdm — Methodfdm(x::Vector; scheme::Symbol = :central)Finite difference method (FDM) applied to x.
Arguments:
x: data vectorscheme: (optional) finite difference method scheme usedbackward: 1st derivative 1st-order backward differenceforward: 1st derivative 1st-order forward differencecentral: 1st derivative 2nd-order central differencebackward2: 1st derivative 2nd-order backward differenceforward2: 1st derivative 2nd-order forward differencefourth: 4th derivative central difference
Returns:
dif: vector of finite differences (length ofx)
MagNav.filter_events! — Methodfilter_events!(flight::Symbol, df_event::DataFrame;
keyword::String = "",
tt_lim::Tuple = ())Filter a DataFrame of in-flight events to only contain relevant events.
Arguments:
flight: flight name (e.g.,:Flt1001)df_event: lookup table (DataFrame) of in-flight events
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
tt | Real | time of event [s] |
event | String | event description |
keyword: (optional) keyword to search within events, case insensitivett_lim: (optional) length-2start & end time limits (inclusive) [s]
Returns:
nothing:df_eventis filtered
MagNav.filter_events — Methodfilter_events(flight::Symbol, df_event::DataFrame;
keyword::String = "",
tt_lim::Tuple = ())Filter a DataFrame of in-flight events to only contain relevant events.
Arguments:
flight: flight name (e.g.,:Flt1001)df_event: lookup table (DataFrame) of in-flight events
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
tt | Real | time of event [s] |
event | String | event description |
keyword: (optional) keyword to search within events, case insensitivett_lim: (optional) length-2start & end time limits (inclusive) [s]
Returns:
df_event: lookup table (DataFrame) of in-flight events, filtered
MagNav.fogm — Methodfogm(sigma, tau, dt, N)First-order Gauss-Markov stochastic process. Represents unmeasureable time-correlated errors.
Arguments:
sigma: FOGM catch-all biastau: FOGM catch-all time constant [s]dt: measurement time step [s]N: number of samples (instances)
Returns:
x: FOGM data
MagNav.get_Axy — Functionget_Axy(lines, df_line::DataFrame,
df_flight::DataFrame, df_map::DataFrame,
features_setup::Vector{Symbol} = [:mag_1_uc,:TL_A_flux_a];
features_no_norm::Vector{Symbol} = Symbol[],
y_type::Symbol = :d,
use_mag::Symbol = :mag_1_uc,
use_mag_c::Symbol = :mag_1_c,
use_vec::Symbol = :flux_a,
terms = [:permanent,:induced,:eddy],
terms_A = [:permanent,:induced,:eddy,:bias],
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
bpf_mag::Bool = false,
reorient_vec::Bool = false,
l_window::Int = -1,
mod_TL::Bool = false,
map_TL::Bool = false,
return_B::Bool = false,
silent::Bool = true)Get "external" Tolles-Lawson A matrix, x data matrix, & y target vector from multiple flight lines, possibly multiple flights. Optionally return Bt & B_dot used to create the "external" Tolles-Lawson A matrix.
Arguments:
lines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line, only used for y_type = :b, :c |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
features_setup: vector of features to includefeatures_no_norm: (optional) vector of features to not normalizey_type: (optional)ytarget type:a= anomaly field #1, compensated tail stinger total field scalar magnetometer measurements:b= anomaly field #2, interpolatedmagnetic anomaly mapvalues:c= aircraft field #1, difference between uncompensated cabin total field scalar magnetometer measurements and interpolatedmagnetic anomaly mapvalues:d= aircraft field #2, difference between uncompensated cabin and compensated tail stinger total field scalar magnetometer measurements:e= BPF'd total field, bandpass filtered uncompensated cabin total field scalar magnetometer measurements
use_mag: (optional) uncompensated scalar magnetometer to use forytarget vector {:mag_1_uc, etc.}, only used fory_type = :c, :d, :euse_mag_c: (optional) compensated scalar magnetometer to use forytarget vector {:mag_1_c, etc.}, only used fory_type = :a, :duse_vec: (optional) vector magnetometer (fluxgate) to use for "external" Tolles-LawsonAmatrix {:flux_a, etc.}terms: (optional) Tolles-Lawson terms to use forAwithinxdata matrix {:permanent,:induced,:eddy,:bias}terms_A: (optional) Tolles-Lawson terms to use for "external" Tolles-LawsonAmatrix {:permanent,:induced,:eddy,:bias}sub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsbpf_mag: (optional) if true, bpf scalar magnetometer measurements inxdata matrixreorient_vec: (optional) if true, align vector magnetometer measurements with body framel_window: (optional) trim data byN % l_window,-1to ignoremod_TL: (optional) if true, create modified "external" Tolles-LawsonAmatrix withuse_magmap_TL: (optional) if true, create map-based "external" Tolles-LawsonAmatrixreturn_B: (optional) if true, also returnBt&B_dotsilent: (optional) if true, no print outs
Returns:
A:NxN_TL"external" Tolles-LawsonAmatrix (N_TLis number of Tolles-Lawson coefficients)x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorno_norm: length-NfBoolean indices of features to not be normalizedfeatures: length-Nffeature vector (including components of TLA, etc.)l_segs: length-N_linesvector of lengths oflines, sum(l_segs) =NBt: ifreturn_B = true, length-Nmagnitude of total field measurements used to createA[nT]B_dot: ifreturn_B = true,Nx3finite differences of total field vector used to createA[nT]
MagNav.get_XYZ — Methodget_XYZ(flight::Symbol, df_flight::DataFrame;
tt_sort::Bool = true,
reorient_vec::Bool = false,
silent::Bool = false)Get XYZ flight data from saved HDF5 file via DataFrame lookup.
Arguments:
flight: flight name (e.g.,:Flt1001)df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
tt_sort: (optional) if true, sort data by time (instead of line)reorient_vec: (optional) if true, align vector magnetometer measurements with body framesilent: (optional) if true, no print outs
Returns:
xyz:XYZflight data struct
MagNav.get_XYZ0 — Functionget_XYZ0(xyz_file::String,
traj_field::Symbol = :traj,
ins_field::Symbol = :ins_data;
info::String = splitpath(xyz_file)[end],
flight = 1,
line = 1,
year = 2023,
doy = 154,
dt = 0.1,
tt_sort::Bool = true,
silent::Bool = false)Get the minimum dataset required for MagNav from saved CSV, HDF5, or MAT file. Not all fields within the XYZ0 flight data struct are required. The minimum data required in the CSV, HDF5, or MAT file includes:
lat,lon,alt(position)mag_1_cORmag_1_uc(scalar magnetometer measurements)
All other fields can be computed/simulated.
If a CSV or HDF5 file is provided, the possible columns/fields in the file are:
| Field | Type | Description |
|---|---|---|
dt | scalar | measurement time step [s] |
tt | vector | time [s] |
lat | vector | latitude [deg] |
lon | vector | longitude [deg] |
alt | vector | altitude [m] |
vn | vector | north velocity [m/s] |
ve | vector | east velocity [m/s] |
vd | vector | down velocity [m/s] |
fn | vector | north specific force [m/s] |
fe | vector | east specific force [m/s] |
fd | vector | down specific force [m/s] |
Cnb | 3x3xN | direction cosine matrix (body to navigation) [-], not valid for CSV file (use roll, pitch, & yaw instead) |
roll | vector | roll [deg] |
pitch | vector | pitch [deg] |
yaw | vector | yaw [deg] |
ins_dt | scalar | INS measurement time step [s] |
ins_tt | vector | INS time [s] |
ins_lat | vector | INS latitude [deg] |
ins_lon | vector | INS longitude [deg] |
ins_alt | vector | INS altitude [m] |
ins_vn | vector | INS north velocity [m/s] |
ins_ve | vector | INS east velocity [m/s] |
ins_vd | vector | INS down velocity [m/s] |
ins_fn | vector | INS north specific force [m/s] |
ins_fe | vector | INS east specific force [m/s] |
ins_fd | vector | INS down specific force [m/s] |
ins_Cnb | 3x3xN | INS direction cosine matrix (body to navigation) [-], not valid for CSV file (use ins_roll, ins_pitch, & ins_yaw instead) |
ins_roll | vector | INS roll [deg] |
ins_pitch | vector | INS pitch [deg] |
ins_yaw | vector | INS yaw [deg] |
ins_P | 17x17xN | INS covariance matrix, only relevant for simulated data, otherwise zeros [-], not valid for CSV file |
flux_a_x | vector | Flux A x-direction magnetic field [nT] |
flux_a_y | vector | Flux A y-direction magnetic field [nT] |
flux_a_z | vector | Flux A z-direction magnetic field [nT] |
flux_a_t | vector | Flux A total magnetic field [nT] |
flight | vector | flight number(s) |
line | vector | line number(s), i.e., segments within flight |
year | vector | year |
doy | vector | day of year |
diurnal | vector | measured diurnal, i.e., temporal variations or space weather effects [nT] |
igrf | vector | International Geomagnetic Reference Field (IGRF), i.e., core field [nT] |
mag_1_c | vector | Mag 1 compensated (clean) scalar magnetometer measurements [nT] |
mag_1_uc | vector | Mag 1 uncompensated (corrupted) scalar magnetometer measurements [nT] |
If a MAT file is provided, the above fields may also be provided, but the non-INS fields should be within the specified traj_field MAT struct and the INS fields should be within the specified ins_field MAT struct and without ins_ prefixes. This is the standard way the MATLAB-companion outputs data.
Arguments:
xyz_file: path/name of flight data CSV, HDF5, or MAT file (.csv,.h5, or.matextension required)traj_field: (optional) trajectory struct field within MAT file to use, not relevant for CSV or HDF5 fileins_field: (optional) INS struct field within MAT file to use,:noneif unavailable, not relevant for CSV or HDF5 filett_sort: (optional) if true, sort data by time (instead of line)silent: (optional) if true, no print outs
If not provided in xyz_file:
info: (optional) flight data informationflight: (optional) flight numberline: (optional) line number, i.e., segment withinflightyear: (optional) yeardoy: (optional) day of yeardt: (optional) measurement time step [s]
Returns:
xyz:XYZ0flight data struct
MagNav.get_XYZ1 — Functionget_XYZ1(xyz_file::String,
traj_field::Symbol = :traj,
ins_field::Symbol = :ins_data;
info::String = splitpath(xyz_file)[end],
flight = 1,
line = 1,
year = 2023,
doy = 154,
dt = 0.1,
tt_sort::Bool = true,
silent::Bool = false)Get the minimum dataset required for MagNav from saved CSV, HDF5, or MAT file. Not all fields within the XYZ1 flight data struct are required. The minimum data required in the CSV, HDF5, or MAT file includes:
lat,lon,alt(position)mag_1_cORmag_1_uc(scalar magnetometer measurements)
All other fields can be computed/simulated, except flux_b, mag_2_c, mag_3_c, mag_2_uc, mag_3_uc, aux_1, aux_2, and aux_3.
If a CSV or HDF5 file is provided, the possible columns/fields in the file are:
| Field | Type | Description |
|---|---|---|
dt | scalar | measurement time step [s] |
tt | vector | time [s] |
lat | vector | latitude [deg] |
lon | vector | longitude [deg] |
alt | vector | altitude [m] |
vn | vector | north velocity [m/s] |
ve | vector | east velocity [m/s] |
vd | vector | down velocity [m/s] |
fn | vector | north specific force [m/s] |
fe | vector | east specific force [m/s] |
fd | vector | down specific force [m/s] |
Cnb | 3x3xN | direction cosine matrix (body to navigation) [-], not valid for CSV file (use roll, pitch, & yaw instead) |
roll | vector | roll [deg] |
pitch | vector | pitch [deg] |
yaw | vector | yaw [deg] |
ins_dt | scalar | INS measurement time step [s] |
ins_tt | vector | INS time [s] |
ins_lat | vector | INS latitude [deg] |
ins_lon | vector | INS longitude [deg] |
ins_alt | vector | INS altitude [m] |
ins_vn | vector | INS north velocity [m/s] |
ins_ve | vector | INS east velocity [m/s] |
ins_vd | vector | INS down velocity [m/s] |
ins_fn | vector | INS north specific force [m/s] |
ins_fe | vector | INS east specific force [m/s] |
ins_fd | vector | INS down specific force [m/s] |
ins_Cnb | 3x3xN | INS direction cosine matrix (body to navigation) [-], not valid for CSV file (use ins_roll, ins_pitch, & ins_yaw instead) |
ins_roll | vector | INS roll [deg] |
ins_pitch | vector | INS pitch [deg] |
ins_yaw | vector | INS yaw [deg] |
ins_P | 17x17xN | INS covariance matrix, only relevant for simulated data, otherwise zeros [-], not valid for CSV file |
flux_a_x | vector | Flux A x-direction magnetic field [nT] |
flux_a_y | vector | Flux A y-direction magnetic field [nT] |
flux_a_z | vector | Flux A z-direction magnetic field [nT] |
flux_a_t | vector | Flux A total magnetic field [nT] |
flux_b_x | vector | Flux B x-direction magnetic field [nT] |
flux_b_y | vector | Flux B y-direction magnetic field [nT] |
flux_b_z | vector | Flux B z-direction magnetic field [nT] |
flux_b_t | vector | Flux B total magnetic field [nT] |
flight | vector | flight number(s) |
line | vector | line number(s), i.e., segments within flight |
year | vector | year |
doy | vector | day of year |
diurnal | vector | measured diurnal, i.e., temporal variations or space weather effects [nT] |
igrf | vector | International Geomagnetic Reference Field (IGRF), i.e., core field [nT] |
mag_1_c | vector | Mag 1 compensated (clean) scalar magnetometer measurements [nT] |
mag_2_c | vector | Mag 2 compensated (clean) scalar magnetometer measurements [nT] |
mag_3_c | vector | Mag 3 compensated (clean) scalar magnetometer measurements [nT] |
mag_1_uc | vector | Mag 1 uncompensated (corrupted) scalar magnetometer measurements [nT] |
mag_2_uc | vector | Mag 2 uncompensated (corrupted) scalar magnetometer measurements [nT] |
mag_3_uc | vector | Mag 3 uncompensated (corrupted) scalar magnetometer measurements [nT] |
aux_1 | vector | flexible-use auxiliary data 1 |
aux_2 | vector | flexible-use auxiliary data 2 |
aux_3 | vector | flexible-use auxiliary data 3 |
If a MAT file is provided, the above fields may also be provided, but the non-INS fields should be within the specified traj_field MAT struct and the INS fields should be within the specified ins_field MAT struct and without ins_ prefixes. This is the standard way the MATLAB-companion outputs data.
Arguments:
xyz_file: path/name of flight data CSV, HDF5, or MAT file (.csv,.h5, or.matextension required)traj_field: (optional) trajectory struct field within MAT file to use, not relevant for CSV or HDF5 fileins_field: (optional) INS struct field within MAT file to use,:noneif unavailable, not relevant for CSV or HDF5 filett_sort: (optional) if true, sort data by time (instead of line)silent: (optional) if true, no print outs
If not provided in xyz_file:
info: (optional) flight data informationflight: (optional) flight numberline: (optional) line number, i.e., segment withinflightyear: (optional) yeardoy: (optional) day of yeardt: (optional) measurement time step [s]
Returns:
xyz:XYZ1flight data struct
MagNav.get_XYZ20 — Methodget_XYZ20(xyz_160_h5::String, xyz_h5::String;
info::String = splitpath(xyz_160_h5)[end] * " & " * splitpath(xyz_h5)[end],
silent::Bool = false)Get 160 Hz (partial) XYZ20 flight data from saved HDF5 file and combine with 10 Hz XYZ20 flight data from another saved HDF5 file. Data is time sorted to ensure data is aligned.
Arguments:
xyz_160_h5: path/name of 160 Hz flight data HDF5 file (.h5extension optional)xyz_h5: path/name of 10 Hz flight data HDF5 file (.h5extension optional)info: (optional) flight data informationsilent: (optional) if true, no print outs
Returns:
xyz:XYZ20flight data struct
MagNav.get_XYZ20 — Methodget_XYZ20(xyz_h5::String;
info::String = splitpath(xyz_h5)[end],
tt_sort::Bool = true,
silent::Bool = false)Get XYZ20 flight data from saved HDF5 file. Based on 2020 SGL data fields.
Arguments:
xyz_h5: path/name of flight data HDF5 file (.h5extension optional)info: (optional) flight data informationtt_sort: (optional) if true, sort data by time (instead of line)silent: (optional) if true, no print outs
Returns:
xyz:XYZ20flight data struct
MagNav.get_XYZ21 — Methodget_XYZ21(xyz_h5::String;
info::String = splitpath(xyz_h5)[end],
tt_sort::Bool = true,
silent::Bool = false)Get XYZ21 flight data from saved HDF5 file. Based on 2021 SGL data fields.
Arguments:
xyz_h5: path/name of flight data HDF5 file (.h5extension optional)info: (optional) flight data informationtt_sort: (optional) if true, sort data by time (instead of line)silent: (optional) if true, no print outs
Returns:
xyz:XYZ21flight data struct
MagNav.get_autocor — Functionget_autocor(x::Vector, dt = 0.1, dt_max = 300.0)Get autocorrelation of data (e.g., actual - expected measurements).
Arguments:
x: data vectordt: (optional) measurement time step [s]dt_max: (optional) maximum time step to evaluate [s]
Returns:
sigma: standard deviationtau: autocorrelation decay to e^-1 ofx[s]
MagNav.get_bpf — Methodget_bpf(; pass1 = 0.1, pass2 = 0.9, fs = 10.0, pole::Int = 4)Create a Butterworth bandpass (or low-pass or high-pass) filter object. Set pass1 = -1 for low-pass filter or pass2 = -1 for high-pass filter.
Arguments:
pass1: (optional) first passband frequency [Hz]pass2: (optional) second passband frequency [Hz]fs: (optional) sampling frequency [Hz]pole: (optional) number of poles for Butterworth filter
Returns:
bpf: filter object
MagNav.get_cached_map — Methodget_cached_map(map_cache::Map_Cache, lat::Real, lon::Real, alt::Real;
silent::Bool = false)Get cached map at specific location.
Arguments:
map_cache:Map_Cachemap cache structlat: latitude [rad]lon: longitude [rad]alt: altitude [m]silent: (optional) if true, no print outs
Returns:
itp_mapS: scalar map interpolation function (f(lat,lon)atalt)
MagNav.get_comp_params — Functionget_comp_params(comp_params_bson::String, silent::Bool = false)Get aeromagnetic compensation parameters from saved BSON file.
Arguments:
comp_params_bson: path/name of aeromagnetic compensation parameters BSON file (.bsonextension optional)silent: (optional) if true, no print outs
Returns:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
MagNav.get_flux — Functionget_flux(flux_file::String,
use_vec::Symbol = :flux_a,
field::Symbol = :traj)Get vector magnetometer data from saved CSV, HDF5, or MAT file.
If a CSV or HDF5 file is provided, the possible columns/fields in the file are:
| Field | Type | Description |
|---|---|---|
use_vec*_x | vector | x-direction magnetic field [nT] |
use_vec*_y | vector | y-direction magnetic field [nT] |
use_vec*_z | vector | z-direction magnetic field [nT] |
use_vec*_t | vector | total magnetic field [nT], optional |
If a MAT file is provided, the above fields may also be provided, but they should be within the specified field MAT struct. This is the standard way the MATLAB-companion outputs data.
Arguments:
flux_file: path/name of vector magnetometer data CSV, HDF5, or MAT file (.csv,.h5, or.matextension required)use_vec: (optional) vector magnetometer (fluxgate) to usefield: (optional) struct field within MAT file to use, not relevant for CSV or HDF5 file
Returns:
flux:MagVvector magnetometer measurement struct
MagNav.get_igrf — Functionget_igrf(xyz::XYZ, ind = trues(xyz.traj.N);
frame::Symbol = :body,
norm_igrf::Bool = false,
check_xyz::Bool = true)Get the IGRF Earth vector in the body or navigation frame given an XYZ flight data struct containing trajectory information, valid indices, a start date in IGRF time (years since 0 CE), and reference frame.
Arguments:
xyz:XYZflight data structind: (optional) selected data indicesframe: (optional) desired reference frame {:body,:nav}norm_igrf: (optional) if true, normalizeigrf_veccheck_xyz: (optional) if true, cross-check withigrffield inxyz
Returns:
igrf_vec: length-Nstacked vector of3IGRF coordinates inframe
MagNav.get_ind — Methodget_ind(xyz::XYZ, lines, df_line::DataFrame;
splits = (1),
l_window::Int = -1)Get BitVector of selected data indices for further analysis via DataFrame lookup.
Arguments:
xyz:XYZflight data structlines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
splits: (optional) data splits, must sum to 1l_window: (optional) trim data byN % l_window,-1to ignore
Returns:
ind: BitVector (or tuple of BitVector) of selected data indices
MagNav.get_ind — Methodget_ind(xyz::XYZ, line::Real, df_line::DataFrame;
splits = (1),
l_window::Int = -1)Get BitVector of indices for further analysis via DataFrame lookup.
Arguments:
xyz:XYZflight data structline: line numberdf_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
splits: (optional) data splits, must sum to 1l_window: (optional) trim data byN % l_window,-1to ignore
Returns:
ind: BitVector (or tuple of BitVector) of selected data indices
MagNav.get_ind — Methodget_ind(xyz::XYZ;
ind = trues(xyz.traj.N),
lines = (),
tt_lim = (),
splits = (1))Get BitVector of indices for further analysis from specified indices (subset), lines, and/or time range. Any or all of these may be used. Defaults to use all indices, lines, and times.
Arguments:
xyz:XYZflight data structind: (optional) selected data indiceslines: (optional) selected line number(s)tt_lim: (optional) end time limit or length-2start & end time limits (inclusive) [s]splits: (optional) data splits, must sum to 1
Returns:
ind: BitVector (or tuple of BitVector) of selected data indices
MagNav.get_ind — Methodget_ind(tt::Vector, line::Vector;
ind = trues(length(tt)),
lines = (),
tt_lim = (),
splits = (1))Get BitVector of indices for further analysis from specified indices (subset), lines, and/or time range. Any or all of these may be used. Defaults to use all indices, lines, and times.
Arguments:
tt: time [s]line: line number(s)ind: (optional) selected data indiceslines: (optional) selected line number(s)tt_lim: (optional) end time limit or length-2start & end time limits (inclusive) [s]splits: (optional) data splits, must sum to 1
Returns:
ind: BitVector (or tuple of BitVector) of selected data indices
MagNav.get_ins — Functionget_ins(xyz::XYZ, ind = trues(xyz.traj.N);
N_zero_ll::Int = 0,
t_zero_ll::Real = 0,
err::Real = 0.0)Get inertial navigation system data at specific indices, possibly zeroed.
Arguments:
xyz:XYZflight data structind: (optional) selected data indicesN_zero_ll: (optional) number of samples (instances) to zero INS lat/lon to truth (xyz.traj)t_zero_ll: (optional) length of time to zero INS lat/lon to truth (xyz.traj), overwritesN_zero_llerr: (optional) additional position error [m]
Returns:
ins:INSinertial navigation system struct atind
MagNav.get_ins — Functionget_ins(ins_file::String, field::Symbol = :ins_data;
dt = 0.1,
tt_sort::Bool = true,
silent::Bool = false)Get inertial navigation system data from saved CSV, HDF5, or MAT file. The only required fields are ins_lat, ins_lon, and ins_alt (position).
If a CSV or HDF5 file is provided, the possible columns/fields in the file are:
| Field | Type | Description |
|---|---|---|
ins_dt | scalar | INS measurement time step [s] |
ins_tt | vector | INS time [s] |
ins_lat | vector | INS latitude [deg] |
ins_lon | vector | INS longitude [deg] |
ins_alt | vector | INS altitude [m] |
ins_vn | vector | INS north velocity [m/s] |
ins_ve | vector | INS east velocity [m/s] |
ins_vd | vector | INS down velocity [m/s] |
ins_fn | vector | INS north specific force [m/s] |
ins_fe | vector | INS east specific force [m/s] |
ins_fd | vector | INS down specific force [m/s] |
ins_Cnb | 3x3xN | INS direction cosine matrix (body to navigation) [-], not valid for CSV file (use ins_roll, ins_pitch, & ins_yaw instead) |
ins_roll | vector | INS roll [deg] |
ins_pitch | vector | INS pitch [deg] |
ins_yaw | vector | INS yaw [deg] |
ins_P | 17x17xN | INS covariance matrix, only relevant for simulated data, otherwise zeros [-], not valid for CSV file |
If a MAT file is provided, the above fields may also be provided, but they should be within the specified field MAT struct and without ins_ prefixes. This is the standard way the MATLAB-companion outputs data.
Arguments:
ins_file: path/name of INS data CSV, HDF5, or MAT file (.csv,.h5, or.matextension required)field: (optional) struct field within MAT file to use, not relevant for CSV or HDF5 filedt: (optional) measurement time step [s], only used if not inins_filesilent: (optional) if true, no print outs
Returns:
ins:INSinertial navigation system struct
MagNav.get_map — Functionget_map(map_name::Symbol, df_map::DataFrame,
map_field::Symbol = :map_data;
map_info::String = "$map_name",
map_units::Symbol = :rad,
file_units::Symbol = :deg,
flip_map::Bool = false)Get map data from saved HDF5 or MAT file or folder containing CSV files via DataFrame lookup. Maps are typically saved in :deg units, while :rad is used internally.
Arguments:
map_name: name of magnetic anomaly mapdf_map: lookup table (DataFrame) of map data HDF5 and/or MAT files and/or folder containing CSV files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) or folder containing CSV files |
map_info: (optional) map information, only used if not inmap_filemap_field: (optional) struct field within MAT file to use, not relevant for CSV or HDF5 filemap_units: (optional) map xx/yy units to use inmap_map{:rad,:deg}file_units: (optional) map xx/yy units used in files withindf_map{:rad,:deg}flip_map: (optional) if true, vertically flip data from map file (possibly useful for CSV map)
Returns:
map_map:Mapmagnetic anomaly map struct
MagNav.get_map — Functionget_map(map_file::String = namad,
map_field::Symbol = :map_data;
map_info::String = splitpath(map_file)[end],
map_units::Symbol = :rad,
file_units::Symbol = :deg,
flip_map::Bool = false)Get map data from saved HDF5 or MAT file or folder containing CSV files. Maps are typically saved in :deg units, while :rad is used internally.
Arguments:
map_file: path/name of map data HDF5 or MAT file (.h5or.matextension required) or folder containing CSV filesmap_field: (optional) struct field within MAT file to use, not relevant for CSV or HDF5 filemap_info: (optional) map information, only used if not inmap_filemap_units: (optional) map xx/yy units to use inmap_map{:rad,:deg}file_units: (optional) map xx/yy units used inmap_file{:rad,:deg}flip_map: (optional) if true, vertically flip data from map file (possibly useful for CSV map)
Returns:
map_map:Mapmagnetic anomaly map struct
MagNav.get_map_val — Functionget_map_val(map_map::Map, path::Path, ind = trues(path.N);
α=200, return_itp::Bool = false)Get scalar magnetic anomaly map values along a flight path. map_map is upward and/or downward continued to alt as necessary.
Arguments:
map_map:Mapmagnetic anomaly map structpath:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indicesα: (optional) regularization parameter for downward continuationreturn_itp: (optional) if true, also returnmap_itp
Returns:
map_val: scalar magnetic anomaly map valuesmap_itp: ifreturn_itp = true, map interpolation function (f(lat,lon)orf(lat,lon,alt))
MagNav.get_map_val — Functionget_map_val(map_map_vec::Vector, path::Path, ind = trues(path.N); α = 200)Get scalar magnetic anomaly map values from multiple maps along a flight path. Each map in map_map_vec is upward and/or downward continued to alt as necessary.
Arguments:
map_map_vec: vector ofMapmagnetic anomaly map structspath:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indicesα: (optional) regularization parameter for downward continuation
Returns:
map_vals: vector of scalar magnetic anomaly map values
MagNav.get_map_val — Methodget_map_val(map_map::Map, lat, lon, alt; α = 200, return_itp::Bool = false)Get scalar magnetic anomaly map values along a flight path. map_map is upward and/or downward continued to alt as necessary (except if drape map).
Arguments:
map_map:Mapmagnetic anomaly map structlat: latitude [rad]lon: longitude [rad]alt: altitude [m]α: (optional) regularization parameter for downward continuationreturn_itp: (optional) if true, also returnitp_map
Returns:
map_val: scalar magnetic anomaly map valuesitp_map: ifreturn_itp = true, map interpolation function (f(lat,lon)orf(lat,lon,alt))
MagNav.get_nn_m — Functionget_nn_m(Nf::Int, Ny::Int = 1;
hidden = [8],
activation::Function = swish,
final_bias::Bool = true,
skip_con::Bool = false,
model_type::Symbol = :m1
l_window::Int = 5,
tf_layer_type::Symbol = :postlayer,
tf_norm_type::Symbol = :batch,
dropout_prob::Real = 0.2,
N_tf_head::Int = 8,
tf_gain::Real = 1.0)Get neural network model. Valid for 0-3 hidden layers, except for model_type = :m3w, :m3tf, which are only valid for 1 or 2 hidden layers.
Arguments:
Nf: length of input (feature) layerNy: (optional) length of output layerhidden: (optional) hidden layers & nodes (e.g.,[8,8]for 2 hidden layers, 8 nodes each)activation: (optional) activation functionrelu= rectified linear unitσ= sigmoid (logistic function)swish= self-gatedtanh= hyperbolic tan- run
plot_activation()for a visual
final_bias: (optional) if true, include final layer biasskip_con: (optional) if true, use skip connections, must have length(hidden) == 1model_type: (optional) aeromagnetic compensation model typel_window: (optional) temporal window length, only used formodel_type = :m3w, :m3tftf_layer_type: (optional) transformer normalization layer before or after skip connection {:prelayer,:postlayer}, only used formodel_type = :m3tftf_norm_type: (optional) normalization for transformer encoder {:batch,:layer,:none}, only used formodel_type = :m3tfdropout_prob: (optional) dropout rate, only used formodel_type = :m3w, :m3tfN_tf_head: (optional) number of attention heads, only used formodel_type = :m3tftf_gain: (optional) weight initialization parameter, only used formodel_type = :m3tf
Returns:
m: neural network model
MagNav.get_optimal_rotation_matrix — Methodget_optimal_rotation_matrix(v1s, v2s)Get the 3 x 3 rotation matrix rotating the directions of v1s into v2s. Uses the Kabsch algorithm.
Reference: https://en.wikipedia.org/wiki/Kabsch_algorithm
Arguments:
v1s:Nx3matrix for first set of 3D pointsv2s:Nx3matrix for second set of 3D points
Returns:
R: 3D matrix rotatating v1s into v2s' directions
MagNav.get_pinson — Methodget_pinson(nx::Int, lat, vn, ve, vd, fn, fe, fd, Cnb;
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
vec_states::Bool = false,
fogm_state::Bool = true,
k1=3e-2, k2=3e-4, k3=1e-6)Get the nx x nx Pinson dynamics matrix. States (errors) are:
| Num | State | Units | Description |
|---|---|---|---|
| 1 | lat | rad | latitude |
| 2 | lon | rad | longitude |
| 3 | alt | m | altitude |
| 4 | vn | m/s | north velocity |
| 5 | ve | m/s | east velocity |
| 6 | vd | m/s | down velocity |
| 7 | tn | rad | north tilt (attitude) |
| 8 | te | rad | east tilt (attitude) |
| 9 | td | rad | down tilt (attitude) |
| 10 | ha | m | barometer aiding altitude |
| 11 | a_hat | m/s^2 | barometer aiding vertical acceleration |
| 12 | ax | m/s^2 | x accelerometer |
| 13 | ay | m/s^2 | y accelerometer |
| 14 | az | m/s^2 | z accelerometer |
| 15 | gx | rad/s | x gyroscope |
| 16 | gy | rad/s | y gyroscope |
| 17 | gz | rad/s | z gyroscope |
| end | S | nT | FOGM catch-all |
Arguments:
nx: total state dimensionlat: latitude [rad]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]baro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]vec_states: (optional) if true, include vector magnetometer statesfogm_state: (optional) if true, include FOGM catch-all bias statek1: (optional) barometer aiding constant [1/s]k2: (optional) barometer aiding constant [1/s^2]k3: (optional) barometer aiding constant [1/s^3]
Returns:
F:nxxnxPinson matrix
MagNav.get_step — Methodget_step(x::AbstractVector)Get the step size (spacing) of elements in x.
MagNav.get_traj — Functionget_traj(xyz::XYZ, ind = trues(xyz.traj.N))Get trajectory data at specific indices.
Arguments:
xyz:XYZflight data structind: (optional) selected data indices
Returns:
traj:Trajtrajectory struct atind
MagNav.get_traj — Functionget_traj(traj_file::String, field::Symbol = :traj;
dt = 0.1,
tt_sort::Bool = true,
silent::Bool = false)Get trajectory data from saved CSV, HDF5, or MAT file. The only required fields are lat, lon, and alt (position).
If a CSV or HDF5 file is provided, the possible columns/fields in the file are:
| Field | Type | Description |
|---|---|---|
dt | scalar | measurement time step [s] |
tt | vector | time [s] |
lat | vector | latitude [deg] |
lon | vector | longitude [deg] |
alt | vector | altitude [m] |
vn | vector | north velocity [m/s] |
ve | vector | east velocity [m/s] |
vd | vector | down velocity [m/s] |
fn | vector | north specific force [m/s] |
fe | vector | east specific force [m/s] |
fd | vector | down specific force [m/s] |
Cnb | 3x3xN | direction cosine matrix (body to navigation) [-], not valid for CSV file (use roll, pitch, & yaw instead) |
roll | vector | roll [deg] |
pitch | vector | pitch [deg] |
yaw | vector | yaw [deg] |
If a MAT file is provided, the above fields may also be provided, but they should be within the specified field MAT struct. This is the standard way the MATLAB-companion outputs data.
Arguments:
traj_file: path/name of trajectory data CSV, HDF5, or MAT file (.csv,.h5, or.matextension required)field: (optional) struct field within MAT file to use, not relevant for CSV or HDF5 filedt: (optional) measurement time step [s], only used if not intraj_filesilent: (optional) if true, no print outs
Returns:
traj:Trajtrajectory struct
MagNav.get_x — Functionget_x(xyz::XYZ, ind = trues(xyz.traj.N),
features_setup::Vector{Symbol} = [:mag_1_uc,:TL_A_flux_a];
features_no_norm::Vector{Symbol} = Symbol[],
terms = [:permanent,:induced,:eddy],
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
bpf_mag::Bool = false)Get x data matrix.
Arguments:
xyz:XYZflight data structind: selected data indicesfeatures_setup: vector of features to includefeatures_no_norm: (optional) vector of features to not normalizeterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}sub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsbpf_mag: (optional) if true, bpf scalar magnetometer measurements
Returns:
x:NxNfdata matrix (Nfis number of features)no_norm: length-NfBoolean indices of features to not be normalizedfeatures: length-Nffeature vector (including components of TLA, etc.)l_segs: length-N_linesvector of lengths oflines, sum(l_segs) =N
MagNav.get_x — Functionget_x(lines, df_line::DataFrame, df_flight::DataFrame,
features_setup::Vector{Symbol} = [:mag_1_uc,:TL_A_flux_a];
features_no_norm::Vector{Symbol} = Symbol[],
terms = [:permanent,:induced,:eddy],
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
bpf_mag::Bool = false,
reorient_vec::Bool = false,
l_window::Int = -1,
silent::Bool = true)Get x data matrix from multiple flight lines, possibly multiple flights.
Arguments:
lines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
features_setup: vector of features to includefeatures_no_norm: (optional) vector of features to not normalizeterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}sub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsbpf_mag: (optional) if true, bpf scalar magnetometer measurementsreorient_vec: (optional) if true, align vector magnetometer measurements with body framel_window: (optional) trim data byN % l_window,-1to ignoresilent: (optional) if true, no print outs
Returns:
x:NxNfdata matrix (Nfis number of features)no_norm: length-NfBoolean indices of features to not be normalizedfeatures: length-Nffeature vector (including components of TLA, etc.)l_segs: length-N_linesvector of lengths oflines, sum(l_segs) =N
MagNav.get_x — Functionget_x(xyz_vec::Vector{XYZ}, ind_vec,
features_setup::Vector{Symbol} = [:mag_1_uc,:TL_A_flux_a];
features_no_norm::Vector{Symbol} = Symbol[],
terms = [:permanent,:induced,:eddy],
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
bpf_mag::Bool = false)Get x data matrix from multiple XYZ flight data structs.
Arguments:
xyz_vec: vector ofXYZflight data structsind_vec: vector of selected data indicesfeatures_setup: vector of features to includefeatures_no_norm: (optional) vector of features to not normalizeterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}sub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsbpf_mag: (optional) if true, bpf scalar magnetometer measurements
Returns:
x:NxNfdata matrix (Nfis number of features)no_norm: length-NfBoolean indices of features to not be normalizedfeatures: length-Nffeature vector (including components of TLA, etc.)l_segs: length-N_linesvector of lengths oflines, sum(l_segs) =N
MagNav.get_y — Functionget_y(xyz::XYZ, ind = trues(xyz.traj.N),
map_val = -1);
y_type::Symbol = :d,
use_mag::Symbol = :mag_1_uc,
use_mag_c::Symbol = :mag_1_c,
sub_diurnal::Bool = false,
sub_igrf::Bool = false)Get y target vector.
Arguments:
xyz:XYZflight data structind: selected data indicesmap_val: (optional) scalar magnetic anomaly map values, only used fory_type = :by_type: (optional)ytarget type:a= anomaly field #1, compensated tail stinger total field scalar magnetometer measurements:b= anomaly field #2, interpolatedmagnetic anomaly mapvalues:c= aircraft field #1, difference between uncompensated cabin total field scalar magnetometer measurements and interpolatedmagnetic anomaly mapvalues:d= aircraft field #2, difference between uncompensated cabin and compensated tail stinger total field scalar magnetometer measurements:e= BPF'd total field, bandpass filtered uncompensated cabin total field scalar magnetometer measurements
use_mag: (optional) uncompensated scalar magnetometer to use forytarget vector {:mag_1_uc, etc.}, only used fory_type = :c, :d, :euse_mag_c: (optional) compensated scalar magnetometer to use forytarget vector {:mag_1_c, etc.}, only used fory_type = :a, :dsub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurements
Returns:
y: length-Ntarget vector
MagNav.get_y — Methodget_y(lines, df_line::DataFrame, df_flight::DataFrame,
df_map::DataFrame;
y_type::Symbol = :d,
use_mag::Symbol = :mag_1_uc,
use_mag_c::Symbol = :mag_1_c,
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
l_window::Int = -1,
silent::Bool = true)Get y target vector from multiple flight lines, possibly multiple flights.
Arguments:
lines: selected line number(s)df_line: lookup table (DataFrame) oflines
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
line | Real | line number, i.e., segments within flight |
t_start | Real | start time of line to use [s] |
t_end | Real | end time of line to use [s] |
map_name | Symbol | (optional) name of magnetic anomaly map relevant to line, only used for y_type = :b, :c |
df_flight: lookup table (DataFrame) of flight data files
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
xyz_type | Symbol | subtype of XYZ to use for flight data {:XYZ0,:XYZ1,:XYZ20,:XYZ21} |
xyz_set | Real | flight dataset number (used to prevent improper mixing of datasets, such as different magnetometer locations) |
xyz_file | String | path/name of flight data CSV, HDF5, or MAT file (.csv, .h5, or .mat extension required) |
df_map: lookup table (DataFrame) of map data files
| Field | Type | Description |
|---|---|---|
map_name | Symbol | name of magnetic anomaly map |
map_file | String | path/name of map data HDF5 or MAT file (.h5 or .mat extension required) |
y_type: (optional)ytarget type:a= anomaly field #1, compensated tail stinger total field scalar magnetometer measurements:b= anomaly field #2, interpolatedmagnetic anomaly mapvalues:c= aircraft field #1, difference between uncompensated cabin total field scalar magnetometer measurements and interpolatedmagnetic anomaly mapvalues:d= aircraft field #2, difference between uncompensated cabin and compensated tail stinger total field scalar magnetometer measurements:e= BPF'd total field, bandpass filtered uncompensated cabin total field scalar magnetometer measurements
use_mag: (optional) uncompensated scalar magnetometer to use forytarget vector {:mag_1_uc, etc.}, only used fory_type = :c, :d, :euse_mag_c: (optional) compensated scalar magnetometer to use forytarget vector {:mag_1_c, etc.}, only used fory_type = :a, :dsub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsl_window: (optional) trim data byN % l_window,-1to ignoresilent: (optional) if true, no print outs
Returns:
y: length-Ntarget vector
MagNav.get_years — Methodget_years(year, doy=0)Get decimal (fractional) year from year and doy (day of year).
Arguments:
year: yeardoy: day of year
Returns:
years: decimal (fractional) year
MagNav.gif_animation_m3 — Functiongif_animation_m3(TL_perm::AbstractMatrix, TL_induced::AbstractMatrix, TL_eddy::AbstractMatrix,
TL_aircraft::AbstractMatrix, B_unit::AbstractMatrix, y_nn::AbstractMatrix,
y::Vector, y_hat::Vector, xyz::XYZ,
filt_lat::Vector = [],
filt_lon::Vector = [];
ind = trues(xyz.traj.N),
tt_lim::Tuple = (0, (xyz.traj(ind).N-1)*xyz.traj.dt/60),
skip_every::Int = 5,
save_plot::Bool = false,
mag_gif::String = "comp_xai.gif")Create a GIF animation of the model 3 components and the true and predicted scalar magnetic field. First run comp_m3_test() to generate the individual model components 3.
Arguments
TL_perm:3xNmatrix of TL permanent vector fieldTL_induced:3xNmatrix of TL induced vector fieldTL_eddy:3xNmatrix of TL eddy current vector fieldTL_aircraft:3xNmatrix of TL aircraft vector fieldB_unit:3xNmatrix of normalized vector magnetometer measurementsy_nn:3xNmatrix of vector neural network correction (for scalar models, in direction ofBt)y: length-Ntarget vectory_hat: length-Nprediction vectorxyz:XYZflight data structfilt_lat: (optional) length-Nfilter output latitude [rad]filt_lon: (optional) length-Nfilter output longitude [rad]ind: (optional) selected data indicestt_lim: (optional) length-2start & end time limits (inclusive) [min]skip_every: (optional) number of time steps to skip between framessave_plot: (optional) if true, saveg1asmag_gifmag_gif: (optional) path/name of magnetic field GIF file to save (.gifextension optional)
Returns
g1: magnetic field GIF animation
Example
gif_animation_m3(TL_perm, TL_induced, TL_eddy, TL_aircraft, B_unit,
y_nn, y, y_hat, xyz, filt_lat, filt_lon; ind=ind, tt_lim=(0.0,10.0),
skip_every=5, save_plot=false, mag_gif="comp_xai.gif")MagNav.gif_ellipse — Functiongif_ellipse(filt_res::FILTres,
filt_out::FILTout,
map_map::Map = mapS_null;
dt = 0.1,
di::Int = 10,
speedup::Int = 60,
conf_units::Symbol = :m,
μ = zeros(eltype(filt_res.P),2),
conf = 0.95,
clip = Inf,
n::Int = 61,
lim = 500,
dpi::Int = 200,
margin::Int = 2,
axis::Bool = true,
plot_eigax::Bool = false,
bg_color::Symbol = :white,
ce_color::Symbol = :black,
map_color::Symbol = :usgs,
clims::Tuple = (),
b_e::AbstractBackend = gr(),
save_plot::Bool = false,
ellipse_gif::String = "conf_ellipse.gif")Create a (position) confidence ellipse GIF animation for a 2 x 2 (x N) covariance matrix.
Arguments:
filt_res:FILTresfilter results structfilt_out:FILToutfilter extracted output structmap_map: (optional)Mapmagnetic anomaly map structdt: (optional) measurement time step [s]di: (optional) GIF measurement interval (e.g.,di = 10uses every 10th measurement)speedup: (optional) GIF speedup (e.g.,speedup = 60is 60x speed)conf_units: (optional) confidence ellipse units {:m,:ft,:deg,:rad}μ: (optional) confidence ellipse center [conf_units]conf: (optional) percentile {0:1}clip: (optional) clipping radius [conf_units]n: (optional) number of confidence ellipse pointslim: (optional) plotx&ylimits (-lim,lim) [conf_units]dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]axis: (optional) if true, show axesplot_eigax: (optional) if true, show major & minor axesbg_color: (optional) background colorce_color: (optional) confidence ellipse colormap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}clims: (optional) length-2map colorbar limits(cmin,cmax)b_e: (optional) plotting backendsave_plot: (optional) if true, saveg1asellipse_gifellipse_gif: (optional) path/name of confidence ellipse GIF file to save (.gifextension optional)
Returns:
g1: confidence ellipse GIF animation
MagNav.gif_ellipse — Functiongif_ellipse(P, lat1 = deg2rad(45);
dt = 0.1,
di::Int = 10,
speedup::Int = 60,
conf_units::Symbol = :m,
μ = zeros(eltype(P),2),
conf = 0.95,
clip = Inf,
n::Int = 61,
lim = 500,
margin::Int = 2,
axis::Bool = true,
plot_eigax::Bool = false,
bg_color::Symbol = :white,
ce_color::Symbol = :black,
b_e::AbstractBackend = gr(),
save_plot::Bool = false,
ellipse_gif::String = "conf_ellipse.gif")Create a (position) confidence ellipse GIF animation for a 2 x 2 (x N) covariance matrix.
Arguments:
P:2x2(xN) covariance matrixlat1: (optional) nominal latitude [rad], only used ifconf_units = :mor:ftdt: (optional) measurement time step [s]di: (optional) GIF measurement interval (e.g.,di = 10uses every 10th measurement)speedup: (optional) GIF speedup (e.g.,speedup = 60is 60x speed)conf_units: (optional) confidence ellipse units {:m,:ft,:deg,:rad}μ: (optional) confidence ellipse center [conf_units]conf: (optional) percentile {0:1}clip: (optional) clipping radius [conf_units]n: (optional) number of confidence ellipse pointslim: (optional) plotx&ylimits (-lim,lim) [conf_units]margin: (optional) margin around plot [mm]axis: (optional) if true, show axesplot_eigax: (optional) if true, show major & minor axesbg_color: (optional) background colorce_color: (optional) confidence ellipse colorb_e: (optional) plotting backendsave_plot: (optional) if true, saveg1asellipse_gifellipse_gif: (optional) path/name of confidence ellipse GIF file to save (.gifextension optional)
Returns:
g1: confidence ellipse GIF animation
MagNav.krr_fit — Functionkrr_fit(x, y, no_norm = falses(size(x,2));
k::Kernel = PolynomialKernel(;degree=1),
λ::Real = 0.5,
norm_type_x::Symbol = :standardize,
norm_type_y::Symbol = :standardize,
data_norms::Tuple = (zeros(1,1),zeros(1,1),[0.0],[0.0]),
l_segs::Vector = [length(y)],
silent::Bool = false)Fit a kernel ridge regression (KRR) model to data.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorno_norm: (optional) length-NfBoolean indices of features to not be normalizedk: (optional) kernelλ: (optional) ridge parameternorm_type_x: (optional) normalization forxdata matrixnorm_type_y: (optional) normalization forytarget vectordata_norms: (optional) length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
model: length-3tuple of KRR-based model, (k, length-Ncoefficients,NxNfdata matrix, normalized)data_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.krr_test — Methodkrr_test(x, y, data_norms::Tuple, model::Tuple;
l_segs::Vector = [length(y)],
silent::Bool = false)Evaluate performance of a kernel ridge regression (KRR) model.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectordata_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)model: length-3tuple of KRR-based model, (k, length-N_traincoefficients,N_trainxNftraining data matrix, normalized)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.linear_fit — Functionlinear_fit(x, y, no_norm = falses(size(x,2));
trim::Int = 0,
λ::Real = 0,
norm_type_x::Symbol = :none,
norm_type_y::Symbol = :none,
data_norms::Tuple = (zeros(1,1),zeros(1,1),[0.0],[0.0]),
l_segs::Vector = [length(y)],
silent::Bool = false)Fit a linear regression model to data.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectorno_norm: (optional) length-NfBoolean indices of features to not be normalizedtrim: (optional) number of elements to trim (e.g., due to bpf)λ: (optional) ridge parameternorm_type_x: (optional) normalization forxdata matrixnorm_type_y: (optional) normalization forytarget vectordata_norms: (optional) length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
model: length-2tuple of linear regression model, (length-Nfcoefficients, bias=0)data_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.linear_test — Methodlinear_test(x_norm, y, y_bias, y_scale, model::Tuple;
l_segs::Vector = [length(y)],
silent::Bool = false)Evaluate performance of a linear model.
Arguments:
x_norm:NxNfnormalized data matrix (Nfis number of features)y: length-Ntarget vectory_bias:ytarget vector bias bias (mean, min, or zero)y_scale:ytarget vector bias scaling factor (std dev, max-min, or one)model: length-2tuple of model, (length-Nfcoefficients, bias)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.linear_test — Methodlinear_test(x, y, data_norms::Tuple, model::Tuple;
l_segs::Vector = [length(y)],
silent::Bool = false)Evaluate performance of a linear model.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectordata_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)model: length-2tuple of model, (length-Nfcoefficients, bias)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nsilent: (optional) if true, no print outs
Returns:
y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) error
MagNav.linreg — Methodlinreg(y, x; λ=0)Linear regression with data matrix.
Arguments:
y: length-Nobserved data vectorx:NxNfinput data matrix (Nfis number of features)λ: (optional) ridge parameter
Returns:
coef: linear regression coefficients
MagNav.linreg — Methodlinreg(y; λ=0)Linear regression to determine best fit line for x = eachindex(y).
Arguments:
y: length-Nobserved data vectorλ: (optional) ridge parameter
Returns:
coef: length-2vector of linear regression coefficients
MagNav.map2kmz — Functionmap2kmz(map_map::Matrix, map_xx::Vector, map_yy::Vector,
map_kmz::String = "map.kmz";
map_units::Symbol = :rad,
plot_alt::Real = 0,
opacity::Real = 0.75,
clims::Tuple = ())Create KMZ file of map for use with Google Earth. Generates an "icon" overlay, thus not suitable for large maps (e.g., > 5 deg x 5 deg).
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinates [rad] or [deg]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg]map_kmz: (optional) path/name of map KMZ file to save (.kmzextension optional)map_units: (optional) map xx/yy units {:rad,:deg}plot_alt: (optional) map altitude in Google Earth [m]opacity: (optional) map opacity {0:1}clims: (optional) length-2map colorbar limits(cmin,cmax)
Returns:
nothing:map_kmzis created
MagNav.map2kmz — Functionmap2kmz(mapS::Union{MapS,MapSd,MapS3D},
map_kmz::String = "map.kmz";
use_mask::Bool = true,
plot_alt::Real = 0,
opacity::Real = 0.75,
clims::Tuple = ())Create KMZ file of map for use with Google Earth. Generates an "icon" overlay, thus not suitable for large maps (e.g., > 5 deg x 5 deg).
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structmap_kmz: (optional) path/name of map KMZ file to save (.kmzextension optional)use_mask: (optional) if true, applymapSmask to mapplot_alt: (optional) map altitude in Google Earth [m]opacity: (optional) map opacity {0:1}clims: (optional) length-2map colorbar limits(cmin,cmax)
Returns:
nothing:map_kmzis created
MagNav.map_border — Methodmap_border(map_map::Matrix, map_xx::Vector, map_yy::Vector;
inner::Bool = true,
sort_border::Bool = false,
return_ind::Bool = false)Get map border from an unfilled map.
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinatesinner: (optional) if true, get inner border, otherwise outer bordersort_border: (optional) if true, sort border data points sequentiallyreturn_ind: (optional) if true, returnind
Returns:
yy: border y-direction (latitude) coordinatesxx: border x-direction (longitude) coordinatesind: ifreturn_ind = true,BitMatrixof border indices withinmap_map
MagNav.map_border — Methodmap_border(mapS::Union{MapS,MapSd,MapS3D};
inner::Bool = true,
sort_border::Bool = false,
return_ind::Bool = false)Get map border from an unfilled map.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structinner: (optional) if true, get inner border, otherwise outer bordersort_border: (optional) if true, sort border data points sequentiallyreturn_ind: (optional) if true, returnind
Returns:
yy: border y-direction (latitude) coordinatesxx: border x-direction (longitude) coordinatesind: ifreturn_ind = true,BitMatrixof border indices withinmap_map
MagNav.map_check — Functionmap_check(map_map::Map, path::Path, ind = trues(path.N))Check if latitude and longitude points are on given map.
Arguments:
map_map:Mapmagnetic anomaly map structpath:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indices
Returns:
bool: if true, allpath[ind] points are onmap_map
MagNav.map_check — Functionmap_check(map_map_vec::Vector, path::Path, ind = trues(path.N))Check if latitude and longitude points are on given maps.
Arguments:
map_map_vec: vector ofMapmagnetic anomaly map structspath:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indices
Returns:
bools: if true, allpath[ind] points are onmap_map_vec[i]
MagNav.map_check — Functionmap_check(map_map::Map, lat, lon, alt = fill(median(map_map.alt),size(lat)))Check if latitude and longitude points are on given map.
Arguments:
map_map:Mapmagnetic anomaly map structlat: latitude [rad]lon: longitude [rad]alt: (optional) altitude [m], only used forMapS3D
Returns:
bool: if true, alllat&lon(&alt) points are onmap_map
MagNav.map_chessboard! — Methodmap_chessboard!(map_map::Matrix, map_alt::Matrix, map_xx::Vector,
map_yy::Vector, alt::Real;
down_cont::Bool = true,
dz = 5,
down_max = 150,
α = 200)Perform the chessboard method, which upward (and possibly downward) continues a map to multiple altitudes to create a 3D map, then vertically interpolates at each horizontal grid point.
Reference: Cordell, Phillips, & Godson, U.S. Geological Survey Potential-Field Software Version 2.0, 1992.
Arguments:
map_map:nyxnx2D gridded target (e.g., magnetic) map data on [m] gridmap_alt:nyxnx2D gridded altitude map data [m] on [m] gridmap_xx:nxmap x-direction (longitude) coordinates [m]map_yy:nymap y-direction (latitude) coordinates [m]alt: final map altitude after upward continuation [m]down_cont: (optional) if true, downward continue if needed, only used ifup_cont = truedz: (optional) upward continuation step size [m]down_max: (optional) maximum downward continuation distance [m]α: (optional) regularization parameter for downward continuation
Returns:
nothing:map_mapis mutated with upward continued map data
MagNav.map_chessboard — Methodmap_chessboard(mapSd::MapSd, alt::Real;
down_cont::Bool = true,
dz = 5,
down_max = 150,
α = 200)Perform the chessboard method, which upward (and possibly downward) continues a map to multiple altitudes to create a 3D map, then vertically interpolates at each horizontal grid point.
Reference: Cordell, Phillips, & Godson, U.S. Geological Survey Potential-Field Software Version 2.0, 1992.
Arguments:
mapSd:MapSdscalar magnetic anomaly map structalt: final map altitude after upward continuation [m]down_cont: (optional) if true, downward continue if neededdz: (optional) upward continuation step size [m]down_max: (optional) maximum downward continuation distance [m]α: (optional) regularization parameter for downward continuation
Returns:
mapS:MapSscalar magnetic anomaly map struct
MagNav.map_combine — Functionmap_combine(mapS::MapS, mapS_fallback::MapS = get_map(namad);
map_info::String = mapS.info,
xx_lim::Tuple = get_lim(mapS.xx,0.1),
yy_lim::Tuple = get_lim(mapS.yy,0.1),
α = 200)Combine two maps at same altitude.
Arguments:
mapS:MapSscalar magnetic anomaly map structmapS_fallback: (optional) fallbackMapSscalar magnetic anomaly map structmap_info: (optional) map informationxx_lim: (optional) length-2x-direction map limits(xx_min,xx_max)yy_lim: (optional) length-2y-direction map limits(yy_min,yy_max)α: (optional) regularization parameter for downward continuation
Returns:
mapS:MapSscalar magnetic anomaly map struct, combined
MagNav.map_combine — Functionmap_combine(mapS_vec::Vector, mapS_fallback::MapS = get_map(namad);
map_info::String = "Combined map",
N_levels::Int = 3,
dx = get_step(mapS_vec[1].xx),
dy = get_step(mapS_vec[1].yy),
xx_lim::Tuple = get_lim(mapS_vec[1].xx,0.5),
yy_lim::Tuple = get_lim(mapS_vec[1].yy,0.5),
α = 200,
use_fallback::Bool = true)Combine maps at different altitudes. Lowest and highest maps are directly used (with resampling & resizing), with intermediate maps determined by N_levels.
Arguments:
mapS_vec: vector ofMapSscalar magnetic anomaly map structsmapS_fallback: (optional) fallbackMapSscalar magnetic anomaly map structmap_info: (optional) map informationN_levels: (optional) number of map altitude levelsdx: (optional) desired x-direction map step sizedy: (optional) desired y-direction map step sizexx_lim: (optional) length-2x-direction map limits(xx_min,xx_max)yy_lim: (optional) length-2y-direction map limits(yy_min,yy_max)α: (optional) regularization parameter for downward continuationuse_fallback: (optional) if true, usemapS_fallbackfor missing map data
Returns:
mapS3D:MapS3D3D (multi-level) scalar magnetic anomaly map struct
MagNav.map_correct_igrf! — Methodmap_correct_igrf!(map_map::Matrix, map_alt,
map_xx::Vector, map_yy::Vector;
sub_igrf_date::Real = get_years(2013,293), # 20-Oct-2013
add_igrf_date::Real = -1,
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad)Correct the International Geomagnetic Reference Field (IGRF), i.e., core field, of a map by subtracting and/or adding the IGRF on specified date(s).
Arguments:
map_map:nyxnx2D gridded map datamap_alt:nyxnx2D gridded altitude map data, single altitude value may be provided [m]map_xx:nxmap x-direction (longitude) coordinates [rad] or [deg] or [m]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg] or [m]sub_igrf_date: (optional) date of IGRF core field to subtract [yr], -1 to ignoreadd_igrf_date: (optional) date of IGRF core field to add [yr], -1 to ignorezone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}
Returns:
nothing:map_mapis mutated with IGRF corrected map data
MagNav.map_correct_igrf! — Methodmap_correct_igrf!(mapS::Union{MapS,MapSd,MapS3D};
sub_igrf_date::Real = get_years(2013,293), # 20-Oct-2013
add_igrf_date::Real = -1,
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad)Correct the International Geomagnetic Reference Field (IGRF), i.e., core field, of a map by subtracting and/or adding the IGRF on specified date(s).
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structsub_igrf_date: (optional) date of IGRF core field to subtract [yr], -1 to ignoreadd_igrf_date: (optional) date of IGRF core field to add [yr], -1 to ignorezone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}
Returns:
nothing:mapfield withinmapSis mutated with IGRF corrected map data
MagNav.map_correct_igrf — Methodmap_correct_igrf(map_map::Matrix, map_alt,
map_xx::Vector, map_yy::Vector;
sub_igrf_date::Real = get_years(2013,293), # 20-Oct-2013
add_igrf_date::Real = -1,
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad)Correct the International Geomagnetic Reference Field (IGRF), i.e., core field, of a map by subtracting and/or adding the IGRF on specified date(s).
Arguments:
map_map:nyxnx2D gridded map datamap_alt:nyxnx2D gridded altitude map data [m]map_xx:nxmap x-direction (longitude) coordinates [rad] or [deg] or [m]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg] or [m]sub_igrf_date: (optional) date of IGRF core field to subtract [yr], -1 to ignoreadd_igrf_date: (optional) date of IGRF core field to add [yr], -1 to ignorezone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}
Returns:
map_map:nyxnx2D gridded map data, IGRF corrected
MagNav.map_correct_igrf — Methodmap_correct_igrf(mapS::Union{MapS,MapSd,MapS3D};
sub_igrf_date::Real = get_years(2013,293), # 20-Oct-2013
add_igrf_date::Real = -1,
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad)Correct the International Geomagnetic Reference Field (IGRF), i.e., core field, of a map by subtracting and/or adding the IGRF on specified date(s).
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structsub_igrf_date: (optional) date of IGRF core field to subtract [yr], -1 to ignoreadd_igrf_date: (optional) date of IGRF core field to add [yr], -1 to ignorezone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}
Returns:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct, IGRF corrected
MagNav.map_expand — Functionmap_expand(map_map::Matrix, pad::Int = 1)Expand a map with padding on each edge to eliminate discontinuities in the discrete Fourier transform. The map is “wrapped around” to make it periodic. Padding expands the map to 7-smooth dimensions, allowing for a faster Fast Fourier Transform algorithm to be used during upward/downward continuation.
Arguments:
map_map:nyxnx2D gridded map datapad: minimum padding (grid cells) along map edges
Returns:
map_map:nyxnx2D gridded map data, expanded (padded)padx: x-direction padding (grid cells) applied on first edgepady: y-direction padding (grid cells) applied on first edge
MagNav.map_fill! — Methodmap_fill!(map_map::Matrix, map_xx::Vector, map_yy::Vector; k::Int = 3)Fill areas that are missing map data.
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinatesk: (optional) number of nearest neighbors for knn
Returns:
nothing:map_mapis mutated with filled map data
MagNav.map_fill! — Methodmap_fill!(mapS::Union{MapS,MapSd,MapS3D}; k::Int = 3)Fill areas that are missing map data.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structk: (optional) number of nearest neighbors for knn
Returns:
nothing:mapfield withinmapSis mutated with filled map data
MagNav.map_fill — Methodmap_fill(map_map::Matrix, map_xx::Vector, map_yy::Vector; k::Int = 3)Fill areas that are missing map data.
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinatesk: (optional) number of nearest neighbors for knn
Returns:
map_map:nyxnx2D gridded map data, filled
MagNav.map_fill — Methodmap_fill(mapS::Union{MapS,MapSd,MapS3D}; k::Int = 3)Fill areas that are missing map data.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structk: (optional) number of nearest neighbors for knn
Returns:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct, filled
MagNav.map_get_gxf — Methodmap_get_gxf(map_gxf::String)Use ArchGDAL to read in map data from GXF file.
Arguments:
map_gxf: path/name of map GXF file (.gxfextension optional)
Returns:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinates
MagNav.map_gxf2h5 — Methodmap_gxf2h5(map_gxf::String, alt::Real;
map_info::String = splitpath(map_gxf)[end],
fill_map::Bool = true,
get_lla::Bool = true,
zone_utm::Int = 18,
is_north::Bool = true,
save_h5::Bool = false,
map_h5::String = "map_data.h5")Convert map data file (with assumed UTM grid) from GXF to HDF5. The order of operations is:
- original map from
map_gxf=> - trim away large areas that are missing map data =>
- fill remaining areas that are missing map data =>
- convert map grid from
UTMtoLLA
Specifically meant for SMALL and LEVEL maps ONLY.
Arguments:
map_gxf: path/name of target (e.g., magnetic) map GXF file (.gxfextension optional)alt: map altitude [m]map_info: (optional) map informationfill_map: (optional) if true, fill areas that are missing map dataget_lla: (optional) if true, convert map grid fromUTMtoLLAzone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispheresave_h5: (optional) if true, savemapStomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
mapS:MapSscalar magnetic anomaly map struct
MagNav.map_gxf2h5 — Methodmap_gxf2h5(map_gxf::String, alt_gxf::String, alt::Real;
map_info::String = splitpath(map_gxf)[end],
pad::Int = 0,
sub_igrf_date::Real = get_years(2013,293),
add_igrf_date::Real = -1,
zone_utm::Int = 18,
is_north::Bool = true,
fill_map::Bool = true,
up_cont::Bool = true,
down_cont::Bool = true,
get_lla::Bool = true,
dz::Real = 5,
down_max::Real = 150,
α::Real = 200,
save_h5::Bool = false,
map_h5::String = "map_data.h5")Convert map data file (with assumed UTM grid) from GXF to HDF5. The order of operations is:
- original map from
map_gxf=> - trim away large areas that are missing map data =>
- subtract and/or add IGRF to map data =>
- fill remaining areas that are missing map data =>
- upward/downward continue to
alt=> - convert map grid from
UTMtoLLA
This can be memory intensive, largely depending on the map size and dz. If up_cont = true, a MapS struct is returned. If up_cont = false, a MapSd struct is returned, which has an included altitude map.
Arguments:
map_gxf: path/name of target (e.g., magnetic) map GXF file (.gxfextension optional)alt_gxf: path/name of altitude map GXF file (.gxfextension optional)alt: final map altitude after upward continuation [m], not used for drape mapmap_info: (optional) map informationpad: (optional) minimum padding (grid cells) along map edgessub_igrf_date: (optional) date of IGRF core field to subtract [yr], -1 to ignoreadd_igrf_date: (optional) date of IGRF core field to add [yr], -1 to ignorezone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispherefill_map: (optional) if true, fill areas that are missing map dataup_cont: (optional) if true, upward/downward continue toaltdown_cont: (optional) if true, downward continue if needed, only used ifup_cont = trueget_lla: (optional) if true, convert map grid fromUTMtoLLAdz: (optional) upward continuation step size [m]down_max: (optional) maximum downward continuation distance [m]α: (optional) regularization parameter for downward continuationsave_h5: (optional) if true, savemapStomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
mapS:MapSorMapSdscalar magnetic anomaly map struct
MagNav.map_interpolate — Functionmap_interpolate(mapS::Union{MapS,MapSd,MapS3D}, type::Symbol = :cubic;
return_vert_deriv::Bool = false)Create map interpolation function, equivalent of griddedInterpolant in MATLAB. Optionally return vertical derivative map interpolation function, which is calculated using finite differences between map and 1 m upward continued map.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structtype: (optional) type of interpolation {:linear,:quad,:cubic}return_vert_deriv: (optional) if true, also returnder_map
Returns:
itp_map: map interpolation function (f(yy,xx)or (f(yy,xx,alt))der_map: ifreturn_vert_deriv = true, vertical derivative map interpolation function (f(yy,xx)or (f(yy,xx,alt))
MagNav.map_interpolate — Functionmap_interpolate(mapV::MapV, dim::Symbol = :X, type::Symbol = :cubic)Create map interpolation function, equivalent of griddedInterpolant in MATLAB.
Arguments:
mapV:MapVvector magnetic anomaly map structdim: map dimension to interpolate {:X,:Y,:Z}type: (optional) type of interpolation {:linear,:quad,:cubic}
Returns:
itp_map: map interpolation function (f(yy,xx))
MagNav.map_interpolate — Methodmap_interpolate(map_map::AbstractArray{T},
map_xx::AbstractVector{T},
map_yy::AbstractVector{T},
type::Symbol = :cubic,
map_alt::AbstractVector{T} = T[0]) where TCreate map interpolation function, equivalent of griddedInterpolant in MATLAB.
Arguments:
map_map:nyxnx(xnz) 2D or 3D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinatestype: (optional) type of interpolation {:linear,:quad,:cubic}map_alt: (optional) map altitude levels
Returns:
itp_map: map interpolation function (f(yy,xx)or (f(yy,xx,alt))
MagNav.map_resample — Methodmap_resample(mapS::MapS, mapS_new::MapS)Resample map with new grid.
Arguments:
mapS:MapSscalar magnetic anomaly map structmapS_new:MapSscalar magnetic anomaly map struct to use for resampling
Returns:
mapS:MapSscalar magnetic anomaly map struct, resampled
MagNav.map_resample — Methodmap_resample(mapS::MapS, map_xx_new::Vector, map_yy_new::Vector)Resample map with new grid.
Arguments:
mapS:MapSscalar magnetic anomaly map structmap_xx_new:nx_newmap x-direction (longitude) coordinates to use for resamplingmap_yy_new:ny_newmap y-direction (latitude) coordinates to use for resampling
Returns:
mapS:MapSscalar magnetic anomaly map struct, resampled
MagNav.map_resample — Methodmap_resample(map_map::Matrix, map_xx::Vector, map_yy::Vector,
map_mask::BitMatrix, map_xx_new::Vector, map_yy_new::Vector)Resample map with new grid.
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinatesmap_yy:nymap y-direction (latitude) coordinatesmap_masknyxnxmask for valid (not filled-in) map datamap_xx_new:nx_newmap x-direction (longitude) coordinates to use for resamplingmap_yy_new:ny_newmap y-direction (latitude) coordinates to use for resampling
Returns:
map_map:ny_newxnx_new2D gridded map data, resampledmap_mask:ny_newxnx_newmask for valid (not filled-in) map data, resampled
MagNav.map_trim — Functionmap_trim(map_map::Matrix,
map_xx::Vector = collect(axes(map_map,2)),
map_yy::Vector = collect(axes(map_map,1)),
pad::Int = 0,
xx_lim::Tuple = (-Inf,Inf),
yy_lim::Tuple = (-Inf,Inf),
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad,
silent::Bool = true)Trim map by removing large areas that are missing map data. Returns indices for the original map that produces the appropriate trimmed map.
Arguments:
map_map:nyxnx2D gridded map datamap_xx: (optional)nxmap x-direction (longitude) coordinates [rad] or [deg] or [m]map_yy: (optional)nymap y-direction (latitude) coordinates [rad] or [deg] or [m]pad: (optional) minimum padding (grid cells) along map edgesxx_lim: (optional) x-direction map limits(xx_min,xx_max)[rad] or [deg] or [m]yy_lim: (optional) y-direction map limits(yy_min,yy_max)[rad] or [deg] or [m]zone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}silent: (optional) if true, no print outs
Returns:
ind_xx:nxtrimmed x-direction map indicesind_yy:nytrimmed y-direction map indices
MagNav.map_trim — Methodmap_trim(map_map::Map, path::Path;
pad::Int = 0,
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad,
silent::Bool = true)Trim map by removing large areas far away from path. Do not use prior to upward continuation, as this causes in edge effect errors. Returns trimmed magnetic anomaly map struct.
Arguments:
map_map:Mapmagnetic anomaly map structpath:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structpad: (optional) minimum padding (grid cells) along map edgeszone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}silent: (optional) if true, no print outs
Returns:
map_map:Mapmagnetic anomaly map struct, trimmed
MagNav.map_trim — Methodmap_trim(map_map::Map;
pad::Int = 0,
xx_lim::Tuple = (-Inf,Inf),
yy_lim::Tuple = (-Inf,Inf),
zone_utm::Int = 18,
is_north::Bool = true,
map_units::Symbol = :rad,
silent::Bool = true)Trim map by removing large areas that are missing map data. Returns trimmed magnetic anomaly map struct.
Arguments:
map_map:Mapmagnetic anomaly map structpad: (optional) minimum padding (grid cells) along map edgesxx_lim: (optional) x-direction map limits(xx_min,xx_max)[rad] or [deg] or [m]yy_lim: (optional) y-direction map limits(yy_min,yy_max)[rad] or [deg] or [m]zone_utm: (optional) UTM zone, only used ifmap_units = :utmis_north: (optional) if true, map is in northern hemisphere, only used ifmap_units = :utmmap_units: (optional) map xx/yy units {:rad,:deg,:utm}silent: (optional) if true, no print outs
Returns:
map_map:Mapmagnetic anomaly map struct, trimmed
MagNav.map_utm2lla! — Methodmap_utm2lla!(map_map::Matrix, map_xx::Vector, map_yy::Vector,
alt, map_mask::BitMatrix;
map_info::String = "Map",
zone_utm::Int = 18,
is_north::Bool = true,
save_h5::Bool = false,
map_h5::String = "map_data.h5")Convert map grid from UTM to LLA.
Arguments:
map_map:nyxnx2D gridded map data onUTMgridmap_xx:nxmap x-direction (longitude) coordinates [m]map_yy:nymap y-direction (latitude) coordinates [m]alt: map altitude(s) or altitude map [m]map_mask:nyxnxmask for valid (not filled-in) map datamap_info: (optional) map informationzone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispheresave_h5: (optional) if true, save map data tomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
nothing:map_map,map_xx,map_yy, &map_mask(&alt) are mutated withLLAgridded map data
MagNav.map_utm2lla! — Methodmap_utm2lla!(mapS::Union{MapS,MapSd,MapS3D};
zone_utm::Int = 18,
is_north::Bool = true,
save_h5::Bool = false
map_h5::String = "map_data.h5")Convert map grid from UTM to LLA.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct onUTMgridzone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispheresave_h5: (optional) if true, savemapStomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
nothing:map,xx,yy, &mask(&alt) fields withinmapSare mutated withLLAgridded map data
MagNav.map_utm2lla — Methodmap_utm2lla(map_map::Matrix, map_xx::Vector, map_yy::Vector,
alt, map_mask::BitMatrix;
map_info::String = "Map",
zone_utm::Int = 18,
is_north::Bool = true,
save_h5::Bool = false,
map_h5::String = "map_data.h5")Convert map grid from UTM to LLA.
Arguments:
map_map:nyxnx2D gridded map data onUTMgridmap_xx:nxmap x-direction (longitude) coordinates [m]map_yy:nymap y-direction (latitude) coordinates [m]alt: map altitude(s) or altitude map [m]map_mask:nyxnxmask for valid (not filled-in) map datamap_info: (optional) map informationzone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispheresave_h5: (optional) if true, save map data tomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
map_map:nyxnx2D gridded map data onLLAgridmap_xx:nxmap x-direction (longitude) coordinates [rad]map_yy:nymap y-direction (latitude) coordinates [rad]map_mask:nyxnxmask for valid (not filled-in) map data onLLAgrid
MagNav.map_utm2lla — Methodmap_utm2lla(mapS::Union{MapS,MapSd,MapS3D};
zone_utm::Int = 18,
is_north::Bool = true,
save_h5::Bool = false
map_h5::String = "map_data.h5")Convert map grid from UTM to LLA.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct onUTMgridzone_utm: (optional) UTM zoneis_north: (optional) if true, map is in northern hemispheresave_h5: (optional) if true, savemapStomap_h5map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)
Returns:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct onLLAgrid
MagNav.mpf — Methodmpf(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas, dt, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
num_part = 1000,
thresh = 0.8,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Rao-Blackwellized (marginalized) particle filter (MPF) for airborne magnetic anomaly navigation. This simplified MPF works only with LINEAR dynamics. This allows the same Kalman filter covariance matrices to be used with each particle, simplifying the filter and reducing the computational load. It is especially suited for map-matching navigation in which there is a highly non-linear, non-Gaussian MEASUREMENT, but NOT non-linear dynamics. The filter also assumes NON-correlated measurements to speed up computation.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancenum_part: (optional) number of particlesthresh: (optional) resampling threshold fraction {0:1}baro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
filt_res:FILTresfilter results struct
MagNav.mpf — Methodmpf(ins::INS, meas, itp_mapS;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
num_part = 1000,
thresh = 0.8,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Rao-Blackwellized (marginalized) particle filter (MPF) for airborne magnetic anomaly navigation. This simplified MPF works only with LINEAR dynamics. This allows the same Kalman filter covariance matrices to be used with each particle, simplifying the filter and reducing the computational load. It is especially suited for map-matching navigation in which there is a highly non-linear, non-Gaussian MEASUREMENT, but NOT non-linear dynamics. The filter also assumes NON-correlated measurements to speed up computation.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancenum_part: (optional) number of particlesthresh: (optional) resampling threshold fraction {0:1}baro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
filt_res:FILTresfilter results struct
MagNav.nekf — Functionnekf(ins::INS, meas, itp_mapS,
x_nn::Matrix = meas[:,:],
m = Dense(1 => 1);
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Measurement noise covariance-adaptive neural extended Kalman filter (nEKF) for airborne magnetic anomaly navigation.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)m: neural network modelP0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
filt_res:FILTresfilter results struct
MagNav.nekf — Functionnekf(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas, dt, itp_mapS,
x_nn::Matrix = meas[:,:],
m = Dense(1 => 1);
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false)Measurement noise covariance-adaptive neural extended Kalman filter (nEKF) for airborne magnetic anomaly navigation.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)m: neural network modelP0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
filt_res:FILTresfilter results struct
MagNav.nekf_train — Methodnekf_train(lat, lon, alt, vn, ve, vd, fn, fe, fd, Cnb, meas, dt,
itp_mapS, x_nn::Matrix, y_nn::Matrix;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
η_adam = 0.1,
epoch_adam::Int = 10,
hidden::Int = 1,
activation::Function = swish,
l_window::Int = 50,
date = get_years(2020,185),
core::Bool = false)Train a measurement noise covariance-adaptive neural extended Kalman filter (nEKF) model for airborne magnetic anomaly navigation.
Arguments:
lat: latitude [rad]lon: longitude [rad]alt: altitude [m]vn: north velocity [m/s]ve: east velocity [m/s]vd: down velocity [m/s]fn: north specific force [m/s^2]fe: east specific force [m/s^2]fd: down specific force [m/s^2]Cnb: direction cosine matrix (body to navigation) [-]meas: scalar magnetometer measurement [nT]dt: measurement time step [s]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)y_nn:ytarget matrix for neural network ([latitude longitude])P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]η_adam: (optional) learning rate for Adam optimizerepoch_adam: (optional) number of epochs for Adam optimizerhidden: (optional) hidden layers & nodes (e.g.,[8,8]for 2 hidden layers, 8 nodes each)activation: (optional) activation functionrelu= rectified linear unitσ= sigmoid (logistic function)swish= self-gatedtanh= hyperbolic tan- run
plot_activation()for a visual
l_window: (optional) temporal window lengthdate: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
m: neural network model
MagNav.nekf_train — Methodnekf_train(ins::INS, meas, itp_mapS, x_nn::Matrix, y_nn::Matrix;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
η_adam = 0.1,
epoch_adam::Int = 10,
hidden::Int = 1,
activation::Function = swish,
l_window::Int = 50,
date = get_years(2020,185),
core::Bool = false)Train a measurement noise covariance-adaptive neural extended Kalman filter (nEKF) model for airborne magnetic anomaly navigation.
Arguments:
ins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x_nn:NxNfdata matrix for neural network (Nfis number of features)y_nn:ytarget matrix for neural network ([latitude longitude])P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]η_adam: (optional) learning rate for Adam optimizerepoch_adam: (optional) number of epochs for Adam optimizerhidden: (optional) hidden layers & nodes (e.g.,[8,8]for 2 hidden layers, 8 nodes each)activation: (optional) activation functionrelu= rectified linear unitσ= sigmoid (logistic function)swish= self-gatedtanh= hyperbolic tan- run
plot_activation()for a visual
l_window: (optional) temporal window lengthdate: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
m: neural network model
MagNav.nekf_train — Methodnekf_train(xyz::XYZ, ind, meas, itp_mapS, x::Matrix;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
η_adam = 0.1,
epoch_adam::Int = 10,
hidden::Int = 1,
activation::Function = swish,
l_window::Int = 50,
date = get_years(2020,185),
core::Bool = false)Train a measurement noise covariance-adaptive neural extended Kalman filter (nEKF) model for airborne magnetic anomaly navigation.
Arguments:
xyz:XYZflight data structind: selected data indicesmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))x:NxNfdata matrix (Nfis number of features)P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancebaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]η_adam: (optional) learning rate for Adam optimizerepoch_adam: (optional) number of epochs for Adam optimizerhidden: (optional) hidden layers & nodes (e.g.,[8,8]for 2 hidden layers, 8 nodes each)activation: (optional) activation functionrelu= rectified linear unitσ= sigmoid (logistic function)swish= self-gatedtanh= hyperbolic tan- run
plot_activation()for a visual
l_window: (optional) temporal window lengthdate: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurement
Returns:
m: neural network modeldata_norms: length-3tuple of data normalizations,(v_scale,x_bias,x_scale)
MagNav.norm_sets — Methodnorm_sets(train, val, test;
norm_type::Symbol = :standardize,
no_norm = falses(size(train,2)))Normalize (or standardize) features (columns) of training, validation, & testing data.
Arguments:
train:N_trainxNftraining dataval:N_valxNfvalidation datatest:N_testxNftesting datanorm_type: (optional) normalization type::standardize= Z-score normalization:normalize= min-max normalization:scale= scale by maximum absolute value, bias = 0:none= scale by 1, bias = 0
no_norm: (optional) length-NfBoolean indices of features to not be normalized
Returns:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalizedval:N_valxNfvalidation data, normalizedtest:N_testxNftesting data, normalized
MagNav.norm_sets — Methodnorm_sets(train, test;
norm_type::Symbol = :standardize,
no_norm = falses(size(train,2)))Normalize (or standardize) features (columns) of training & testing data.
Arguments:
train:N_trainxNftraining datatest:N_testxNftesting datanorm_type: (optional) normalization type::standardize= Z-score normalization:normalize= min-max normalization:scale= scale by maximum absolute value, bias = 0:none= scale by 1, bias = 0
no_norm: (optional) length-NfBoolean indices of features to not be normalized
Returns:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalizedtest:N_testxNftesting data, normalized
MagNav.norm_sets — Methodnorm_sets(train;
norm_type::Symbol = :standardize,
no_norm = falses(size(train,2)))Normalize (or standardize) features (columns) of training data.
Arguments:
train:N_trainxNftraining datanorm_type: (optional) normalization type::standardize= Z-score normalization:normalize= min-max normalization:scale= scale by maximum absolute value, bias = 0:none= scale by 1, bias = 0
no_norm: (optional) length-NfBoolean indices of features to not be normalized
Returns:
train_bias:1xNftraining data biases (means, mins, or zeros)train_scale:1xNftraining data scaling factors (std devs, maxs-mins, or ones)train:N_trainxNftraining data, normalized
MagNav.ottawa_area_maps — Functionottawa_area_maps(f::Union{String,Symbol} = "")Magnetic anomaly maps near Ottawa, Ontario, Canada, contains:
Eastern_395.h5: Eastern Ontario at 395 m HAEEastern_drape.h5: Eastern Ontario on drapeRenfrew_395.h5: Renfrew at 395 m HAERenfrew_555.h5: Renfrew at 555 m HAERenfrew_drape.h5: Renfrew on drapeHighAlt_5181.h5: High Altitude mini-survey (within Renfrew) at 5181 m HAEPerth_800.h5: Perth mini-survey (within Eastern Ontario) at 800 m HAE
NOTE: Missing map data within each map has been filled in (using k-nearest neighbors) so that the maps are fully filled. Care must be taken to not navigate in the filled-in areas, as this is not real data and only done for more accurate and consistent upward continuation of the maps. Use the map_check function with the desired map and flight path data to check if the map may be used without navigating into filled-in (artificial) areas.
Arguments:
f: (optional) name of data file (.h5extension optional)
Returns:
p: path of folder orfdata file
MagNav.ottawa_area_maps_gxf — Functionottawa_area_maps_gxf(f::Union{String,Symbol} = "")GXF versions of small magnetic anomaly maps near Ottawa, Ontario, Canada, contains:
HighAlt_Mag.gxf: High Altitude mini-survey (within Renfrew) at 5181 m HAEPerth_Mag.gxf: Perth mini-survey (within Eastern Ontario) at 800 m HAE
Arguments:
f: (optional) name of data file (_Mag&.gxfextension optional)
Returns:
p: path of folder orfdata file
MagNav.path2kml — Functionpath2kml(path::Path,
path_kml::String = "path.kml";
width::Int = 3,
color1::String = "",
color2::String = "00ffffff",
points::Bool = false)Create KML file of flight path for use with Google Earth.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structpath_kml: (optional) path/name of flight path KML file to save (.kmlextension optional)width: (optional) line widthcolor1: (optional) path colorcolor2: (optional) below-path colorpoints: (optional) if true, create points instead of line
Returns:
nothing:path_kmlis created
MagNav.path2kml — Functionpath2kml(lat::Vector, lon::Vector, alt::Vector,
path_kml::String = "path.kml";
path_units::Symbol = :rad,
width::Int = 3,
color1::String = "ff000000",
color2::String = "80000000",
points::Bool = false)Create KML file of flight path for use with Google Earth.
Arguments:
lat: latitude [rad] or [deg]lon: longitude [rad] or [deg]alt: altitude [m]path_kml: (optional) path/name of flight path KML file to save (.kmlextension optional)path_units: (optional)lat/lonunits {:rad,:deg}width: (optional) line widthcolor1: (optional) path colorcolor2: (optional) below-path colorpoints: (optional) if true, create points instead of line
Returns:
nothing:path_kmlis created
MagNav.plot_activation — Functionplot_activation(activation = [:relu,:σ,:swish,:tanh];
plot_deriv::Bool = false,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "act_func.png")Plot activation function(s) or their derivative(s).
Arguments:
activation: activation function(s) to plotrelu= rectified linear unitσ= sigmoid (logistic function)swish= self-gatedtanh= hyperbolic tan
plot_deriv: (optional) if true, plot activation function(s) derivative(s)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of activation function(s) or their derivative(s)
MagNav.plot_autocor — Functionplot_autocor(x::Vector, dt = 0.1, dt_max = 300.0;
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "autocor.png")Plot autocorrelation of data (e.g., actual - expected measurements). Prints out σ = standard deviation & τ = autocorrelation decay to e^-1 of x.
Arguments:
x: data vectordt: (optional) measurement time step [s]dt_max: (optional) maximum time step to evaluate [s]show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of autocorrelation ofx
MagNav.plot_basic — Functionplot_basic(tt::Vector, y::Vector, ind = trues(length(tt));
lab::String = "",
xlab::String = "time [min]",
ylab::String = "",
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "data_vs_time.png")Plot data vs time.
Arguments:
tt: length-Ntime vector [s]y: length-Ndata vectorind: (optional) selected data indiceslab: (optional) data (legend) labelxlab: (optional) x-axis labelylab: (optional) y-axis labelshow_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot ofyvstt
MagNav.plot_correlation — Functionplot_correlation(x::Vector, y::Vector,
xfeature::Symbol = :feature_1,
yfeature::Symbol = :feature_2;
lim::Real = 0,
dpi::Int = 200,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "$xfeature-$yfeature.png",
silent::Bool = true)Plot the correlation between 2 features.
Arguments:
x: x-axis datay: y-axis dataxfeature: x-axis feature nameyfeature: y-axis feature namelim: (optional) only plot if Pearson correlation coefficient >limdpi: (optional) dots per inch (image resolution)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)silent: (optional) if true, no print outs
Returns:
p1: plot ofyfeaturevsxfeaturecorrelation
MagNav.plot_correlation — Functionplot_correlation(xyz::XYZ,
xfeature::Symbol = :mag_1_c,
yfeature::Symbol = :mag_1_uc,
ind = trues(xyz.traj.N);
lim::Real = 0,
dpi::Int = 200,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "$xfeature-$yfeature.png",
silent::Bool = true)Plot the correlation between 2 features.
Arguments:
xyz:XYZflight data structxfeature: x-axis feature nameyfeature: y-axis feature nameind: (optional) selected data indiceslim: (optional) only plot if Pearson correlation coefficient >limdpi: (optional) dots per inch (image resolution)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)silent: (optional) if true, no print outs
Returns:
p1: plot ofyfeaturevsxfeaturecorrelation
MagNav.plot_correlation_matrix — Functionplot_correlation_matrix(xyz::XYZ, ind = trues(xyz.traj.N),
features_setup::Vector{Symbol} = [:mag_1_uc,:TL_A_flux_a];
terms = [:permanent],
sub_diurnal::Bool = false,
sub_igrf::Bool = false,
bpf_mag::Bool = false,
dpi::Int = 200,
Nmax::Int = 1000,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "correlation_matrix.png")Plot the correlation matrix for 2-5 features.
Arguments:
xyz:XYZflight data structind: selected data indicesfeatures_setup: vector of features to includeterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}sub_diurnal: (optional) if true, subtract diurnal from scalar magnetometer measurementssub_igrf: (optional) if true, subtract IGRF from scalar magnetometer measurementsbpf_mag: (optional) if true, bpf scalar magnetometer measurementsdpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of correlation matrix betweenfeatures(created fromfeatures_setup)
MagNav.plot_correlation_matrix — Methodplot_correlation_matrix(x::AbstractMatrix, features::Vector{Symbol};
dpi::Int = 200,
Nmax::Int = 1000,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "correlation_matrix.png")Plot the correlation matrix for 2-5 features.
Arguments:
x:NxNfdata matrix (Nfis number of features)features: length-Nffeature vector (including components of TLA, etc.)dpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of correlation matrix betweenfeatures
MagNav.plot_events! — Functionplot_events!(p1::Plot, t::Real, lab::String = "";
legend::Symbol = :outertopright)Plot in-flight event on an existing plot.
Arguments:
p1: plot (i.e., time series of data)t: time of in-flight eventlab: (optional) in-flight event (legend) labellegend: (optional) legend position (e.g.,:topleft,:outertopright)
Returns:
nothing: in-flight event is plotted onp1
MagNav.plot_events! — Methodplot_events!(p1::Plot, flight::Symbol, df_event::DataFrame;
keyword::String = "",
show_lab::Bool = true,
t0::Real = 0,
t_units::Symbol = :sec,
legend::Symbol = :outertopright)Plot in-flight event(s) on an existing plot.
Arguments:
p1: plot (i.e., time series of data)flight: flight name (e.g.,:Flt1001)df_event: lookup table (DataFrame) of in-flight events
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
tt | Real | time of event [s] |
event | String | event description |
keyword: (optional) keyword to search within events, case insensitiveshow_lab: (optional) if true, show in-flight event (legend) label(s)t0: (optional) time offset [t_units]t_units: (optional) time units {:sec,:min}legend: (optional) legend position (e.g.,:topleft,:outertopright)
Returns:
nothing: in-flight events are plotted onp1
MagNav.plot_events — Methodplot_events(p1::Plot, flight::Symbol, df_event::DataFrame;
keyword::String = "",
show_lab::Bool = true,
t0::Real = 0,
t_units::Symbol = :sec,
legend::Symbol = :outertopright)Plot in-flight event(s) on an existing plot.
Arguments:
p1: plot (i.e., time series of data)flight: flight name (e.g.,:Flt1001)df_event: lookup table (DataFrame) of in-flight events
| Field | Type | Description |
|---|---|---|
flight | Symbol | flight name (e.g., :Flt1001) |
tt | Real | time of event [s] |
event | String | event description |
keyword: (optional) keyword to search within events, case insensitiveshow_lab: (optional) if true, show in-flight event (legend) label(s)t0: (optional) time offset [t_units]t_units: (optional) time units {:sec,:min}legend: (optional) legend position (e.g.,:topleft,:outertopright)
Returns:
p2:p1with in-flight events
MagNav.plot_filt! — Methodplot_filt!(p1::Plot, traj::Traj, ins::INS, filt_out::FILTout;
dpi::Int = 200,
Nmax::Int = 5000,
show_plot::Bool = true,
save_plot::Bool = false)Plot flights paths on an existing plot.
Arguments:
p1: plot (i.e., map)traj:Trajtrajectory structins:INSinertial navigation system structfilt_out:FILToutfilter extracted output structdpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotssave_plot: (optional) if true, save plots with default file names
Returns:
nothing: flight paths are plotted onp1
MagNav.plot_filt — Methodplot_filt(traj::Traj, ins::INS, filt_out::FILTout;
dpi::Int = 200,
Nmax::Int = 5000,
plot_vel::Bool = false,
show_plot::Bool = true,
save_plot::Bool = false)Plot flights paths and latitudes & longitudes vs time.
Arguments:
traj:Trajtrajectory structins:INSinertial navigation system structfilt_out:FILToutfilter extracted output structdpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedplot_vel: (optional) if true, plot velocitiesshow_plot: (optional) if true, show plotssave_plot: (optional) if true, save plots with default file names
Returns:
p1: flight pathsp2: latitudes vs timep3: longitudes vs timep4: ifplot_vel = true, north velocities vs timep5: ifplot_vel = true, east velocities vs time
MagNav.plot_filt — Methodplot_filt(p1::Plot, traj::Traj, ins::INS, filt_out::FILTout;
dpi::Int = 200,
Nmax::Int = 5000,
plot_vel::Bool = false,
show_plot::Bool = true,
save_plot::Bool = false)Plot flights paths on an existing plot and latitudes & longitudes vs time.
Arguments:
p1: plot (i.e., map)traj:Trajtrajectory structins:INSinertial navigation system structfilt_out:FILToutfilter extracted output structdpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedplot_vel: (optional) if true, plot velocitiesshow_plot: (optional) if true, show plotssave_plot: (optional) if true, save plots with default file names
Returns:
p2:p1with flight pathsp3: latitudes vs timep4: longitudes vs timep5: ifplot_vel = true, north velocities vs timep6: ifplot_vel = true, east velocities vs time
MagNav.plot_filt_err — Methodplot_filt_err(traj::Traj, filt_out::FILTout, crlb_out::CRLBout;
dpi::Int = 200,
Nmax::Int = 5000,
plot_vel::Bool = false,
show_plot::Bool = true,
save_plot::Bool = false)Plot northing & easting errors vs time.
Arguments:
traj:Trajtrajectory structfilt_out:FILToutfilter extracted output structcrlb_out:CRLBoutCramér–Rao lower bound extracted output structdpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plottedplot_vel: (optional) if true, plot velocity errorsshow_plot: (optional) if true, show plotssave_plot: (optional) if true, save plots with default file names
Returns:
p1: northing errors vs timep2: easting errors vs timep3: ifplot_vel = true, north velocity errors vs timep4: ifplot_vel = true, east velocity errors vs time
MagNav.plot_frequency — Methodplot_frequency(xyz::XYZ;
ind = trues(xyz.traj.N),
field::Symbol = :mag_1_uc,
freq_type::Symbol = :PSD,
detrend_data::Bool = true,
window::Function = hamming,
dpi::Int = 200,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "PSD.png")Plot frequency data, either Welch power spectral density (PSD) or spectrogram.
Arguments:
xyz:XYZflight data structind: (optional) selected data indicesfield: (optional) data field inxyzto plotfreq_type: (optional) frequency plot type {:PSD,:psd,:spectrogram,:spec}detrend_data: (optional) if true, detrend plot datawindow: (optional) type of window useddpi: (optional) dots per inch (image resolution)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of Welch power spectral density (PSD) or spectrogram
MagNav.plot_mag — Methodplot_mag(xyz::XYZ;
ind = trues(xyz.traj.N),
detrend_data::Bool = false,
use_mags::Vector{Symbol} = [:all_mags],
vec_terms::Vector{Symbol} = [:all],
ylim::Tuple = (),
dpi::Int = 200,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "scalar_mags.png")Plot scalar or vector (fluxgate) magnetometer data from a given flight test.
Arguments:
xyz:XYZflight data structind: (optional) selected data indicesdetrend_data: (optional) if true, detrend plot datause_mags: (optional) scalar or vector (fluxgate) magnetometers to plot {:all_mags,:comp_magsor:mag_1_c,:mag_1_uc,:flux_a, etc.}:all_mags= all provided scalar magnetometer fields (e.g.,:mag_1_c,:mag_1_uc, etc.):comp_mags= provided compensation(s) between:mag_1_uc&:mag_1_c, etc.
vec_terms: (optional) vector magnetometer (fluxgate) terms to plot {:allor:x,:y,:z,:t}ylim: (optional) length-2plotylimits (ymin,ymax) [nT]dpi: (optional) dots per inch (image resolution)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of scalar or vector (fluxgate) magnetometer data
MagNav.plot_mag_c — Methodplot_mag_c(xyz::XYZ,xyz_comp::XYZ;
ind = trues(xyz.traj.N),
ind_comp = trues(xyz_comp.traj.N),
detrend_data::Bool = true,
λ = 0.025,
terms = [:permanent,:induced,:eddy],
pass1 = 0.1,
pass2 = 0.9,
fs = 10.0,
use_mags::Vector{Symbol} = [:all_mags],
use_vec::Symbol = :flux_a,
plot_diff::Bool = false,
plot_mag_1_uc::Bool = true,
plot_mag_1_c::Bool = true,
ylim = (),
dpi::Int = 200,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "scalar_mags_comp.png")Plot compensated magnetometer(s) data from a given flight test. Assumes Mag 1 (i.e., :mag_1_uc & :mag_1_c) is the best magnetometer (i.e., stinger).
Arguments:
xyz:XYZflight data structxyz_comp:XYZflight data struct to use for compensationind: (optional) selected data indicesind_comp: (optional) selected data indices to use for compensationdetrend_data: (optional) if true, detrend plot dataλ: (optional) ridge parameterterms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}pass1: (optional) filter first passband frequency [Hz]pass2: (optional) filter second passband frequency [Hz]fs: (optional) filter sampling frequency [Hz]use_mags: (optional) scalar or vector (fluxgate) magnetometers to plot {:all_mags, or:mag_1_uc, etc.}:all_mags= all uncompensated scalar magnetometer fields (e.g.,:mag_1_uc, etc.)
use_vec: (optional) vector magnetometer (fluxgate) to use for Tolles-LawsonAmatrix {:flux_a, etc.}plot_diff: (optional) if true, plot difference betweenprovidedcompensated data & compensated magsas performed hereplot_mag_1_uc: (optional) if true, plot mag1uc (uncompensated mag_1)plot_mag_1_c: (optional) if true, plot mag1c (compensated mag_1)ylim: (optional) length-2plotylimits (ymin,ymax) [nT]dpi: (optional) dots per inch (image resolution)show_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: plot of compensated magnetometer data
MagNav.plot_mag_map — Methodplot_mag_map(path::Path, mag, itp_mapS;
lab::String = "magnetometer",
order::Symbol = :magmap,
dpi::Int = 200,
Nmax::Int = 5000,
detrend_data::Bool = true,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "mag_vs_map.png")Plot scalar magnetometer measurements vs map values.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structmag: scalar magnetometer measurements [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))lab: (optional) magnetometer data (legend) labelorder: (optional) plotting order {:magmap,:mapmag}dpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plotteddetrend_data: (optional) if true, detrend plot datashow_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: scalar magnetometer measurements vs map values
MagNav.plot_mag_map_err — Methodplot_mag_map_err(path::Path, mag, itp_mapS;
lab::String = "",
dpi::Int = 200,
Nmax::Int = 5000,
detrend_data::Bool = true,
show_plot::Bool = true,
save_plot::Bool = false,
plot_png::String = "mag_map_err.png")Plot scalar magnetometer measurement vs map value errors.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structmag: scalar magnetometer measurements [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))lab: (optional) data (legend) labeldpi: (optional) dots per inch (image resolution)Nmax: (optional) maximum number of data points plotteddetrend_data: (optional) if true, detrend plot datashow_plot: (optional) if true, showp1save_plot: (optional) if true, savep1asplot_pngplot_png: (optional) plot file name to save (.pngextension optional)
Returns:
p1: scalar magnetometer measurement vs map value errors
MagNav.plot_map — Functionplot_map(map_map::Matrix,
map_xx::Vector = [],
map_yy::Vector = [];
clims::Tuple = (),
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 6*dpi,
legend::Bool = true,
axis::Bool = true,
map_color::Symbol = :usgs,
bg_color::Symbol = :white,
map_units::Symbol = :rad,
plot_units::Symbol = :deg,
b_e::AbstractBackend = gr())Plot map.
Arguments:
map_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinates [rad] or [deg]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg]clims: (optional) length-2colorbar limits(cmin,cmax)dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plotted (per axis)legend: (optional) if true, show legendaxis: (optional) if true, show axesmap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}bg_color: (optional) background colormap_units: (optional) map xx/yy units {:rad,:deg}plot_units: (optional) plot xx/yy units {:rad,:deg,:m}b_e: (optional) plotting backend
Returns:
p1: plot of map
MagNav.plot_map! — Functionplot_map!(p1::Plot, map_map::Matrix,
map_xx::Vector = [],
map_yy::Vector = [];
clims::Tuple = (),
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 6*dpi,
legend::Bool = true,
axis::Bool = true,
map_color::Symbol = :usgs,
bg_color::Symbol = :white,
map_units::Symbol = :rad,
plot_units::Symbol = :deg,
b_e::AbstractBackend = gr())Plot map on an existing plot.
Arguments:
p1: plotmap_map:nyxnx2D gridded map datamap_xx:nxmap x-direction (longitude) coordinates [rad] or [deg]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg]clims: (optional) length-2colorbar limits(cmin,cmax)dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plotted (per axis)legend: (optional) if true, show legendaxis: (optional) if true, show axesmap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}bg_color: (optional) background colormap_units: (optional) map xx/yy units {:rad,:deg}plot_units: (optional) plot xx/yy units {:rad,:deg,:m}b_e: (optional) plotting backend
Returns:
nothing: map is plotted onp1
MagNav.plot_map! — Methodplot_map!(p1::Plot, p2::Plot, p3::Plot, mapV::MapV;
use_mask::Bool = true,
clims::Tuple = (),
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 6*dpi,
legend::Bool = true,
axis::Bool = true,
map_color::Symbol = :usgs,
bg_color::Symbol = :white,
map_units::Symbol = :rad,
plot_units::Symbol = :deg
b_e::AbstractBackend = gr())Plot map on an existing plot.
Arguments:
p1: plotp2: plotp3: plotmapV:MapVvector magnetic anomaly map structuse_mask: (optional) if true, applymapVmask to mapclims: (optional) length-2colorbar limits(cmin,cmax)dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plotted (per axis)legend: (optional) if true, show legendaxis: (optional) if true, show axesmap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}bg_color: (optional) background colormap_units: (optional) map xx/yy units {:rad,:deg}plot_units: (optional) plot xx/yy units {:rad,:deg,:m}b_e: (optional) plotting backend
Returns:
nothing:mapXis plotted onp1nothing:mapYis plotted onp2nothing:mapZis plotted onp3
MagNav.plot_map! — Methodplot_map!(p1::Plot, mapS::Union{MapS,MapSd,MapS3D};
use_mask::Bool = true,
clims::Tuple = (),
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 6*dpi,
legend::Bool = true,
axis::Bool = true,
map_color::Symbol = :usgs,
bg_color::Symbol = :white,
map_units::Symbol = :rad,
plot_units::Symbol = :deg
b_e::AbstractBackend = gr())Plot map on an existing plot.
Arguments:
p1: plotmapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map structuse_mask: (optional) if true, applymapSmask to mapclims: (optional) length-2colorbar limits(cmin,cmax)dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plotted (per axis)legend: (optional) if true, show legendaxis: (optional) if true, show axesmap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}bg_color: (optional) background colormap_units: (optional) map xx/yy units {:rad,:deg}plot_units: (optional) plot xx/yy units {:rad,:deg,:m}b_e: (optional) plotting backend
Returns:
nothing: map is plotted onp1
MagNav.plot_map — Methodplot_map(map_map::Map;
use_mask::Bool = true,
clims::Tuple = (),
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 6*dpi,
legend::Bool = true,
axis::Bool = true,
map_color::Symbol = :usgs,
bg_color::Symbol = :white,
map_units::Symbol = :rad,
plot_units::Symbol = :deg,
b_e::AbstractBackend = gr())Plot map.
Arguments:
map_map:Mapmagnetic anomaly map structuse_mask: (optional) if true, applymap_mapmask to mapclims: (optional) length-2colorbar limits(cmin,cmax)dpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plotted (per axis)legend: (optional) if true, show legendaxis: (optional) if true, show axesmap_color: (optional) filled contour color scheme {:usgs,:gray,:gray1,:gray2,:plasma,:magma}bg_color: (optional) background colormap_units: (optional) map xx/yy units {:rad,:deg}plot_units: (optional) plot xx/yy units {:rad,:deg,:m}b_e: (optional) plotting backend
Returns:
p1: plot of map (ifmap_map isa MapV,mapX)p2: ifmap_map isa MapV,mapYp3: ifmap_map isa MapV,mapZ
MagNav.plot_path — Functionplot_path(path::Path, ind = trues(path.N);
lab::String = "",
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = true,
path_color::Symbol = :ignore)Plot flight path.
Arguments:
path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indiceslab: (optional) data (legend) labeldpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
p1: plot of flight path
MagNav.plot_path — Functionplot_path(p1::Plot, path::Path, ind = trues(path.N);
lab::String = "",
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = false,
path_color::Symbol = :ignore)Plot flight path on an existing plot.
Arguments:
p1: plot (i.e., map)path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indiceslab: (optional) data (legend) labelNmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
p2:p1with flight path
MagNav.plot_path! — Functionplot_path!(p1::Plot, path::Path, ind = trues(path.N);
lab::String = "",
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = false,
path_color::Symbol = :ignore)Plot flight path on an existing plot.
Arguments:
p1: plot (i.e., map)path:Pathstruct, i.e.,Trajtrajectory struct,INSinertial navigation system struct, orFILToutfilter extracted output structind: (optional) selected data indiceslab: (optional) data (legend) labelNmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
nothing: flight path is plotted onp1
MagNav.plot_path! — Methodplot_path!(p1::Plot, lat, lon;
lab::String = "",
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = false,
path_color::Symbol = :ignore)Plot flight path on an existing plot.
Arguments:
p1: plot (i.e., map)lat: latitude [rad]lon: longitude [rad]lab: (optional) data (legend) labelNmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
nothing: flight path is plotted onp1
MagNav.plot_path — Methodplot_path(lat, lon;
lab::String = "",
dpi::Int = 200,
margin::Int = 2,
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = true,
path_color::Symbol = :ignore)Plot flight path.
Arguments:
lat: latitude [rad]lon: longitude [rad]lab: (optional) data (legend) labeldpi: (optional) dots per inch (image resolution)margin: (optional) margin around plot [mm]Nmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
p1: plot of flight path
MagNav.plot_path — Methodplot_path(p1::Plot, lat, lon;
lab::String = "",
Nmax::Int = 5000,
show_plot::Bool = true,
zoom_plot::Bool = false,
path_color::Symbol = :ignore)Plot flight path on an existing plot.
Arguments:
p1: plot (i.e., map)lat: latitude [rad]lon: longitude [rad]lab: (optional) data (legend) labelNmax: (optional) maximum number of data points plottedshow_plot: (optional) if true, show plotzoom_plot: (optional) if true, zoom plot onto flight pathpath_color: (optional) path color {:ignore,:black,:gray,:red,:orange,:yellow,:green,:cyan,:blue,:purple}
Returns:
p2:p1with flight path
MagNav.plot_shapley — Functionplot_shapley(df_shap, baseline_shap,
range_shap::UnitRange = UnitRange(axes(df_shap,1));
title::String = "features $range_shap",
dpi::Int = 200)Plot horizontal bar graph of feature importance (Shapley effects).
Arguments:
df_shap: DataFrame of Shapley effects
| Field | Type | Description |
|---|---|---|
feature_name | Symbol | feature name |
mean_effect | Real | mean Shapley effect |
baseline_shap: intercept of Shapley effectsrange_shap: (optional) range of Shapley effects to plot (limit to length ~20)dpi: (optional) dots per inch (image resolution)title: (optional) plot title
Returns:
p1: plot of Shapley effects
MagNav.plsr_fit — Functionplsr_fit(x, y, k::Int = size(x,2), no_norm = falses(size(x,2));
data_norms::Tuple = (zeros(1,1),zeros(1,1),[0.0],[0.0]),
l_segs::Vector = [length(y)],
return_set::Bool = false,
silent::Bool = false)Fit a multi-input, multi-output (MIMO) partial least squares regression (PLSR) model to data with a specified output dimension. PLSR is a type of regularized linear regression where the number of components controls the strength of the regularization.
Arguments:
x:NxNfdata matrix (Nfis number of features)y: length-Ntarget vectork: (optional) number of componentsno_norm: (optional) length-NfBoolean indices of features to not be normalizeddata_norms: (optional) length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)l_segs: (optional) length-N_linesvector of lengths oflines, sum(l_segs) =Nreturn_set: (optional) if true, returncoef_setinstead of other outputssilent: (optional) if true, no print outs
Returns:
model: length-2tuple of PLSR-based model, (length-Nfcoefficients, bias=0)data_norms: length-4tuple of data normalizations,(x_bias,x_scale,y_bias,y_scale)y_hat: length-Nprediction vectorerr: length-Nmean-corrected (per line) errorcoef_set: ifreturn_set = true, set of coefficients (sizeNfxNyxk)
MagNav.predict_rnn_full — Methodpredict_rnn_full(m, x)Apply model m to data matrix x.
Arguments:
m: recurrent neural network modelx:NxNfdata matrix (Nfis number of features)
Returns:
y_hat: length-Nprediction vector
MagNav.predict_rnn_windowed — Methodpredict_rnn_windowed(m, x, l_window::Int)Apply model m to data matrix x with sliding window of length-l_window.
Arguments:
m: recurrent neural network modelx:NxNfdata matrix (Nfis number of features)l_window: temporal window length
Returns:
y_hat: length-Nprediction vector
MagNav.project_body_field_to_2d_igrf — Methodproject_body_field_to_2d_igrf(vec_body, igrf_nav, Cnb)Project a body frame vector onto a 2D plane defined by the direction of the IGRF and a tangent vector to the Earth ellipsoid, which is computed by taking the cross product of the IGRF with the upward direction. Returns a 2D vector whose components describe the amount of the body field that is in alignment with the Earth field and an orthogonal direction to the Earth field (roughly to the east).
Arguments:
vec_body: vector in body frame (e.g., aircraft induced field)igrf_nav: IGRF unit vector in navigation frameCnb:3x3xNdirection cosine matrix (body to navigation) [-]
Returns:
v_out: 2D vector whose components illustrate projection onto the Earth field and an orthogonal component
MagNav.psd — Methodpsd(map_map::Matrix, dx, dy)Power spectral density of a potential field (i.e., magnetic anomaly field) map. Uses the Fast Fourier Transform to determine the spectral energy distribution across the radial wavenumbers (spatial frequencies) in the Fourier transform.
Arguments:
map_map:nyxnx2D gridded map datadx: x-direction map step size [m]dy: y-direction map step size [m]
Returns:
map_psd:nyxnxpower spectral density of 2D gridded map datakx:nyxnxx-direction radial wavenumberky:nyxnxy-direction radial wavenumber
MagNav.psd — Methodpsd(mapS::Union{MapS,MapSd,MapS3D})Power spectral density of a potential field (i.e., magnetic anomaly field) map. Uses the Fast Fourier Transform to determine the spectral energy distribution across the radial wavenumbers (spatial frequencies) in the Fourier transform.
Arguments:
mapS:MapS,MapSd, orMapS3Dscalar magnetic anomaly map struct
Returns:
map_psd:nyxnxpower spectral density of 2D gridded map datakx:nyxnxx-direction radial wavenumberky:nyxnxy-direction radial wavenumber
MagNav.run_filt — Functionrun_filt(traj::Traj, ins::INS, meas, itp_mapS, filt_type::Symbol = :ekf;
P0 = create_P0(),
Qd = create_Qd(),
R = 1.0,
num_part = 1000,
thresh = 0.8,
baro_tau = 3600.0,
acc_tau = 3600.0,
gyro_tau = 3600.0,
fogm_tau = 600.0,
date = get_years(2020,185),
core::Bool = false,
map_alt = 0,
x_nn = nothing,
m = nothing,
y_norms = nothing,
terms = [:permanent,:induced,:eddy,:bias],
flux::MagV = MagV([0.0],[0.0],[0.0],[0.0]),
x0_TL = ones(eltype(P0),19),
extract::Bool = true,
run_crlb::Bool = true)Run navigation filter and optionally compute Cramér–Rao lower bound (CRLB).
Arguments:
traj:Trajtrajectory structins:INSinertial navigation system structmeas: scalar magnetometer measurement [nT]itp_mapS: scalar map interpolation function (f(lat,lon)orf(lat,lon,alt))filt_type: (optional) filter type {:ekf,:mpf}P0: (optional) initial covariance matrixQd: (optional) discrete time process/system noise matrixR: (optional) measurement (white) noise variancenum_part: (optional) number of particles, only used forfilt_type = :mpfthresh: (optional) resampling threshold fraction {0:1}, only used forfilt_type = :mpfbaro_tau: (optional) barometer time constant [s]acc_tau: (optional) accelerometer time constant [s]gyro_tau: (optional) gyroscope time constant [s]fogm_tau: (optional) FOGM catch-all time constant [s]date: (optional) measurement date (decimal year) for IGRF [yr]core: (optional) if true, include core magnetic field in measurementmap_alt: (optional) map altitude [m]x_nn: (optional)NxNfdata matrix for neural network (Nfis number of features)m: (optional) neural network modely_norms: (optional) tuple ofynormalizations, i.e.,(y_bias,y_scale)terms: (optional) Tolles-Lawson terms to use {:permanent,:induced,:eddy,:bias}flux: (optional)MagVvector magnetometer measurement structx0_TL: (optional) initial Tolles-Lawson coefficient statesextract: (optional) if true, extract output structsrun_crlb: (optional) if true, compute the Cramér–Rao lower bound (CRLB)
Returns:
- if
extract = true&run_crlb = truecrlb_out:CRLBoutCramér–Rao lower bound extracted output structins_out:INSoutinertial navigation system extracted output structfilt_out:FILToutfilter extracted output struct
- if
extract = true&run_crlb = falsefilt_out:FILToutfilter extracted output struct
- if
extract = false&run_crlb = truefilt_res:FILTresfilter results structcrlb_P: Cramér–Rao lower bound non-linear covariance matrix
- if
extract = false&run_crlb = falsefilt_res:FILTresfilter results struct
MagNav.run_filt — Methodrun_filt(traj::Traj, ins::INS, meas, itp_mapS,
filt_type::Vector{Symbol}; ...)Run multiple filter models and print results (nothing returned).
Arguments:
filt_type: multiple filter types, e.g., [:ekf,:ekf_online_nn]
MagNav.save_comp_params — Functionsave_comp_params(comp_params::CompParams,
comp_params_bson::String = "comp_params.bson")Save aeromagnetic compensation parameters to BSON file.
Arguments:
comp_params:CompParamsaeromagnetic compensation parameters struct, either:NNCompParams: neural network-based aeromagnetic compensation parameters structLinCompParams: linear aeromagnetic compensation parameters struct
comp_params_bson: (optional) path/name of aeromagnetic compensation parameters BSON file to save (.bsonextension optional)
Returns:
nothing:comp_params_bsonis created
MagNav.save_map — Functionsave_map(map_map, map_xx, map_yy, map_alt, map_h5::String = "map_data.h5";
map_info::String = "Map",
map_mask::BitArray = falses(1,1),
map_border::Matrix = zeros(eltype(map_alt),1,1),
map_units::Symbol = :rad,
file_units::Symbol = :deg)Save map data to HDF5 file. Maps are typically saved in :deg units, while :rad is used internally.
Arguments:
map_map:nyxnx(xnz) 2D or 3D gridded map datamap_xx:nxmap x-direction (longitude) coordinates [rad] or [deg]map_yy:nymap y-direction (latitude) coordinates [rad] or [deg]map_alt: map altitude(s) ornyxnx2D gridded altitude map data [m]map_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)map_info: (optional) map informationmap_mask: (optional)nyxnx(xnz) mask for valid (not filled-in) map datamap_border: (optional) [xx yy] border for valid (not filled-in) map data [rad] or [deg]map_units: (optional) map xx/yy units used inmap_xx&map_yy{:rad,:deg}file_units: (optional) map xx/yy units to use inmap_h5{:rad,:deg}
Returns:
nothing:map_h5is created
MagNav.save_map — Functionsave_map(map_map::Map, map_h5::String = "map_data.h5";
map_border::Matrix = zeros(eltype(map_map.alt),1,1),
map_units::Symbol = :rad,
file_units::Symbol = :deg)Save map data to HDF5 file. Maps are typically saved in :deg units, while :rad is used internally.
Arguments:
map_map:Mapmagnetic anomaly map structmap_h5: (optional) path/name of map data HDF5 file to save (.h5extension optional)map_border: (optional) [xx yy] border for valid (not filled-in) map data [rad] or [deg]map_units: (optional) map xx/yy units used inmap_map{:rad,:deg}file_units: (optional) map xx/yy units to use inmap_h5{:rad,:deg}
Returns:
nothing:map_h5is created
MagNav.sgl_2020_train — Functionsgl_2020_train(f::Union{String,Symbol} = "")Flight data from the 2020 SGL flight data collection - training portion. Collected from 20-Jun-2020 to 07-Jul-2020 near Ottawa, Ontario, Canada by Sander Geophysics Ltd. (SGL) using a Cessna Grand Caravan. Contains:
Flt1002_train.h5Flt1003_train.h5Flt1004_train.h5Flt1005_train.h5Flt1006_train.h5Flt1007_train.h5
Arguments:
f: (optional) name of data file (_train&.h5extension optional)
Returns:
p: path of folder orfdata file
MagNav.sgl_2021_train — Functionsgl_2021_train(f::Union{String,Symbol} = "")Flight data from the 2021 SGL flight data collection - training portion. Collected from 13-Dec-2021 to 05-Jan-2022 near Ottawa, Ontario, Canada by Sander Geophysics Ltd. (SGL) using a Cessna Grand Caravan. Contains:
Flt2001_train.h5Flt2002_train.h5Flt2004_train.h5Flt2005_train.h5Flt2006_train.h5Flt2007_train.h5Flt2008_train.h5Flt2015_train.h5Flt2016_train.h5Flt2017_train.h5
Arguments:
f: (optional) name of data file (_train&.h5extension optional)
Returns:
p: path of folder orfdata file
MagNav.sparse_group_lasso — Functionsparse_group_lasso(m::Chain, α=1)Get the sparse group Lasso term for sparse-input regularization, which is the combined L1 & L2 norm of the first-layer neural network weights corresponding to each input feature.
Reference: Feng & Simon, Sparse-Input Neural Networks for High-dimensional Nonparametric Regression and Classification, 2017 (pg. 4).
Arguments:
m: neural network modelα: (optional) Lasso (α=0) vs group Lasso (α=1) balancing parameter {0:1}
Returns:
w_norm: sparse group Lasso term
MagNav.sparse_group_lasso — Functionsparse_group_lasso(weights::Params, α=1)Get the sparse group Lasso term for sparse-input regularization, which is the combined L1 & L2 norm of the first-layer neural network weights corresponding to each input feature.
Reference: Feng & Simon, Sparse-Input Neural Networks for High-dimensional Nonparametric Regression and Classification, 2017 (pg. 4).
Arguments:
weights: neural network model weightsα: (optional) Lasso (α=0) vs group Lasso (α=1) balancing parameter {0:1}
Returns:
w_norm: sparse group Lasso term
MagNav.upward_fft — Methodupward_fft(map_map::Map, alt; expand::Bool = true, α = 0)Upward continuation of a potential field (i.e., magnetic anomaly field) map. Uses the Fast Fourier Transform (FFT) to convert the map to the frequency domain, applies an upward continuation filter, and uses the inverse FFT to convert the map back to the spatial domain. Optionally expands the map temporarily with periodic padding. Downward continuation may be performed to a limited degree as well, but be careful, as this is generally unstable and amplifies high frequencies (i.e., noise).
Reference: Blakely, Potential Theory in Gravity and Magnetic Applications, 2009, Chapter 12 & Appendix B (pg. 315-317 & 402).
Arguments:
map_map:Mapmagnetic anomaly map structalt: target upward continuation altitude(s) [m]expand: (optional) if true, expand map temporarily to reduce edge effectsα: (optional) regularization parameter for downward continuation
Returns:
map_map:Mapmagnetic anomaly map struct, upward/downward continued (MapSwithaltvector =>MapS3D)
MagNav.upward_fft — Methodupward_fft(map_map::Matrix, dx, dy, dz; expand::Bool = true, α = 0)Upward continuation of a potential field (i.e., magnetic anomaly field) map. Uses the Fast Fourier Transform (FFT) to convert the map to the frequency domain, applies an upward continuation filter, and uses the inverse FFT to convert the map back to the spatial domain. Optionally expands the map temporarily with periodic padding. Downward continuation may be performed to a limited degree as well, but be careful, as this is generally unstable and amplifies high frequencies (i.e., noise).
Reference: Blakely, Potential Theory in Gravity and Magnetic Applications, 2009, Chapter 12 & Appendix B (pg. 315-317 & 402).
Arguments:
map_map:nyxnx2D gridded map datadx: x-direction map step size [m]dy: y-direction map step size [m]dz:nzz-direction upward/downward continuation distance(s) [m]expand: (optional) if true, expand map temporarily to reduce edge effectsα: (optional) regularization parameter for downward continuation
Returns:
map_map:nyxnx(xnz) 2D or 3D gridded map data, upward/downward continued
MagNav.vector_fft — Methodvector_fft(map_map::Matrix, dx, dy, D, I)Get potential field (i.e., magnetic anomaly field) map vector components using declination and inclination.
Arguments:
map_map:nyxnx2D gridded map datadx: x-direction map step size [m]dy: y-direction map step size [m]D: map declination (Earth core field) [rad]I: map inclination (Earth core field) [rad]
Returns:
Bx, By, Bz: map vector components
MagNav.xyz2h5 — Methodxyz2h5(data::Array, xyz_h5::String, flight::Symbol;
tt_sort::Bool = true,
lines::Vector = [()],
lines_type::Symbol = :exclude)MagNav.xyz2h5 — Methodxyz2h5(xyz_xyz::String, xyz_h5::String, flight::Symbol;
lines::Vector = [()],
lines_type::Symbol = :exclude,
tt_sort::Bool = true,
downsample_160::Bool = true,
return_data::Bool = false)Convert SGL flight data file from .xyz to HDF5.
- Valid for SGL flights:
:Flt1001:Flt1002:Flt1003:Flt1004_1005:Flt1004:Flt1005:Flt1006:Flt1007:Flt1008:Flt1009:Flt1001_160Hz:Flt1002_160Hz:Flt2001_2017
May take 1+ hr for 1+ GB files. For reference, a 1.23 GB file took 46.8 min to process using a 64 GB MacBook Pro.
Arguments:
xyz_xyz: path/name of flight data .xyz file (.xyzextension optional)xyz_h5: path/name of flight data HDF5 file to save (.h5extension optional)flight: flight name (e.g.,:Flt1001)lines: (optional) selected line number(s) to ONLY include or exclude, must be a vector of 3-element (line,start_time,end_time) tuple(s)lines_type: (optional) whether to ONLY:include(i.e., to generate testing data) or:exclude(i.e., to generate training data)linestt_sort: (optional) if true, sort data by time (instead of line)downsample_160: (optional) if true, downsample 160 Hz data to 10 Hz (only for 160 Hz data files)return_data: (optional) if true, returndatainstead of creatingxyz_h5
Returns:
data: ifreturn_data = true, internal data matrix