geodezyx.reffram package
Submodules
geodezyx.reffram.geometry module
@author: psakic
This sub-module of geodezyx.reffram contains functions for low-level geometry operations.
it can be imported directly with: from geodezyx import reffram
The GeodeZYX Toolbox is a software for simple but useful functions for Geodesy and Geophysics under the GNU LGPL v3 License
Copyright (C) 2019 Pierre Sakic et al. (IPGP, sakic@ipgp.fr) GitHub repository : https://github.com/GeodeZYX/geodezyx-toolbox
- geodezyx.reffram.geometry.BL_from_points(listpointin)
From a list of 2-D or 3-dD points, returns the a matrix with distance between each points
- Parameters:
listpointin (list or numpy.array) – List of N 2D or 3D points [[x1,y1,z1] … [xn , yn , zn]]
- Returns:
BL – matrix with distances between each points
- Return type:
numpy.array
- geodezyx.reffram.geometry.R2_calc(y_obs, y_fit, with_r2_bis=False)
- geodezyx.reffram.geometry.R2_from_a_line_regress(Xobs, Yobs, a, b)
- geodezyx.reffram.geometry.calc_pos_speed_itrf(x0, y0, z0, t0, vx, vy, vz, t)
just a wrapper of itrf_speed_calc for legacy reasons
- geodezyx.reffram.geometry.circle_draw(xc, yc, R, N)
- geodezyx.reffram.geometry.estimated_autocorrelation(x)
http://stackoverflow.com/q/14297012/190597 http://en.wikipedia.org/wiki/Autocorrelation#Estimation
- geodezyx.reffram.geometry.group_consecutives(vals, step=1)
Return list of consecutive lists of numbers from vals (number list).
- geodezyx.reffram.geometry.guess_seq_len(seq)
- geodezyx.reffram.geometry.helmert_trans_apply(Xin, SevenParam_in, legacy_mode=False)
Apply an Helmert transformation (7-parameters) to a set of points
- Parameters:
Xin (list or np.array) – input point set list of N (x,y,z) points, or an (N,3)-shaped numpy array.
SevenParam_in (7 element list or array) – 7 Helmert params. : x,y,z translations, x,y,z rotations, scale.
legacy_mode (bool, optional) – Use a non-optimized and slow computation approach (but same result). This option should be removed in the Future. The default is False.
- Returns:
Xout – output transformed points. Same type as the input.
- Return type:
list/array of N (x,y,z) points
- geodezyx.reffram.geometry.helmert_trans_estim(X1list, X2list, Weights=[])
estimates 7 parameters of a 3D Helmert transformation between a set of points X1 and a set of points X2 (compute transformation X1 => X2)
- Parameters:
X2list (X1list &) – Input point set list of N (x,y,z) points, or an (N,3)-shaped numpy array
Weights (list of N Weights,) – or an numpy array of shape (N,)
- Returns:
HParam – 7 Helmert params. : x,y,z translations, x,y,z rotations, scale translations are given in meters, rotations in arcsec, scale in unitless ratio
A – Design matrix
l – Differences X2 - X1 (before transformation !!!)
References
https://elib.uni-stuttgart.de/bitstream/11682/9661/1/BscThesis_GaoYueqing.pdf
- geodezyx.reffram.geometry.helmert_trans_estim_minimisation(X1in, X2in, HParam_apri=array([0., 0., 0., 0., 0., 0., 0.]), L1norm=True, tol=1e-09, full_output=False, method='Powell')
estimates 7 parameters of a 3D Helmert transformation between a set of points X1 and a set of points X2 (compute transformation X1 => X2) using a Minimization approach (and not a Least Square inversion)
- Parameters:
X2in (X1in &) – Input point set list of N (x,y,z) points, or an (N,3)-shaped numpy array
HParam_apri (list) – list of 7 values, The Apriori for the Helmert parameter
L1norm (bool) – Use the L1-norm as a criteria, use the quadratic sum instead if False
tol (float) – tolerence for the convergence
full_output (bool) – return only the result if True, return the scipy optimize result if False
method (str, optional) – minimization method. see scipy.optimize.minimize for details The default is “Powell”.
- Returns:
7 Helmert params. : x,y,z translations, x,y,z rotations, scale translations are given in meters, rotations in arcsec, scale in unitless ratio
- Return type:
Res
- geodezyx.reffram.geometry.helmert_trans_estim_minimisation_scalar(X1, X2, HParam_opti_apriori, L1norm=True, itera=2)
Estimates the Helmert parameters but based on a minimisation approach between X1 and X2 (as suggested in the IGS combination software)
NOT STABLE AVOID THE USE (it does the optimization of the 1st parameter, and then the optimization of the 2nd one etc, etc…)
- Parameters:
X2 (X1 &) – Input point sets
HParam_opti_apriori (list of 7 values) – The Apriori for the Helmert parameter
L1norm (bool) – Use the L1-norm as a criteria, use the quadratic sum instead if False
itera (int, optional) – number of iterations. The default is 2.
- Returns:
DESCRIPTION.
- Return type:
TYPE
- geodezyx.reffram.geometry.helmert_trans_legacy(Xa, params='itrf2008_2_etrf2000', invert=True, workepoc=2009.0)
This function is highly unstable and should be used only for legacy reason
- Parameters:
Xa (TYPE) – DESCRIPTION.
params (TYPE, optional) – DESCRIPTION. The default is ‘itrf2008_2_etrf2000’.
invert (TYPE, optional) – DESCRIPTION. The default is True.
workepoc (TYPE, optional) – DESCRIPTION. The default is 2009..
- Returns:
Xb – DESCRIPTION.
- Return type:
TYPE
Notes
NB 1 : http://etrs89.ensg.ign.fr/memo-V8.pdf NB 2 : Transformation inverse : * -1 pour les paramètres cf https://en.wikipedia.org/wiki/Helmert_transformation
optimisé pour RGF93 => ITRF2008 (d’ou le invert = True & workepoc = 2009.)
NB3 : Attention losque l’on compare avec le convertisseur EUREF avec une vitesse parce que elle aussi est modifiée dans la conversion ETRS => ITRS. Conclusion : le faire en 2 étapes: ETRS => ITRS dans la même epoc ITRS epoc 1 => ITRS epoc 2
NE MARCHE PAS PARCE BESOIN DU RATE(TAUX) DES PARAMS D’HELMERT !!!!!! (160923)
if manual then params is a tuple params = (t1,t2,t3,dab,r1,r2,r3)
- class geodezyx.reffram.geometry.interp1d_ang(T, A, angtype='deg', kind='linear', bounds_error=False)
Bases:
object
- geodezyx.reffram.geometry.itrf_helmert_get_parameters(TRF_Input_name, TRF_Ext_name, verbose=False, convert=True)
Get the Helmert parameters for a transformation I/ETRFxx <=> I/ETRFxx
- Parameters:
TRF_Input_name (str) – Name of the initial Reference Frame. Can be ITRFyyyy, ITRFyy, ETRFyyyy, ETRFyy where yyyy or yy is the TRF realization
TRF_Ext_name (str) – Name of the wished Reference Frame. Can be ITRFyyyy, ITRFyy, ETRFyyyy, ETRFyy where yyyy or yy is the TRF realization
verbose (bool, optional) – Print the parameters. The default is True.
convert (bool, optional) – Gives directly useful values in the good units for a proper transformation
- Returns:
T (3-Array) – Translation.
Tdot (3-Array) – Translation rate.
D (float) – Scale factor.
Ddot (float) – Scale factor rate.
R (3-Array) – Rotation.
Rdot (3-Array) – Rotation rte.
epoch_ref_Xe (float) – Reference epoch of TRF_Ext_name.
Notes
Based on the values of http://etrs89.ensg.ign.fr/pub/EUREF-TN-1.pdf
Warning: a validation unsing the EUREF converter is highly recomended https://www.epncb.oma.be/_productsservices/coord_trans/
- geodezyx.reffram.geometry.itrf_helmert_trans(Xi, epoch_Xi, T, Tdot, D, Ddot, R, Rdot, epoch_ref_Xe)
Do the Helmert transformation I/ETRFxx <=> I/ETRFxx
The seven Helmert Parameters are generated with
itrf_helmert_get_parameters
- Parameters:
Xi (3-Array, Nx3-Array OR 6-Array, Nx6-Array) – Points in the initial Reference Frame. 3 columns if positions only 6 columns if positions+velocities
epoch_Xi (float) – epoch of the initial Reference Frame points (decimal year).
T (3-Array) – Translation parameter.
Tdot (3-Array) – Translation rate parameter.
D (float) – Scale factor parameter.
Ddot (float) – Scale factor rate parameter.
R (3-Array) – Rotation parameter.
Rdot (3-Array) – Rotation rate parameter.
epoch_ref_Xe (float) – Reference epoch of the destination TRF.
- Returns:
Xe – Points in the destination Reference Frame.
- Return type:
3-Array or Nx3-Array
Notes
** Warning ** This function does not the velocity shift from one epoch to another
The output coordinates will be provided at the same epoch as the input ones
Use
itrf_speed_calc
function to perform this potential step before (with the initial velocities) or after (with the transformed velocities) calling the present function.Based on the theory and values of EUREF Technical Note 1: Relationship and Transformation between the International and the European Terrestrial Reference Systems, Z. Altamimi, 2018 http://etrs89.ensg.ign.fr/pub/EUREF-TN-1.pdf
We recommend to confirm the values with the official EUREF converter http://www.epncb.oma.be/_productsservices/coord_trans/index.php
** Notes for the French geodetic system ** By definition, since the 2021/01/04 RGF93(v2b) = ETRF2000@2019.0
References
https://geodesie.ign.fr/index.php?page=rgf93
https://geodesie.ign.fr/contenu/fichiers/RGF93v2b-RAF18b.pdf
https://geodesie.ign.fr/contenu/fichiers/rgf93v2b_information_cnig.pdf
- geodezyx.reffram.geometry.itrf_psd_fundamuntal_formula(t, A_l, t_l, tau_l, A_e, t_e, tau_e)
Get the Post Seismic Deformation of a station
- Parameters:
t (float) – epoch in decimal years.
A_l (float) – Amplitude of the logarithmic term.
t_l (float) – Earthquake time(date) corresponding to logarithmic term
tau_l (float) – Relaxation time of the logarithmic term.
A_e (float) – Amplitude of the exponential term.
t_e (float) – Earthquake time(date) corresponding to the exponential term.
tau_e (float) – Relaxation time of the exponential term..
- Returns:
dL – the total sum of PSD corrections.
- Return type:
TYPE
References
http://itrf.ensg.ign.fr/ITRF_solutions/2014/doc/ITRF2014-PSD-model-eqs-IGN.pdf
- geodezyx.reffram.geometry.itrf_speed_calc(x0, y0, z0, t0, vx, vy, vz, t)
Translate a position to a given epoch using point velocity
- Parameters:
x0 (float) – coordinates at the reference epoch (m).
y0 (float) – coordinates at the reference epoch (m).
z0 (float) – coordinates at the reference epoch (m).
t0 (float) – reference epoch (decimal year).
vx (float) – speed of the point (m/yr).
vy (float) – speed of the point (m/yr).
vz (float) – speed of the point (m/yr).
t (float) – output epoch.
- Returns:
xout,yout,zout – coordinates of the point @ the ref. epoch (m).
- Return type:
floats
- geodezyx.reffram.geometry.mat_poids(Sinp, Ninp, fuvinp=1)
discontinued
- geodezyx.reffram.geometry.points_circle_border(Npts, r, r_sigma, az_type_normal=True, main_dir=3.14159, dir_range=3.14159, seed=None)
- geodezyx.reffram.geometry.project_point_on_plan(N, M, A)
Project a point on a Plan
- Parameters:
N (3-Array (Vector)) – Vector describing the plan.
M (3-Array (Vector)) – Point of the plan.
A (3-Array (Vector)) – Point we want to project.
- Returns:
P – Projected point.
- Return type:
3-Array (Vector)
Note
It can be an ambiguity on the sign
References
https://fr.wikipedia.org/wiki/Distance_d%27un_point_%C3%A0_un_plan https://www.mathematex.fr/viewtopic.php?t=923
- geodezyx.reffram.geometry.randn_bool(N, true_ratio=0.5, RandGene=None)
- geodezyx.reffram.geometry.random_walk_in_a_circle(x0, y0, xc, yc, R, N, step_size, param=1, polar=True, uniform_or_normal='n', rand_seed=-1)
generate random walk in a circle
- Parameters:
x0 (float) – x start of the random walk.
y0 (float) – y start of the random walk.
xc (float) – x center of the circle.
yc (float) – y center of the circle.
R (float) – radius of the circle.
N (int) – number of points.
step_size (float) – size of the step.
param (int, optional) – control the random. The default is 1.
polar (TYPE, optional) – DESCRIPTION. The default is True.
uniform_or_normal ('u' or 'n', optional) – uniform of normal walk. The default is ‘n’.
rand_seed (int, optional) – control the random. The default is -1.
- Returns:
X (np.array) – X-coordinates of the random walk.
Y (np.array) – Y-coordinates of the random walk.
Xcircle (np.array) – X-coordinates of the circle.
Ycircle (np.array) – Y-coordinates of the circle.
Exemple
——-
for un in (‘u’,’n’)
for pol in range(2) – X,Y , Xcircle , Ycircle = random_walk_in_a_circle(10,10,0,0,50,10000,polar = pol,uniform_or_normal=un)
plt.figure() plt.plot(Xcircle,Ycircle) plt.plot(X,Y) plt.axis(‘equal’) plt.suptitle(un + str(pol))
- geodezyx.reffram.geometry.randomwalk_normal(N=100, d=2, moy=0, sigma=1)
d = dimension
- geodezyx.reffram.geometry.randomwalk_uniform(N=100, d=2, bound=0.5)
d = dimension bound = contraint of the random walk
- geodezyx.reffram.geometry.rotate_points(alphal, betal, gammal, pointlin, Rtype='R1', xyzreftuple=([1, 0, 0], [0, 1, 0], [0, 0, 1]), angtype='deg', fullout=False)
- R1 = Rz(g) * Ry(b) * Rx(a)
si les RPY sont donnés dans le NED alors les positions résultantes sont dans le NED
- R2 = matrice RPY2ENU
si les RPY sont donnés dans le NED alors les résultantes sont DANS LE ENU pas besoin de rotation NED2ENU
Grewal et al. 2007
- Entrée :
Angles n = A liste de listes de P * [ points ]
- Sortie :
liste de listes [ [ xA ] [ xA ] … xP [ xA ] ]
- geodezyx.reffram.geometry.rotmat2(theta, angtype='deg')
- geodezyx.reffram.geometry.rotmat3(alpha, beta, gamma, xyzreftuple=([1, 0, 0], [0, 1, 0], [0, 0, 1]), angtype='deg')
- geodezyx.reffram.geometry.savage_buford_formula(Vs, X, d)
X : distance à la faille , un iterable pour toutes le profil, un nombre pour la longeur max
d : profondeur de la faille retourne X , et Vdeform(X)
X et d doivent être dans la même unité, Vs pas forcément
- geodezyx.reffram.geometry.semi_major_axis_from_mean_motion(n)
source : https://space.stackexchange.com/questions/18289/how-to-get-semi-major-axis-from-tle
- geodezyx.reffram.geometry.unwrap180(anglist, angtype='deg')
- geodezyx.reffram.geometry.wrap360(anglist, angtype='deg')
- geodezyx.reffram.geometry.wrapTo180(lonin)
wrapTo180 Wrap angle in degrees to [-180 180]
lonWrapped = wrapTo180(LON) wraps angles in LON, in degrees, to the interval [-180 180] such that 180 maps to 180 and -180 maps to -180. (In general, odd, positive multiples of 180 map to 180 and odd, negative multiples of 180 map to -180.)
See also wrapTo360, wrapTo2Pi, wrapToPi.
- geodezyx.reffram.geometry.wrapTo2Pi(lon)
wrapTo2Pi Wrap angle in radians to [0 2*pi]
lambdaWrapped = wrapTo2Pi(LAMBDA) wraps angles in LAMBDA, in radians, to the interval [0 2*pi] such that zero maps to zero and 2*pi maps to 2*pi. (In general, positive multiples of 2*pi map to 2*pi and negative multiples of 2*pi map to zero.)
See also wrapToPi, wrapTo180, wrapTo360.
- geodezyx.reffram.geometry.wrapTo360(lonin)
wrapTo360 Wrap angle in degrees to [0 360]
lonWrapped = wrapTo360(LON) wraps angles in LON, in degrees, to the interval [0 360] such that zero maps to zero and 360 maps to 360. (In general, positive multiples of 360 map to 360 and negative multiples of 360 map to zero.)
See also wrapTo180, wrapToPi, wrapTo2Pi.
- geodezyx.reffram.geometry.wrapToPi(lon)
wrapToPi Wrap angle in radians to [-pi pi]
lambdaWrapped = wrapToPi(LAMBDA) wraps angles in LAMBDA, in radians, to the interval [-pi pi] such that pi maps to pi and -pi maps to -pi. (In general, odd, positive multiples of pi map to pi and odd, negative multiples of pi map to -pi.)
See also wrapTo2Pi, wrapTo180, wrapTo360.
geodezyx.reffram.gnss_products module
@author: psakic
This sub-module of geodezyx.reffram contains functions for operations related to GNSS-products
it can be imported directly with: from geodezyx import reffram
The GeodeZYX Toolbox is a software for simple but useful functions for Geodesy and Geophysics under the GNU LGPL v3 License
Copyright (C) 2019 Pierre Sakic et al. (IPGP, sakic@ipgp.fr) GitHub repository : https://github.com/GeodeZYX/geodezyx-toolbox
- geodezyx.reffram.gnss_products.ClkDF_common_epoch_finder(ClkDFa_in, ClkDFb_in, return_index=False, supplementary_sort=False, order=['name', 'epoch'])
Find common sats/station and epochs in to Clock DF, and output the corresponding Clock DFs
Is an adapted version of OrbDF_common_epoch_finder
- geodezyx.reffram.gnss_products.ClkDF_common_epoch_finder_multi(ClkDF_list_in, return_index=False, supplementary_sort=False, order=['name', 'epoch'])
Find common sats/station and epochs in to Clock DF, and output the corresponding Clock DFs
Is is the multi version of ClkDF_common_epoch_finder
- geodezyx.reffram.gnss_products.ClkDF_filter(ClkDF_in, typ=('AS', 'AR'), name=None, ac=None, epoch_strt=datetime.datetime(1980, 1, 1, 0, 0), epoch_end=datetime.datetime(2099, 1, 1, 0, 0), name_regex=False)
Filter the content of a Clock DataFrame
- Parameters:
ClkDF_in (DataFrame) – Input Clock DataFrame (a concatenation of DF generated by files_rw.read_clk.
typ (iterable of str, optional) – List of the types of clocks: AS (satellite) or AR (receiver). The default is (“AS”,”AR”).
name (iterable of str, optional) – List of wished satellites/stations. Can be a regex (see also name_regex) The default is None.
ac (iterable of str, optional) – List of wished ACs. The default is None.
epoch_strt (datetime, optional) – Start epoch. The default is dt.datetime(1980,1,1).
epoch_end (datetime, optional) – End epoch (not included). The default is dt.datetime(2099,1,1).
name_regex (bool, optional) – the given names as ‘name’ arguments are regular expressions Some useful regex are given bellow The default is False
- Returns:
Output Clock DataFrame.
- Return type:
Clock DataFrame
Notes
‘^E[0-9]{2}’: Galileo Satellites ‘^G[0-9]{2}’: GPS Satellites
- geodezyx.reffram.gnss_products.ClkDF_filter2(ClkDF_in, typ=('AS', 'AR'), name=None, ac=None, epoch_strt=datetime.datetime(1980, 1, 1, 0, 0), epoch_end=datetime.datetime(2099, 1, 1, 0, 0), name_regex=False)
attempt for a faster version of ClkDF_filter, but the original is faster
- geodezyx.reffram.gnss_products.ClkDF_reg_2_multidx(ClkDFin, index_order=['name', 'epoch'])
From an regular Clock DF generated by read_clk(), set some columns (typically [“name”,”epoch”]) as indexes The outputed DF is then a Multi-index DF
It an adapted version of OrbDF_reg_2_multidx
- geodezyx.reffram.gnss_products.DFOrb_velocity_calc(DFOrb_in, drop_nan=False)
Compute the velocity of satellites from a DFOrb dataframe (differentiate the position)
- Parameters:
DFOrb_in (Pandas DataFrame) – an Orbit DataFrame.
drop_nan (bool, optional) – Remove the nan values (the first ones, since it is a numerical differentiation). The default is False. It is recommended to keep the NaN values, thus the output will keep same size and index as input
- Returns:
df_vel – an Orbit DataFrame with velocities (vx, vy ,vz columns).
- Return type:
Pandas DataFrame
- geodezyx.reffram.gnss_products.OrbDF_common_epoch_finder(OrbDFa_in, OrbDFb_in, return_index=False, supplementary_sort=False, order=['prn', 'epoch'], skip_reg2multidx_OrbDFa=False, skip_reg2multidx_OrbDFb=False)
This function finds common satellites and epochs in two Orbit DataFrames and outputs the corresponding Orbit DataFrames.
- Parameters:
OrbDFa_in (DataFrame) – The first input Orbit DataFrame.
OrbDFb_in (DataFrame) – The second input Orbit DataFrame.
return_index (bool, optional) – If True, the function also returns the common index. Default is False.
supplementary_sort (bool, optional) – If True, an additional sort is performed. This is useful for multi GNSS where the output DataFrame may not be well sorted. Default is False.
order (list of str, optional) – The order of the index for the multi-index DataFrame. Default is [“prn”,”epoch”].
skip_reg2multidx_OrbDFa (bool, optional) – If True, skips the conversion of the first input DataFrame to a multi-index DataFrame. Default is False. The inputs are assumed to be already in multi-index format to optimize execution speed. (For advanced use only)
skip_reg2multidx_OrbDFb (bool, optional) – If True, skips the conversion of the second input DataFrame to a multi-index DataFrame. Default is False. The inputs are assumed to be already in multi-index format to optimize execution speed. (For advanced use only)
- Returns:
OrbDFa_out (DataFrame) – The first output Orbit DataFrame with common satellites and epochs.
OrbDFb_out (DataFrame) – The second output Orbit DataFrame with common satellites and epochs.
Iinter (Index, optional) – The common index. Only returned if return_index is True.
Note
designed for orbits/sp3 first with sat and epoch as order parmeter, but can be used also for instance for snx files with STAT and epoch as order parmeter
- geodezyx.reffram.gnss_products.OrbDF_const_sv_columns_maker(OrbDFin, inplace=True)
(re)generate the const and sv columns from the sat one
- geodezyx.reffram.gnss_products.OrbDF_crf2trf(DForb_inp, DF_EOP_inp, time_scale_inp='gps', inv_trf2crf=False)
Convert an Orbit DataFrame from Celetrial Reference Frame to Terrestrial Reference Frame (.
Requires EOP to work. Cf. note below.
- Parameters:
DForb_inp (DataFrame) – Input Orbit DataFrame in Celetrial Reference Frame.
DF_EOP_inp (DataFrame) – EOP DataFrame (C04 format).
time_scale_inp (str, optional) – The time scale used in. manage ‘utc’, ‘tai’ and ‘gps’. The default is “gps”.
inv_trf2crf (bool, optional) – Provide the inverse transformation TRF => CRF. The default is False.
- Returns:
DForb_out – Output Orbit DataFrame in Terrestrial Reference Frame. (or Celestrial if inv_trf2crf is True)
- Return type:
DataFrame
Note
The EOP can be obtained from the IERS C04 products. e.g. https://datacenter.iers.org/data/latestVersion/224_EOP_C04_14.62-NOW.IAU2000A224.txt To get them as a Compatible DataFrame, use the function files_rw.read_eop_C04()
- geodezyx.reffram.gnss_products.OrbDF_lagrange_interpolate(DForb_in, Titrp, n=10, append_to_input_DF=False, plot=False)
High level function to interpolate an orbit DataFrame
- Parameters:
DForb_in (DataFrame) – an Orbit DataFrame.
Titrp (iterable of datetime) – Epochs of the wished points.
n (int, optional) – degree of the polynom. Better if even. The default is 10.
append_to_input_DF (bool, optional) – append the interpolated DF to the input DF. The default is False.
plot (bool, optional) – Plot the values. For debug only. The default is False.
- Returns:
DForb_out (DataFrame) – Interpolated orbits.
Tips
—-
Use conv.dt_range to generate the wished epochs range
- geodezyx.reffram.gnss_products.OrbDF_multidx_2_reg(OrbDFin, index_order=['prn', 'epoch'])
Convert a Multi-index formatted OrbitDF to his original form
- geodezyx.reffram.gnss_products.OrbDF_reg_2_multidx(OrbDFin, index_order=['prn', 'epoch'])
From an regular Orbit DF generated by read_sp3(), set some columns (typically [“prn”,”epoch”]) as indexes The outputed DF is then a Multi-index DF
- geodezyx.reffram.gnss_products.beta_angle_calc(DFOrb_in, calc_beta_sun_ra_dec=True, calc_beta_sun_eclip_long=True, beta_rad2deg=True)
Compute beta angle for GNSS satellite’s orbits stored in an orbit DataFrame
- Parameters:
DFOrb_in (Pandas DataFrame) – an Orbit DataFrame (ECEF frame).
calc_beta_sun_ra_dec (bool, optional) – compute beta angle with Sun’s right ascension and declination. The default is True.
calc_beta_sun_eclip_long (bool, optional) – compute beta angle with Sun’s Ecliptic longitude. The default is True.
beta_rad2deg (bool, optional) – convert beta angle to degrees. The default is True.
- Returns:
df_out (Pandas DataFrame) – DFOrb_in with a new ‘beta’ column.
df_wrk (Pandas DataFrame) – intermediate values dataframe for debug. here, coordinates are in ECI frame
- geodezyx.reffram.gnss_products.beta_sun_eclip_long(sun_ecl_long, sat_o_lan, sat_i, earth_i)
Compute beta angle based on Sun’s Ecliptic longitude Angles are in radians
- Parameters:
sun_ecl_long (float) – Sun’s Ecliptic longitude.
sat_i (float) – Satellite’s inclination.
sat_o_lan (float) – Satellite’s Longitude of the ascending node.
earth_i (float) – Earth’s inclination.
- Returns:
beta (float) – beta angle (in radians).
Source
——
Calculation of the Eclipse Factor for Elliptical Satellite Orbits NASA (1962)
https (//ntrs.nasa.gov/citations/19630000622)
Computation of Eclipse Time for Low-Earth Orbiting Small Satellites Sumanth R M (2019)
https (//commons.erau.edu/ijaaa/vol6/iss5/15/)
Simpler formula – https://en.wikipedia.org/wiki/Beta_angle
Note
you can use
pyorbital.astronomy.sun_ecliptic_longitude()
to computesun_ecl_long
If so, Sun’s right ascension and declination computation polynoms are based on: Astronomical algorithms, Jean Meeus (1st edition, 1991)
- geodezyx.reffram.gnss_products.beta_sun_ra_dec(sun_dec, sun_ra, sat_i, sat_o_lan)
Compute beta angle based on Sun’s right ascension and declination Angles are in radians
- Parameters:
sun_dec (float) – Sun’s declination.
sun_ra (float) – Sun’s right ascension.
sat_i (float) – Satellite’s inclination.
sat_o_lan (float) – Satellite’s Longitude of the ascending node.
- Returns:
beta (float) – beta angle (in radians).
Source
——
`Satellites (de Kepler au GPS`, Michel Capderou (2012))
simpler formula here – https://www.fxsolver.com/browse/formulas/Beta+Angle
Note
you can use
pyorbital.astronomy.sun_ra_dec()
to computesun_dec
andsun_ra
If so, Sun’s right ascension and declination computation polynoms are based on: Astronomical algorithms, Jean Meeus (1st edition, 1991)
- geodezyx.reffram.gnss_products.compar_clk_plot(Diff_sat_all_df_in, save_plot=False, save_plot_dir='', save_plot_name='auto', save_plot_name_suffix=None, save_plot_ext=('.pdf', '.png', '.svg'), yaxis_limit=None, yaxis_label_unit='psec', col_name='name', bias_Col_name='bias')
General description
- Parameters:
Diff_sat_all_df_in (DataFrame) – a DataFrame produced by compar_clk
yaxis_limit (3-tuple iterable or 2-element tuple) – force the y axis limits. must look like (ymin,ymax) to set all th axis at the same limits
col_name (Normally the name of the column with the sat names) –
bias_Col_name (The column with the clk values) –
Default ('bias') –
- Returns:
the Figure and the 3 Axes if no save is asked
export path (str) if save is asked
but plot a plot anyway
- geodezyx.reffram.gnss_products.compar_clock(DFclk_inp_1, DFclk_inp_2, col_name='name', bias_Col_name='bias')
Compares 2 GNSS clock bias DataFrames (from .clk), to a statistics table (with compar_clock_table)
- Parameters:
DFclk_inp_2 (DFclk_inp_1 &) – Clock DataFrame provided by files_rw.read_clk()
- Returns:
DFclk_diff – Clock bias difference DataFrame
- Return type:
DataFrame
- geodezyx.reffram.gnss_products.compar_clock_table(DFclk_diff_in, col_name='name', bias_Col_name='bias_diff')
Generate a table with statistical indicators for a clock comparison (RMS mean, standard dev, …)
- Parameters:
DFclk_diff_in (DataFrame) – Clock bias difference DataFrame (from compar_clock)
- Returns:
DFcompar_out – Statistical results of the comparison.
- Return type:
DataFrame
- geodezyx.reffram.gnss_products.compar_orbit(Data_inp_1, Data_inp_2, step_data=900, sats_used_list=['G'], name1='', name2='', use_name_1_2_for_table_name=False, RTNoutput=True, convert_ECEF_ECI=True, clean_null_values=True, conv_coef=1000, return_satNull=False)
Compares 2 GNSS orbits files (SP3), and gives a summary plot and a statistics table
- Parameters:
Data_inp_2 (Data_inp_1 &) – contains the orbits or path (string) to the sp3
step_data (int) – per default data sampling
sats_used_list (list of str) – used constellation or satellite : G E R C … E01 , G02 … Individuals satellites are prioritary on whole constellations e.g. [‘G’,”E04”]
RTNoutput (bool) – select output, Radial Transverse Normal or XYZ
convert_ECEF_ECI (bool) – convert sp3 ECEF => ECI (Terrestrial => Celestrial) must be True in operational to avoid artifacts.
name2 (name1 &) – optional custom names for the 2 orbits
use_name_1_2_for_table_name (bool) – False : use name 1 and 2 for table name, use datafile instead
clean_null_values (bool or str) – if True or “all” remove sat position in all X,Y,Z values are null (0.000000) if “any”, remove sat position if X or Y or Z is null if False, keep everything
conv_coef (int) – conversion coefficient, km to m 10**3, km to mm 10**6
- Returns:
Diff_sat_all (Pandas DataFrame)
contains differences b/w Data_inp_1 & Data_inp_2
in Radial Transverse Normal OR XYZ frame –
- Attributes of Diff_sat_all :
Diff_sat_all.name : title of the table
Note
clean_null_values if useful (and necessary) only if convert_ECEF_ECI = False if convert_ECEF_ECI = True, the cleaning will be done by a side effect trick : the convertion ECEF => ECI will generate NaN for a zero-valued position But, nevertheless, activating clean_null_values = True is better This Note is in fact usefull if you want to see bad positions on a plot => Then convert_ECEF_ECI = False and clean_null_values = False
References
“Coordinate Systems”, ASEN 3200 1/24/06 George H. Born
- geodezyx.reffram.gnss_products.compar_orbit_frontend(DataDF1, DataDF2, ac1, ac2, sats_used_list=['G'])
- geodezyx.reffram.gnss_products.compar_orbit_plot(Diff_sat_all_df_in, save_plot=False, save_plot_dir='', save_plot_name='auto', save_plot_name_suffix=None, save_plot_ext=('.pdf', '.png', '.svg'), yaxis_limit=None, yaxis_label_unit='m')
General description
- Parameters:
Diff_sat_all_df_in (DataFrame) – a DataFrame produced by compar_orbit
yaxis_limit (3-tuple iterable or 2-element tuple) – force the y axis limits. must look like [(ymin_r,ymax_r),(ymin_t,ymax_t),(ymin_n,ymax_n)] to control all the axis independely OR (ymin,ymax) to set all th axis at the same limits
- Returns:
the Figure and the 3 Axes if no save is asked
export path (str) if save is asked
but plot a plot anyway
- geodezyx.reffram.gnss_products.compar_orbit_table(Diff_sat_all_df_in, RMS_style='natural', light_tab=False)
Generate a table with statistical indicators for an orbit comparison (RMS mean, standard dev, …)
- Parameters:
Diff_sat_all_df_in (Pandas DataFrame) – a DataFrame produced by compar_orbit
RMS_style (str) – ‘natural’: use the natural definition of the RMS ‘GRGS’: RMS calc based on the GRGS definition of the RMS (OV help) is actually the standard deviation ‘kouba’: RMS as defined in Kouba et al. 1994, p75 using the degree of freedom (3*Nobs - 7)
light_tab (bool) – produce a table with only RMS, with min/max/arithmetic instead
- Returns:
Compar_tab_out – Statistical results of the comparison
- Return type:
DataFrame
Note
you can pretty print the output DataFrame using tabular module here is a template:
>>> from tabulate import tabulate >>> print(tabulate(ComparTable,headers="keys",floatfmt=".4f"))
- geodezyx.reffram.gnss_products.compar_sinex(snx1, snx2, stat_select=None, invert_select=False, out_means_summary=True, out_meta=True, out_dataframe=True, manu_wwwwd=None)
- geodezyx.reffram.gnss_products.eop_interpotate(DF_EOP, Epochs_intrp, eop_params=['x', 'y'])
Interopolate the EOP provided in a C04-like DataFrame
- Parameters:
DF_EOP (DataFrame) – Input EOP DataFrame (C04 format). Can be generated by files_rw.read_eop_C04
Epochs_intrp (datetime of list of datetimes) – Wished epochs for the interpolation.
eop_params (list of str, optional) – Wished EOP parameter to be interpolated. The default is [“x”,”y”].
- Returns:
OUT – Interpolated parameters. Series if onely one epoch is provided, DF_EOP elsewere
- Return type:
DataFrame or Series
- geodezyx.reffram.gnss_products.get_block_svn(sat_in, svn_prn_equiv_DF)
Get the equivalence SVN block type
- Parameters:
sat_in (str) – Satellite SVN.
svn_prn_equiv_DF (DataFrame) – Equivalence table generated by svn_prn_equiv_DF.
- Return type:
str with the block name
- geodezyx.reffram.gnss_products.stats_slr(DFin, grpby_keys=['sat'], threshold=0.5)
computes statistics for SLR Residuals
- Parameters:
DFin (Pandas DataFrame) – Input residual Dataframe from read_pdm_res_slr.
grpby_keys (list of str, optional) – The default is [‘sat’]. per day, per solution, per satellite: [‘day’,’sol’,’sat’] per day, per solution, per station: [‘day’,’sol’,’sta’] per day, per solution, per satellite, per station: [‘day’,’sol’,’sta’,’sat’]
threshold (float) – apply a Threshold
- Returns:
DD – return the mean, the rms and the std.
- Return type:
Output statistics DataFrame
- geodezyx.reffram.gnss_products.svn_prn_equiv(sat_in, date_in, svn_prn_equiv_DF, mode='svn2prn', full_output=False)
Get the equivalence SVN <> PRN for a given epoch
- Parameters:
sat_in (str) – Satellite “ID”, SVN or PRN.
date_in (datetime) – wished epoch.
svn_prn_equiv_DF (DataFrame) – Equivalence table generated by svn_prn_equiv_DF.
mode (str, optional) – prn2svn: PRN > SVN svn2prn: SVN > PRN. The default is “svn2prn”.
full_output (bool, optional) – get the complete Equivalence table row. The default is False.
- Return type:
str or DataFrame
- geodezyx.reffram.gnss_products.svn_prn_equiv_DF(path_meta_snx)
generate a SVN <> PRN equivalent DataFrame
- Parameters:
path_meta_snx (str) – path of the MGEX metadata sinex. last version avaiable here http://mgex.igs.org/IGS_MGEX_Metadata.php
- Returns:
DFfin – SVN <> PRN equivalent DataFrame.
- Return type:
Pandas DataFrame
geodezyx.reffram.kepler_gzyx module
@author: psakic
- This sub-module of geodezyx.reffram contains functions related to
Kepler’s element operations.
it can be imported directly with: from geodezyx import reffram
The GeodeZYX Toolbox is a software for simple but useful functions for Geodesy and Geophysics under the GNU LGPL v3 License
Copyright (C) 2019 Pierre Sakic et al. (IPGP, sakic@ipgp.fr) GitHub repository : https://github.com/GeodeZYX/geodezyx-toolbox
- geodezyx.reffram.kepler_gzyx.ECI_2_kepler_elts(pos, vel, rad2deg=True, mu=398600441880000.0)
Convert ECI coordinates > Kepler’s elements
- Parameters:
pos (array) – Position in m. Can be a 3-element (x,y,z) array (simple point) or a (N,3)-shaped numpy array
vel (array) – Velocity in m/s. Can be a 3-element (vx,vy,vz) array (simple point) or a (N,3)-shaped numpy array
rad2deg (bool, optional) – convert Keplerian’s angles (anomalies) to deg. The default is True.
mu (float, optional) – Standard gravitational parameter. The default is 3.9860044188e14.
- Returns:
a,ecc,i,o_peri,o_lan,m – Kepler’s elements:
a : semi-major axis
ecc : orbit eccentricity
i : orbit inclination
o_peri : argument of periapsis ω
o_lan : longitude of the ascending node Ω
m : mean anomaly m
- Return type:
floats
Note
Be sure your input is in ECI (SP3 are in ECEF for instance).
If neeeded, use
conv.ECEF2ECI()
(gross results) orreffram.OrbDF_crf2trf(inv_trf2crf=True)
(precise results)You can get the velocity with
reffram.DFOrb_velocity_calc()
- geodezyx.reffram.kepler_gzyx.ECI_2_kepler_elts_mono(P, V, rad2deg=True, mu=398600441880000.0)
Kept for debug purposes, use
ECI_2_kepler_elts(pos, vel)
insteadConvert ECI coordinates > Kepler’s elements
- Parameters:
P (array) – Position.
V (array) – Velocity.
rad2deg (bool, optional) – canvert Keplerian’s angles (anomalies) to deg. The default is True.
mu (float, optional) – Standard gravitational parameter. The default is 3.9860044188e14.
- Returns:
a,ecc,i,omega_periarg,omega_LAN,m – Kepler’s elements.
- Return type:
floats
- geodezyx.reffram.kepler_gzyx.extrapolate_orbit_kepler(P, V, t, t0, mu=398600441880000.0)
- geodezyx.reffram.kepler_gzyx.extrapolate_sp3_DataFrame(DFsp3, step=900, n_step=9, backward=True, forward=True, until_backward=None, until_forward=None, return_all=True)
Extrapolate the positions in a SP3 based on the first/last epochs’ position
- Parameters:
DFsp3 (DataFrame) – Input Orbit DataFrame (i.e. generated by files_rw.read_sp3).
step (int, optional) – step between two epochs. The default is 900.
n_step (int, optional) – number of epochs to interpolate. The default is 9.
forward (backward and) – extrapolate for the day before/after. The default is True.
return_all (bool, optional) – if True, returns the input DF enhanced with the extrapolated values if False, returns only the extrapolated values. The default is True.
until_backward (until_backward &) – epoch until then the extrapolation has to be done. Override n_step
- Returns:
DForb_out – Orbit DataFrame with extrapolated values (see return_all).
- Return type:
DataFrame
geodezyx.reffram.quaternions module
@author: psakic
This sub-module of geodezyx.reffram contains functions for quaternion operations.
it can be imported directly with: from geodezyx import reffram
The GeodeZYX Toolbox is a software for simple but useful functions for Geodesy and Geophysics under the GNU LGPL v3 License
Copyright (C) 2019 Pierre Sakic et al. (IPGP, sakic@ipgp.fr) GitHub repository : https://github.com/GeodeZYX/geodezyx-toolbox
- geodezyx.reffram.quaternions.interp_quat(Tlis, Quatlis, t)
Interpolate quaterinon
- Parameters:
Tlis (float) – list of time, in POSIX time
Quatlis – list of quaternion (produced by quaterion fct)
t – the instant for the interpolation
- Returns:
an interpoled Quaternion
- geodezyx.reffram.quaternions.interp_quat_multi(Tlis, Quatlis, tlis)
- geodezyx.reffram.quaternions.quatern_multi(rlist, plist, hlist, angtype='rad')
wrapper of quaternion, with lists/arrays of angles
- geodezyx.reffram.quaternions.quaternion(r, p, h, angtype='rad', euler_sequence='zyx')
- Parameters:
roll (xaxis) pitch (yaxis) heading (zaxis) –
angtype – ‘rad’ or ‘deg’
- Returns:
the associate quaternion
Note
the input args are in XYZ but the rotation is in the order ZYX
Need of the cgkit quaternion module http://cgkit.sourceforge.net/doc2/quat.html
- geodezyx.reffram.quaternions.rotate_quat(ptin, quatin)
- geodezyx.reffram.quaternions.rotate_quat2(ptin, quatin)
DISCONTINUED plus long que rotate_quat : benchmarking 5.53084397316 vs 0.703444957733 pour 100000
- geodezyx.reffram.quaternions.rotate_quat3(ptin, quatin)
DISCONTINUED use the old quaterinons toolbox
- geodezyx.reffram.quaternions.rotate_quat_multi(pts_list_in, quats_list_in)
Multi version of rotate_quat
- Parameters:
pts_list_in – a list of Points object
quats_list_in – a list of quaternions
- Returns:
a list of list of rotated points