Title: | S3 Classes and Methods for Tidy Functional Data |
---|---|
Description: | Defines S3 vector data types for vectors of functional data (grid-based, spline-based or functional principal components-based) with all arithmetic and summary methods, derivation, integration and smoothing, plotting, data import and export, and data wrangling, such as re-evaluating, subsetting, sub-assigning, zooming into sub-domains, or extracting functional features like minima/maxima and their locations. The implementation allows including such vectors in data frames for joint analysis of functional and scalar variables. |
Authors: | Fabian Scheipl [aut, cre] , Jeff Goldsmith [aut], Julia Wrobel [ctb] , Maximilian Mücke [ctb] , Sebastian Fischer [ctb] , Trevor Hastie [ctb] (softImpute author), Rahul Mazumder [ctb] (softImpute author), Chen Meng [ctb] (mogsa author) |
Maintainer: | Fabian Scheipl <[email protected]> |
License: | AGPL (>= 3) |
Version: | 0.3.5 |
Built: | 2024-11-20 18:28:44 UTC |
Source: | https://github.com/tidyfun/tf |
Various converters to turn tfb
- or tfd
-vectors into data.frames or
matrices, or even an actual R function.
## S3 method for class 'tf' as.data.frame(x, row.names = NULL, optional = FALSE, unnest = FALSE, ...) ## S3 method for class 'tf' as.matrix(x, arg, interpolate = FALSE, ...) ## S3 method for class 'tf' as.function(x, ...)
## S3 method for class 'tf' as.data.frame(x, row.names = NULL, optional = FALSE, unnest = FALSE, ...) ## S3 method for class 'tf' as.matrix(x, arg, interpolate = FALSE, ...) ## S3 method for class 'tf' as.function(x, ...)
x |
a |
row.names |
|
optional |
not used |
unnest |
if |
... |
additional arguments to be passed to or from methods. |
arg |
a vector of argument values / evaluation points for |
interpolate |
should functions be evaluated (i.e., inter-/extrapolated)
for values in |
for as.data.frame.tf
: if unnest
is FALSE
(default), a
one-column data.frame
with a tf
-column containing x
. if unnest
is
TRUE
, a 3-column data frame with columns id
for the (unique) names of
x
or a numeric identifier, arg
and value
, with each row containing
one function evaluation at the original arg
-values.
for as.matrix.tf
: a matrix with one row per function and one
column per arg
.
for as.function.tf
: an R function with argument arg
that
evaluates x
on arg
and returns the list of function values
See above.
ensure_list(x)
ensure_list(x)
x |
any input |
x
turned into a list.
Other tidyfun developer tools:
prep_plotting_arg()
,
unique_id()
Compute (truncated) orthonormal eigenfunctions and scores for (partially missing) data on a common (potentially non-equidistant) grid.
fpc_wsvd(data, arg, pve = 0.995) ## S3 method for class 'matrix' fpc_wsvd(data, arg, pve = 0.995) ## S3 method for class 'data.frame' fpc_wsvd(data, arg, pve = 0.995)
fpc_wsvd(data, arg, pve = 0.995) ## S3 method for class 'matrix' fpc_wsvd(data, arg, pve = 0.995) ## S3 method for class 'data.frame' fpc_wsvd(data, arg, pve = 0.995)
data |
numeric matrix of function evaluations (each row is one curve, no NAs) |
arg |
numeric vector of argument values |
pve |
percentage of variance explained |
Performs a weighted SVD with trapezoidal quadrature weights s.t. returned
vectors represent (evaluations of)
orthonormal eigenfunctions , not eigenvectors
, specifically:
given quadrature weights
, not
;
not
,
c.f.
mogsa::wsvd()
.
For incomplete data, this uses an adaptation of softImpute::softImpute()
,
see references. Note that will not work well for data on a common grid if more
than a few percent of data points are missing, and it breaks down completely
for truly irregular data with no/few common timepoints, even if observed very
densely. For such data, either re-evaluate on a common grid first or use more
advanced FPCA approaches like refund::fpca_sc()
,
see last example for tfb_fpc()
a list with entries
mu
estimated mean function (numeric vector)
efunctions
estimated FPCs (numeric matrix, columns represent FPCs)
scores
estimated FPC scores (one row per observed curve)
npc
how many FPCs were returned for the given pve
(integer)
scoring_function
a function that returns FPC scores for new data
and given eigenfunctions, see tf:::.fpc_wsvd_scores
for an example.
Trevor Hastie, Rahul Mazumder, Cheng Meng, Fabian Scheipl
code adapted from / inspired by mogsa::wsvd()
by Cheng Meng
and softImpute::softImpute()
by Trevor Hastie and Rahul Mazumder.
Meng C (2023).
mogsa
: Multiple omics data integrative clustering and gene set analysis.
https://bioconductor.org/packages/mogsa.
Mazumder, Rahul, Hastie, Trevor, Tibshirani, Robert (2010). “Spectral regularization algorithms for learning large incomplete matrices.” The Journal of Machine Learning Research, 11, 2287-2322.
Hastie T, Mazumder R (2021).
softImpute
: Matrix Completion via Iterative Soft-Thresholded SVD.
R package version 1.4-1, https://CRAN.R-project.org/package=softImpute.
Other tfb-class:
tfb
,
tfb_fpc()
,
tfb_spline()
Other tfb_fpc-class:
tfb_fpc()
tf
in a vectorThese functions extract (user-specified) function-wise summary statistics
from each entry in a tf
-vector. To summarize a vector of functions at each
argument value, see ?tfsummaries
. Note that these will tend to yield lots
of NA
s for irregular tfd
unless you set a tf_evaluator()
-function
that does inter- and extrapolation for them beforehand.
tf_fwise(x, .f, arg = tf_arg(x), ...) tf_fmax(x, arg = tf_arg(x), na.rm = FALSE) tf_fmin(x, arg = tf_arg(x), na.rm = FALSE) tf_fmedian(x, arg = tf_arg(x), na.rm = FALSE) tf_frange(x, arg = tf_arg(x), na.rm = FALSE, finite = FALSE) tf_fmean(x, arg = tf_arg(x)) tf_fvar(x, arg = tf_arg(x)) tf_fsd(x, arg = tf_arg(x)) tf_crosscov(x, y, arg = tf_arg(x)) tf_crosscor(x, y, arg = tf_arg(x))
tf_fwise(x, .f, arg = tf_arg(x), ...) tf_fmax(x, arg = tf_arg(x), na.rm = FALSE) tf_fmin(x, arg = tf_arg(x), na.rm = FALSE) tf_fmedian(x, arg = tf_arg(x), na.rm = FALSE) tf_frange(x, arg = tf_arg(x), na.rm = FALSE, finite = FALSE) tf_fmean(x, arg = tf_arg(x)) tf_fvar(x, arg = tf_arg(x)) tf_fsd(x, arg = tf_arg(x)) tf_crosscov(x, y, arg = tf_arg(x)) tf_crosscor(x, y, arg = tf_arg(x))
x |
a |
.f |
a function or formula that is applied to each entry of |
arg |
defaults to standard argument values of |
... |
additional arguments for |
na.rm |
a logical indicating whether missing values should be removed. |
finite |
logical, indicating if all non-finite elements should be omitted. |
y |
a |
tf_fwise
turns x
into a list of data.frames with columns arg
and values
internally, so the function/formula in .f
gets a data.frame
.x
with these columns, see examples below or source code for tf_fmin()
,
tf_fmax()
, etc
a list (or vector) of the same length as x
with the respective
summaries
tf_fwise()
: User-specified function-wise summary statistics
tf_fmax()
: maximal value of each function
tf_fmin()
: minimal value of each function
tf_fmedian()
: median value of each function
tf_frange()
: range of values of each function
tf_fmean()
: mean of each function:
tf_fvar()
: variance of each function:
tf_fsd()
: standard deviation of each function:
tf_crosscov()
: cross-covariances between two functional vectors:
tf_crosscor()
: cross-correlation between two functional vectors:
tf_crosscov(x, y) / (tf_fsd(x) * tf_fsd(y))
Other tidyfun summary functions:
tfsummaries
x <- tf_rgp(3) layout(t(1:3)) plot(x, col = 1:3) # each function's values to [0,1]: x_clamp <- (x - tf_fmin(x)) / (tf_fmax(x) - tf_fmin(x)) plot(x_clamp, col = 1:3) # standardize each function to have mean / integral 0 and sd 1: x_std <- (x - tf_fmean(x)) / tf_fsd(x) tf_fvar(x_std) == c(1, 1, 1) plot(x_std, col = 1:3) # Custom functions: # 80%tiles of each function's values: tf_fwise(x, \(.x) quantile(.x$value, 0.8)) |> unlist() # minimal value of each function for t > 0.5 tf_fwise(x, \(.x) min(.x$value[.x$arg > 0.5])) |> unlist() tf_crosscor(x, -x) tf_crosscov(x, x) == tf_fvar(x)
x <- tf_rgp(3) layout(t(1:3)) plot(x, col = 1:3) # each function's values to [0,1]: x_clamp <- (x - tf_fmin(x)) / (tf_fmax(x) - tf_fmin(x)) plot(x_clamp, col = 1:3) # standardize each function to have mean / integral 0 and sd 1: x_std <- (x - tf_fmean(x)) / tf_fsd(x) tf_fvar(x_std) == c(1, 1, 1) plot(x_std, col = 1:3) # Custom functions: # 80%tiles of each function's values: tf_fwise(x, \(.x) quantile(.x$value, 0.8)) |> unlist() # minimal value of each function for t > 0.5 tf_fwise(x, \(.x) min(.x$value[.x$arg > 0.5])) |> unlist() tf_crosscor(x, -x) tf_crosscov(x, x) == tf_fvar(x)
in_range
and its infix-equivalent %inr%
return TRUE
for all values in
the numeric vector f
that are within the range of values in r
.
in_range(f, r) f %inr% r
in_range(f, r) f %inr% r
f |
a numeric vector |
r |
numeric vector used to specify a range, only the minimum and maximum
of |
a logical
vector of the same length as f
Other tidyfun utility functions:
tf_arg()
,
tf_zoom()
base
plots for tf
sSome base
functions for displaying functional data in
spaghetti- (i.e., line plots) and lasagna- (i.e., heat map) flavors.
## S3 method for class 'tf' plot( x, y, n_grid = 50, points = is_irreg(x), type = c("spaghetti", "lasagna"), alpha = min(1, max(0.05, 2/length(x))), ... ) ## S3 method for class 'tf' lines(x, arg, n_grid = 50, alpha = min(1, max(0.05, 2/length(x))), ...) ## S3 method for class 'tf' points( x, arg, n_grid = NA, alpha = min(1, max(0.05, 2/length(x))), interpolate = FALSE, ... )
## S3 method for class 'tf' plot( x, y, n_grid = 50, points = is_irreg(x), type = c("spaghetti", "lasagna"), alpha = min(1, max(0.05, 2/length(x))), ... ) ## S3 method for class 'tf' lines(x, arg, n_grid = 50, alpha = min(1, max(0.05, 2/length(x))), ...) ## S3 method for class 'tf' points( x, arg, n_grid = NA, alpha = min(1, max(0.05, 2/length(x))), interpolate = FALSE, ... )
x |
an |
y |
(optional) numeric vector to be used as |
n_grid |
minimal size of equidistant grid used for plotting, defaults to 50. See details. |
points |
should the original evaluation points be marked by points?
Defaults to |
type |
"spaghetti": line plots, "lasagna": heat maps. |
alpha |
alpha-value (see |
... |
additional arguments for |
arg |
evaluation grid (vector) |
interpolate |
should functions be evaluated (i.e., inter-/extrapolated) for arg for which no original data is available? Only relevant for tfd, defaults to FALSE |
If no second argument y
is given, evaluation points (arg
) for the functions
are given by the union of the tf
's arg
and an equidistant grid
over its domain with n_grid
points. If you want to only see the original
data for tfd
-objects without inter-/extrapolation, use n_grid < 1
or
n_grid = NA
.
the plotted tf
-object, invisibly.
Swihart, J B, Caffo, Brian, James, D B, Strand, Matthew, Schwartz, S B, Punjabi, M N (2010). “Lasagna plots: a saucy alternative to spaghetti plots.” Epidemiology (Cambridge, Mass.), 21(5), 621–625.
(internal function exported for re-use in upstream packages)
prep_plotting_arg(f, n_grid)
prep_plotting_arg(f, n_grid)
f |
a |
n_grid |
length of evaluation grid |
a semi-regular grid rounded down to appropriate resolution
Other tidyfun developer tools:
ensure_list()
,
unique_id()
Print/format tf
-objects.
## S3 method for class 'tf' print(x, n = 5, ...) ## S3 method for class 'tfd_reg' print(x, n = 5, ...) ## S3 method for class 'tfd_irreg' print(x, n = 5, ...) ## S3 method for class 'tfb' print(x, n = 5, ...) ## S3 method for class 'tf' format( x, digits = 2, nsmall = 0, width = options()$width, n = 5, prefix = TRUE, ... )
## S3 method for class 'tf' print(x, n = 5, ...) ## S3 method for class 'tfd_reg' print(x, n = 5, ...) ## S3 method for class 'tfd_irreg' print(x, n = 5, ...) ## S3 method for class 'tfb' print(x, n = 5, ...) ## S3 method for class 'tf' format( x, digits = 2, nsmall = 0, width = options()$width, n = 5, prefix = TRUE, ... )
x |
any R object (conceptually); typically numeric. |
n |
how many elements of |
... |
further arguments passed to or from other methods. |
digits |
a positive integer indicating how many significant digits
are to be used for
numeric and complex |
nsmall |
the minimum number of digits to the right of the decimal
point in formatting real/complex numbers in non-scientific formats.
Allowed values are |
width |
|
prefix |
used internally. |
prints out x
and returns it invisibly
tfd
-objectsThese are the currently available evaluator
-functions for tfd
-objects,
which control how the entries are inter-/extrapolated to previously unseen
arg
-values. They all are merely wrappers around zoo::na.fill()
,
zoo::na.approx()
, etc... Note that these are not meant to be called directly –
they are internal functions used by tf_evaluate.tfd()
to do its thing.
The list:
tf_approx_linear
for linear interpolation without extrapolation (i.e.,
zoo::na.approx()
with na.rm = FALSE
) – this is the default,
tf_approx_spline
for cubic spline interpolation, (i.e., zoo::na.spline()
with na.rm = FALSE
),
tf_approx_none
in order to not inter-/extrapolate ever (i.e., zoo::na.fill()
with fill = NA
)
tf_approx_fill_extend
for linear interpolation and constant extrapolation
(i.e., zoo::na.fill()
with fill = "extend"
)
tf_approx_locf
for "last observation carried forward" (i.e.,
zoo::na.locf()
with na.rm = FALSE
and
tf_approx_nocb
for "next observation carried backward" (i.e.,
zoo::na.locf()
with na.rm = FALSE, fromLast = TRUE
).
For implementing your own, see source code of tf:::zoo_wrapper
.
tf_approx_linear(x, arg, evaluations) tf_approx_spline(x, arg, evaluations) tf_approx_none(x, arg, evaluations) tf_approx_fill_extend(x, arg, evaluations) tf_approx_locf(x, arg, evaluations) tf_approx_nocb(x, arg, evaluations)
tf_approx_linear(x, arg, evaluations) tf_approx_spline(x, arg, evaluations) tf_approx_none(x, arg, evaluations) tf_approx_fill_extend(x, arg, evaluations) tf_approx_locf(x, arg, evaluations) tf_approx_nocb(x, arg, evaluations)
x |
new |
arg |
the |
evaluations |
the function values at |
a vector of values of the function defined by the given
=
(arg, evaluations)
-tuples at new argument values x
.
tfd
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
Other tidyfun inter/extrapolation functions:
tf_evaluate()
,
tf_interpolate()
tf
-objectsA bunch of methods & utilities that do what they say: get or set the
respective attributes of a tf
-object.
tf_arg(f) tf_evaluations(f) tf_count(f) tf_domain(f) tf_domain(x) <- value tf_evaluator(f) tf_evaluator(x) <- value tf_basis(f, as_tfd = FALSE) tf_arg(x) <- value ## S3 replacement method for class 'tfd_irreg' tf_arg(x) <- value ## S3 replacement method for class 'tfd_reg' tf_arg(x) <- value ## S3 replacement method for class 'tfb' tf_arg(x) <- value ## S3 method for class 'tfb' coef(object, ...) ## S3 method for class 'tf' rev(x) ## S3 method for class 'tf' is.na(x) ## S3 method for class 'tfd_irreg' is.na(x) is_tf(x) is_tfd(x) is_reg(x) is_tfd_reg(x) is_irreg(x) is_tfd_irreg(x) is_tfb(x) is_tfb_spline(x) is_tfb_fpc(x)
tf_arg(f) tf_evaluations(f) tf_count(f) tf_domain(f) tf_domain(x) <- value tf_evaluator(f) tf_evaluator(x) <- value tf_basis(f, as_tfd = FALSE) tf_arg(x) <- value ## S3 replacement method for class 'tfd_irreg' tf_arg(x) <- value ## S3 replacement method for class 'tfd_reg' tf_arg(x) <- value ## S3 replacement method for class 'tfb' tf_arg(x) <- value ## S3 method for class 'tfb' coef(object, ...) ## S3 method for class 'tf' rev(x) ## S3 method for class 'tf' is.na(x) ## S3 method for class 'tfd_irreg' is.na(x) is_tf(x) is_tfd(x) is_reg(x) is_tfd_reg(x) is_irreg(x) is_tfd_irreg(x) is_tfb(x) is_tfb_spline(x) is_tfb_fpc(x)
f |
an |
x |
an |
value |
for |
as_tfd |
should the basis be returned as a |
object |
as usual |
... |
dots |
either the respective attribute or, for setters (assignment functions), the input object with modified properties.
Other tidyfun utility functions:
in_range()
,
tf_zoom()
Data depths for functional data. Currently implemented:
Modified Band-2 Depth
Modified Epigraph Index.
tf_depth(x, arg, depth = "MBD", na.rm = TRUE, ...) ## S3 method for class 'matrix' tf_depth(x, arg, depth = c("MBD", "MEI"), na.rm = TRUE, ...) ## S3 method for class 'tf' tf_depth(x, arg, depth = "MBD", na.rm = TRUE, ...)
tf_depth(x, arg, depth = "MBD", na.rm = TRUE, ...) ## S3 method for class 'matrix' tf_depth(x, arg, depth = c("MBD", "MEI"), na.rm = TRUE, ...) ## S3 method for class 'tf' tf_depth(x, arg, depth = "MBD", na.rm = TRUE, ...)
x |
|
arg |
grid of evaluation points |
depth |
currently available: "MBD", i.e. modified 2-band depth, and "MEI" |
na.rm |
TRUE remove missing observations? |
... |
further arguments handed to the function computing the respective tf_depth. |
Roughly, modified band depth computes the centrality
of a function (scale of 0 (extreme) to 0.5 (central)), while the epigraph index computes how often it is above
other functions (scale of 0 (lowest) to 1 (highest)).
The two are closely related – for functions that never cross other functions,
MBD is .
vector of tf_depth values
Sun, Ying, Genton, G M, Nychka, W D (2012). “Exact fast computation of band depth for large functional datasets: How quickly can one million curves be ranked?” Stat, 1(1), 68–74.
López-Pintado, Sara, Romo, Juan (2009). “On the concept of depth for functional data.” Journal of the American Statistical Association, 104(486), 718–734.
López-Pintado, Sara, Romo, Juan (2011). “A half-region depth for functional data.” Computational Statistics & Data Analysis, 55(4), 1679–1695.
Derivatives of tf
-objects use finite differences of the evaluations for
tfd
and finite differences of the basis functions for tfb
.
tf_derive(f, arg, order = 1, ...) ## S3 method for class 'matrix' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfd' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfb_spline' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfb_fpc' tf_derive(f, arg, order = 1, ...)
tf_derive(f, arg, order = 1, ...) ## S3 method for class 'matrix' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfd' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfb_spline' tf_derive(f, arg, order = 1, ...) ## S3 method for class 'tfb_fpc' tf_derive(f, arg, order = 1, ...)
f |
a |
arg |
grid to use for the finite differences.
Not the |
order |
order of differentiation. Maximal value for |
... |
not used |
The derivatives of tfd
objects use centered finite differences, e.g. for
first derivatives , so the domains of differentiated
tfd
will
shrink (slightly) at both ends. Unless the tfd
has a rather fine and
regular grid, representing the data in a suitable basis representation with
tfb()
and then computing the derivatives or integrals of those is usually
preferable.
Note that, for some spline bases like "cr"
or "tp"
which always begin/end
linearly, computing second derivatives will produce artefacts at the outer
limits of the functions' domain due to these boundary constraints. Basis
"bs"
does not have this problem for sufficiently high orders, but tends to
yield slightly less stable fits.
a tf
(with slightly different arg
or basis
for the
derivatives, see Details)
tf_derive(matrix)
: row-wise finite differences
tf_derive(tfd)
: derivatives by finite differencing.
tf_derive(tfb_spline)
: derivatives by finite differencing.
tf_derive(tfb_fpc)
: derivatives by finite differencing.
Other tidyfun calculus functions:
tf_integrate()
tf
-vectors for given argument valuesAlso used internally by the [
-operator for tf
data (see ?tfbrackets
) to
evaluate object
, see examples.
tf_evaluate(object, arg, ...) ## Default S3 method: tf_evaluate(object, arg, ...) ## S3 method for class 'tfd' tf_evaluate(object, arg, evaluator = tf_evaluator(object), ...) ## S3 method for class 'tfb' tf_evaluate(object, arg, ...)
tf_evaluate(object, arg, ...) ## Default S3 method: tf_evaluate(object, arg, ...) ## S3 method for class 'tfd' tf_evaluate(object, arg, evaluator = tf_evaluator(object), ...) ## S3 method for class 'tfb' tf_evaluate(object, arg, ...)
object |
a |
arg |
optional evaluation grid (vector or list of vectors).
Defaults to |
... |
not used |
evaluator |
optional. The function to use for inter/extrapolating the
|
A list of numeric vectors containing the function
evaluations on arg
.
Other tidyfun inter/extrapolation functions:
tf_approx_linear()
,
tf_interpolate()
f <- tf_rgp(3, arg = seq(0, 1, length.out = 11)) tf_evaluate(f) |> str() tf_evaluate(f, arg = 0.5) |> str() # equivalent, as matrix: f[, 0.5] new_grid <- seq(0, 1, length.out = 6) tf_evaluate(f, arg = new_grid) |> str() # equivalent, as matrix: f[, new_grid]
f <- tf_rgp(3, arg = seq(0, 1, length.out = 11)) tf_evaluate(f) |> str() tf_evaluate(f, arg = 0.5) |> str() # equivalent, as matrix: f[, 0.5] new_grid <- seq(0, 1, length.out = 6) tf_evaluate(f, arg = new_grid) |> str() # equivalent, as matrix: f[, new_grid]
Integrals of tf
-objects are computed by simple quadrature (trapezoid rule).
By default the scalar definite integral
is returned (option
definite = TRUE
),
alternatively for definite = FALSE
the anti-derivative on
[lower, upper]
, e.g. a tfd
or tfb
object representing , for
[lower, upper]
, is returned.
tf_integrate(f, arg, lower, upper, ...) ## Default S3 method: tf_integrate(f, arg, lower, upper, ...) ## S3 method for class 'tfd' tf_integrate( f, arg, lower = tf_domain(f)[1], upper = tf_domain(f)[2], definite = TRUE, ... ) ## S3 method for class 'tfb' tf_integrate( f, arg, lower = tf_domain(f)[1], upper = tf_domain(f)[2], definite = TRUE, ... )
tf_integrate(f, arg, lower, upper, ...) ## Default S3 method: tf_integrate(f, arg, lower, upper, ...) ## S3 method for class 'tfd' tf_integrate( f, arg, lower = tf_domain(f)[1], upper = tf_domain(f)[2], definite = TRUE, ... ) ## S3 method for class 'tfb' tf_integrate( f, arg, lower = tf_domain(f)[1], upper = tf_domain(f)[2], definite = TRUE, ... )
f |
a |
arg |
(optional) grid to use for the quadrature. |
lower |
lower limits of the integration range. For |
upper |
upper limits of the integration range (but see |
... |
not used |
definite |
should the definite integral be returned (default) or the antiderivative. See Description. |
For definite = TRUE
, the definite integrals of the functions in
f
. For definite = FALSE
and tf
-inputs, a tf
object containing their
anti-derivatives
Other tidyfun calculus functions:
tf_derive()
tf
-objects on a new grid of argument values.Change the internal representation of a tf
-object so that it
uses a different grid of argument values (arg
). Useful for
thinning out dense grids to make data smaller
filling out sparse grids to make derivatives/integrals and locating extrema or zero crossings more accurate (... if the interpolation works well ...)
making irregular functional data into (more) regular data.
For tfd
-objects, this is just syntactic sugar for tfd(object, arg = arg)
.
To inter/extrapolate more reliably and avoid NA
s, call
tf_interpolate
with evaluator = tf_approx_fill_extend
.
For tfb
-objects, this re-evaluates basis functions on the new grid which can
speed up subsequent computations if they all use that grid.
NB: To reliably impute very irregular data on a regular, common grid,
you'll be better off doing FPCA-based imputation or other model-based
approaches in most cases.
tf_interpolate(object, arg, ...) ## S3 method for class 'tfb' tf_interpolate(object, arg, ...) ## S3 method for class 'tfd' tf_interpolate(object, arg, ...)
tf_interpolate(object, arg, ...) ## S3 method for class 'tfb' tf_interpolate(object, arg, ...) ## S3 method for class 'tfd' tf_interpolate(object, arg, ...)
object |
an object inheriting from |
arg |
a vector of argument values on which to evaluate the functions in
|
... |
additional arguments handed over to |
a tfd
or tfb
object on the new grid given by arg
tf_rebase()
, which is more general.
Other tidyfun inter/extrapolation functions:
tf_approx_linear()
,
tf_evaluate()
# thinning out a densely observed tfd dense <- tf_rgp(10, arg = seq(0, 1, length.out = 1001)) less_dense <- tf_interpolate(dense, arg = seq(0, 1, length.out = 101)) dense less_dense # filling out sparse data (use a suitable evaluator-function!) sparse <- tf_rgp(10, arg = seq(0, 5, length.out = 11)) plot(sparse, points = TRUE) # change evaluator for better interpolation tfd(sparse, evaluator = tf_approx_spline) |> tf_interpolate(arg = seq(0, 5, length.out = 201)) |> lines(col = 2, lty = 2) set.seed(1860) sparse_irregular <- tf_rgp(5) |> tf_sparsify(0.5) |> tf_jiggle() tf_interpolate(sparse_irregular, arg = seq(0, 1, length.out = 51))
# thinning out a densely observed tfd dense <- tf_rgp(10, arg = seq(0, 1, length.out = 1001)) less_dense <- tf_interpolate(dense, arg = seq(0, 1, length.out = 101)) dense less_dense # filling out sparse data (use a suitable evaluator-function!) sparse <- tf_rgp(10, arg = seq(0, 5, length.out = 11)) plot(sparse, points = TRUE) # change evaluator for better interpolation tfd(sparse, evaluator = tf_approx_spline) |> tf_interpolate(arg = seq(0, 5, length.out = 201)) |> lines(col = 2, lty = 2) set.seed(1860) sparse_irregular <- tf_rgp(5) |> tf_sparsify(0.5) |> tf_jiggle() tf_interpolate(sparse_irregular, arg = seq(0, 1, length.out = 51))
tf
(more) irregularRandomly create some irregular functional data from regular ones.
jiggle it by randomly moving around its arg
-values. Only for tfd
.
sparsify it by setting (100*dropout
)% of its values to NA
.
tf_jiggle(f, amount = 0.4, ...) tf_sparsify(f, dropout = 0.5, ...)
tf_jiggle(f, amount = 0.4, ...) tf_sparsify(f, dropout = 0.5, ...)
f |
a |
amount |
how far away from original grid points can the new grid points lie, at most (relative to original distance to neighboring grid points). Defaults to at most 40% (0.4) of the original grid distances. Must be lower than 0.5 |
... |
additional args for the returned |
dropout |
how many values of |
an (irregular) tfd
object
Other tidyfun RNG functions:
tf_rgp()
Other tidyfun RNG functions:
tf_rgp()
tf
-objectApply the representation of one tf
-object to another; i.e. re-express it in
the other's basis, on its grid, etc.
Useful for making different functional data objects compatible so they can
be combined, compared or computed with.
tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...) ## S3 method for class 'tfd' tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...) ## S3 method for class 'tfb' tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...)
tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...) ## S3 method for class 'tfd' tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...) ## S3 method for class 'tfb' tf_rebase(object, basis_from, arg = tf_arg(basis_from), ...)
object |
a |
basis_from |
the |
arg |
optional new |
... |
forwarded to the |
This uses double dispatch (S3) internally, so the methods defined below are
themselves generics for methods tf_rebase.tfd.tfd
,
tf_rebase.tfd.tfb_spline
, tf_rebase.tfd.tfb_fpc
, tf_rebase.tfb.tfd
,
tf_rebase.tfb.tfb
that dispatch on object_from
.
a tf
-vector containing the data of object
in the same representation
as basis_from
(potentially modified by the arguments given in ...
).
tf_rebase(tfd)
: re-express a tfd
-vector in the same representation as
some other tf
-vector
tf_rebase(tfb)
: re-express a tfb
-vector in the same representation as
some other tf
-vector.
Generates n
realizations of a zero-mean Gaussian process. The function also
accepts user-defined covariance functions (without "nugget" effect, see
cov
), The implemented defaults with scale
parameter ,
order
and
nugget
effect variance are:
squared exponential covariance .
Wiener process covariance ,
Matèrn process
covariance
tf_rgp( n, arg = 51L, cov = c("squareexp", "wiener", "matern"), scale = diff(range(arg))/10, nugget = scale/200, order = 1.5 )
tf_rgp( n, arg = 51L, cov = c("squareexp", "wiener", "matern"), scale = diff(range(arg))/10, nugget = scale/200, order = 1.5 )
n |
how many realizations to draw |
arg |
vector of evaluation points ( |
cov |
type of covariance function to use. Implemented defaults are
|
scale |
scale parameter (see Description). Defaults to the width of the domain divided by 10. |
nugget |
nugget effect for additional white noise / unstructured
variability. Defaults to |
order |
order of the Matèrn covariance (if used, must be >0), defaults
to 1.5. The higher, the smoother the process. Evaluation of the covariance
function becomes numerically unstable for large (>20) |
an tfd
-vector of length n
Other tidyfun RNG functions:
tf_jiggle()
tf
objectsApply running means or medians, lowess
or Savitzky-Golay
filtering to smooth functional data. This does nothing for tfb
-objects,
which should be smoothed by using a smaller basis / stronger penalty.
tf_smooth(x, ...) ## S3 method for class 'tfb' tf_smooth(x, verbose = TRUE, ...) ## S3 method for class 'tfd' tf_smooth( x, method = c("lowess", "rollmean", "rollmedian", "savgol"), verbose = TRUE, ... )
tf_smooth(x, ...) ## S3 method for class 'tfb' tf_smooth(x, verbose = TRUE, ...) ## S3 method for class 'tfd' tf_smooth( x, method = c("lowess", "rollmean", "rollmedian", "savgol"), verbose = TRUE, ... )
x |
a |
... |
arguments for the respective |
verbose |
give lots of diagnostic messages? Defaults to TRUE |
method |
one of "lowess" (see |
tf_smooth.tfd
overrides/automatically sets some defaults of the
used methods:
lowess
uses a span parameter of f
= 0.15 (instead of 0.75)
by default.
rollmean
/median
use a window size of k
= $<$number of
grid points$>$/20 (i.e., the nearest odd integer to that) and sets fill= "extend"
(i.e., constant extrapolation to replace missing values at the
extremes of the domain) by default. Use fill= NA
for zoo
's default
behavior of shortening the smoothed series.
savgol
uses a window size of k
= $<$number of
grid points$>$/10 (i.e., the nearest odd integer to that).
a smoothed version of the input. For some methods/options, the smoothed functions may be shorter than the original ones (at both ends).
library(zoo) library(pracma) f <- tf_sparsify(tf_jiggle(tf_rgp(4, 201, nugget = 0.05))) f_lowess <- tf_smooth(f, "lowess") # these methods ignore the distances between arg-values: f_mean <- tf_smooth(f, "rollmean") f_median <- tf_smooth(f, "rollmean", k = 31) f_sg <- tf_smooth(f, "savgol", fl = 31) layout(t(1:4)) plot(f, points = FALSE, main = "original") plot(f_lowess, points = FALSE, col = "blue", main = "lowess (default,\n span 0.9 in red)" ) lines(tf_smooth(f, "lowess", f = 0.9), col = "red", alpha = 0.2) plot(f_mean, points = FALSE, col = "blue", main = "rolling means &\n medians (red)" ) lines(f_median, col = "red", alpha = 0.2) # note constant extrapolation at both ends! plot(f, points = FALSE, main = "orginal and\n savgol (red)") lines(f_sg, col = "red")
library(zoo) library(pracma) f <- tf_sparsify(tf_jiggle(tf_rgp(4, 201, nugget = 0.05))) f_lowess <- tf_smooth(f, "lowess") # these methods ignore the distances between arg-values: f_mean <- tf_smooth(f, "rollmean") f_median <- tf_smooth(f, "rollmean", k = 31) f_sg <- tf_smooth(f, "savgol", fl = 31) layout(t(1:4)) plot(f, points = FALSE, main = "original") plot(f_lowess, points = FALSE, col = "blue", main = "lowess (default,\n span 0.9 in red)" ) lines(tf_smooth(f, "lowess", f = 0.9), col = "red", alpha = 0.2) plot(f_mean, points = FALSE, col = "blue", main = "rolling means &\n medians (red)" ) lines(f_median, col = "red", alpha = 0.2) # note constant extrapolation at both ends! plot(f, points = FALSE, main = "orginal and\n savgol (red)") lines(f_sg, col = "red")
tf_where
allows to define a logical expression about the function values
and returns the argument values for which that condition is true.tf_anywhere
is syntactic sugar for tf_where
with return = "any"
to
get a logical flag for each function if the condition is TRUE
anywhere,
see below.
tf_where(f, cond, return = c("all", "first", "last", "range", "any"), arg) tf_anywhere(f, cond, arg)
tf_where(f, cond, return = c("all", "first", "last", "range", "any"), arg) tf_anywhere(f, cond, arg)
f |
a |
cond |
a logical expression about |
return |
for each entry in |
arg |
optional |
Entries in f
that do not fulfill cond
anywhere yield numeric(0)
.cond
is evaluated as a base::subset()
-statement on a data.frame
containing a single entry in f
with columns arg
and value
, so most
of the usual dplyr
tricks are available as well, see examples.
Any cond
ition evaluates to NA
on NA
-entries in f
.
depends on return
:
return = "any"
, i.e, anywhere
:
a logical vector of the same length as f
.
return = "all"
: a list of vectors of the same length as f
, with
empty vectors for the functions that never fulfill the cond
ition.
return = "range"
: a data frame with columns "begin" and "end".
else, a numeric vector of the same length as f
with NA
for entries of
f
that nowhere fulfill the cond
ition.
lin <- 1:4 * tfd(seq(-1, 1, length.out = 11), seq(-1, 1, length.out = 11)) tf_where(lin, value %inr% c(-1, 0.5)) tf_where(lin, value %inr% c(-1, 0.5), "range") a <- 1 tf_where(lin, value > a, "first") tf_where(lin, value < a, "last") tf_where(lin, value > 2, "any") tf_anywhere(lin, value > 2) set.seed(4353) f <- tf_rgp(5, 11) plot(f, pch = as.character(1:5), points = TRUE) tf_where(f, value == max(value)) # where is the function increasing/decreasing? tf_where(f, value > dplyr::lag(value, 1, value[1])) tf_where(f, value < dplyr::lead(value, 1, tail(value, 1))) # where are the (interior) extreme points (sign changes of `diff(value)`)? tf_where( f, sign(c(diff(value)[1], diff(value))) != sign(c(diff(value), tail(diff(value), 1))) ) # where in its second half is the function positive? tf_where(f, arg > 0.5 & value > 0) # does the function ever exceed? tf_anywhere(f, value > 1)
lin <- 1:4 * tfd(seq(-1, 1, length.out = 11), seq(-1, 1, length.out = 11)) tf_where(lin, value %inr% c(-1, 0.5)) tf_where(lin, value %inr% c(-1, 0.5), "range") a <- 1 tf_where(lin, value > a, "first") tf_where(lin, value < a, "last") tf_where(lin, value > 2, "any") tf_anywhere(lin, value > 2) set.seed(4353) f <- tf_rgp(5, 11) plot(f, pch = as.character(1:5), points = TRUE) tf_where(f, value == max(value)) # where is the function increasing/decreasing? tf_where(f, value > dplyr::lag(value, 1, value[1])) tf_where(f, value < dplyr::lead(value, 1, tail(value, 1))) # where are the (interior) extreme points (sign changes of `diff(value)`)? tf_where( f, sign(c(diff(value)[1], diff(value))) != sign(c(diff(value), tail(diff(value), 1))) ) # where in its second half is the function positive? tf_where(f, arg > 0.5 & value > 0) # does the function ever exceed? tf_anywhere(f, value > 1)
These are used to redefine or restrict the domain
of tf
objects.
tf_zoom(f, begin, end, ...) ## S3 method for class 'tfd' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...) ## S3 method for class 'tfb' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...) ## S3 method for class 'tfb_fpc' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...)
tf_zoom(f, begin, end, ...) ## S3 method for class 'tfd' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...) ## S3 method for class 'tfb' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...) ## S3 method for class 'tfb_fpc' tf_zoom(f, begin = tf_domain(f)[1], end = tf_domain(f)[2], ...)
f |
a |
begin |
numeric vector of length 1 or |
end |
numeric vector of length 1 or |
... |
not used |
an object like f
on a new domain (potentially).
Note that regular functional data and functions in basis representation will
be turned into irregular tfd
-objects iff begin
or end
are not scalar.
Other tidyfun utility functions:
in_range()
,
tf_arg()
x <- tf_rgp(10) plot(x) tf_zoom(x, 0.5, 0.9) tf_zoom(x, 0.5, 0.9) |> lines(col = "red") tf_zoom(x, seq(0, 0.5, length.out = 10), seq(0.5, 1, length.out = 10)) |> lines(col = "blue", lty = 3)
x <- tf_rgp(10) plot(x) tf_zoom(x, 0.5, 0.9) tf_zoom(x, 0.5, 0.9) |> lines(col = "red") tf_zoom(x, seq(0, 0.5, length.out = 10), seq(0.5, 1, length.out = 10)) |> lines(col = "blue", lty = 3)
Various constructors for tfb
-vectors from different kinds of inputs.
tfb( data = data_frame(.name_repair = "minimal"), basis = c("spline", "fpc", "wavelet"), ... ) tfb_wavelet(data, ...) as.tfb(data, basis = c("spline", "fpc"), ...)
tfb( data = data_frame(.name_repair = "minimal"), basis = c("spline", "fpc", "wavelet"), ... ) tfb_wavelet(data, ...) as.tfb(data, basis = c("spline", "fpc"), ...)
data |
a |
basis |
either " |
... |
further arguments for |
tfb
is a wrapper for functions that set up spline-, principal component- or
wavelet-based representations of functional data. For all three, the input
data are represented as weighted sums of a set of common basis
functions
identical for all observations and
weight or coefficient vectors
estimated
for each observation:
. Depending on
the value of
basis
, the basis functions will either be
spline
functions or the first few estimated eigenfunctions of the covariance
operator of the (
fpc
) or wavelets (wavelet
).
See tfb_spline()
for more details on spline basis representation (the
default). See tfb_fpc()
for using an functional principal component
representation with an orthonormal basis estimated from the data instead.
a tfb
-object (or a data.frame
/matrix
for the conversion
functions, obviously.)
Other tfb-class:
fpc_wsvd()
,
tfb_fpc()
,
tfb_spline()
Other tfb-class:
fpc_wsvd()
,
tfb_fpc()
,
tfb_spline()
These functions perform a (functional) principal component analysis (FPCA) of
the input data and return an tfb_fpc
tf
-object that uses the empirical
eigenfunctions as basis functions for representing the data. The default
("method = fpc_wsvd
") uses a (truncated) weighted SVD for complete
data on a common grid and a nuclear-norm regularized (truncated) weighted SVD
for partially missing data on a common grid, see fpc_wsvd()
.
The latter is likely to break down for high PVE and/or high amounts of
missingness.
tfb_fpc(data, ...) ## S3 method for class 'data.frame' tfb_fpc( data, id = 1, arg = 2, value = 3, domain = NULL, method = fpc_wsvd, ... ) ## S3 method for class 'matrix' tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...) ## S3 method for class 'numeric' tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...) ## S3 method for class 'tf' tfb_fpc(data, arg = NULL, method = fpc_wsvd, ...) ## Default S3 method: tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...)
tfb_fpc(data, ...) ## S3 method for class 'data.frame' tfb_fpc( data, id = 1, arg = 2, value = 3, domain = NULL, method = fpc_wsvd, ... ) ## S3 method for class 'matrix' tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...) ## S3 method for class 'numeric' tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...) ## S3 method for class 'tf' tfb_fpc(data, arg = NULL, method = fpc_wsvd, ...) ## Default S3 method: tfb_fpc(data, arg = NULL, domain = NULL, method = fpc_wsvd, ...)
data |
a |
... |
arguments to the |
id |
The name or number of the column defining which data belong to which function. |
arg |
|
value |
The name or number of the column containing the function evaluations. |
domain |
range of the |
method |
the function to use that computes eigenfunctions and scores.
Defaults to |
For the FPC basis, any factorization method that accepts a data.frame
with
columns id
, arg
, value
containing the functional data and returns a
list with eigenfunctions and FPC scores structured like the return object
of fpc_wsvd()
can be used for the 'method“ argument, see example below.
Note that the mean function, with a fixed "score" of 1 for all functions,
is used as the first basis function for all FPC bases.
an object of class tfb_fpc
, inheriting from tfb
.
The basis used by tfb_fpc
is a tfd
-vector containing the estimated
mean and eigenfunctions.
tfb_fpc(default)
: convert tfb
: default method, returning prototype when
data is NULL
fpc_wsvd()
for FPCA options.
Other tfb-class:
fpc_wsvd()
,
tfb
,
tfb_spline()
Other tfb_fpc-class:
fpc_wsvd()
set.seed(13121) x <- tf_rgp(25, nugget = .02) x_pc <- tfb_fpc(x, pve = .9) x_pc plot(x, lwd = 3) lines(x_pc, col = 2, lty = 2) x_pc_full <- tfb_fpc(x, pve = .995) x_pc_full lines(x_pc_full, col = 3, lty = 2) # partially missing data on common grid: x_mis <- x |> tf_sparsify(dropout = .05) x_pc_mis <- tfb_fpc(x_mis, pve = .9) x_pc_mis plot(x_mis, lwd = 3) lines(x_pc_mis, col = 4, lty = 2) # extract FPC basis -- # first "eigenvector" in black is (always) the mean function x_pc |> tf_basis(as_tfd = TRUE) |> plot(col = 1:5) # Apply FPCA for sparse, irregular data using refund::fpca.sc: set.seed(99290) # create small, sparse, irregular data: x_irreg <- x[1:8] |> tf_jiggle() |> tf_sparsify(dropout = 0.3) plot(x_irreg) x_df <- x_irreg |> as.data.frame(unnest = TRUE) # wrap refund::fpca_sc for use as FPCA method in tfb_fpc -- # 1. define scoring function (simple weighted LS fit) fpca_scores <- function(data_matrix, efunctions, mean, weights) { w_mat <- matrix(weights, ncol = length(weights), nrow = nrow(data_matrix), byrow = TRUE) w_mat[is.na(data_matrix)] <- 0 data_matrix[is.na(data_matrix)] <- 0 data_wc <- t((t(data_matrix) - mean) * sqrt(t(w_mat))) t(qr.coef(qr(efunctions), t(data_wc) / sqrt(weights))) } # 2. define wrapper for fpca_sc: fpca_sc_wrapper <- function(data, arg, pve = 0.995, ...) { data_mat <- tfd(data) |> as.matrix(interpolate = TRUE) fpca <- refund::fpca.sc( Y = data_mat, argvals = attr(data_mat, "arg"), pve = pve, ... ) c(fpca[c("mu", "efunctions", "scores", "npc")], scoring_function = fpca_scores) } x_pc <- tfb_fpc(x_df, method = fpca_sc_wrapper) lines(x_pc, col = 2, lty = 2)
set.seed(13121) x <- tf_rgp(25, nugget = .02) x_pc <- tfb_fpc(x, pve = .9) x_pc plot(x, lwd = 3) lines(x_pc, col = 2, lty = 2) x_pc_full <- tfb_fpc(x, pve = .995) x_pc_full lines(x_pc_full, col = 3, lty = 2) # partially missing data on common grid: x_mis <- x |> tf_sparsify(dropout = .05) x_pc_mis <- tfb_fpc(x_mis, pve = .9) x_pc_mis plot(x_mis, lwd = 3) lines(x_pc_mis, col = 4, lty = 2) # extract FPC basis -- # first "eigenvector" in black is (always) the mean function x_pc |> tf_basis(as_tfd = TRUE) |> plot(col = 1:5) # Apply FPCA for sparse, irregular data using refund::fpca.sc: set.seed(99290) # create small, sparse, irregular data: x_irreg <- x[1:8] |> tf_jiggle() |> tf_sparsify(dropout = 0.3) plot(x_irreg) x_df <- x_irreg |> as.data.frame(unnest = TRUE) # wrap refund::fpca_sc for use as FPCA method in tfb_fpc -- # 1. define scoring function (simple weighted LS fit) fpca_scores <- function(data_matrix, efunctions, mean, weights) { w_mat <- matrix(weights, ncol = length(weights), nrow = nrow(data_matrix), byrow = TRUE) w_mat[is.na(data_matrix)] <- 0 data_matrix[is.na(data_matrix)] <- 0 data_wc <- t((t(data_matrix) - mean) * sqrt(t(w_mat))) t(qr.coef(qr(efunctions), t(data_wc) / sqrt(weights))) } # 2. define wrapper for fpca_sc: fpca_sc_wrapper <- function(data, arg, pve = 0.995, ...) { data_mat <- tfd(data) |> as.matrix(interpolate = TRUE) fpca <- refund::fpca.sc( Y = data_mat, argvals = attr(data_mat, "arg"), pve = pve, ... ) c(fpca[c("mu", "efunctions", "scores", "npc")], scoring_function = fpca_scores) } x_pc <- tfb_fpc(x_df, method = fpca_sc_wrapper) lines(x_pc, col = 2, lty = 2)
Represent curves as a weighted sum of spline basis functions.
tfb_spline(data, ...) ## S3 method for class 'data.frame' tfb_spline( data, id = 1, arg = 2, value = 3, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'matrix' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'numeric' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'list' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'tfd' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'tfb' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## Default S3 method: tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... )
tfb_spline(data, ...) ## S3 method for class 'data.frame' tfb_spline( data, id = 1, arg = 2, value = 3, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'matrix' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'numeric' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'list' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'tfd' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## S3 method for class 'tfb' tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... ) ## Default S3 method: tfb_spline( data, arg = NULL, domain = NULL, penalized = TRUE, global = FALSE, verbose = TRUE, ... )
data |
a |
... |
arguments to the calls to |
id |
The name or number of the column defining which data belong to which function. |
arg |
|
value |
The name or number of the column containing the function evaluations. |
domain |
range of the |
penalized |
|
global |
Defaults to |
verbose |
|
The basis to be used is set up via a call to mgcv::s()
and all the spline
bases discussed in mgcv::smooth.terms()
are available, in principle.
Depending on the value of the penalized
- and global
-flags, the
coefficient vectors for each observation are then estimated via fitting a GAM
(separately for each observation, if !global
) via mgcv::magic()
(least
square error, the default) or mgcv::gam()
(if a family
argument was
supplied) or unpenalized least squares / maximum likelihood.
After the "smoothed" representation is computed, the amount of smoothing that
was performed is reported in terms of the "percentage of variability
preserved", which is the variance (or the explained deviance, in the general
case if family
was specified) of the smoothed function values divided by the variance of the original
values (the null deviance, in the general case). Reporting can be switched off
with verbose = FALSE
.
The ...
arguments supplies arguments to both the
spline basis (via mgcv::s()
) and the estimation (via
mgcv::magic()
or mgcv::gam()
), the most important arguments are:
k
: how many basis functions should the spline basis use, default is 25.
bs
: which type of spline basis should be used, the default is cubic
regression splines (bs = "cr"
)
family
argument: use this if minimizing squared errors is not
a reasonable criterion for the representation accuracy (see
mgcv::family.mgcv()
for what's available) and/or if function values are
restricted to be e.g. positive (family = Gamma()/tw()/...
), in
(
family = betar()
), etc.
sp
: numeric value for the smoothness penalty weight, for manually
setting the amount of smoothing for all curves, see mgcv::s()
. This
(drastically) reduces computation time. Defaults to -1
, i.e., automatic
optimization of sp
using mgcv::magic()
(LS fits) or mgcv::gam()
(GLM),
source code in R/tfb-spline-utils.R
.
If global == TRUE
, this uses a small subset of curves (10%
of curves,
at least 5, at most 100; non-random sample using every j-th curve in the
data) on which smoothing parameters per curve are estimated and then takes
the mean of the log smoothing parameter of those as sp
for all curves. This
is much faster than optimizing for each curve on large data sets. For very
sparse or noisy curves, estimating a common smoothing parameter based on the
data for all curves simultaneously is likely to yield better results, this is
not what's implemented here.
a tfb
-object
tfb_spline(data.frame)
: convert data frames
tfb_spline(matrix)
: convert matrices
tfb_spline(numeric)
: convert matrices
tfb_spline(list)
: convert lists
tfb_spline(tfd)
: convert tfd
(raw functional data)
tfb_spline(tfb)
: convert tfb
: modify basis representation, smoothing.
tfb_spline(default)
: convert tfb
: default method, returning prototype
when data is missing
mgcv::smooth.terms()
for spline basis options.
Other tfb-class:
fpc_wsvd()
,
tfb
,
tfb_fpc()
tf
vectorsThese functions access, subset, replace and evaluate tf
objects.
For more information on creating tf
objects and converting them to/from
list
, data.frame
or matrix
, see tfd()
and tfb()
. See Details.
## S3 method for class 'tf' x[i, j, interpolate = TRUE, matrix = TRUE] ## S3 replacement method for class 'tf' x[i] <- value
## S3 method for class 'tf' x[i, j, interpolate = TRUE, matrix = TRUE] ## S3 replacement method for class 'tf' x[i] <- value
x |
an |
i |
index of the observations ( |
j |
The |
interpolate |
should functions be evaluated (i.e., inter-/extrapolated)
for values in |
matrix |
should the result be returned as a |
value |
|
Note that these break certain (terrible) R conventions for vector-like objects:
no argument recycling,
no indexing with NA
,
no indexing with names not present in x
,
no indexing with integers > length(x)
All of the above will trigger errors.
If j
is missing, a subset of the functions in x
as given by
i
.
If j
is given and matrix == TRUE
, a numeric matrix of function
evaluations in which each row represents one function and each column
represents one argval
as given in argument j
, with an attribute
arg
=j
and row- and column-names derived from x[i]
and j
.
If
j
is given and matrix == FALSE
, a list of tbl_df
s with columns
arg
= j
and value
= evaluations at j
for each observation in
i
.
x <- 1:3 * tfd(data = 0:10, arg = 0:10) plot(x) # this operator's 2nd argument is quite overloaded -- you can: # 1. simply extract elements from the vector if no second arg is given: x[1] x[c(TRUE, FALSE, FALSE)] x[-(2:3)] # 2. use the second argument and optional additional arguments to # extract specific function evaluations in a number of formats: x[1:2, c(4.5, 9)] # returns a matrix of function evaluations x[1:2, c(4.5, 9), interpolate = FALSE] # NA for arg-values not in the original data x[-3, seq(1, 9, by = 2), matrix = FALSE] # list of data.frames for each function # in order to evaluate a set of observed functions on a new grid and # save them as a functional data vector again, use `tfd` or `tfb` instead: tfd(x, arg = seq(0, 10, by = 0.01))
x <- 1:3 * tfd(data = 0:10, arg = 0:10) plot(x) # this operator's 2nd argument is quite overloaded -- you can: # 1. simply extract elements from the vector if no second arg is given: x[1] x[c(TRUE, FALSE, FALSE)] x[-(2:3)] # 2. use the second argument and optional additional arguments to # extract specific function evaluations in a number of formats: x[1:2, c(4.5, 9)] # returns a matrix of function evaluations x[1:2, c(4.5, 9), interpolate = FALSE] # NA for arg-values not in the original data x[-3, seq(1, 9, by = 2), matrix = FALSE] # list of data.frames for each function # in order to evaluate a set of observed functions on a new grid and # save them as a functional data vector again, use `tfd` or `tfb` instead: tfd(x, arg = seq(0, 10, by = 0.01))
Various constructor methods for tfd
-objects.
tfd.matrix
accepts a numeric matrix with one function per
row (!). If arg
is not provided, it tries to guess arg
from the
column names and falls back on 1:ncol(data)
if that fails.
tfd.data.frame
uses the first 3 columns of data
for
function information by default: (id
, arg
, value
)
tfd.list
accepts a list of vectors of identical lengths
containing evaluations or a list of 2-column matrices/data.frames with
arg
in the first and evaluations in the second column
tfd.default
returns class prototype when argument to tfd() is
NULL or not a recognised class
tfd(data, ...) ## S3 method for class 'matrix' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'numeric' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'data.frame' tfd( data, id = 1, arg = 2, value = 3, domain = NULL, evaluator = tf_approx_linear, ... ) ## S3 method for class 'list' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'tf' tfd(data, arg = NULL, domain = NULL, evaluator = NULL, ...) ## Default S3 method: tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) as.tfd(data, ...) as.tfd_irreg(data, ...)
tfd(data, ...) ## S3 method for class 'matrix' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'numeric' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'data.frame' tfd( data, id = 1, arg = 2, value = 3, domain = NULL, evaluator = tf_approx_linear, ... ) ## S3 method for class 'list' tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) ## S3 method for class 'tf' tfd(data, arg = NULL, domain = NULL, evaluator = NULL, ...) ## Default S3 method: tfd(data, arg = NULL, domain = NULL, evaluator = tf_approx_linear, ...) as.tfd(data, ...) as.tfd_irreg(data, ...)
data |
a |
... |
not used in |
arg |
|
domain |
range of the |
evaluator |
a function accepting arguments |
id |
The name or number of the column defining which data belong to which function. |
value |
The name or number of the column containing the function evaluations. |
evaluator
: must be the (quoted or bare) name of a
function with signature function(x, arg, evaluations)
that returns
the functions' (approximated/interpolated) values at locations x
based on
the function evaluations
available at locations arg
.
Available evaluator
-functions:
tf_approx_linear
for linear interpolation without extrapolation (i.e.,
zoo::na.approx()
with na.rm = FALSE
) – this is the default,
tf_approx_spline
for cubic spline interpolation, (i.e., zoo::na.spline()
with na.rm = FALSE
),
tf_approx_fill_extend
for linear interpolation and constant extrapolation
(i.e., zoo::na.fill()
with fill = "extend"
)
tf_approx_locf
for "last observation carried forward" (i.e.,
zoo::na.locf()
with na.rm = FALSE
and
tf_approx_nocb
for "next observation carried backward" (i.e.,
zoo::na.locf()
with na.rm = FALSE, fromLast = TRUE
).
See tf:::zoo_wrapper
and tf:::tf_approx_linear
, which is simply
zoo_wrapper(zoo::na.tf_approx, na.rm = FALSE)
, for examples of
implementations of this.
an tfd
-object (or a data.frame
/matrix
for the conversion
functions, obviously.)
# turn irregular to regular tfd by evaluating on a common grid: f <- c( tf_rgp(1, arg = seq(0, 1, length.out = 11)), tf_rgp(1, arg = seq(0, 1, length.out = 21)) ) tfd(f, arg = seq(0, 1, length.out = 21)) set.seed(1213) f <- tf_rgp(3, arg = seq(0, 1, length.out = 51)) |> tf_sparsify(0.9) # does not yield regular data because linear extrapolation yields NAs # outside observed range: tfd(f, arg = seq(0, 1, length.out = 101)) # this "works" (but may not yield sensible values..!!) for # e.g. constant extrapolation: tfd(f, evaluator = tf_approx_fill_extend, arg = seq(0, 1, length.out = 101)) plot(f, col = 2) tfd(f, arg = seq(0, 1, length.out = 151), evaluator = tf_approx_fill_extend ) |> lines()
# turn irregular to regular tfd by evaluating on a common grid: f <- c( tf_rgp(1, arg = seq(0, 1, length.out = 11)), tf_rgp(1, arg = seq(0, 1, length.out = 21)) ) tfd(f, arg = seq(0, 1, length.out = 21)) set.seed(1213) f <- tf_rgp(3, arg = seq(0, 1, length.out = 51)) |> tf_sparsify(0.9) # does not yield regular data because linear extrapolation yields NAs # outside observed range: tfd(f, arg = seq(0, 1, length.out = 101)) # this "works" (but may not yield sensible values..!!) for # e.g. constant extrapolation: tfd(f, evaluator = tf_approx_fill_extend, arg = seq(0, 1, length.out = 101)) plot(f, col = 2) tfd(f, arg = seq(0, 1, length.out = 151), evaluator = tf_approx_fill_extend ) |> lines()
tf
These methods and operators mostly work arg
-value-wise on tf
objects, see
vctrs::vec_arith()
etc. for implementation details.
## S3 method for class 'tfd' e1 == e2 ## S3 method for class 'tfd' e1 != e2 ## S3 method for class 'tfb' e1 == e2 ## S3 method for class 'tfb' e1 != e2 ## S3 method for class 'tfd' vec_arith(op, x, y, ...) ## S3 method for class 'tfb' vec_arith(op, x, y, ...) ## S3 method for class 'tfd' Math(x, ...) ## S3 method for class 'tfb' Math(x, ...) ## S3 method for class 'tf' Summary(...) ## S3 method for class 'tfd' cummax(...) ## S3 method for class 'tfd' cummin(...) ## S3 method for class 'tfd' cumsum(...) ## S3 method for class 'tfd' cumprod(...) ## S3 method for class 'tfb' cummax(...) ## S3 method for class 'tfb' cummin(...) ## S3 method for class 'tfb' cumsum(...) ## S3 method for class 'tfb' cumprod(...)
## S3 method for class 'tfd' e1 == e2 ## S3 method for class 'tfd' e1 != e2 ## S3 method for class 'tfb' e1 == e2 ## S3 method for class 'tfb' e1 != e2 ## S3 method for class 'tfd' vec_arith(op, x, y, ...) ## S3 method for class 'tfb' vec_arith(op, x, y, ...) ## S3 method for class 'tfd' Math(x, ...) ## S3 method for class 'tfb' Math(x, ...) ## S3 method for class 'tf' Summary(...) ## S3 method for class 'tfd' cummax(...) ## S3 method for class 'tfd' cummin(...) ## S3 method for class 'tfd' cumsum(...) ## S3 method for class 'tfd' cumprod(...) ## S3 method for class 'tfb' cummax(...) ## S3 method for class 'tfb' cummin(...) ## S3 method for class 'tfb' cumsum(...) ## S3 method for class 'tfb' cumprod(...)
e1 |
an |
e2 |
an |
op |
An arithmetic operator as a string |
x |
a |
y |
a |
... |
|
See examples below. Equality checks of functional objects are even more iffy
than usual for computer math and not very reliable. Note that max
and min
are not guaranteed to be maximal/minimal over the entire domain, only on the
evaluation grid used for computation. With the exception of addition and
multiplication, operations on tfb
-objects first evaluate the data on their
arg
, perform computations on these evaluations and then convert back to an
tfb
- object, so a loss of precision should be expected – especially so for
small spline bases and/or very wiggly data.
a tf
- or logical
vector with the computed result
tf_fwise()
for scalar summaries of each function in a tf
-vector
set.seed(1859) f <- tf_rgp(4) 2 * f == f + f sum(f) == f[1] + f[2] + f[3] + f[4] log(exp(f)) == f plot(f, points = FALSE) lines(range(f), col = 2, lty = 2) f2 <- tf_rgp(5) |> exp() |> tfb(k = 25) layout(t(1:3)) plot(f2, col = gray.colors(5)) plot(cummin(f2), col = gray.colors(5)) plot(cumsum(f2), col = gray.colors(5)) # ?tf_integrate for integrals, ?tf_fwise for scalar summaries of each function
set.seed(1859) f <- tf_rgp(4) 2 * f == f + f sum(f) == f[1] + f[2] + f[3] + f[4] log(exp(f)) == f plot(f, points = FALSE) lines(range(f), col = 2, lty = 2) f2 <- tf_rgp(5) |> exp() |> tfb(k = 25) layout(t(1:3)) plot(f2, col = gray.colors(5)) plot(cummin(f2), col = gray.colors(5)) plot(cumsum(f2), col = gray.colors(5)) # ?tf_integrate for integrals, ?tf_fwise for scalar summaries of each function
tf
objects across argument valuesThese will return a tf
object containing the respective functional
statistic. See tf_fwise()
for scalar summaries
(e.g. tf_fmean
for means, tf_fmax
for max. values) of each entry
in a tf
-vector.
## S3 method for class 'tf' mean(x, ...) ## S3 method for class 'tf' median(x, na.rm = FALSE, depth = c("MBD", "pointwise"), ...) sd(x, na.rm = FALSE) ## Default S3 method: sd(x, na.rm = FALSE) ## S3 method for class 'tf' sd(x, na.rm = FALSE) var(x, y = NULL, na.rm = FALSE, use) ## Default S3 method: var(x, y = NULL, na.rm = FALSE, use) ## S3 method for class 'tf' var(x, y = NULL, na.rm = FALSE, use) ## S3 method for class 'tf' summary(object, ...)
## S3 method for class 'tf' mean(x, ...) ## S3 method for class 'tf' median(x, na.rm = FALSE, depth = c("MBD", "pointwise"), ...) sd(x, na.rm = FALSE) ## Default S3 method: sd(x, na.rm = FALSE) ## S3 method for class 'tf' sd(x, na.rm = FALSE) var(x, y = NULL, na.rm = FALSE, use) ## Default S3 method: var(x, y = NULL, na.rm = FALSE, use) ## S3 method for class 'tf' var(x, y = NULL, na.rm = FALSE, use) ## S3 method for class 'tf' summary(object, ...)
x |
a |
... |
optional additional arguments. |
na.rm |
logical. Should missing values be removed? |
depth |
method used to determine the most central element in |
y |
|
use |
an optional character string giving a
method for computing covariances in the presence
of missing values. This must be (an abbreviation of) one of the strings
|
object |
a |
a tf
object with the computed result.summary.tf
returns a tf
-vector with the mean function, the
variance function, the functional median, and the functional range
(i.e., the pointwise min/max) of the central half of the functions,
as defined by tf_depth()
.
Other tidyfun summary functions:
functionwise
See above.
unique_id(x)
unique_id(x)
x |
any input |
x
turned into a list.
Other tidyfun developer tools:
ensure_list()
,
prep_plotting_arg()
vctrs
methods for tf
objectsThese functions are the extensions that allow tf
vectors
to work with vctrs
.
## S3 method for class 'tfd_reg.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfd_irreg' vec_ptype2(x, y, ...)
## S3 method for class 'tfd_reg.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfd_irreg.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfb_spline' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfb_fpc' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfb_spline.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfd_reg' vec_cast(x, to, ...) ## S3 method for class 'tfb_fpc.tfd_irreg' vec_cast(x, to, ...) ## S3 method for class 'tfd_reg.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_reg.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfd_irreg.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_spline.tfd_irreg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfb_spline' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfb_fpc' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfd_reg' vec_ptype2(x, y, ...) ## S3 method for class 'tfb_fpc.tfd_irreg' vec_ptype2(x, y, ...)
x |
Vectors to cast. |
to |
Type to cast to. If |
... |
For |
y |
Vectors to cast. |
Notes on vec_cast
:
Use tf_rebase()
to change the representations of tf
-vectors,
these methods are only for internal use –
automatic/implicit casting of tf
objects is tricky
because it's hard to determine automatically whether such an operation would
lose precision (different bases with different expressivity? different
argument grids?), and it's not generally clear which instances of which
tf
-subclasses should be considered the "richer" objects.
Rules for casting:
If the casted object's domain
would not contain the entire original domain
,
no casting is possible (would lose data).
Every cast that evaluates (basis) functions on different arg
values is a lossy cast,
since it might lose precision (vctrs::maybe_lossy_cast
).
As long as the casted object's domain
contains the entire original domain
:
every tfd_reg
, tfd_irreg
or tfb
can always be cast into an equivalent
tfd_irreg
(which may also change its evaluator
and domain
).
every tfd_reg
can always be cast to tfd_reg
(which may change its evaluator
and domain
)
every tfb
can be cast losslessly to tfd
(regular or irregular,
note it's lossless only on the original arg
-grid)
Any cast of a tfd
into tfb
is potentially lossy (because we don't know how expressive the chosen basis is)
Only tfb
with identical bases and domains can be cast into one another losslessly
for vec_cast
: the casted tf
-vector, for vec_ptype2
: the common prototype
vctrs::vec_cast()
, vctrs::vec_ptype2()