LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zhbevd_2stage.f
Go to the documentation of this file.
1*> \brief <b> ZHBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* @precisions fortran z -> s d c
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> Download ZHBEVD_2STAGE + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhbevd_2stage.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhbevd_2stage.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhbevd_2stage.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZHBEVD_2STAGE( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ,
22* WORK, LWORK, RWORK, LRWORK, IWORK,
23* LIWORK, INFO )
24*
25* IMPLICIT NONE
26*
27* .. Scalar Arguments ..
28* CHARACTER JOBZ, UPLO
29* INTEGER INFO, KD, LDAB, LDZ, LIWORK, LRWORK, LWORK, N
30* ..
31* .. Array Arguments ..
32* INTEGER IWORK( * )
33* DOUBLE PRECISION RWORK( * ), W( * )
34* COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> ZHBEVD_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
44*> a complex Hermitian band matrix A using the 2stage technique for
45*> the reduction to tridiagonal. If eigenvectors are desired, it
46*> uses a divide and conquer algorithm.
47*>
48*> \endverbatim
49*
50* Arguments:
51* ==========
52*
53*> \param[in] JOBZ
54*> \verbatim
55*> JOBZ is CHARACTER*1
56*> = 'N': Compute eigenvalues only;
57*> = 'V': Compute eigenvalues and eigenvectors.
58*> Not available in this release.
59*> \endverbatim
60*>
61*> \param[in] UPLO
62*> \verbatim
63*> UPLO is CHARACTER*1
64*> = 'U': Upper triangle of A is stored;
65*> = 'L': Lower triangle of A is stored.
66*> \endverbatim
67*>
68*> \param[in] N
69*> \verbatim
70*> N is INTEGER
71*> The order of the matrix A. N >= 0.
72*> \endverbatim
73*>
74*> \param[in] KD
75*> \verbatim
76*> KD is INTEGER
77*> The number of superdiagonals of the matrix A if UPLO = 'U',
78*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
79*> \endverbatim
80*>
81*> \param[in,out] AB
82*> \verbatim
83*> AB is COMPLEX*16 array, dimension (LDAB, N)
84*> On entry, the upper or lower triangle of the Hermitian band
85*> matrix A, stored in the first KD+1 rows of the array. The
86*> j-th column of A is stored in the j-th column of the array AB
87*> as follows:
88*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
89*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
90*>
91*> On exit, AB is overwritten by values generated during the
92*> reduction to tridiagonal form. If UPLO = 'U', the first
93*> superdiagonal and the diagonal of the tridiagonal matrix T
94*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
95*> the diagonal and first subdiagonal of T are returned in the
96*> first two rows of AB.
97*> \endverbatim
98*>
99*> \param[in] LDAB
100*> \verbatim
101*> LDAB is INTEGER
102*> The leading dimension of the array AB. LDAB >= KD + 1.
103*> \endverbatim
104*>
105*> \param[out] W
106*> \verbatim
107*> W is DOUBLE PRECISION array, dimension (N)
108*> If INFO = 0, the eigenvalues in ascending order.
109*> \endverbatim
110*>
111*> \param[out] Z
112*> \verbatim
113*> Z is COMPLEX*16 array, dimension (LDZ, N)
114*> If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
115*> eigenvectors of the matrix A, with the i-th column of Z
116*> holding the eigenvector associated with W(i).
117*> If JOBZ = 'N', then Z is not referenced.
118*> \endverbatim
119*>
120*> \param[in] LDZ
121*> \verbatim
122*> LDZ is INTEGER
123*> The leading dimension of the array Z. LDZ >= 1, and if
124*> JOBZ = 'V', LDZ >= max(1,N).
125*> \endverbatim
126*>
127*> \param[out] WORK
128*> \verbatim
129*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
130*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
131*> \endverbatim
132*>
133*> \param[in] LWORK
134*> \verbatim
135*> LWORK is INTEGER
136*> The length of the array WORK. LWORK >= 1, when N <= 1;
137*> otherwise
138*> If JOBZ = 'N' and N > 1, LWORK must be queried.
139*> LWORK = MAX(1, dimension) where
140*> dimension = (2KD+1)*N + KD*NTHREADS
141*> where KD is the size of the band.
142*> NTHREADS is the number of threads used when
143*> openMP compilation is enabled, otherwise =1.
144*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.
145*>
146*> If LWORK = -1, then a workspace query is assumed; the routine
147*> only calculates the optimal sizes of the WORK, RWORK and
148*> IWORK arrays, returns these values as the first entries of
149*> the WORK, RWORK and IWORK arrays, and no error message
150*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
151*> \endverbatim
152*>
153*> \param[out] RWORK
154*> \verbatim
155*> RWORK is DOUBLE PRECISION array,
156*> dimension (LRWORK)
157*> On exit, if INFO = 0, RWORK(1) returns the optimal LRWORK.
158*> \endverbatim
159*>
160*> \param[in] LRWORK
161*> \verbatim
162*> LRWORK is INTEGER
163*> The dimension of array RWORK.
164*> If N <= 1, LRWORK must be at least 1.
165*> If JOBZ = 'N' and N > 1, LRWORK must be at least N.
166*> If JOBZ = 'V' and N > 1, LRWORK must be at least
167*> 1 + 5*N + 2*N**2.
168*>
169*> If LRWORK = -1, then a workspace query is assumed; the
170*> routine only calculates the optimal sizes of the WORK, RWORK
171*> and IWORK arrays, returns these values as the first entries
172*> of the WORK, RWORK and IWORK arrays, and no error message
173*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
174*> \endverbatim
175*>
176*> \param[out] IWORK
177*> \verbatim
178*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
179*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
180*> \endverbatim
181*>
182*> \param[in] LIWORK
183*> \verbatim
184*> LIWORK is INTEGER
185*> The dimension of array IWORK.
186*> If JOBZ = 'N' or N <= 1, LIWORK must be at least 1.
187*> If JOBZ = 'V' and N > 1, LIWORK must be at least 3 + 5*N .
188*>
189*> If LIWORK = -1, then a workspace query is assumed; the
190*> routine only calculates the optimal sizes of the WORK, RWORK
191*> and IWORK arrays, returns these values as the first entries
192*> of the WORK, RWORK and IWORK arrays, and no error message
193*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
194*> \endverbatim
195*>
196*> \param[out] INFO
197*> \verbatim
198*> INFO is INTEGER
199*> = 0: successful exit.
200*> < 0: if INFO = -i, the i-th argument had an illegal value.
201*> > 0: if INFO = i, the algorithm failed to converge; i
202*> off-diagonal elements of an intermediate tridiagonal
203*> form did not converge to zero.
204*> \endverbatim
205*
206* Authors:
207* ========
208*
209*> \author Univ. of Tennessee
210*> \author Univ. of California Berkeley
211*> \author Univ. of Colorado Denver
212*> \author NAG Ltd.
213*
214*> \ingroup hbevd_2stage
215*
216*> \par Further Details:
217* =====================
218*>
219*> \verbatim
220*>
221*> All details about the 2stage techniques are available in:
222*>
223*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
224*> Parallel reduction to condensed forms for symmetric eigenvalue problems
225*> using aggregated fine-grained and memory-aware kernels. In Proceedings
226*> of 2011 International Conference for High Performance Computing,
227*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
228*> Article 8 , 11 pages.
229*> http://doi.acm.org/10.1145/2063384.2063394
230*>
231*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
232*> An improved parallel singular value algorithm and its implementation
233*> for multicore hardware, In Proceedings of 2013 International Conference
234*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
235*> Denver, Colorado, USA, 2013.
236*> Article 90, 12 pages.
237*> http://doi.acm.org/10.1145/2503210.2503292
238*>
239*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
240*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
241*> calculations based on fine-grained memory aware tasks.
242*> International Journal of High Performance Computing Applications.
243*> Volume 28 Issue 2, Pages 196-209, May 2014.
244*> http://hpc.sagepub.com/content/28/2/196
245*>
246*> \endverbatim
247*
248* =====================================================================
249 SUBROUTINE zhbevd_2stage( JOBZ, UPLO, N, KD, AB, LDAB, W, Z,
250 $ LDZ,
251 $ WORK, LWORK, RWORK, LRWORK, IWORK,
252 $ LIWORK, INFO )
253*
254 IMPLICIT NONE
255*
256* -- LAPACK driver routine --
257* -- LAPACK is a software package provided by Univ. of Tennessee, --
258* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
259*
260* .. Scalar Arguments ..
261 CHARACTER JOBZ, UPLO
262 INTEGER INFO, KD, LDAB, LDZ, LIWORK, LRWORK, LWORK, N
263* ..
264* .. Array Arguments ..
265 INTEGER IWORK( * )
266 DOUBLE PRECISION RWORK( * ), W( * )
267 COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
268* ..
269*
270* =====================================================================
271*
272* .. Parameters ..
273 DOUBLE PRECISION ZERO, ONE
274 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
275 COMPLEX*16 CZERO, CONE
276 parameter( czero = ( 0.0d0, 0.0d0 ),
277 $ cone = ( 1.0d0, 0.0d0 ) )
278* ..
279* .. Local Scalars ..
280 LOGICAL LOWER, LQUERY, WANTZ
281 INTEGER IINFO, IMAX, INDE, INDWK2, INDRWK, ISCALE,
282 $ LLWORK, INDWK, LHTRD, LWTRD, IB, INDHOUS,
283 $ liwmin, llrwk, llwk2, lrwmin, lwmin
284 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
285 $ SMLNUM
286* ..
287* .. External Functions ..
288 LOGICAL LSAME
289 INTEGER ILAENV2STAGE
290 DOUBLE PRECISION DLAMCH, ZLANHB
291 EXTERNAL lsame, dlamch, zlanhb, ilaenv2stage
292* ..
293* .. External Subroutines ..
294 EXTERNAL dscal, dsterf, xerbla, zgemm,
295 $ zlacpy,
297* ..
298* .. Intrinsic Functions ..
299 INTRINSIC dble, sqrt
300* ..
301* .. Executable Statements ..
302*
303* Test the input parameters.
304*
305 wantz = lsame( jobz, 'V' )
306 lower = lsame( uplo, 'L' )
307 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 .OR. lrwork.EQ.-1 )
308*
309 info = 0
310 IF( n.LE.1 ) THEN
311 lwmin = 1
312 lrwmin = 1
313 liwmin = 1
314 ELSE
315 ib = ilaenv2stage( 2, 'ZHETRD_HB2ST', jobz, n, kd, -1,
316 $ -1 )
317 lhtrd = ilaenv2stage( 3, 'ZHETRD_HB2ST', jobz, n, kd, ib,
318 $ -1 )
319 lwtrd = ilaenv2stage( 4, 'ZHETRD_HB2ST', jobz, n, kd, ib,
320 $ -1 )
321 IF( wantz ) THEN
322 lwmin = 2*n**2
323 lrwmin = 1 + 5*n + 2*n**2
324 liwmin = 3 + 5*n
325 ELSE
326 lwmin = max( n, lhtrd + lwtrd )
327 lrwmin = n
328 liwmin = 1
329 END IF
330 END IF
331 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
332 info = -1
333 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
334 info = -2
335 ELSE IF( n.LT.0 ) THEN
336 info = -3
337 ELSE IF( kd.LT.0 ) THEN
338 info = -4
339 ELSE IF( ldab.LT.kd+1 ) THEN
340 info = -6
341 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
342 info = -9
343 END IF
344*
345 IF( info.EQ.0 ) THEN
346 work( 1 ) = lwmin
347 rwork( 1 ) = real( lrwmin )
348 iwork( 1 ) = liwmin
349*
350 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
351 info = -11
352 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
353 info = -13
354 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
355 info = -15
356 END IF
357 END IF
358*
359 IF( info.NE.0 ) THEN
360 CALL xerbla( 'ZHBEVD_2STAGE', -info )
361 RETURN
362 ELSE IF( lquery ) THEN
363 RETURN
364 END IF
365*
366* Quick return if possible
367*
368 IF( n.EQ.0 )
369 $ RETURN
370*
371 IF( n.EQ.1 ) THEN
372 w( 1 ) = dble( ab( 1, 1 ) )
373 IF( wantz )
374 $ z( 1, 1 ) = cone
375 RETURN
376 END IF
377*
378* Get machine constants.
379*
380 safmin = dlamch( 'Safe minimum' )
381 eps = dlamch( 'Precision' )
382 smlnum = safmin / eps
383 bignum = one / smlnum
384 rmin = sqrt( smlnum )
385 rmax = sqrt( bignum )
386*
387* Scale matrix to allowable range, if necessary.
388*
389 anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
390 iscale = 0
391 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
392 iscale = 1
393 sigma = rmin / anrm
394 ELSE IF( anrm.GT.rmax ) THEN
395 iscale = 1
396 sigma = rmax / anrm
397 END IF
398 IF( iscale.EQ.1 ) THEN
399 IF( lower ) THEN
400 CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab,
401 $ info )
402 ELSE
403 CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab,
404 $ info )
405 END IF
406 END IF
407*
408* Call ZHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
409*
410 inde = 1
411 indrwk = inde + n
412 llrwk = lrwork - indrwk + 1
413 indhous = 1
414 indwk = indhous + lhtrd
415 llwork = lwork - indwk + 1
416 indwk2 = indwk + n*n
417 llwk2 = lwork - indwk2 + 1
418*
419 CALL zhetrd_hb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
420 $ rwork( inde ), work( indhous ), lhtrd,
421 $ work( indwk ), llwork, iinfo )
422*
423* For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
424*
425 IF( .NOT.wantz ) THEN
426 CALL dsterf( n, w, rwork( inde ), info )
427 ELSE
428 CALL zstedc( 'I', n, w, rwork( inde ), work, n,
429 $ work( indwk2 ),
430 $ llwk2, rwork( indrwk ), llrwk, iwork, liwork,
431 $ info )
432 CALL zgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
433 $ work( indwk2 ), n )
434 CALL zlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
435 END IF
436*
437* If matrix was scaled, then rescale eigenvalues appropriately.
438*
439 IF( iscale.EQ.1 ) THEN
440 IF( info.EQ.0 ) THEN
441 imax = n
442 ELSE
443 imax = info - 1
444 END IF
445 CALL dscal( imax, one / sigma, w, 1 )
446 END IF
447*
448 work( 1 ) = lwmin
449 rwork( 1 ) = real( lrwmin )
450 iwork( 1 ) = liwmin
451 RETURN
452*
453* End of ZHBEVD_2STAGE
454*
455 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zhbevd_2stage(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine zhetrd_hb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
ZHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:142
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zstedc(compz, n, d, e, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
ZSTEDC
Definition zstedc.f:204
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84