fitFMM: Fitting FMM models

Description Usage Arguments Details Value References Examples

View source: R/fitFMM.R

Description

fitFMM() is used to fit FMM models. The only required argument to fit FMM models is the input data. By default it is assumed that time points, corresponding to a single time period, are equally spaced from 0 to .

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
fitFMM(
  vData,
  nPeriods = 1,
  timePoints = NULL,
  nback = 1,
  betaOmegaRestrictions = 1:nback,
  maxiter = nback,
  stopFunction = alwaysFalse,
  lengthAlphaGrid = 48,
  lengthOmegaGrid = 24,
  numReps = 3,
  showProgress = TRUE,
  showTime = FALSE,
  parallelize = FALSE,
  restrExactSolution = FALSE
)

Arguments

vData

A numeric vector containing the data to be fitted a FMM model.

nPeriods

A numeric value specifying the number of periods at which vData is observed.

timePoints

A numeric vector containing the time points at which each data of one single period is observed. The default value is NULL, in which case they are equally spaced in range [0, 2π]. It must be between 0 and .

nback

Number of FMM components to be fitted. Its default value is 1.

betaOmegaRestrictions

An integer vector of length nback indicating which FMM waves are constrained to have equal beta and omega parameters. For example, c(1,1,1,2,2) indicates that beta1=beta2=beta3 and beta4=beta5 as well as omega1=omega2=omega3 and omega4=omega5. In brief, some waves are restricted to have the same shape. Its default value is the sequence 1:nback to fit the FMM model without restrictions on shape parameters (beta and omega).

maxiter

Maximum number of iterations for the backfitting algorithm. By default, it is set at nback.

stopFunction

Function to check the convergence criterion for the backfitting algorithm (see Details).

lengthAlphaGrid

Precision of the grid of alpha in the search of the best model. If it is increased, more possible values of alpha will be considered, resulting in an increasing in the computation time too. By default, it is set to 48 possible values of alpha, equally spaced between 0 and .

lengthOmegaGrid

Precision of the grid of omega in the search of the best model. If it is increased, more possible values of omega will be considered, resulting in an increasing in the computation time too. By default it is set to 24 possible values of omega, equally spaced between 0 and 1 in a logarithmic way.

numReps

Number of times (alpha, omega) parameters are refined.

showProgress

TRUE to display a progress indicator on the console.

showTime

TRUE to display execution time on the console.

parallelize

TRUE to use parallelized procedure to fit FMM model. Its default value is FALSE. When it is TRUE, the number of cores to be used is equal to 12, or if the machine has less, the number of cores - 1.

restrExactSolution

FALSE to use an aproximated algorithm to fit the FMM model with restrictions (default). If TRUE is specified, the exact solution is computed.

Details

Data will be collected over nPeriods periods. When nPeriods > 1 the fitting is carried out by averaging the data collected at each time point across all considered periods. The model is fitting to summarized data. timePoints is a n-length numeric vector where n is the number of different time points per period.

The stopFunction argument can either be the functions alwaysFalse or R2 included in the package or user-defined functions that have the same arguments. The included functions serve for the following:

Value

An S4 object of class 'FMM' with information about the fitted model. The object contains the following slots:

@timePoints

The time points as specified by the input argument. It is a numeric vector containing the time points at which each data of one single period is observed.

@data

The data as specified by the input argument. It is a numeric vector containing the data to be fitted a FMM model. Data could be collected over multiple periods.

@summarizedData

When the data has more than one period, a numeric vector containing data averaging the data at each time point across all considered periods.

@nPeriods

A numeric value containing the number of periods in data as specified by the input argument.

@fittedValues

A numeric vector of the fitted values by the FMM model.

@M

A numeric value of the estimated intercept parameter M.

@A

A numeric value or vector of the estimated FMM wave amplitude parameter(s) A.

@alpha

A numeric value or vector of the estimated FMM wave phase translation parameter(s) α.

@beta

A numeric value or vector of the estimated FMM wave skewness parameter(s) β.

@omega

A numeric value or vector of the estimated FMM wave kurtosis parameter(s) ω.

@SSE

A numeric value of the sum of the residual squares values.

@R2

A numeric vector specifying the explained variance by each of the fitted FMM components.

@nIter

A numeric value specifying the number of iterations of the fitting algorithm.

References

Rueda C, Larriba Y, Peddada SD (2019). Frequency Modulated Moebius Model Accurately Predicts Rhythmic Signals in Biological and Physical Sciences. Scientific reports, 9 (1), 18701. https://www.nature.com/articles/s41598-019-54569-1

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
# A monocomponent FMM model is fitted.
FMM_data <- generateFMM(2, 3, 1.5, 2.3, 0.1,
                        from = 0, to = 2*pi, length.out = 100,
                        outvalues = TRUE, sigmaNoise = 0.3, plot = FALSE)
fit <- fitFMM(FMM_data$y, lengthAlphaGrid = 10, lengthOmegaGrid = 10)
summary(fit)

# To see the differences between number of repetitions.
FMM_data <- generateFMM(2, 1, 1.5, 1.1 ,0.14 , outvalues = TRUE, sigmaNoise = 0.15, plot=TRUE)
fit1 <- fitFMM(FMM_data$y, lengthAlphaGrid = 6, lengthOmegaGrid = 3, numReps = 1)
fit2 <- fitFMM(FMM_data$y, lengthAlphaGrid = 6, lengthOmegaGrid = 3, numReps = 6,
               showProgress = FALSE)  # suppress progress messages
getSSE(fit1)
getSSE(fit2)

# Finer resolution grid.
fit3 <- fitFMM(FMM_data$y, lengthAlphaGrid = 10, lengthOmegaGrid = 5, numReps = 1)
getSSE(fit3)

# Two component FMM model with beta and omega restricted
restFMM2w_data <- generateFMM(M = 3, A = c(7, 4), alpha = c(0.5, 5), beta = c(rep(3, 2)),
                              omega = rep(0.05, 2), from = 0, to = 2*pi, length.out = 100,
                              sigmaNoise = 0.3, plot = FALSE)
fit2w.rest <- fitFMM(restFMM2w_data$y, nback = 2, maxiter = 1, numReps = 1,
                     lengthAlphaGrid = 15, lengthOmegaGrid = 10,
                     betaOmegaRestrictions = c(1, 1))
plotFMM(fit2w.rest, components = TRUE)

FMM documentation built on Dec. 17, 2021, 5:06 p.m.