next up previous contents index
Next: 5. DSP Function Summary Up: UDI 2.1A Unified DSP Previous: 3.2 Examples

4. The DSP Functions

 acort(a, b, s)
  UDI_object a; /* Input Vector */
  UDI_object b; /* Output Correlation vector */
  short s; /* Flag for the choice of the estimator */


acort() performs an auto-correlation of n=length(b), using one of the three time-domain estimators specified by the flag s, on the data vector a of length m=length(a). The lags are from 0 to n-1.

for(j=0; j<m; j++) b[j] = 0.0;

for(h=0; h<=n; h++)
{
for(j=0; j< (m-h); j++)
b[k] += (a[k] * a[k + h]) / scale_factor;
}

The estimators chosen by s differ by the value of the scale factor:

Estimator flag s scale factor
Energy 0 1.0
Bartlett 1 1.0/(m-h)
Blackmann-tukey 2 1.0/m

   


The Energy estimator gives an estimation which has the dimension of an energy. The two others, Blackmann-tukey and Bartlett, give an estimation which has the dimension of a power. In statistics terms, the Bartlett estimator is an unbiased one, and the Energy and Blackmann-tukey estimators are biased, but the Energy and the Blackmann-tukey estimators have less variances for the high order of lags.

If you are not familiar with these concepts, use the Energy estimator.


 autoyw(a, b, c, d)
  UDI_object a; /* Real Autocorrelation coefficient vector */
  UDI_object b; /* Real output scalar (variance estimation) */
  UDI_object c; /* Real regression coefficient vector */
  short d; /* */


autoyw() calculates the regression coefficient vector c for an AR model of a stochastic process of order n-1, where n is the size of the input autocorrelation coefficient vector a. This is the solution for the Yule-Walker equations when the autocorrelation matrix is a Toeplitz one. The algorithm returns the variance of the error of the estimation in the scalar b. Note that the vectors a and c are of size n. autoyw() does not support increments.


 bkman(a, m)
  UDI_object a; /* Output vector */
  short m; /* Flag for the desired Blackman window */


bkman() computes one of the six Blackman windows specified by the flag m, placing the result in the vector a according to the equation:

a[k] = a0 + a1 * cos(2*PI*k/n) + a2 * cos(4*PI*k/n) + a3 * cos(6*PI*k/n) for k = 0 to n-1

The coefficients a0, a1, a2 and a3 of the desired Blackman window are listed in the following table:

m window a0 a1 a2 a3
0 Exact blackman 0.42659 -0.49656 0.07685 0.0
1 Blackman 0.42 -0.50 0.08 0.0
2 Blackman-Harris 3 (-67 dB) 0.42323 -0.49755 0.07922 0.0
3 Blackman-Harris 3 (-61 dB) 0.44959 -0.49364 0.05677 0.0
4 Blackman-Harris 4 (-92 dB) 0.35875 -0.48829 0.14128 -0.01168
5 Blackman-Harris 4 (-74 dB) 0.40217 -0.49703 0.09392 -0.00183

         


 cdotpr(a, b, c)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Second complex input vector */
  UDI_object c; /* Complex dot-product output scalar */


cdotpr() computes the complex dot-product of the two vectors a and b, placing the result in the scalar described by c. A C-like expression of the complex dot product is :

c = 0.0;

for(k = 0; k < length(a); k++)
c += a[k] * b[k];
for k = 0 to n-1


 cfft(a, b, f)
  UDI_object a; /* Complex input vector to be transformed */
  UDI_object b; /* Complex weights vector */
  short f; /* Flag for direct/inverse FFT */


cfft() computes an in-place FFT of the complex input vector a if the flag f is 1, and an in-place inverse FFT if f is -1. It assumes that the length of a is a power of two. Before you use cfft() for the first time, you have to create the FFT weight vector b with fftwts().


 cholesky(r, a, t, d)
  UDI_object r; /* Square input matrix */
  UDI_object a; /* First output matrix */
  UDI_object t; /* Second output matrix */
  UDI_object d; /* Output vector */


