PDSYTRD reduces a real symmetric matrix A to symmetric tridiagonal form T by an orthogonal similarity transformation:
where A represents the global real symmetric submatrix Aia:ia+n-1, ja:ja+n-1.
PZHETRD reduces a complex Hermitian matrix A to symmetric tridiagonal form T by a unitary similarity transformation:
where A represents the global complex Hermitian submatrix Aia:ia+n-1, ja:ja+n-1
If n = 0, no computation is performed and the subroutine returns after doing some parameter checking.
A,tau, work | d,e | Subroutine |
Long-precision real | Long-precision real | PDSYTRD |
Long-precision complex | Long-precision real | PZHETRD |
Fortran | CALL PDSYTRD | PZHETRD (uplo, n, a, ia, ja, desc_a, d, e, tau, work, lwork, info) |
C and C++ | pdsytrd | pzhetrd (uplo, n, a, ia, ja, desc_a, d, e, tau, work, lwork, info); |
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 104. 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.
If lwork = 0, work is ignored.
If lwork <> 0, work is the work area used by this subroutine, where:
Scope: local
Specified as: an area of storage containing numbers of data type indicated in Table 104.
Scope:
Specified as: a fullword integer; where:
lwork >= max(nb(np+1), 3nb)
where:
See Function, for more information.
Scope: local
Returned as: an LLD_A by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 104. Details about the square block-cyclic data distribution of global matrix A are stored in desc_a.
This identifies the first element of the local array D. This subroutine computes the location of the first element of the local subarray used, based on ja, desc_a, p, q, myrow, and mycol; therefore, the leading 1 by LOCq(ja+n-1) part of the local array D must contain the local pieces of the leading 1 by ja+n-1 part of the global matrix D.
A copy of the vector d, with a block size of NB_A and global index ja, is returned to each row of the process grid. The process column over which the first column of d is distributed is CSRC_A.
Scope: local
Returned as: a 1 by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 104.
If uplo = 'U', then eja = 0 and eja+1:ja+n-1 contains the superdiagonal elements of the tridiagonal matrix T.
If uplo = 'L', then eja:ja+n-2 contains the subdiagonal elements of the tridiagonal matrix T, and eja+n-1 = 0.
This identifies the first element of the local array E. This subroutine computes the location of the first element of the local subarray used, based on ja, desc_a, p, q, myrow, and mycol; therefore, the leading 1 by LOCq(ja+n-1) part of the local array E must contain the local pieces of the leading 1 by ja+n-1 part of the global matrix E.
A copy of the vector e, with a block size of NB_A and global index ja, is returned to each row of the process grid. The process column over which the first column of E is distributed is CSRC_A.
Scope: local
Returned as: a 1 by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 104.
If uplo = 'U', then tauja is zero and tauja+1:ja+n-1 contains the scalar factors of the elementary reflectors.
If uplo = 'L', then tauja:ja+n-2 contains the scalar factors of the elementary reflectors and tauja+n-1 is zero.
This identifies the first element of the local array tau. This subroutine computes the location of the first element of the local subarray used, based on ja, desc_a, p, q, myrow, and mycol; therefore, the leading 1 by LOCq(ja+n-1) part of the local array tau must contain the local pieces of the leading 1 by ja+n-1 part of the global matrix tau.
A copy of the vector tau, with a block size of NB_A and global index ja, is returned to each row of the process grid. The process column over which the first column of tau is distributed is CSRC_A.
Scope: local
Returned as: a 1 by (at least) LOCq(N_A) array, containing numbers of the data type indicated in Table 104.
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, where:
If lwork >= 1 or lwork = -1, then work1 is set to the minimum lwork value and contains numbers of the data type indicated in Table 104. Except for work1, the contents of work are overwritten on return.
Scope: global
Returned as: a fullword integer; info = 0.
PDSYTRD reduces a real symmetric matrix A to symmetric tridiagonal form T by an orthogonal similarity transformation:
where:
where:
If uplo = 'U', then the following example shows the contents of A on return with n = 5 and ia = ja = 1:
where:
where:
If uplo = 'L', then the following example shows the contents of A on return with n = 5 and ia = ja = 1:
where:
PZHETRD reduces a complex Hermitian matrix A to symmetric tridiagonal form T by a unitary similarity transformation:
where:
where:
If uplo = 'U', then the following example shows the contents of A on return with n = 5 and ia = ja = 1:
where:
where:
If uplo = 'L', then the following example shows the contents of A on return with n = 5 and ia = ja = 1:
where:
None
Stage 5: If n <> 0:
In all cases:
where:
Each of the following global input arguments are checked to determine whether its value differs from the value specified on process P00:
Also:
This example shows the reduction of a real symmetric matrix of order 4 to symmetric tridiagonal form, using a 2 × 2 process grid.
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) UPLO N A IA JA DESC_A D E TAU WORK LWORK INFO | | | | | | | | | | | | CALL PDSYTRD( 'U' , 4 , A , 1 , 1 , DESC_A , D , E , TAU , WORK , 0 , INFO )
| DESC_A |
---|---|
DTYPE_ | 1 |
CTXT_ | icontxt(IITOO10) |
M_ | 4 |
N_ | 4 |
MB_ | 1 |
NB_ | 1 |
RSRC_ | 0 |
CSRC_ | 0 |
LLD_ | See below(EPSST10) |
Notes: |
Global real symmetric matrix A of order 4 with block sizes 1 × 1:
B,D 0 1 2 3 * * 0 | 5.0 | 4.0 | 1.0 | 1.0 | | ------|--------|--------|------ | 1 | . | 5.0 | 1.0 | 1.0 | | ------|--------|--------|------ | 2 | . | . | 4.0 | 2.0 | | ------|--------|--------|------ | 3 | . | . | . | 4.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 | 5.0 1.0 | 4.0 1.0 | . 4.0 | . 2.0 -----|------------|------------ 1 | . 1.0 | 5.0 1.0 | . . | . 4.0
Output:
Global real symmetric matrix A of order 4 with block sizes 1 × 1:
B,D 0 1 2 3 * * 0 | 1.00 | 0.00 | 0.41 | 0.22 | | -------|---------|---------|------- | 1 | . | 6.00 | 2.83 | 0.22 | | -------|---------|---------|------- | 2 | . | . | 7.00 | -2.45 | | -------|---------|---------|------- | 3 | . | . | . | 4.00 | * *
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.00 0.41 | 0.00 0.22 | . 7.00 | . -2.45 -----|--------------|-------------- 1 | . 2.83 | 6.00 0.22 | . . | . 4.00
Global row vector D of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | 1.00 | 6.00 | 7.00 | 4.00 | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 -----| ------- |----- | P10 | P11
Local arrays for D:
p,q | 0 | 1 -----|--------------|-------------- 0 | 1.00 7.00 | 6.00 4.00 -----|--------------|-------------- 1 | 1.00 7.00 | 6.00 4.00
Global row vector E of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | 0.00 | 0.00 | 2.83 | -2.45 | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 -----| ------- |----- | P10 | P11
Local arrays for E:
p,q | 0 | 1 -----|--------------|-------------- 0 | 0.00 2.83 | 0.00 -2.45 -----|--------------|-------------- 1 | 0.00 2.83 | 0.00 -2.45
Global row vector tau of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | 0.00 | 0.00 | 1.71 | 1.82 | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 -----| ------- |----- | P10 | P11
Local arrays for tau:
p,q | 0 | 1 -----|--------------|-------------- 0 | 0.00 1.71 | 0.00 1.82 -----|--------------|-------------- 1 | 0.00 1.71 | 0.00 1.82
The value of info is 0 on all processes.
This example shows the reduction of a complex Hermitian matrix of order 4 to symmetric tridiagonal form, using a 2 × 2 process grid.
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) UPLO N A IA JA DESC_A D E TAU WORK LWORK INFO | | | | | | | | | | | | CALL PZHETRD( 'L' , 4 , A , 1 , 1 , DESC_A , D , E , TAU , WORK , 0 , INFO )
| DESC_A |
---|---|
DTYPE_ | 1 |
CTXT_ | icontxt(IITOOT2) |
M_ | 4 |
N_ | 4 |
MB_ | 1 |
NB_ | 1 |
RSRC_ | 0 |
CSRC_ | 0 |
LLD_ | See below(EPSSTL2) |
Notes: |
Global complex Hermitian matrix A of order 4 with block sizes 1 × 1:
B,D 0 1 2 3 * * 0 | (5.0, 0.0) | . | . | . | |------------|------------|------------|------------| 1 | (4.0, 1.0) | (5.0, 0.0) | . | . | |------------|------------|------------|------------| 2 | (1.0, 2.0) | (1.0, 0.0) | (4.0, 0.0) | . | |------------|------------|------------|------------| 3 | (2.0, 3.0) | (3.0, 2.0) | (5.0, 1.0) | (4.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 | (5.0, . ) . | . . | (1.0, 2.0) (4.0, . ) | (1.0, 0.0) . -----|------------------------|------------------------ 1 | (4.0, 1.0) . | (5.0, . ) . | (2.3, 3.0) (5.0, 1.0) | (3.0, 2.0) (4.0, . )
Output:
Global complex Hermitian matrix A of order 4 with block sizes 1 × 1:
B,D 0 1 2 3 * * 0 | ( 5.00, 0.00) | . | . | . | |---------------|---------------|--------------|---------------| 1 | (-5.92, 0.00) | (10.09, 0.00) | . | . | |---------------|---------------|--------------|---------------| 2 | ( 0.12, 0.19) | ( 2.36, 0.00) | ( 4.16, 0.0) | . | |---------------|---------------|--------------|---------------| 3 | ( 0.23, 0.28) | ( 0.14, 0.19) | ( 1.62, 0.00)| (-1.25, 0.00) | * *
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 | ( 5.00, 0.00) . | . . | ( 0.12, 0.19) ( 4.16, 0.00) | ( 2.36, 0.00) . -----|------------------------------|----------------------------- 1 | (-5.92, 0.00) . | (10.09, 0.00) . | ( 0.23, 0.28) ( 1.62, 0.00) | ( 0.14, 0.19) (-1.25, 0.00)
Global row vector D of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | 5.00 | 10.09 | 4.16 | -1.25 | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 -----| ------- |----- | P10 | P11
Local arrays for D:
p,q | 0 | 1 -----|--------------|-------------- 0 | 5.00 4.16 | 10.09 -1.25 -----|--------------|-------------- 1 | 5.00 4.16 | 10.09 -1.25
Global row vector E of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | -5.92 | 2.36 | 1.62 | 0.00 | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 | | -----| ------- |----- | P10 | P11 | |
Local arrays for E:
p,q | 0 | 1 -----|--------------|-------------- 0 | -5.92 1.62 | 2.36 0.00 -----|--------------|-------------- 1 | -5.92 1.62 | 2.36 0.00
Global row vector tau of length 4 with block size 1:
B,D 0 1 2 3 * * 0 | (1.68, 0.17) | (1.87, 0.21) | (1.96, -0.27) | (0.00, 0.00) | * *
The following is the 2 × 2 process grid:
B,D | 0 2 | 1 3 -----| ------- |----- | P00 | P01 -----| ------- |----- | P10 | P11
Local arrays for tau:
p,q | 0 | 1 -----|----------------------------|-------------------------- 0 | (1.68, 0.17) (1.96, -0.27) | (1.87, 0.21) (0.00, 0.00) -----|----------------------------|-------------------------- 1 | (1.68, 0.17) (1.96, -0.27) | (1.87, 0.21) (0.00, 0.00)
The value of info is 0 on all processes.