IBM Books

Parallel Engineering and Scientific Subroutine Library for AIX Version 2 Release 3: Guide and Reference

PDSYGVX and PZHEGVX--Selected Eigenvalues and, Optionally, the Eigenvectors of a Real Symmetric or Complex Hermitian Positive Definite Generalized Eigenproblem

|These subroutines compute selected eigenvalues and, optionally, the |eigenvectors of a real symmetric or complex Hermitian positive definite |generalized eigenproblem: |

|In the formulas above: |

|A represents the global real symmetric or complex Hermitian |submatrix Aia:ia+n-1, |ja:ja+n-1
|B represents the global real symmetric or complex Hermitian |positive definite submatrix |Bib:ib+n-1, |jb:jb+n-1 |

Eigenvalues and eigenvectors can be selected by specifying a range of values or a range of indices for the desired eigenvalues.

If n = 0, no computation is performed and the subroutine returns after doing some parameter checking.

See references [13], [24], [25], and [26].

Table 103. Data Types

vl, vu, abstol, w, orfac, rwork, gap A, B, Z, work iwork, ifail, iclustr Subroutine
Long-precision real Long-precision real Integer PDSYGVX
Long-precision real Long-precision complex Integer PZHEGVX

Syntax

Fortran CALL PDSYGVX (ibtype, jobz, range, uplo, n, a, ia, ja, desc_a, b, ib, jb, desc_b, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, desc_z, work, lwork, iwork, liwork, ifail, iclustr, gap, info)
CALL PZHEGVX (ibtype, jobz, range, uplo, n, a, ia, ja, desc_a, b, ib, jb, desc_b, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, desc_z, work, lwork, rwork, lrwork, iwork, liwork, ifail, iclustr, gap, info)
C and C++ pdsygvx (ibtype, jobz, range, uplo, n, a, ia, ja, desc_a, b, ib, jb, desc_b, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, desc_z, work, lwork, iwork, liwork, ifail, iclustr, gap, info);
pzhegvx (ibtype, jobz, range, uplo, n, a, ia, ja, desc_a, b, ib, jb, desc_b, vl, vu, il, iu, abstol, m, nz, w, orfac, z, iz, jz, desc_z, work, lwork, rwork, lrwork, iwork, liwork, ifail, iclustr, gap, info);

On Entry

ibtype
specifies the problem type, where:

If ibtype = 1, the problem is Ax = lambdaBx

If ibtype = 2, the problem is ABx = lambdax

If ibtype = 3, the problem is BAx = lambdax

Scope: global

Specified as: a fullword integer; ibtype = 1, 2, or 3.

jobz
indicates the type of computation to be performed, where:

If jobz = 'N', eigenvalues only are computed.

If jobz = 'V', eigenvalues and eigenvectors are computed.

Scope: global

Specified as: a single character; jobz = 'N' or 'V'.

range
indicates which eigenvalues to compute, where:

If range = 'A', all eigenvalues are to be found.

If range = 'V', all eigenvalues in the interval [vl, vu] are to be found.

If range = 'I', the il-th through iu-th eigenvalues are to be found.

Scope: global

Specified as: a single character; range = 'A', 'V', or 'I'.

uplo
indicates whether the upper or lower triangular part of the global submatrices A and B are referenced, where:

If uplo = 'U', the upper triangular part is referenced.

If uplo = 'L', the lower triangular part is referenced.

Scope: global

Specified as: a single character; uplo = 'U' or 'L'.

n
is the order of submatrices A and B used in the computation.

Scope: global

Specified as: a fullword integer; n >= 0.

a
is the local part of the global real symmetric |or complex Hermitian matrix A. This identifies the first element of the local array A. This subroutine computes the location of the first element of the local subarray used, based on ia, ja, desc_a, p, q, myrow, and mycol; therefore, the leading LOCp(ia+n-1) by LOCq(ja+n-1) part of the local array A must contain the local pieces of the leading ia+n-1 by ja+n-1 part of the global matrix, and:

Scope: local

Specified as: an LLD_A by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 103. Details about the square block-cyclic data distribution of global matrix A are stored in desc_a.

ia
is the row index of the global matrix A, identifying the first row of the submatrix A.

Scope: global

Specified as: a fullword integer; 1 <= ia <= M_A and ia+n-1 <= M_A.

ja
is the column index of the global matrix A, identifying the first column of the submatrix A.

Scope: global

Specified as: a fullword integer; 1 <= ja <= N_A and ja+n-1 <= N_A.

desc_a
is the array descriptor for global matrix A, described in the following table:
desc_a Name Description Limits Scope
1 DTYPE_A Descriptor type DTYPE_A=1 Global
2 CTXT_A BLACS context Valid value, as returned by BLACS_GRIDINIT or BLACS_GRIDMAP Global
3 M_A Number of rows in the global matrix If n = 0: M_A >= 0

Otherwise: M_A >= 1

Global
4 N_A Number of columns in the global matrix If n = 0: N_A >= 0

Otherwise: N_A >= 1

Global
5 MB_A Row block size MB_A >= 1 Global
6 NB_A Column block size NB_A >= 1 Global
7 RSRC_A The process row of the p × q grid over which the first row of the global matrix is distributed 0 <= RSRC_A < p Global
8 CSRC_A The process column of the p × q grid over which the first column of the global matrix is distributed 0 <= CSRC_A < q Global
9 LLD_A The leading dimension of the local array LLD_A >= max(1,LOCp(M_A)) Local

Specified as: an array of (at least) length 9, containing fullword integers.

b
is the local part of the global real symmetric |or complex Hermitian positive definite matrix B. This identifies the first element of the local array B. This subroutine computes the location of the first element of the local subarray used, based on ib, jb, desc_b, p, q, myrow, and mycol; therefore, the leading LOCp(ib+n-1) by LOCq(jb+n-1) part of the local array B must contain the local pieces of the leading ib+n-1 by jb+n-1 part of the global matrix, and:

Scope: local

Specified as: an LLD_B by (at least) LOCq(N_B) array, containing numbers of the data type indicated in Table 105. Details about the square block-cyclic data distribution of global matrix B are stored in desc_b.

ib
is the row index of the global matrix B, identifying the first row of the submatrix B.

Scope: global

Specified as: a fullword integer; 1 <= ib <= M_B and ib+n-1 <= M_B.

jb
is the column index of the global matrix B, identifying the first column of the submatrix B.

Scope: global

Specified as: a fullword integer; 1 <= jb <= N_B and jb+n-1 <= N_B.

desc_b
is the array descriptor for global matrix B, described in the following table:
desc_b Name Description Limits Scope
1 DTYPE_B Descriptor type DTYPE_B=1 Global
2 CTXT_B BLACS context Valid value, as returned by BLACS_GRIDINIT or BLACS_GRIDMAP Global
3 M_B Number of rows in the global matrix If n = 0: M_B >= 0

Otherwise: M_B >= 1

Global
4 N_B Number of columns in the global matrix If n = 0: N_B >= 0

Otherwise: N_B >= 1

Global
5 MB_B Row block size MB_B >= 1 Global
6 NB_B Column block size NB_B >= 1 Global
7 RSRC_B The process row of the p × q grid over which the first row of the global matrix is distributed 0 <= RSRC_B < p Global
8 CSRC_B The process column of the p × q grid over which the first column of the global matrix is distributed 0 <= CSRC_B < q Global
9 LLD_B The leading dimension of the local array LLD_B >= max(1,LOCp(M_B)) Local

Specified as: an array of (at least) length 9, containing fullword integers.

vl
has the following meaning:

If range = 'V', it is the lower bound of the interval to be searched for eigenvalues.

If range <> 'V', this argument is ignored.

Scope: global

Specified as: a number of the data type indicated in Table 103. If range = 'V', vl < vu.

vu
has the following meaning:

If range = 'V', it is the upper bound of the interval to be searched for eigenvalues.

If range <> 'V', this argument is ignored.

Scope: global

Specified as: a number of the data type indicated in Table 103. If range = 'V', vl < vu.

il
has the following meaning:

If range = 'I', it is the index (from smallest to largest) of the smallest eigenvalue to be returned.

If range <> 'I', this argument is ignored.

Scope: global

Specified as: a fullword integer; il >= 1.

iu
has the following meaning:

If range = 'I', it is the index (from smallest to largest) of the largest eigenvalue to be returned.

If range <> 'I', this argument is ignored.

Scope: global

Specified as: a fullword integer; min(il, n) <= iu <= n.

abstol
is the absolute tolerance for the eigenvalues. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a, b] of width less than or equal to:
abstol+epsilon(max(|a|, |b|))

where epsilon is the machine precision. If abstol is less than or equal to zero, then epsilon(norm(T)) is used in its place, where norm(T) is the 1-norm of the tridiagonal matrix obtained by reducing the standard form of the generalized problem to tridiagonal form. For most problems, this is the appropriate level of accuracy to request.

For certain strongly graded matrices, greater accuracy can be obtained in very small eigenvalues by setting abstol to a very small positive number. However, if abstol is less than:



Math Graphic

where unfl is the underflow threshold, then:



Math Graphic

is used in its place.

Eigenvalues are computed most accurately when abstol is set to twice the underflow threshold--that is, (2)(unfl).

If jobz = 'V', then setting abstol to unfl, the underflow threshold, yields the most orthogonal eigenvectors.



Math Graphic

Scope: global

Specified as: a number of the data type indicated in Table 103.

m
See On Return.

nz
See On Return.

w
See On Return.

orfac
specifies which eigenvectors should be reorthogonalized. Eigenvectors that correspond to eigenvalues which are within:
ortol = (orfac)(norm(R))

of each other (where norm(R) is the 1-norm of the standard form of the generalized eigenproblem) are to be reorthogonalized.

However, if the workspace is insufficient (see lwork and |lrwork), ortol may be decreased until all eigenvectors to be reorthogonalized can be stored in one process.

If orfac is zero, no reorthogonalization is done.

If orfac is less than zero, a default value of 10-3 is used.

Scope: global

Specified as: a number of the data type indicated in Table 103.

z
See On Return.

iz
is the row index of the global matrix Z, identifying the first row of the submatrix Z.

Scope: global

Specified as: a fullword integer; 1 <= iz <= M_Z and iz+n-1 <= M_Z.

jz
is the column index of the global matrix Z, identifying the first column of the submatrix Z.

Scope: global

Specified as: a fullword integer; 1 <= jz <= N_Z and jz+n-1 <= N_Z.

desc_z
is the array descriptor for global matrix Z, described in the following table:
desc_z Name Description Limits Scope
1 DTYPE_Z Descriptor type DTYPE_Z=1 Global
2 CTXT_Z BLACS context Valid value, as returned by BLACS_GRIDINIT or BLACS_GRIDMAP Global
3 M_Z Number of rows in the global matrix If n = 0: M_Z >= 0

Otherwise: M_Z >= 1

Global
4 N_Z Number of columns in the global matrix If n = 0: N_Z >= 0

Otherwise: N_Z >= 1

Global
5 MB_Z Row block size MB_Z >= 1 Global
6 NB_Z Column block size NB_Z >= 1 Global
7 RSRC_Z The process row of the p × q grid over which the first row of the global matrix is distributed 0 <= RSRC_Z < p Global
8 CSRC_Z The process column of the p × q grid over which the first column of the global matrix is distributed 0 <= CSRC_Z < q Global
9 LLD_Z The leading dimension of the local array LLD_Z >= max(1,LOCp(M_Z)) Local

Specified as: an array of (at least) length 9, containing fullword integers. |

|work
|has the following meaning:

|If lwork = 0, work is ignored.

|If lwork <> 0, work is a work area used by this |subroutine, where: |

|Scope: local

|Specified as: an area of storage containing numbers of data type |indicated in |Table 103.

| lwork
|is the number of elements in array WORK.

|Scope: |

|Specified as: a fullword integer; where: |

|

|rwork
|has the following meaning:

|If lrwork = 0, rwork is ignored.

|If lrwork <> 0, rwork is a work area used by |this subroutine, where: |

|Scope: local

|Specified as: an area of storage containing numbers of data type |indicated in Table 103.

| lrwork
|is the number of elements in array RWORK.

|Scope: |

|Specified as: a fullword integer; where: |

iwork
has the following meaning:

If liwork = 0, iwork is ignored.

If liwork <> 0, iwork is a work area used by this subroutine, where:

Scope: local

Specified as: an area of storage containing numbers of data type indicated in Table 103.

liwork
is the number of elements in array IWORK.

Scope:

Specified as: a fullword integer; where:

ifail
See On Return.

iclustr
See On Return.

gap
See On Return.

info
See On Return.

On Return

a
a is the local part of the |global matrix A, where:

If uplo = 'U', the upper triangle and diagonal of submatrix Aia:ia+n-1, ja:ja+n-1 are overwritten; that is, the original input is not preserved.

If uplo = 'L', the lower triangle and diagonal of submatrix Aia:ia+n-1, ja:ja+n-1 are overwritten; that is, the original input is not preserved.

Scope: local

Returned as: an LLD_A by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 103.

b
is the local part of the global |real symmetric |or complex Hermitian matrix B, containing the results of the Cholesky |factorization.

Scope: local

Specified as: an LLD_B by (at least) LOCq(N_B) array, containing numbers of the data type indicated in Table 103. Details about the square block-cyclic data distribution of global matrix B are stored in desc_b.

m
is the number of eigenvalues found.

Scope: global

Returned as: a fullword integer; 0 <= m <= n.

nz
has the following meaning:

If jobz <> 'V', then nz is ignored.

If jobz = 'V', then nz is the number of eigenvectors computed--that is, the number of columns of Z used in the computation. On output, nz = m unless you provide insufficient space. To get all the eigenvectors requested, you must supply both sufficient space to hold the eigenvectors in Z and sufficient workspace to compute them (see lwork and |lrwork).

If range = 'A' or 'I', this subroutine does not perform any computations if the work space supplied is insufficient. In this case, an input-argument error is issued and your job is terminated. For range = 'V', the number of requested eigenvectors is unknown until the eigenvalues are found. In this case, the subroutine computes as many eigenvectors as space allows. Then, if nz <> m, a computational error message is issued.

Scope: global

Returned as: a fullword integer; 0 <= nz <= m.

w
On normal exit (see info), it is the vector w, containing the selected eigenvalues in ascending order in the first m elements of w.

Scope: global

Returned as: a one-dimensional array of (at least) length n, containing numbers of the data type indicated in Table 103.

z
has the following meaning:

If jobz = 'N', then z is ignored.

If jobz = 'V' and there is a normal exit (see info), then this is the updated local part of the global matrix Z, where columns jz to jz+m-1 of the global matrix Z contain the orthonormal eigenvectors, corresponding to the selected eigenvalues. If an eigenvector fails to converge, then the corresponding column of the global matrix Z contains the last approximation to the eigenvector, and the index of the eigenvector is returned in ifail.

This identifies the first element of the local array Z. This subroutine computes the location of the first element of the local subarray used, based on iz, jz, desc_z, p, q, myrow, and mycol; therefore, the leading LOCp(iz+n-1) by LOCq(jz+n-1) part of the local array Z must contain the local pieces of the leading iz+n-1 by jz+n-1 part of the global matrix Z.

Scope: local

Returned as: an LLD_Z by (at least) LOCq(N_Z) array, containing numbers of the data type indicated in Table 103.

| work
|is the work area used by this subroutine if lwork  <>  |0,where:

|If lwork  <>  0 and lwork  <>  -1, |its size is (at least) of length lwork.

|If lwork  =  -1, its size is (at least) of length |1.

|Scope: local

|Returned as: an area of storage, containing numbers of the data type |indicated in Table 103, where:

|If lwork  >=  1 or lwork  =  -1, |then: |

| rwork
|is the work area used by this subroutine if lrwork  <>  |0,where:

|If lrwork  <>  0 and lrwork  <>  -1, |its size is (at least) of length lrwork.

|If lrwork  =  -1, its size is (at least) of length |1.

|Scope: local

|Returned as: an area of storage, containing numbers of the data type |indicated in Table 103, where:

|If lrwork  >=  1 or lrwork  =  -1, |then: |

iwork
is the work area used by this subroutine if liwork  <>  0, where:

If liwork  <>  0 and liwork  <>  -1, then its size is (at least) of length liwork.

If liwork  =  -1, then its size is (at least) of length 1.

Scope: local

Returned as: an area of storage, where:

If liwork  >=  1 or liwork  =  -1, then iwork1 is set to the minimum liwork value and contains numbers of the data type indicated in Table 103. Except for iwork1, the contents of iwork are overwritten on return.

ifail
has the following meaning:

If B is not positive definite (see info), then ifail1 is set to the order of the smallest minor which is not positive definite.

If B is positive definite and if jobz = 'V', it is vector ifail, where:

Scope: global

Returned as: a one-dimensional array of (at least) length n, containing fullword integers; 0 <= ifaili <= n.

iclustr
has the following meaning:

If jobz = 'N', then iclustr is ignored.

If jobz = 'V', it is vector iclustr, containing the indices of the eigenvectors corresponding to a cluster of eigenvalues that could not be reorthogonalized due to insufficient workspace. Eigenvectors corresponding to clusters of eigenvalues indexed iclustr2i-1 to iclustr2i could not be reorthogonalized due to lack of workspace. Hence, the eigenvectors corresponding to these clusters may not be orthogonal.

iclustr is a zero-terminated vector; that is, the last element of iclustr is set to zero. Assuming that k is the number of clusters, then:

iclustr2k <> 0 and iclustr2k+1 = 0

