Title: | Time-Frequency analysis of 1-D signals |
---|---|
Description: | Rwave is a library of R functions which provide an environment for the Time-Frequency analysis of 1-D signals (and especially for the wavelet and Gabor transforms of noisy signals). It was originally written for Splus by Rene Carmona, Bruno Torresani, and Wen L. Hwang, first at the University of California at Irvine and then at Princeton University. Credit should also be given to Andrea Wang whose functions on the dyadic wavelet transform are included. Rwave is based on the book: "Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S", by Rene Carmona, Wen L. Hwang and Bruno Torresani, Academic Press, 1998. This package is no longer actively maintained. A C++ rewrite of core functionality is in progress. If you'd like to participate, please contact Christian Gunning. |
Authors: | S original by Rene Carmona <[email protected]> and Bruno Torresani <[email protected]>; R port by Brandon Whitcher <[email protected]> |
Maintainer: | Brandon Whitcher <[email protected]> and Christian Gunning <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.25-2 |
Built: | 2024-11-21 05:20:11 UTC |
Source: | https://github.com/r-forge/rwave |
Transient signal.
data(A0)
data(A0)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(A4)
data(A4)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Add zeros to the end of the data if necessary so that its length is a power of 2. It returns the data with zeros added if nessary and the length of the adjusted data.
adjust.length(inputdata)
adjust.length(inputdata)
inputdata |
either a text file or an S object containing data. |
Zero-padded 1D array.
See discussions in the text of “Practical Time-Frequency Analysis”.
Pixel from amber camara.
data(amber7)
data(amber7)
A vector containing 7000 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Pixel from amber camara.
data(amber8)
data(amber8)
A vector containing 7000 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Pixel from amber camara.
data(amber9)
data(amber9)
A vector containing 7000 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(B0)
data(B0)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(B4)
data(B4)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Acoustic returns from natural underwater clutter.
data(back1.000)
data(back1.000)
A vector containing 7936 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Acoustic returns from ...
data(back1.180)
data(back1.180)
A vector containing 7936 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Acoustic returns from an underwater metallic object.
data(back1.220)
data(back1.220)
A vector containing 7936 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(C0)
data(C0)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(C4)
data(C4)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Chains the ridge estimates produced by the function crc
.
cfamily(ccridge, bstep=1, nbchain=100, ptile=0.05)
cfamily(ccridge, bstep=1, nbchain=100, ptile=0.05)
ccridge |
unchained ridge set as the output of the function |
bstep |
maximal length for a gap in a ridge. |
nbchain |
maximal number of chains produced by the function. |
ptile |
relative threshold for the ridges. |
crc
returns a measure in time-frequency (or time-scale)
space. cfamily
turns it into a series of one-dimensional
objects (ridges). The measure is first thresholded, with a relative
threshold value set to the input parameter ptile. During the chaining
procedure, gaps within a given ridge are allowed and filled in. The
maximal length of such gaps is the input parameter bstep.
Returns the results of the chaining algorithm
ordered map |
image containing the ridges (displayed with different colors) |
chain |
2D array containing the chained ridges, according to the chain data
structure |
nbchain |
number of chains produced by the algorithm |
See discussion in text of “Practical Time-Frequency Analysis”.
crc
for the ridge estimation, and crcrec
,
gcrcrec
and scrcrec
for corresponding
reconstruction functions.
Computes the continuous Gabor transform with Gaussian window.
cgt(input, nvoice, freqstep=(1/nvoice), scale=1, plot=TRUE)
cgt(input, nvoice, freqstep=(1/nvoice), scale=1, plot=TRUE)
input |
input signal (possibly complex-valued). |
nvoice |
number of frequencies for which gabor transform is to be computed. |
freqstep |
Sampling rate for the frequency axis. |
scale |
Size parameter for the window. |
plot |
logical variable set to TRUE to display the modulus of the continuous gabor transform on the graphic device. |
The output contains the (complex) values of the gabor transform of the input signal. The format of the output is a 2D array (signal\_size x nb\_scales).
continuous (complex) gabor transform (2D array).
freqstep must be less than 1/nvoice to avoid aliasing. freqstep=1/nvoice corresponds to the Nyquist limit.
See discussion in text of “Practical Time-Frequency Analysis”.
cwt
, cwtp
, DOG
for continuous wavelet transforms.
cwtsquiz
for synchrosqueezed wavelet transform.
Chen's chirp.
data(ch)
data(ch)
A vector containing 15,000 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Stop when is larger than the signal size.
check.maxresoln(maxresoln, np)
check.maxresoln(maxresoln, np)
maxresoln |
number of decomposition scales. |
np |
signal size. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Sets to zero the phase of time-frequency transform when modulus is below a certain value.
cleanph(tfrep, thresh=0.01, plot=TRUE)
cleanph(tfrep, thresh=0.01, plot=TRUE)
tfrep |
continuous time-frequency transform (2D array) |
thresh |
(relative) threshold. |
plot |
if set to TRUE, displays the maxima of cwt on the graphic device. |
thresholded phase (2D array)
See discussion in text of “Practical Time-Frequency Analysis”.
Dolphin click data.
data(click)
data(click)
A vector containing 2499 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Estimate a (single) ridge from a time-frequency representation, using the corona method.
corona(tfrep, guess, tfspec=numeric(dim(tfrep)[2]), subrate=1, temprate=3, mu=1, lambda=2 * mu, iteration=1000000, seed=-7, stagnant=20000, costsub=1, plot=TRUE)
corona(tfrep, guess, tfspec=numeric(dim(tfrep)[2]), subrate=1, temprate=3, mu=1, lambda=2 * mu, iteration=1000000, seed=-7, stagnant=20000, costsub=1, plot=TRUE)
tfrep |
Time-Frequency representation (real valued). |
guess |
Initial guess for the algorithm. |
tfspec |
Estimate for the contribution of the noise to modulus. |
subrate |
Subsampling rate for ridge estimation. |
temprate |
Initial value of temperature parameter. |
mu |
Coefficient of the ridge's second derivative in cost function. |
lambda |
Coefficient of the ridge's derivative in cost function. |
iteration |
Maximal number of moves. |
seed |
Initialization of random number generator. |
stagnant |
Maximum number of stationary iterations before stopping. |
costsub |
Subsampling of cost function in output. |
plot |
When set(default), some results will be shown on the display. |
To accelerate convergence, it is useful to preprocess modulus before
running annealing method. Such a preprocessing (smoothing and
subsampling of modulus) is implemented in corona
. The
parameter subrate specifies the subsampling rate.
Returns the estimated ridge and the cost function.
ridge |
1D array (of same length as the signal) containing the ridge. |
cost |
1D array containing the cost function. |
The returned cost may be a large array, which is time consuming. The argument costsub allows subsampling the cost function.
See discussion in text of “Practical Time-Frequency Analysis”.
Estimate a ridge using the modified corona method (modified cost function).
coronoid(tfrep, guess, tfspec=numeric(dim(tfrep)[2]), subrate=1, temprate=3, mu=1, lambda=2 * mu, iteration=1000000, seed=-7, stagnant=20000, costsub=1, plot=TRUE)
coronoid(tfrep, guess, tfspec=numeric(dim(tfrep)[2]), subrate=1, temprate=3, mu=1, lambda=2 * mu, iteration=1000000, seed=-7, stagnant=20000, costsub=1, plot=TRUE)
tfrep |
Estimate for the contribution of the noise to modulus. |
guess |
Initial guess for the algorithm. |
tfspec |
Estimate for the contribution of the noise to modulus. |
subrate |
Subsampling rate for ridge estimation. |
temprate |
Initial value of temperature parameter. |
mu |
Coefficient of the ridge's derivative in cost function. |
lambda |
Coefficient of the ridge's second derivative in cost function. |
iteration |
Maximal number of moves. |
seed |
Initialization of random number generator. |
stagnant |
Maximum number of stationary iterations before stopping. |
costsub |
Subsampling of cost function in output. |
plot |
When set(default), some results will be shown on the display. |
To accelerate convergence, it is useful to preprocess modulus before
running annealing method. Such a preprocessing (smoothing and
subsampling of modulus) is implemented in coronoid
. The
parameter subrate specifies the subsampling rate.
Returns the estimated ridge and the cost function.
ridge |
1D array (of same length as the signal) containing the ridge. |
cost |
1D array containing the cost function. |
The returned cost may be a large array. The argument costsub allows subsampling the cost function.
See discussion in text of “Practical Time-Frequency Analysis”.
Uses the "crazy climber algorithm" to detect ridges in the modulus of a continuous wavelet or a Gabor transform.
crc(tfrep, tfspec=numeric(dim(tfrep)[2]), bstep=3, iteration=10000, rate=0.001, seed=-7, nbclimb=10, flag.int=TRUE, chain=TRUE, flag.temp=FALSE)
crc(tfrep, tfspec=numeric(dim(tfrep)[2]), bstep=3, iteration=10000, rate=0.001, seed=-7, nbclimb=10, flag.int=TRUE, chain=TRUE, flag.temp=FALSE)
tfrep |
modulus of the (wavelet or Gabor) transform. |
tfspec |
numeric vector which gives, for each value of the scale or frequency the expected size of the noise contribution. |
bstep |
stepsize for random walk of the climbers. |
iteration |
number of iterations. |
rate |
initial value of the temperature. |
seed |
initial value of the random number generator. |
nbclimb |
number of crazy climbers. |
flag.int |
if set to TRUE, the weighted occupation measure is computed. |
chain |
if set to TRUE, chaining of the ridges is done. |
flag.temp |
if set to TRUE: constant temperature. |
Returns a 2D array called beemap containing the (weighted or unweighted) occupation measure (integrated with respect to time)
See discussion in text of “Practical Time-Frequency Analysis”.
corona
, icm
, coronoid
,
snake
, snakoid
for ridge estimation,
cfamily
for chaining and
crcrec
,gcrcrec
,scrcrec
for
reconstruction.
Reconstructs a real valued signal from the output of crc
(wavelet case) by minimizing an appropriate quadratic form.
crcrec(siginput, inputwt, beemap, noct, nvoice, compr, minnbnodes=2, w0=2 * pi, bstep=5, ptile=0.01, epsilon=0, fast=FALSE, para=5, real=FALSE, plot=2)
crcrec(siginput, inputwt, beemap, noct, nvoice, compr, minnbnodes=2, w0=2 * pi, bstep=5, ptile=0.01, epsilon=0, fast=FALSE, para=5, real=FALSE, plot=2)
siginput |
original signal. |
inputwt |
wavelet transform. |
beemap |
occupation measure, output of |
noct |
number of octaves. |
nvoice |
number of voices per octave. |
compr |
compression rate for sampling the ridges. |
minnbnodes |
minimal number of points per ridge. |
w0 |
center frequency of the wavelet. |
bstep |
size (in the time direction) of the steps for chaining. |
ptile |
relative threshold of occupation measure. |
epsilon |
constant in front of the smoothness term in penalty function. |
fast |
if set to TRUE, uses trapezoidal rule to evaluate $Q_2$. |
para |
scale parameter for extrapolating the ridges. |
real |
if set to TRUE, uses only real constraints. |
plot |
1: displays signal,components,and reconstruction one after another. 2: displays signal, components and reconstruction. |
When ptile is high, boundary effects may appeare. para controls extrapolation of the ridge.
Returns a structure containing the following elements:
rec |
reconstructed signal. |
ordered |
image of the ridges (with different colors). |
comp |
2D array containing the signals reconstructed from ridges. |
displays a family of chained ridges, output of cfamily
.
crfview(beemap, twod=TRUE)
crfview(beemap, twod=TRUE)
beemap |
Family of chained ridges, output of |
twod |
If set to T, displays the ridges as an image. If set to F, displays as a series of curves. |
See discussions in the text of “Practical Time-Frequency Analysis”.
crc
,cfamily
for crazy climbers and corresponding chaining algorithms.
Computes the continuous wavelet transform with for the (complex-valued) Morlet wavelet.
cwt(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
cwt(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
input |
input signal (possibly complex-valued) |
noctave |
number of powers of 2 for the scale variable |
nvoice |
number of scales in each octave (i.e. between two consecutive powers of 2). |
w0 |
central frequency of the wavelet. |
twoD |
logical variable set to T to organize the output as a 2D array (signal\_size x nb\_scales), otherwise, the output is a 3D array (signal\_size x noctave x nvoice). |
plot |
if set to T, display the modulus of the continuous wavelet transform on the graphic device. |
The output contains the (complex) values of the wavelet transform of the input signal. The format of the output can be
2D array (signal\_size x nb\_scales)
3D array (signal\_size x noctave x nvoice)
Since Morlet's wavelet is not strictly speaking a wavelet (it is not of vanishing integral), artifacts may occur for certain signals.
continuous (complex) wavelet transform
See discussions in the text of “Practical Time-Frequency Analysis”.
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12)
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12)
Converts the output (modulus or argument) of cwtpolar to a 2D array and displays on the graphic device.
cwtimage(input)
cwtimage(input)
input |
3D array containing a continuous wavelet transform |
The output contains the (complex) values of the wavelet transform of the input signal. The format of the output can be
2D array (signal\_size x nb\_scales)
3D array (signal\_size x noctave x nvoice)
2D array continuous (complex) wavelet transform
See discussions in the text of “Practical Time-Frequency Analysis”.
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12, twoD=FALSE, plot=FALSE) retPolar <- cwtpolar(retChirp) retImageMod <- cwtimage(retPolar$modulus) retImageArg <- cwtimage(retPolar$argument)
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12, twoD=FALSE, plot=FALSE) retPolar <- cwtpolar(retChirp) retImageMod <- cwtimage(retPolar$modulus) retImageArg <- cwtimage(retPolar$argument)
Computes the continuous wavelet transform with (complex-valued) Morlet wavelet and its phase derivative.
cwtp(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
cwtp(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
input |
input signal (possibly complex-valued) |
noctave |
number of powers of 2 for the scale variable |
nvoice |
number of scales in each octave (i.e., between two consecutive powers of 2). |
w0 |
central frequency of the wavelet. |
twoD |
logical variable set to |
plot |
if set to |
list containing the continuous (complex) wavelet transform and the phase derivative
wt |
array of complex numbers for the values of the continuous wavelet transform. |
f |
array of the same dimensions containing the values of the derivative of the phase of the continuous wavelet transform. |
See discussions in the text of “Practical Time-Frequency Analysis”.
cgt
, cwt
, cwtTh
,
DOG
for wavelet transform, and gabor
for
continuous Gabor transform.
## discards imaginary part with error, ## c code does not account for Im(input) x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) chirp <- chirp + 1i * sin(2*pi * (x + 0.004 * (x-256)^2 ) / 16) retChirp <- cwtp(chirp, noctave=5, nvoice=12)
## discards imaginary part with error, ## c code does not account for Im(input) x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) chirp <- chirp + 1i * sin(2*pi * (x + 0.004 * (x-256)^2 ) / 16) retChirp <- cwtp(chirp, noctave=5, nvoice=12)
Converts one of the possible outputs of the function cwt
to modulus and phase.
cwtpolar(cwt, threshold=0)
cwtpolar(cwt, threshold=0)
cwt |
3D array containing the values of a continuous wavelet
transform in the format (signal size |
threshold |
value of a level for the absolute value of the modulus below which
the value of the argument of the output is set to |
The output contains the (complex) values of the wavelet transform of the input signal. The format of the output can be
2D array (signal size nb\_scales)
3D array (signal size noctave
nvoice)
Modulus and Argument of the values of the continuous wavelet transform
output1 |
3D array giving the values (in the same format as the input) of the modulus of the input. |
output2 |
3D array giving the values of the argument of the input. |
See discussions in the text of “Practical Time-Frequency Analysis”.
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12, twoD=FALSE, plot=FALSE) retPolar <- cwtpolar(retChirp)
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwt(chirp, noctave=5, nvoice=12, twoD=FALSE, plot=FALSE) retPolar <- cwtpolar(retChirp)
Computes the synchrosqueezed continuous wavelet transform with the (complex-valued) Morlet wavelet.
cwtsquiz(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
cwtsquiz(input, noctave, nvoice=1, w0=2 * pi, twoD=TRUE, plot=TRUE)
input |
input signal (possibly complex-valued) |
noctave |
number of powers of 2 for the scale variable |
nvoice |
number of scales in each octave (i.e. between two consecutive powers of 2). |
w0 |
central frequency of the wavelet. |
twoD |
logical variable set to T to organize the output as a 2D array
(signal size |
plot |
logical variable set to T to T to display the modulus of the squeezed wavelet transform on the graphic device. |
The output contains the (complex) values of the squeezed wavelet transform of the input signal. The format of the output can be
2D array (signal size nb scales),
3D array (signal size noctave
nvoice).
synchrosqueezed continuous (complex) wavelet transform
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute the continuous wavelet transform with (complex-valued) Cauchy's wavelet.
cwtTh(input, noctave, nvoice=1, moments, twoD=TRUE, plot=TRUE)
cwtTh(input, noctave, nvoice=1, moments, twoD=TRUE, plot=TRUE)
input |
input signal (possibly complex-valued). |
noctave |
number of powers of 2 for the scale variable. |
nvoice |
number of scales in each octave (i.e. between two consecutive powers of 2). |
moments |
number of vanishing moments. |
twoD |
logical variable set to |
plot |
if set to |
The output contains the (complex) values of the wavelet transform of the input signal. The format of the output can be
2D array (signal size nb scales)
3D array (signal size noctave
nvoice)
tmp |
continuous (complex) wavelet transform. |
See discussions in the text of “Practical Time-Frequency Analysis”.
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwtTh(chirp, noctave=5, nvoice=12, moments=20)
x <- 1:512 chirp <- sin(2*pi * (x + 0.002 * (x-256)^2 ) / 16) retChirp <- cwtTh(chirp, noctave=5, nvoice=12, moments=20)
Transient signal.
data(D0)
data(D0)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Transient signal.
data(D4)
data(D4)
A vector containing 1024 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Computes the continuous wavelet transform with for (complex-valued) derivative of Gaussian wavelets.
DOG(input, noctave, nvoice=1, moments, twoD=TRUE, plot=TRUE)
DOG(input, noctave, nvoice=1, moments, twoD=TRUE, plot=TRUE)
input |
input signal (possibly complex-valued). |
noctave |
number of powers of 2 for the scale variable. |
moments |
number of vanishing moments of the wavelet (order of the derivative). |
nvoice |
number of scales in each octave (i.e. between two consecutive powers of 2) |
twoD |
logical variable set to T to organize the output as a 2D array (signal\_size x nb\_scales), otherwise, the output is a 3D array (signal\_size x noctave x nvoice) |
plot |
if set to T, display the modulus of the continuous wavelet transform on the graphic device |
The output contains the (complex) values of the wavelet transform of the input signal. The format of the output can be
2D array (signal\_size x nb\_scales)
3D array (signal\_size x noctave x nvoice)
continuous (complex) wavelet transform
See discussions in the text of “Practical Time-Frequency Analysis”.
Invert the dyadic wavelet transform.
dwinverse(wt, filtername="Gaussian1")
dwinverse(wt, filtername="Gaussian1")
wt |
dyadic wavelet transform |
filtername |
filters used. ("Gaussian1" stands for the filters corresponds to those of Mallat and Zhong's wavlet. And "Haar" stands for the filters of Haar basis. |
Reconstructed signal
See discussions in the text of “Practical Time-Frequency Analysis”.
Successive beat-to-beat intervals for a normal patient.
data(Ekg)
data(Ekg)
A vector containing 16,042 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Compute the local extrema of the dyadic wavelet transform modulus.
ext(wt, scale=FALSE, plot=TRUE)
ext(wt, scale=FALSE, plot=TRUE)
wt |
dyadic wavelet transform. |
scale |
flag indicating if the extrema at each resolution will be plotted at the same scale. |
plot |
if set to TRUE, displays the transform on the graphics device. |
Structure containing:
original |
original signal. |
extrema |
extrema representation. |
Sf |
coarse resolution of signal. |
maxresoln |
number of decomposition scales. |
np |
size of signal. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Computes the cost from the sample of points on the estimated ridge and the matrix used in the reconstruction of the original signal, using simple trapezoidal rule for integrals.
fastgkernel(node, phinode, freqstep, scale, x.inc=1, x.min=node[1], x.max=node[length(node)], plot=FALSE)
fastgkernel(node, phinode, freqstep, scale, x.inc=1, x.min=node[1], x.max=node[length(node)], plot=FALSE)
node |
values of the variable b for the nodes of the ridge |
phinode |
values of the frequency variable |
freqstep |
sampling rate for the frequency axis |
scale |
size of the window |
x.inc |
step unit for the computation of the kernel. |
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
plot |
if set to TRUE, displays the modulus of the matrix of |
Uses trapezoidal rule (instead of Romberg's method) to evaluate the kernel.
matrix of the kernel.
See discussions in the text of “Time-Frequency Analysis”.
gkernel
, fastkernel
, rkernel
,
zerokernel
.
Computes the cost from the sample of points on the estimated ridge and the matrix used in the reconstruction of the original signal, using simple trapezoidal rule for integrals.
fastkernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
fastkernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
node |
values of the variable b for the nodes of the ridge. |
phinode |
values of the scale variable a for the nodes of the ridge. |
nvoice |
number of scales within 1 octave. |
x.inc |
step unit for the computation of the kernel |
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
w0 |
central frequency of the wavelet |
plot |
if set to TRUE, displays the modulus of the matrix of |
Uses trapezoidal rule (instead of Romberg's method) to evaluate the kernel.
matrix of the kernel.
See discussions in the text of “Practical Time-Frequency Analysis”.
kernel
, rkernel
, gkernel
,
zerokernel
.
# The function is currently defined as function(node, phinode, nvoice, x.inc = 1, x.min = node[1], x.max = node[length(node)], w0 = 2 * pi, plot = F) { ######################################################################### # fastkernel: # ----------- # Same as kernel, except that the kernel is computed # using Riemann sums instead of Romberg integration. # # Input: # ------ # node: values of the variable b for the nodes of the ridge # phinode: values of the scale variable a for the nodes of the ridge # nvoice: number of scales within 1 octave # x.inc: step unit for the computation of the kernel # x.min: minimal value of x for the computation of Q2 # x.max: maximal value of x for the computation of Q2 # w0: central frequency of the wavelet # plot: if set to TRUE, displays the modulus of the matrix of Q2 # # Output: # ------- # ker: matrix of the Q2 kernel # ######################################################################### lng <- as.integer((x.max - x.min)/x.inc) + 1 nbnode <- length(node) b.start <- x.min - 50 b.end <- x.max + 50 ker.r <- matrix(0, lng, lng) ker.i <- matrix(0, lng, lng) dim(ker.r) <- c(lng * lng, 1) dim(ker.i) <- c(lng * lng, 1) phinode <- 2 * 2^(phinode/nvoice) z <- .C(fastkernel, ker.r = as.double(ker.r), ker.i = as.double(ker.i), as.integer(x.min), as.integer(x.max), as.integer(x.inc), as.integer(lng), as.double(node), as.double(phinode), as.integer(nbnode), as.double(w0), as.double(b.start), as.double(b.end)) ker.r <- z$ker.r ker.i <- z$ker.i dim(ker.r) <- c(lng, lng) dim(ker.i) <- c(lng, lng) ker <- matrix(0, lng, lng) i <- sqrt(as.complex(-1)) ker <- ker.r + i * ker.i if(plot == T) { par(mfrow = c(1, 1)) image(Mod(ker)) title("Matrix of the reconstructing kernel (modulus)") } ker }
# The function is currently defined as function(node, phinode, nvoice, x.inc = 1, x.min = node[1], x.max = node[length(node)], w0 = 2 * pi, plot = F) { ######################################################################### # fastkernel: # ----------- # Same as kernel, except that the kernel is computed # using Riemann sums instead of Romberg integration. # # Input: # ------ # node: values of the variable b for the nodes of the ridge # phinode: values of the scale variable a for the nodes of the ridge # nvoice: number of scales within 1 octave # x.inc: step unit for the computation of the kernel # x.min: minimal value of x for the computation of Q2 # x.max: maximal value of x for the computation of Q2 # w0: central frequency of the wavelet # plot: if set to TRUE, displays the modulus of the matrix of Q2 # # Output: # ------- # ker: matrix of the Q2 kernel # ######################################################################### lng <- as.integer((x.max - x.min)/x.inc) + 1 nbnode <- length(node) b.start <- x.min - 50 b.end <- x.max + 50 ker.r <- matrix(0, lng, lng) ker.i <- matrix(0, lng, lng) dim(ker.r) <- c(lng * lng, 1) dim(ker.i) <- c(lng * lng, 1) phinode <- 2 * 2^(phinode/nvoice) z <- .C(fastkernel, ker.r = as.double(ker.r), ker.i = as.double(ker.i), as.integer(x.min), as.integer(x.max), as.integer(x.inc), as.integer(lng), as.double(node), as.double(phinode), as.integer(nbnode), as.double(w0), as.double(b.start), as.double(b.end)) ker.r <- z$ker.r ker.i <- z$ker.i dim(ker.r) <- c(lng, lng) dim(ker.i) <- c(lng, lng) ker <- matrix(0, lng, lng) i <- sqrt(as.complex(-1)) ker <- ker.r + i * ker.i if(plot == T) { par(mfrow = c(1, 1)) image(Mod(ker)) title("Matrix of the reconstructing kernel (modulus)") } ker }
Generates a Gabor for given location and frequency.
gabor(sigsize, location, frequency, scale)
gabor(sigsize, location, frequency, scale)
sigsize |
length of the Gabor function. |
location |
position of the Gabor function. |
frequency |
frequency of the Gabor function. |
scale |
size parameter for the Gabor function. |
complex 1D array of size sigsize.
See discussions in the text of “Practical Time-Frequency Analysis”.
Reconstructs a real-valued signal from ridges found by crazy climbers on a Gabor transform.
gcrcrec(siginput, inputgt, beemap, nvoice, freqstep, scale, compr, bstep=5, ptile=0.01, epsilon=0, fast=TRUE, para=5, minnbnodes=3, hflag=FALSE, real=FALSE, plot=2)
gcrcrec(siginput, inputgt, beemap, nvoice, freqstep, scale, compr, bstep=5, ptile=0.01, epsilon=0, fast=TRUE, para=5, minnbnodes=3, hflag=FALSE, real=FALSE, plot=2)
siginput |
original signal. |
inputgt |
Gabor transform. |
beemap |
occupation measure, output of |
nvoice |
number of frequencies. |
freqstep |
sampling step for frequency axis. |
scale |
size of windows. |
compr |
compression rate to be applied to the ridges. |
bstep |
size (in the time direction) of the steps for chaining. |
ptile |
threshold of ridge |
epsilon |
constant in front of the smoothness term in penalty function. |
fast |
if set to TRUE, uses trapezoidal rule to evaluate |
para |
scale parameter for extrapolating the ridges. |
minnbnodes |
minimal number of points per ridge. |
hflag |
if set to FALSE, uses the identity as first term
in the kernel. If not, uses |
real |
if set to |
plot |
|
When ptile
is high, boundary effects may appear. para
controls extrapolation of the ridge.
Returns a structure containing the following elements:
rec |
reconstructed signal. |
ordered |
image of the ridges (with different colors). |
comp |
2D array containing the signals reconstructed from ridges. |
See discussions in the text of “Practical Time-Frequency Analysis”.
crc
, cfamily
, crcrec
,
scrcrec
.
Computes the cost from the sample of points on the estimated ridge and the matrix used in the reconstruction of the original signal.
gkernel(node, phinode, freqstep, scale, x.inc=1, x.min=node[1], x.max=node[length(node)], plot=FALSE)
gkernel(node, phinode, freqstep, scale, x.inc=1, x.min=node[1], x.max=node[length(node)], plot=FALSE)
node |
values of the variable b for the nodes of the ridge. |
phinode |
values of the scale variable a for the nodes of the ridge. |
freqstep |
sampling rate for the frequency axis. |
scale |
size of the window. |
x.inc |
step unit for the computation of the kernel. |
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
plot |
if set to TRUE, displays the modulus of the matrix of |
matrix of the kernel
See discussions in the text of “Time-Frequency Analysis”.
fastgkernel
, kernel
, rkernel
,
fastkernel
, zerokernel
.
Reconstructs signal from a “regularly sampled” ridge, in the Gabor case.
gregrec(siginput, gtinput, phi, nbnodes, nvoice, freqstep, scale, epsilon=0, fast=FALSE, plot=FALSE, para=0, hflag=FALSE, real=FALSE, check=FALSE)
gregrec(siginput, gtinput, phi, nbnodes, nvoice, freqstep, scale, epsilon=0, fast=FALSE, plot=FALSE, para=0, hflag=FALSE, real=FALSE, check=FALSE)
siginput |
input signal. |
gtinput |
Gabor transform, output of |
phi |
unsampled ridge. |
nbnodes |
number of nodes used for the reconstruction. |
nvoice |
number of different scales per octave |
freqstep |
sampling rate for the frequency axis |
scale |
size parameter for the Gabor function. |
epsilon |
coefficient of the |
fast |
if set to T, the kernel is computed using trapezoidal rule. |
plot |
if set to TRUE, displays original and reconstructed signals |
para |
scale parameter for extrapolating the ridges. |
hflag |
if set to TRUE, uses |
real |
if set to TRUE, uses only real constraints on the transform. |
check |
if set to TRUE, computes |
Returns a list containing:
sol |
reconstruction from a ridge. |
A |
<gaborlets,dualgaborlets> matrix. |
lam |
coefficients of dual wavelets in reconstructed signal. |
dualwave |
array containing the dual wavelets. |
gaborets |
array containing the wavelets on sampled ridge. |
solskel |
Gabor transform of sol, restricted to the ridge. |
inputskel |
Gabor transform of signal, restricted to the ridge. |
Q2 |
second part of the reconstruction kernel. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Reconstructs signal from sample of a ridge, in the Gabor case.
gridrec(gtinput, node, phinode, nvoice, freqstep, scale, Qinv, epsilon, np, real=FALSE, check=FALSE)
gridrec(gtinput, node, phinode, nvoice, freqstep, scale, Qinv, epsilon, np, real=FALSE, check=FALSE)
gtinput |
Gabor transform, output of |
node |
time coordinates of the ridge samples. |
phinode |
frequency coordinates of the ridge samples. |
nvoice |
number of different frequencies. |
freqstep |
sampling rate for the frequency axis. |
scale |
scale of the window. |
Qinv |
inverse of the matrix |
epsilon |
coefficient of the |
np |
number of samples of the reconstructed signal. |
real |
if set to TRUE, uses only constraints on the real part of the transform. |
check |
if set to TRUE, computes |
Returns a list containing the reconstructed signal and the chained ridges.
sol |
reconstruction from a ridge. |
A |
<gaborlets,dualgaborlets> matrix. |
lam |
coefficients of dual gaborlets in reconstructed signal. |
dualwave |
array containing the dual gaborlets. |
gaborets |
array of gaborlets located on the ridge samples. |
solskel |
Gabor transform of sol, restricted to the ridge. |
inputskel |
Gabor transform of signal, restricted to the ridge. |
See discussions in the text of “Practical Time-Frequency Analysis”.
sridrec
, gregrec
, regrec
,
regrec2
.
Generate a sampled identity matrix.
gsampleOne(node, scale, np)
gsampleOne(node, scale, np)
node |
location of the reconstruction gabor functions. |
scale |
scale of the gabor functions. |
np |
size of the reconstructed signal. |
diagonal of the “sampled” term (1D vector)
See discussions in the text of “Time-Frequency Analysis”.
Generation of Gabor functions located on the ridge.
gwave(bridge, omegaridge, nvoice, freqstep, scale, np, N)
gwave(bridge, omegaridge, nvoice, freqstep, scale, np, N)
bridge |
time coordinates of the ridge samples |
omegaridge |
frequency coordinates of the ridge samples |
nvoice |
number of different scales per octave |
freqstep |
sampling rate for the frequency axis |
scale |
scale of the window |
np |
size of the reconstruction kernel |
N |
number of complex constraints |
Array of Gabor functions located on the ridge samples
See discussions in the text of "Time-Frequency Analysis”.
Generation of the real parts of gabor functions located on a ridge.
(modification of gwave
.)
gwave2(bridge, omegaridge, nvoice, freqstep, scale, np, N)
gwave2(bridge, omegaridge, nvoice, freqstep, scale, np, N)
bridge |
time coordinates of the ridge samples |
omegaridge |
frequency coordinates of the ridge samples |
nvoice |
number of different scales per octave |
freqstep |
sampling rate for the frequency axis |
scale |
scale of the window |
np |
size of the reconstruction kernel |
N |
number of complex constraints |
Array of real Gabor functions located on the ridge samples
See discussions in the text of “Time-Frequency Analysis”.
Example of speech signal.
data(HOWAREYOU)
data(HOWAREYOU)
A vector containing 5151 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Estimates Hurst exponent from a wavelet transform.
hurst.est(wspec, range, nvoice, plot=TRUE)
hurst.est(wspec, range, nvoice, plot=TRUE)
wspec |
wavelet spectrum (output of |
range |
range of scales from which estimate the exponent. |
nvoice |
number of scales per octave of the wavelet transform. |
plot |
if set to |
complex 1D array of size sigsize.
See discussions in the text of “Practical Time-Frequency Analysis”.
# White Noise Hurst Exponent: The plots on the top row of Figure 6.8 # were produced by the folling S-commands. These make use of the two # functions Hurst.est (estimation of Hurst exponent from CWT) and # wspec.pl (display wavelet spectrum). # Compare the periodogram and the wavelet spectral estimate. wnoise <- rnorm(8192) plot.ts(wnoise) spwnoise <- fft(wnoise) spwnoise <- Mod(spwnoise) spwnoise <- spwnoise*spwnoise plot(spwnoise[1:4096], log="xy", type="l") lswnoise <- lsfit(log10(1:4096), log10(spwnoise[1:4096])) abline(lswnoise$coef) cwtwnoise <- DOG(wnoise, 10, 5, 1, plot=FALSE) mcwtwnoise <- Mod(cwtwnoise) mcwtwnoise <- mcwtwnoise*mcwtwnoise wspwnoise <- tfmean(mcwtwnoise, plot=FALSE) wspec.pl(wspwnoise, 5) hurst.est(wspwnoise, 1:50, 5)
# White Noise Hurst Exponent: The plots on the top row of Figure 6.8 # were produced by the folling S-commands. These make use of the two # functions Hurst.est (estimation of Hurst exponent from CWT) and # wspec.pl (display wavelet spectrum). # Compare the periodogram and the wavelet spectral estimate. wnoise <- rnorm(8192) plot.ts(wnoise) spwnoise <- fft(wnoise) spwnoise <- Mod(spwnoise) spwnoise <- spwnoise*spwnoise plot(spwnoise[1:4096], log="xy", type="l") lswnoise <- lsfit(log10(1:4096), log10(spwnoise[1:4096])) abline(lswnoise$coef) cwtwnoise <- DOG(wnoise, 10, 5, 1, plot=FALSE) mcwtwnoise <- Mod(cwtwnoise) mcwtwnoise <- mcwtwnoise*mcwtwnoise wspwnoise <- tfmean(mcwtwnoise, plot=FALSE) wspec.pl(wspwnoise, 5) hurst.est(wspwnoise, 1:50, 5)
Estimate a (single) ridge from a time-frequency representation, using the ICM minimization method.
icm(modulus, guess, tfspec=numeric(dim(modulus)[2]), subrate=1, mu=1, lambda=2 * mu, iteration=100)
icm(modulus, guess, tfspec=numeric(dim(modulus)[2]), subrate=1, mu=1, lambda=2 * mu, iteration=100)
modulus |
Time-Frequency representation (real valued). |
guess |
Initial guess for the algorithm. |
tfspec |
Estimate for the contribution of the noise to modulus. |
subrate |
Subsampling rate for ridge estimation. |
mu |
Coefficient of the ridge's second derivative in cost function. |
lambda |
Coefficient of the ridge's derivative in cost function. |
iteration |
Maximal number of moves. |
To accelerate convergence, it is useful to preprocess modulus before
running annealing method. Such a preprocessing (smoothing and
subsampling of modulus) is implemented in icm
. The
parameter subrate specifies the subsampling rate.
Returns the estimated ridge and the cost function.
ridge |
1D array (of same length as the signal) containing the ridge. |
cost |
1D array containing the cost function. |
See discussions in the text of “Practical Time-Frequency Analysis”.
corona
, coronoid
, and snake
,
snakoid
.
Computes the cost from the sample of points on the estimated ridge and the matrix used in the reconstruction of the original signal
kernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
kernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
node |
values of the variable b for the nodes of the ridge. |
phinode |
values of the scale variable a for the nodes of the ridge. |
nvoice |
number of scales within 1 octave. |
x.inc |
step unit for the computation of the kernel. |
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
w0 |
central frequency of the wavelet. |
plot |
if set to TRUE, displays the modulus of the matrix of |
The kernel is evaluated using Romberg's method.
matrix of the kernel
See discussions in the text of "Time-Frequency Analysis".
Trimming of dyadic wavelet transform local extrema, using bootstrapping.
mbtrim(extrema, scale=FALSE, prct=0.95)
mbtrim(extrema, scale=FALSE, prct=0.95)
extrema |
dyadic wavelet transform extrema (output of |
scale |
when set, the wavelet transform at each scale will be plotted with the same scale. |
prct |
percentage critical value used for thresholding |
The distribution of extrema of dyadic wavelet transform at each scale is generated by bootstrap method, and the 95% critical value is used for thresholding the extrema of the signal.
Structure containing
original |
original signal. |
extrema |
trimmed extrema representation. |
Sf |
coarse resolution of signal. |
maxresoln |
number of decomposition scales. |
np |
size of signal. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Trimming of dyadic wavelet transform local extrema, assuming normal distribution.
mntrim(extrema, scale=FALSE, prct=0.95)
mntrim(extrema, scale=FALSE, prct=0.95)
extrema |
dyadic wavelet transform extrema (output of |
scale |
when set, the wavelet transform at each scale will be plotted with the same scale. |
prct |
percentage critical value used for thresholding |
The distribution of extrema of dyadic wavelet transform at each scale is generated by simulation, assuming a normal distribution, and the 95% critical value is used for thresholding the extrema of the signal.
Structure containing
original |
original signal. |
extrema |
trimmed extrema representation. |
Sf |
coarse resolution of signal. |
maxresoln |
number of decomposition scales. |
np |
size of signal. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Computes a Morlet wavelet at the point of the time-scale plane given in the input
morlet(sigsize, location, scale, w0=2 * pi)
morlet(sigsize, location, scale, w0=2 * pi)
sigsize |
length of the output. |
location |
time location of the wavelet. |
scale |
scale of the wavelet. |
w0 |
central frequency of the wavelet. |
The details of this construction (including the definition formulas) are given in the text.
Returns the values of the complex Morlet wavelet at the point of the time-scale plane given in the input
See discussions in the text of “Practical Time-Frequency Analysis”.
Generates the Morlet wavelets at the sample points of the ridge.
morwave(bridge, aridge, nvoice, np, N, w0=2 * pi)
morwave(bridge, aridge, nvoice, np, N, w0=2 * pi)
bridge |
time coordinates of the ridge samples. |
aridge |
scale coordinates of the ridge samples. |
nvoice |
number of different scales per octave. |
np |
number of samples in the input signal. |
N |
size of reconstructed signal. |
w0 |
central frequency of the wavelet. |
Returns the Morlet wavelets at the samples of the time-scale plane given in the input: complex array of Morlet wavelets located on the ridge samples
See discussions in the text of “Time-Frequency Analysis”.
Generates the real parts of the Morlet wavelets at the sample points of a ridge
morwave2(bridge, aridge, nvoice, np, N, w0=2 * pi)
morwave2(bridge, aridge, nvoice, np, N, w0=2 * pi)
bridge |
time coordinates of the ridge samples. |
aridge |
scale coordinates of the ridge samples. |
nvoice |
number of different scales per octave. |
np |
number of samples in the input signal. |
N |
size of reconstructed signal. |
w0 |
central frequency of the wavelet. |
Returns the real parts of the Morlet wavelets at the samples of the time-scale plane given in the input: array of Morlet wavelets located on the ridge samples
See discussions in the text of “Time-Frequency Analysis”.
Reconstruct from dyadic wavelet transform modulus extrema. The reconstructed signal preserves locations and values at extrema.
mrecons(extrema, filtername="Gaussian1", readflag=FALSE)
mrecons(extrema, filtername="Gaussian1", readflag=FALSE)
extrema |
the extrema representation. |
filtername |
filter used for dyadic wavelet transform. |
readflag |
if set to T, read reconstruction kernel from precomputed file. This is not supported in the current package, and will cause an error. |
The reconstruction involves only the wavelet coefficients, without taking care of the coarse scale component. The latter may be added a posteriori.
Structure containing
f |
the reconstructed signal. |
g |
reconstructed signal plus mean of original signal. |
h |
reconstructed signal plus coarse scale component of original signal. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Dyadic wavelet transform, with Mallat's wavelet. The reconstructed signal preserves locations and values at extrema.
mw(inputdata, maxresoln, filtername="Gaussian1", scale=FALSE, plot=TRUE)
mw(inputdata, maxresoln, filtername="Gaussian1", scale=FALSE, plot=TRUE)
inputdata |
either a text file or an R object containing data. |
maxresoln |
number of decomposition scales. |
filtername |
name of filter (either Gaussian1 for Mallat and Zhong's wavelet or Haar wavelet). |
scale |
when set, the wavelet transform at each scale is plotted with the same scale. |
plot |
indicate if the wavelet transform at each scale will be plotted. |
The decomposition goes from resolution 1 to the given maximum resolution.
Structure containing
original |
original signal. |
Wf |
dyadic wavelet transform of signal. |
Sf |
multiresolution of signal. |
maxresoln |
number of decomposition scales. |
np |
size of signal. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Noisy gravitational wave.
data(noisywave)
data(noisywave)
A vector containing 8192 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Splits the graphics device into prescrivbed number of windows.
npl(nbrow)
npl(nbrow)
nbrow |
number of plots. |
Plot extrema of dyadic wavelet transform.
plotResult(result, original, maxresoln, scale=FALSE, yaxtype="s")
plotResult(result, original, maxresoln, scale=FALSE, yaxtype="s")
result |
result. |
original |
input signal. |
maxresoln |
number of decomposition scales. |
scale |
when set, the extrema at each scale is plotted withe the same scale. |
yaxtype |
y axis type (see R manual). |
See discussions in the text of “Time-Frequency Analysis”.
Plot dyadic wavelet transform.
plotwt(original, psi, phi, maxresoln, scale=FALSE, yaxtype="s")
plotwt(original, psi, phi, maxresoln, scale=FALSE, yaxtype="s")
original |
input signal. |
psi |
dyadic wavelet transform. |
phi |
scaling function transform at last resolution. |
maxresoln |
number of decomposition scales. |
scale |
when set, the wavelet transform at each scale is plotted with the same scale. |
yaxtype |
axis type (see R manual). |
See discussions in the text of “Time-Frequency Analysis”.
plotResult
, epl
, wpl
.
Pure gravitational wave.
data(purwave)
data(purwave)
A vector containing 8192 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Reconstructs signal from a “regularly sampled” ridge, in the wavelet case.
regrec(siginput, cwtinput, phi, compr, noct, nvoice, epsilon=0, w0=2 * pi, fast=FALSE, plot=FALSE, para=0, hflag=FALSE, check=FALSE, minnbnodes=2, real=FALSE)
regrec(siginput, cwtinput, phi, compr, noct, nvoice, epsilon=0, w0=2 * pi, fast=FALSE, plot=FALSE, para=0, hflag=FALSE, check=FALSE, minnbnodes=2, real=FALSE)
siginput |
input signal. |
cwtinput |
wavelet transform, output of |
phi |
unsampled ridge. |
compr |
subsampling rate for the wavelet coefficients (at scale 1) |
noct |
number of octaves (powers of 2) |
nvoice |
number of different scales per octave |
epsilon |
coefficient of the |
w0 |
central frequency of Morlet wavelet |
fast |
if set to TRUE, the kernel is computed using trapezoidal rule. |
plot |
if set to TRUE, displays original and reconstructed signals |
para |
scale parameter for extrapolating the ridges. |
hflag |
if set to TRUE, uses |
check |
if set to TRUE, computes |
minnbnodes |
minimum number of nodes for the reconstruction. |
real |
if set to TRUE, uses only real constraints on the transform. |
Returns a list containing:
sol |
reconstruction from a ridge. |
A |
<wavelets,dualwavelets> matrix. |
lam |
coefficients of dual wavelets in reconstructed signal. |
dualwave |
array containing the dual wavelets. |
morvelets |
array containing the wavelets on sampled ridge. |
solskel |
wavelet transform of sol, restricted to the ridge. |
inputskel |
wavelet transform of signal, restricted to the ridge. |
Q2 |
second part of the reconstruction kernel. |
nbnodes |
number of nodes used for the reconstruction. |
See discussions in the text of “Practical Time-Frequency Analysis”.
regrec2
, ridrec
, gregrec
,
gridrec
.
Reconstructs signal from a “regularly sampled” ridge, in the wavelet case, from a precomputed kernel.
regrec2(siginput, cwtinput, phi, nbnodes, noct, nvoice, Q2, epsilon=0.5, w0=2 * pi, plot=FALSE)
regrec2(siginput, cwtinput, phi, nbnodes, noct, nvoice, Q2, epsilon=0.5, w0=2 * pi, plot=FALSE)
siginput |
input signal. |
cwtinput |
wavelet transform, output of |
phi |
unsampled ridge. |
nbnodes |
number of samples on the ridge |
noct |
number of octaves (powers of 2) |
nvoice |
number of different scales per octave |
Q2 |
second term of the reconstruction kernel |
epsilon |
coefficient of the |
w0 |
central frequency of Morlet wavelet |
plot |
if set to TRUE, displays original and reconstructed signals |
The computation of the kernel may be time consuming. This function avoids recomputing it if it was computed already.
Returns a list containing:
sol |
reconstruction from a ridge. |
A |
<wavelets,dualwavelets> matrix. |
lam |
coefficients of dual wavelets in reconstructed signal. |
dualwave |
array containing the dual wavelets. |
morvelets |
array containing the wavelets on sampled ridge. |
solskel |
wavelet transform of sol, restricted to the ridge. |
inputskel |
wavelet transform of signal, restricted to the ridge. |
nbnodes |
number of nodes used for the reconstruction. |
See discussions in the text of “Practical Time-Frequency Analysis”.
regrec
, gregrec
, ridrec
,
sridrec
.
Given a ridge phi (for the Gabor transform), returns a (regularly) subsampled version of length nbnodes.
RidgeSampling(phi, nbnodes)
RidgeSampling(phi, nbnodes)
phi |
ridge (1D array). |
nbnodes |
number of samples. |
Gabor ridges are sampled uniformly.
Returns a list containing the discrete values of the ridge.
node |
time coordinates of the ridge samples. |
phinode |
frequency coordinates of the ridge samples. |
See discussions in the text of "Time-Frequency Analysis”.
Reconstructs signal from sample of a ridge, in the wavelet case.
ridrec(cwtinput, node, phinode, noct, nvoice, Qinv, epsilon, np, w0=2 * pi, check=FALSE, real=FALSE)
ridrec(cwtinput, node, phinode, noct, nvoice, Qinv, epsilon, np, w0=2 * pi, check=FALSE, real=FALSE)
cwtinput |
wavelet transform, output of |
node |
time coordinates of the ridge samples. |
phinode |
scale coordinates of the ridge samples. |
noct |
number of octaves (powers of 2). |
nvoice |
number of different scales per octave. |
Qinv |
inverse of the matrix |
epsilon |
coefficient of the |
np |
number of samples of the reconstructed signal. |
w0 |
central frequency of Morlet wavelet. |
check |
if set to TRUE, computes |
real |
if set to TRUE, uses only constraints on the real part of the transform. |
Returns a list containing the reconstructed signal and the chained ridges.
sol |
reconstruction from a ridge |
A |
<wavelets,dualwavelets> matrix |
lam |
coefficients of dual wavelets in reconstructed signal. |
dualwave |
array containing the dual wavelets. |
morvelets |
array of morlet wavelets located on the ridge samples. |
solskel |
wavelet transform of sol, restricted to the ridge |
inputskel |
wavelet transform of signal, restricted to the ridge |
See discussions in the text of “Practical Time-Frequency Analysis”.
Computes the cost from the sample of points on the estimated ridge
and the matrix used in the reconstruction of the original signal,
in the case of real constraints. Modification of the function
kernel
.
rkernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
rkernel(node, phinode, nvoice, x.inc=1, x.min=node[1], x.max=node[length(node)], w0=2 * pi, plot=FALSE)
node |
values of the variable b for the nodes of the ridge. |
phinode |
values of the scale variable a for the nodes of the ridge. |
nvoice |
number of scales within 1 octave. |
x.inc |
step unit for the computation of the kernel. |
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
w0 |
central frequency of the wavelet. |
plot |
if set to TRUE, displays the modulus of the matrix of |
Uses Romberg's method for computing the kernel.
matrix of the kernel
See discussions in the text of "Time-Frequency Analysis".
kernel
, fastkernel
, gkernel
,
zerokernel
.
Reconstructs signal from ridges obtained by crc
,
using the restriction of the transform to the ridge.
scrcrec(siginput, tfinput, beemap, bstep=5, ptile=0.01, plot=2)
scrcrec(siginput, tfinput, beemap, bstep=5, ptile=0.01, plot=2)
siginput |
input signal. |
tfinput |
|
beemap |
output of crazy climber algorithm |
bstep |
used for the chaining (see |
ptile |
threshold on the measure beemap (see |
plot |
1: displays signal,components, and reconstruction one after another. |
Returns a list containing the reconstructed signal and the chained ridges.
rec |
reconstructed signal |
ordered |
image of the ridges (with different colors) |
comp |
2D array containing the signals reconstructed from ridges |
See discussions in the text of “Practical Time-Frequency Analysis”.
crc
,cfamily
for crazy climbers method,
crcrec
for reconstruction methods.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(sig_W_tilda.1)
data(sig_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(sig_W_tilda.1)
data(sig_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(sig_W_tilda.1)
data(sig_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(sig_W_tilda.1)
data(sig_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(sig_W_tilda.1)
data(sig_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
The package maintainer believes this file was read or written by a C function (signal_W_tilda) called from mrecons, and was a precomputed kernel. All C function disk reads and writes have been disabled, but the files are preserved for historical purposes.
data(signal_W_tilda.1)
data(signal_W_tilda.1)
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
mrecons
.
Computes the reconstructed signal from the ridge, given the inverse of the matrix Q.
skeleton(cwtinput, Qinv, morvelets, bridge, aridge, N)
skeleton(cwtinput, Qinv, morvelets, bridge, aridge, N)
cwtinput |
continuous wavelet transform (as the output of cwt) |
Qinv |
inverse of the reconstruction kernel (2D array) |
morvelets |
array of Morlet wavelets located at the ridge samples |
bridge |
time coordinates of the ridge samples |
aridge |
scale coordinates of the ridge samples |
N |
size of reconstructed signal |
Returns a list of the elements of the reconstruction of a signal from sample points of a ridge
sol |
reconstruction from a ridge |
A |
matrix of the inner products |
lam |
coefficients of dual wavelets in reconstructed signal.
They are the Lagrange multipliers |
dualwave |
array containing the dual wavelets. |
See discussions in the text of “Practical Time-Frequency Analysis”.
skeleton2
, zeroskeleton
,
zeroskeleton2
.
Computes the reconstructed signal from the ridge in the case of real constraints.
skeleton2(cwtinput, Qinv, morvelets, bridge, aridge, N)
skeleton2(cwtinput, Qinv, morvelets, bridge, aridge, N)
cwtinput |
continuous wavelet transform (as the output of cwt). |
Qinv |
inverse of the reconstruction kernel (2D array). |
morvelets |
array of Morlet wavelets located at the ridge samples. |
bridge |
time coordinates of the ridge samples. |
aridge |
scale coordinates of the ridge samples. |
N |
size of reconstructed signal. |
Returns a list of the elements of the reconstruction of a signal from sample points of a ridge
sol |
reconstruction from a ridge. |
A |
matrix of the inner products. |
lam |
coefficients of dual wavelets in reconstructed signal. They
are the Lagrange multipliers |
dualwave |
array containing the dual wavelets. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Smooth a time series by averaging window.
smoothts(ts, windowsize)
smoothts(ts, windowsize)
ts |
Time series. |
windowsize |
Length of smoothing window. |
Smoothed time series (1D array).
See discussions in the text of “Time-Frequency Analysis”.
smooth the wavelet (or Gabor) transform in the time direction.
smoothwt(modulus, subrate, flag=FALSE)
smoothwt(modulus, subrate, flag=FALSE)
modulus |
Time-Frequency representation (real valued). |
subrate |
Length of smoothing window. |
flag |
If set to TRUE, subsample the representation. |
2D array containing the smoothed transform.
See discussions in the text of “Time-Frequency Analysis”.
corona
, coronoid
, snake
,
snakoid
.
Estimate a ridge from a time-frequency representation, using the snake method.
snake(tfrep, guessA, guessB, snakesize=length(guessB), tfspec=numeric(dim(modulus)[2]), subrate=1, temprate=3, muA=1, muB=muA, lambdaB=2 * muB, lambdaA=2 * muA, iteration=1000000, seed=-7, costsub=1, stagnant=20000, plot=TRUE)
snake(tfrep, guessA, guessB, snakesize=length(guessB), tfspec=numeric(dim(modulus)[2]), subrate=1, temprate=3, muA=1, muB=muA, lambdaB=2 * muB, lambdaA=2 * muA, iteration=1000000, seed=-7, costsub=1, stagnant=20000, plot=TRUE)
tfrep |
Time-Frequency representation (real valued). |
guessA |
Initial guess for the algorithm (frequency variable). |
guessB |
Initial guess for the algorithm (time variable). |
snakesize |
the length of the initial guess of time variable. |
tfspec |
Estimate for the contribution of the noise to modulus. |
subrate |
Subsampling rate for ridge estimation. |
temprate |
Initial value of temperature parameter. |
muA |
Coefficient of the ridge's derivative in cost function (frequency component). |
muB |
Coefficient of the ridge's derivative in cost function (time component). |
lambdaB |
Coefficient of the ridge's second derivative in cost function (time component). |
lambdaA |
Coefficient of the ridge's second derivative in cost function (frequency component). |
iteration |
Maximal number of moves. |
seed |
Initialization of random number generator. |
costsub |
Subsampling of cost function in output. |
stagnant |
maximum number of steps without move (for the stopping criterion) |
plot |
when set (by default), certain results will be displayed |
Returns a structure containing:
ridge |
1D array (of same length as the signal) containing the ridge. |
cost |
1D array containing the cost function. |
See discussions in the text of “Practical Time-Frequency Analysis”.
corona
, coronoid
, icm
,
snakoid
.
Restrict time-frequency transform to a snake.
snakeview(modulus, snake)
snakeview(modulus, snake)
modulus |
Time-Frequency representation (real valued). |
snake |
Time and frequency components of a snake. |
Recall that a snake is a (two components) R structure.
2D array containing the restriction of the transform modulus to the snake.
See discussions in the text of “Time-Frequency Analysis”.
Estimate a ridge from a time-frequency representation, using the modified snake method (modified cost function).
snakoid(modulus, guessA, guessB, snakesize=length(guessB), tfspec=numeric(dim(modulus)[2]), subrate=1, temprate=3, muA=1, muB=muA, lambdaB=2 * muB, lambdaA=2 * muA, iteration=1000000, seed=-7, costsub=1, stagnant=20000, plot=TRUE)
snakoid(modulus, guessA, guessB, snakesize=length(guessB), tfspec=numeric(dim(modulus)[2]), subrate=1, temprate=3, muA=1, muB=muA, lambdaB=2 * muB, lambdaA=2 * muA, iteration=1000000, seed=-7, costsub=1, stagnant=20000, plot=TRUE)
modulus |
Time-Frequency representation (real valued). |
guessA |
Initial guess for the algorithm (frequency variable). |
guessB |
Initial guess for the algorithm (time variable). |
snakesize |
The length of the first guess of time variable. |
tfspec |
Estimate for the contribution of srthe noise to modulus. |
subrate |
Subsampling rate for ridge estimation. |
temprate |
Initial value of temperature parameter. |
muA |
Coefficient of the ridge's derivative in cost function (frequency component). |
muB |
Coefficient of the ridge's derivative in cost function (time component). |
lambdaB |
Coefficient of the ridge's second derivative in cost function (time component). |
lambdaA |
Coefficient of the ridge's second derivative in cost function (frequency component). |
iteration |
Maximal number of moves. |
seed |
Initialization of random number generator. |
costsub |
Subsampling of cost function in output. |
stagnant |
Maximum number of stationary iterations before stopping. |
plot |
when set(default), some results will be displayed |
Returns a structure containing:
ridge |
1D array (of same length as the signal) containing the ridge. |
cost |
1D array containing the cost function. |
plot |
when set(default), some results will be displayed. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Simple reconstruction of a real valued signal from a ridge, by restriction of the transform to the ridge.
sridrec(tfinput, ridge)
sridrec(tfinput, ridge)
tfinput |
time-frequency representation. |
ridge |
ridge (1D array). |
(real) reconstructed signal (1D array)
See discussions in the text of “Practical Time-Frequency Analysis”.
Computes singular value decomposition of a matrix.
SVD(a)
SVD(a)
a |
input matrix. |
R interface for Numerical Recipes singular value decomposition routine.
a structure containing the 3 matrices of the singular value decomposition of the input.
See discussions in the text of “Time-Frequency Analysis”.
Computes the maxima (for each fixed value of the time variable) of the modulus of a continuous wavelet transform.
tfgmax(input, plot=TRUE)
tfgmax(input, plot=TRUE)
input |
wavelet transform (as the output of the function |
plot |
if set to TRUE, displays the values of the energy as a function of the scale. |
output |
values of the maxima (1D array) |
pos |
positions of the maxima (1D array) |
See discussions in the text of “Practical Time-Frequency Analysis”.
Computes the local maxima (for each fixed value of the time variable) of the modulus of a time-frequency transform.
tflmax(input, plot=TRUE)
tflmax(input, plot=TRUE)
input |
time-frequency transform (real 2D array). |
plot |
if set to T, displays the local maxima on the graphic device. |
values of the maxima (2D array).
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute the mean of time-frequency representation frequency by frequency.
tfmean(input, plot=TRUE)
tfmean(input, plot=TRUE)
input |
|
plot |
if set to T, displays the values of the energy as a function of the scale (or frequency). |
1D array containing the noise estimate.
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute a percentile of time-frequency representation frequency by frequency.
tfpct(input, percent=0.8, plot=TRUE)
tfpct(input, percent=0.8, plot=TRUE)
input |
|
percent |
percentile to be retained. |
plot |
if set to T, displays the values of the energy as a function of the scale (or frequency). |
1D array containing the noise estimate.
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute the variance of time-frequency representation frequency by frequency.
tfvar(input, plot=TRUE)
tfvar(input, plot=TRUE)
input |
|
plot |
if set to T, displays the values of the energy as a function of the scale (or frequency). |
1D array containing the noise estimate.
See discussions in the text of “Practical Time-Frequency Analysis”.
Numerous functions were not documented in the original Swave help files.
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute DOG wavelet transform at one scale.
vDOG(input, scale, moments)
vDOG(input, scale, moments)
input |
Input signal (1D array). |
scale |
Scale at which the wavelet transform is to be computed. |
moments |
number of vanishing moments. |
1D (complex) array containing wavelet transform at one scale.
See discussions in the text of “Practical Time-Frequency Analysis”.
Generate Gabor functions at specified positions on a ridge.
vecgabor(sigsize, nbnodes, location, frequency, scale)
vecgabor(sigsize, nbnodes, location, frequency, scale)
sigsize |
Signal size. |
nbnodes |
Number of wavelets to be generated. |
location |
b coordinates of the ridge samples (1D array of length nbnodes). |
frequency |
frequency coordinates of the ridge samples (1D array of length nbnodes). |
scale |
size parameter for the Gabor functions. |
size parameter for the Gabor functions.
Generate Morlet wavelets at specified positions on a ridge.
vecmorlet(sigsize, nbnodes, bridge, aridge, w0=2 * pi)
vecmorlet(sigsize, nbnodes, bridge, aridge, w0=2 * pi)
sigsize |
Signal size. |
nbnodes |
Number of wavelets to be generated. |
bridge |
b coordinates of the ridge samples (1D array of length nbnodes). |
aridge |
a coordinates of the ridge samples (1D array of length nbnodes). |
w0 |
Center frequency of the wavelet. |
2D (complex) array containing wavelets located at the specific points.
Compute Gabor transform for fixed frequency.
vgt(input, frequency, scale, plot=FALSE)
vgt(input, frequency, scale, plot=FALSE)
input |
Input signal (1D array). |
frequency |
frequency at which the Gabor transform is to be computed. |
scale |
frequency at which the Gabor transform is to be computed. |
plot |
if set to TRUE, plots the real part of cgt on the graphic device. |
1D (complex) array containing Gabor transform at specified frequency.
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute Morlet's wavelet transform at one scale.
vwt(input, scale, w0=2 * pi)
vwt(input, scale, w0=2 * pi)
input |
Input signal (1D array). |
scale |
Scale at which the wavelet transform is to be computed. |
w0 |
Center frequency of the wavelet. |
1D (complex) array containing wavelet transform at one scale.
See discussions in the text of “Practical Time-Frequency Analysis”.
Given a ridge (for the wavelet transform), returns a
(appropriately) subsampled version with a given subsampling rate.
wRidgeSampling(phi, compr, nvoice)
wRidgeSampling(phi, compr, nvoice)
phi |
ridge (1D array). |
compr |
subsampling rate for the ridge. |
nvoice |
number of voices per octave. |
To account for the variable sizes of wavelets, the sampling rate of a wavelet ridge is not uniform, and is proportional to the scale.
Returns a list containing the discrete values of the ridge.
node |
time coordinates of the ridge samples. |
phinode |
scale coordinates of the ridge samples. |
nbnode |
number of nodes of the ridge samples. |
Displays normalized log of wavelet spectrum.
wspec.pl(wspec, nvoice)
wspec.pl(wspec, nvoice)
wspec |
wavelet spectrum. |
nvoice |
number of voices. |
See discussions in the text of “Practical Time-Frequency Analysis”.
Compute the Wigner-Ville transform, without any smoothing.
WV(input, nvoice, freqstep = (1/nvoice), plot = TRUE)
WV(input, nvoice, freqstep = (1/nvoice), plot = TRUE)
input |
input signal (possibly complex-valued) |
nvoice |
number of frequency bands |
freqstep |
sampling rate for the frequency axis |
plot |
if set to TRUE, displays the modulus of CWT on the graphic device. |
(complex) Wigner-Ville transform.
See discussions in the text of “Practical Time-Frequency Analysis”.
Logarithms of the prices of a contract of Japanese yen.
data(YN)
data(YN)
A vector containing 500 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Daily differences of YN
.
data(YNdiff)
data(YNdiff)
A vector containing 499 observations.
See discussions in the text of “Practical Time-Frequency Analysis”.
Carmona, R. A., W. L. Hwang and B Torresani (1998) Practical Time-Frequency Analysis: Gabor and Wavelet Transforms with an Implementation in S, Academic Press, San Diego.
Generate a zero kernel for reconstruction from ridges.
zerokernel(x.inc=1, x.min, x.max)
zerokernel(x.inc=1, x.min, x.max)
x.min |
minimal value of x for the computation of |
x.max |
maximal value of x for the computation of |
x.inc |
step unit for the computation of the kernel. |
matrix of the kernel
kernel
, fastkernel
, gkernel
,
gkernel
.
Computes the the reconstructed signal from the ridge when the epsilon parameter is set to zero
zeroskeleton(cwtinput, Qinv, morvelets, bridge, aridge, N)
zeroskeleton(cwtinput, Qinv, morvelets, bridge, aridge, N)
cwtinput |
continuous wavelet transform (as the output of |
Qinv |
inverse of the reconstruction kernel (2D array). |
morvelets |
array of Morlet wavelets located at the ridge samples. |
bridge |
time coordinates of the ridge samples. |
aridge |
scale coordinates of the ridge samples. |
N |
size of reconstructed signal. |
The details of this reconstruction are the same as for the function skeleton. They can be found in the text
Returns a list of the elements of the reconstruction of a signal from sample points of a ridge
sol |
reconstruction from a ridge. |
A |
matrix of the inner products. |
lam |
coefficients of dual wavelets in reconstructed signal. They are the
Lagrange multipliers |
dualwave |
array containing the dual wavelets. |
See discussions in the text of “Practical Time-Frequency Analysis”.
skeleton
, skeleton2
,
zeroskeleton2
.
Computes the the reconstructed signal from the ridge when the epsilon parameter is set to zero, in the case of real constraints.
zeroskeleton2(cwtinput, Qinv, morvelets, bridge, aridge, N)
zeroskeleton2(cwtinput, Qinv, morvelets, bridge, aridge, N)
cwtinput |
continuous wavelet transform (output of |
Qinv |
inverse of the reconstruction kernel (2D array). |
morvelets |
array of Morlet wavelets located at the ridge samples. |
bridge |
time coordinates of the ridge samples. |
aridge |
scale coordinates of the ridge samples. |
N |
size of reconstructed signal. |
The details of this reconstruction are the same as for the function skeleton. They can be found in the text
Returns a list of the elements of the reconstruction of a signal from sample points of a ridge
sol |
reconstruction from a ridge. |
A |
matrix of the inner products. |
lam |
coefficients of dual wavelets in reconstructed signal. They
are the Lagrange multipliers |
dualwave |
array containing the dual wavelets. |
See discussions in the text of “Practical Time-Frequency Analysis”.
skeleton
, skeleton2
,
zeroskeleton
.