1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
|
SUBROUTINE ZHEGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
*
* -- LAPACK routine (version 3.3.1) --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
* -- April 2011 --
*
* .. Scalar Arguments ..
CHARACTER UPLO
INTEGER INFO, ITYPE, LDA, LDB, N
* ..
* .. Array Arguments ..
COMPLEX*16 A( LDA, * ), B( LDB, * )
* ..
*
* Purpose
* =======
*
* ZHEGST reduces a complex Hermitian-definite generalized
* eigenproblem to standard form.
*
* If ITYPE = 1, the problem is A*x = lambda*B*x,
* and A is overwritten by inv(U**H)*A*inv(U) or inv(L)*A*inv(L**H)
*
* If ITYPE = 2 or 3, the problem is A*B*x = lambda*x or
* B*A*x = lambda*x, and A is overwritten by U*A*U**H or L**H*A*L.
*
* B must have been previously factorized as U**H*U or L*L**H by ZPOTRF.
*
* Arguments
* =========
*
* ITYPE (input) INTEGER
* = 1: compute inv(U**H)*A*inv(U) or inv(L)*A*inv(L**H);
* = 2 or 3: compute U*A*U**H or L**H*A*L.
*
* UPLO (input) CHARACTER*1
* = 'U': Upper triangle of A is stored and B is factored as
* U**H*U;
* = 'L': Lower triangle of A is stored and B is factored as
* L*L**H.
*
* N (input) INTEGER
* The order of the matrices A and B. N >= 0.
*
* A (input/output) COMPLEX*16 array, dimension (LDA,N)
* On entry, the Hermitian matrix A. If UPLO = 'U', the leading
* N-by-N upper triangular part of A contains the upper
* triangular part of the matrix A, and the strictly lower
* triangular part of A is not referenced. If UPLO = 'L', the
* leading N-by-N lower triangular part of A contains the lower
* triangular part of the matrix A, and the strictly upper
* triangular part of A is not referenced.
*
* On exit, if INFO = 0, the transformed matrix, stored in the
* same format as A.
*
* LDA (input) INTEGER
* The leading dimension of the array A. LDA >= max(1,N).
*
* B (input) COMPLEX*16 array, dimension (LDB,N)
* The triangular factor from the Cholesky factorization of B,
* as returned by ZPOTRF.
*
* LDB (input) INTEGER
* The leading dimension of the array B. LDB >= max(1,N).
*
* INFO (output) INTEGER
* = 0: successful exit
* < 0: if INFO = -i, the i-th argument had an illegal value
*
* =====================================================================
*
* .. Parameters ..
DOUBLE PRECISION ONE
PARAMETER ( ONE = 1.0D+0 )
COMPLEX*16 CONE, HALF
PARAMETER ( CONE = ( 1.0D+0, 0.0D+0 ),
$ HALF = ( 0.5D+0, 0.0D+0 ) )
* ..
* .. Local Scalars ..
LOGICAL UPPER
INTEGER K, KB, NB
* ..
* .. External Subroutines ..
EXTERNAL XERBLA, ZHEGS2, ZHEMM, ZHER2K, ZTRMM, ZTRSM
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MIN
* ..
* .. External Functions ..
LOGICAL LSAME
INTEGER ILAENV
EXTERNAL LSAME, ILAENV
* ..
* .. Executable Statements ..
*
* Test the input parameters.
*
INFO = 0
UPPER = LSAME( UPLO, 'U' )
IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
INFO = -1
ELSE IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
INFO = -2
ELSE IF( N.LT.0 ) THEN
INFO = -3
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -5
ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
INFO = -7
END IF
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'ZHEGST', -INFO )
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
* Determine the block size for this environment.
*
NB = ILAENV( 1, 'ZHEGST', UPLO, N, -1, -1, -1 )
*
IF( NB.LE.1 .OR. NB.GE.N ) THEN
*
* Use unblocked code
*
CALL ZHEGS2( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
ELSE
*
* Use blocked code
*
IF( ITYPE.EQ.1 ) THEN
IF( UPPER ) THEN
*
* Compute inv(U**H)*A*inv(U)
*
DO 10 K = 1, N, NB
KB = MIN( N-K+1, NB )
*
* Update the upper triangle of A(k:n,k:n)
*
CALL ZHEGS2( ITYPE, UPLO, KB, A( K, K ), LDA,
$ B( K, K ), LDB, INFO )
IF( K+KB.LE.N ) THEN
CALL ZTRSM( 'Left', UPLO, 'Conjugate transpose',
$ 'Non-unit', KB, N-K-KB+1, CONE,
$ B( K, K ), LDB, A( K, K+KB ), LDA )
CALL ZHEMM( 'Left', UPLO, KB, N-K-KB+1, -HALF,
$ A( K, K ), LDA, B( K, K+KB ), LDB,
$ CONE, A( K, K+KB ), LDA )
CALL ZHER2K( UPLO, 'Conjugate transpose', N-K-KB+1,
$ KB, -CONE, A( K, K+KB ), LDA,
$ B( K, K+KB ), LDB, ONE,
$ A( K+KB, K+KB ), LDA )
CALL ZHEMM( 'Left', UPLO, KB, N-K-KB+1, -HALF,
$ A( K, K ), LDA, B( K, K+KB ), LDB,
$ CONE, A( K, K+KB ), LDA )
CALL ZTRSM( 'Right', UPLO, 'No transpose',
$ 'Non-unit', KB, N-K-KB+1, CONE,
$ B( K+KB, K+KB ), LDB, A( K, K+KB ),
$ LDA )
END IF
10 CONTINUE
ELSE
*
* Compute inv(L)*A*inv(L**H)
*
DO 20 K = 1, N, NB
KB = MIN( N-K+1, NB )
*
* Update the lower triangle of A(k:n,k:n)
*
CALL ZHEGS2( ITYPE, UPLO, KB, A( K, K ), LDA,
$ B( K, K ), LDB, INFO )
IF( K+KB.LE.N ) THEN
CALL ZTRSM( 'Right', UPLO, 'Conjugate transpose',
$ 'Non-unit', N-K-KB+1, KB, CONE,
$ B( K, K ), LDB, A( K+KB, K ), LDA )
CALL ZHEMM( 'Right', UPLO, N-K-KB+1, KB, -HALF,
$ A( K, K ), LDA, B( K+KB, K ), LDB,
$ CONE, A( K+KB, K ), LDA )
CALL ZHER2K( UPLO, 'No transpose', N-K-KB+1, KB,
$ -CONE, A( K+KB, K ), LDA,
$ B( K+KB, K ), LDB, ONE,
$ A( K+KB, K+KB ), LDA )
CALL ZHEMM( 'Right', UPLO, N-K-KB+1, KB, -HALF,
$ A( K, K ), LDA, B( K+KB, K ), LDB,
$ CONE, A( K+KB, K ), LDA )
CALL ZTRSM( 'Left', UPLO, 'No transpose',
$ 'Non-unit', N-K-KB+1, KB, CONE,
$ B( K+KB, K+KB ), LDB, A( K+KB, K ),
$ LDA )
END IF
20 CONTINUE
END IF
ELSE
IF( UPPER ) THEN
*
* Compute U*A*U**H
*
DO 30 K = 1, N, NB
KB = MIN( N-K+1, NB )
*
* Update the upper triangle of A(1:k+kb-1,1:k+kb-1)
*
CALL ZTRMM( 'Left', UPLO, 'No transpose', 'Non-unit',
$ K-1, KB, CONE, B, LDB, A( 1, K ), LDA )
CALL ZHEMM( 'Right', UPLO, K-1, KB, HALF, A( K, K ),
$ LDA, B( 1, K ), LDB, CONE, A( 1, K ),
$ LDA )
CALL ZHER2K( UPLO, 'No transpose', K-1, KB, CONE,
$ A( 1, K ), LDA, B( 1, K ), LDB, ONE, A,
$ LDA )
CALL ZHEMM( 'Right', UPLO, K-1, KB, HALF, A( K, K ),
$ LDA, B( 1, K ), LDB, CONE, A( 1, K ),
$ LDA )
CALL ZTRMM( 'Right', UPLO, 'Conjugate transpose',
$ 'Non-unit', K-1, KB, CONE, B( K, K ), LDB,
$ A( 1, K ), LDA )
CALL ZHEGS2( ITYPE, UPLO, KB, A( K, K ), LDA,
$ B( K, K ), LDB, INFO )
30 CONTINUE
ELSE
*
* Compute L**H*A*L
*
DO 40 K = 1, N, NB
KB = MIN( N-K+1, NB )
*
* Update the lower triangle of A(1:k+kb-1,1:k+kb-1)
*
CALL ZTRMM( 'Right', UPLO, 'No transpose', 'Non-unit',
$ KB, K-1, CONE, B, LDB, A( K, 1 ), LDA )
CALL ZHEMM( 'Left', UPLO, KB, K-1, HALF, A( K, K ),
$ LDA, B( K, 1 ), LDB, CONE, A( K, 1 ),
$ LDA )
CALL ZHER2K( UPLO, 'Conjugate transpose', K-1, KB,
$ CONE, A( K, 1 ), LDA, B( K, 1 ), LDB,
$ ONE, A, LDA )
CALL ZHEMM( 'Left', UPLO, KB, K-1, HALF, A( K, K ),
$ LDA, B( K, 1 ), LDB, CONE, A( K, 1 ),
$ LDA )
CALL ZTRMM( 'Left', UPLO, 'Conjugate transpose',
$ 'Non-unit', KB, K-1, CONE, B( K, K ), LDB,
$ A( K, 1 ), LDA )
CALL ZHEGS2( ITYPE, UPLO, KB, A( K, K ), LDA,
$ B( K, K ), LDB, INFO )
40 CONTINUE
END IF
END IF
END IF
RETURN
*
* End of ZHEGST
*
END
|