LFSDH

more...

more...

Solves a complex Hermitian positive definite system of linear equations given the RH R factorization of the coefficient matrix.

Required Arguments

FACT — Complex N by N matrix containing the factorization of the coefficient matrix A as output from routine LFCDH/DLFCDH or LFTDH/DLFTDH. (Input)

B — Complex vector of length N containing the right-hand side of the linear system. (Input)

X — Complex vector of length N containing the solution to the linear system. (Output)
If B is not needed, B and X can share the same storage locations.

Optional Arguments

N — Number of equations. (Input)
Default: N = size (FACT,2).

LDFACT — Leading dimension of FACT exactly as specified in the dimension statement of the calling program. (Input)
Default: LDFACT = size (FACT,1).

FORTRAN 90 Interface

Generic: CALL LFSDH (FACT, B, X [, …])

Specific: The specific interface names are S_LFSDH and D_LFSDH.

FORTRAN 77 Interface

Single: CALL LFSDH (N, FACT, LDFACT, B, X)

Double: The double precision name is DLFSDH.

ScaLAPACK Interface

Generic: CALL LFSDH (FACT0, B0, X0 [, …])

Specific: The specific interface names are S_LFSDH and D_LFSDH.

See the ScaLAPACK Usage Notes below for a description of the arguments for distributed computing.

Description

Routine LFSDH computes the solution for a system of linear algebraic equations having a complex Hermitian positive definite coefficient matrix. To compute the solution, the coefficient matrix must first undergo an RHR factorization. This may be done by calling either LFCDH or LFTDH. R is an upper triangular matrix.

The solution to Ax = b is found by solving the triangular systems RH y = b and Rx = y.

LFSDH and LFIDH both solve a linear system given its RHR factorization. LFIDH generally takes more time and produces a more accurate answer than LFSDH. Each iteration of the iterative refinement algorithm used by LFIDH calls LFSDH.

The underlying code is based on either LINPACK , LAPACK, or ScaLAPACK code depending upon which supporting libraries are used during linking. For a detailed explanation see Using ScaLAPACK, LAPACK, LINPACK, and EISPACK in the Introduction section of this manual.

Comments

Informational error

 

Type

Code

Description

4

1

The input matrix is singular.

ScaLAPACK Usage Notes

The arguments which differ from the standard version of this routine are:

FACT0MXLDA by MXCOL complex local matrix containing the local portions of the distributed matrix FACT as output from routine LFCDH/DLFCDH or LFTDH/DLFTDH. FACT contains the factorization of the matrix A. (Input)

B0 — Complex local vector of length MXLDA containing the local portions of the distributed vector B. B contains the right-hand side of the linear system. (Input)

X0 — Complex local vector of length MXLDA containing the local portions of the distributed vector X. X contains the solution to the linear system. (Output)
If B is not needed, B and X can share the same storage locations.

All other arguments are global and are the same as described for the standard version of the routine. In the argument descriptions above, MXLDA and MXCOL can be obtained through a call to SCALAPACK_GETDIM (see Utilities) after a call to SCALAPACK_SETUP (Utilities) has been made. See the ScaLAPACK Example below.

Examples

Example

A set of linear systems is solved successively. LFTDH is called to factor the coefficient matrix. LFSDH is called to compute the four solutions for the four right-hand sides. In this case, the coefficient matrix is assumed to be well-conditioned and correctly scaled. Otherwise, it would be better to call LFCDH to perform the factorization, and LFIDH to compute the solutions.

 

USE LFSDH_INT

USE LFTDH_INT

USE WRCRN_INT

! Declare variables

INTEGER LDA, LDFACT, N

PARAMETER (LDA=5, LDFACT=5, N=5)

COMPLEX A(LDA,LDA), B(N,3), FACT(LDFACT,LDFACT), X(N,3)

 

! Set values for A and B

!

! A = ( 2.0+0.0i -1.0+1.0i 0.0+0.0i 0.0+0.0i 0.0+0.0i )

! ( 4.0+0.0i 1.0+2.0i 0.0+0.0i 0.0+0.0i )

! ( 10.0+0.0i 0.0+4.0i 0.0+0.0i )

! ( 6.0+0.0i 1.0+1.0i )

! ( 9.0+0.0i )

!

! B = ( 3.0+3.0i 4.0+0.0i 29.0-9.0i )

! ( 5.0-5.0i 15.0-10.0i -36.0-17.0i )

! ( 5.0+4.0i -12.0-56.0i -15.0-24.0i )

! ( 9.0+7.0i -12.0+10.0i -23.0-15.0i )

! (-22.0+1.0i 3.0-1.0i -23.0-28.0i )

 

DATA A /(2.0,0.0), 4*(0.0,0.0), (-1.0,1.0), (4.0,0.0),&

4*(0.0,0.0), (1.0,2.0), (10.0,0.0), 4*(0.0,0.0),&

(0.0,4.0), (6.0,0.0), 4*(0.0,0.0), (1.0,1.0), (9.0,0.0)/

DATA B /(3.0,3.0), (5.0,-5.0), (5.0,4.0), (9.0,7.0), (-22.0,1.0),&

(4.0,0.0), (15.0,-10.0), (-12.0,-56.0), (-12.0,10.0),&

(3.0,-1.0), (29.0,-9.0), (-36.0,-17.0), (-15.0,-24.0),&

(-23.0,-15.0), (-23.0,-28.0)/

 

! Factor the matrix A

CALL LFTDH (A, FACT)

! Compute the solutions

DO 10 I=1, 3

CALL LFSDH (FACT, B(:,I), X(:,I))

10 CONTINUE

! Print solutions

CALL WRCRN (’X’, X)

!

END

Output

 

X

1 2 3

