├── +matmap3d ├── private │ ├── EOP-All.h5 │ ├── mustBeEqualSize.m │ ├── mustBeSizeOrEmpty.m │ ├── Frac.m │ ├── PoleMatrix.m │ ├── R_y.m │ ├── gast.m │ ├── R_x.m │ ├── R_z.m │ ├── GHAMatrix.m │ ├── MeanObliquity.m │ ├── NutMatrix.m │ ├── gmst.m │ ├── EqnEquinox.m │ ├── timediff.m │ ├── PrecMatrix.m │ ├── Mjday.m │ ├── load_eop.m │ ├── IERS.m │ ├── SAT_const.m │ └── NutAngles.m ├── R3.m ├── wgs84Ellipsoid.m ├── get_radius_normal.m ├── aer2ned.m ├── eci2aer.m ├── eci2ecef_naive.m ├── ecef2eci_naive.m ├── greenwichsrt.m ├── aer2enu.m ├── enu2uvw.m ├── ecef2ned.m ├── enu2ecefv.m ├── ecef2enuv.m ├── ecef2enu.m ├── aer2eci.m ├── geodetic2ned.m ├── enu2ecef.m ├── enu2geodetic.m ├── ecef2aer.m ├── enu2aer.m ├── geodetic2enu.m ├── aer2geodetic.m ├── geodetic2ecef.m ├── geodetic2aer.m ├── aer2ecef.m ├── ecef2geodetic.m ├── eci2ecef.m ├── @referenceEllipsoid │ └── referenceEllipsoid.m ├── ecef2eci.m ├── lookAtSpheroid.m ├── vreckon.m └── vdist.m ├── .gitignore ├── .github └── workflows │ ├── composite-install-matlab │ └── action.yml │ ├── composite-buildtool │ └── action.yml │ ├── ci.yml │ └── publish.yml ├── octave └── juliandate.m ├── codemeta.json ├── LICENSE ├── private └── publish_gen_index_html.m ├── README.md └── test └── TestUnit.m /+matmap3d/private/EOP-All.h5: -------------------------------------------------------------------------------- https://raw.githubusercontent.com/geospace-code/matmap3d/HEAD/+matmap3d/private/EOP-All.h5 -------------------------------------------------------------------------------- /.gitignore: -------------------------------------------------------------------------------- 1 | docs/ 2 | resources/ 3 | *.asv 4 | *.m~ 5 | code-coverage.xml 6 | .buildtool/ 7 | *.sarif 8 | *.xml 9 | -------------------------------------------------------------------------------- /+matmap3d/R3.m: -------------------------------------------------------------------------------- 1 | function A = R3(x) 2 | %% R3 rotation matrix for ECI 3 | A = [cos(x), sin(x), 0; 4 | -sin(x), cos(x), 0; 5 | 0, 0, 1]; 6 | end 7 | -------------------------------------------------------------------------------- /+matmap3d/private/mustBeEqualSize.m: -------------------------------------------------------------------------------- 1 | function mustBeEqualSize(a,b) 2 | 3 | sa = size(a); 4 | sb = size(b); 5 | if ~isequal(sa, sb) 6 | throwAsCaller(MException('MATLAB:validators:mustBeEqualSize', 'Size of inputs [%s] [%s] must equal each other', num2str(sa), num2str(sb))); 7 | end 8 | 9 | end 10 | -------------------------------------------------------------------------------- /.github/workflows/composite-install-matlab/action.yml: -------------------------------------------------------------------------------- 1 | name: "install-matlab" 2 | 3 | runs: 4 | 5 | using: "composite" 6 | 7 | steps: 8 | 9 | - name: Install MATLAB 10 | uses: matlab-actions/setup-matlab@v2 11 | with: 12 | release: ${{ matrix.release }} 13 | cache: true 14 | -------------------------------------------------------------------------------- /+matmap3d/private/mustBeSizeOrEmpty.m: -------------------------------------------------------------------------------- 1 | function mustBeSizeOrEmpty(a, s) 2 | 3 | if ~isempty(a) 4 | if isvector(a) && sum(s ~= 1) == 1 5 | ok = numel(a) == prod(s); 6 | else 7 | ok = isequal(size(a), s); 8 | end 9 | if ~ok 10 | throwAsCaller(MException('MATLAB:validators:mustBeEqualSize','Input must be size or empty')) 11 | end 12 | end 13 | 14 | end 15 | -------------------------------------------------------------------------------- /+matmap3d/private/Frac.m: -------------------------------------------------------------------------------- 1 | %% FRAC Fractional part of a number (y=x-[x]) 2 | % 3 | % Last modified: 2018/01/27 Meysam Mahooti 4 | % 5 | % Meysam Mahooti (2025). 6 | % ECI2ECEF & ECEF2ECI Transformations 7 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 8 | % MATLAB Central File Exchange. 9 | 10 | function [res] = Frac(x) 11 | 12 | res = x-floor(x); 13 | 14 | end 15 | -------------------------------------------------------------------------------- /.github/workflows/composite-buildtool/action.yml: -------------------------------------------------------------------------------- 1 | name: "matlab-buildtool" 2 | 3 | runs: 4 | 5 | using: "composite" 6 | 7 | steps: 8 | 9 | - name: Check Matlab Syntax 10 | if: ${{ matrix.release >= 'R2023b' }} 11 | uses: matlab-actions/run-build@v2 12 | with: 13 | tasks: check 14 | 15 | - name: Test Matlab Code 16 | if: ${{ matrix.release >= 'R2023b' }} 17 | uses: matlab-actions/run-build@v2 18 | with: 19 | tasks: test 20 | -------------------------------------------------------------------------------- /+matmap3d/wgs84Ellipsoid.m: -------------------------------------------------------------------------------- 1 | %% WGS84ELLIPSOID - generate a WGS84 referenceEllipsoid 2 | % 3 | %%% Inputs 4 | % 5 | % * lengthUnit: currently not implemented, for compatibility either nargin=0 or 'm' or 'meter' is accepted without conversion 6 | % 7 | %%% Outputs 8 | % 9 | % * E: referenceEllipsoid 10 | function E = wgs84Ellipsoid(lengthUnit) 11 | if nargin < 1 12 | lengthUnit = 'm'; 13 | end 14 | 15 | E = matmap3d.referenceEllipsoid('wgs84', lengthUnit); 16 | 17 | end 18 | -------------------------------------------------------------------------------- /+matmap3d/get_radius_normal.m: -------------------------------------------------------------------------------- 1 | %% GET_RADIUS_NORMAL normal along the prime vertical section ellipsoidal radius of curvature 2 | % 3 | %%% Inputs 4 | % * lat: geodetic latitude in Radians 5 | % * ell: referenceEllipsoid 6 | % 7 | %%% Outputs 8 | % * N: normal along the prime vertical section ellipsoidal radius of curvature, at a given geodetic latitude. 9 | function N = get_radius_normal(lat, E) 10 | 11 | if isempty(E) 12 | E = matmap3d.wgs84Ellipsoid(); 13 | end 14 | 15 | N = E.SemimajorAxis^2 ./ sqrt( E.SemimajorAxis^2 .* cos(lat).^2 + E.SemiminorAxis^2 .* sin(lat).^2 ); 16 | end 17 | -------------------------------------------------------------------------------- /+matmap3d/private/PoleMatrix.m: -------------------------------------------------------------------------------- 1 | %% POLEMATRIX Transformation from pseudo Earth-fixed to Earth-fixed coordinates for a given date 2 | % 3 | % Input: 4 | % Pole coordinte(xp,yp) 5 | % 6 | % Output: 7 | % PoleMat Pole matrix 8 | % 9 | % Last modified: 2018/01/27 Meysam Mahooti 10 | % 11 | % Meysam Mahooti (2025). 12 | % ECI2ECEF & ECEF2ECI Transformations 13 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 14 | % MATLAB Central File Exchange. 15 | 16 | function PoleMat = PoleMatrix(xp,yp) 17 | 18 | PoleMat = R_y(-xp) * R_x(-yp); 19 | 20 | end 21 | -------------------------------------------------------------------------------- /+matmap3d/private/R_y.m: -------------------------------------------------------------------------------- 1 | %% R_Y rotation matrix around y-axis 2 | % Input: 3 | % angle angle of rotation [rad] 4 | % 5 | % Output: 6 | % rotmat rotation matrix 7 | % 8 | % Last modified: 2018/01/27 Meysam Mahooti 9 | % 10 | % Meysam Mahooti (2025). 11 | % ECI2ECEF & ECEF2ECI Transformations 12 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 13 | % MATLAB Central File Exchange. 14 | 15 | function rotmat = R_y(angle) 16 | 17 | C = cos(angle); 18 | S = sin(angle); 19 | rotmat = [ C 0 -S; 20 | 0 1 0; 21 | S 0 C ]; 22 | 23 | end 24 | -------------------------------------------------------------------------------- /+matmap3d/private/gast.m: -------------------------------------------------------------------------------- 1 | %% GAST Greenwich Apparent Sidereal Time 2 | % 3 | % Inputs: 4 | % Mjd_UT1 Modified Julian Date UT1 5 | % Mjd_TT Modified Julian Date TT 6 | % 7 | % Output: 8 | % gstime GAST in [rad] 9 | % 10 | % Last modified: 2022/06/16 Meysam Mahooti 11 | % 12 | % Meysam Mahooti (2025). 13 | % ECI2ECEF & ECEF2ECI Transformations 14 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 15 | % MATLAB Central File Exchange. 16 | 17 | function gstime = gast(Mjd_UT1,Mjd_TT) 18 | 19 | gstime = mod(gmst(Mjd_UT1) + EqnEquinox(Mjd_TT), 2*pi); 20 | 21 | end 22 | -------------------------------------------------------------------------------- /.github/workflows/ci.yml: -------------------------------------------------------------------------------- 1 | name: ci 2 | 3 | on: 4 | push: 5 | paths: 6 | - "**.m" 7 | - ".github/workflows/ci.yml" 8 | 9 | # avoid wasted runs 10 | concurrency: 11 | group: ${{ github.workflow }}-${{ github.ref }} 12 | cancel-in-progress: true 13 | 14 | jobs: 15 | 16 | matlab: 17 | runs-on: ubuntu-latest 18 | timeout-minutes: 10 19 | 20 | strategy: 21 | matrix: 22 | release: [R2023b, R2025b] 23 | 24 | steps: 25 | - uses: actions/checkout@v5 26 | 27 | - uses: ./.github/workflows/composite-install-matlab 28 | 29 | - uses: ./.github/workflows/composite-buildtool 30 | -------------------------------------------------------------------------------- /+matmap3d/private/R_x.m: -------------------------------------------------------------------------------- 1 | %% R_X rotation matrix around x-axis 2 | % Input: 3 | % angle angle of rotation [rad] 4 | % 5 | % Output: 6 | % rotmat rotation matrix 7 | % 8 | % Last modified: 2018/01/27 Meysam Mahooti 9 | % 10 | % Meysam Mahooti (2025). 11 | % ECI2ECEF & ECEF2ECI Transformations 12 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 13 | % MATLAB Central File Exchange. 14 | 15 | function [rotmat] = R_x(angle) 16 | 17 | C = cos(angle); 18 | S = sin(angle); 19 | rotmat = [ 1, 0, 0; 20 | 0, C, S; 21 | 0, -S, C ]; 22 | 23 | end 24 | -------------------------------------------------------------------------------- /+matmap3d/private/R_z.m: -------------------------------------------------------------------------------- 1 | % R_Z rotation matrix around z-axis 2 | % Input: 3 | % angle angle of rotation [rad] 4 | % 5 | % Output: 6 | % rotmat rotation matrix 7 | % 8 | % Last modified: 2018/01/27 Meysam Mahooti 9 | % 10 | % Meysam Mahooti (2025). 11 | % ECI2ECEF & ECEF2ECI Transformations 12 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 13 | % MATLAB Central File Exchange. 14 | 15 | function [rotmat] = R_z(angle) 16 | 17 | C = cos(angle); 18 | S = sin(angle); 19 | rotmat = [ C, S, 0.0; 20 | -S, C, 0.0; 21 | 0.0, 0.0, 1.0 ]; 22 | 23 | end 24 | -------------------------------------------------------------------------------- /octave/juliandate.m: -------------------------------------------------------------------------------- 1 | %% juliandate datetime to Julian time 2 | % 3 | % from D.Vallado Fundamentals of Astrodynamics and Applications p.187 4 | % and J. Meeus Astronomical Algorithms 1991 Eqn. 7.1 pg. 61 5 | % 6 | % parameters: 7 | % t: datetime vector 8 | % 9 | % results: 10 | % jd: julian date (days since Jan 1, 4173 BCE 11 | 12 | function jd = juliandate(t) 13 | 14 | [y, mon, d, h, m, s] = datevec(t); 15 | 16 | i = mon < 3; 17 | y(i) = y(i) - 1; 18 | mon(i) = mon(i) + 12; 19 | 20 | A = fix(y / 100.); 21 | B = 2 - A + fix(A / 4.); 22 | C = ((s / 60. + m) / 60. + h) / 24.; 23 | 24 | jd = fix(365.25 * (y + 4716)) + fix(30.6001 * (mon + 1)) + d + B - 1524.5 + C; 25 | 26 | end 27 | -------------------------------------------------------------------------------- /+matmap3d/aer2ned.m: -------------------------------------------------------------------------------- 1 | %% AER2NED convert azimuth, elevation, range to NED coordinates 2 | % 3 | %%% Inputs 4 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 5 | % * az: azimuth clockwise from local north 6 | % * el: elevation angle above local horizon 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% Outputs 10 | % * north, east, down: coordinates of points (meters) 11 | function [north, east, down] = aer2ned(az, el, slantRange, angleUnit) 12 | if nargin < 4 13 | angleUnit = 'd'; 14 | end 15 | 16 | [east, north, up] = matmap3d.aer2enu(az, el, slantRange, angleUnit); 17 | 18 | down = -up; 19 | 20 | end 21 | -------------------------------------------------------------------------------- /+matmap3d/private/GHAMatrix.m: -------------------------------------------------------------------------------- 1 | %% GHAMATRIX Transformation from true equator and equinox to Earth equator and Greenwich meridian system 2 | % 3 | % Inputs: 4 | % Mjd_UT1 Modified Julian Date UT1 5 | % Mjd_TT Modified Julian Date TT 6 | % 7 | % Output: 8 | % GHAmat Greenwich Hour Angle matrix 9 | % 10 | % Last modified: 2022/06/16 Meysam Mahooti 11 | % 12 | % Meysam Mahooti (2025). 13 | % ECI2ECEF & ECEF2ECI Transformations 14 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 15 | % MATLAB Central File Exchange. 16 | 17 | function GHAmat = GHAMatrix(Mjd_UT1,Mjd_TT) 18 | 19 | GHAmat = R_z(gast(Mjd_UT1,Mjd_TT)); 20 | 21 | end 22 | -------------------------------------------------------------------------------- /+matmap3d/private/MeanObliquity.m: -------------------------------------------------------------------------------- 1 | %% MEANOBLIQUITY Computes the mean obliquity of the ecliptic 2 | % 3 | % Input: 4 | % Mjd_TT Modified Julian Date (Terrestrial Time) 5 | % 6 | % Output: 7 | % MOblq Mean obliquity of the ecliptic 8 | % 9 | % Last modified: 2018/01/27 Meysam Mahooti 10 | % 11 | % Meysam Mahooti (2025). 12 | % ECI2ECEF & ECEF2ECI Transformations 13 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 14 | % MATLAB Central File Exchange. 15 | 16 | function MOblq = MeanObliquity(Mjd_TT) 17 | 18 | T = (Mjd_TT - 51544.5)/36525; 19 | 20 | MOblq = (pi/180)*(23.43929111-(46.8150+(0.00059-0.001813*T)*T)*T/3600); 21 | 22 | end 23 | -------------------------------------------------------------------------------- /+matmap3d/eci2aer.m: -------------------------------------------------------------------------------- 1 | %% ECI2AER convert ECI to AER (azimuth, elevation, slant range) 2 | % 3 | % parameters: 4 | % utc: datetime UTC 5 | % x0, y0, z0: ECI x, y, z 6 | % lat, lon, alt: latitude, longitude, altiude of observer (degrees, meters) 7 | % 8 | % outputs: 9 | % az,el,rng: Azimuth (degrees), Elevation (degrees), Slant Range (meters) 10 | function [az, el, rng] = eci2aer(utc, x0, y0, z0, lat, lon, alt) 11 | 12 | az = nan(like=x0); 13 | el = nan(like=x0); 14 | rng = nan(like=x0); 15 | 16 | for i = 1:numel(x0) 17 | r_ecef = matmap3d.eci2ecef(utc, [x0(i); y0(i); z0(i)]); 18 | 19 | [az(i), el(i), rng(i)] = matmap3d.ecef2aer(r_ecef(1), r_ecef(2), r_ecef(3), lat, lon, alt); 20 | end 21 | 22 | end 23 | -------------------------------------------------------------------------------- /+matmap3d/eci2ecef_naive.m: -------------------------------------------------------------------------------- 1 | %% ECI2ECEF_NAIVE rotate ECI coordinates to ECEF 2 | % because this doesn't account for nutation, etc. error is generally significant 3 | % 4 | % x_eci, y_eci, z_eci: eci position vectors 5 | % utc: Matlab datetime UTC 6 | % 7 | % x,y,z: ECEF position (meters) 8 | 9 | function [x,y,z] = eci2ecef_naive(utc, x_eci, y_eci, z_eci) 10 | 11 | % Greenwich hour angles (radians) 12 | gst = matmap3d.greenwichsrt(juliandate(utc)); 13 | 14 | x = nan(like=x_eci); 15 | y = nan(like=y_eci); 16 | z = nan(like=z_eci); 17 | 18 | for j = 1:length(x) 19 | ecef = matmap3d.R3(gst(j)) * [x_eci(j), y_eci(j), z_eci(j)].'; 20 | x(j) = ecef(1); 21 | y(j) = ecef(2); 22 | z(j) = ecef(3); 23 | end 24 | end 25 | -------------------------------------------------------------------------------- /+matmap3d/ecef2eci_naive.m: -------------------------------------------------------------------------------- 1 | %% ECEF2ECI_NAIVE rotate ECEF coordinates to ECI 2 | % because this doesn't account for nutation, etc. error is generally significant 3 | % 4 | %%% Inputs 5 | % x0, y0, z0: ECEF position (meters) 6 | % utc: time UTC 7 | %%% Outputs 8 | % * x,y,z: ECI position (meters) 9 | 10 | function [x,y,z] = ecef2eci_naive(utc, x0, y0, z0) 11 | 12 | % Greenwich hour angles (radians) 13 | gst = matmap3d.greenwichsrt(juliandate(utc)); 14 | 15 | % Convert into ECEF 16 | x = nan(like=gst); 17 | y = nan(like=gst); 18 | z = nan(like=gst); 19 | 20 | for j = 1:length(x) 21 | eci = matmap3d.R3(gst(j)).' * [x0(j), y0(j), z0(j)].'; 22 | x(j) = eci(1); 23 | y(j) = eci(2); 24 | z(j) = eci(3); 25 | end 26 | end 27 | -------------------------------------------------------------------------------- /+matmap3d/greenwichsrt.m: -------------------------------------------------------------------------------- 1 | %% GREENWICHSRT Compute greenwich sidereal time from Julian date 2 | % compute greenwich sidereal time from D. Vallado 4th edition 3 | % 4 | %%% Inputs 5 | % Jdate: Julian days from Jan 1, 4713 BCE from juliantime(utc) or juliandate(utc) 6 | %%% Outputs 7 | % gst: greenwich sidereal time [0, 2pi) 8 | function gst = greenwichsrt(Jdate) 9 | 10 | % Vallado Eq. 3-42 p. 184, Seidelmann 3.311-1 11 | tUT1 = (Jdate - 2451545) / 36525; 12 | % Eqn. 3-47 p. 188 13 | gmst_sec = 67310.54841 + (876600 * 3600 + 8640184.812866) * tUT1 + 0.093104 * tUT1 .^ 2 - 6.2e-6 * tUT1 .^ 3; 14 | % 1/86400 and 2*pi implied by units of radians 15 | tau = 2*pi; 16 | gst = mod(gmst_sec * tau / 86400.0, tau); 17 | 18 | end 19 | -------------------------------------------------------------------------------- /+matmap3d/aer2enu.m: -------------------------------------------------------------------------------- 1 | %% AER2ENU convert azimuth, elevation, range to ENU coordinates 2 | % 3 | %%% Inputs 4 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 5 | % * az: azimuth clockwise from local north 6 | % * el: elevation angle above local horizon 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% Outputs 10 | % * e,n,u: East, North, Up coordinates of test points (meters) 11 | function [e, n, u] = aer2enu (az, el, slantRange, angleUnit) 12 | if nargin < 4 13 | angleUnit = 'd'; 14 | end 15 | 16 | if strncmp(angleUnit, 'd', 1) 17 | az = deg2rad(az); 18 | el = deg2rad(el); 19 | end 20 | 21 | u = slantRange .* sin(el); 22 | r = slantRange .* cos(el); 23 | e = r .* sin(az); 24 | n = r .* cos(az); 25 | 26 | end 27 | -------------------------------------------------------------------------------- /+matmap3d/enu2uvw.m: -------------------------------------------------------------------------------- 1 | %% ENU2UVW convert from ENU to UVW coordinates 2 | % 3 | %%% Inputs 4 | % * e,n,up: East, North, Up coordinates of point(s) (meters) 5 | % * lat0,lon0: geodetic coordinates of observer/reference point (degrees) 6 | % * angleUnit: string for angular units. Default 'd': degrees 7 | % 8 | %%% outputs 9 | % * u,v,w: coordinates of test point(s) (meters) 10 | function [u,v,w] = enu2uvw(east, north, up, lat0, lon0, angleUnit) 11 | if nargin < 6 12 | angleUnit = 'd'; 13 | end 14 | 15 | if strncmp(angleUnit, 'd', 1) 16 | lat0 = deg2rad(lat0); 17 | lon0 = deg2rad(lon0); 18 | end 19 | 20 | t = cos(lat0) * up - sin(lat0) * north; 21 | w = sin(lat0) * up + cos(lat0) * north; 22 | 23 | u = cos(lon0) * t - sin(lon0) * east; 24 | v = sin(lon0) * t + cos(lon0) * east; 25 | 26 | end 27 | -------------------------------------------------------------------------------- /+matmap3d/ecef2ned.m: -------------------------------------------------------------------------------- 1 | %% ECEF2NED Convert ECEF coordinates to NED 2 | % 3 | %%% Inputs 4 | % * x,y,z: Earth Centered Earth Fixed (ECEF) coordinates of test point (meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% outputs 10 | % * North,East,Down: coordinates of points (meters) 11 | function [north, east, down] = ecef2ned(x, y, z, lat0, lon0, alt0, spheroid, angleUnit) 12 | if nargin < 7 || isempty(spheroid) 13 | spheroid = matmap3d.wgs84Ellipsoid(); 14 | end 15 | if nargin < 8 16 | angleUnit = 'd'; 17 | end 18 | 19 | [east, north, up] = matmap3d.ecef2enu(x,y,z,lat0,lon0,alt0,spheroid,angleUnit); 20 | 21 | down = -up; 22 | 23 | end 24 | -------------------------------------------------------------------------------- /+matmap3d/enu2ecefv.m: -------------------------------------------------------------------------------- 1 | %% ENU2ECEFV convert from ENU to ECEF coordinates 2 | % 3 | %%% Inputs 4 | % * e,n,u: East, North, Up coordinates of point(s) (meters) 5 | % * lat0,lon0: geodetic coordinates of observer/reference point (degrees) 6 | % * angleUnit: string for angular units. Default 'd': degrees 7 | % 8 | %%% outputs 9 | % * u,v,w: coordinates of test point(s) (meters) 10 | 11 | function [u, v, w] = enu2ecefv(east, north, up, lat0, lon0, angleUnit) 12 | if nargin < 6 13 | angleUnit = 'd'; 14 | end 15 | 16 | if strncmp(angleUnit, 'd', 1) 17 | lat0 = deg2rad(lat0); 18 | lon0 = deg2rad(lon0); 19 | end 20 | 21 | t = cos(lat0) .* up - sin(lat0) .* north; 22 | w = sin(lat0) .* up + cos(lat0) .* north; 23 | 24 | u = cos(lon0) .* t - sin(lon0) .* east; 25 | v = sin(lon0) .* t + cos(lon0) .* east; 26 | 27 | end 28 | -------------------------------------------------------------------------------- /+matmap3d/private/NutMatrix.m: -------------------------------------------------------------------------------- 1 | %% NUTMATRIX Transformation from mean to true equator and equinox 2 | % 3 | % Input: 4 | % Mjd_TT Modified Julian Date (Terrestrial Time) 5 | % 6 | % Output: 7 | % NutMat Nutation matrix 8 | % 9 | % Last modified: 2018/01/27 Meysam Mahooti 10 | % 11 | % Meysam Mahooti (2025). 12 | % ECI2ECEF & ECEF2ECI Transformations 13 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 14 | % MATLAB Central File Exchange. 15 | 16 | function NutMat = NutMatrix(Mjd_TT) 17 | 18 | % Mean obliquity of the ecliptic 19 | ep = MeanObliquity(Mjd_TT); 20 | 21 | % Nutation in longitude and obliquity 22 | [dpsi, deps] = NutAngles(Mjd_TT); 23 | 24 | % Transformation from mean to true equator and equinox 25 | NutMat = R_x(-ep-deps)*R_z(-dpsi)*R_x(ep); 26 | 27 | end 28 | -------------------------------------------------------------------------------- /+matmap3d/ecef2enuv.m: -------------------------------------------------------------------------------- 1 | %% ECEF2ENUV convert *vector projection* UVW to ENU 2 | % 3 | %%% Inputs 4 | % * u,v,w: meters 5 | % * lat0,lon0: geodetic latitude and longitude (degrees) 6 | % * angleUnit: string for angular system. Default 'd' degrees 7 | % 8 | %%% Outputs 9 | % * e,n,Up: East, North, Up vector 10 | 11 | function [e, n, Up] = ecef2enuv(u, v, w, lat0, lon0, angleUnit) 12 | if nargin < 6 13 | angleUnit = 'd'; 14 | end 15 | 16 | if strncmp(angleUnit, 'd', 1) 17 | lat0 = deg2rad(lat0); 18 | lon0 = deg2rad(lon0); 19 | end 20 | 21 | t = cos(lon0) .* u + sin(lon0) .* v; 22 | e = -sin(lon0) .* u + cos(lon0) .* v; 23 | Up = cos(lat0) .* t + sin(lat0) .* w; 24 | n = -sin(lat0) .* t + cos(lat0) .* w; 25 | 26 | % 1mm precision 27 | if abs(e) < 1e-3, e=0; end 28 | if abs(n) < 1e-3, n=0; end 29 | if abs(Up) < 1e-3, Up=0; end 30 | end 31 | -------------------------------------------------------------------------------- /+matmap3d/private/gmst.m: -------------------------------------------------------------------------------- 1 | %% GMST Greenwich Mean Sidereal Time 2 | % 3 | % Input: 4 | % Mjd_UT1 Modified Julian Date UT1 5 | % 6 | % Output: 7 | % gmstime GMST in [rad] 8 | % 9 | % Last modified: 2018/01/27 Meysam Mahooti 10 | % 11 | % Meysam Mahooti (2025). 12 | % ECI2ECEF & ECEF2ECI Transformations 13 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 14 | % MATLAB Central File Exchange. 15 | 16 | function gmstime = gmst(Mjd_UT1) 17 | 18 | Mjd_0 = floor(Mjd_UT1); 19 | UT1 = 86400*(Mjd_UT1-Mjd_0); % [s] 20 | T_0 = (Mjd_0 - 51544.5)/36525; 21 | T = (Mjd_UT1- 51544.5)/36525; 22 | 23 | gmst = 24110.54841 + 8640184.812866*T_0 + 1.002737909350795*UT1... 24 | + (0.093104-6.2e-6*T)*T*T; % [s] 25 | 26 | gmstime = 2*pi*Frac(gmst/86400); % [rad], 0..2pi 27 | 28 | end 29 | -------------------------------------------------------------------------------- /+matmap3d/ecef2enu.m: -------------------------------------------------------------------------------- 1 | 2 | %% ECEF2ENU convert ECEF to ENU 3 | % 4 | %%% Inputs 5 | % * x,y,z: Earth Centered Earth Fixed (ECEF) coordinates of test point (meters) 6 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 7 | % * spheroid: referenceEllipsoid 8 | % * angleUnit: string for angular units. Default 'd': degrees 9 | % 10 | %%% outputs 11 | % * East, North, Up coordinates of test points (meters) 12 | 13 | function [east, north, up] = ecef2enu(x, y, z, lat0, lon0, alt0, spheroid, angleUnit) 14 | if nargin < 7 || isempty(spheroid) 15 | spheroid = matmap3d.wgs84Ellipsoid(); 16 | end 17 | if nargin < 8 18 | angleUnit = 'd'; 19 | end 20 | 21 | [x0, y0, z0] = matmap3d.geodetic2ecef(spheroid, lat0, lon0, alt0, angleUnit); 22 | [east, north, up] = matmap3d.ecef2enuv(x - x0, y - y0, z - z0, lat0, lon0, angleUnit); 23 | 24 | end 25 | -------------------------------------------------------------------------------- /+matmap3d/aer2eci.m: -------------------------------------------------------------------------------- 1 | %% AER2ECI convert AER (azimuth, elevation, slant range) to ECI 2 | % 3 | %%% Inputs 4 | % * utc: datetime UTC 5 | % * az,el,rng: (degrees, meters) 6 | % * lat, lon, alt: latitude, longitude, altiude of observer (degrees, meters) 7 | % 8 | %%% Outputs 9 | % * x, y, z: ECI x, y, z 10 | function [x, y, z] = aer2eci(utc, az, el, rng, lat, lon, alt, spheroid, angleUnit) 11 | if nargin < 8 || isempty(spheroid) 12 | spheroid = matmap3d.wgs84Ellipsoid(); 13 | end 14 | if nargin < 9 15 | angleUnit = 'd'; 16 | end 17 | 18 | [x1, y1, z1] = matmap3d.aer2ecef(az, el, rng, lat, lon, alt, spheroid, angleUnit); 19 | 20 | x = nan(like=az); 21 | y = nan(like=az); 22 | z = nan(like=az); 23 | 24 | for i = 1:numel(x) 25 | r_eci = matmap3d.ecef2eci(utc, [x1(i); y1(i); z1(i)]); 26 | 27 | x(i) = r_eci(1); 28 | y(i) = r_eci(2); 29 | z(i) = r_eci(3); 30 | end 31 | 32 | end 33 | -------------------------------------------------------------------------------- /+matmap3d/geodetic2ned.m: -------------------------------------------------------------------------------- 1 | %% GEODETIC2NED convert from geodetic to NED coordinates 2 | % 3 | %%% Inputs 4 | % * lat,lon, alt: ellipsoid geodetic coordinates of point under test (degrees, degrees, meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | %% outputs 9 | % * n,e,d: North, East, Down coordinates of test points (meters) 10 | 11 | function [north, east, down] = geodetic2ned(lat, lon, alt, lat0, lon0, alt0, spheroid, angleUnit) 12 | if nargin < 7 || isempty(spheroid) 13 | spheroid = matmap3d.wgs84Ellipsoid(); 14 | end 15 | if nargin < 8 || isempty(angleUnit) 16 | angleUnit = 'd'; 17 | end 18 | 19 | [east, north, up] = matmap3d.geodetic2enu(lat, lon, alt, lat0, lon0, alt0, spheroid, angleUnit); 20 | 21 | down = -up; 22 | 23 | end 24 | -------------------------------------------------------------------------------- /+matmap3d/enu2ecef.m: -------------------------------------------------------------------------------- 1 | %% ENU2ECEF convert from ENU to ECEF coordiantes 2 | % 3 | %%% Inputs 4 | % * east, north, up: coordinates of test points (meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% outputs 10 | % * x,y,z: Earth Centered Earth Fixed (ECEF) coordinates of test point (meters) 11 | function [x, y, z] = enu2ecef(east, north, up, lat0, lon0, alt0, spheroid, angleUnit) 12 | if nargin < 7 || isempty(spheroid) 13 | spheroid = matmap3d.wgs84Ellipsoid(); 14 | end 15 | if nargin < 8 16 | angleUnit = 'd'; 17 | end 18 | 19 | [x0, y0, z0] = matmap3d.geodetic2ecef(spheroid, lat0, lon0, alt0, angleUnit); 20 | [dx, dy, dz] = matmap3d.enu2uvw(east, north, up, lat0, lon0, angleUnit); 21 | 22 | x = x0 + dx; 23 | y = y0 + dy; 24 | z = z0 + dz; 25 | end 26 | -------------------------------------------------------------------------------- /+matmap3d/enu2geodetic.m: -------------------------------------------------------------------------------- 1 | %% ENU2GEODETIC convert from ENU to geodetic coordinates 2 | % 3 | %%% Inputs 4 | % * east,north,up: ENU coordinates of point(s) (meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees, otherwise Radians 8 | % 9 | %%% outputs 10 | % * lat,lon,alt: geodetic coordinates of test points (degrees,degrees,meters) 11 | 12 | function [lat, lon, alt] = enu2geodetic(east, north, up, lat0, lon0, alt0, spheroid, angleUnit) 13 | if nargin < 7 || isempty(spheroid) 14 | spheroid = matmap3d.referenceEllipsoid(); 15 | end 16 | if nargin < 8 17 | angleUnit = 'd'; 18 | end 19 | 20 | [x, y, z] = matmap3d.enu2ecef(east, north, up, lat0, lon0, alt0, spheroid, angleUnit); 21 | [lat, lon, alt] = matmap3d.ecef2geodetic(spheroid, x, y, z, angleUnit); 22 | 23 | end 24 | -------------------------------------------------------------------------------- /+matmap3d/private/EqnEquinox.m: -------------------------------------------------------------------------------- 1 | %% EQNEQUINOX Computation of the equation of the equinoxes 2 | % 3 | % Input: 4 | % Mjd_TT Modified Julian Date (Terrestrial Time) 5 | % 6 | % Output: 7 | % EqE Equation of the equinoxes 8 | % 9 | % Notes: 10 | % The equation of the equinoxes dpsi*cos(eps) is the right ascension of 11 | % the mean equinox referred to the true equator and equinox and is equal 12 | % to the difference between apparent and mean sidereal time. 13 | % 14 | % Last modified: 2018/01/27 Meysam Mahooti 15 | % 16 | % Meysam Mahooti (2025). 17 | % ECI2ECEF & ECEF2ECI Transformations 18 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 19 | % MATLAB Central File Exchange. 20 | 21 | function EqE = EqnEquinox (Mjd_TT) 22 | 23 | % Nutation in longitude and obliquity 24 | dpsi = NutAngles (Mjd_TT); 25 | 26 | % Equation of the equinoxes 27 | EqE = dpsi * cos(MeanObliquity(Mjd_TT)); 28 | 29 | end 30 | -------------------------------------------------------------------------------- /+matmap3d/ecef2aer.m: -------------------------------------------------------------------------------- 1 | %% ECEF2AER convert ECEF of target to azimuth, elevation, slant range from observer 2 | % 3 | %%% Inputs 4 | % * x,y,z: Earth Centered Earth Fixed (ECEF) coordinates of test point (meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% Outputs 10 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 11 | % * az: azimuth clockwise from local north 12 | % * el: elevation angle above local horizon 13 | function [az, el, slantRange] = ecef2aer(x, y, z, lat0, lon0, alt0, spheroid, angleUnit) 14 | if nargin < 7 || isempty(spheroid) 15 | spheroid = matmap3d.wgs84Ellipsoid(); 16 | end 17 | if nargin < 8 18 | angleUnit = 'd'; 19 | end 20 | 21 | [e, n, u] = matmap3d.ecef2enu(x, y, z, lat0, lon0, alt0, spheroid, angleUnit); 22 | [az,el,slantRange] = matmap3d.enu2aer(e, n, u, angleUnit); 23 | 24 | end 25 | -------------------------------------------------------------------------------- /+matmap3d/enu2aer.m: -------------------------------------------------------------------------------- 1 | %% ENU2AER convert ENU to azimuth, elevation, slant range 2 | % 3 | %%% Inputs 4 | % * e,n,u: East, North, Up coordinates of test points (meters) 5 | % * angleUnit: string for angular units. Default 'd': degrees 6 | % 7 | %%% outputs 8 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 9 | % * az: azimuth clockwise from local north 10 | % * el: elevation angle above local horizon 11 | function [az, elev, slantRange] = enu2aer(east, north, up, angleUnit) 12 | if nargin < 4 13 | angleUnit = 'd'; 14 | end 15 | 16 | if abs(east) < 1e-3, east = 0.; end % singularity, 1 mm precision 17 | if abs(north) < 1e-3, north = 0.; end % singularity, 1 mm precision 18 | if abs(up) < 1e-3, up = 0.; end % singularity, 1 mm precision 19 | 20 | r = hypot(east, north); 21 | slantRange = hypot(r,up); 22 | % radians 23 | elev = atan2(up,r); 24 | az = mod(atan2(east, north), 2*pi); 25 | 26 | if strncmp(angleUnit, 'd', 1) 27 | elev = rad2deg(elev); 28 | az = rad2deg(az); 29 | end 30 | 31 | end 32 | -------------------------------------------------------------------------------- /+matmap3d/geodetic2enu.m: -------------------------------------------------------------------------------- 1 | %% GEODETIC2ENU convert from geodetic to ENU coordinates 2 | % 3 | %%% Inputs 4 | % * lat,lon, alt: ellipsoid geodetic coordinates of point under test (degrees, degrees, meters) 5 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 6 | % * spheroid: referenceEllipsoid 7 | % * angleUnit: string for angular units. Default 'd': degrees 8 | % 9 | %%% outputs 10 | % * east,north,up: coordinates of points (meters) 11 | function [east, north, up] = geodetic2enu(lat, lon, alt, lat0, lon0, alt0, spheroid, angleUnit) 12 | if nargin < 7 || isempty(spheroid) 13 | spheroid = matmap3d.wgs84Ellipsoid(); 14 | end 15 | if nargin < 8 || isempty(angleUnit) 16 | angleUnit = 'd'; 17 | end 18 | 19 | [x1,y1,z1] = matmap3d.geodetic2ecef(spheroid, lat,lon,alt,angleUnit); 20 | [x2,y2,z2] = matmap3d.geodetic2ecef(spheroid, lat0,lon0,alt0,angleUnit); 21 | 22 | dx = x1-x2; 23 | dy = y1-y2; 24 | dz = z1-z2; 25 | 26 | [east, north, up] = matmap3d.ecef2enuv(dx, dy, dz, lat0, lon0, angleUnit); 27 | 28 | end 29 | -------------------------------------------------------------------------------- /+matmap3d/aer2geodetic.m: -------------------------------------------------------------------------------- 1 | %% aer2geodetic convert azimuth, elevation, range of target from observer to geodetic coordiantes 2 | % 3 | %%% Inputs 4 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 5 | % * az: azimuth clockwise from local north 6 | % * el: elevation angle above local horizon 7 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 8 | % * spheroid: referenceEllipsoid 9 | % * angleUnit: string for angular units. Default 'd': degrees 10 | % 11 | %%% Outputs 12 | % * lat1,lon1,alt1: geodetic coordinates of test points (degrees,degrees,meters) 13 | function [lat1, lon1, alt1] = aer2geodetic(az, el, slantRange, lat0, lon0, alt0, spheroid, angleUnit) 14 | if nargin < 7 || isempty(spheroid) 15 | spheroid = matmap3d.wgs84Ellipsoid(); 16 | end 17 | if nargin < 8 18 | angleUnit = 'd'; 19 | end 20 | 21 | [x, y, z] = matmap3d.aer2ecef(az, el, slantRange, lat0, lon0, alt0, spheroid, angleUnit); 22 | 23 | [lat1, lon1, alt1] = matmap3d.ecef2geodetic(spheroid, x, y, z, angleUnit); 24 | 25 | end 26 | -------------------------------------------------------------------------------- /+matmap3d/geodetic2ecef.m: -------------------------------------------------------------------------------- 1 | %% GEODETIC2ECEF convert from geodetic to ECEF coordiantes 2 | % 3 | %%% Inputs 4 | % * lat,lon, alt: ellipsoid geodetic coordinates of point(s) (degrees, degrees, meters) 5 | % * spheroid: referenceEllipsoid 6 | % * angleUnit: string for angular units. Default 'd': degrees 7 | % 8 | %%% outputs 9 | % * x,y,z: ECEF coordinates of test point(s) (meters) 10 | function [x,y,z] = geodetic2ecef(spheroid, lat, lon, alt, angleUnit) 11 | if nargin < 5 12 | angleUnit = 'd'; 13 | end 14 | 15 | %% defaults 16 | if isempty(spheroid) 17 | spheroid = matmap3d.wgs84Ellipsoid(); 18 | end 19 | 20 | if strncmp(angleUnit, 'd', 1) 21 | lat = deg2rad(lat); 22 | lon = deg2rad(lon); 23 | end 24 | 25 | % Radius of curvature of the prime vertical section 26 | N = matmap3d.get_radius_normal(lat, spheroid); 27 | % Compute cartesian (geocentric) coordinates given (curvilinear) geodetic coordinates. 28 | 29 | x = (N + alt) .* cos(lat) .* cos(lon); 30 | y = (N + alt) .* cos(lat) .* sin(lon); 31 | z = (N .* (spheroid.SemiminorAxis / spheroid.SemimajorAxis)^2 + alt) .* sin(lat); 32 | 33 | end 34 | -------------------------------------------------------------------------------- /+matmap3d/private/timediff.m: -------------------------------------------------------------------------------- 1 | %% TIMEDIFF Time differences [s] 2 | % 3 | % Last modified: 2018/01/27 M. Mahooti 4 | % 5 | % Meysam Mahooti (2025). 6 | % ECI2ECEF & ECEF2ECI Transformations 7 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 8 | % MATLAB Central File Exchange. 9 | 10 | function [UT1_TAI,UTC_GPS,UT1_GPS,TT_UTC,GPS_UTC] = timediff(UT1_UTC,TAI_UTC) 11 | 12 | TT_TAI = +32.184; % TT-TAI time difference [s] 13 | 14 | GPS_TAI = -19.0; % GPS-TAI time difference [s] 15 | 16 | % TT_GPS = TT_TAI-GPS_TAI; % TT-GPS time difference [s] 17 | 18 | % TAI_GPS = -GPS_TAI; % TAI-GPS time difference [s] 19 | 20 | UT1_TAI = UT1_UTC-TAI_UTC; % UT1-TAI time difference [s] 21 | 22 | UTC_TAI = -TAI_UTC; % UTC-TAI time difference [s] 23 | 24 | UTC_GPS = UTC_TAI-GPS_TAI; % UTC_GPS time difference [s] 25 | 26 | UT1_GPS = UT1_TAI-GPS_TAI; % UT1-GPS time difference [s] 27 | 28 | TT_UTC = TT_TAI-UTC_TAI; % TT-UTC time difference [s] 29 | 30 | GPS_UTC = GPS_TAI-UTC_TAI; % GPS-UTC time difference [s] 31 | 32 | end 33 | -------------------------------------------------------------------------------- /codemeta.json: -------------------------------------------------------------------------------- 1 | { 2 | "@context": "https://doi.org/10.5063/schema/codemeta-2.0", 3 | "@type": "SoftwareSourceCode", 4 | "license": "https://spdx.org/licenses/BSD-3-Clause", 5 | "codeRepository": "https://github.com/geospace-code/matmap3d", 6 | "contIntegration": "https://github.com/geospace-code/matmap3d/actions", 7 | "downloadUrl": "https://github.com/geospace-code/matmap3d/releases", 8 | "issueTracker": "https://github.com/geospace-code/matmap3d/issues", 9 | "name": "matmap3d", 10 | "identifier": "10.5281/zenodo.3966173", 11 | "description": "Matlab / GNU Octave coordinate conversions for geospace ecef enu eci. Similar to Python PyMap3D.", 12 | "applicationCategory": "compuation", 13 | "developmentStatus": "active", 14 | "keywords": [ 15 | "geodesy" 16 | ], 17 | "programmingLanguage": [ 18 | "Matlab" 19 | ], 20 | "author": [ 21 | { 22 | "@type": "Person", 23 | "@id": "https://orcid.org/0000-0002-1637-6526", 24 | "givenName": "Michael", 25 | "familyName": "Hirsch" 26 | } 27 | ] 28 | } 29 | -------------------------------------------------------------------------------- /+matmap3d/private/PrecMatrix.m: -------------------------------------------------------------------------------- 1 | %% PRECMATRIX Precession transformation of equatorial coordinates 2 | % 3 | % Inputs: 4 | % Mjd_1 Epoch given (Modified Julian Date TT) 5 | % MjD_2 Epoch to precess to (Modified Julian Date TT) 6 | % 7 | % Output: 8 | % PrecMat Precession transformation matrix 9 | % 10 | % Last modified: 2018/01/27 Meysam Mahooti 11 | % 12 | % Meysam Mahooti (2025). 13 | % ECI2ECEF & ECEF2ECI Transformations 14 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 15 | % MATLAB Central File Exchange. 16 | 17 | function PrecMat = PrecMatrix(Mjd_1, Mjd_2) 18 | 19 | T = (Mjd_1-51544.5)/36525; 20 | dT = (Mjd_2-Mjd_1)/36525; 21 | 22 | Arcs = 3600*180/pi; 23 | 24 | % Precession angles 25 | zeta = ( (2306.2181+(1.39656-0.000139*T)*T)+ ... 26 | ((0.30188-0.000344*T)+0.017998*dT)*dT )*dT/Arcs; 27 | z = zeta + ( (0.79280+0.000411*T)+0.000205*dT)*dT*dT/Arcs; 28 | theta = ( (2004.3109-(0.85330+0.000217*T)*T)- ... 29 | ((0.42665+0.000217*T)+0.041833*dT)*dT )*dT/Arcs; 30 | 31 | % Precession matrix 32 | PrecMat = R_z(-z) * R_y(theta) * R_z(-zeta); 33 | 34 | end 35 | -------------------------------------------------------------------------------- /+matmap3d/geodetic2aer.m: -------------------------------------------------------------------------------- 1 | 2 | %% GEODETIC2AER converts geodetic coordinates to azimuth, elevation, slant range 3 | % from an observer's perspective, convert target coordinates to azimuth, elevation, slant range. 4 | % 5 | %%% Inputs 6 | % * lat,lon, alt: ellipsoid geodetic coordinates of point under test (degrees, degrees, meters) 7 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 8 | % * spheroid: referenceEllipsoid 9 | % * angleUnit: string for angular units. Default 'd': degrees, otherwise Radians 10 | % 11 | %%% Outputs 12 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 13 | % * az: azimuth clockwise from local north 14 | % * el: elevation angle above local horizon 15 | function [az, el, slantRange] = geodetic2aer(lat, lon, alt, lat0, lon0, alt0, spheroid, angleUnit) 16 | if nargin < 7 || isempty(spheroid) 17 | spheroid = matmap3d.wgs84Ellipsoid(); 18 | end 19 | if nargin < 8 20 | angleUnit = 'd'; 21 | end 22 | 23 | [e, n, u] = matmap3d.geodetic2enu(lat, lon, alt, lat0, lon0, alt0, spheroid, angleUnit); 24 | [az, el, slantRange] = matmap3d.enu2aer(e, n, u, angleUnit); 25 | 26 | end 27 | -------------------------------------------------------------------------------- /+matmap3d/aer2ecef.m: -------------------------------------------------------------------------------- 1 | %% AER2ECEF convert azimuth, elevation, range to target from observer to ECEF coordinates 2 | % 3 | %%% Inputs 4 | % 5 | % * az, el, slantrange: look angles and distance to point under test (degrees, degrees, meters) 6 | % * az: azimuth clockwise from local north 7 | % * el: elevation angle above local horizon 8 | % * lat0, lon0, alt0: ellipsoid geodetic coordinates of observer/reference (degrees, degrees, meters) 9 | % * spheroid: referenceEllipsoid 10 | % * angleUnit: string for angular units. Default 'd': degrees 11 | % 12 | %%% outputs 13 | % 14 | % * x,y,z: Earth Centered Earth Fixed (ECEF) coordinates of test point (meters) 15 | 16 | function [x,y,z] = aer2ecef(az, el, slantRange, lat0, lon0, alt0, spheroid, angleUnit) 17 | if nargin < 7 || isempty(spheroid) 18 | spheroid = matmap3d.wgs84Ellipsoid(); 19 | end 20 | if nargin < 8 21 | angleUnit = 'd'; 22 | end 23 | 24 | % Origin of the local system in geocentric coordinates. 25 | [x0, y0, z0] = matmap3d.geodetic2ecef(spheroid, lat0, lon0, alt0, angleUnit); 26 | % Convert Local Spherical AER to ENU 27 | [e, n, u] = matmap3d.aer2enu(az, el, slantRange, angleUnit); 28 | % Rotating ENU to ECEF 29 | [dx, dy, dz] = matmap3d.enu2uvw(e, n, u, lat0, lon0, angleUnit); 30 | % Origin + offset from origin equals position in ECEF 31 | x = x0 + dx; 32 | y = y0 + dy; 33 | z = z0 + dz; 34 | 35 | end 36 | -------------------------------------------------------------------------------- /+matmap3d/private/Mjday.m: -------------------------------------------------------------------------------- 1 | %% MJDAY Modified Julian Date from calendar date and time 2 | % 3 | % Inputs: 4 | % Year Calendar date components 5 | % Month 6 | % Day 7 | % Hour Time components 8 | % Min 9 | % Sec 10 | % Output: 11 | % Modified Julian Date 12 | % 13 | % Last modified: 2022/09/24 Meysam Mahooti 14 | % 15 | % Reference: 16 | % Montenbruck O., Gill E.; Satellite Orbits: Models, Methods and 17 | % Applications; Springer Verlag, Heidelberg; Corrected 3rd Printing (2005). 18 | % 19 | % Meysam Mahooti (2025). 20 | % ECI2ECEF & ECEF2ECI Transformations 21 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 22 | % MATLAB Central File Exchange. 23 | 24 | function Mjd = Mjday(utc) 25 | arguments 26 | utc (1,1) datetime 27 | end 28 | 29 | Month = utc.Month; 30 | Year = utc.Year; 31 | 32 | if (Month<=2) 33 | Month = Month+12; 34 | Year = Year-1; 35 | end 36 | 37 | if ( (10000*Year + 100*Month + utc.Day) <= 15821004 ) 38 | b = (-2 + fix((Year+4716)/4) - 1179); % Julian calendar 39 | else 40 | b = fix(Year/400)-fix(Year/100)+fix(Year/4); % Gregorian calendar 41 | end 42 | 43 | MjdMidnight = 365*Year - 679004 + b + fix(30.6001*(Month+1)) + utc.Day; 44 | FracOfDay = (utc.Hour + utc.Minute/60+ utc.Second/3600)/24; 45 | 46 | Mjd = MjdMidnight + FracOfDay; 47 | 48 | 49 | end 50 | -------------------------------------------------------------------------------- /LICENSE: -------------------------------------------------------------------------------- 1 | BSD 2-Clause License 2 | 3 | Copyright (c) 2014-2022, SciVision, Inc. 4 | Copyright (c) 2013, Felipe Geremia Nievinski 5 | 6 | All rights reserved. 7 | 8 | Redistribution and use in source and binary forms, with or without 9 | modification, are permitted provided that the following conditions are met: 10 | 11 | * Redistributions of source code must retain the above copyright notice, this 12 | list of conditions and the following disclaimer. 13 | 14 | * Redistributions in binary form must reproduce the above copyright notice, 15 | this list of conditions and the following disclaimer in the documentation 16 | and/or other materials provided with the distribution. 17 | 18 | THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" 19 | AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 20 | IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE 21 | DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE 22 | FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL 23 | DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR 24 | SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER 25 | CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, 26 | OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE 27 | OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 28 | -------------------------------------------------------------------------------- /+matmap3d/private/load_eop.m: -------------------------------------------------------------------------------- 1 | %% LOAD_EOP: Reads Earth Orientation Parameters (EOP) from a file. 2 | 3 | function eopdata = load_eop(filename) 4 | arguments 5 | filename {mustBeFile} = fullfile(fileparts(mfilename('fullpath')), 'EOP-All.txt') 6 | end 7 | 8 | fid = fopen(filename, 'r'); 9 | % ---------------------------------------------------------------------------------------------------- 10 | % | Date MJD x y UT1-UTC LOD dPsi dEpsilon dX dY DAT 11 | % |(0h UTC) " " s s " " " " s 12 | % ---------------------------------------------------------------------------------------------------- 13 | while ~feof(fid) 14 | tline = fgetl(fid); 15 | k = strfind(tline,'NUM_OBSERVED_POINTS'); 16 | if (k ~= 0) 17 | numrecsobs = str2num(tline(21:end)); %#ok 18 | 19 | fgetl(fid); 20 | for i=1:numrecsobs 21 | eopdata(:,i) = fscanf(fid,'%i %d %d %i %f %f %f %f %f %f %f %f %i',[13 1]); %#ok 22 | end 23 | for i=1:4 24 | tline = fgetl(fid); 25 | end 26 | numrecspred = str2num(tline(22:end)); %#ok 27 | 28 | fgetl(fid); 29 | 30 | for i=numrecsobs+1:numrecsobs+numrecspred 31 | eopdata(:,i) = fscanf(fid,'%i %d %d %i %f %f %f %f %f %f %f %f %i',[13 1]); 32 | end 33 | break 34 | end 35 | end 36 | 37 | fclose(fid); 38 | 39 | end 40 | -------------------------------------------------------------------------------- /+matmap3d/ecef2geodetic.m: -------------------------------------------------------------------------------- 1 | %% ECEF2GEODETIC convert ECEF to geodetic coordinates 2 | % 3 | %%% Inputs 4 | % * x,y,z: ECEF coordinates of test point(s) (meters) 5 | % * spheroid: referenceEllipsoid 6 | % * angleUnit: string for angular units. Default 'd': degrees 7 | % 8 | %%% Outputs 9 | % * lat,lon, alt: ellipsoid geodetic coordinates of point(s) (degrees, degrees, meters) 10 | % 11 | % based on: 12 | % You, Rey-Jer. (2000). Transformation of Cartesian to Geodetic Coordinates without Iterations. 13 | % Journal of Surveying Engineering. doi: 10.1061/(ASCE)0733-9453 14 | 15 | function [lat,lon,alt] = ecef2geodetic(spheroid, x, y, z, angleUnit) 16 | if nargin < 5 17 | angleUnit = 'd'; 18 | end 19 | 20 | if isempty(spheroid) 21 | spheroid = matmap3d.wgs84Ellipsoid(); 22 | end 23 | 24 | a = spheroid.SemimajorAxis; 25 | b = spheroid.SemiminorAxis; 26 | 27 | r = sqrt(x.^2 + y.^2 + z.^2); 28 | 29 | E = sqrt(a.^2 - b.^2); 30 | 31 | % eqn. 4a 32 | u = sqrt(0.5 * (r.^2 - E.^2) + 0.5 * sqrt((r.^2 - E.^2).^2 + 4 * E.^2 .* z.^2)); 33 | 34 | Q = hypot(x, y); 35 | 36 | huE = hypot(u, E); 37 | 38 | % eqn. 4b 39 | Beta = atan(huE ./ u .* z ./ hypot(x, y)); 40 | 41 | % eqn. 13 42 | eps = ((b * u - a * huE + E.^2) .* sin(Beta)) ./ (a * huE ./ cos(Beta) - E.^2 .* cos(Beta)); 43 | 44 | Beta = Beta + eps; 45 | % final output 46 | lat = atan(a / b * tan(Beta)); 47 | 48 | lon = atan2(y, x); 49 | 50 | % eqn. 7 51 | alt = hypot(z - b * sin(Beta), Q - a * cos(Beta)); 52 | 53 | % inside ellipsoid? 54 | inside = (x.^2 ./ a.^2) + (y.^2 ./ a.^2) + (z.^2 ./ b.^2) < 1; 55 | alt(inside) = -alt(inside); 56 | 57 | 58 | if strncmp(angleUnit, 'd', 1) 59 | lat = rad2deg(lat); 60 | lon = rad2deg(lon); 61 | end 62 | 63 | end 64 | -------------------------------------------------------------------------------- /.github/workflows/publish.yml: -------------------------------------------------------------------------------- 1 | # Simple workflow for deploying static content to GitHub Pages 2 | name: Deploy static content to Pages 3 | 4 | on: 5 | # Runs on pushes targeting the default branch 6 | push: 7 | branches: ["main"] 8 | 9 | # Allows you to run this workflow manually from the Actions tab 10 | workflow_dispatch: 11 | 12 | # Sets permissions of the GITHUB_TOKEN to allow deployment to GitHub Pages 13 | permissions: 14 | contents: read 15 | pages: write 16 | id-token: write 17 | 18 | # Allow only one concurrent deployment, skipping runs queued between the run in-progress and latest queued. 19 | # However, do NOT cancel in-progress runs as we want to allow these production deployments to complete. 20 | concurrency: 21 | group: "pages" 22 | cancel-in-progress: false 23 | 24 | 25 | jobs: 26 | 27 | deploy: 28 | 29 | strategy: 30 | matrix: 31 | release: ["latest"] 32 | 33 | environment: 34 | name: github-pages 35 | url: ${{ steps.deployment.outputs.page_url }} 36 | 37 | runs-on: ubuntu-latest 38 | 39 | steps: 40 | - uses: actions/checkout@v5 41 | 42 | - name: Setup Pages 43 | uses: actions/configure-pages@v5 44 | 45 | - uses: ./.github/workflows/composite-install-matlab 46 | 47 | - name: Run Matlab buildtool 48 | uses: matlab-actions/run-build@v2 49 | with: 50 | tasks: publish 51 | # don't run "test" task as then even "clean" will leave Matlab in cached state where 52 | # it won't publish due to thinking a function is a Mex file. Error is 53 | # "Error using publish Only MATLAB code can be published" 54 | 55 | - name: Upload artifact 56 | uses: actions/upload-pages-artifact@v4 57 | with: 58 | path: 'docs/' 59 | 60 | - name: Deploy to GitHub Pages 61 | id: deployment 62 | uses: actions/deploy-pages@v4 63 | -------------------------------------------------------------------------------- /+matmap3d/eci2ecef.m: -------------------------------------------------------------------------------- 1 | %% ECI2ECEF Transforms Earth Centered Inertial (ECI) coordinates to Earth 2 | % Centered Earth Fixed (ECEF) coordinates 3 | % 4 | % original by: Meysam Mahooti 5 | % 6 | 7 | function [r_ecef, v_ecef] = eci2ecef(utc, r_eci, v_eci, eopdata) 8 | arguments 9 | utc (1,1) datetime 10 | r_eci (3,1) {mustBeReal} 11 | v_eci {mustBeReal} = [] 12 | eopdata (13, :) {mustBeReal} = h5read(fullfile(fileparts(mfilename('fullpath')), "private/EOP-All.h5"), '/eop') 13 | end 14 | 15 | mjd_utc = Mjday(utc); 16 | 17 | [x_pole,y_pole,UT1_UTC,LOD,~,~,~,~,TAI_UTC] = IERS(eopdata, mjd_utc,'l'); 18 | [~,~,~,TT_UTC] = timediff(UT1_UTC,TAI_UTC); 19 | 20 | MJD_UT1 = mjd_utc + UT1_UTC/86400; 21 | MJD_TT = mjd_utc + TT_UTC/86400; 22 | 23 | % ICRS to ITRS transformation matrix and its derivative 24 | P = PrecMatrix(51544.5, MJD_TT); % IAU 1976 Precession 25 | N = NutMatrix(MJD_TT); % IAU 1980 Nutation 26 | Theta = GHAMatrix(MJD_UT1,MJD_TT); % Earth rotation 27 | Po = PoleMatrix(x_pole,y_pole); % Polar motion 28 | 29 | S = zeros(3); 30 | S(1,2) = 1; 31 | S(2,1) = -1; 32 | Omega = 15.04106717866910/3600*pi/180 - 0.843994809*1e-9*LOD; % [rad/s] IERS 33 | dTheta = Omega*S*Theta; % Derivative of Earth rotation matrix [1/s] 34 | U = Po*Theta*N*P; % ICRS to ITRS transformation 35 | dU = Po*dTheta*N*P; % Derivative [1/s] % Derivative [1/s] 36 | 37 | % Transformation from ICRS to ITRS 38 | r_ecef = U * r_eci; 39 | if nargout > 1 40 | mustBeEqualSize(r_eci, v_eci) 41 | v_ecef = U * v_eci + dU * r_eci; 42 | end 43 | 44 | end 45 | 46 | %!test 47 | %! pkg load tablicious 48 | %! pkg load hdf5oct 49 | %! eopdata = h5read("private/EOP-All.h5", '/eop'); 50 | %! utc = datetime(2019, 1, 4, 12,0,0); 51 | %! eci = [-2981784; 5207055; 3161595]; 52 | %! r_ecef = eci2ecef(utc, eci, [], eopdata); 53 | %! assert(abs(r_ecef - [-5762654.65677142; -1682688.09235503; 3156027.98313692]) < 2e-5 * max(abs(r_ecef))) 54 | -------------------------------------------------------------------------------- /+matmap3d/@referenceEllipsoid/referenceEllipsoid.m: -------------------------------------------------------------------------------- 1 | %% referenceEllipsoid Select available ellipsoid 2 | % 3 | %%% inputs 4 | % * name: string of model name. Default: 'wgs84' 5 | %%% outputs 6 | % * E: referenceEllipsoid 7 | 8 | classdef referenceEllipsoid 9 | properties 10 | Code 11 | Name 12 | LengthUnit 13 | SemimajorAxis 14 | Flattening 15 | SemiminorAxis 16 | Eccentricity 17 | MeanRadius 18 | Volume 19 | end 20 | 21 | methods 22 | function e = referenceEllipsoid(name, lengthUnit) 23 | arguments 24 | name {mustBeTextScalar} = 'wgs84' 25 | lengthUnit {mustBeTextScalar} = 'm' 26 | end 27 | 28 | assert(any(strcmp(lengthUnit, {'m', 'meter', 'meters'}))) 29 | 30 | switch name 31 | case 'wgs84' 32 | % WGS-84 ellipsoid parameters. 33 | % 34 | % 35 | e.Code = 7030; 36 | e.Name = 'World Geodetic System 1984'; 37 | e.LengthUnit = 'meter'; 38 | e.SemimajorAxis = 6378137.0; 39 | e.Flattening = 1/298.2572235630; 40 | e.SemiminorAxis = e.SemimajorAxis * (1 - e.Flattening); 41 | e.Eccentricity = get_eccentricity(e); 42 | e.MeanRadius = meanradius(e); 43 | e.Volume = spheroidvolume(e); 44 | case 'grs80' 45 | % GRS-80 ellipsoid parameters 46 | % (accessed 2018-01-22) 47 | e.Code = 7019; 48 | e.Name = 'Geodetic Reference System 1980'; 49 | e.LengthUnit = 'meter'; 50 | e.SemimajorAxis = 6378137.0; 51 | e.Flattening = 1/298.257222100882711243; 52 | e.SemiminorAxis = e.SemimajorAxis * (1 - e.Flattening); 53 | e.Eccentricity = get_eccentricity(e); 54 | e.MeanRadius = meanradius(e); 55 | e.Volume = spheroidvolume(e); 56 | otherwise, error('%s not yet implemented', name) 57 | end 58 | end 59 | 60 | function v = spheroidvolume(E) 61 | v = 4*pi/3 * E.SemimajorAxis^2 * E.SemiminorAxis; 62 | 63 | assert(v>=0) 64 | end 65 | 66 | function r = meanradius(E) 67 | r = (2*E.SemimajorAxis + E.SemiminorAxis) / 3; 68 | 69 | assert(r>=0) 70 | end 71 | 72 | function ecc = get_eccentricity(E) 73 | ecc = sqrt ( (E.SemimajorAxis^2 - E.SemiminorAxis^2) / (E.SemimajorAxis^2)); 74 | end 75 | 76 | end 77 | end 78 | -------------------------------------------------------------------------------- /+matmap3d/ecef2eci.m: -------------------------------------------------------------------------------- 1 | %% ECEF2ECI Transforms Earth Centered Earth Fixed (ECEF) coordinates to 2 | % Earth Centered Inertial (ECI) coordinates 3 | % 4 | %%% Inputs 5 | % * utc: datetime UTC 6 | % * r_ecef: 3x1 vector of ECEF position [meters] 7 | % * v_ecef: 3x1 vector if ECEF velocity [m/s] (optional) 8 | %%% Outputs 9 | % * r_eci: 3x1 vector of ECI position [meters] 10 | % * v_eci; 3x1 vector of ECI velocity [m/s] (optional) 11 | % 12 | % Original by: Meysam Mahooti 13 | 14 | function [r_eci, v_eci] = ecef2eci(utc, r_ecef, v_ecef, eopdata) 15 | arguments 16 | utc (1,1) datetime 17 | r_ecef (3,1) {mustBeReal} 18 | v_ecef {mustBeReal} = [] 19 | eopdata (13, :) {mustBeReal} = h5read(fullfile(fileparts(mfilename('fullpath')), "private/EOP-All.h5"), '/eop') 20 | end 21 | 22 | mjd_utc = Mjday(utc); 23 | 24 | [x_pole,y_pole,UT1_UTC,LOD,~,~,~,~,TAI_UTC] = IERS(eopdata, mjd_utc, 'l'); 25 | [~,~,~,TT_UTC] = timediff(UT1_UTC,TAI_UTC); 26 | 27 | MJD_UT1 = mjd_utc + UT1_UTC/86400; 28 | MJD_TT = mjd_utc + TT_UTC/86400; 29 | 30 | % ICRS to ITRS transformation matrix and its derivative 31 | P = PrecMatrix(51544.5, MJD_TT); % IAU 1976 Precession 32 | N = NutMatrix(MJD_TT); % IAU 1980 Nutation 33 | Theta = GHAMatrix(MJD_UT1,MJD_TT); % Earth rotation 34 | Po = PoleMatrix(x_pole,y_pole); % Polar motion 35 | U = Po*Theta*N*P; % ICRS to ITRS transformation 36 | 37 | % Transformation from WGS to ICRS 38 | r_eci = U'*r_ecef; 39 | if nargout > 1 40 | mustBeEqualSize(r_ecef, v_ecef) 41 | S = zeros(3); 42 | S(1,2) = 1; 43 | S(2,1) = -1; 44 | Omega = 15.04106717866910/3600*pi/180 - 0.843994809*1e-9*LOD; % [rad/s] IERS 45 | dTheta = Omega*S*Theta; % Derivative of Earth rotation matrix [1/s] 46 | dU = Po*dTheta*N*P; % Derivative [1/s] 47 | v_eci = U'*v_ecef + dU'*r_ecef; 48 | end 49 | 50 | end 51 | 52 | %!test 53 | %! pkg load tablicious 54 | %! pkg load hdf5oct 55 | %! eopdata = h5read("private/EOP-All.h5", '/eop'); 56 | %! utc = datetime(2019, 1, 4, 12,0,0); 57 | %! ecef = [-5762640; -1682738; 3156028]; 58 | %! r_eci = ecef2eci(utc, ecef, [], eopdata); 59 | %! assert(abs(r_eci - [-2981829.07728415; 5207029.04470791; 3161595.0981204]) < 1e-5 * max(abs(r_eci))) 60 | -------------------------------------------------------------------------------- /+matmap3d/private/IERS.m: -------------------------------------------------------------------------------- 1 | %% IERS Management of IERS time and polar motion data 2 | % 3 | % Last modified: 2018/02/01 Meysam Mahooti 4 | % 5 | % Meysam Mahooti (2025). 6 | % ECI2ECEF & ECEF2ECI Transformations 7 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 8 | % MATLAB Central File Exchange. 9 | 10 | function [x_pole,y_pole,UT1_UTC,LOD,dpsi,deps,dx_pole,dy_pole,TAI_UTC] = IERS(eop,Mjd_UTC,interp) 11 | arguments 12 | eop (13,:) {mustBeReal} 13 | Mjd_UTC (1,1) {mustBeReal, mustBePositive, mustBeFinite} 14 | interp {mustBeTextScalar} = 'n' 15 | end 16 | 17 | Arcs = 3600*180/pi; 18 | 19 | mjd = floor(Mjd_UTC); 20 | i = find(mjd==eop(4,:),1,'first'); 21 | 22 | switch interp 23 | case 'l' 24 | % linear interpolation 25 | preeop = eop(:,i); 26 | nexteop = eop(:,i+1); 27 | fixf = Mjd_UTC-floor(Mjd_UTC); 28 | % Setting of IERS Earth rotation parameters 29 | % (UT1-UTC [s], TAI-UTC [s], x ["], y ["]) 30 | x_pole = preeop(5)+(nexteop(5)-preeop(5))*fixf; 31 | y_pole = preeop(6)+(nexteop(6)-preeop(6))*fixf; 32 | UT1_UTC = preeop(7)+(nexteop(7)-preeop(7))*fixf; 33 | LOD = preeop(8)+(nexteop(8)-preeop(8))*fixf; 34 | dpsi = preeop(9)+(nexteop(9)-preeop(9))*fixf; 35 | deps = preeop(10)+(nexteop(10)-preeop(10))*fixf; 36 | dx_pole = preeop(11)+(nexteop(11)-preeop(11))*fixf; 37 | dy_pole = preeop(12)+(nexteop(12)-preeop(12))*fixf; 38 | TAI_UTC = preeop(13); 39 | 40 | x_pole = x_pole/Arcs; % Pole coordinate [rad] 41 | y_pole = y_pole/Arcs; % Pole coordinate [rad] 42 | dpsi = dpsi/Arcs; 43 | deps = deps/Arcs; 44 | dx_pole = dx_pole/Arcs; % Pole coordinate [rad] 45 | dy_pole = dy_pole/Arcs; % Pole coordinate [rad] 46 | case 'n' 47 | eop = eop(:,i); 48 | % Setting of IERS Earth rotation parameters 49 | % (UT1-UTC [s], TAI-UTC [s], x ["], y ["]) 50 | x_pole = eop(5)/Arcs; % Pole coordinate [rad] 51 | y_pole = eop(6)/Arcs; % Pole coordinate [rad] 52 | UT1_UTC = eop(7); % UT1-UTC time difference [s] 53 | LOD = eop(8); % Length of day [s] 54 | dpsi = eop(9)/Arcs; 55 | deps = eop(10)/Arcs; 56 | dx_pole = eop(11)/Arcs; % Pole coordinate [rad] 57 | dy_pole = eop(12)/Arcs; % Pole coordinate [rad] 58 | TAI_UTC = eop(13); % TAI-UTC time difference [s] 59 | end 60 | 61 | end 62 | -------------------------------------------------------------------------------- /+matmap3d/lookAtSpheroid.m: -------------------------------------------------------------------------------- 1 | %% LOOKATSPHEROID 2 | % Calculates line-of-sight intersection with Earth (or other ellipsoid) surface from above surface ./ orbit 3 | % 4 | %%% Inputs 5 | % * lat0, lon0: latitude and longitude of starting point 6 | % * h0: altitude of starting point in meters 7 | % * az: azimuth angle of line-of-sight, clockwise from North 8 | % * tilt: tilt angle of line-of-sight with respect to local vertical (nadir = 0) 9 | % 10 | %%% Outputs 11 | % * lat, lon: latitude and longitude where the line-of-sight intersects with the Earth ellipsoid 12 | % * d: slant range in meters from the starting point to the intersect point 13 | % 14 | % Values will be NaN if the line of sight does not intersect. 15 | % 16 | % Algorithm based on: 17 | % https://medium.com/@stephenhartzell/satellite-line-of-sight-intersection-with-earth-d786b4a6a9b6 18 | % Stephen Hartzell 19 | function [lat, lon, d] = lookAtSpheroid(lat0, lon0, h0, az, tilt, spheroid, angleUnit) 20 | if nargin < 6 || isempty(spheroid) 21 | spheroid = matmap3d.wgs84Ellipsoid(); 22 | end 23 | if nargin < 7 || isempty(angleUnit) 24 | angleUnit = 'd'; 25 | end 26 | 27 | if strncmp(angleUnit, 'd', 1) 28 | el = tilt - 90; 29 | else 30 | el = tilt - pi/2; 31 | end 32 | 33 | a = spheroid.SemimajorAxis; 34 | b = a; 35 | c = spheroid.SemiminorAxis; 36 | 37 | [e, n, u] = matmap3d.aer2enu(az, el, 1., angleUnit); % fixed 1 km slant range 38 | [u, v, w] = matmap3d.enu2uvw(e, n, u, lat0, lon0, angleUnit); 39 | [x, y, z] = matmap3d.geodetic2ecef([], lat0, lon0, h0, angleUnit); 40 | 41 | value = -a.^2 .* b.^2 .* w .* z - a.^2 .* c.^2 .* v .* y - b.^2 .* c.^2 .* u .* x; 42 | radical = a.^2 .* b.^2 .* w.^2 + a.^2 .* c.^2 .* v.^2 - a.^2 .* v.^2 .* z.^2 + 2 .* a.^2 .* v .* w .* y .* z - ... 43 | a.^2 .* w.^2 .* y.^2 + b.^2 .* c.^2 .* u.^2 - b.^2 .* u.^2 .* z.^2 + 2 .* b.^2 .* u .* w .* x .* z - ... 44 | b.^2 .* w.^2 .* x.^2 - c.^2 .* u.^2 .* y.^2 + 2 .* c.^2 .* u .* v .* x .* y - c.^2 .* v.^2 .* x.^2; 45 | 46 | magnitude = a.^2 .* b.^2 .* w.^2 + a.^2 .* c.^2 .* v.^2 + b.^2 .* c.^2 .* u.^2; 47 | 48 | % Return nan if radical < 0 or d < 0 because LOS vector does not point towards Earth 49 | d = (value - a .* b .* c .* sqrt(radical)) ./ magnitude; 50 | d(radical < 0 | d < 0) = nan; % separate line 51 | 52 | % altitude should be zero 53 | [lat, lon] = matmap3d.ecef2geodetic([], x + d .* u, y + d .* v, z + d .* w, angleUnit); 54 | 55 | end 56 | -------------------------------------------------------------------------------- /+matmap3d/private/SAT_const.m: -------------------------------------------------------------------------------- 1 | %% SAT_CONST Definition of astronomical and mathematical constants 2 | % 3 | % Last modified: 2022/06/03 Meysam Mahooti 4 | % 5 | % Meysam Mahooti (2025). 6 | % ECI2ECEF & ECEF2ECI Transformations 7 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 8 | % MATLAB Central File Exchange. 9 | 10 | function const = SAT_const() 11 | % Mathematical constants 12 | const.pi2 = 2*pi; % 2pi 13 | const.Rad = pi/180; % Radians per degree 14 | const.Deg = 180/pi; % Degrees per radian 15 | const.Secs = 86400; % Seconds per day 16 | const.Arcs = 3600*180/pi; % Arcseconds per radian 17 | const.TURNAS = 1296000.0; % Arcseconds in a full circle 18 | const.DAS2R = 4.848136811095359935899141e-6; % Arcseconds to radians 19 | 20 | % General 21 | const.MJD_J2000 = 51544.5; % Modified Julian Date of J2000 22 | const.T_B1950 = -0.500002108; % Epoch B1950 23 | const.c_light = 299792457.999999984; % Speed of light [m/s]; DE440 24 | const.AU = 149597870699.999988; % Astronomical unit [m]; DE440 25 | 26 | % Physical parameters of the Earth, Sun and Moon 27 | 28 | % Equatorial radius and flattening 29 | const.R_Earth = 6378.137e3; % Earth's radius [m]; WGS-84 30 | const.f_Earth = 1/298.257223563; % Flattening; WGS-84 31 | const.R_Sun = 696000.0e3; % Sun's radius [m]; DE440 32 | const.R_Moon = 1738.0e3; % Moon's radius [m]; DE440 33 | 34 | % Earth rotation (derivative of GMST at J2000; differs from inertial period by precession) 35 | const.omega_Earth = 15.04106717866910/3600*const.Rad; % [rad/s]; WGS-84 36 | 37 | % Gravitational coefficients 38 | const.GM_Earth = 398600.4415e9; % [m^3/s^2]; GGM03C & GGM03S 39 | const.GM_Sun = 132712440041.279419e9; % [m^3/s^2]; DE440 40 | const.GM_Moon = const.GM_Earth/81.3005682214972154; % [m^3/s^2]; DE440 41 | const.GM_Mercury = 22031.868551e9; % [m^3/s^2]; DE440 42 | const.GM_Venus = 324858.592000e9; % [m^3/s^2]; DE440 43 | const.GM_Mars = 42828.375816e9; % [m^3/s^2]; DE440 44 | const.GM_Jupiter = 126712764.100000e9; % [m^3/s^2]; DE440 45 | const.GM_Saturn = 37940584.841800e9; % [m^3/s^2]; DE440 46 | const.GM_Uranus = 5794556.400000e9; % [m^3/s^2]; DE440 47 | const.GM_Neptune = 6836527.100580e9; % [m^3/s^2]; DE440 48 | const.GM_Pluto = 975.500000e9; % [m^3/s^2]; DE440 49 | 50 | % Solar radiation pressure at 1 AU 51 | const.P_Sol = 1367/const.c_light; % [N/m^2] (1367 W/m^2); IERS 96 52 | 53 | end 54 | -------------------------------------------------------------------------------- /private/publish_gen_index_html.m: -------------------------------------------------------------------------------- 1 | %% PUBLISH_GEN_INDEX_HTML generate index.html for package docs 2 | % publish (generate) docs from Matlab project 3 | % called from buildfile.m 4 | % 5 | % Ref: 6 | % * https://www.mathworks.com/help/matlab/ref/publish.html 7 | % * https://www.mathworks.com/help/matlab/matlab_prog/marking-up-matlab-comments-for-publishing.html 8 | % 9 | % for package code -- assumes no classes and depth == 1 10 | % 11 | % if *.mex* files are present, publish fails 12 | 13 | function publish_gen_index_html(pkg_name, tagline, project_url, outdir) 14 | arguments 15 | pkg_name {mustBeTextScalar} 16 | tagline {mustBeTextScalar} 17 | project_url {mustBeTextScalar} 18 | outdir {mustBeTextScalar} 19 | end 20 | 21 | pkg = what("+" + pkg_name); 22 | % "+" avoids picking up cwd of same name 23 | assert(~isempty(pkg), pkg_name + " is not detected as a Matlab directory or package") 24 | 25 | %% Git info 26 | repo = gitrepo(pkg.path); 27 | git_txt = "Git branch / commit: " + repo.CurrentBranch.Name + " " + repo.LastCommit.ID{1}(1:8); 28 | 29 | %% generate docs 30 | readme = fullfile(outdir, "index.html"); 31 | 32 | if ~isfolder(outdir) 33 | mkdir(outdir); 34 | end 35 | 36 | txt = ["", ... 37 | "",... 38 | '', ... 39 | '', ... 40 | '', ... 41 | "" + pkg_name + " API", ... 42 | "", ... 43 | "", ... 44 | "

