CHOL

Computes the Cholesky factorization of a positive-definite, symmetric or self-adjoint matrix.

Function Return Value

Matrix containing the Cholesky factorization of A.  The factor is upper triangular, RTR = A.  (Output)

Required Argument

A — Matrix to be factored. This argument must be a rank-2 or rank-3 array that contains a positive-definite, symmetric or self-adjoint matrix.  It may be real, double, complex, double complex. (Input)
For rank-3 arrays each rank-2 array, (for fixed third subscript), is a positive-definite, symmetric or self-adjoint matrix. In this case, the output is a rank-3 array of Cholesky factors for the individual problems. 

Optional Arguments, Packaged Options

This function uses lin_sol_self (See Chapter 1, “Linear Systems”), using the appropriate options to obtain the Cholesky factorization.

The option and derived type names are given in the following tables:

Option Names for CHOL

Option Value

Use_lin_sol_gen_only

4

Use_lin_sol_lsq_only

5

 

Name of Unallocated Option Array
to Use for Setting Options

Use

Derived Type

?_chol_options(:)

Use when setting options for calls hereafter.

?_options

?_chol_options_once(:)

Use when setting options for next call only.

?_options

For a description on how to use these options, see Matrix Optional Data Changes.  See lin_sol_self located in Chapter 1, “Linear Systems”  for the specific options for this routine.

FORTRAN 90 Interface

CHOL(A)

Description

Computes the Cholesky factorization of a positive-definite, symmetric or self-adjoint matrix, A. The factor is upper triangular, RTR = A.

Examples

Dense Matrix Example  (operator_ex06.f90)

 

      use linear_operators

 

      implicit none

 

! This is the equivalent of Example 2 for LIN_SOL_SELF using operators
! and functions.

 

      integer, parameter :: m=64, n=32

      real(kind(1e0)) :: one=1e0, zero=0e0, err

      real(kind(1e0)) A(n,n), b(n), C(m,n), d(m), cov(n,n), x(n)

           

! Generate a random rectangular matrix and right-hand side.

      C = rand(C); d=rand(d)

 

! Form the normal equations for the rectangular system.

      A = C .tx. C; b = C .tx. d

      COV = .i. CHOL(A); COV = COV .xt. COV

 

! Compute the least-squares solution.

       x = C .ix. d

 

! Compare with solution obtained using the inverse matrix.

      err = norm(x - (COV .x. b))/norm(cov)

 

! Scale the inverse to obtain the sample covariance matrix.

      COV = sum((d - (C .x. x))**2)/(m-n) * COV

! Check the results.

      if (err <= sqrt(epsilon(one))) then

         write (*,*) 'Example 2 for LIN_SOL_SELF (operators) is correct.'

      end if

 

      end 

Parallel Example  (parallel_ex06.f90)

 

      use linear_operators

      use mpi_setup_int

 

      implicit none

 

! This is the equivalent of Parallel Example 6 for box data types, operators ! and functions.

 

      integer, parameter :: m=64, n=32, nr=4

      real(kind(1e0)) :: one=1e0, zero=0e0, err(nr)

      real(kind(1e0)), dimension(m,n,nr) :: C, d(m,1,nr)

      real(kind(1e0)), dimension(n,n,nr) :: A, cov

      real(kind(1e0)), dimension(n,1,nr) :: b, x

 

! Setup for MPI:

      mp_nprocs=mp_setup()

          

! Generate a random rectangular matrix and right-hand side.

      if(mp_rank == 0) then

         C = rand(C); d=rand(d)

      endif

 

! Form the normal equations for the rectangular system.

      A = C .tx. C; b = C .tx. d

      COV = .i. CHOL(A); COV = COV .xt. COV

 

! Compute the least-squares solution.

       x = C .ix. d

 

! Compare with solution obtained using the inverse matrix.

      err = norm(x - (COV .x. b))/norm(cov)

 

! Check the results.

      if (ALL(err <= sqrt(epsilon(one))) .and. mp_rank == 0) &

         write (*,*) 'Parallel Example 6 is correct.'

     

! See to any eror messages and quit MPI

      mp_nprocs=mp_setup('Final')

 

      end


Visual Numerics, Inc.
Visual Numerics - Developers of IMSL and PV-WAVE
http://www.vni.com/
PHONE: 713.784.3131
FAX:713.781.9260