Chapter 12: Random Number Generation > random_orthogonal_matrix

random_orthogonal_matrix

Generates a pseudorandom orthogonal matrix or a correlation matrix.

Synopsis

#include <imsls.h>

float *imsls_f_random_orthogonal_matrix (int n, ..., 0)

The type double function is imsls_d_random_orthogonal_matrix.

Required Arguments

int n   (Input)
The order of the matrix to be generated.

Return Value

n by n random orthogonal matrix. To release this space, use free.

Synopsis with Optional Arguments

#include <imsls.h>

float *imsls_f_random_orthogonal_matrix (int n,
IMSLS_EIGENVALUES, float *eignevalues[],
IMSLS_A_MATRIX, float *a,
IMSLS_A_COL_DIM, int a_col_dim,
IMSLS_RETURN_USER, float r[],
 0)

Optional Arguments

IMSLS_EIGENVALUES, float *eigenvalues   (Input)
A vector of length n containing the eigenvalues of the correlation matrix to be generated.   The elements of eigenvalues must be positive, they must sum to n, and they cannot all be equal.

IMSLS_A_MATRIX, float *a   (Input)
n by n random orthogonal matrix.   A random correlation matrix is generated using the orthogonal matrix input in a.  The option IMSLS_EIGENVALUES must also be supplied if IMSLS_A_MATRIX is used.

IMSLS_A_COL_DIM, int a_col_dim   (Input)
Column dimension of the matrix a.
Default: a_col_dim = n 

IMSLS_RETURN_USER, float r[]   (Output)
User-supplied array of length n × n containing the random correlation matrix.

Description

Function imsls_f_random_orthogonal_matrix generates a pseudorandom orthogonal matrix from the invariant Haar measure. For each column, a random vector from a uniform distribution on a hypersphere is selected and then is projected onto the orthogonal complement of the columns already formed. The method is described by Heiberger (1978). (See also Tanner and Thisted 1982.)

If the optional argument IMSLS_EIGENVALUES is used, a correlation matrix is formed by applying a sequence of planar rotations to the matrix AT DA, where
D = diag(
eigenvalues[0], …, eigenvalues [n-1]), so as to yield ones along the diagonal. The planar rotations are applied in such an order that in the two by two matrix that determines the rotation, one diagonal element is less than 1.0 and one is greater than 1.0. This method is discussed by Bendel and Mickey (1978) and by Lin and Bendel (1985).

The distribution of the correlation matrices produced by this method is not known. Bendel and Mickey (1978) and Johnson and Welch (1980) discuss the distribution.

For larger matrices, rounding can become severe; and the double precision results may differ significantly from single precision results.

Example

In this example, imsls_f_random_orthogonal_matrix is used to generate a 4 by 4 pseudorandom correlation matrix with eigenvalues in the ratio 1:2:3:4.

 

#include <stdio.h>

#include <imsls.h>

 

int main()

{

  int   i, n = 4;

  float *a, *cor;

  float ev[] = {1., 2., 3., 4.};

  for (i=0;i<4;i++) ev[i] = 4.*ev[i]/10.;

 

  imsls_random_seed_set(123457);

  a = imsls_f_random_orthogonal_matrix(n, 0);

  imsls_f_write_matrix("Random orthogonal matrix",

                 4, 4, (float*)a, 0);

 

  cor = imsls_f_random_orthogonal_matrix(n,

                           IMSLS_EIGENVALUES, ev,

                           IMSLS_A_MATRIX, a,

                           0);

  imsls_f_write_matrix("Random correlation matrix",

                 4, 4, (float*)cor, 0);

 

}

 Output

 

            Random orthogonal matrix

            1           2           3           4

1     -0.8804     -0.2417      0.4065     -0.0351

2      0.3088     -0.3002      0.5520      0.7141

3     -0.3500      0.5256     -0.3874      0.6717

4     -0.0841     -0.7584     -0.6165      0.1941

 

            Random correlation matrix

            1           2           3           4

1       1.000      -0.236      -0.326      -0.110

2      -0.236       1.000       0.191      -0.017

3      -0.326       0.191       1.000      -0.435

4      -0.110      -0.017      -0.435       1.000

 


RW_logo.jpg
Contact Support