" + pkg_name + " API

", ... 45 | "

" + tagline + "

", ... 46 | "

" + git_txt + "

", ... 47 | "

Project URL: " + project_url + "

", ... 48 | "

API Reference

"]; 49 | fid = fopen(readme, 'w'); 50 | fprintf(fid, join(txt, "\n")); 51 | 52 | for sub = pkg.m.' 53 | 54 | s = sub{1}; 55 | [~, name] = fileparts(s); 56 | 57 | doc_fn = publish(pkg_name + "." + name, evalCode=false, outputDir=outdir); 58 | disp(doc_fn) 59 | 60 | % inject summary for each function into Readme.md 61 | help_txt = split(string(help(pkg_name + "." + name)), newline); 62 | words = split(strip(help_txt(1)), " "); 63 | 64 | % error if no docstring 65 | fname = words(1); 66 | assert(lower(fname) == lower(name), "fname %s does not match name %s \nis there a docstring at the top of the .m file?", fname, name) 67 | 68 | line = "" + fname + " "; 69 | if(length(words) > 1) 70 | line = line + join(words(2:end)); 71 | end 72 | 73 | req = ""; 74 | 75 | if contains(help_txt(2), "requires:") || contains(help_txt(2), "optional:") 76 | req = "(" + strip(help_txt(2), " ") + ")"; 77 | end 78 | 79 | fprintf(fid, line + " " + req + "
\n"); 80 | 81 | end 82 | 83 | fprintf(fid, " "); 84 | 85 | fclose(fid); 86 | 87 | end 88 | -------------------------------------------------------------------------------- /README.md: -------------------------------------------------------------------------------- 1 | # MatMap3d 2 | 3 | [![DOI](https://zenodo.org/badge/144219717.svg)](https://zenodo.org/badge/latestdoi/144219717) 4 | [![View matmap3d on File Exchange](https://www.mathworks.com/matlabcentral/images/matlab-file-exchange.svg)](https://www.mathworks.com/matlabcentral/fileexchange/68480-matmap3d) 5 | [![ci](https://github.com/geospace-code/matmap3d/actions/workflows/ci.yml/badge.svg)](https://github.com/geospace-code/matmap3d/actions/workflows/ci.yml) 6 | 7 | Matlab coordinate conversions for geospace ecef enu eci. 8 | Similar to Python 9 | [PyMap3D](https://github.com/geospace-code/pymap3d). 10 | 11 | ## Usage 12 | 13 | MatMap3D is setup as a 14 | [Matlab package](https://www.mathworks.com/help/matlab/matlab_oop/scoping-classes-with-packages.html), 15 | which means `import matmap3d` statements allow scoped use of this code. 16 | 17 | ```matlab 18 | 19 | [x,y,z] = matmap3d.geodetic2ecef([],lat,lon,alt) 20 | 21 | [az,el,range] = matmap3d.geodetic2aer(lat, lon, alt, observer_lat, observer_lon, observer_alt) 22 | ``` 23 | 24 | Optionally, run self-tests: 25 | 26 | ```matlab 27 | buildtool check test 28 | ``` 29 | 30 | ### Functions 31 | 32 | Popular mapping & aerospace toolbox functions ported to Matlab, where the source coordinate system (before the "2") is converted to the desired coordinate system: 33 | 34 | Abbreviations: 35 | 36 | * [AER: Azimuth, Elevation, Range](https://en.wikipedia.org/wiki/Spherical_coordinate_system) 37 | * [ECEF: Earth-centered, Earth-fixed](https://en.wikipedia.org/wiki/ECEF) 38 | * [ECI: Earth-centered Inertial](https://en.wikipedia.org/wiki/Earth-centered_inertial) 39 | * [ENU: East North Up](https://en.wikipedia.org/wiki/Axes_conventions#Ground_reference_frames:_ENU_and_NED) 40 | * [NED: North East Down](https://en.wikipedia.org/wiki/North_east_down) 41 | * [radec: right ascension, declination](https://en.wikipedia.org/wiki/Right_ascension) 42 | 43 | ### Caveats 44 | 45 | * Atmospheric effects neglected in all functions not invoking AstroPy. 46 | Would need to update code to add these input parameters (just start a GitHub Issue to request). 47 | * Planetary perturbations and nutation etc. not fully considered. 48 | 49 | These functions present a similar API of a subset of functions in the Mathworks Matlab: 50 | 51 | * [Mapping Toolbox](https://www.mathworks.com/products/mapping.html) 52 | * [Aerospace Blockset](https://www.mathworks.com/products/aerospace-blockset.html) 53 | 54 | ## Notes 55 | 56 | Python PyMap3d has more conversions. 57 | PyMap3d can be accessed from Matlab by commands like: 58 | 59 | ```matlab 60 | lla = py.pymap3d.geodetic2ecef(x,y,z) 61 | ``` 62 | 63 | ### GNU Octave 64 | 65 | GNU Octave users should consider the 66 | [Octave Mapping Toolbox](https://gnu-octave.github.io/packages/mapping/), 67 | which added similar functions in 68 | [version 1.4.2](https://sourceforge.net/p/octave/mapping/ci/Release-1.4.2/tree/NEWS). 69 | -------------------------------------------------------------------------------- /+matmap3d/vreckon.m: -------------------------------------------------------------------------------- 1 | % VRECKON - Using the WGS-84 Earth ellipsoid, travel a given distance along 2 | % a given azimuth starting at a given initial point, and return 3 | % the endpoint within a few millimeters of accuracy, using 4 | % Vincenty's algorithm. 5 | % 6 | % USAGE: 7 | % [lat2,lon2] = vreckon(lat1, lon1, s, a12) 8 | % 9 | % VARIABLES: 10 | % lat1 = inital latitude (degrees) 11 | % lon1 = initial longitude (degrees) 12 | % s = distance (meters) 13 | % a12 = intial azimuth (degrees) 14 | % lat2, lon2 = second point (degrees) 15 | % a21 = reverse azimuth (degrees), at final point facing back toward the 16 | % intial point 17 | % 18 | % Original algorithm source: 19 | % T. Vincenty, "Direct and Inverse Solutions of Geodesics on the Ellipsoid 20 | % with Application of Nested Equations", Survey Review, vol. 23, no. 176, 21 | % April 1975, pp 88-93. 22 | % Available at: http://www.ngs.noaa.gov/PUBS_LIB/inverse.pdf 23 | % 24 | % Notes: 25 | % (1) The Vincenty reckoning algorithm was transcribed verbatim into 26 | % JavaScript by Chris Veness. It was modified and translated to Matlab 27 | % by Michael Kleder. Mr. Veness's website is: 28 | % http://www.movable-type.co.uk/scripts/latlong-vincenty-direct.html 29 | % (2) Error correcting code, polar error corrections, WGS84 ellipsoid 30 | % parameters, testing, and comments by Michael Kleder. 31 | % (3) By convention, when starting at a pole, the longitude of the initial 32 | % point (otherwise meaningless) determines the longitude line along 33 | % which to traverse, and hence the longitude of the final point. 34 | % (4) The convention noted in (3) above creates a discrepancy with VDIST 35 | % when the the intial or final point is at a pole. In the VDIST 36 | % function, when traversing from a pole, the azimuth is 0 when 37 | % heading away from the south pole and 180 when heading away from the 38 | % north pole. In contrast, this VRECKON function uses the azimuth as 39 | % noted in (3) above when traversing away form a pole. 40 | % (5) In testing, where the traversal subtends no more than 178 degrees, 41 | % this function correctly inverts the VDIST function to within 0.2 42 | % millimeters of distance, 5e-10 degrees of forward azimuth, 43 | % and 5e-10 degrees of reverse azimuth. Precision reduces as test 44 | % points approach antipodal because the precision of VDIST is reduced 45 | % for nearly antipodal points. (A warning is given by VDIST.) 46 | % (6) Tested but no warranty. Use at your own risk. 47 | % (7) Ver 1.0, Michael Kleder, November 2007 48 | function [lat2,lon2,a21] = vreckon(lat1,lon1,s,a12) 49 | 50 | %% compute 51 | a = 6378137; % semimajor axis 52 | b = 6356752.31424518; % semiminor axis 53 | f = 1/298.257223563; % flattening coefficient 54 | lat1 = lat1 * .1745329251994329577e-1; % intial latitude in radians 55 | lon1 = lon1 * .1745329251994329577e-1; % intial longitude in radians 56 | % correct for errors at exact poles by adjusting 0.6 millimeters: 57 | kidx = abs(pi/2-abs(lat1)) < 1e-10; 58 | if any(kidx) 59 | lat1(kidx) = sign(lat1(kidx))*(pi/2-(1e-10)); 60 | end 61 | alpha1 = a12 * .1745329251994329577e-1; % inital azimuth in radians 62 | sinAlpha1 = sin(alpha1); 63 | cosAlpha1 = cos(alpha1); 64 | tanU1 = (1-f) * tan(lat1); 65 | cosU1 = 1 / sqrt(1 + tanU1*tanU1); 66 | sinU1 = tanU1*cosU1; 67 | sigma1 = atan2(tanU1, cosAlpha1); 68 | sinAlpha = cosU1 * sinAlpha1; 69 | cosSqAlpha = 1 - sinAlpha*sinAlpha; 70 | uSq = cosSqAlpha * (a*a - b*b) / (b*b); 71 | A = 1 + uSq/16384*(4096+uSq*(-768+uSq*(320-175*uSq))); 72 | B = uSq/1024 * (256+uSq*(-128+uSq*(74-47*uSq))); 73 | sigma = s / (b*A); 74 | sigmaP = 2*pi; 75 | while (abs(sigma-sigmaP) > 1e-12) 76 | cos2SigmaM = cos(2*sigma1 + sigma); 77 | sinSigma = sin(sigma); 78 | cosSigma = cos(sigma); 79 | deltaSigma = B*sinSigma*(cos2SigmaM+B/4*(cosSigma*(-1+... 80 | 2*cos2SigmaM*cos2SigmaM)-... 81 | B/6*cos2SigmaM*(-3+4*sinSigma*sinSigma)*(-3+... 82 | 4*cos2SigmaM*cos2SigmaM))); 83 | sigmaP = sigma; 84 | sigma = s / (b*A) + deltaSigma; 85 | end 86 | tmp = sinU1*sinSigma - cosU1*cosSigma*cosAlpha1; 87 | lat2 = atan2(sinU1*cosSigma + cosU1*sinSigma*cosAlpha1,... 88 | (1-f)*sqrt(sinAlpha*sinAlpha + tmp*tmp)); 89 | lambda = atan2(sinSigma*sinAlpha1, cosU1*cosSigma - ... 90 | sinU1*sinSigma*cosAlpha1); 91 | C = f/16*cosSqAlpha*(4+f*(4-3*cosSqAlpha)); 92 | L = lambda - (1-C) * f * sinAlpha * (sigma + C*sinSigma*(cos2SigmaM+... 93 | C*cosSigma*(-1+2*cos2SigmaM*cos2SigmaM))); 94 | lon2 = lon1 + L; 95 | % output degrees 96 | lat2 = rad2deg(lat2); 97 | lon2 = rad2deg(lon2); 98 | lon2 = mod(lon2,360); % follow [0,360] convention 99 | if nargout > 2 100 | a21 = atan2(sinAlpha, -tmp); 101 | a21 = 180 + rad2deg(a21); % note direction reversal 102 | a21=mod(a21,360); 103 | end 104 | 105 | end 106 | -------------------------------------------------------------------------------- /+matmap3d/private/NutAngles.m: -------------------------------------------------------------------------------- 1 | %% NUTANGLES Nutation in longitude and obliquity 2 | % 3 | % Input: 4 | % Mjd_TT Modified Julian Date (Terrestrial Time) 5 | % 6 | % Output: 7 | % dpsi, deps 8 | % 9 | % Last modified: 2018/01/27 Meysam Mahooti 10 | % 11 | % Meysam Mahooti (2025). 12 | % ECI2ECEF & ECEF2ECI Transformations 13 | % https://www.mathworks.com/matlabcentral/fileexchange/61957-eci2ecef-ecef2eci-transformations) 14 | % MATLAB Central File Exchange. 15 | 16 | function [dpsi, deps] = NutAngles(Mjd_TT) 17 | 18 | T = (Mjd_TT - 51544.5)/36525; 19 | T2 = T*T; 20 | T3 = T2*T; 21 | rev = 360*3600; % arcsec/revolution 22 | 23 | N_coeff = 106; 24 | C = [ 25 | % l l' F D Om dpsi *T deps *T 26 | 0, 0, 0, 0, 1,-1719960,-1742, 920250, 89 ; % 1 27 | 0, 0, 0, 0, 2, 20620, 2, -8950, 5 ; % 2 28 | -2, 0, 2, 0, 1, 460, 0, -240, 0 ; % 3 29 | 2, 0,-2, 0, 0, 110, 0, 0, 0 ; % 4 30 | -2, 0, 2, 0, 2, -30, 0, 10, 0 ; % 5 31 | 1,-1, 0,-1, 0, -30, 0, 0, 0 ; % 6 32 | 0,-2, 2,-2, 1, -20, 0, 10, 0 ; % 7 33 | 2, 0,-2, 0, 1, 10, 0, 0, 0 ; % 8 34 | 0, 0, 2,-2, 2, -131870, -16, 57360, -31 ; % 9 35 | 0, 1, 0, 0, 0, 14260, -34, 540, -1 ; % 10 36 | 0, 1, 2,-2, 2, -5170, 12, 2240, -6 ; % 11 37 | 0,-1, 2,-2, 2, 2170, -5, -950, 3 ; % 12 38 | 0, 0, 2,-2, 1, 1290, 1, -700, 0 ; % 13 39 | 2, 0, 0,-2, 0, 480, 0, 10, 0 ; % 14 40 | 0, 0, 2,-2, 0, -220, 0, 0, 0 ; % 15 41 | 0, 2, 0, 0, 0, 170, -1, 0, 0 ; % 16 42 | 0, 1, 0, 0, 1, -150, 0, 90, 0 ; % 17 43 | 0, 2, 2,-2, 2, -160, 1, 70, 0 ; % 18 44 | 0,-1, 0, 0, 1, -120, 0, 60, 0 ; % 19 45 | -2, 0, 0, 2, 1, -60, 0, 30, 0 ; % 20 46 | 0,-1, 2,-2, 1, -50, 0, 30, 0 ; % 21 47 | 2, 0, 0,-2, 1, 40, 0, -20, 0 ; % 22 48 | 0, 1, 2,-2, 1, 40, 0, -20, 0 ; % 23 49 | 1, 0, 0,-1, 0, -40, 0, 0, 0 ; % 24 50 | 2, 1, 0,-2, 0, 10, 0, 0, 0 ; % 25 51 | 0, 0,-2, 2, 1, 10, 0, 0, 0 ; % 26 52 | 0, 1,-2, 2, 0, -10, 0, 0, 0 ; % 27 53 | 0, 1, 0, 0, 2, 10, 0, 0, 0 ; % 28 54 | -1, 0, 0, 1, 1, 10, 0, 0, 0 ; % 29 55 | 0, 1, 2,-2, 0, -10, 0, 0, 0 ; % 30 56 | 0, 0, 2, 0, 2, -22740, -2, 9770, -5 ; % 31 57 | 1, 0, 0, 0, 0, 7120, 1, -70, 0 ; % 32 58 | 0, 0, 2, 0, 1, -3860, -4, 2000, 0 ; % 33 59 | 1, 0, 2, 0, 2, -3010, 0, 1290, -1 ; % 34 60 | 1, 0, 0,-2, 0, -1580, 0, -10, 0 ; % 35 61 | -1, 0, 2, 0, 2, 1230, 0, -530, 0 ; % 36 62 | 0, 0, 0, 2, 0, 630, 0, -20, 0 ; % 37 63 | 1, 0, 0, 0, 1, 630, 1, -330, 0 ; % 38 64 | -1, 0, 0, 0, 1, -580, -1, 320, 0 ; % 39 65 | -1, 0, 2, 2, 2, -590, 0, 260, 0 ; % 40 66 | 1, 0, 2, 0, 1, -510, 0, 270, 0 ; % 41 67 | 0, 0, 2, 2, 2, -380, 0, 160, 0 ; % 42 68 | 2, 0, 0, 0, 0, 290, 0, -10, 0 ; % 43 69 | 1, 0, 2,-2, 2, 290, 0, -120, 0 ; % 44 70 | 2, 0, 2, 0, 2, -310, 0, 130, 0 ; % 45 71 | 0, 0, 2, 0, 0, 260, 0, -10, 0 ; % 46 72 | -1, 0, 2, 0, 1, 210, 0, -100, 0 ; % 47 73 | -1, 0, 0, 2, 1, 160, 0, -80, 0 ; % 48 74 | 1, 0, 0,-2, 1, -130, 0, 70, 0 ; % 49 75 | -1, 0, 2, 2, 1, -100, 0, 50, 0 ; % 50 76 | 1, 1, 0,-2, 0, -70, 0, 0, 0 ; % 51 77 | 0, 1, 2, 0, 2, 70, 0, -30, 0 ; % 52 78 | 0,-1, 2, 0, 2, -70, 0, 30, 0 ; % 53 79 | 1, 0, 2, 2, 2, -80, 0, 30, 0 ; % 54 80 | 1, 0, 0, 2, 0, 60, 0, 0, 0 ; % 55 81 | 2, 0, 2,-2, 2, 60, 0, -30, 0 ; % 56 82 | 0, 0, 0, 2, 1, -60, 0, 30, 0 ; % 57 83 | 0, 0, 2, 2, 1, -70, 0, 30, 0 ; % 58 84 | 1, 0, 2,-2, 1, 60, 0, -30, 0 ; % 59 85 | 0, 0, 0,-2, 1, -50, 0, 30, 0 ; % 60 86 | 1,-1, 0, 0, 0, 50, 0, 0, 0 ; % 61 87 | 2, 0, 2, 0, 1, -50, 0, 30, 0 ; % 62 88 | 0, 1, 0,-2, 0, -40, 0, 0, 0 ; % 63 89 | 1, 0,-2, 0, 0, 40, 0, 0, 0 ; % 64 90 | 0, 0, 0, 1, 0, -40, 0, 0, 0 ; % 65 91 | 1, 1, 0, 0, 0, -30, 0, 0, 0 ; % 66 92 | 1, 0, 2, 0, 0, 30, 0, 0, 0 ; % 67 93 | 1,-1, 2, 0, 2, -30, 0, 10, 0 ; % 68 94 | -1,-1, 2, 2, 2, -30, 0, 10, 0 ; % 69 95 | -2, 0, 0, 0, 1, -20, 0, 10, 0 ; % 70 96 | 3, 0, 2, 0, 2, -30, 0, 10, 0 ; % 71 97 | 0,-1, 2, 2, 2, -30, 0, 10, 0 ; % 72 98 | 1, 1, 2, 0, 2, 20, 0, -10, 0 ; % 73 99 | -1, 0, 2,-2, 1, -20, 0, 10, 0 ; % 74 100 | 2, 0, 0, 0, 1, 20, 0, -10, 0 ; % 75 101 | 1, 0, 0, 0, 2, -20, 0, 10, 0 ; % 76 102 | 3, 0, 0, 0, 0, 20, 0, 0, 0 ; % 77 103 | 0, 0, 2, 1, 2, 20, 0, -10, 0 ; % 78 104 | -1, 0, 0, 0, 2, 10, 0, -10, 0 ; % 79 105 | 1, 0, 0,-4, 0, -10, 0, 0, 0 ; % 80 106 | -2, 0, 2, 2, 2, 10, 0, -10, 0 ; % 81 107 | -1, 0, 2, 4, 2, -20, 0, 10, 0 ; % 82 108 | 2, 0, 0,-4, 0, -10, 0, 0, 0 ; % 83 109 | 1, 1, 2,-2, 2, 10, 0, -10, 0 ; % 84 110 | 1, 0, 2, 2, 1, -10, 0, 10, 0 ; % 85 111 | -2, 0, 2, 4, 2, -10, 0, 10, 0 ; % 86 112 | -1, 0, 4, 0, 2, 10, 0, 0, 0 ; % 87 113 | 1,-1, 0,-2, 0, 10, 0, 0, 0 ; % 88 114 | 2, 0, 2,-2, 1, 10, 0, -10, 0 ; % 89 115 | 2, 0, 2, 2, 2, -10, 0, 0, 0 ; % 90 116 | 1, 0, 0, 2, 1, -10, 0, 0, 0 ; % 91 117 | 0, 0, 4,-2, 2, 10, 0, 0, 0 ; % 92 118 | 3, 0, 2,-2, 2, 10, 0, 0, 0 ; % 93 119 | 1, 0, 2,-2, 0, -10, 0, 0, 0 ; % 94 120 | 0, 1, 2, 0, 1, 10, 0, 0, 0 ; % 95 121 | -1,-1, 0, 2, 1, 10, 0, 0, 0 ; % 96 122 | 0, 0,-2, 0, 1, -10, 0, 0, 0 ; % 97 123 | 0, 0, 2,-1, 2, -10, 0, 0, 0 ; % 98 124 | 0, 1, 0, 2, 0, -10, 0, 0, 0 ; % 99 125 | 1, 0,-2,-2, 0, -10, 0, 0, 0 ; % 100 126 | 0,-1, 2, 0, 1, -10, 0, 0, 0 ; % 101 127 | 1, 1, 0,-2, 1, -10, 0, 0, 0 ; % 102 128 | 1, 0,-2, 2, 0, -10, 0, 0, 0 ; % 103 129 | 2, 0, 0, 2, 0, 10, 0, 0, 0 ; % 104 130 | 0, 0, 2, 4, 2, -10, 0, 0, 0 ; % 105 131 | 0, 1, 0, 1, 0, 10, 0, 0, 0 % 106 132 | ]; 133 | 134 | % Mean arguments of luni-solar motion 135 | 136 | % l mean anomaly of the Moon 137 | % l' mean anomaly of the Sun 138 | % F mean argument of latitude 139 | % D mean longitude elongation of the Moon from the Sun 140 | % Om mean longitude of the ascending node 141 | l = mod ( 485866.733 + (1325.0*rev + 715922.633)*T ... 142 | + 31.310*T2 + 0.064*T3, rev ); 143 | lp = mod ( 1287099.804 + ( 99.0*rev + 1292581.224)*T ... 144 | - 0.577*T2 - 0.012*T3, rev ); 145 | F = mod ( 335778.877 + (1342.0*rev + 295263.137)*T ... 146 | - 13.257*T2 + 0.011*T3, rev ); 147 | D = mod ( 1072261.307 + (1236.0*rev + 1105601.328)*T ... 148 | - 6.891*T2 + 0.019*T3, rev ); 149 | Om = mod ( 450160.280 - ( 5.0*rev + 482890.539)*T ... 150 | + 7.455*T2 + 0.008*T3, rev ); 151 | 152 | % Nutation in longitude and obliquity [rad] 153 | dpsi = 0; 154 | deps = 0; 155 | Arcs = 3600*180/pi; 156 | 157 | for i=1:N_coeff 158 | arg = ( C(i,1)*l+C(i,2)*lp+C(i,3)*F+C(i,4)*D+C(i,5)*Om )/ Arcs; 159 | dpsi = dpsi + ( C(i,6)+C(i,7)*T ) * sin(arg); 160 | deps = deps + ( C(i,8)+C(i,9)*T ) * cos(arg); 161 | end 162 | 163 | dpsi = 1e-5 * dpsi / Arcs; 164 | deps = 1e-5 * deps / Arcs; 165 | 166 | end 167 | -------------------------------------------------------------------------------- /+matmap3d/vdist.m: -------------------------------------------------------------------------------- 1 | %% VDIST - Using the WGS-84 Earth ellipsoid, compute the distance between two points 2 | % 3 | % within a few millimeters of accuracy, compute forward 4 | % azimuth, and compute backward azimuth, all using a vectorized 5 | % version of Vincenty's algorithm. 6 | % 7 | % s = vdist(lat1,lon1,lat2,lon2) 8 | % [s,a12] = vdist(lat1,lon1,lat2,lon2) 9 | % [s,a12,a21] = vdist(lat1,lon1,lat2,lon2) 10 | % 11 | % * s = distance in meters (inputs may be scalars, vectors, or matrices) 12 | % * a12 = azimuth in degrees from first point to second point (forward) 13 | % * a21 = azimuth in degrees from second point to first point (backward) 14 | % (Azimuths are in degrees clockwise from north.) 15 | % * lat1 = GEODETIC latitude of first point (degrees) 16 | % * lon1 = longitude of first point (degrees) 17 | % * lat2, lon2 = second point (degrees) 18 | % 19 | % from https://www.mathworks.com/matlabcentral/fileexchange/8607-vectorized-geodetic-distance-and-azimuth-on-the-wgs84-earth-ellipsoid 20 | 21 | function varargout = vdist(lat1,lon1,lat2,lon2) 22 | % 23 | %%% Original algorithm source: 24 | % 25 | % T. Vincenty, "Direct and Inverse Solutions of Geodesics on the Ellipsoid with Application of Nested Equations", 26 | % Survey Review, vol. 23, no. 176, April 1975, pp 88-93. 27 | % 28 | % Notes: 29 | % 30 | % # lat1,lon1,lat2,lon2 can be any (identical) size/shape. Outputs will have the same size and shape. 31 | % # Error correcting code, convergence failure traps, antipodal corrections, polar error corrections, WGS84 ellipsoid parameters, testing, and comments: Michael Kleder, 2004. 32 | % # Azimuth implementation (including quadrant abiguity resolution) and code vectorization, Michael Kleder, Sep 2005. 33 | % # Vectorization is convergence sensitive; that is, quantities which have already converged to within tolerance are not recomputed during subsequent iterations (while other quantities are still converging). 34 | % # Vincenty describes his distance algorithm as precise to within 0.01 millimeters, subject to the ellipsoidal model. 35 | % # For distance calculations, essentially antipodal points are 36 | % treated as exactly antipodal, potentially reducing accuracy 37 | % slightly. 38 | % # Distance failures for points exactly at the poles are 39 | % eliminated by moving the points by 0.6 millimeters. 40 | % # The Vincenty distance algorithm was transcribed verbatim by 41 | % Peter Cederholm, August 12, 2003. It was modified and 42 | % translated to English by Michael Kleder. 43 | % Mr. Cederholm's website is http://www.plan.aau.dk/~pce/ 44 | % # Distances agree with the Mapping Toolbox, version 2.2 (R14SP3) 45 | % with a max relative difference of about 5e-9, except when the 46 | % two points are nearly antipodal, and except when one point is 47 | % near the equator and the two longitudes are nearly 180 degrees 48 | % apart. This function (vdist) is more accurate in such cases. 49 | % For example, note this difference (as of this writing): 50 | % >>vdist(0.2,305,15,125) 51 | % 18322827.0131551 52 | % >>distance(0.2,305,15,125,[6378137 0.08181919]) 53 | % 0 54 | % # Azimuths FROM the north pole (either forward starting at the 55 | % north pole or backward when ending at the north pole) are set 56 | % to 180 degrees by convention. Azimuths FROM the south pole are 57 | % set to 0 degrees by convention. 58 | % # Azimuths agree with the Mapping Toolbox, version 2.2 (R14SP3) 59 | % to within about a hundred-thousandth of a degree, except when 60 | % traversing to or from a pole, where the convention for this 61 | % function is described in (10), and except in the cases noted 62 | % above in (9). 63 | % # No warranties; use at your own risk. 64 | 65 | keepsize = size(lat1); 66 | % Supply WGS84 earth ellipsoid axis lengths in meters: 67 | a = 6378137; % definitionally 68 | b = 6356752.31424518; % computed from WGS84 earth flattening coefficient 69 | % preserve true input latitudes: 70 | lat1tr = lat1; 71 | lat2tr = lat2; 72 | % convert inputs in degrees to radians: 73 | lat1 = lat1 * 0.0174532925199433; 74 | lon1 = lon1 * 0.0174532925199433; 75 | lat2 = lat2 * 0.0174532925199433; 76 | lon2 = lon2 * 0.0174532925199433; 77 | % correct for errors at exact poles by adjusting 0.6 millimeters: 78 | kidx = abs(pi/2-abs(lat1)) < 1e-10; 79 | if any(kidx) 80 | lat1(kidx) = sign(lat1(kidx))*(pi/2-(1e-10)); 81 | end 82 | kidx = abs(pi/2-abs(lat2)) < 1e-10; 83 | if any(kidx) 84 | lat2(kidx) = sign(lat2(kidx))*(pi/2-(1e-10)); 85 | end 86 | f = (a-b)/a; 87 | U1 = atan((1-f)*tan(lat1)); 88 | U2 = atan((1-f)*tan(lat2)); 89 | lon1 = mod(lon1,2*pi); 90 | lon2 = mod(lon2,2*pi); 91 | L = abs(lon2-lon1); 92 | kidx = L > pi; 93 | if any(kidx) 94 | L(kidx) = 2*pi - L(kidx); 95 | end 96 | lambda = L; 97 | lambdaold = 0*lat1; 98 | itercount = 0; 99 | notdone = logical(1+0*lat1); 100 | alpha = 0*lat1; 101 | sigma = 0*lat1; 102 | cos2sigmam = 0*lat1; 103 | C = 0*lat1; 104 | warninggiven = false; 105 | while any(notdone) % force at least one execution 106 | %disp(['lambda(21752) = ' num2str(lambda(21752),20)]); 107 | itercount = itercount+1; 108 | if itercount > 50 109 | if ~warninggiven 110 | warning(['Essentially antipodal points encountered. ' ... 111 | 'Precision may be reduced slightly.']); 112 | end 113 | lambda(notdone) = pi; 114 | break 115 | end 116 | lambdaold(notdone) = lambda(notdone); 117 | 118 | sinsigma(notdone) = sqrt((cos(U2(notdone)).*sin(lambda(notdone)))... 119 | .^2+(cos(U1(notdone)).*sin(U2(notdone))-sin(U1(notdone)).*... 120 | cos(U2(notdone)).*cos(lambda(notdone))).^2); %#ok<*AGROW> 121 | 122 | cossigma(notdone) = sin(U1(notdone)).*sin(U2(notdone))+... 123 | cos(U1(notdone)).*cos(U2(notdone)).*cos(lambda(notdone)); 124 | 125 | % eliminate rare imaginary portions at limit of numerical precision: 126 | sinsigma(notdone)=real(sinsigma(notdone)); 127 | cossigma(notdone)=real(cossigma(notdone)); 128 | sigma(notdone) = atan2(sinsigma(notdone),cossigma(notdone)); 129 | alpha(notdone) = asin(cos(U1(notdone)).*cos(U2(notdone)).*... 130 | sin(lambda(notdone))./sin(sigma(notdone))); 131 | cos2sigmam(notdone) = cos(sigma(notdone))-2*sin(U1(notdone)).*... 132 | sin(U2(notdone))./cos(alpha(notdone)).^2; 133 | C(notdone) = f/16*cos(alpha(notdone)).^2.*(4+f*(4-3*... 134 | cos(alpha(notdone)).^2)); 135 | lambda(notdone) = L(notdone)+(1-C(notdone)).*f.*sin(alpha(notdone))... 136 | .*(sigma(notdone)+C(notdone).*sin(sigma(notdone)).*... 137 | (cos2sigmam(notdone)+C(notdone).*cos(sigma(notdone)).*... 138 | (-1+2.*cos2sigmam(notdone).^2))); 139 | %disp(['then, lambda(21752) = ' num2str(lambda(21752),20)]); 140 | % correct for convergence failure in the case of essentially antipodal 141 | % points 142 | if any(lambda(notdone) > pi) 143 | warning('Essentially antipodal points encountered. Precision may be reduced slightly.') 144 | warninggiven = true; 145 | lambdaold(lambda>pi) = pi; 146 | lambda(lambda>pi) = pi; 147 | end 148 | notdone = abs(lambda-lambdaold) > 1e-12; 149 | end 150 | u2 = cos(alpha).^2.*(a^2-b^2)/b^2; 151 | A = 1+u2./16384.*(4096+u2.*(-768+u2.*(320-175.*u2))); 152 | B = u2./1024.*(256+u2.*(-128+u2.*(74-47.*u2))); 153 | deltasigma = B.*sin(sigma).*(cos2sigmam+B./4.*(cos(sigma).*(-1+2.*... 154 | cos2sigmam.^2)-B./6.*cos2sigmam.*(-3+4.*sin(sigma).^2).*(-3+4*... 155 | cos2sigmam.^2))); 156 | varargout{1} = b.*A.*(sigma-deltasigma); 157 | if nargout > 1 158 | % From point #1 to point #2 159 | % correct sign of lambda for azimuth calcs: 160 | lambda = abs(lambda); 161 | kidx=sign(sin(lon2-lon1)) .* sign(sin(lambda)) < 0; 162 | lambda(kidx) = -lambda(kidx); 163 | numer = cos(U2).*sin(lambda); 164 | denom = cos(U1).*sin(U2)-sin(U1).*cos(U2).*cos(lambda); 165 | a12 = atan2(numer,denom); 166 | kidx = a12<0; 167 | a12(kidx)=a12(kidx)+2*pi; 168 | % from poles: 169 | a12(lat1tr <= -90) = 0; 170 | a12(lat1tr >= 90 ) = pi; 171 | varargout{2} = rad2deg(a12); % to degrees 172 | end 173 | if nargout > 2 174 | % From point #2 to point #1 175 | % correct sign of lambda for azimuth calcs: 176 | lambda = abs(lambda); 177 | kidx=sign(sin(lon1-lon2)) .* sign(sin(lambda)) < 0; 178 | lambda(kidx)=-lambda(kidx); 179 | numer = cos(U1).*sin(lambda); 180 | denom = sin(U1).*cos(U2)-cos(U1).*sin(U2).*cos(lambda); 181 | a21 = atan2(numer,denom); 182 | kidx=a21<0; 183 | a21(kidx)= a21(kidx)+2*pi; 184 | % backwards from poles: 185 | a21(lat2tr >= 90) = pi; 186 | a21(lat2tr <= -90) = 0; 187 | varargout{3} = reshape(rad2deg(a21),keepsize); % to degrees 188 | end 189 | 190 | end 191 | -------------------------------------------------------------------------------- /test/TestUnit.m: -------------------------------------------------------------------------------- 1 | classdef (SharedTestFixtures={ matlab.unittest.fixtures.PathFixture(fileparts(fileparts(mfilename('fullpath'))))}) ... 2 | TestUnit < matlab.unittest.TestCase 3 | 4 | properties 5 | atol = 1e-9 6 | rtol = 1e-6 7 | angleUnit='d' 8 | 9 | lat = 42 10 | lon = -82 11 | alt = 200 12 | x0 = 660.6752518e3 13 | y0 = -4700.9486832e3 14 | z0 = 4245.7376622e3 15 | lat1 = 42.0026 16 | lon1 = -81.9978 17 | alt1 = 1.1397e3 18 | atol_dist = 1e-3 % 1 mm 19 | 20 | er = 186.277521 21 | nr = 286.84222 22 | ur = 939.69262 23 | az = 33 24 | el = 70 25 | srange = 1e3 26 | 27 | xl = 6.609301927610815e+5 28 | yl = -4.701424222957011e6 29 | zl = 4.246579604632881e+06 % aer2ecef 30 | 31 | a90 = 90 32 | E 33 | end 34 | 35 | 36 | methods (TestMethodSetup) 37 | function setupEllipsoid(tc) 38 | tc.E = matmap3d.wgs84Ellipsoid(); 39 | end 40 | end 41 | 42 | 43 | methods(Test) 44 | 45 | function test_ellipsoid(tc) 46 | tc.verifyClass(tc.E, 'matmap3d.referenceEllipsoid') 47 | end 48 | 49 | function test_geodetic2ecef(tc) 50 | 51 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, tc.lat, tc.lon, tc.alt, tc.angleUnit); 52 | tc.verifyEqual([x,y,z], [tc.x0, tc.y0, tc.z0], AbsTol=tc.atol, RelTol=tc.rtol) 53 | 54 | [x,y,z] = matmap3d.geodetic2ecef([], 0,0,-1); 55 | tc.verifyEqual([x,y,z], [tc.E.SemimajorAxis-1,0,0], AbsTol=tc.atol, RelTol=tc.rtol) 56 | 57 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, 0,90,-1); 58 | tc.verifyEqual([x,y,z], [0, tc.E.SemimajorAxis-1,0], AbsTol=tc.atol, RelTol=tc.rtol) 59 | 60 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, 0,-90,-1); 61 | tc.verifyEqual([x,y,z], [0, -tc.E.SemimajorAxis+1,0], AbsTol=tc.atol, RelTol=tc.rtol) 62 | 63 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, 90,0,-1); 64 | tc.verifyEqual([x,y,z], [0, 0, tc.E.SemiminorAxis-1], AbsTol=tc.atol, RelTol=tc.rtol) 65 | 66 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, 90,15,-1); 67 | tc.verifyEqual([x,y,z], [0,0, tc.E.SemiminorAxis-1], AbsTol=tc.atol, RelTol=tc.rtol) 68 | 69 | [x,y,z] = matmap3d.geodetic2ecef(tc.E, -90,0,-1); 70 | tc.verifyEqual([x,y,z], [0,0, -tc.E.SemiminorAxis+1], AbsTol=tc.atol, RelTol=tc.rtol) 71 | 72 | end 73 | 74 | function test_ecef2geodetic(tc) 75 | 76 | ea = tc.E.SemimajorAxis; 77 | eb = tc.E.SemiminorAxis; 78 | 79 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, tc.x0, tc.y0, tc.z0, tc.angleUnit); 80 | tc.verifyEqual([lt, ln, at], [tc.lat, tc.lon, tc.alt], AbsTol=tc.atol, RelTol=tc.rtol) 81 | 82 | [lt, ln, at] = matmap3d.ecef2geodetic([], ea-1, 0, 0); 83 | tc.verifyEqual([lt, ln, at], [0, 0, -1], AbsTol=tc.atol, RelTol=tc.rtol) 84 | 85 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, 0, ea-1, 0); 86 | tc.verifyEqual([lt, ln, at], [0, 90, -1], AbsTol=tc.atol, RelTol=tc.rtol) 87 | 88 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, 0, 0, eb-1); 89 | tc.verifyEqual([lt, ln, at], [90, 0, -1], AbsTol=tc.atol, RelTol=tc.rtol) 90 | 91 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, 0, 0, -eb+1); 92 | tc.verifyEqual([lt, ln, at], [-90, 0, -1], AbsTol=tc.atol, RelTol=tc.rtol) 93 | 94 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, -ea+1, 0, 0); 95 | tc.verifyEqual([lt, ln, at], [0, 180, -1], AbsTol=tc.atol, RelTol=tc.rtol) 96 | 97 | [lt, ln, at] = matmap3d.ecef2geodetic(tc.E, (ea-1000)/sqrt(2), (ea-1000)/sqrt(2), 0); 98 | tc.verifyEqual([lt,ln,at], [0,45,-1000], AbsTol=tc.atol, RelTol=tc.rtol) 99 | 100 | end 101 | 102 | function test_enu2aer(tc) 103 | 104 | [a, e, r] = matmap3d.enu2aer(tc.er, tc.nr, tc.ur, tc.angleUnit); 105 | tc.verifyEqual([a,e,r], [tc.az, tc.el, tc.srange] , AbsTol=tc.atol, RelTol=tc.rtol) 106 | 107 | [a, e, r] = matmap3d.enu2aer(1, 0, 0, tc.angleUnit); 108 | tc.verifyEqual([a,e,r], [tc.a90, 0, 1], AbsTol=tc.atol, RelTol=tc.rtol) 109 | end 110 | 111 | function test_aer2enu(tc) 112 | 113 | [e,n,u] = matmap3d.aer2enu(tc.az, tc.el, tc.srange, tc.angleUnit); 114 | tc.verifyEqual([e,n,u], [tc.er, tc.nr, tc.ur], AbsTol=tc.atol, RelTol=tc.rtol) 115 | 116 | [n1,e1,d] = matmap3d.aer2ned(tc.az, tc.el, tc.srange, tc.angleUnit); 117 | tc.verifyEqual([e,n,u], [e1,n1,-d]) 118 | 119 | [a,e,r] = matmap3d.enu2aer(e,n,u, tc.angleUnit); 120 | tc.verifyEqual([a,e,r], [tc.az, tc.el, tc.srange], AbsTol=tc.atol, RelTol=tc.rtol) 121 | end 122 | 123 | function test_ecef2aer(tc) 124 | 125 | [a, e, r] = matmap3d.ecef2aer(tc.xl, tc.yl, tc.zl, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 126 | % round-trip 127 | tc.verifyEqual([a,e,r], [tc.az, tc.el, tc.srange], AbsTol=tc.atol, RelTol=tc.rtol) 128 | 129 | % singularity check 130 | [a, e, r] = matmap3d.ecef2aer(tc.E.SemimajorAxis-1, 0, 0, 0,0,0, tc.E, tc.angleUnit); 131 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 132 | 133 | [a, e, r] = matmap3d.ecef2aer(-tc.E.SemimajorAxis+1, 0, 0, 0, 2*tc.a90,0, tc.E, tc.angleUnit); 134 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 135 | 136 | [a, e, r] = matmap3d.ecef2aer(0, tc.E.SemimajorAxis-1, 0,0, tc.a90,0, tc.E, tc.angleUnit); 137 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 138 | 139 | [a, e, r] = matmap3d.ecef2aer(0, -tc.E.SemimajorAxis+1, 0,0, -tc.a90,0, tc.E, tc.angleUnit); 140 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 141 | 142 | [a, e, r] = matmap3d.ecef2aer(0, 0, tc.E.SemiminorAxis-1, tc.a90, 0, 0, tc.E, tc.angleUnit); 143 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 144 | 145 | [a, e, r] = matmap3d.ecef2aer(0, 0, -tc.E.SemiminorAxis+1,-tc.a90,0,0, tc.E, tc.angleUnit); 146 | tc.verifyEqual([a,e,r], [0, -tc.a90, 1], AbsTol=tc.atol, RelTol=tc.rtol) 147 | 148 | [a, e, r] = matmap3d.ecef2aer((tc.E.SemimajorAxis-1000)/sqrt(2), (tc.E.SemimajorAxis-1000)/sqrt(2), 0, 0, 45, 0); 149 | tc.verifyEqual([a,e,r],[0,-90,1000], AbsTol=tc.atol, RelTol=tc.rtol) 150 | 151 | [x,y,z] = matmap3d.aer2ecef(tc.az, tc.el, tc.srange, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 152 | tc.verifyEqual([x,y,z], [tc.xl, tc.yl, tc.zl], AbsTol=tc.atol, RelTol=tc.rtol) 153 | 154 | [a,e,r] = matmap3d.ecef2aer(x,y,z, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 155 | tc.verifyEqual([a,e,r], [tc.az, tc.el, tc.srange], AbsTol=tc.atol, RelTol=tc.rtol) 156 | end 157 | 158 | function test_geodetic2aer(tc) 159 | 160 | [lt,ln,at] = matmap3d.aer2geodetic(tc.az, tc.el, tc.srange, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 161 | tc.verifyEqual([lt,ln,at], [tc.lat1, tc.lon1, tc.alt1], AbsTol=2*tc.atol_dist) 162 | 163 | [a, e, r] = matmap3d.geodetic2aer(lt,ln,at, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); % round-trip 164 | tc.verifyEqual([a,e,r], [tc.az, tc.el, tc.srange], AbsTol=tc.atol, RelTol=tc.rtol) 165 | end 166 | 167 | 168 | function test_geodetic2enu(tc) 169 | 170 | [e, n, u] = matmap3d.geodetic2enu(tc.lat, tc.lon, tc.alt-1, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 171 | tc.verifyEqual([e,n,u], [0,0,-1], AbsTol=tc.atol, RelTol=tc.rtol) 172 | 173 | [lt, ln, at] = matmap3d.enu2geodetic(e,n,u,tc.lat,tc.lon,tc.alt, tc.E, tc.angleUnit); % round-trip 174 | tc.verifyEqual([lt, ln, at],[tc.lat, tc.lon, tc.alt-1], AbsTol=tc.atol, RelTol=tc.rtol) 175 | 176 | [n1, e1, d] = matmap3d.geodetic2ned(tc.lat, tc.lon, tc.alt-1, tc.lat, tc.lon, tc.alt, tc.E, tc.angleUnit); 177 | tc.verifyEqual([e,n,u],[e1,n1,-d], AbsTol=tc.atol, RelTol=tc.rtol) 178 | end 179 | 180 | 181 | function test_enu2ecef(tc) 182 | 183 | [x, y, z] = matmap3d.enu2ecef(tc.er, tc.nr, tc.ur, tc.lat,tc.lon,tc.alt, tc.E, tc.angleUnit); 184 | tc.verifyEqual([x,y,z],[tc.xl, tc.yl, tc.zl], AbsTol=tc.atol, RelTol=tc.rtol) 185 | 186 | [e,n,u] = matmap3d.ecef2enu(x,y,z,tc.lat,tc.lon,tc.alt, tc.E, tc.angleUnit); % round-trip 187 | tc.verifyEqual([e,n,u],[tc.er, tc.nr, tc.ur], AbsTol=tc.atol, RelTol=tc.rtol) 188 | 189 | [n1, e1, d] = matmap3d.ecef2ned(x,y,z,tc.lat,tc.lon,tc.alt, tc.E, tc.angleUnit); 190 | tc.verifyEqual([e,n,u],[e1,n1,-d]) 191 | 192 | end 193 | 194 | 195 | function test_enu_vector(tc) 196 | uvw = [27.9799, -1.0990, -15.7723]; 197 | lat0 = 17.4114; 198 | lon0 = 78.270; 199 | [e,n,u] = matmap3d.ecef2enuv(uvw(1), uvw(2), uvw(3), lat0, lon0); 200 | [u,v,w] = matmap3d.enu2ecefv(e,n,u, lat0, lon0); 201 | tc.verifyEqual([u,v,w], uvw, AbsTol=1e-3) 202 | end 203 | 204 | 205 | 206 | function test_lookAtSpheroid(tc) 207 | az5 = [0., 10., 125.]; 208 | tilt = [30, 45, 90]; 209 | 210 | [lat5, lon5, rng5] = matmap3d.lookAtSpheroid(tc.lat, tc.lon, tc.alt, az5, 0.); 211 | tc.verifyEqual(lat5, repmat(tc.lat, 1, 3), AbsTol=tc.atol, RelTol=tc.rtol) 212 | tc.verifyEqual(lon5, repmat(tc.lon, 1, 3), AbsTol=tc.atol, RelTol=tc.rtol) 213 | tc.verifyEqual(rng5, repmat(tc.alt, 1, 3), AbsTol=tc.atol, RelTol=tc.rtol) 214 | 215 | 216 | [lat5, lon5, rng5] = matmap3d.lookAtSpheroid(tc.lat, tc.lon, tc.alt, az5, tilt); 217 | 218 | truth = [42.00103959, tc.lon, 230.9413173; 219 | 42.00177328, -81.9995808, 282.84715651; 220 | nan, nan, nan]; 221 | 222 | tc.verifyEqual([lat5, lon5, rng5], truth(:).', AbsTol=tc.atol, RelTol=tc.rtol) 223 | end 224 | 225 | 226 | function test_eci2ecef(tc) 227 | utc = datetime(2019, 1, 4, 12,0,0); 228 | eci = [-2981784; 5207055; 3161595]; 229 | r_ecef = matmap3d.eci2ecef(utc, eci); 230 | tc.verifyEqual(r_ecef, [-5762654.65677142; -1682688.09235503; 3156027.98313692], RelTol=2e-5) 231 | end 232 | 233 | 234 | function test_naive(tc) 235 | utc = datetime(2019, 1, 4, 12,0,0); 236 | [x,y,z] = matmap3d.ecef2eci_naive(utc, tc.x0, tc.y0, tc.z0); 237 | [x1,y1,z1] = matmap3d.eci2ecef_naive(utc, x,y,z); 238 | tc.verifyEqual([x1,y1,z1], [tc.x0,tc.y0,tc.z0], RelTol=1e-9) 239 | end 240 | 241 | 242 | function test_ecef2eci(tc) 243 | ecef = [-5762640; -1682738; 3156028]; 244 | utc = datetime(2019, 1, 4, 12,0,0); 245 | r_eci = matmap3d.ecef2eci(utc, ecef); 246 | tc.verifyEqual(r_eci, [-2981829.07728415; 5207029.04470791; 3161595.0981204], RelTol=1e-5) 247 | end 248 | 249 | 250 | function test_ecef2eci_null(tc) 251 | [x, y, z] = matmap3d.geodetic2ecef(tc.E, 0, 0, 0); 252 | t = datetime(2000, 1, 1, 12, 0, 0, TimeZone='UTCLeapSeconds'); 253 | r_eci = matmap3d.ecef2eci(t, [x;y;z]); 254 | tc.verifyEqual(r_eci, [1158174.72525987; -6272101.9503871; -143.138407305876], RelTol=tc.rtol) 255 | end 256 | 257 | 258 | function test_eci2aer(tc) 259 | eci = [-3.8454e8, -0.5099e8, -0.3255e8]; 260 | utc = datetime(1969, 7, 20, 21, 17, 40); 261 | lla = [28.4, -80.5, 2.7]; 262 | 263 | [a, e, r] = matmap3d.eci2aer(utc, eci(1), eci(2), eci(3), lla(1), lla(2), lla(3)); 264 | tc.verifyEqual([a, e, r], [162.548042074738, 55.1223823017527, 384013992.914642], RelTol=tc.rtol) 265 | end 266 | 267 | function test_aer2eci(tc) 268 | aer = [162.55, 55.12, 384013940.9]; 269 | lla = [28.4, -80.5, 2.7]; 270 | utc = datetime(1969, 7, 20, 21, 17, 40); 271 | 272 | [x,y,z] = matmap3d.aer2eci(utc, aer(1), aer(2), aer(3), lla(1), lla(2), lla(3)); 273 | tc.verifyEqual([x, y, z], [-384538755.067354, -50986804.9565394, -32567306.0200869], RelTol=2e5) 274 | end 275 | 276 | end 277 | end 278 | --------------------------------------------------------------------------------