LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dsbevd_2stage()

subroutine dsbevd_2stage ( character  jobz,
character  uplo,
integer  n,
integer  kd,
double precision, dimension( ldab, * )  ab,
integer  ldab,
double precision, dimension( * )  w,
double precision, dimension( ldz, * )  z,
integer  ldz,
double precision, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

DSBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

Download DSBEVD_2STAGE + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 DSBEVD_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
 a real symmetric band matrix A using the 2stage technique for
 the reduction to tridiagonal. If eigenvectors are desired, it uses
 a divide and conquer algorithm.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
                  Not available in this release.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 0.
[in,out]AB
          AB is DOUBLE PRECISION array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the symmetric band
          matrix A, stored in the first KD+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).

          On exit, AB is overwritten by values generated during the
          reduction to tridiagonal form.  If UPLO = 'U', the first
          superdiagonal and the diagonal of the tridiagonal matrix T
          are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
          the diagonal and first subdiagonal of T are returned in the
          first two rows of AB.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD + 1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is DOUBLE PRECISION array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
          eigenvectors of the matrix A, with the i-th column of Z
          holding the eigenvector associated with W(i).
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension LWORK
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK. LWORK >= 1, when N <= 1;
          otherwise
          If JOBZ = 'N' and N > 1, LWORK must be queried.
                                   LWORK = MAX(1, dimension) where
                                   dimension = (2KD+1)*N + KD*NTHREADS + N
                                   where KD is the size of the band.
                                   NTHREADS is the number of threads used when
                                   openMP compilation is enabled, otherwise =1.
          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If JOBZ  = 'N' or N <= 1, LIWORK must be at least 1.
          If JOBZ  = 'V' and N > 2, LIWORK must be at least 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the algorithm failed to converge; i
                off-diagonal elements of an intermediate tridiagonal
                form did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  All details about the 2stage techniques are available in:

  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
  Parallel reduction to condensed forms for symmetric eigenvalue problems
  using aggregated fine-grained and memory-aware kernels. In Proceedings
  of 2011 International Conference for High Performance Computing,
  Networking, Storage and Analysis (SC '11), New York, NY, USA,
  Article 8 , 11 pages.
  http://doi.acm.org/10.1145/2063384.2063394

  A. Haidar, J. Kurzak, P. Luszczek, 2013.
  An improved parallel singular value algorithm and its implementation
  for multicore hardware, In Proceedings of 2013 International Conference
  for High Performance Computing, Networking, Storage and Analysis (SC '13).
  Denver, Colorado, USA, 2013.
  Article 90, 12 pages.
  http://doi.acm.org/10.1145/2503210.2503292

  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
  A novel hybrid CPU-GPU generalized eigensolver for electronic structure
  calculations based on fine-grained memory aware tasks.
  International Journal of High Performance Computing Applications.
  Volume 28 Issue 2, Pages 196-209, May 2014.
  http://hpc.sagepub.com/content/28/2/196

Definition at line 226 of file dsbevd_2stage.f.

228*
229 IMPLICIT NONE
230*
231* -- LAPACK driver routine --
232* -- LAPACK is a software package provided by Univ. of Tennessee, --
233* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
234*
235* .. Scalar Arguments ..
236 CHARACTER JOBZ, UPLO
237 INTEGER INFO, KD, LDAB, LDZ, LIWORK, LWORK, N
238* ..
239* .. Array Arguments ..
240 INTEGER IWORK( * )
241 DOUBLE PRECISION AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
242* ..
243*
244* =====================================================================
245*
246* .. Parameters ..
247 DOUBLE PRECISION ZERO, ONE
248 parameter( zero = 0.0d+0, one = 1.0d+0 )
249* ..
250* .. Local Scalars ..
251 LOGICAL LOWER, LQUERY, WANTZ
252 INTEGER IINFO, INDE, INDWK2, INDWRK, ISCALE, LIWMIN,
253 $ LLWORK, LWMIN, LHTRD, LWTRD, IB, INDHOUS,
254 $ LLWRK2
255 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
256 $ SMLNUM
257* ..
258* .. External Functions ..
259 LOGICAL LSAME
260 INTEGER ILAENV2STAGE
261 DOUBLE PRECISION DLAMCH, DLANSB
262 EXTERNAL lsame, dlamch, dlansb, ilaenv2stage
263* ..
264* .. External Subroutines ..
265 EXTERNAL dgemm, dlacpy, dlascl, dscal, 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, -1 )
285 lhtrd = ilaenv2stage( 3, 'DSYTRD_SB2ST', jobz, n, kd, ib, -1 )
286 lwtrd = ilaenv2stage( 4, 'DSYTRD_SB2ST', jobz, n, kd, ib, -1 )
287 IF( wantz ) THEN
288 liwmin = 3 + 5*n
289 lwmin = 1 + 5*n + 2*n**2
290 ELSE
291 liwmin = 1
292 lwmin = max( 2*n, n+lhtrd+lwtrd )
293 END IF
294 END IF
295 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
296 info = -1
297 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
298 info = -2
299 ELSE IF( n.LT.0 ) THEN
300 info = -3
301 ELSE IF( kd.LT.0 ) THEN
302 info = -4
303 ELSE IF( ldab.LT.kd+1 ) THEN
304 info = -6
305 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
306 info = -9
307 END IF
308*
309 IF( info.EQ.0 ) THEN
310 work( 1 ) = lwmin
311 iwork( 1 ) = liwmin
312*
313 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
314 info = -11
315 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
316 info = -13
317 END IF
318 END IF
319*
320 IF( info.NE.0 ) THEN
321 CALL xerbla( 'DSBEVD_2STAGE', -info )
322 RETURN
323 ELSE IF( lquery ) THEN
324 RETURN
325 END IF
326*
327* Quick return if possible
328*
329 IF( n.EQ.0 )
330 $ RETURN
331*
332 IF( n.EQ.1 ) THEN
333 w( 1 ) = ab( 1, 1 )
334 IF( wantz )
335 $ z( 1, 1 ) = one
336 RETURN
337 END IF
338*
339* Get machine constants.
340*
341 safmin = dlamch( 'Safe minimum' )
342 eps = dlamch( 'Precision' )
343 smlnum = safmin / eps
344 bignum = one / smlnum
345 rmin = sqrt( smlnum )
346 rmax = sqrt( bignum )
347*
348* Scale matrix to allowable range, if necessary.
349*
350 anrm = dlansb( 'M', uplo, n, kd, ab, ldab, work )
351 iscale = 0
352 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
353 iscale = 1
354 sigma = rmin / anrm
355 ELSE IF( anrm.GT.rmax ) THEN
356 iscale = 1
357 sigma = rmax / anrm
358 END IF
359 IF( iscale.EQ.1 ) THEN
360 IF( lower ) THEN
361 CALL dlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
362 ELSE
363 CALL dlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
364 END IF
365 END IF
366*
367* Call DSYTRD_SB2ST to reduce band symmetric matrix to tridiagonal form.
368*
369 inde = 1
370 indhous = inde + n
371 indwrk = indhous + lhtrd
372 llwork = lwork - indwrk + 1
373 indwk2 = indwrk + n*n
374 llwrk2 = lwork - indwk2 + 1
375*
376 CALL dsytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
377 $ work( inde ), work( indhous ), lhtrd,
378 $ work( indwrk ), llwork, iinfo )
379*
380* For eigenvalues only, call DSTERF. For eigenvectors, call SSTEDC.
381*
382 IF( .NOT.wantz ) THEN
383 CALL dsterf( n, w, work( inde ), info )
384 ELSE
385 CALL dstedc( 'I', n, w, work( inde ), work( indwrk ), n,
386 $ work( indwk2 ), llwrk2, iwork, liwork, info )
387 CALL dgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ), n,
388 $ zero, work( indwk2 ), n )
389 CALL dlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
390 END IF
391*
392* If matrix was scaled, then rescale eigenvalues appropriately.
393*
394 IF( iscale.EQ.1 )
395 $ CALL dscal( n, one / sigma, w, 1 )
396*
397 work( 1 ) = lwmin
398 iwork( 1 ) = liwmin
399 RETURN
400*
401* End of DSBEVD_2STAGE
402*
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 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
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlansb(norm, uplo, n, k, ab, ldab, work)
DLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansb.f:129
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:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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:182
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
Here is the call graph for this function:
Here is the caller graph for this function: