| Title: | Estimate Functional and Stochastic Parameters of Linear Models with Correlated Residuals |
|---|---|
| Description: | Implements the Generalized Method of Wavelet Moments with Exogenous Inputs estimator (GMWMX) presented in Cucci, D. A., Voirol, L., Kermarrec, G., Montillet, J. P., and Guerrier, S. (2023) <doi:10.1007/s00190-023-01702-8>. The GMWMX estimator allows to estimate functional and stochastic parameters of linear models with correlated residuals. The 'gmwmx' package provides functions to estimate, compare and analyze models, utilities to load and work with Global Navigation Satellite System (GNSS) data as well as methods to compare results with the Maximum Likelihood Estimator (MLE) implemented in Hector. |
| Authors: | Davide Antonio Cucci [aut], Lionel Voirol [aut, cre], Stéphane Guerrier [aut], Jean-Philippe Montillet [ctb], Gaël Kermarrec [ctb] |
| Maintainer: | Lionel Voirol <[email protected]> |
| License: | AGPL-3 |
| Version: | 1.0.3 |
| Built: | 2026-05-31 09:13:49 UTC |
| Source: | https://github.com/cran/gmwmx |
Data from station COLA of the Plate Boundary Observatory
colacola
A gnssts object of the East position (dE) of the COLA station
https://data.unavco.org/archive/gnss/products/position/COLA/COLA.pbo.igs14.pos
gnsstsmodel objects.Compare graphically two gnsstsmodel objects.
compare_fits(fit_1, fit_2, main = NULL, y_unit = "mm", x_unit = "days")compare_fits(fit_1, fit_2, main = NULL, y_unit = "mm", x_unit = "days")
fit_1 |
A |
fit_2 |
A |
main |
A |
y_unit |
A |
x_unit |
A |
No return value. Produce a plot comparing two estimated models.
## Not run: data(cola) fit_gmwmx_1 = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") fit_gmwmx_2 = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+powerlaw") compare_fits(fit_gmwmx_1, fit_gmwmx_2) ## End(Not run)## Not run: data(cola) fit_gmwmx_1 = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") fit_gmwmx_2 = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+powerlaw") compare_fits(fit_gmwmx_1, fit_gmwmx_2) ## End(Not run)
Define matrix A of the functional model
create_A_matrix(t_nogap, jumps, n_seasonal)create_A_matrix(t_nogap, jumps, n_seasonal)
t_nogap |
A |
jumps |
A |
n_seasonal |
An |
Matrix A in order to compute the functional component of the model in a linear fashion
n= 10*365 jump_vec <- c(200, 300, 500) nbr_sin = 2 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) head(A) A <- create_A_matrix(1:n, jumps = NULL, n_seasonal = nbr_sin) head(A)n= 10*365 jump_vec <- c(200, 300, 500) nbr_sin = 2 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) head(A) A <- create_A_matrix(1:n, jumps = NULL, n_seasonal = nbr_sin) head(A)
Create a gnssts object
create.gnssts(t, y, jumps = NULL, sampling_period = 1)create.gnssts(t, y, jumps = NULL, sampling_period = 1)
t |
A |
y |
A |
jumps |
A |
sampling_period |
An |
A gnssts object.
phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj)phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj)
Estimate a stochastic model in a two-steps procedure using the GMWMX estimator.
estimate_gmwmx( x, theta_0, n_seasonal = 1, model_string, method = "L-BFGS-B", maxit = 1e+06, ci = FALSE, k_iter = 1 )estimate_gmwmx( x, theta_0, n_seasonal = 1, model_string, method = "L-BFGS-B", maxit = 1e+06, ci = FALSE, k_iter = 1 )
x |
A |
theta_0 |
A |
n_seasonal |
An |
model_string |
A |
method |
A |
maxit |
An |
ci |
A |
k_iter |
An |
A gnsstsmodel object.
## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") ## End(Not run)## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") ## End(Not run)
Estimate a stochastic model based on the MLE and the Hector implementation.
estimate_hector( x, n_seasonal = 1, model_string, likelihood_method = "AmmarGrag", cleanup = TRUE )estimate_hector( x, n_seasonal = 1, model_string, likelihood_method = "AmmarGrag", cleanup = TRUE )
x |
A |
n_seasonal |
An |
model_string |
A |
likelihood_method |
A |
cleanup |
A |
A gnsstsmodel object.
## Not run: cola = PBO_get_station(station_name = "COLA", column = "dE", time_range = c(51130, 52000)) fit_mle = estimate_hector(x = cola, n_seasonal = 1, model_string = "wn+matern") ## End(Not run)## Not run: cola = PBO_get_station(station_name = "COLA", column = "dE", time_range = c(51130, 52000)) fit_mle = estimate_hector(x = cola, n_seasonal = 1, model_string = "wn+matern") ## End(Not run)
Extract offsets for a PBO station
PBO_get_offsets(station_name)PBO_get_offsets(station_name)
station_name |
A |
A vector specifying the offsets of a PBO station.
## Not run: pbo_cola_offsets = PBO_get_offsets(station_name = "COLA") pbo_cola_offsets ## End(Not run)## Not run: pbo_cola_offsets = PBO_get_offsets(station_name = "COLA") pbo_cola_offsets ## End(Not run)
Load station data from PBO
PBO_get_station(station_name, column, time_range = c(-Inf, Inf), scale = 1)PBO_get_station(station_name, column, time_range = c(-Inf, Inf), scale = 1)
station_name |
A |
column |
A |
time_range |
A |
scale |
A |
A gnssts object that contains the data associated with the specified PBO station.
## Not run: pbo_cola_data = PBO_get_station("COLA", column="dE") str(pbo_cola_data) ## End(Not run)## Not run: pbo_cola_data = PBO_get_station("COLA", column="dE") str(pbo_cola_data) ## End(Not run)
gnsstsmodel object.Plotting method for a gnsstsmodel object.
## S3 method for class 'gnsstsmodel' plot( x, main = NULL, y_unit = "mm", x_unit = "days", legend_position = "bottomright", legend_position_wv = "bottomleft", ... )## S3 method for class 'gnsstsmodel' plot( x, main = NULL, y_unit = "mm", x_unit = "days", legend_position = "bottomright", legend_position_wv = "bottomleft", ... )
x |
A |
main |
A |
y_unit |
A |
x_unit |
A |
legend_position |
A |
legend_position_wv |
A |
... |
Additional graphical parameters. |
No return value. Plot a gnsstsmodel object.
## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") plot(fit_gmwmx) ## End(Not run)## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") plot(fit_gmwmx) ## End(Not run)
gnsstsmodel object.Print method for a gnsstsmodel object.
## S3 method for class 'gnsstsmodel' print(x, ...)## S3 method for class 'gnsstsmodel' print(x, ...)
x |
A |
... |
Additional graphical parameters. |
No return value. Print a gnsstsmodel object.
## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") print(fit_gmwmx) ## End(Not run)## Not run: data(cola) fit_gmwmx = estimate_gmwmx(x = cola, theta_0 = c(0.1,0.1,0.1,0.1), n_seasonal = 1, model_string = "wn+matern") print(fit_gmwmx) ## End(Not run)
Read a gnssts object
read.gnssts(filename, format = "mom")read.gnssts(filename, format = "mom")
filename |
A |
format |
A |
Return a gnssts object.
phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj) ## Not run: write.gnssts(x = gnssts_obj, filename = "test.mom") gnssts_obj <-read.gnssts(filename = "test.mom", format = "mom") ## End(Not run)phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj) ## Not run: write.gnssts(x = gnssts_obj, filename = "test.mom") gnssts_obj <-read.gnssts(filename = "test.mom", format = "mom") ## End(Not run)
Remove outliers from a gnssts object using Hector
remove_outliers_hector(x, n_seasonal, IQ_factor = 3, cleanup = TRUE)remove_outliers_hector(x, n_seasonal, IQ_factor = 3, cleanup = TRUE)
x |
A |
n_seasonal |
An |
IQ_factor |
A |
cleanup |
An |
A gnssts object.
phase = 0.45 amplitude = 2.5 sigma2_wn = 15 bias = 0 trend = 5/365.25 cosU = amplitude*cos(phase) sinU = amplitude*sin(phase) n= 2*365 # define time at which there are jumps jump_vec = c(100, 200) jump_height = c(10, 20) # generate residuals eps = rnorm(n = n, sd = sqrt(sigma2_wn)) # add trend, gaps and sin A = create_A_matrix(1:length(eps), jump_vec, n_seasonal = 1) # define beta x_0 = c(bias, trend, jump_height, cosU, sinU) # create time series yy = A %*% x_0 + eps plot(yy, type="l") n_outliers = 30 set.seed(123) id_outliers=sample(150:350, size = n_outliers) val_outliers = rnorm(n = n_outliers, mean = max(yy)+10, sd = 5) yy[id_outliers] = val_outliers plot(yy, type="l") # save signal in temp gnssts_obj = create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) ## Not run: clean_yy = remove_outliers_hector(x=gnssts_obj, n_seasonal = 1) plot(clean_yy$t, clean_yy$y, type="l") ## End(Not run)phase = 0.45 amplitude = 2.5 sigma2_wn = 15 bias = 0 trend = 5/365.25 cosU = amplitude*cos(phase) sinU = amplitude*sin(phase) n= 2*365 # define time at which there are jumps jump_vec = c(100, 200) jump_height = c(10, 20) # generate residuals eps = rnorm(n = n, sd = sqrt(sigma2_wn)) # add trend, gaps and sin A = create_A_matrix(1:length(eps), jump_vec, n_seasonal = 1) # define beta x_0 = c(bias, trend, jump_height, cosU, sinU) # create time series yy = A %*% x_0 + eps plot(yy, type="l") n_outliers = 30 set.seed(123) id_outliers=sample(150:350, size = n_outliers) val_outliers = rnorm(n = n_outliers, mean = max(yy)+10, sd = 5) yy[id_outliers] = val_outliers plot(yy, type="l") # save signal in temp gnssts_obj = create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) ## Not run: clean_yy = remove_outliers_hector(x=gnssts_obj, n_seasonal = 1) plot(clean_yy$t, clean_yy$y, type="l") ## End(Not run)
Write a gnssts object
write.gnssts(x, filename, format = "mom")write.gnssts(x, filename, format = "mom")
x |
A |
filename |
A |
format |
A |
No return value. Write a gnssts object in a .mom file by default.
phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj) ## Not run: write.gnssts(x = gnssts_obj, filename = "test.mom") ## End(Not run)phase <- 0.45 amplitude <- 2.5 sigma2_wn <- 15 bias <- 0 trend <- 5 / 365.25 cosU <- amplitude * cos(phase) sinU <- amplitude * sin(phase) year <- 5 n <- year * 365 jump_vec <- c(200, 300, 500) jump_height <- c(10, 15, 20) nbr_sin <- 1 A <- create_A_matrix(1:n, jump_vec, n_seasonal = nbr_sin) x_0 <- c(bias, trend, cosU, sinU, jump_height) eps <- rnorm(n = n, sd = sqrt(sigma2_wn)) yy <- A %*% x_0 + eps gnssts_obj <- create.gnssts(t = 1:length(yy), y = yy, jumps = jump_vec) str(gnssts_obj) ## Not run: write.gnssts(x = gnssts_obj, filename = "test.mom") ## End(Not run)