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
261
262
263
264
265
266
267
|
*> \brief \b ZLAGHE
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
* Definition:
* ===========
*
* SUBROUTINE ZLAGHE( N, K, D, A, LDA, ISEED, WORK, INFO )
*
* .. Scalar Arguments ..
* INTEGER INFO, K, LDA, N
* ..
* .. Array Arguments ..
* INTEGER ISEED( 4 )
* DOUBLE PRECISION D( * )
* COMPLEX*16 A( LDA, * ), WORK( * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> ZLAGHE generates a complex hermitian matrix A, by pre- and post-
*> multiplying a real diagonal matrix D with a random unitary matrix:
*> A = U*D*U'. The semi-bandwidth may then be reduced to k by additional
*> unitary transformations.
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrix A. N >= 0.
*> \endverbatim
*>
*> \param[in] K
*> \verbatim
*> K is INTEGER
*> The number of nonzero subdiagonals within the band of A.
*> 0 <= K <= N-1.
*> \endverbatim
*>
*> \param[in] D
*> \verbatim
*> D is DOUBLE PRECISION array, dimension (N)
*> The diagonal elements of the diagonal matrix D.
*> \endverbatim
*>
*> \param[out] A
*> \verbatim
*> A is COMPLEX*16 array, dimension (LDA,N)
*> The generated n by n hermitian matrix A (the full matrix is
*> stored).
*> \endverbatim
*>
*> \param[in] LDA
*> \verbatim
*> LDA is INTEGER
*> The leading dimension of the array A. LDA >= N.
*> \endverbatim
*>
*> \param[in,out] ISEED
*> \verbatim
*> ISEED is INTEGER array, dimension (4)
*> On entry, the seed of the random number generator; the array
*> elements must be between 0 and 4095, and ISEED(4) must be
*> odd.
*> On exit, the seed is updated.
*> \endverbatim
*>
*> \param[out] WORK
*> \verbatim
*> WORK is COMPLEX*16 array, dimension (2*N)
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> = 0: successful exit
*> < 0: if INFO = -i, the i-th argument had an illegal value
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \date November 2011
*
*> \ingroup complex16_matgen
*
* =====================================================================
SUBROUTINE ZLAGHE( N, K, D, A, LDA, ISEED, WORK, INFO )
*
* -- LAPACK auxiliary routine (version 3.4.0) --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
* November 2011
*
* .. Scalar Arguments ..
INTEGER INFO, K, LDA, N
* ..
* .. Array Arguments ..
INTEGER ISEED( 4 )
DOUBLE PRECISION D( * )
COMPLEX*16 A( LDA, * ), WORK( * )
* ..
*
* =====================================================================
*
* .. Parameters ..
COMPLEX*16 ZERO, ONE, HALF
PARAMETER ( ZERO = ( 0.0D+0, 0.0D+0 ),
$ ONE = ( 1.0D+0, 0.0D+0 ),
$ HALF = ( 0.5D+0, 0.0D+0 ) )
* ..
* .. Local Scalars ..
INTEGER I, J
DOUBLE PRECISION WN
COMPLEX*16 ALPHA, TAU, WA, WB
* ..
* .. External Subroutines ..
EXTERNAL XERBLA, ZAXPY, ZGEMV, ZGERC, ZHEMV, ZHER2,
$ ZLARNV, ZSCAL
* ..
* .. External Functions ..
DOUBLE PRECISION DZNRM2
COMPLEX*16 ZDOTC
EXTERNAL DZNRM2, ZDOTC
* ..
* .. Intrinsic Functions ..
INTRINSIC ABS, DBLE, DCONJG, MAX
* ..
* .. Executable Statements ..
*
* Test the input arguments
*
INFO = 0
IF( N.LT.0 ) THEN
INFO = -1
ELSE IF( K.LT.0 .OR. K.GT.N-1 ) THEN
INFO = -2
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -5
END IF
IF( INFO.LT.0 ) THEN
CALL XERBLA( 'ZLAGHE', -INFO )
RETURN
END IF
*
* initialize lower triangle of A to diagonal matrix
*
DO 20 J = 1, N
DO 10 I = J + 1, N
A( I, J ) = ZERO
10 CONTINUE
20 CONTINUE
DO 30 I = 1, N
A( I, I ) = D( I )
30 CONTINUE
*
* Generate lower triangle of hermitian matrix
*
DO 40 I = N - 1, 1, -1
*
* generate random reflection
*
CALL ZLARNV( 3, ISEED, N-I+1, WORK )
WN = DZNRM2( N-I+1, WORK, 1 )
WA = ( WN / ABS( WORK( 1 ) ) )*WORK( 1 )
IF( WN.EQ.ZERO ) THEN
TAU = ZERO
ELSE
WB = WORK( 1 ) + WA
CALL ZSCAL( N-I, ONE / WB, WORK( 2 ), 1 )
WORK( 1 ) = ONE
TAU = DBLE( WB / WA )
END IF
*
* apply random reflection to A(i:n,i:n) from the left
* and the right
*
* compute y := tau * A * u
*
CALL ZHEMV( 'Lower', N-I+1, TAU, A( I, I ), LDA, WORK, 1, ZERO,
$ WORK( N+1 ), 1 )
*
* compute v := y - 1/2 * tau * ( y, u ) * u
*
ALPHA = -HALF*TAU*ZDOTC( N-I+1, WORK( N+1 ), 1, WORK, 1 )
CALL ZAXPY( N-I+1, ALPHA, WORK, 1, WORK( N+1 ), 1 )
*
* apply the transformation as a rank-2 update to A(i:n,i:n)
*
CALL ZHER2( 'Lower', N-I+1, -ONE, WORK, 1, WORK( N+1 ), 1,
$ A( I, I ), LDA )
40 CONTINUE
*
* Reduce number of subdiagonals to K
*
DO 60 I = 1, N - 1 - K
*
* generate reflection to annihilate A(k+i+1:n,i)
*
WN = DZNRM2( N-K-I+1, A( K+I, I ), 1 )
WA = ( WN / ABS( A( K+I, I ) ) )*A( K+I, I )
IF( WN.EQ.ZERO ) THEN
TAU = ZERO
ELSE
WB = A( K+I, I ) + WA
CALL ZSCAL( N-K-I, ONE / WB, A( K+I+1, I ), 1 )
A( K+I, I ) = ONE
TAU = DBLE( WB / WA )
END IF
*
* apply reflection to A(k+i:n,i+1:k+i-1) from the left
*
CALL ZGEMV( 'Conjugate transpose', N-K-I+1, K-1, ONE,
$ A( K+I, I+1 ), LDA, A( K+I, I ), 1, ZERO, WORK, 1 )
CALL ZGERC( N-K-I+1, K-1, -TAU, A( K+I, I ), 1, WORK, 1,
$ A( K+I, I+1 ), LDA )
*
* apply reflection to A(k+i:n,k+i:n) from the left and the right
*
* compute y := tau * A * u
*
CALL ZHEMV( 'Lower', N-K-I+1, TAU, A( K+I, K+I ), LDA,
$ A( K+I, I ), 1, ZERO, WORK, 1 )
*
* compute v := y - 1/2 * tau * ( y, u ) * u
*
ALPHA = -HALF*TAU*ZDOTC( N-K-I+1, WORK, 1, A( K+I, I ), 1 )
CALL ZAXPY( N-K-I+1, ALPHA, A( K+I, I ), 1, WORK, 1 )
*
* apply hermitian rank-2 update to A(k+i:n,k+i:n)
*
CALL ZHER2( 'Lower', N-K-I+1, -ONE, A( K+I, I ), 1, WORK, 1,
$ A( K+I, K+I ), LDA )
*
A( K+I, I ) = -WA
DO 50 J = K + I + 1, N
A( J, I ) = ZERO
50 CONTINUE
60 CONTINUE
*
* Store full hermitian matrix
*
DO 80 J = 1, N
DO 70 I = J + 1, N
A( J, I ) = DCONJG( A( I, J ) )
70 CONTINUE
80 CONTINUE
RETURN
*
* End of ZLAGHE
*
END
|