LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlasda ( integer  ICOMPQ,
integer  SMLSIZ,
integer  N,
integer  SQRE,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
double precision, dimension( ldu, * )  U,
integer  LDU,
double precision, dimension( ldu, * )  VT,
integer, dimension( * )  K,
double precision, dimension( ldu, * )  DIFL,
double precision, dimension( ldu, * )  DIFR,
double precision, dimension( ldu, * )  Z,
double precision, dimension( ldu, * )  POLES,
integer, dimension( * )  GIVPTR,
integer, dimension( ldgcol, * )  GIVCOL,
integer  LDGCOL,
integer, dimension( ldgcol, * )  PERM,
double precision, dimension( ldu, * )  GIVNUM,
double precision, dimension( * )  C,
double precision, dimension( * )  S,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

DLASDA computes the singular value decomposition (SVD) of a real upper bidiagonal matrix with diagonal d and off-diagonal e. Used by sbdsdc.

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

Purpose:
 Using a divide and conquer approach, DLASDA computes the singular
 value decomposition (SVD) of a real upper bidiagonal N-by-M matrix
 B with diagonal D and offdiagonal E, where M = N + SQRE. The
 algorithm computes the singular values in the SVD B = U * S * VT.
 The orthogonal matrices U and VT are optionally computed in
 compact form.

 A related subroutine, DLASD0, computes the singular values and
 the singular vectors in explicit form.
Parameters
[in]ICOMPQ
          ICOMPQ is INTEGER
         Specifies whether singular vectors are to be computed
         in compact form, as follows
         = 0: Compute singular values only.
         = 1: Compute singular vectors of upper bidiagonal
              matrix in compact form.
[in]SMLSIZ
          SMLSIZ is INTEGER
         The maximum size of the subproblems at the bottom of the
         computation tree.
[in]N
          N is INTEGER
         The row dimension of the upper bidiagonal matrix. This is
         also the dimension of the main diagonal array D.
[in]SQRE
          SQRE is INTEGER
         Specifies the column dimension of the bidiagonal matrix.
         = 0: The bidiagonal matrix has column dimension M = N;
         = 1: The bidiagonal matrix has column dimension M = N + 1.
[in,out]D
          D is DOUBLE PRECISION array, dimension ( N )
         On entry D contains the main diagonal of the bidiagonal
         matrix. On exit D, if INFO = 0, contains its singular values.
[in]E
          E is DOUBLE PRECISION array, dimension ( M-1 )
         Contains the subdiagonal entries of the bidiagonal matrix.
         On exit, E has been destroyed.
[out]U
          U is DOUBLE PRECISION array,
         dimension ( LDU, SMLSIZ ) if ICOMPQ = 1, and not referenced
         if ICOMPQ = 0. If ICOMPQ = 1, on exit, U contains the left
         singular vector matrices of all subproblems at the bottom
         level.
[in]LDU
          LDU is INTEGER, LDU = > N.
         The leading dimension of arrays U, VT, DIFL, DIFR, POLES,
         GIVNUM, and Z.
[out]VT
          VT is DOUBLE PRECISION array,
         dimension ( LDU, SMLSIZ+1 ) if ICOMPQ = 1, and not referenced
         if ICOMPQ = 0. If ICOMPQ = 1, on exit, VT**T contains the right
         singular vector matrices of all subproblems at the bottom
         level.
[out]K
          K is INTEGER array,
         dimension ( N ) if ICOMPQ = 1 and dimension 1 if ICOMPQ = 0.
         If ICOMPQ = 1, on exit, K(I) is the dimension of the I-th
         secular equation on the computation tree.
[out]DIFL
          DIFL is DOUBLE PRECISION array, dimension ( LDU, NLVL ),
         where NLVL = floor(log_2 (N/SMLSIZ))).