1 ( 1.00, 0.00) ( 3.00, -1.00) ( 11.00, -1.00)

2 ( 1.00, -2.00) ( 2.00, 0.00) ( -7.00, 0.00)

3 ( 2.00, 0.00) ( -1.00, -6.00) ( -2.00, -3.00)

4 ( 2.00, 3.00) ( 2.00, 1.00) ( -2.00, -3.00)

5 ( -3.00, 0.00) ( 0.00, 0.00) ( -2.00, -3.00)

ScaLAPACK Example

The same set of linear systems as in in the preceding example is solved successively as a distributed computing example. LFTDH is called to factor the matrix. LFSDH is called to compute the four solutions for the four right-hand sides. In this case, the coefficient matrix is assumed to be well-conditioned and correctly scaled. Otherwise, it would be better to call LFCDH to perform the factorization, and LFIDH to compute the solutions.

SCALAPACK_MAP and SCALAPACK_UNMAP are IMSL utility routines (see Utilities) used to map and unmap arrays to and from the processor grid. They are used here for brevity. DESCINIT is a ScaLAPACK tools routine which initializes the descriptors for the local arrays.

 

USE MPI_SETUP_INT

USE LFTDH_INT

USE LFSDH_INT

USE WRCRN_INT

USE SCALAPACK_SUPPORT

IMPLICIT NONE

INCLUDE ‘mpif.h’

! Declare variables

INTEGER J, LDA, N, DESCA(9), DESCX(9)

INTEGER INFO, MXCOL, MXLDA

COMPLEX, ALLOCATABLE :: A(:,:), B(:,:), B0(:), X(:,:)

COMPLEX, ALLOCATABLE :: A0(:,:), FACT0(:,:), X0(:)

PARAMETER (LDA=5, N=5)

! Set up for MPI

MP_NPROCS = MP_SETUP()

IF(MP_RANK .EQ. 0) THEN

ALLOCATE (A(LDA,N), B(LDA,3), X(LDA,3))

! Set values for A and B

A(1,:) = (/(2.0, 0.0),(-1.0, 1.0),( 0.0, 0.0),(0.0, 0.0),(0.0, 0.0)/)

A(2,:) = (/(0.0, 0.0),( 4.0, 0.0),( 1.0, 2.0),(0.0, 0.0),(0.0, 0.0)/)

A(3,:) = (/(0.0, 0.0),( 0.0, 0.0),(10.0, 0.0),(0.0, 4.0),(0.0, 0.0)/)

A(4,:) = (/(0.0, 0.0),( 0.0, 0.0),( 0.0, 0.0),(6.0, 0.0),(1.0, 1.0)/)

A(5,:) = (/(0.0, 0.0),( 0.0, 0.0),( 0.0, 0.0),(0.0, 0.0),(9.0, 0.0)/)

!

B(1,:) = (/(3.0, 3.0), ( 4.0, 0.0), ( 29.0, -9.0)/)

B(2,:) = (/(5.0, -5.0), ( 15.0,-10.0), (-36.0,-17.0)/)

B(3,:) = (/(5.0, 4.0), (-12.0,-56.0), (-15.0,-24.0)/)

B(4,:) = (/(9.0, 7.0), (-12.0, 10.0), (-23.0,-15.0)/)

B(5,:) = (/(-22.0,1.0), ( 3.0, -1.0), (-23.0,-28.0)/)

ENDIF

! Set up a 1D processor grid and define

! its context ID, MP_ICTXT

CALL SCALAPACK_SETUP(N, N, .TRUE., .TRUE.)

! Get the array descriptor entities MXLDA,

! and MXCOL

CALL SCALAPACK_GETDIM(N, N, MP_MB, MP_NB, MXLDA, MXCOL)

! Set up the array descriptors

CALL DESCINIT(DESCA, N, N, MP_MB, MP_NB, 0, 0, MP_ICTXT, MXLDA, INFO)

CALL DESCINIT(DESCX, N, 1, MP_MB, 1, 0, 0, MP_ICTXT, MXLDA, INFO)

! Allocate space for the local arrays

ALLOCATE(A0(MXLDA,MXCOL), X0(MXLDA),FACT0(MXLDA,MXCOL), &

B0(MXLDA))

! Map input arrays to the processor grid

CALL SCALAPACK_MAP(A, DESCA, A0)

! Factor the matrix A

CALL LFTDH (A0, FACT0)

! Compute the solutions

DO 10 J=1, 3

CALL SCALAPACK_MAP(B(:,J), DESCX, B0)

CALL LFSDH (FACT0, B0, X0)

! Unmap the results from the distributed

! array back to a non-distributed array

CALL SCALAPACK_UNMAP(X0, DESCX, X(:,J))

10 CONTINUE

! Print the results.

! After the unmap, only Rank=0 has the full

! array.

IF(MP_RANK .EQ. 0) CALL WRCRN (’X’, X)

IF (MP_RANK .EQ. 0) DEALLOCATE(A, B, X)

DEALLOCATE(A0, B0, FACT0, X0)

! Exit ScaLAPACK usage

CALL SCALAPACK_EXIT(MP_ICTXT)

! Shut down MPI

MP_NPROCS = MP_SETUP(‘FINAL’)

END

Output

 

X

1 2 3

1 ( 1.00, 0.00) ( 3.00, -1.00) ( 11.00, -1.00)

2 ( 1.00, -2.00) ( 2.00, 0.00) ( -7.00, 0.00)

3 ( 2.00, 0.00) ( -1.00, -6.00) ( -2.00, -3.00)

4 ( 2.00, 3.00) ( 2.00, 1.00) ( -2.00, -3.00)

5 ( -3.00, 0.00) ( 0.00, 0.00) ( -2.00, -3.00)