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