[out]DIFR
          DIFR is DOUBLE PRECISION array,
                  dimension ( LDU, 2 * NLVL ) if ICOMPQ = 1 and
                  dimension ( N ) if ICOMPQ = 0.
         If ICOMPQ = 1, on exit, DIFL(1:N, I) and DIFR(1:N, 2 * I - 1)
         record distances between singular values on the I-th
         level and singular values on the (I -1)-th level, and
         DIFR(1:N, 2 * I ) contains the normalizing factors for
         the right singular vector matrix. See DLASD8 for details.
[out]Z
          Z is DOUBLE PRECISION array,
                  dimension ( LDU, NLVL ) if ICOMPQ = 1 and
                  dimension ( N ) if ICOMPQ = 0.
         The first K elements of Z(1, I) contain the components of
         the deflation-adjusted updating row vector for subproblems
         on the I-th level.
[out]POLES
          POLES is DOUBLE PRECISION array,
         dimension ( LDU, 2 * NLVL ) if ICOMPQ = 1, and not referenced
         if ICOMPQ = 0. If ICOMPQ = 1, on exit, POLES(1, 2*I - 1) and
         POLES(1, 2*I) contain  the new and old singular values
         involved in the secular equations on the I-th level.
[out]GIVPTR
          GIVPTR is INTEGER array,
         dimension ( N ) if ICOMPQ = 1, and not referenced if
         ICOMPQ = 0. If ICOMPQ = 1, on exit, GIVPTR( I ) records
         the number of Givens rotations performed on the I-th
         problem on the computation tree.
[out]GIVCOL
          GIVCOL is INTEGER array,
         dimension ( LDGCOL, 2 * NLVL ) if ICOMPQ = 1, and not
         referenced if ICOMPQ = 0. If ICOMPQ = 1, on exit, for each I,
         GIVCOL(1, 2 *I - 1) and GIVCOL(1, 2 *I) record the locations
         of Givens rotations performed on the I-th level on the
         computation tree.
[in]LDGCOL
          LDGCOL is INTEGER, LDGCOL = > N.
         The leading dimension of arrays GIVCOL and PERM.
[out]PERM
          PERM is INTEGER array,
         dimension ( LDGCOL, NLVL ) if ICOMPQ = 1, and not referenced
         if ICOMPQ = 0. If ICOMPQ = 1, on exit, PERM(1, I) records
         permutations done on the I-th level of the computation tree.
[out]GIVNUM
          GIVNUM is DOUBLE PRECISION array,
         dimension ( LDU,  2 * NLVL ) if ICOMPQ = 1, and not
         referenced if ICOMPQ = 0. If ICOMPQ = 1, on exit, for each I,
         GIVNUM(1, 2 *I - 1) and GIVNUM(1, 2 *I) record the C- and S-
         values of Givens rotations performed on the I-th level on
         the computation tree.
[out]C
          C is DOUBLE PRECISION array,
         dimension ( N ) if ICOMPQ = 1, and dimension 1 if ICOMPQ = 0.
         If ICOMPQ = 1 and the I-th subproblem is not square, on exit,
         C( I ) contains the C-value of a Givens rotation related to
         the right null space of the I-th subproblem.
[out]S
          S is DOUBLE PRECISION array, dimension ( N ) if
         ICOMPQ = 1, and dimension 1 if ICOMPQ = 0. If ICOMPQ = 1
         and the I-th subproblem is not square, on exit, S( I )
         contains the S-value of a Givens rotation related to
         the right null space of the I-th subproblem.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
         (6 * N + (SMLSIZ + 1)*(SMLSIZ + 1)).
[out]IWORK
          IWORK is INTEGER array.
         Dimension must be at least (7 * N).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = 1, a singular value did not converge
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Contributors:
Ming Gu and Huan Ren, Computer Science Division, University of California at Berkeley, USA

Definition at line 276 of file dlasda.f.

276 *
277 * -- LAPACK auxiliary routine (version 3.4.2) --
278 * -- LAPACK is a software package provided by Univ. of Tennessee, --
279 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
280 * September 2012
281 *
282 * .. Scalar Arguments ..
283  INTEGER icompq, info, ldgcol, ldu, n, smlsiz, sqre
284 * ..
285 * .. Array Arguments ..
286  INTEGER givcol( ldgcol, * ), givptr( * ), iwork( * ),
287  $ k( * ), perm( ldgcol, * )
288  DOUBLE PRECISION c( * ), d( * ), difl( ldu, * ), difr( ldu, * ),
289  $ e( * ), givnum( ldu, * ), poles( ldu, * ),
290  $ s( * ), u( ldu, * ), vt( ldu, * ), work( * ),
291  $ z( ldu, * )
292 * ..
293 *
294 * =====================================================================
295 *
296 * .. Parameters ..
297  DOUBLE PRECISION zero, one
298  parameter ( zero = 0.0d+0, one = 1.0d+0 )
299 * ..
300 * .. Local Scalars ..
301  INTEGER i, i1, ic, idxq, idxqi, im1, inode, itemp, iwk,
302  $ j, lf, ll, lvl, lvl2, m, ncc, nd, ndb1, ndiml,
303  $ ndimr, nl, nlf, nlp1, nlvl, nr, nrf, nrp1, nru,
304  $ nwork1, nwork2, smlszp, sqrei, vf, vfi, vl, vli
305  DOUBLE PRECISION alpha, beta
306 * ..
307 * .. External Subroutines ..
308  EXTERNAL dcopy, dlasd6, dlasdq, dlasdt, dlaset, xerbla
309 * ..
310 * .. Executable Statements ..
311 *
312 * Test the input parameters.
313 *
314  info = 0
315 *
316  IF( ( icompq.LT.0 ) .OR. ( icompq.GT.1 ) ) THEN
317  info = -1
318  ELSE IF( smlsiz.LT.3 ) THEN
319  info = -2
320  ELSE IF( n.LT.0 ) THEN
321  info = -3
322  ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
323  info = -4
324  ELSE IF( ldu.LT.( n+sqre ) ) THEN
325  info = -8
326  ELSE IF( ldgcol.LT.n ) THEN
327  info = -17
328  END IF
329  IF( info.NE.0 ) THEN
330  CALL xerbla( 'DLASDA', -info )
331  RETURN
332  END IF
333 *
334  m = n + sqre
335 *
336 * If the input matrix is too small, call DLASDQ to find the SVD.
337 *
338  IF( n.LE.smlsiz ) THEN
339  IF( icompq.EQ.0 ) THEN
340  CALL dlasdq( 'U', sqre, n, 0, 0, 0, d, e, vt, ldu, u, ldu,
341  $ u, ldu, work, info )
342  ELSE
343  CALL dlasdq( 'U', sqre, n, m, n, 0, d, e, vt, ldu, u, ldu,
344  $ u, ldu, work, info )
345  END IF
346  RETURN
347  END IF
348 *
349 * Book-keeping and set up the computation tree.
350 *
351  inode = 1
352  ndiml = inode + n
353  ndimr = ndiml + n
354  idxq = ndimr + n
355  iwk = idxq + n
356 *
357  ncc = 0
358  nru = 0
359 *
360  smlszp = smlsiz + 1
361  vf = 1
362  vl = vf + m
363  nwork1 = vl + m
364  nwork2 = nwork1 + smlszp*smlszp
365 *
366  CALL dlasdt( n, nlvl, nd, iwork( inode ), iwork( ndiml ),
367  $ iwork( ndimr ), smlsiz )
368 *
369 * for the nodes on bottom level of the tree, solve
370 * their subproblems by DLASDQ.
371 *
372  ndb1 = ( nd+1 ) / 2
373  DO 30 i = ndb1, nd
374 *
375 * IC : center row of each node
376 * NL : number of rows of left subproblem
377 * NR : number of rows of right subproblem
378 * NLF: starting row of the left subproblem
379 * NRF: starting row of the right subproblem
380 *
381  i1 = i - 1
382  ic = iwork( inode+i1 )
383  nl = iwork( ndiml+i1 )
384  nlp1 = nl + 1
385  nr = iwork( ndimr+i1 )
386  nlf = ic - nl
387  nrf = ic + 1
388  idxqi = idxq + nlf - 2
389  vfi = vf + nlf - 1
390  vli = vl + nlf - 1
391  sqrei = 1
392  IF( icompq.EQ.0 ) THEN
393  CALL dlaset( 'A', nlp1, nlp1, zero, one, work( nwork1 ),
394  $ smlszp )
395  CALL dlasdq( 'U', sqrei, nl, nlp1, nru, ncc, d( nlf ),
396  $ e( nlf ), work( nwork1 ), smlszp,
397  $ work( nwork2 ), nl, work( nwork2 ), nl,
398  $ work( nwork2 ), info )
399  itemp = nwork1 + nl*smlszp
400  CALL dcopy( nlp1, work( nwork1 ), 1, work( vfi ), 1 )
401  CALL dcopy( nlp1, work( itemp ), 1, work( vli ), 1 )
402  ELSE
403  CALL dlaset( 'A', nl, nl, zero, one, u( nlf, 1 ), ldu )
404  CALL dlaset( 'A', nlp1, nlp1, zero, one, vt( nlf, 1 ), ldu )
405  CALL dlasdq( 'U', sqrei, nl, nlp1, nl, ncc, d( nlf ),
406  $ e( nlf ), vt( nlf, 1 ), ldu, u( nlf, 1 ), ldu,
407  $ u( nlf, 1 ), ldu, work( nwork1 ), info )
408  CALL dcopy( nlp1, vt( nlf, 1 ), 1, work( vfi ), 1 )
409  CALL dcopy( nlp1, vt( nlf, nlp1 ), 1, work( vli ), 1 )
410  END IF
411  IF( info.NE.0 ) THEN
412  RETURN
413  END IF
414  DO 10 j = 1, nl
415  iwork( idxqi+j ) = j
416  10 CONTINUE
417  IF( ( i.EQ.nd ) .AND. ( sqre.EQ.0 ) ) THEN
418  sqrei = 0
419  ELSE
420  sqrei = 1
421  END IF
422  idxqi = idxqi + nlp1
423  vfi = vfi + nlp1
424  vli = vli + nlp1
425  nrp1 = nr + sqrei
426  IF( icompq.EQ.0 ) THEN
427  CALL dlaset( 'A', nrp1, nrp1, zero, one, work( nwork1 ),
428  $ smlszp )
429  CALL dlasdq( 'U', sqrei, nr, nrp1, nru, ncc, d( nrf ),
430  $ e( nrf ), work( nwork1 ), smlszp,
431  $ work( nwork2 ), nr, work( nwork2 ), nr,
432  $ work( nwork2 ), info )
433  itemp = nwork1 + ( nrp1-1 )*smlszp
434  CALL dcopy( nrp1, work( nwork1 ), 1, work( vfi ), 1 )
435  CALL dcopy( nrp1, work( itemp ), 1, work( vli ), 1 )
436  ELSE
437  CALL dlaset( 'A', nr, nr, zero, one, u( nrf, 1 ), ldu )
438  CALL dlaset( 'A', nrp1, nrp1, zero, one, vt( nrf, 1 ), ldu )
439  CALL dlasdq( 'U', sqrei, nr, nrp1, nr, ncc, d( nrf ),
440  $ e( nrf ), vt( nrf, 1 ), ldu, u( nrf, 1 ), ldu,
441  $ u( nrf, 1 ), ldu, work( nwork1 ), info )
442  CALL dcopy( nrp1, vt( nrf, 1 ), 1, work( vfi ), 1 )
443  CALL dcopy( nrp1, vt( nrf, nrp1 ), 1, work( vli ), 1 )
444  END IF
445  IF( info.NE.0 ) THEN
446  RETURN
447  END IF
448  DO 20 j = 1, nr
449  iwork( idxqi+j ) = j
450  20 CONTINUE
451  30 CONTINUE
452 *
453 * Now conquer each subproblem bottom-up.
454 *
455  j = 2**nlvl
456  DO 50 lvl = nlvl, 1, -1
457  lvl2 = lvl*2 - 1
458 *
459 * Find the first node LF and last node LL on
460 * the current level LVL.
461 *
462  IF( lvl.EQ.1 ) THEN
463  lf = 1
464  ll = 1
465  ELSE
466  lf = 2**( lvl-1 )
467  ll = 2*lf - 1
468  END IF
469  DO 40 i = lf, ll
470  im1 = i - 1
471  ic = iwork( inode+im1 )
472  nl = iwork( ndiml+im1 )
473  nr = iwork( ndimr+im1 )
474  nlf = ic - nl
475  nrf = ic + 1
476  IF( i.EQ.ll ) THEN
477  sqrei = sqre
478  ELSE
479  sqrei = 1
480  END IF
481  vfi = vf + nlf - 1
482  vli = vl + nlf - 1
483  idxqi = idxq + nlf - 1
484  alpha = d( ic )
485  beta = e( ic )
486  IF( icompq.EQ.0 ) THEN
487  CALL dlasd6( icompq, nl, nr, sqrei, d( nlf ),
488  $ work( vfi ), work( vli ), alpha, beta,
489  $ iwork( idxqi ), perm, givptr( 1 ), givcol,
490  $ ldgcol, givnum, ldu, poles, difl, difr, z,
491  $ k( 1 ), c( 1 ), s( 1 ), work( nwork1 ),
492  $ iwork( iwk ), info )
493  ELSE
494  j = j - 1
495  CALL dlasd6( icompq, nl, nr, sqrei, d( nlf ),
496  $ work( vfi ), work( vli ), alpha, beta,
497  $ iwork( idxqi ), perm( nlf, lvl ),
498  $ givptr( j ), givcol( nlf, lvl2 ), ldgcol,
499  $ givnum( nlf, lvl2 ), ldu,
500  $ poles( nlf, lvl2 ), difl( nlf, lvl ),
501  $ difr( nlf, lvl2 ), z( nlf, lvl ), k( j ),
502  $ c( j ), s( j ), work( nwork1 ),
503  $ iwork( iwk ), info )
504  END IF
505  IF( info.NE.0 ) THEN
506  RETURN
507  END IF
508  40 CONTINUE
509  50 CONTINUE
510 *
511  RETURN
512 *
513 * End of DLASDA
514 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dlasdt(N, LVL, ND, INODE, NDIML, NDIMR, MSUB)
DLASDT creates a tree of subproblems for bidiagonal divide and conquer. Used by sbdsdc.
Definition: dlasdt.f:107
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlasdq(UPLO, SQRE, N, NCVT, NRU, NCC, D, E, VT, LDVT, U, LDU, C, LDC, WORK, INFO)
DLASDQ computes the SVD of a real bidiagonal matrix with diagonal d and off-diagonal e...
Definition: dlasdq.f:213
subroutine dlasd6(ICOMPQ, NL, NR, SQRE, D, VF, VL, ALPHA, BETA, IDXQ, PERM, GIVPTR, GIVCOL, LDGCOL, GIVNUM, LDGNUM, POLES, DIFL, DIFR, Z, K, C, S, WORK, IWORK, INFO)
DLASD6 computes the SVD of an updated upper bidiagonal matrix obtained by merging two smaller ones by...
Definition: dlasd6.f:315

Here is the call graph for this function:

Here is the caller graph for this function: