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

◆ dlalsa()

subroutine dlalsa ( integer  icompq,
integer  smlsiz,
integer  n,
integer  nrhs,
double precision, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( ldbx, * )  bx,
integer  ldbx,
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 
)

DLALSA computes the SVD of the coefficient matrix in compact form. Used by sgelsd.

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

Purpose:
 DLALSA is an intermediate step in solving the least squares problem
 by computing the SVD of the coefficient matrix in compact form (The
 singular vectors are computed as products of simple orthogonal
 matrices.).

 If ICOMPQ = 0, DLALSA applies the inverse of the left singular vector
 matrix of an upper bidiagonal matrix to the right hand side; and if
 ICOMPQ = 1, DLALSA applies the right singular vector matrix to the
 right hand side. The singular vector matrices were generated in
 compact form by DLALSA.
Parameters
[in]ICOMPQ
          ICOMPQ is INTEGER
         Specifies whether the left or the right singular vector
         matrix is involved.
         = 0: Left singular vector matrix
         = 1: Right singular vector matrix
[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 and column dimensions of the upper bidiagonal matrix.
[in]NRHS
          NRHS is INTEGER
         The number of columns of B and BX. NRHS must be at least 1.
[in,out]B
          B is DOUBLE PRECISION array, dimension ( LDB, NRHS )
         On input, B contains the right hand sides of the least
         squares problem in rows 1 through M.
         On output, B contains the solution X in rows 1 through N.
[in]LDB
          LDB is INTEGER
         The leading dimension of B in the calling subprogram.
         LDB must be at least max(1,MAX( M, N ) ).
[out]BX
          BX is DOUBLE PRECISION array, dimension ( LDBX, NRHS )
         On exit, the result of applying the left or right singular
         vector matrix to B.
[in]LDBX
          LDBX is INTEGER
         The leading dimension of BX.
[in]U
          U is DOUBLE PRECISION array, dimension ( LDU, SMLSIZ ).
         On entry, 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.
[in]VT
          VT is DOUBLE PRECISION array, dimension ( LDU, SMLSIZ+1 ).
         On entry, VT**T contains the right singular vector matrices of
         all subproblems at the bottom level.
[in]K
          K is INTEGER array, dimension ( N ).
[in]DIFL
          DIFL is DOUBLE PRECISION array, dimension ( LDU, NLVL ).
         where NLVL = INT(log_2 (N/(SMLSIZ+1))) + 1.
[in]DIFR
          DIFR is DOUBLE PRECISION array, dimension ( LDU, 2 * NLVL ).
         On entry, DIFL(*, I) and DIFR(*, 2 * I -1) record
         distances between singular values on the I-th level and
         singular values on the (I -1)-th level, and DIFR(*, 2 * I)
         record the normalizing factors of the right singular vectors
         matrices of subproblems on I-th level.
[in]Z
          Z is DOUBLE PRECISION array, dimension ( LDU, NLVL ).
         On entry, Z(1, I) contains the components of the deflation-
         adjusted updating row vector for subproblems on the I-th
         level.
[in]POLES
          POLES is DOUBLE PRECISION array, dimension ( LDU, 2 * NLVL ).
         On entry, POLES(*, 2 * I -1: 2 * I) contains the new and old
         singular values involved in the secular equations on the I-th
         level.
[in]GIVPTR
          GIVPTR is INTEGER array, dimension ( N ).
         On entry, GIVPTR( I ) records the number of Givens
         rotations performed on the I-th problem on the computation
         tree.
[in]GIVCOL
          GIVCOL is INTEGER array, dimension ( LDGCOL, 2 * NLVL ).
         On entry, for each I, GIVCOL(*, 2 * I - 1: 2 * I) records 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.
[in]PERM
          PERM is INTEGER array, dimension ( LDGCOL, NLVL ).
         On entry, PERM(*, I) records permutations done on the I-th
         level of the computation tree.
[in]GIVNUM
          GIVNUM is DOUBLE PRECISION array, dimension ( LDU, 2 * NLVL ).
         On entry, GIVNUM(*, 2 *I -1 : 2 * I) records the C- and S-
         values of Givens rotations performed on the I-th level on the
         computation tree.
[in]C
          C is DOUBLE PRECISION array, dimension ( N ).
         On entry, if the I-th subproblem is not square,
         C( I ) contains the C-value of a Givens rotation related to
         the right null space of the I-th subproblem.
[in]S
          S is DOUBLE PRECISION array, dimension ( N ).
         On entry, if the I-th subproblem is not square,
         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 (N)
[out]IWORK
          IWORK is INTEGER array, dimension (3*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Ming Gu and Ren-Cang Li, Computer Science Division, University of California at Berkeley, USA
Osni Marques, LBNL/NERSC, USA

Definition at line 263 of file dlalsa.f.

267*
268* -- LAPACK computational routine --
269* -- LAPACK is a software package provided by Univ. of Tennessee, --
270* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
271*
272* .. Scalar Arguments ..
273 INTEGER ICOMPQ, INFO, LDB, LDBX, LDGCOL, LDU, N, NRHS,
274 $ SMLSIZ
275* ..
276* .. Array Arguments ..
277 INTEGER GIVCOL( LDGCOL, * ), GIVPTR( * ), IWORK( * ),
278 $ K( * ), PERM( LDGCOL, * )
279 DOUBLE PRECISION B( LDB, * ), BX( LDBX, * ), C( * ),
280 $ DIFL( LDU, * ), DIFR( LDU, * ),
281 $ GIVNUM( LDU, * ), POLES( LDU, * ), S( * ),
282 $ U( LDU, * ), VT( LDU, * ), WORK( * ),
283 $ Z( LDU, * )
284* ..
285*
286* =====================================================================
287*
288* .. Parameters ..
289 DOUBLE PRECISION ZERO, ONE
290 parameter( zero = 0.0d0, one = 1.0d0 )
291* ..
292* .. Local Scalars ..
293 INTEGER I, I1, IC, IM1, INODE, J, LF, LL, LVL, LVL2,
294 $ ND, NDB1, NDIML, NDIMR, NL, NLF, NLP1, NLVL,
295 $ NR, NRF, NRP1, SQRE
296* ..
297* .. External Subroutines ..
298 EXTERNAL dcopy, dgemm, dlals0, dlasdt, xerbla
299* ..
300* .. Executable Statements ..
301*
302* Test the input parameters.
303*
304 info = 0
305*
306 IF( ( icompq.LT.0 ) .OR. ( icompq.GT.1 ) ) THEN
307 info = -1
308 ELSE IF( smlsiz.LT.3 ) THEN
309 info = -2
310 ELSE IF( n.LT.smlsiz ) THEN
311 info = -3
312 ELSE IF( nrhs.LT.1 ) THEN
313 info = -4
314 ELSE IF( ldb.LT.n ) THEN
315 info = -6
316 ELSE IF( ldbx.LT.n ) THEN
317 info = -8
318 ELSE IF( ldu.LT.n ) THEN
319 info = -10
320 ELSE IF( ldgcol.LT.n ) THEN
321 info = -19
322 END IF
323 IF( info.NE.0 ) THEN
324 CALL xerbla( 'DLALSA', -info )
325 RETURN
326 END IF
327*
328* Book-keeping and setting up the computation tree.
329*
330 inode = 1
331 ndiml = inode + n
332 ndimr = ndiml + n
333*
334 CALL dlasdt( n, nlvl, nd, iwork( inode ), iwork( ndiml ),
335 $ iwork( ndimr ), smlsiz )
336*
337* The following code applies back the left singular vector factors.
338* For applying back the right singular vector factors, go to 50.
339*
340 IF( icompq.EQ.1 ) THEN
341 GO TO 50
342 END IF
343*
344* The nodes on the bottom level of the tree were solved
345* by DLASDQ. The corresponding left and right singular vector
346* matrices are in explicit form. First apply back the left
347* singular vector matrices.
348*
349 ndb1 = ( nd+1 ) / 2
350 DO 10 i = ndb1, nd
351*
352* IC : center row of each node
353* NL : number of rows of left subproblem
354* NR : number of rows of right subproblem
355* NLF: starting row of the left subproblem
356* NRF: starting row of the right subproblem
357*
358 i1 = i - 1
359 ic = iwork( inode+i1 )
360 nl = iwork( ndiml+i1 )
361 nr = iwork( ndimr+i1 )
362 nlf = ic - nl
363 nrf = ic + 1
364 CALL dgemm( 'T', 'N', nl, nrhs, nl, one, u( nlf, 1 ), ldu,
365 $ b( nlf, 1 ), ldb, zero, bx( nlf, 1 ), ldbx )
366 CALL dgemm( 'T', 'N', nr, nrhs, nr, one, u( nrf, 1 ), ldu,
367 $ b( nrf, 1 ), ldb, zero, bx( nrf, 1 ), ldbx )
368 10 CONTINUE
369*
370* Next copy the rows of B that correspond to unchanged rows
371* in the bidiagonal matrix to BX.
372*
373 DO 20 i = 1, nd
374 ic = iwork( inode+i-1 )
375 CALL dcopy( nrhs, b( ic, 1 ), ldb, bx( ic, 1 ), ldbx )
376 20 CONTINUE
377*
378* Finally go through the left singular vector matrices of all
379* the other subproblems bottom-up on the tree.
380*
381 j = 2**nlvl
382 sqre = 0
383*
384 DO 40 lvl = nlvl, 1, -1
385 lvl2 = 2*lvl - 1
386*
387* find the first node LF and last node LL on
388* the current level LVL
389*
390 IF( lvl.EQ.1 ) THEN
391 lf = 1
392 ll = 1
393 ELSE
394 lf = 2**( lvl-1 )
395 ll = 2*lf - 1
396 END IF
397 DO 30 i = lf, ll
398 im1 = i - 1
399 ic = iwork( inode+im1 )
400 nl = iwork( ndiml+im1 )
401 nr = iwork( ndimr+im1 )
402 nlf = ic - nl
403 nrf = ic + 1
404 j = j - 1
405 CALL dlals0( icompq, nl, nr, sqre, nrhs, bx( nlf, 1 ), ldbx,
406 $ b( nlf, 1 ), ldb, perm( nlf, lvl ),
407 $ givptr( j ), givcol( nlf, lvl2 ), ldgcol,
408 $ givnum( nlf, lvl2 ), ldu, poles( nlf, lvl2 ),
409 $ difl( nlf, lvl ), difr( nlf, lvl2 ),
410 $ z( nlf, lvl ), k( j ), c( j ), s( j ), work,
411 $ info )
412 30 CONTINUE
413 40 CONTINUE
414 GO TO 90
415*
416* ICOMPQ = 1: applying back the right singular vector factors.
417*
418 50 CONTINUE
419*
420* First now go through the right singular vector matrices of all
421* the tree nodes top-down.
422*
423 j = 0
424 DO 70 lvl = 1, nlvl
425 lvl2 = 2*lvl - 1
426*
427* Find the first node LF and last node LL on
428* the current level LVL.
429*
430 IF( lvl.EQ.1 ) THEN
431 lf = 1
432 ll = 1
433 ELSE
434 lf = 2**( lvl-1 )
435 ll = 2*lf - 1
436 END IF
437 DO 60 i = ll, lf, -1
438 im1 = i - 1
439 ic = iwork( inode+im1 )
440 nl = iwork( ndiml+im1 )
441 nr = iwork( ndimr+im1 )
442 nlf = ic - nl
443 nrf = ic + 1
444 IF( i.EQ.ll ) THEN
445 sqre = 0
446 ELSE
447 sqre = 1
448 END IF
449 j = j + 1
450 CALL dlals0( icompq, nl, nr, sqre, nrhs, b( nlf, 1 ), ldb,
451 $ bx( nlf, 1 ), ldbx, perm( nlf, lvl ),
452 $ givptr( j ), givcol( nlf, lvl2 ), ldgcol,
453 $ givnum( nlf, lvl2 ), ldu, poles( nlf, lvl2 ),
454 $ difl( nlf, lvl ), difr( nlf, lvl2 ),
455 $ z( nlf, lvl ), k( j ), c( j ), s( j ), work,
456 $ info )
457 60 CONTINUE
458 70 CONTINUE
459*
460* The nodes on the bottom level of the tree were solved
461* by DLASDQ. The corresponding right singular vector
462* matrices are in explicit form. Apply them back.
463*
464 ndb1 = ( nd+1 ) / 2
465 DO 80 i = ndb1, nd
466 i1 = i - 1
467 ic = iwork( inode+i1 )
468 nl = iwork( ndiml+i1 )
469 nr = iwork( ndimr+i1 )
470 nlp1 = nl + 1
471 IF( i.EQ.nd ) THEN
472 nrp1 = nr
473 ELSE
474 nrp1 = nr + 1
475 END IF
476 nlf = ic - nl
477 nrf = ic + 1
478 CALL dgemm( 'T', 'N', nlp1, nrhs, nlp1, one, vt( nlf, 1 ), ldu,
479 $ b( nlf, 1 ), ldb, zero, bx( nlf, 1 ), ldbx )
480 CALL dgemm( 'T', 'N', nrp1, nrhs, nrp1, one, vt( nrf, 1 ), ldu,
481 $ b( nrf, 1 ), ldb, zero, bx( nrf, 1 ), ldbx )
482 80 CONTINUE
483*
484 90 CONTINUE
485*
486 RETURN
487*
488* End of DLALSA
489*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dlals0(icompq, nl, nr, sqre, nrhs, b, ldb, bx, ldbx, perm, givptr, givcol, ldgcol, givnum, ldgnum, poles, difl, difr, z, k, c, s, work, info)
DLALS0 applies back multiplying factors in solving the least squares problem using divide and conquer...
Definition dlals0.f:268
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:105
Here is the call graph for this function:
Here is the caller graph for this function: