|These subroutines compute selected eigenvalues and, optionally, the |eigenvectors of a real symmetric or complex Hermitian positive definite |generalized eigenproblem: |
|In the formulas above: |
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].
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 |
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); |
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.
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'.
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'.
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'.
Scope: global
Specified as: a fullword integer; n >= 0.
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.
Scope: global
Specified as: a fullword integer; 1 <= ia <= M_A and ia+n-1 <= M_A.
Scope: global
Specified as: a fullword integer; 1 <= ja <= N_A and ja+n-1 <= N_A.
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.
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.
Scope: global
Specified as: a fullword integer; 1 <= ib <= M_B and ib+n-1 <= M_B.
Scope: global
Specified as: a fullword integer; 1 <= jb <= N_B and jb+n-1 <= N_B.
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.
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.
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.
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.
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.
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:
where unfl is the underflow threshold, then:
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.
Scope: global
Specified as: a number of the data type indicated in Table 103.
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.
Scope: global
Specified as: a fullword integer; 1 <= iz <= M_Z and iz+n-1 <= M_Z.
Scope: global
Specified as: a fullword integer; 1 <= jz <= N_Z and jz+n-1 <= N_Z.
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. |
|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.
|Scope: |
|Specified as: a fullword integer; where: |
|For PDSYGVX, lwork >= |3(n+ja-1)+2n +max(5nn, |(nb)(np+1))
|For PZHEGVX, lwork >= n+ja-1 |+max(3nb, nb(np+1))
|For PDSYGVX, lwork >= |3(n+ja-1)+2n +max(5nn, |(np0)(mq0)+2(nb)(nb)) + |iceil(neig, (nprow)(npcol))(nn)
|For PZHEGVX, lwork >= n+ja-1 |+(np0+mq0+nb)(nb)
|where: |
|For PDSYGVX, the computed eigenvectors may not be orthogonal if the minimum |workspace is supplied and ortol is too small; therefore, if you |want to guarantee orthogonality (at the cost of potentially compromising |performance), you should add the following to lwork:
|(clustersize-1)(n)
|where clustersize is the number of eigenvalues in the largest |cluster, where a cluster is defined as a set of close eigenvalues: |
|
|When lwork is too small: |
|For the relationship between workspace, orthogonality, and performance, see |Notes and Coding Rules 18. |
|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.
|Scope: |
|Specified as: a fullword integer; where: |
|lrwork >= 2(n+ja-1) |+2n +5nn
|lrwork >= 2(n+ja-1) |+2n +max(5nn, (np0)(mq0)) + |iceil(neig, (nprow)(npcol))(nn)
|where: |
|For PZHEGVX, the computed eigenvectors may not be orthogonal if the minimum |workspace is supplied and ortol is too small; therefore, if you |want to guarantee orthogonality (at the cost of potentially compromising |performance), you should add the following to lrwork:
|(clustersize-1)(n)
|where clustersize is the number of eigenvalues in the largest |cluster, where a cluster is defined as a set of close eigenvalues: |
|
|When lrwork is too small: |
|For the relationship between workspace, orthogonality, and performance, see |Notes and Coding Rules 18. |
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.
Scope:
Specified as: a fullword integer; where:
liwork >= max(isizestein, isizestebz)+2n
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.
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.
Scope: global
Returned as: a fullword integer; 0 <= m <= n.
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.
Scope: global
Returned as: a one-dimensional array of (at least) length n, containing numbers of the data type indicated in Table 103.
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.
|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: |
|For PDSYGVX, work1 is set to the minimum lwork |needed to compute all eigenvectors, but not necessarily sufficient to |guarantee orthogonality of the eigenvectors.
|For PZHEGVX, work1 is set to the minimum lwork |needed to compute all eigenvectors.
|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: |
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.
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.
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:
Scope: global
Returned as: a one-dimensional array of (at least) length 2(nprow)(npcol), containing fullword integers; 0 <= iclustri <= n.
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.
If info = 0, then no input-argument errors or computational errors occurred. This indicates a normal exit.
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.
where:
where:
where:
|Let clustersize be the number of eigenvalues in the largest |cluster, where a cluster is defined as a set of close eigenvalues: |
|
|
|
|then providing enough space to compute all the eigenvectors orthogonally |causes serious degradation in performance. In the limit |(clustersize = n-1), performance may be no |better than using one process.
|
|
|then reorthogonalizing all eigenvectors increases the total execution time |by a factor of 2 or more.
|
|
|then execution time grows as the square of the cluster size, assuming all |other factors remain equal and there is enough workspace. Less |workspace means less reorthogonalization, but faster execution. |
|For PDSYGVX, see the description of work ***. For PZHEGVX, see the description of rwork ***.
This subroutine computes selected eigenvalues and, optionally, the eigenvectors of a real symmetric |or complex Hermitian positive definite generalized eigenproblem:
In the formulas above:
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:
If jobz = 'V':
If n <> 0:
If n <> 0 and jobz = 'V':
In all cases:
If jobz = 'V':
If jobz = 'V':
Each of the following global input arguments are checked to determine whether its value differs from the value specified on process P00:
Also:
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':
If range = 'I':
If jobz = 'V':
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:
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:
|
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:
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:
Global vector iclustr of length 8 ( = 2(nprow)(npcol)) is the same on all processes:
Global vector gap of length 4 ( = (nprow)(npcol)) is the same on all processes:
The value of info is 0 on all processes.
|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:
|
| 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:
|
|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: |
|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: |
|Global vector iclustr of length 8 |( = 2(nprow)(npcol)) is the same on all |processes: |
|Global vector gap of length 4 |( = (nprow)(npcol)) is the same on all |processes: |
|The value of info is 0 on all processes.