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

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

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

Purpose:
 Using a divide and conquer approach, SLASDA 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, SLASD0, 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 REAL 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 REAL array, dimension ( M-1 )
         Contains the subdiagonal entries of the bidiagonal matrix.
         On exit, E has been destroyed.
[out]U
          U is REAL 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 REAL 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 REAL array, dimension ( LDU, NLVL ),
         where NLVL = floor(log_2 (N/SMLSIZ))).
[out]DIFR
          DIFR is REAL 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 SLASD8 for details.
[out]Z
          Z is REAL 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 REAL 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 REAL 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 REAL 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 REAL 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 REAL array, dimension
         (6 * N + (SMLSIZ + 1)*(SMLSIZ + 1)).
[out]IWORK
          IWORK is INTEGER array, dimension (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 275 of file slasda.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: