DET

 


   more...

Computes the determinant of a rectangular matrix.

Function Return Value

Determinant of matrix A. This is a scalar for the case where A is rank 2. It is a rank-1 array of determinant values for the case where A is rank 3. (Output)

Required Argument

A — Matrix for which the determinant is to be computed. This argument must be a rank-2 or rank-3 array that contains a rectangular matrix. It may be real, double, complex, double complex. (Input)

For rank-3 arrays, each rank-2 array (for fixed third subscript) is a separate matrix. In this case, the output is a rank-1 array of determinant values for each problem.

Optional Arguments, Packaged Options

This function uses LIN_SOL_LSQ (see Chapter 1, “Linear Systems”) to compute the QR decomposition of A, and the logarithmic value of det(A), which is exponentiated for the result.

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

Option Names for DET

Option Value

?_det_for_lin_sol_lsq

1

 

Name of Unallocated Option Array to Use for Setting Options

Use

Derived Type

?_det_options(:)

Use when setting options for calls hereafter.

?_options

?_det_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_LSQ in Chapter 1, “Linear Systems” for the specific options for these routines.

FORTRAN 90 Interface

DET (A)

Description

Computes the determinant of a rectangular matrix, A. The evaluation is based on the QR decomposition,

 

and k = rank(A). Thus det(A) = s × det(R) where s = det(Q× det(P) = ±1.

Even well-conditioned matrices can have determinants with values that have very large or very tiny magnitudes. The values may overflow or underflow. For this class of problems, the use of the logarithmic representation of the determinant found in LIN_SOL_GEN or LIN_SOL_LSQ is required.

Examples

Dense Matrix Example (operator_ex02.f90)

 

use linear_operators

implicit none

 

! This is Example 2 for LIN_SOL_GEN using operators and functions.

 

integer, parameter :: n=32

real(kind(1e0)) :: one=1e0, err, det_A, det_i

real(kind(1e0)), dimension(n,n) :: A, inv

 

! Generate a random matrix.

A = rand(A)

! Compute the matrix inverse and its determinant.

inv = .i.A; det_A = det(A)

! Compute the determinant for the inverse matrix.

det_i = det(inv)

! Check the quality of both left and right inverses.

err = (norm(EYE(n)-(A .x. inv))+norm(EYE(n)-(inv.x.A)))/cond(A)

if (err <= sqrt(epsilon(one)) .and. abs(det_A*det_i - one) <= &

sqrt(epsilon(one))) &

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

end

Parallel Example (parallel_ex02.f90)

 

use linear_operators

use mpi_setup_int

 

implicit none

 

! This is the equivalent of Parallel Example 2 for .i. and det() with box

! data types, operators and functions.

 

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

integer J

real(kind(1e0)) :: one=1e0

real(kind(1e0)), dimension(nr) :: err, det_A, det_i

real(kind(1e0)), dimension(n,n,nr) :: A, inv, R, S

 

! Setup for MPI.

MP_NPROCS=MP_SETUP()

! Generate a random matrix.

A = rand(A)

! Compute the matrix inverse and its determinant.

inv = .i.A; det_A = det(A)

! Compute the determinant for the inverse matrix.

det_i = det(inv)

! Check the quality of both left and right inverses.

DO J=1,nr; R(:,:,J)=EYE(N); END DO

 

S=R; R=R-(A .x. inv); S=S-(inv .x. A)

err = (norm(R)+norm(S))/cond(A)

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

abs(det_A*det_i - one) <= sqrt(epsilon(one)))&

.and. MP_RANK == 0) &

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

 

! See to any error messages and quit MPI.

MP_NPROCS=MP_SETUP('Final')

 

end