cholesky() performs a Cholesky decomposition of a positive definite matrix R = A * D * T, where A is lower triangular with diagonal elements set to one. T is the transpose of A, and D is a diagonal matrix. In fact D is returned as a vector.


 convt(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


convt() performs a convolution of the two real input vectors a and , placing the result in the real output vector c. Note that the length of output vector c must be equal to the sum of the lengths of input vectors a and b minus one.


 cpvmul(a, b, c, f)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Second complex input vector */
  UDI_object c; /* Complex output vector */
  short f; /* Conjugation flag (1=normal, -1=conjugate) */


cpvmul() multiplies each corresponding element of the two polar coordinate complex input vectors a and b, placing the result in polar coordinate complex output vector c. The vectors are considered as complex vectors each containing pairs of data representing amplitude and phase. The corresponding amplitude values are multiplied, while the phase values are either added (f=1), or subtracted (f=-1 indicating subtraction of the phase component of b).

real(c[k]) = real(a[k]) * real(b[k])
imag(c[k]) = imag(a[k]) + f * imag(b[k])
for k = 0 to length(a)


 cpvsmml(a, b, c, d)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Real scalar input */
  UDI_object c; /* Second complex input vector */
  UDI_object d; /* Complex output vector */


cpvsmml() multiplies each element of the polar coordinate complex vector a by the real scalar b and by the corresponding element of the second polar coordinate complex vector c, placing the result in the polar coordinate complex vector d according to the equation :

real(d[k]) = b * real(a[k]) * real(c[k])
imag(d[k]) = imag(a[k]) + imag(c[k])
for k = 0 to length(a)


 creorgo(a, b, c)
  UDI_object a; /* Complex input vector to be sorted */
  UDI_object b; /* Sorting addresses (REAL_LONG) */
  UDI_object c; /* Complex output vector */


creorgo() sorts the complex elements of a by putting them in the complex output vector c at a place specified by the sorting address stored in b. A C-like expression of the process is:

c[b[k]] = a[k] for k = 0 to length(a)-1


 crergoa(a, b, c)
  UDI_object a; /* Complex input vector to be sorted */
  UDI_object b; /* Sorting addresses (REAL_LONG) */
  UDI_object c; /* Complex output vector */


crergoa() sorts the complex elements of a by adding them to the element of the complex output vector c which is at the place specified by the corresponding address stored in b and finally puts the result at the same place in the complex output vector c. A C-like expression of the process is:

c[b[k]] += a[k] for k = 0 to length(a)-1


 cvadd(a, b, c)
  UDI_object a; /* First complex input vector to be added */
  UDI_object b; /* Second complex input vector to be added */
  UDI_object c; /* Output vector */


cvadd() adds the two complex input vectors a and b, placing the result in the complex vector c.


 cvcomb(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Complex output vector */


cvcomb() combines two real vectors to create a complex vector c. a is taken as the real part, and b as the imaginary part.


 cvconj(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


cvconj() conjugates each element of the complex input vector a, placing the result in the complex output vector b.


 cvcsml(a, b, c)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex input scalar */
  UDI_object c; /* Complex output vector */


cvcsml() multiplies each element of the complex input vector a by the complex scalar b, placing the result in the complex output vector c.


 cvdiv(a, b, c)
  UDI_object a; /* Complex denominators input vector */
  UDI_object b; /* Complex numerators input vector */
  UDI_object c; /* Complex output vector */


cvdiv() divides each element of the complex vector b by the corresponding element of the complex vector a, placing the result in the complex vector c.


 cvexp(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Complex output vector */


cvexp() computes the complex exponential value of the phase-elements (expressed in radians) of the real vector a, placing the result in the complex vector b. A C-like expression of the process is :

b[k] = exp(i * a[k]) for k = 0 to length(a)-1


 cvmags(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Real squared magnitude output vector */


cvmags() computes the squared magnitudes of the elements of the vector a, placing the result in the real vector b according to the equation:

b[k] = real(a[k]) * real(a[k]) + imag(a[k]) * imag(a[k])
for k = 0 to length(a)-1


 cvmov(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


cvmov() copies the elements of the complex input vector a into the complex output vector b.


 cvmul(a, b, c, m)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Second complex input vector */
  UDI_object c; /* Complex output vector */
  short m; /* Conjugate flag */


cvmul() multiplies the corresponding complex elements of the vectors a and b, placing the complex result in the vector c according to the equation :

c[k] = a[k] * b[k] /* if m = 1 */
c[k] = conj(a[k]) * b[k] /* if m = -1 */
for k = 0 to length(a)-1

The flag m indicates the conjugation of the first vector a by multiplying the imaginary part of its elements, so be careful that m takes the value 1 or -1... The C-version takes care of this.


 cvneg(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


cvneg() evaluates the negative value of the complex elements of the vector a, placing the result in the complex vector b.


 cvinv(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


cvinv() computes the reciprocal of each element of the complex vector a, placing the result in the complex vector b according to the equation :

b[k] = 1.0 / a[k] for k = 0 to length(a)-1


 cvrev(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


cvrev() reverses the order of the elements of complex vector a, placing the result in the complex output vector b according to the equation :

b[k] = a[n - k - 1] for k = 0 to n-1 (where n = length(a) )


 cvcsma(a, b, c, d)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Complex scalar input */
  UDI_object c; /* Second complex input vector */
  UDI_object d; /* Complex output vector */


cvcsma() multiplies each element of the complex vector a by the complex scalar b, then adds the corresponding elements of the second complex vector c, placing the result in the complex vector d according to the equation :

d[k] = a[k] * b + c[k] for k = 0 to length(a)-1


 cvsmul(a, b, c)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Real scalar input */
  UDI_object c; /* Complex output vector */


cvsmul() multiplies each element of the complex input vector a by the real scalar b, placing the result in the complex output vector c. A C-like expression of this scalar multiply is :

c[k] = a[k] * b for k = 0 to length(a)-1


 cvsub(a, b, c)
  UDI_object a; /* First complex input vector */
  UDI_object b; /* Second complex input vector */
  UDI_object c; /* Complex output vector */


cvsub() substracts each complex element of the complex vector b from each corresponding element of the complex vector a, placing the result into the complex output vector c according to the equation :

c[k] = a[k] - b[k] for k = 0 to length(a)-1


 dotpr(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output scalar */


dotpr() computes the real dot-product of the two vectors a and b, placing the result in the scalar described by c. The dot-product is computed by summing the element-wise product of a and b.


 fftwts(a, n)
  UDI_object a; /* Complex vector (roots of unity) */
  long n; /* FFT size */


fftwts() computes the weight values needed to run an FFT algorithm of size n. The values are placed in the complex vector a. A C-like expression for the weight values is :

a[k] = exp(2 * PI * i * k / n) for k = 0 to n-1 with PI = 3.1415926...

Be sure to call fftwts the first time before calling any of the UDI FFT functions on a vector of a given size n.


 hamm(a)
  UDI_object a; /* Output vector */


hamm() computes a Hamming window and puts it in the vector a according to the equation:

a[k] = 0.54 - 0.46 * cos(2*PI*k/n) for k = 0 to n=length(a)-1


 hann(a)
  UDI_object a; /* Output vector */


hann() computes the Hanning window of the length n and puts it in the vector a according to the equation:

a[k] = 0.5 - 0.5 * cos(2*PI*k/n) for k = 0 to n=length(a)-1


 hist(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar (lower threshold) */
  UDI_object c; /* Real input scalar (upper threshold) */
  UDI_object d; /* Output vector (REAL_LONG) */


hist() performs a histogram on the input vector a. It adds 1.0 to the "bin" of the histogram output vector d which corresponds to the value of the vector a quantized according to the size of d, the upper threshold value c and the lower thresold value b. A C-like expression for the process is:

for(k = 0; k < length(a); k++)
{
if(a[k] < b) d[0] += 1.0;
else {
if(a[k] > c) d[m - 1] += 1.0;
else {
j = (long)(length(d) * (a[k] - b) / (c - b));
d[j] += 1.0;
}
}
}

The histogram output vector is not initialized to zero when hist is called, so you can accumulate histograms or initialize it each time it's called with vclr().


 ilevins(a, k)
  UDI_object a; /* Input AR filter coefficient vector */
  UDI_object k; /* Output reflection coefficient vector */


ilevins() performs the recursive step-down levinson algorithm on the AR-filter coefficients contained in real input vector a, placing the corresponding reflection coefficients in the real output vector k. length(k) must be greater than or equal to length(a)+1 and the first element of a must be set to 1.0. If any reflection coefficient is found to be more than one in magnitude a warning message is displayed on the standard error output.


 levins(k, a)
  UDI_object k; /* Input reflection coefficient vector */
  UDI_object a; /* Output AR filter coefficient vector */


levins() performs the recursive levinson algorithm on the reflection coefficients in real input vector k, placing the corresponding AR-filter coefficients in real output vector a. length(k) must be greater than or equal to length(a)+1. The reflection coefficients must be of magnitude less than one so that the corresponding AR-filter will be stable. The first element of a will be set to 1.0.


 lveq(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


lveq() compares the corresponding elements of each of the two vectors a and b, and sets the corresponding element of c to 1.0 if the element of a is equal to the element of b and to 0.0 if not. A C-like expression of the process is :

c[k] = (a[k] == b[k]) ? 1.0 : 0.0 for k = 0 to length(a)-1


 lveqzs(a, b)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Real output vector */


lveqzs() compares each element of the vector a with zero and sets the corresponding element of b to 1.0 if it equals zero and leaves it unchanged if not. A C-like expression of the process is :

b[k] = (a[k] == 0.0 ) ? 1.0 : a[k] for k = 0 to length(a)-1


 lvge(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


lvge() compares the corresponding elements of each of the two vectors a and b, and sets the corresponding element of c to 1.0 if the element of a is greater or equal to the element of b and to 0.0 if not. A C-like expression of the process is :

c[k] = (a[k] >= b[k]) ? 1.0 : 0.0 for k = 0 to length(a)-1


 lvgt(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


lvgt() compares the corresponding elements of each of the two vectors a and b, and sets the corresponding element of c to 1.0 if the element of a is greater than the element of b and to 0.0 if not. A C-like expression of the process is :

c[k] = (a[k] > b[k]) ? 1.0 : 0.0 for k = 0 to length(a)-1


 lvne(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


lvne() compares the corresponding elements of each of the two vectors a and b, and sets the corresponding element of c to 0.0 if the element of a is equal to the element of b and sets it to 1.0 if not. A C-like expression of the process is :

c[k] = (a[k] != b[k]) ? 1.0 : 0.0 for k = 0 to length(a)-1


 lvnot(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


lvnot() make a logical not test on each element of the vector a and sets the corresponding element of b to 1.0 if the element of a is equal to 0.0. Otherwise, it's set to 0.0. A C-like expression of the process is :

b[k] = (a[k] == 0.0) ? 1.0 : 0.0 for k = 0 to length(a)-1


 lvsmth(a, b, c, m)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar (threshold) */
  UDI_object c; /* Real output vector */
  long m; /* Cutoff number of a 0- or 1-domain. */


lvsmth() smooths the real input vector a by first thresholding each element by real input scalar b. Elements less than or equal to the threshold receive a value of 0.0, otherwise they receive a value of 1.0. Then contiguous runs of 0.0 less than or equal to the cutoff number m are changed to 1.0. In the final pass, contiguous runs of 1.0 less than or equal to m are changed to 0.0.


 matinv(a, b, c)
  UDI_object a; /* Input matrix to be inverted */
  UDI_object b; /* Real input scalar (threshold) */
  UDI_object c; /* Inverted matrix */


matinv() computes the inverse matrix of the input matrix a using the JORDAN algorithm and places the result in the output matrix c. The real input scalar, b, is a threshold used to determine if the matrix is singular. The default threshold is 0.0001 which is used if Krf_0 is passed as the second argument.


 maxmgv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


maxmgv() finds the maximum magnitude value of the elements of the real vector a and sets the real scalar b to this value. An expression for the process is :

b = MAX(
| a[k] | for k = 0 to length(a)-1)

 maxmgvi(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar (max value) */
  UDI_object c; /* Long output scalar (index) */


maxmgv() finds the maximum magnitude value of the elements of the real vector a and sets real scalar b to this value and scalar c to its location. An expression for the process is :

b = MAX(
| a[k] | for k = 0 to length(a)-1)
c = m for which a[m] = MAX(
| a[k] | for k = 0 to length(a)-1)

 maxv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


maxv() finds the maximum value of the elements of the real vector a and sets real scalar b to this value. An expression for the process is :

b = MAX(a[k] for k = 0 to length(a)-1)


 maxvi(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar (max value) */
  UDI_object c; /* Long output scalar (index) */


maxvi() finds the maximum value of the elements of the real vector a and sets real scalar b to this value and scalar c to its location. An expression for the process is :

b = MAX(a[k] for k = 0 to length(a)-1)
c = m for which a[m] = MAX(a[k] for k = 0 to length(a)-1)


 meamgv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


meamgv() computes the mean value of the magnitudes of the elements of the real vector a and sets real scalar b to this value. A C-like expression of the process is :

b = 0.0;

b +=
| a[k] | / length(a) for k = 0 to length(a)-1

 meanv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


meanv() computes the mean value of the elements of the real vector a and sets real scalar b to this value. A C-like expression of the process is :

b = 0.0;

b += a[k] / length(a) for k = 0 to length(a)-1


 measqv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


measqv() computes the mean value of the square of the elements of the real vector a and sets real scalar b to this value. A C-like expression of the process is :

b = 0.0;

b += (a[k] * a[k]) / length(a) for k = 0 to length(a)-1


 minmgv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


minmgv() finds the minimum magnitude value of the elements of the real vector a and sets real scalar b to this value. An expression for the process is :

b = MIN(
| a[k] | for k = 0 to length(a)-1 )

 minmgvi(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar (min value) */
  UDI_object c; /* Long output scalar (index) */


minmgvi() finds the minimum magnitude value of the elements of the real vector a and sets real scalar b to this value and scalar c to its location. An expression for the process is :

b = MIN(
| a[k] | for k = 0 to length(a)-1 )
c = m for which a[m] = MIN(
| a[k] | for k = 0 to length(a)-1 )

 minv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


minv() finds the minimum value of the elements of the real vector a and sets real scalar b to this value. An expression for the process is :

b = MIN(a[k] for k = 0 to length(a)-1 )


 minvi(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar (min value) */
  UDI_object c; /* Long output scalar (index) */


minvi() finds the minimum value of the elements of the real vector a and sets real scalar b to this value and scalar c to its location. An expression for the process is :

b = MIN(a[k] for k = 0 to length(a)-1 )
c = m with a[m] = MIN(a[k] for k = 0 to length(a)-1)


 mmul(a, b, c)
  UDI_object a; /* First input matrix */
  UDI_object b; /* Second input matrix */
  UDI_object c; /* Output product matrix */


mmul() multiplies the input matrix a by the input matrix b and places the result in the output matrix c.


 mtrans(a, c)
  UDI_object a; /* Input matrix to be transposed */
  UDI_object c; /* Output transposed matrix */


mtrans() transposes the input matrix a and places the result in the output matrix c.


 mvmul(a, b, c)
  UDI_object a; /* Real input matrix */
  UDI_object b; /* Real input vector */
  UDI_object c; /* Real output product vector */


mvmul() computes the product of real input matrix a and real input vector b, placing the result in real output vector c.


 polar(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


polar() converts the rectangular coordinates contained in the vector a and fills b with the corresponding polar coordinates. An expression of the process is :

real(b[k]) = sqrt(real(a[k])*real(a[k]) + imag(a[k])*imag(a[k]))
imag(b[k]) = arctan(imag(a[k]) / real(a[k]))
for k = 0 to length(a)


 rblsq(a, b, c, p)
  UDI_object a; /* Real input vector of data observation */
  UDI_object b; /* Input matrix of least-square design */
  UDI_object c; /* Real output vector of parameters */
  short p; /* Flag of matrix inversion */


rblsq() computes the optimal bloc least-square solution C of the following linear problem

a = Bc + noise

where a is an observation vector, B is a design matrix of the least-square problem, and c is the desired vector of parameters. Flag p allows one to choose the matrix inversion which is involved in the procedure. At this time, this flag has no effect. Be sure, that the size of a is greater than the size of c.


 rect(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Complex output vector */


rect() converts the polar coordinates contained in the vector a and fills b with the corresponding rectangular coordinates. An expression of the process is :

real(b[k]) = real(a[k]) * cos(imag(a[k]))
imag(b[k]) = real(a[k]) * sin(imag(a[k]))
for k = 0 to length(a)



 reorgo(a, b, c)
  UDI_object a; /* Real input vector to be sorted */
  UDI_object b; /* Sorting addresses (REAL_LONG) */
  UDI_object c; /* Real output vector */


reorgo() sorts the elements of the vector a by putting them in the output vector c at a place specified by the sorting address stored in b. A C-like expression of the process is:

c[b[k]] = a[k] for k = 0 to length(a)-1


 rfft(a, b, m)
  UDI_object a; /* Real vector */
  UDI_object b; /* Complex vector */
  short m; /* Flag for Normal/Inverse (1 / -1) */


rfft() computes an out of place forward or inverse FFT. If the flag m equals 1, rfft() assumes that the vector a contains n=length(a) real data and computes the n/2+1 significant complex values of the transform and stores them in complex vector b (length(b) = length(a)/2) in a packed form (the last real value is stored in the first imaginary position). If the flag m equals -1 (inverse), rfft() assumes that the vector b contains n/2+1 complex data which correspond to an FFT-transform (in packed form) and computes the n real values of the inverse FFT placing the result in vector a. rfft() assumes that n (the size of a) is a power of two.


 rfftsig(a, b, w, m)
  UDI_object a; /* Real vector */
  UDI_object b; /* Complex vector */
  UDI_object w; /* Weights vector (obtained by fftwts) */
  short m; /* Flag for Normal/Inverse (1 / -1) */


rfftsig() is identical to rfft(), except that only the n/2 complex values are used in b (length(b) = length(a)/2). The last complex value is considered to be zero, since in signal processing applications there should be extremely little energy in this bin.


 rmsqv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


rmsqv() computes the square root of the mean value of the squares of the elements of the real vector a and sets real scalar b to this value. An expression for the process is :

b = sqrt( sum(a[k] * a[k]) / length(a) for k = 0 to length(a)-1 )


 sve(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


sve() computes the sum of the elements of the real vector a and sets real scalar b to this value.


 svemg(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


svemg() computes the sum of the absolute values of the elements of the real vector a and sets real scalar b to this value.


 svesq(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


svesq() computes the sum of the squares of the elements of the real vector a and sets real scalar b to this value.


 svs(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */


svs() computes the sum of the signed squares of the elements of the real vector a and sets real scalar b to this value. A C-like expression of the process is :

b = 0.0;
b += (a[k] * abs(a[k])) for k = 0 to length(a)-1


 trian(a)
  UDI_object a; /* Output vector */


trian() computes a triangular window and puts it in the real output vector a. The first element of the vector a will contain 0.0, with successive elements linearly increasing to a peak value of 1.0 in the middle of the vector. Then, the elements linearly decrease, with the last element equal to the second element.


 trianginf(a, b, x)
  UDI_object a; /* inferior triangular matrix */
  UDI_object b; /* column vector */
  UDI_object x; /* output column vector */


Resolution of the matrix equation A * x = b, where A is a lower (inferior) triangular matrix. Can be used with Cholesky decomposition.


 triangsup(a, b, x)
  UDI_object a; /* superior triangular matrix */
  UDI_object b; /* column vector */
  UDI_object x; /* output column vector */


Resolution of the matrix equation A * x = b, where A is an upper (superior) triangular matrix. Can be used with Cholesky decomposition.


 csort(a, b, c)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Sorting address vector (long) */
  UDI_object c; /* Complex output vector */


csort() selects the elements of the complex vector a by the indices stored in the vector b and puts them in the complex output vector c. A C-like expression of the process is:

c[k] = a[b[k]] for k = 0 to length(a)-1


 sort(a, b, c)
  UDI_object a; /* Real input vector to be sorted */
  UDI_object b; /* Sorting address vector (long) */
  UDI_object c; /* Real output vector */


sort() selects the elements of the real vector a by the indices stored in the vector b and puts them in the real output vector c. A C-like expression of the process is:

c[k] = a[b[k]] for k = 0 to length(a)-1


 vabs(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vabs() computes the absolute value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = fabs(a[k]) for k = 0 to length(a)-1


 vacos(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vacos() computes the arc-cosine value of the elements of real vector a, placing the result in real output vector b. A C-like expression of the process is :

b[k] = acos(a[k]) for k = 0 to length(a)-1


 vacosh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vacosh() computes the hyperbolic arc-cosine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = acosh(a[k]) for k = 0 to length(a)-1


 vadd(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Output vector */


vadd() adds two real vectors a and b, placing the result in the real vector c according to the equation :

c[k] = a[k] + b[k] for k = 0 to length(a)-1


 vaddosb(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Output vector */


vaddosb() computes the quotient of the sum of real vectors a and b divided by the difference of a and b. The result is put in real vector c. A C-like expression of the process is :

c[k] = (a[k] + b[k]) / (a[k] - b[k]) for k = 0 to length(a)-1


 valog(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


valog() computes the antilog value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = pow(10.0, a[k]) for k = 0 to length(a)-1


 vam(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Third real input vector */
  UDI_object d; /* Output vector */


vam() adds the corresponding elements of the real vectors a and b, multiplies the sums by the corresponding elements of c, and then puts the result in the vector d. A C-like expression of the process is :

d[k] = (a[k] + b[k]) * c[k] for k = 0 to length(a)-1


 vasin(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vasin() computes the arc-sine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = asin(a[k]) for k = 0 to length(a)-1


 vasinh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vasinh() computes the hyberbolic arc-sine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = asinh(a[k]) for k = 0 to length(a)-1


 vatan(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vatan() computes the arc-tangent value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = atan(a[k]) for k = 0 to length(a)-1


 vatan2(a, b, c)
  UDI_object a; /* Real denominators input vector */
  UDI_object b; /* Real numerators input vector */
  UDI_object c; /* Real output vector */


vatan2() computes the arc-tangent of the quotients of the corresponding elements of the vectors b and a, while taking into account their signs. The result is placed in the vector c. A C-like expression of the process is :

c[k] = atan2(b[k], a[k]) for k = 0 to length(a)-1


 vatanh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vatanh() computes the hyperbolic arc-tangent value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = atanh(a[k]) for k = 0 to length(a)-1


 vclip(a, b, c, d)
  UDI_object a; /* Real input vector to be clipped */
  UDI_object b; /* Real scalar (lower threshold) */
  UDI_object c; /* Real scalar (upper threshold) */
  UDI_object d; /* Real output data */


vclip() compares each element of the real vector a with the upper and lower thresholds. If the element lies between the lower threshold b and the upper threshold c then it remains unchanged. If the element is greater than the upper threshold c, then c is taken as the value. If the element is less than the lower threshold b, then b is taken as the value. The result is then stored in the vector d. A C-like expression of the process is :

if (a[k]) < b) d[k] = b;
else if (a[k]) > c) d[k] = c;
else d[k] = a[k];
for k = 0 to length(a)-1


 vclr(a)
  UDI_object a; /* Real vector to be cleared */


vclr() sets all the elements of real vector a to 0.0. A C-like expression of the process is :

a[k] = 0.0 for k = 0 to length(a)-1


 vcos(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vcos() computes the cosine value of the elements of real vector a, placing the result in real vector b. A C-like expression of the process is :

b[k] = cos(a[k]) for k = 0 to length(a)-1


 vcosh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vcosh() computes the hyperbolic cosine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = cosh(a[k]) for k = 0 to length(a)-1


 vcvmul(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Complex input vector */
  UDI_object c; /* Complex output vector */


vcvmul() computes the element-wise vector product between real input vector a and complex input vector b, placing the result in complex output vector c.


 vdbsafe(a, b, s)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */
  short s; /* Flag for choosing the kind of decibel. */


vdbsafe() converts each element of the input real vector a into decibels, placing the result in real vector b. A C-like expression of the process is :

b[k] = (a[k] > linear_limit) ? scale * log10(a[k]) : log_limit for k = 0 to length(a)-1

where scale, linear_limit and log_limit are chosen by the flag s according to the nature of the input data and of the output data.

Decibel conversions Flag s scale linear limit log limit
Ampl to ampl 0 10.0 4.539989e-5 -100
Ampl to pow 1 20.0 4.539989e-5 -200
Pow to pow 2 10.0 2.06115e-9 -200
Pow to ampl 3 5.0 2.06115e-9 -100

       


 vdist(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vdist() computes the square root of the sum of the squares of each element of the two real vectors a and b, placing the result in the real vector c according to the equation :

c[k] = sqrt(a[k] * a[k] + b[k] * b[k]) for k = 0 to length(a)-1


 vdiv(a, b, c)
  UDI_object a; /* Real input vector (denominators) */
  UDI_object b; /* Real input vector (numerators) */
  UDI_object c; /* Real output vector */


vdiv() divides each element of the real vector b by the corresponding element of the real vector a, placing the result in the real vector c according to the equation :

c[k] = b[k] / a[k] for k = 0 to length(a)-1


 vexp(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vexp() computes the exponential value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = exp(a[k]) for k = 0 to length(a)-1


 vexpdev(a, m, p)
  UDI_object a; /* Real output vector */
  long m; /* Algorithm type */
  long p; /* Seed value */


vexpdev() generates a vector a of real random numbers distributed according to a Gaussian distribution with mean zero and variance one. Algorithm type m should be between 0 and 3 inclusive, and represents the algorithm used to calculate a uniformly distributed real number which is the basis for the final gaussian distributed result. p represents a seed value for the random number generator, and must be different to insure a different result each time that vexpdev() is called.


 vfill(a, b)
  UDI_object a; /* Real input scalar */
  UDI_object b; /* Real output vector */


vfill() sets all the elements of real vector b to the real scalar value a. A C-like expression of the process is :

b[k] = a for k = 0 to length(a)-1


 vfiliir1(a, b, m, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar (AR coefficient) */
  UDI_object m; /* Real scalar (filter memory) */
  UDI_object c; /* Real output vector */


vfilt1() performs a first order autoregressive filtering with unity gain on the input signal contained in real vector a, placing the filtered signal in the real output vector c. The real input scalar b is the coefficient of the AR filter, and the real scalar m is the memory of the filter (the c[-1] element) to insure the continuity when successive blocks of signal are processed. This scalar should be set to zero the first time the filter is called. When vfilt1() returns, the scalar contain the last output value from the filter.
A C-like expression of the process is :

c[k] = a[k] - (c[k - 1] * b) for k = 0 to length(a)-1


 vfiliir2(a, b, c, m, n, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar (AR coefficient 1) */
  UDI_object c; /* Real input scalar (AR coefficient 2) */
  UDI_object m; /* Real scalar (filter memory -1) */
  UDI_object n; /* Real scalar (filter memory -2) */
  UDI_object d; /* Real output vector */


vfilt2() performs a second-order autoregressive filtering with unity gain on the real input signal contained in a, placing the real filtered signal in the real output vector d. The real input scalars b and c are, respectively, the first and the second coefficients of the AR filter. The real scalars m and n are the memory of the filter (the d[-1] and d[-2] values, respectively). The memory values are used to insure the continuity of the filtering when successive blocks of signal are processed. When the filter is first called, these scalars should be set to zero. When vfilt2() returns, these scalars contain the last two output values from the filter.
A C-like expression of the process is :

d[k] = a[k] - (d[k - 1] * b) - (d[k - 2] * c) for k = 0 to length(a)-1


 vshort(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Output short integer vector */


vshort() converts all the elements of the real vector a to short integers, placing the result in the vector b.


 vlong(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Output long integer vector */


vlong() converts all the elements of the real vector a to long integers, placing the result in the vector b.


 vfloat(a, b)
  UDI_object a; /* Short integer input vector */
  UDI_object b; /* Real output vector */


vfloat() converts all the elements of the vector a to float, placing the result in the real vector b. A C-like expression of the process is :

b[k] = (float)a[k] for k = 0 to length(a)-1


 vlfloat(a, b)
  UDI_object a; /* Long integer input vector */
  UDI_object b; /* Real output vector */


vlfloat() converts all the elements of a to float, placing the result in b. A C-like expression of the process is :

b[k] = float(a[k]) for k = 0 to length(a)-1


 vfrac(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vfrac() extracts the fractional part of each element of the real vector a, placing the result in b. A C-like expression of the process is :

b[k] = a[k] - (float)((long)a[k]) for k = 0 to length(a)-1


 vgeomea(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output geometric mean vector */


vgeomea() computes the geometric mean of each corresponding elements of the real vectors a and b, placing the result in the real vector c.

A C-like expression of the process is :

c[k] = sqrt(a[k] * b[k]) for k = 0 to length(a)-1


 vgeoser(a, b, c)
  UDI_object a; /* Real input scalar (initial value) */
  UDI_object b; /* Real input scalar (coefficient) */
  UDI_object c; /* Real output vector */


vgeoser() generates a geometric series of initial value a and of multiplicative coefficient b and stores this series in the real vector c. A C-like expression of the process is :

c[k] = a * pow(b,k) for k = 0 to length(a)-1
If you want to initalize this series with k = 1 instead of 0, then a and b must be the same value.


 vharmea(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output harmonic mean vector */


vharmea() computes the harmonic mean of each corresponding element of the real vectors a and b, placing the result in the real vector c.

A C-like expression of the process is :

c[k] = (a[k] * b[k]) / (a[k] + b[k]) for k = 0 to length(a)-1


 vicexpr(a, b, c, d)
  UDI_object a; /* Real input vector (numerator) */
  UDI_object b; /* Real input vector (Part of denominator) */
  UDI_object c; /* Real input vector (Part of denominator) */
  UDI_object d; /* Real output vector */


vicexpr() computes an intercorrelation-index like expression between the real vectors a, b, and c, placing the result in the real vector d.

A C-like expression of the process is :

d[k] = a[k] / sqrt(b[k] * c[k]) for k = 0 to length(a)-1


 vimag(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Real output vector */


vimag() extracts the imaginary part of the complex elements of a, placing the result in the real vector b.


 vinter(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real interpolation vector */
  UDI_object c; /* Real input vector */
  UDI_object d; /* Real output vector */


vinter() computes the real linear interpolation vector d of the two real vectors a and c with the real interpolation vector b. An expression for the process is :

d[k] = (1.0 - b[k])*a[k] + b[k]*c[k] for k = 0 to length(a)-1

The elements of b should be between 0.0 and 1.0. A value of 0.0 means that d = a. A value of 1.0 means that d = c


 vinv(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vinv() computes the reciprocal value of the real elements of a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = 1.0 / a[k] for k = 0 to length(a)-1


 vinvsafe(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vinv() computes the reciprocal value of the real elements of a, placing the result in the real vector b. If the element is zero, then a very large value is returned as its reciprocal, instead of your program crashing. A C-like expression of the process is :

b[k] = 1.0 / a[k] for k = 0 to length(a)-1


 vj0(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vj0() computes the zero order Bessel function of the first kind for each element of the real vector a, placing the result in the vector b. A C-like expression of the process is :

b[k] = j0(a[k]) for k = 0 to length(a)-1


 vj1(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vj1() computes the first order Bessel function of the first kind for each element of the real vector a, placing the result in the vector b. A C-like expression of the process is :

b[k] = j1(a[k]) for k = 0 to length(a)-1


 vjn(a, b, p)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */
  long p; /* Order of the Bessel Function */


vjn() computes the Bessel function of the first kind and of the order p for each element of the real vector a, placing the result in the vector b. A C-like expression of the process is :

b[k] = jn(p, a[k]) for k = 0 to length(a)-1


 vlicomb(a, b, c, d, e)
  UDI_object a; /* First real input vector */
  UDI_object b; /* First real input scalar */
  UDI_object c; /* Second real input scalar */
  UDI_object d; /* Second real input vector */
  UDI_object e; /* Real output vector */


vlicomb() computes a linear combination of the first input vector a multiplied by the scalar b with the second vector d multiplied by the scalar c, and finally puts the result in the output vector d. An expression for the process is :

e[k] = (a[k] * b) + (d[k] * c) for k = 0 to length(a)-1


 vlim(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* First real input scalar */
  UDI_object c; /* Second real input scalar */
  UDI_object d; /* Real output vector */


vlim() compares each element of the real input vector a with the input scalar b. If it's greater or equal to b, then the result is input scalar c. Otherwise, the result is -c.


 vlisnl(a, b, n)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */
  long n; /* Order of the filter */


vlisnl() performs an n-order median filtering of the real input vector a, placing the result in the real output vector b. The filter order n should be an odd number greater than one, and represents the number of samples for which to calculate the median value. This median value result is then placed in the output vector at the position representing the center of the neighborhood for which it was calculated. For example, if n is five, then each median value is calculated for the two elements on the left, the two elements on the right, and the center element of the input vector. The first and last (n-1)/2 elements of the input vector remain unchanged in the output vector.


 vln(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vln() computes the natural logarithm value of the real elements of a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = log(a[k]) for k = 0 to length(a)-1


 vlog(a, b)
  UDI_object a; /* Real input data */
  UDI_object b; /* Real output data */


vlog() computes the base 10 logarithm value of the real elements of a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = log10(a[k]) for k = 0 to length(a)-1


 vma(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Third real input vector */
  UDI_object d; /* Real output vector */


vma() multiplies each of the corresponding elements of the two real vectors a and b, then adds the corresponding element of the real vector c and finally puts the result in the output vector d. An expression for the process is :

d[k] = (a[k] * b[k]) + c[k] for k = 0 to length(a)-1


 vmax(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vmax() compares each corresponding element of the real vectors a and b, placing the maximum value in the real output vector c. An expression for the process is :

c[k] = MAX(a[k], b[k]) for k = 0 to length(a)-1


 vmaxmg(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vmaxmg() compares the absolute value of each corresponding element of the real vectors a and b, placing the maximum value in the real output vector c. An expression for the process is :

c[k] = MAX(fabs(a[k]), fabs(b[k])) for k = 0 to length(a)-1


 vmin(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vmin() compares each corresponding element of the real vectors a and b, placing the minimum value in the real output vector c. An expression for the process is :

c[k] = MIN(a[k], b[k]) for k = 0 to length(a)-1


 vminmg(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vminmg() compares the absolute value of each corresponding element of the real vectors a and b, placing the minimum value in the real output vector c. An expression for the process is :

c[k] = MIN(fabs(a[k]), fabs(b[k])) for k = 0 to length(a)-1


 vmmul(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input matrix */
  UDI_object c; /* Real output product vector */


vmmul() multiplies the real input row vector a by the real input matrix b, placing the result in the real output vector c.


 vmov(a, b)
  UDI_object a; /* Real source vector */
  UDI_object b; /* Real destination vector */


vmov() copies the elements of the vector a into b. An expression for the process is :

b[k] = a[k] for k = 0 to length(a)-1

Be careful! if b and a share the same physical location, data may be destroyed in the process.


 vmsa(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Second real input vector */
  UDI_object d; /* Real output vector */


vmsa() multiplies each of the corresponding elements of the two real vectors a and c, then adds the real scalar b and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (a[k] * c[k]) + b for k = 0 to length(a)-1


 vmsb(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Third real input vector */
  UDI_object d; /* Real output vector */


vmsb() multiplies each of the corresponding elements of the two real vectors a and b, then subtracts the corresponding element of the real vector c and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (a[k] * b[k]) - c[k] for k = 0 to length(a)-1


 vmul(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vmul() multiplies the corresponding elements of the two real vectors a and b, placing the result in the real vector c according to the equation:

c[k] = a[k] * b[k] for k = 0 to length(a)-1


 vneg(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vneg() evaluates the negative value of the elements of real vector a, placing the result in the real vector b. An expression for the process is :

b[k] = - a[k] for k = 0 to length(a)-1


 vpeaks(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output scalar */
  UDI_object c; /* Real output vector */
  long d; /* Minimum peak separation */


vpeaks() find peaks in the real vector a, with each peak separated by at least the number of points given by d. The output vector c contains the indices of the peaks found, and the number of peaks is placed in the output scalar b. c has to be of length at least length(a)/d and d must be between 1 and length(a)/2 - 1.

                          *
                *       *   *           *
          *   *   *   *       *   *   *   *   *     signal
        *   *       *           *   *       *   *

          +     +         +       +     +     +     d=1, 6 found
                +         +             +           d=2, 3 found
                          +                         d=3, 1 found
                                                    d=4, 0 found


 vpoly(a, b, c)
  UDI_object a; /* Real input vector of coefficients */
  UDI_object b; /* Real input vector */
  UDI_object c; /* Real output vector */


vpoly() computes the value of a polynomial function. This polynomial function is described by the real vector a of its coefficients, with the lower order coefficients appearing first. The polynomial is evaluated for each of the elements of the real vector b and the result is stored in the real vector c. A C-like expression of the process is:

for(k=0; k<length(b); k++)
{
for(q=0; q<length(a); q++)
c[k] += a[q] * pow(b[k], q)
}


 vramp(a, b, c)
  UDI_object a; /* Real input scalar (initial value) */
  UDI_object b; /* Real input scalar (increment) */
  UDI_object c; /* Real output vector */


vramp() generates a linear ramp of initial value a and of increment b and stores this ramp in the vector c. A C-like expression of the process is :

c[k] = a + b * k for k = 0 to length(a)-1


 vran0(a, m)
  UDI_object a; /* Real random output vector */
  long m; /* Seed value */


vran0() generates a vector a of real random numbers distributed according to a uniform distribution on the interval [0,1], using the algorithm of Bays and Duhram. The seed value m must be changed each time vran0() is called to insure unique results.


 vreal(a, b)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Real output vector */


vreal() extracts the real part of the complex elements of the vector a, placing the result in the real vector b.


 vrelerr(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vrelerr() calculates the relative error between the supposed exact values stored in real vector a and those contained in real vector b. The result is placed in the real data vector c according to the equation :

c[k] = (b[k] - a[k]) / a[k] for k = 0 to length(a)-1

Notice that vrelerr() performs (b - a)/a and NOT (a - b)/a.


 vresamp(a, b, c, d, e)
  UDI_object a; /* Complex input vector */
  UDI_object b; /* Real input scalar (multiplicative coefficient) */
  UDI_object c; /* Real input scalar (constant coefficient) */
  UDI_object d; /* Real output scalar (number of elements) */
  UDI_object e; /* Complex output vector */


vresamp() considers real vector a as a series of time, frequency pairs, representing an instantaneous frequency value at a specific time, and resamples (interpolates) these points, placing the result in the real output vector e. The rate at which this resampling takes place is determined by the real input scalars b and c. b is a multiplicative coefficient of the instantaneous period. This time increment is then added to constant-time increment c, determining the overall time increment from one resampled point to the next. Real output scalar d will contain the number of resampled points.


 vrev(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vrev() reverses the real vector a, placing the result in the real output vector b according to the equation :

b[k] = a[n - k - 1] for k = 0 to n-1 (where n = length(a) )

Notice that vrev() can reverse a vector in place if the input and output buffers are the same.


 vrot(a, b, n)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */
  short n; /* Shift amount */


vrot() rotates the input vector a to the left by n places, placing the result in output vector b.


 vsadd(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vsadd() adds the scalar value b to each element of the real vector a, placing the result in the real output vector c. A C-like expression of this scalar addition is :

c[k] = a[k] + b for k = 0 to length(a)-1


 vsasm(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* First real scalar input */
  UDI_object c; /* Second real scalar input */
  UDI_object d; /* Real output vector */


vsasm() adds to each element of the real data vector a the real scalar value b, and then multiplies the result by the real scalar value c and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (a[k] + b) * c for k = 0 to length(a)-1


 vsboadd(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vsboadd() computes a-b/a+b for each corresponding element of real input vectors a and b, placing the result in real output vector c. No checking is performed to insure that the denominator is non-zero!


 vsboads(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vsboads() computes (a-b)/(a+b) for each corresponding element of real input vectors a and b, placing the result in real output vector c. If the denominator is zero, then the result is considered as 1.0.


 vsbm(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Third real input vector */
  UDI_object d; /* Real output vector */


vsbm() subtracts each element of real vector a from the corresponding element of real vector b, then multiplies the result by the corresponding element of the real vector c and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (b[k] - a[k]) * c[k] for k = 0 to length(a)-1


 vsdiv(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vsdiv() divides each element of the real vector a by the real scalar b, placing the result in the real output vector c. A C-like expression of this process is :

c[k] = a[k] / b for k = 0 to length(a)-1


 vsin(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vsin() computes the sine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = sin(a[k]) for k = 0 to length(a)-1


 vsinh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vsinh() computes the hyperbolic-sine value of the elements of real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = sinh(a[k]) for k = 0 to length(a)-1


 vsinter(a, s, c, d)
  UDI_object a; /* Real input vector */
  UDI_object s; /* Real interpolation scalar */
  UDI_object c; /* Real input vector */
  UDI_object d; /* Real output vector */


vsinter() computes the real linear interpolation vector d of the two real vectors a and c with the real interpolation scalar s. An expression for the process is :

d[k] = (1.0 - s)*a[k] + s*c[k] for k = 0 to length(a)-1

The elements of s should be between 0 and 1. A value of 0 means that d = a. A value of 1 means that d = c


 vsma(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Real input scalar buffer */
  UDI_object c; /* Second real input vector */
  UDI_object d; /* Output vector */


vsma() multiplies each element of the real vector a by the real scalar b, then adds the corresponding element of the real vector c and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (a[k] * b) + c[k] for k = 0 to length(a)-1


 vsmax(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vsmax() compares each element of the real vector a with the real scalar value b, placing the maximum in the real output vector c. An expression for the process is :

c[k] = MAX(a[k], b) for k = 0 to length(a)-1


 vsmin(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vsmin() compares each element of the real vector a with the real scalar value b, placing the minimum in the real output vector c. An expression for the process is :

c[k] = MIN(a[k], b) for k = 0 to length(a)-1


 vsmmul(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Second real input vector */
  UDI_object d; /* Real output vector */


vsmmul() multiplies each element of the real vector a by the corresponding element of the real vector c and by the real scalar value b, placing the result in the real output vector d. A C-like expression of this process is :

d[k] = b * a[k] * c[k] for k = 0 to length(a)-1


 vsmsa(a, b, c, d)
  UDI_object a; /* Real input vector */
  UDI_object b; /* First real input scalar */
  UDI_object c; /* Second real input scalar */
  UDI_object d; /* Real output vector */


vsmsa() multiplies each element of the real vector a by the real scalar value b, then adds the real scalar value c and finally puts the result in the real output vector d. An expression for the process is : d[k] = (a[k] * b) + c for k = 0 to length(a)-1


 vsmsb(a, b, c, d)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Real scalar input */
  UDI_object c; /* Second real input vector */
  UDI_object d; /* Real output vector */


vsmsb() multiplies each element of the real vector a by the real scalar b, then subtracts the corresponding element of the real vector c and finally puts the result in the real output vector d. An expression for the process is :

d[k] = (a[k] * b) - c[k] for k = 0 to length(a)-1


 vsmul(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vsmul() multiplies each element of the real vector a by the real scalar value b, placing the result in the output vector c. A C-like expression of this process is :

c[k] = a[k] * b for k = 0 to length(a)-1


 vsq(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vsq() computes the square of each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = a[k] * a[k] for k = 0 to length(a)-1


 vsqicex(a, b, c, d)
  UDI_object a; /* Real input vector (numerator) */
  UDI_object b; /* Real input vector (part of denominator) */
  UDI_object c; /* Real input vector (part of denominator) */
  UDI_object d; /* Real output vector */


vsqicex() computes the square value of an intercorrelation-index like expression between the real vectors a, b and c, placing the result in the real vector d.

A C-like expression of the process is :

d[k] = (a[k] * a[k]) / (b[k] * c[k]) for k = 0 to length(a)-1


 vsqrt(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vsqrt() computes the square-root of each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = sqrt(a[k]) for k = 0 to length(a)-1


 vssq(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vssq() computes the signed-square of each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = a[k] * fabs(a[k]) for k = 0 to length(a)-1


 vssub(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Real output vector */


vssub() subtracts each element of the real vector a by the real scalar b, placing the result in the output vector c. A C-like expression of this process is :

c[k] = a[k] - b for k = 0 to length(a)-1


 vsub(a, b, c)
  UDI_object a; /* First real input vector */
  UDI_object b; /* Second real input vector */
  UDI_object c; /* Real output vector */


vsub() subtracts the real vector a from the real vector b, placing the result in the real vector c according to the equation :

c[k] = b[k] - a[k] for k = 0 to length(a)-1

Notice that vsub performs (b - a) and NOT (a - b).


 vswap(a, b)
  UDI_object a; /* First real vector */
  UDI_object b; /* Second real vector */


vsawp() exchanges each element of the real vector a with each corresponding element of the real vector b. An expression for the process is :

a[k] <-> b[k] for k = 0 to length(a)-1


 vtan(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vtan() computes the tangent of each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = tan(a[k]) for k = 0 to length(a)-1


 vtanh(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vtanh() computes the hyperbolic tangent of each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = tanh(a[k]) for k = 0 to length(a)-1


 vthres(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar (threshold) */
  UDI_object c; /* Real output vector */


vthres() thresholds each element of the real vector a by comparing it with the real scalar threshold value b. If the element is less than the threshold then the result is 0.0, otherwise the result is the unchanged element. The result is put in the real vector c. A C-like expression of the process is :

c[k] = (a[k] >= b ) ? a[k] : 0.0 for k = 0 to length(a)-1


 vwrap(a, b, c)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real input scalar */
  UDI_object c; /* Output vector */


vwrap() computes the value modulo b of each element of the real input vector a, placing result in the real output vector c. A C-like expression of the process is :

c[k] = a[k] - ((int)( a[k] / b )) * b for k = 0 to length(a)-1


 vy0(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vy0() computes the zero order Bessel function of the second kind for each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = y0(a[k]) for k = 0 to length(a)-1


 vy1(a, b)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */


vy1() computes the first order Bessel function of the second kind for each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = y1(a[k]) for k = 0 to length(a)-1


 vyn(a, b, p)
  UDI_object a; /* Real input vector */
  UDI_object b; /* Real output vector */
  long p; /* Bessel function order */


vyn() computes the p-th order Bessel function of the second kind for each element of the real vector a, placing the result in the real vector b. A C-like expression of the process is :

b[k] = yn(p, a[k]) for k = 0 to length(a)-1



next up previous contents index
Next: 5. DSP Function Summary Up: UDI 2.1A Unified DSP Previous: 3.2 Examples
Diemo Schwarz
1999-03-04