Gts methods¶
Gts format methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- read_pos(tsdir='.', tsfile=None, xyz=True, verbose=False)¶
Read GAMIT/GLOBK PBO pos file in a directory and actually loads the time series
- Parameters
tsdir – directory of pos file
tsfile – pos file to be loaded
xyz – reads xyz sx sy sz corr_xy corr_xz corr_yz columns
verbose – verbose mode
- Note
Since a pos file includes (almost) all the information, data, code, X0,Y0,Z0,t0 will be populated
- Note
If tsfile=None, then read_pos will look for a file named CODE*.pos
- write_pos(idir, add_key='', force=None, verbose=False)¶
Write a time series in GAMIT/GLOBK PBO pos format
- Parameters
idir – output directory
add_key – if not blank then the output pos file will be CODE_add_key.pos, CODE.pos otherwise.
force – set force to ‘data’ or ‘data_xyz’ to force pos to be written from .data or .data_xyz
- :note1:default behaviour (force = None)
if data and data_xyz are not None, then print them independently if there are data only, then uses X0,Y0,Z0 to write data_xyz if there are data_xyz only, recreate data and write it
- read_cats_file(idir='.', ifile=None, gmt=True, verbose=False)¶
Read cats files in a directory and actually loads the time series
- write_cats(idir, offsets_dates=None, add_key='')¶
Writes a file for a cats analysis if offsets_dates is not None then offsets are added at the beginning of the file
- force_daily(in_place=False)¶
force a time series to be daily with dates at 12:00:00 of every day
- get_unr(site, verbose=False)¶
Get a time series from http://geodesy.unr.edu/gps_timeseries/txyz/IGS14/ in PYACS
- Parameters
site – 4-letters code
verbose – verbose mode
- read_kenv(ifile, date_type='jd')¶
Read kenv file (magnet) format for time series
- read_mb_file(idir='.', ifile=None, gmt=True, verbose=False)¶
Read GAMIT/GLOBK mb_files in a directory and actually loads the time series
- read_pride(tsdir='.', tsfile=None, xyz=True, verbose=False)¶
Read PRIDE-PPPAR kinematic result file :param tsdir: directory of pride-pppar kinematic files :param tsfile: pride-pppar kinematic file to be loaded :param verbose: verbose mode :return Nothing: :note: If file=None, then read_pride will look for a files named kin_*code
- read_pride_pos(tsdir='.', tsfile=None, verbose=False)¶
Read PRIDE-PPPAR static result file
- Parameters
tsdir – directory of pride-pppar pos static files
tsfile – pride-pppar pos static file to be loaded
verbose – verbose mode
:note:If file=None, then read_pride will look for a files named pos_*code
- read_tdp(idir='.', ifile=None, gmt=True, verbose=False)¶
Read tdp (Gipsy kinematics provided by Cedric Twardzik 17/04/2018) format for time series
- to_pandas_df(data_xyz=False, uncertainty=False, round=False)¶
Converts a pyacs Gts to a pandas dataframe
- Returns
pandas DataFrame
- Note
uncertainties are not imported.
- read_track_NEU(tsdir='.', tsfile=None, leap_sec=0.0, eq_time=None, verbose=False)¶
Read a GAMIT/GLOBK Track output file generated with the option out_type NEU in this case dates are seconds by default the seconds are with respect to the first epoch of measurements If option leap_sec is provided with a value > 0.0, then GPS time is corrected for the difference between GPTS time and UTC If eq_time is provided, it is assumed to be UTC. Expected format is YYYY:MM:MD:HH:MM:SS.S
Gts primitive methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- cdata(data=False, data_xyz=False, tol=0.001, verbose=False)¶
Check data/data_xyz attributes
- Parameters
data – boolean, if True, data attribute will be checked
data_xyz – boolean, if True, data_xyz attribute will be checked
tol – tolerance in days for two dates to be considered as the same (default 0.001 of day)
verbose – boolean, verbose mode
:return : boolean, True if everything is OK, False otherwise
:note : in future, this routine should also whether .data and .data_xyz value are consistent
- copy(data=True, data_xyz=True, loutliers=True)¶
makes a (deep) copy of the time series.
By default, all attributes are also copied, including .data, .data_xyz, loutliers etc.
Default behaviour can be modified for the following attribute:
- Parameters
data – can be set to None or a 2D numpy array of shape (n,10)
data_xyz – can be set to None or a 2D numpy array of shape (n,10)
loutliers – False will not copy the loutliers atrribute
- differentiate()¶
differentiate the current time series :return: the differentiated time series as a new Gts object :note : differentiation is made on .data. .data_xyz is set to None.
- extract_periods(lperiod, in_place=False, verbose=False, no_reset=False)¶
extract periods of a Gts
- Parameters
lperiod – a list [start_date,end_date] or a list of periods e.g. periods=[[2000.1,2003.5],[2009.3,2010.8]]
in_place – if True, will make change in place, if False, return s a new time series
- Note 1
X0,Y0,Z0 attributes will be changed if necessary
- Note 2
handles both .data and .data_xyz
- exclude_periods(lperiod, in_place=False, verbose=False)¶
exclude periods of a Gts
- Parameters
lperiod – a list [start_date,end_date] or a list of periods e.g. periods=[[2000.1,2003.5],[2009.3,2010.8]]
in_place – if True, will make change in place, if False, return s a new time series
- Note 1
X0,Y0,Z0 attributes will be changed if necessary
- Note 2
handles both .data and .data_xyz
- extract_dates(dates, tol=0.05, in_place=False, verbose=True)¶
Returns a time series extracted for a given list of dates
- Parameters
dates – dates either as a list or 1D numpy array of decimal dates
tol – date tolerance in days to assert that two dates are equal (default 0.05 day)
in_place – if True, will make change in place, if False, return s a new time series
verbose – boolean, verbose mode
- substract_ts(ts, tol=0.05, verbose=True)¶
substract the ts provided as argument to the current time series
- Parameters
ts – time series to be substracted as a Gts instance
tol – date tolerance to decide whether two dates are identical in both time series. default = 1/4 day
verbose – verbose mode
:return : new Gts
- substract_ts_daily(ts, verbose=True)¶
substract the ts provided as argument to the current time series
- Parameters
ts – time series to be substracted as a Gts instance
verbose – verbose mode
:return : new Gts
- Note
this method assumes daily time series
- add_offsets_dates(offsets_dates, in_place=False)¶
add_offsets_dates to a time series if in_place = True then replace the current time series
- remove_velocity(vel_neu, in_place=False)¶
remove velocity from a time series vel_neu is a 1D array of any arbitrary length, but with the velocities (NEU) to be removed in the first 3 columns if in_place = True then replace the current time series
- set_zero_at_date(date, offset=None, in_place=False)¶
make a translation of a time series, setting to 0 at a given date if the provided date does not exist, uses the next date available
- Parameters
date – date in decimal year
offset – an offset (in mm) to be added. Could be a float, a list or 1D numpy array with 3 elements
- decimate(time_step=30.0, dates=[], method='median', verbose=False)¶
decimate a time series
- Parameters
time_step – time step in days
dates – list of dates where point are forced to be written regardless time_step
method – method used to be used to calculated the position. choose among [‘median’,’mean’,’exact’]
verbose – verbose mode
:return : new Gts
- displacement(sdate=None, edate=None, window=None, method='median', speriod=[], eperiod=[], rounding='day', verbose=True)¶
Calculates displacements between two dates or two periods
- Parameters
sdate – start date in decimal year
edate – start date in decimal year
window – time window in days for searching available dates
method – method to calculate the position. ‘median’ or ‘mean’. default is ‘median’.
speriod – period for calculating the start position
eperiod – period for calculating the end position
rounding – rounding for dates. Choose among ‘day’,’hour’,’minute’ or ‘second’. default is ‘day’.
verbose – verbose mode
- Returns
displacement as np.array([dn,de,du,sdn,sde,sdu])
- add_obs(date, NEUSNSESUCNECNUCEU, in_place=False, check=True, verbose=False)¶
Adds observation(s) as DN,DE,DU to a time series
- Parameters
date – date in decimal year. float, a list or 1D numpy array
NEUSNSESUCNECNUCEU – value to be added in the Gts, provided as a list, a 1D numpy array or a 2D numpy array. requires at least NEU: North, East, UP values optional: SN, SE, SU, CNE, CNU, CEU: standard deviations and correlation coefficient between North, East and Up components. If not provided, SN=SE=SU=0.001 (1 mm) and CNE=CNU=CEU=0
in_place – boolean, if True add_obs to the current Gts, if False, returns a new Gts
check – check time order and duplicate dates
verbose – verbose mode
- Returns
new Gts or the modified Gts if in_place
- Note
if it exists, .data_xyz will be set to None for consistency.
- add_obs_xyz(date, XYZSXSYSZCXYCXZCYZ, in_place=False, check=True, neu=True, verbose=False)¶
Adds observation(s) as XYZ to a time series
- Parameters
date – date in decimal year. float, a list or 1D numpy array
XYZSXSYSZCXYCXZCYZ – value to be added in the Gts, provided as a list, a 1D numpy array or a 2D numpy array. requires at least X,Y,Z. Optional: SX, SY, SZ, CXY, CXZ, CYZ: standard deviations and correlation coefficients. If not provided, SX=SY=SZ=0.001 (1 mm) and CXY=CXZ=CYZ=0
in_place – boolean, if True add_obs to the current Gts, if False, returns a new Gts
check – check time order , duplicate dates and re-generate NEU time series (.data)
neu – regenerate .data from the updated .data_xyz
verbose – verbose mode
:return : new Gts or the modified Gts if in_place :note 1: by default .data will be updated from .data_xyz, and X0,Y0,Z0 will be updated. :note 2:
- xyz2neu(corr=False, ref_xyz=None, verbose=False)¶
populates neu (data) using xyz (data_xyz) lon, lat and h will also be set.
- Parameters
corr – if True, then standard deviation and correlations will also be calculated
ref_xyz – [X,Y,Z] corresponding to the 0 of the local NEU frame. If not provided, the first position is used as a reference
verbose – verbose mode
- Note
this method is always in place
- neu2xyz(corr=False, verbose=False)¶
populates .data_xyz from .data requires X0,Y0,Z0 attributes to be set
- Parameters
corr – if True, then standard deviation and correlations will also be calculated
verbose – verbose mode
- reorder(verbose=False)¶
reorder data and/or data_xyz by increasing dates always in place
- Parameters
verbose – verbose mode
- extract_ndates_before_date(date, n, verbose=False)¶
Extract n values before a given date If n values are not available, returns all available values before date .data is set to None if no value at all is available
- Parameters
date – date in decimal year
n – number of observations to be extracted
- Returns
a new Gts
- extract_ndates_after_date(date, n, verbose=False)¶
Extract n values after a given date If n values are not available, returns all available values after date .data is set to None if no value at all is available
- Parameters
date – date in decimal year
n – number of observations to be extracted
- Returns
a new Gts
- extract_ndates_around_date(date, n)¶
Extract n values before and n values after a given date If n values are not available, returns all available values .data is set to None if no value at all is available
- Parameters
date – date in decimal year
n – number of observations to be extracted
- Returns
a new Gts
- get_coseismic(eq_date, window_days=5, sample_after=1, method='median', in_place=False)¶
Get coseismic displacement at a given date. Coseismic displacement is estimated as the position difference between the median of window_days before the earthquake date and the median of sample_after samples after the earthquake date.
note: only median method implemented
- correct_duplicated_dates(action='correct', tol=0.1, in_place=False, verbose=False)¶
Check or remove duplicated dates in a time series
- Parameters
action – ‘correct’ (default) or ‘check’
tol – tolerance for two dates to be considered as the same (default = 0.1 day)
in_place – boolean, if True,
verbose – verbose mode
- rotate(angle, in_place=False)¶
rotates the axis by an angle
- Parameters
angle – angle in decimal degrees clockwise
if in_place = True then replace the current time series
- insert_gts_data(gts, in_place=False, verbose=False)¶
insert data (and/or) .data_xyz of a gts into the current gts
- Parameters
gts – time series to be inserted
in_place – boolean, if True add_obs to the current Gts, if False, returns a new Gts
verbose – verbose mode
:return : new Gts or the modified Gts if in_place
- insert_ts(ts, rounding='day', data='xyz', overlap=True)¶
- Parameters
ts – Gts to be inserted
rounding – data rounding, used to decide whether an entry should be replaced. Choose among [‘second’,’minute’,’hour’,’day]
data – Gts attribute to be updated. ‘xyz’ for .data_xyz or None for .data
overlap – if True, update occurs only on dates. If False, then ts overwrites the current Gts over the ts period
- Returns
a new gts
- Note
The returned gts will have .data or .data_xyz will be set to None according to data argument
- find_large_uncertainty(sigma_thresold=10, verbose=True, lcomponent='NE')¶
Find dates with large uncertainty and flag them as outliers.
- Parameters
sigma_threshold – value (mm) for a date to be flagged.
verbose – verbose mode
lcomponent – list of components to be checked. default = ‘NE’
- split_gap(gap=10, verbose=False)¶
- Parameters
gap – gap in number of days to split the time series
verbose – verbose mode
- Returns
a list a gts split from the original
- interpolate(date='day', kind='linear', gap=10, in_place=False, verbose=False)¶
- Parameters
self – Gts instance
date – ‘day’ will perform daily interpolation, alternatively date is a 1D numpy array with either datetime or decimal year
method – scipy.interpolate.interp1d kind argument
gap – maximum gap for interpolation
in_place – boolean.
verbose – verbose mode
- Returns
Gts model methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- detrend(method='L2', in_place=False, periods=[], exclude_periods=[])¶
detrends a time series and save velocity estimates in velocity attribute
:param periods : periods used to estimate the velocity :param exclude_periods : periods to be excluded for the velocity estimate :param in_place : if True then replace the current time series :return : the detrended time series :note : outliers from Gts.outliers are omitted in the estimation and offsets given Gts.offsets_dates are estimated simultaneously
- detrend_seasonal(method='L2', in_place=False, periods=None, exclude_periods=None)¶
estimates a trend + annual + semi-annual terms in a time series and removes them velocity, annual and semi-annual attributes are saved in Gts.velocity, Gts.annual, Gts.semi_annual
:param periods : periods used for estimation :param exclude_periods : periods to be excluded from estimation :param in_place : if True then replace the current time series :return : the detrended time series :note : outliers from Gts.outliers are ommitted in the estimation and offsets given Gts.offsets_dates are estimated simultaneously
- detrend_annual(method='L2', in_place=False, periods=None, exclude_periods=None)¶
estimates a trend + annual terms in a time series and removes them velocity and annual attribute are saved in Gts.velocity & Gts.annual
:param periods : periods used for estimation :param exclude_periods : periods to be excluded from estimation :param in_place : if True then replace the current time series :return : the detrended time series :note : outliers from Gts.outliers are ommitted in the estimation and offsets given Gts.offsets_dates are estimated simultaneously
- detrend_median(delta_day=None, in_place=False, periods=[], exclude_periods=[], verbose=False, auto=False)¶
Calculates a velocity using the median of pair of displacements exactly separated by one year, inspired from MIDAS If the time series has less than a year of data, then the time series is kept untouched. :param delta_day: if None, it is one year, if 0 then it is the relax mode for campaign data,
any integer is the time delta (in days) used to compute velocity.
- Parameters
in_place – boolean, if True, in_place, if False, returns a new Gts instance (default)
periods – periods (list of lists) to be included for trend calculation
exclude_periods – periods (list of lists) to be excluded for trend calculation
verbose – verbose mode
auto – if True, then start will delta_day=None, if it fails or found less than 100 pairs then use delta_day=0, if fails then use regular detrend
- Note
returns None if time series is shorter than 1 year
- detrend_seasonal_median(wl=11, in_place=False, verbose=False)¶
Calculates a velocity using the median of pair of displacements exactly separated by one year, inspired from MIDAS and then removes repeating yearly signal If the time series has less than three years of data, then the time series is kept untouched.
- frame(frame=None, in_place=False, verbose=False)¶
Rotates a time series according to an Euler pole Returns a new Gts instance
- make_model(option='detrend', method='L2', loutlier=None, in_place=False)¶
Estimate linear model parameters using least squares input: data: Gts format option are: ‘detrend’/’detrend_annual’/’detrend_seasonal’ output: new Gts object: time series is now the residuals wrt to the model and its associated values (vel, annual, semi-annual etc)
- mmodel()¶
Generates a modeled time series from the parameters read in self
- remove_pole(pole, pole_type='euler', in_place=False, verbose=True)¶
remove velocity predicted by an Euler pole or a rotation rate vector from a time series pole is a 1D array with 3 values requires self.lon & self.lat attributes to have been filled before if in_place = True then replace the current time series
- add_vel_sigma(in_place=False, b_fn=4, verbose=True)¶
calculates realistic sigma on velocity components assuming white & flicker using eq (19) & (23) from Williams (J. of Geodesy, 2003) b_fn is the value for flicker noise, taken as 4 mm/yr^1/4 model can be detrend, detrend_annual, detrend_seasonal if in_place = True then replace the current time series
- trajectory(model_type, offset_dates=[], eq_dates=[], H_fix={}, H_constraints={}, H_bounds={}, component='NEU', verbose=False)¶
Calculates the parameters of a (non-linear) trajectory model for a Geodetic Time Series. The trajectory model is:
y(t) =
trend : trend_cst + trend * ( t - t0 ) +
annual: a_annual * cos( 2*pi + phi_annual ) +
semi-annual: a_semi_annual * cos( 2*pi + phi_semi_annual ) +
offset : Heaviside( t - t_offset_i ) * offset_i +
post-seismic_deformation as decaying log (psd_log): psd_eq_i * np.log( 1 + Heaviside( t - eq_i )/tau_i )
- Parameters
model_type – string made of the key-word the parameters to be estimated.
Key-word parameters are
‘trend’,’annual’,’semi-annual’,’seasonal’,’offset’,’psd_log’.
‘trend-seasonal-offset-psd_log’ will do the full trajectory model.
- Parameters
offset_dates – a list of offset_dates in decimal year
eq_dates – a list of earthquake dates for which post-seismic deformation (psd_log) will be estimated
H_fix – a dictionary including the name of the parameter to be hold fixed and the value.
For instance to impose the co-seismic offset (North-East-Up) and relaxation time of 100 days for the first earthquake use:
H_fix = { ‘psd_log_offset_00’:[10., 15., 0.] , ‘psd_log_tau_00’:[100., 100., 100.]}
- Parameters
H_constraints – a dictionary including the name of the parameter to be constrained.
For instance to impose a 50 days constraints around 500 days on the relaxation time of the second earthquake for all NEU components use: H_fix = { ‘psd_log_tau_01’:[[500.,50], [500.,50] , [500.,50]]}
- Parameters
H_bounds – a dictionary including the bounds.
For instance to impose a relaxation time for the third earthquake to be in the range of 2 to 3 years, for all NEU components use: H_bounds = { ‘psd_log_tau_02’:[[2*365.,3*365.], [[2*365.,3*365.] , [[2*365.,3*365.]]}
- Parameters
component – string , component for which the trajectory model will be estimated.
verbose – verbose mode
- Note
Unlike most pyacs.gts functions, trajectory returns 4 elements: the results as a dictionary, the model Gts,
the residual Gts and a Gts with model predictions at every day.
Gts filter methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- l1_trend(vlambda, gap=10, in_place=False, verbose=True, component='NEU')¶
return a piecewise linear filtered Gts
- Parameters
vlambda – weight of regularization
gap – gap in days to split the time series before applying the filter.default is 10.
in_place – if True then replace the current time series
verbose – boolean, verbose mode
component – string. Default ‘NEU’
- Returns
the filtered time series
- Note
if there are less than 4 points in a segment, return an L1 estimated trend
- minimum_component(mask_period=[], p=1, fcut=None, Q=None, in_place=False, verbose=True)¶
Minimum component filtering for Gts. Minimum component filtering is useful for determining the background component of a signal in the presence of spikes :param mask_periods: periods (list or list of lists) which should be ignored for smoothing :param p: integer (optional). polynomial degree to be used for the fit (default = 1) :param fcut: float (optional). the cutoff frequency for the low-pass filter. Default value is f_nyq / sqrt(N) :param Q: float (optional). the strength of the low-pass filter. Larger Q means a steeper cutoff. default value is 0.1 * fcut :param in_place: if True then replace the current time series :param verbose: boolean, verbose mode :return: the filtered time series :note: This code follows the procedure explained in the book “Practical Statistics for Astronomers” by Wall & Jenkins book, as well as in Wall, J, A&A 122:371, 1997
- savitzky_golay(in_place=False, verbose=True, window_length=15, polyorder=3, deriv=0, delta=1.0, mode='interp', cval=0.0)¶
returns a filtered time series using scipy.signal.savgol_filter
See documentation for the filter parameters. http://https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.savgol_filter.html#scipy.signal.savgol_filter
- Parameters
in_place – if True then replace the current time series
verbose – boolean, verbose mode
- Returns
the filtered time series
- smooth(window_len=11, window='hanning', in_place=False, verbose=False, component='NEU')¶
smooth a time series
- spline(smoothing=1, degree=5, date=None)¶
- Parameters
smoothing – Positive smoothing factor used to choose the number of knots. Number of knots will be increased
- until the smoothing condition is satisfied:
sum((w[i] * (y[i]-spl(x[i])))**2, axis=0) <= s
- Parameters
degree – Degree of the smoothing spline. Must be <= 5. Default is k=3, a cubic spline.
date – 1D array of interpolation dates in decimal year, or ‘day’ for every day. defualt None will interpolate
at data date only. :return: new gts instance
- vondrak(fc, in_place=False, verbose=True, component='NEU')¶
returned a filtered Gts using a Vondrak filter
- Parameters
fc – cutoff frequence in cycle per year
in_place – if True then replace the current time series
verbose – boolean, verbose mode
- Returns
the filtered time series
- wiener(in_place=False, verbose=True, my_size=15, noise=None)¶
returns a filtered time series using scipy.signal.wiener
See documentation for the filter parameters.
- Parameters
in_place – if True then replace the current time series
verbose – boolean, verbose mode
- Returns
the filtered time series
Gts outliers methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- remove_outliers(periods=None, in_place=False)¶
removes outliers provided in self.outliers return a new Gts without the outliers if in_place = True then self has the outliers removed as well (in _place)
- find_outliers_sliding_window(threshold=3, in_place=False, verbose=True, periods=[[]], excluded_periods=[[]], component='NE', window_len=15, automatic=True)¶
Find outliers using sliding windows
- find_outliers_vondrak(threshold=10, fc=2.0, in_place=False, verbose=True, periods=[[]], excluded_periods=[[]], component='NE')¶
Find outliers using a Vondrak filter
Gts plot methods¶
- class pyacs.gts.Gts.Gts(code=None, lat=None, lon=None, h=None, X0=None, Y0=None, Z0=None, t0=None, data=None, data_xyz=None, data_corr_neu=None, data_corr_xyz=None, offsets_dates=[], offsets_values=None, outliers=[], annual=None, semi_annual=None, velocity=None, ifile=None, log=None, metadata=None)[source]¶
- plot(title=None, loffset=True, loutliers=True, verbose=False, date=[], yaxis=None, min_yaxis=None, yupaxis=None, xticks_minor_locator=1, lcomponent=['N', 'E', 'U'], error_scale=1.0, lperiod=[[]], lvline=[], save_dir_plots='.', save=None, show=True, unit='mm', date_unit='cal', date_ref=0.0, center=True, superimposed=None, lcolor=['r', 'g', 'c', 'm', 'y', 'k', 'b'], label=None, legend=False, set_zero_at_date=None, grid=True, plot_size=None, info=[], xlabel_fmt=None, **kwargs)¶
Create a plot of a North-East-Up time series and related info (offsets, outliers) using Matplotlib
Coordinates of the time series are assumed to be in meters default plots units will be mm; Use unit=’m’ to get meters instead
- Parameters
title – string to be added to the site name as a plot title
loffset – boolean print a dash vertical line at offset dates
loutliers – boolean print outliers
verbose – boolean verbose mode
date – [sdate,edate] start and end date for plots sdate and edate in decimal years if date_unit is either ‘decyear’ or ‘cal’, or in days if date_unit is ‘days’
yaxis – [min_y,max_y] min and max value for the yaxis if not provided automatically adjusted
yupaxis – same as yaxis but applies to the up component only
xticks_minor_locator – where xticks_minor_locator will be placed. Float when date_unit is ‘decyear’ or ‘days’, a string ‘%Y’,’%m’,’%d’ is date_unit is ‘cal’.
lcomponent – list of components to be plotted (default =[‘N’,’E’,’U’])
error_scale – scaling factor for error bars (default = 1.0, 0 means no error bar)
lperiod – list of periods to be drawn in background (color=light salmon)
lvline – list of dates where vertical lines will be drawn in background (color=green)
save_dir_plots – directory used for saving the plots
save – name, save the plot into name, if simply True an automatic name is given
show – boolean, is True show the plot
unit – ‘m’,’cm’,’mm’, default=’mm’
date_unit – ‘decyear’ or ‘cal’ or ‘days’, default=’decyear’
date_ref – reference date, default=0.0
center – boolean, if True the y_axis is centered around the mean value for the plotted period
superimposed – if a gts is provided, it is superimposed to the master, default=None
lcolor – color list used for the superimposed time series, default=[‘r’,’g’,’c’,’m’,’y’,’k’,’b’]
label – label for superimposed time series to be displayed in legend, default=None
legend – boolean. Set true to display label for superimposed time series, default=False
set_zero_at_date – date at which the master and superimposed gts will be equal (default=None). date can also be a list with [date,offset_north,offset_east,offset_up]
plot_size – plot size as a tuple. Default, best guess.
grid – boolean
info – title to appear in time series subplots
**kwargs –
any argument to be passed to matplotlib.pyplot.errorbar
- Note
The list of kwargs are:
{ ‘linewidth’ : 0, ‘marker’ : marker_main_symbol , ‘markersize’ : marker_main_size, ‘markerfacecolor’ : marker_main_color, ‘markeredgecolor’ : marker_main_color, ‘markeredgewidth’ : marker_main_colorlw, ‘ecolor’ : error_bar_color, ‘elinewidth’ : error_bar_linewidth, ‘capsize’ : error_bar_capsize }