Scope: global

Returned as: a one-dimensional array of (at least) length 2(nprow)(npcol), containing fullword integers; 0 <= iclustri <= n.

gap
has the following meaning:

If jobz = 'N', then gap is ignored.

If jobz = 'V', it is vector gap, containing the gap between the eigenvalues whose eigenvectors could not be reorthogonalized. The values in this vector correspond to the clusters indicated by iclustr. As a result, the dot product between the eigenvectors corresponding to the i-th cluster may be as high as (C)(n)/gapi, where C is a small constant.

Scope: global

Returned as: a one-dimensional array of (at least) length (nprow)(npcol), containing numbers of the data type indicated in Table 103.

info
has the following meaning:

If info = 0, then no input-argument errors or computational errors occurred. This indicates a normal exit.

Note:
One use of info in ScaLAPACK is to identify whether input-argument errors occurred. Because Parallel ESSL terminates the application if input-argument errors occur, the setting of info is irrelevant for these errors.

If info > 0, then one or more of the following computational errors occurred and the appropriate error messages were issued, indicating an error exit, where:

Scope: global

Returned as: a fullword integer; info >= 0.

Notes and Coding Rules
  1. This subroutine accepts lowercase letters for the jobz, range, and uplo arguments.
  2. In your C program, argument info must be passed by reference.
  3. A, B, Z, w, ifail, iclustr, gap, work, |rwork, and iwork must have no common elements; otherwise, results are unpredictable.
  4. The NUMROC utility subroutine can be used to determine the values of LOCp(M_) and LOCq(N_) used in the argument descriptions above. For details, see Determining the Number of Rows and Columns in Your Local Arrays and NUMROC--Compute the Number of Rows or Columns of a Block-Cyclically Distributed Matrix Contained in a Process.
  5. The global matrix A must be distributed using a square block-cyclic distribution; that is, MB_A = NB_A.
  6. The global matrix A must be aligned on a block boundary; that is:
  7. In the process grid, the process row containing the first row of the submatrix A must also contain the first row of the submatrix B; that is: iarow = ibrow

    where:

  8. In the process grid, the process column containing the first column of the submatrix A must also contain the first column of the submatrix B; that is: iacol = ibcol

    where:

  9. The block row offset of the global matrix A must be equal to the block row offset of the global matrix B; that is:
  10. The block column offset of the global matrix A must be equal to the block column offset of the global matrix B; that is:
  11. The following values must be equal:
  12. If jobz = 'V', then also follow these rules:
  13. Eigenvectors associated with tightly clustered eigenvalues may not be orthogonal.
  14. Eigenvectors that are on different processes are not reorthogonalized. For details, see the lwork and |lrwork argument.
  15. If lwork = -1 on any process, it must equal -1 on all processes. That is, if a subset of the processes specifies -1 for the work area size, they must all specify -1.
  16. |If lrwork = -1 on any process, it must equal |-1 on all processes. That is, if a subset of the processes |specifies -1 for the work area size, they must all specify |-1.
  17. If liwork = -1 on any process, it must equal -1 on all processes. That is, if a subset of the processes specifies -1 for the work area size, they must all specify -1.
  18. |The following describes the relationship between workspace, |orthogonality, and performance.

    |Let clustersize be the number of eigenvalues in the largest |cluster, where a cluster is defined as a set of close eigenvalues: |

    |{wk, ..., |wk+iclustrsz-1 | |wj+1 <= |wj+orfac(2)(norm(A))} |

    |

    |For PDSYGVX, see the description of work ***. For PZHEGVX, see the description of rwork ***.

Function

This subroutine computes selected eigenvalues and, optionally, the eigenvectors of a real symmetric |or complex Hermitian positive definite generalized eigenproblem:

In the formulas above:

|A represents the global submatrix |Aia:ia+n-1, |ja:ja+n-1
|B represents the global submatrix |Bib:ib+n-1, |jb:jb+n-1

Eigenvalues and eigenvectors can be selected by specifying a range of values or a range of indices for the eigenvalues. The computation involves the following steps:

  1. Compute the Cholesky factorization of B using PDPOTRF |or PZPOTRF.
  2. Reduce the real symmetric |or complex Hermitian positive definite generalized eigenproblem to standard form using PDSYGST |or PZHEGST.
  3. Compute the requested eigenvalues and, optionally, the eigenvectors of the standard form using PDSYEVX |or PZHEEVX.
  4. Backtransform the eigenvectors to obtain the eigenvectors of the original problem.

Error Conditions

Computational Errors
Note:
For more details, see output argument info.
  1. The matrix B is not positive definite. The order of the smallest minor which is not positive definite is stored in ifail1.
  2. Bisection failed to converge for some eigenvalues. The eigenvalues may not be as accurate as the absolute and relative tolerances.
  3. The number of eigenvalues computed does not match the number of eigenvalues requested.
  4. No eigenvalues were computed, because the Gershgorin interval initially used was incorrect.
  5. Some eigenvectors failed to converge. The indices are stored in ifail.
  6. Eigenvectors corresponding to one or more clusters of eigenvalues could not be reorthogonalized because of insufficient workspace. The indices of the clusters are stored in iclustr.
  7. All the eigenvectors between vl and vu could not be computed due to insufficient workspace. The number of eigenvectors computed is returned in nz.

Resource Errors
  1. (lwork = 0, |lrwork = 0, or liwork = 0) and unable to allocate work space

Input-Argument and Miscellaneous Errors

Stage 1 

  1. DTYPE_A is invalid.
  2. DTYPE_B is invalid.
  3. DTYPE_Z is invalid and jobz = 'V'

Stage 2 

  1. CTXT_A is invalid.

Stage 3 

  1. This subroutine has been called from outside the process grid.

Stage 4 

  1. ibtype <> 1, 2, or 3
  2. jobz <> 'N' or 'V'
  3. range <> 'A', 'V', or 'I'
  4. uplo <> 'U' or 'L'
  5. n < 0
  6. M_A < 0 and n = 0; M_A < 1 otherwise
  7. N_A < 0 and n = 0; N_A < 1 otherwise
  8. MB_A < 1
  9. NB_A < 1
  10. RSRC_A < 0 or RSRC_A >= p
  11. CSRC_A < 0 or CSRC_A >= q
  12. ia < 1
  13. ja < 1
  14. M_B < 0 and n = 0; M_B < 1 otherwise
  15. N_B < 0 and n = 0; N_B < 1 otherwise
  16. MB_B < 1
  17. NB_B < 1
  18. RSRC_B < 0 or RSRC_B >= p
  19. CSRC_B < 0 or CSRC_B >= q
  20. ib < 1
  21. jb < 1
  22. CTXT_A <> CTXT_B

    If jobz = 'V':

  23. M_Z < 0 and n = 0; M_Z < 1 otherwise
  24. N_Z < 0 and n = 0; N_Z < 1 otherwise
  25. MB_Z < 1
  26. NB_Z < 1
  27. RSRC_Z < 0 or RSRC_Z >= p
  28. CSRC_Z < 0 or CSRC_Z >= q
  29. iz < 1
  30. jz < 1
  31. CTXT_A <> CTXT_Z

Stage 5 

  1. vu <= vl and range = 'V' and n <> 0
  2. il < 1 and range = 'I'
  3. (iu < min(n, il) or iu > n) and range = 'I'

    If n <> 0:

  4. ia > M_A
  5. ja > N_A
  6. ia+n-1 > M_A
  7. ja+n-1 > N_A
  8. ib > M_B
  9. jb > N_B
  10. ib+n-1 > M_B
  11. jb+n-1 > N_B

    If n <> 0 and jobz = 'V':

  12. iz > M_Z
  13. jz > N_Z
  14. iz+n-1 > M_Z
  15. jz+n-1 > N_Z

    In all cases:

  16. MB_A <> NB_A
  17. mod(ia-1, MB_A)  <>  0
  18. mod(ja-1, NB_A)  <>  0
  19. MB_A <> MB_B
  20. NB_A <> NB_B
  21. mod(ib-1, MB_B) <> mod(ia-1, MB_A)
  22. mod(jb-1, NB_B) <> mod(ja-1, NB_A)
  23. In the process grid, the process row containing the first row of the submatrix A does not contain the first row of the submatrix B; that is, iarow <> ibrow, where:
    iarow = mod(RSRC_A + (ia-1)/MB_A, p)
    ibrow = mod(RSRC_B + (ib-1)/MB_B, p)
  24. In the process grid, the process column containing the first column of the submatrix A does not contain the first column of the submatrix B; that is, iacol <> ibcol, where:
    iacol = mod(CSRC_A + (ja-1)/NB_A, q)
    ibcol = mod(CSRC_B + (jb-1)/NB_B, q)

    If jobz = 'V':

  25. M_A <> M_Z
  26. MB_A <> MB_Z
  27. NB_A <> NB_Z
  28. mod(iz-1, MB_Z)  <>  mod(ia-1, MB_A)
  29. In the process grid, the process row containing the first row of the submatrix A does not contain the first row of the submatrix Z; that is, iarow  <>  izrow, where:
    iarow = mod(RSRC_A + (ia-1)/MB_A, p)
    izrow = mod(RSRC_Z + (iz-1)/MB_Z, p)
  30. RSRC_A <> RSRC_Z
  31. CSRC_A <> CSRC_Z

Stage 6 

  1. LLD_A < max(1, LOCp(M_A))
  2. LLD_B < max(1, LOCp(M_B))
  3. lwork <> 0, lwork <> -1, and lwork < (minimum value). (For the minimum value, see the lwork argument description.)
  4. |lrwork <> 0, lrwork <> -1, and |lrwork < (minimum value). (For |the minimum value, see the lrwork argument description.)
  5. liwork <> 0, liwork <> -1, and liwork < (minimum value). (For the minimum value, see the liwork argument description.)

    If jobz = 'V':

  6. LLD_Z < max(1, LOCp(M_Z))

Stage 7 

    Each of the following global input arguments are checked to determine whether its value differs from the value specified on process P00:

  1. ibtype differs.
  2. jobz differs.
  3. range differs.
  4. uplo differs.
  5. n differs.
  6. ia differs.
  7. ja differs.
  8. DTYPE_A differs.
  9. M_A differs.
  10. N_A differs.
  11. MB_A differs.
  12. NB_A differs.
  13. RSRC_A differs.
  14. CSRC_A differs.
  15. ib differs.
  16. jb differs.
  17. DTYPE_B differs.
  18. M_B differs.
  19. N_B differs.
  20. MB_B differs.
  21. NB_B differs.
  22. RSRC_B differs.
  23. CSRC_B differs.
  24. ABSTOL differs.

    Also:

  25. lwork = -1 on a subset of processes.
  26. |lrwork = -1 on a subset of |processes.
  27. liwork = -1 on a subset of processes.

Stage 8 

    Each of the following global input arguments are checked to determine whether its value differs from the value specified on process P00:

    If range = 'V':

  1. vl differs.
  2. vu differs.

    If range = 'I':

  3. il differs.
  4. iu differs.

    If jobz = 'V':

  5. iz differs.
  6. jz differs.
  7. DTYPE_Z differs.
  8. M_Z differs.
  9. N_Z differs.
  10. MB_Z differs.
  11. NB_Z differs.
  12. RSRC_Z differs.
  13. CSRC_Z differs.
  14. ORFAC differs.

Example 1

This example shows how to find all the eigenvalues and eigenvectors of a |real symmetric positive definite generalized eigenproblem using a 2 × 2 process grid.

Notes:

  1. Because range = 'A', arguments vl, vu, il, and iu are not referenced.

  2. Because lwork = 0 and liwork = 0, PDSYGVX dynamically allocates the work areas used by this subroutine.

Call Statements and Input


 ORDER = 'R'
 NPROW = 2
 NPCOL = 2
 CALL BLACS_GET(0, 0, ICONTXT)
 CALL BLACS_GRIDINIT(ICONTXT, ORDER, NPROW, NPCOL)
 CALL BLACS_GRIDINFO(ICONTXT, NPROW, NPCOL, MYROW, MYCOL)
 
               IBTYPE JOBZ RANGE UPLO  N  A IA JA  DESC_A  B IB JB  DESC_B   VL     VU   IL IU  ABSTOL
                  |    |     |     |   |  |  |  |    |     |  |  |    |       |      |    |  |    |   
 CALL PDSYGVX(    1,  'V',  'A',  'L', 4, A, 1, 1, DESC_A, B, 1, 1, DESC_B, 0.0D0, 0.0D0, 0, 0, -1.0, 
 
                 M  NZ  W   ORFAC  Z IZ JZ  DESC_Z  WORK LWORK IWORK LIWORK IFAIL  ICLUSTR  GAP  INFO
+                |   |  |     |    |  |  |    |      |     |     |      |     |       |      |    |
                 M, NZ, W, -1.0D0, Z, 1, 1, DESC_Z, WORK , 0 , IWORK ,  0 , IFAIL, ICLUSTR, GAP, INFO)             


DESC_A DESC_B DESC_Z
DTYPE_ 1 1 1
CTXT_ icontxt(IITOOTB) icontxt(IITOOTB) icontxt(IITOOTB)
M_ 4 4 4
N_ 4 4 4
MB_ 1 1 1
NB_ 1 1 1
RSRC_ 0 0 0
CSRC_ 0 0 0
LLD_ See below(EPSSTLA) See below(EPSSTLA) See below(EPSSTLA)

Notes:

  1. icontxt is the output of the BLACS_GRIDINIT call.

  2. Each process should set the LLD_ as follows:
    LLD_A = MAX(1,NUMROC(M_A, MB_A, MYROW, RSRC_A, NPROW))
    LLD_B = MAX(1,NUMROC(M_B, MB_B, MYROW, RSRC_B, NPROW))
    LLD_Z = MAX(1,NUMROC(M_Z, MB_Z, MYROW, RSRC_Z, NPROW))
    

    In this example, LLD_A = LLD_B = LLD_Z = 2 on all processes.

Global |real symmetric matrix A of order 4, stored in lower storage mode, with block sizes 1 × 1:

B,D     0        1        2        3
     *                                 *
 0   | -1.0  |    .   |    .   |    .  |
     | ------|--------|--------|------ |
 1   |  1.0  |   1.0  |    .   |    .  |
     | ------|--------|--------|------ |
 2   | -1.0  |  -1.0  |   1.0  |    .  |
     | ------|--------|--------|------ |
 3   |  1.0  |   1.0  |  -1.0  |   1.0 |
     *                                 *

The following is the 2 × 2 process grid:

B,D  |   0 2   | 1 3 
-----| ------- |-----
0    |   P00   |  P01
2    |         |
-----| ------- |-----
1    |   P10   |  P11
3    |         |

Local arrays for A:

p,q  |     0      |      1
-----|------------|------------
 0   | -1.0   .   |    .    . 
     | -1.0  1.0  |  -1.0   . 
-----|------------|------------
 1   |  1.0   .   |   1.0   . 
     |  1.0 -1.0  |   1.0  1.0

Global |real symmetric positive definite matrix B of order 4, stored in lower storage mode, with block sizes 1 × 1:

B,D      0         1         2         3
     *                                     *
 0   |  2.0   |    .    |    .    |    .   |
     | -------|---------|---------|------- |
 1   |  1.0   |   2.0   |    .    |    .   |
     | -------|---------|---------|------- |
 2   |  0.0   |   1.0   |   2.0   |    .   |
     | -------|---------|---------|------- |
 3   |  0.0   |   0.0   |   1.0   |   2.0  |
     *                                     *

The following is the 2 × 2 process grid:

B,D  |   0 2   | 1 3 
-----| ------- |-----
0    |   P00   |  P01
2    |         |
-----| ------- |-----
1    |   P10   |  P11
3    |         |

Local arrays for B:

p,q  |      0       |       1
-----|--------------|--------------
 0   |  2.0    .    |    .     .  
     |  0.0   2.0   |   1.0    .  
-----|--------------|--------------
 1   |  1.0    .    |   2.0    .  
     |  0.0   1.0   |   0.0   2.0 

Output:

The lower triangle, including the diagonal, of the global |real symmetric matrix A is overwritten; i.e., the original input is not preserved.

On all processes, m = 4 and nz = 4.

Global vector w of length 4 is the same on all processes:

w = (-1.5526, 0.0000, 0.0000, 5.1526)

Global |real symmetric positive definite matrix B of order 4, stored in lower storage mode, with block sizes 1 × 1:

B,D       0           1           2           3
     *                                             *
 0   |  1.4142  |     .     |     .     |     .    |
     | ---------|-----------|-----------|--------- |
 1   |  0.7071  |   1.2247  |     .     |     .    |
     | ---------|-----------|-----------|--------- |
 2   |  0.0000  |   0.8165  |   1.1547  |     .    |
     | ---------|-----------|-----------|--------- |
 3   |  0.0000  |   0.0000  |   0.8660  |  1.1180  |
     *                                             *

The following is the 2 × 2 process grid:

B,D  |   0 2   | 1 3 
-----| ------- |-----
0    |   P00   |  P01
2    |         |
-----| ------- |-----
1    |   P10   |  P11
3    |         |

Local arrays for B:

p,q  |        0          |         1
-----|-------------------|-------------------
 0   |  1.4142    .      |    .        .     
     |  0.0000   1.1547  |   0.8165    .     
-----|-------------------|-------------------
 1   |  0.7071    .      |   1.2247    .     
     |  0.0000   0.8660  |   0.0000   1.1180 

Global general matrix Z of order 4, with block sizes 1 × 1:

B,D       0           1           2           3
     *                                             *
 0   |  0.8842  |   0.0000  |   0.0000  |  0.1349  |
     | ---------|-----------|-----------|--------- |
 1   | -0.5619  |   0.3634  |  -0.4098  | -0.7643  |
     | ---------|-----------|-----------|--------- |
 2   |  0.3072  |   0.4266  |   0.1343  |  0.9517  |
     | ---------|-----------|-----------|--------- |
 3   | -0.1198  |   0.0631  |   0.5441  | -0.6969  |
     *                                             *

The following is the 2 × 2 process grid:

B,D  |   0 2   | 1 3 
-----| ------- |-----
0    |   P00   |  P01
2    |         |
-----| ------- |-----
1    |   P10   |  P11
3    |         |

Local arrays for Z:

p,q  |        0          |         1
-----|-------------------|-------------------
 0   |  0.8842   0.0000  |   0.0000   0.1349 
     |  0.3072   0.1343  |   0.4266   0.9517 
-----|-------------------|-------------------
 1   | -0.5619  -0.4098  |   0.3634  -0.7643 
     | -0.1198   0.5441  |   0.0631  -0.6969 

Global vector ifail of length 4 is the same on all processes:

ifail = (0, 0, 0, 0)

Global vector iclustr of length 8 ( = 2(nprow)(npcol)) is the same on all processes:

iclustr = (0, 0, 0, 0, 0, 0, 0, 0)

Global vector gap of length 4 ( = (nprow)(npcol)) is the same on all processes:

gap = (-1.0, -1.0, -1.0, -1.0)

The value of info is 0 on all processes.

|Example 2

|This example shows how to find all the eigenvalues and eigenvectors of a |complex Hermitian positive definite generalized eigenproblem using a |2 × 2 process grid.

|Notes:

  1. |Because range = 'A', arguments vl, |vu, il, and iu are not referenced.

  2. |Because lwork = 0, lrwork = 0, and |liwork = 0, PZHEGVX dynamically allocates the work areas used |by this subroutine. |

|Call Statements and Input


|

| ORDER = 'R'
| NPROW = 2
| NPCOL = 2
| CALL BLACS_GET(0, 0, ICONTXT)
| CALL BLACS_GRIDINIT(ICONTXT, ORDER, NPROW, NPCOL)
| CALL BLACS_GRIDINFO(ICONTXT, NPROW, NPCOL, MYROW, MYCOL)
| 
|               IBTYPE JOBZ RANGE UPLO  N  A IA JA  DESC_A  B IB JB  DESC_B   VL     VU   IL IU  ABSTOL
|                  |    |     |     |   |  |  |  |    |     |  |  |    |       |      |    |  |    |   
| CALL PZHEGVX(    1,  'V',  'A',  'L', 4, A, 1, 1, DESC_A, B, 1, 1, DESC_B, 0.0D0, 0.0D0, 0, 0, -1.0, 
| 
|                 M  NZ  W   ORFAC  Z IZ JZ  DESC_Z  WORK LWORK  RWORK  LRWORK  IWORK LIWORK IFAIL  ICLUSTR
|+                |   |  |     |    |  |  |    |      |     |      |       |      |      |     |       |   
|                 M, NZ, W, -1.0D0, Z, 1, 1, DESC_Z, WORK,  0,   RWORK,    0,   IWORK,   0,  IFAIL, ICLUSTR,   
| 
|                 GAP  INFO
|+                 |    |
|                 GAP, INFO)   
|          

|

DESC_A DESC_B DESC_Z
DTYPE_ 1 1 1
CTXT_ icontxt(IITOOT9) icontxt(IITOOT9) icontxt(IITOOT9)
M_ 4 4 4
N_ 4 4 4
MB_ 1 1 1
NB_ 1 1 1
RSRC_ 0 0 0
CSRC_ 0 0 0
LLD_ See below(EPSSTL9) See below(EPSSTL9) See below(EPSSTL9)

Notes:

  1. icontxt is the output of the BLACS_GRIDINIT call.

  2. Each process should set the LLD_ as follows:
    LLD_A = MAX(1,NUMROC(M_A, MB_A, MYROW, RSRC_A, NPROW))
    LLD_B = MAX(1,NUMROC(M_B, MB_B, MYROW, RSRC_B, NPROW))
    LLD_Z = MAX(1,NUMROC(M_Z, MB_Z, MYROW, RSRC_Z, NPROW))
    

    In this example, LLD_A = LLD_B = LLD_Z = 2 on all processes.

|Global complex Hermitian matrix A of order 4, stored in lower |storage mode, with block sizes 1 × 1:

|B,D          0               1               2                3
|     *                                                               *
| 0   | ( 1.0,  0.0) |        .      |        .       |        .      |
|     | -------------|---------------|----------------|-------------- |
| 1   | ( 5.0, -2.0) | ( 10.0,  0.0) |        .       |        .      |
|     | -------------|---------------|----------------|-------------- |
| 2   | ( 7.0,  4.0) | ( 15.0,  6.0) | ( 20.0,  0.0)  |        .      |
|     | -------------|---------------|----------------|-------------- |
| 3   | ( 9.0, -6.0) | ( 20.0, -4.0) | ( 25.0, -9.0)  | ( 30.0,  0.0) |
|     *                                                               *

|The following is the 2 × 2 process grid:

|B,D  |   0 2   | 1 3 
|-----| ------- |-----
|0    |   P00   |  P01
|2    |         |
|-----| ------- |-----
|1    |   P10   |  P11
|3    |         |

|Local arrays for A:

|p,q  |              0              |              1          
|-----|-----------------------------|----------------------------
| 0   | ( 1.0,   . )         .      |        .              .  
|     | ( 7.0,  4.0)  ( 20.0,  . )  | ( 15.0,  6.0)         .  
|-----|-----------------------------|----------------------------
| 1   | ( 5.0, -2.0)         .      | ( 10.0,   . )         .   
|     | ( 9.0,  6.0)  ( 25.0, -9.0) | ( 20.0, -4.0)  ( 30.0,  . )
| 

|Global complex Hermitian positive definite matrix B of order 4, |stored in lower storage mode, with block sizes 1 × 1:

|B,D          0             1              2               3
|     *                                                          *
| 0   |( 9.0,  0.0) |       .      |       .      |        .     |
|     |-------------|--------------|--------------|--------------|
| 1   |( 3.0, -3.0) |( 18.0,  0.0) |       .      |        .     |
|     |-------------|--------------|--------------|--------------|
| 2   |( 3.0,  3.0) |(  8.0,  6.0) |( 27.0,  0.0) |        .     |
|     |-------------|--------------|--------------|--------------|
| 3   |( 3.0, -3.0) |(  8.0, -6.0) |( 12.0, -9.0) | ( 61.0, 0.0) |
|     *                                                          *

|The following is the 2 × 2 process grid:

|B,D  |   0 2   | 1 3 
|-----| ------- |-----
|0    |   P00   |  P01
|2    |         |
|-----| ------- |-----
|1    |   P10   |  P11
|3    |         |

|Local arrays for B:

|p,q  |              0              |              1          
|-----|-----------------------------|----------------------------
| 0   | ( 9.0,   . )         .      |        .              .  
|     | ( 3.0,  3.0)  ( 27.0,  . )  | (  8.0,  6.0)         .  
|-----|-----------------------------|----------------------------
| 1   | ( 3.0, -3.0)         .      | ( 18.0,   . )         .   
|     | ( 3.0, -3.0)  ( 12.0, -9.0) | (  8.0, -6.0)  ( 61.0,  . )
| 

|Output:

|The lower triangle, including the diagonal, of the global complex Hermitian |matrix A is overwritten; i.e., the original |input is not preserved.

|On all processes, m = 4 and |nz = 4.

|Global vector w of length 4 is the same on all processes: |

|w = (-0.7969, -0.1152, |-0.1669, -1.2304) |

|Global complex Hermitian positive definite matrix B of order 4, |stored in lower storage mode, with block sizes 1 × 1:

|B,D          0             1               2                 3
|     *                                                              *
| 0   | ( 3.0, 0.0) |       .     |         .        |        .      |
|     | ------------|-------------|------------------|-------------- |
| 1   | ( 1.0,-1.0) | ( 4.0, 0.0) |         .        |        .      |
|     | ------------|-------------|------------------|-------------- |
| 2   | ( 1.0, 1.0) | ( 2.0, 1.0) | (4.4721, 0.0   ) |        .      |
|     | ------------|-------------|------------------|-------------- |
| 3   | ( 1.0,-1.0) | ( 1.5,-1.5) | (2.3479, -.5590) | (6.9767, 0.0) |
|     *                                                              *

|The following is the 2 × 2 process grid:

|B,D  |   0 2   | 1 3 
|-----| ------- |-----
|0    |   P00   |  P01
|2    |         |
|-----| ------- |-----
|1    |   P10   |  P11
|3    |         |

|Local arrays for B:

|p,q  |                 0               |              1          
|-----|---------------------------------|------------------------------
| 0   | ( 3.0,  0.0)           .        |       .               .      
|     | ( 1.0,  1.0)  ( 4.4721, 0.0)    | ( 2.0, 1.0)           .      
|-----|---------------------------------|----------------------------- 
| 1   | ( 1.0, -1.0)           .        | ( 4.0, 0.0)           .      
|     | ( 1.0, -1.0)  ( 2.3479, -.5590) | ( 1.5,-1.5)  ( 6.9767, 0.0 ) 
| 

|Global general matrix Z of order 4, with block sizes |1 × 1:

|B,D       0           1           2           3
|     *                                                                           *
| 0   | ( .2784, -.0218) | (-.1300, -.0562) | ( .1842, -.0064) | (-.0692,  .0118) |
|     | -----------------|------------------|------------------|----------------- |
| 1   | ( .0601,  .1506) | ( .1844,  .1016) | (-.0064, -.0516) | (-.0948,  .0008) |
|     | -----------------|------------------|------------------|----------------- |
| 2   | ( .0000, -.1675) | ( .0004, -.0082) | (-.0591,  .1230) | (-.0946,  .0178) |
|     | -----------------|------------------|------------------|----------------- |
| 3   | (-.0406,  .0627) | (-.0726, -.0362) | (-.0414, -.0635) | (-.0513, -.0023) |
|     *                                                                           *

|The following is the 2 × 2 process grid:

|B,D  |   0 2   | 1 3 
|-----| ------- |-----
|0    |   P00   |  P01
|2    |         |
|-----| ------- |-----
|1    |   P10   |  P11
|3    |         |

|Local arrays for Z:

|p,q  |                0                 |                1          
|-----|----------------------------------|--------------------------------
| 0   | ( .2784,-.0218)  ( .1842,-.0064) | (-.1300,-.0562)  (-.0692, .0118) 
|     | ( .0000,-.1675)  (-.0591, .1230) | ( .0004,-.0082)  (-.0946, .0178) 
|-----|----------------------------------|--------------------------------- 
| 1   | ( .0601, .1506)  (-.0064,-.0516) | ( .1844, .1016)  (-.0948, .0008) 
|     | (-.0406, .0627)  (-.0414,-.0635) | (-.0726,-.0362)  (-.0513,-.0023) 

|Global vector ifail of length 4 is the same on all |processes: |

|ifail = (0, 0, 0, 0) |

|Global vector iclustr of length 8 |( = 2(nprow)(npcol)) is the same on all |processes: |

|iclustr = (0, 0, 0, 0, 0, 0, 0, 0) |

|Global vector gap of length 4 |( = (nprow)(npcol)) is the same on all |processes: |

|gap = (-1.0, -1.0, |-1.0, -1.0) |

|The value of info is 0 on all processes.


[ Top of Page | Previous Page | Next Page | Table of Contents | Index ]