Title: | Discrete Trigonometric Transforms |
---|---|
Description: | This package provides functions for 1D and 2D Discrete Cosine Transform (DCT), Discrete Sine Transform (DST) and Discrete Hartley Transform (DHT). |
Authors: | Lukasz Komsta |
Maintainer: | Lukasz Komsta <[email protected]> |
License: | GPL (>= 2) |
Version: | 0.1-2 |
Built: | 2025-01-15 03:53:32 UTC |
Source: | https://github.com/cran/dtt |
This package provides functions for 1D and 2D Discrete Cosine Transform (DCT), Discrete Sine Transform (DST) and Discrete Hartley Transform (DHT).
Package: | dtt |
Type: | Package |
Version: | 0.1-1 |
Date: | 2007-02-25 |
License: | GPL version 2 or newer. |
Maintainer: Lukasz Komsta <[email protected]>
Performs univariate discrete sine, cosine or Hartley transform.
dtt(x, type = c("dct", "dst", "dht"), variant = 2, inverted = FALSE) dct(x, variant = 2, inverted = FALSE) dst(x, variant = 2, inverted = FALSE) dht(x, inverted = FALSE)
dtt(x, type = c("dct", "dst", "dht"), variant = 2, inverted = FALSE) dct(x, variant = 2, inverted = FALSE) dst(x, variant = 2, inverted = FALSE) dht(x, inverted = FALSE)
x |
a vector or matrix to be transformed |
type |
type of transform. Default "dct" is discrete cosine, "dst" is discrete sine and "dht" is discrete Hartley |
variant |
a transformation variant - 1...4 for DCT-I...DCT-IV or DST-I...DST-IV. Default is DCT-II or DST-II. Ignored when type = "dht" |
inverted |
if the inverted transform should be performed? |
This function transforms a vector of real numbers into a vector of its DCT, DST or DHT components, of the same length.
If the x is a matrix, the transform goes by rows (each row of a result is a transform of corresponding row in x).
The dct, dst and dht functions are simple wrappers for choosing the type by function name.
A transformed vector.
Lukasz Komsta
1. N. Ahmed, T. Natarajan, and K. R. Rao, "Discrete Cosine Transform", IEEE Trans. Computers, 90-93, Jan 1974. 2. S. A. Martucci, "Symmetric convolution and the discrete sine and cosine transforms", IEEE Trans. Sig. Processing SP-42, 1038-1051 (1994). 3. R. V. L. Hartley, "A more symmetrical Fourier analysis applied to transmission problems," Proc. IRE 30, 144-150 (1942).
x=seq(0,20,length=200) y=x*sin(x)+cos(x)+5*cos(10*x)+rnorm(200,sd=0.1) plot(y) z=dct(y); z[151:200]=0; lines(dct(z,inverted=TRUE),col=2); z=dct(y); z[21:200]=0; lines(dct(z,inverted=TRUE),col=4);
x=seq(0,20,length=200) y=x*sin(x)+cos(x)+5*cos(10*x)+rnorm(200,sd=0.1) plot(y) z=dct(y); z[151:200]=0; lines(dct(z,inverted=TRUE),col=2); z=dct(y); z[21:200]=0; lines(dct(z,inverted=TRUE),col=4);
Performs multivariate (2D) discrete sine, cosine or Hartley transform.
mvdtt(x, type = c("dct", "dst", "dht"), variant = 2, inverted = FALSE) mvdct(x, variant = 2, inverted = FALSE) mvdst(x, variant = 2, inverted = FALSE) mvdht(x, inverted = FALSE)
mvdtt(x, type = c("dct", "dst", "dht"), variant = 2, inverted = FALSE) mvdct(x, variant = 2, inverted = FALSE) mvdst(x, variant = 2, inverted = FALSE) mvdht(x, inverted = FALSE)
x |
a matrix to be transformed |
type |
type of transform. Default "dct" is discrete cosine, "dst" is discrete sine and "dht" is discrete Hartley |
variant |
a transformation variant - 1...4 for DCT-I...DCT-IV or DST-I...DST-IV. Default is DCT-II or DST-II. Ignored when type = "dht" |
inverted |
if the inverted transform should be performed? |
This function transforms a matrix of real numbers into a matrix of its DCT, DST or DHT components, of the same dimensions. It is done by so-called row-matrix algorithm.
The mvdct, mvdst and mvdht functions are simple wrappers for choosing the type by function name.
A transformed matrix.
Lukasz Komsta
1. N. Ahmed, T. Natarajan, and K. R. Rao, "Discrete Cosine Transform", IEEE Trans. Computers, 90-93, Jan 1974. 2. S. A. Martucci, "Symmetric convolution and the discrete sine and cosine transforms", IEEE Trans. Sig. Processing SP-42, 1038-1051 (1994). 3. R. V. L. Hartley, "A more symmetrical Fourier analysis applied to transmission problems," Proc. IRE 30, 144-150 (1942).
x = rnorm(100); dim(x) = c(10,10); x mvdct(x) mvdct(mvdct(x),inverted=TRUE)
x = rnorm(100); dim(x) = c(10,10); x mvdct(x) mvdct(mvdct(x),inverted=TRUE)