LAPACK  3.5.0
LAPACK: Linear Algebra PACKage
dlals0.f File Reference

Go to the source code of this file.

Functions/Subroutines

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 SVD approach. Used by sgelsd. More...
 

Function/Subroutine Documentation

subroutine dlals0 ( integer  ICOMPQ,
integer  NL,
integer  NR,
integer  SQRE,
integer  NRHS,
double precision, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( ldbx, * )  BX,
integer  LDBX,
integer, dimension( * )  PERM,
integer  GIVPTR,
integer, dimension( ldgcol, * )  GIVCOL,
integer  LDGCOL,
double precision, dimension( ldgnum, * )  GIVNUM,
integer  LDGNUM,
double precision, dimension( ldgnum, * )  POLES,
double precision, dimension( * )  DIFL,
double precision, dimension( ldgnum, * )  DIFR,
double precision, dimension( * )  Z,
integer  K,
double precision  C,
double precision  S,
double precision, dimension( * )  WORK,
integer  INFO 
)

DLALS0 applies back multiplying factors in solving the least squares problem using divide and conquer SVD approach. Used by sgelsd.

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

Purpose:
 DLALS0 applies back the multiplying factors of either the left or the
 right singular vector matrix of a diagonal matrix appended by a row
 to the right hand side matrix B in solving the least squares problem
 using the divide-and-conquer SVD approach.

 For the left singular vector matrix, three types of orthogonal
 matrices are involved:

 (1L) Givens rotations: the number of such rotations is GIVPTR; the
      pairs of columns/rows they were applied to are stored in GIVCOL;
      and the C- and S-values of these rotations are stored in GIVNUM.

 (2L) Permutation. The (NL+1)-st row of B is to be moved to the first
      row, and for J=2:N, PERM(J)-th row of B is to be moved to the
      J-th row.

 (3L) The left singular vector matrix of the remaining matrix.

 For the right singular vector matrix, four types of orthogonal
 matrices are involved:

 (1R) The right singular vector matrix of the remaining matrix.

 (2R) If SQRE = 1, one extra Givens rotation to generate the right
      null space.

 (3R) The inverse transformation of (2L).

 (4R) The inverse transformation of (1L).
Parameters
[in]ICOMPQ
          ICOMPQ is INTEGER
         Specifies whether singular vectors are to be computed in
         factored form:
         = 0: Left singular vector matrix.
         = 1: Right singular vector matrix.
[in]NL
          NL is INTEGER
         The row dimension of the upper block. NL >= 1.
[in]NR
          NR is INTEGER
         The row dimension of the lower block. NR >= 1.
[in]SQRE
          SQRE is INTEGER
         = 0: the lower block is an NR-by-NR square matrix.
         = 1: the lower block is an NR-by-(NR+1) rectangular matrix.

         The bidiagonal matrix has row dimension N = NL + NR + 1,
         and column dimension M = N + SQRE.
[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. LDB must be at least
         max(1,MAX( M, N ) ).
[out]BX
          BX is DOUBLE PRECISION array, dimension ( LDBX, NRHS )
[in]LDBX
          LDBX is INTEGER
         The leading dimension of BX.
[in]PERM
          PERM is INTEGER array, dimension ( N )
         The permutations (from deflation and sorting) applied
         to the two blocks.
[in]GIVPTR
          GIVPTR is INTEGER
         The number of Givens rotations which took place in this
         subproblem.
[in]GIVCOL
          GIVCOL is INTEGER array, dimension ( LDGCOL, 2 )
         Each pair of numbers indicates a pair of rows/columns
         involved in a Givens rotation.
[in]LDGCOL
          LDGCOL is INTEGER
         The leading dimension of GIVCOL, must be at least N.
[in]GIVNUM
          GIVNUM is DOUBLE PRECISION array, dimension ( LDGNUM, 2 )
         Each number indicates the C or S value used in the
         corresponding Givens rotation.
[in]LDGNUM
          LDGNUM is INTEGER
         The leading dimension of arrays DIFR, POLES and
         GIVNUM, must be at least K.
[in]POLES
          POLES is DOUBLE PRECISION array, dimension ( LDGNUM, 2 )
         On entry, POLES(1:K, 1) contains the new singular
         values obtained from solving the secular equation, and
         POLES(1:K, 2) is an array containing the poles in the secular
         equation.
[in]DIFL
          DIFL is DOUBLE PRECISION array, dimension ( K ).
         On entry, DIFL(I) is the distance between I-th updated
         (undeflated) singular value and the I-th (undeflated) old
         singular value.
[in]DIFR
          DIFR is DOUBLE PRECISION array, dimension ( LDGNUM, 2 ).
         On entry, DIFR(I, 1) contains the distances between I-th
         updated (undeflated) singular value and the I+1-th
         (undeflated) old singular value. And DIFR(I, 2) is the
         normalizing factor for the I-th right singular vector.
[in]Z
          Z is DOUBLE PRECISION array, dimension ( K )
         Contain the components of the deflation-adjusted updating row
         vector.
[in]K
          K is INTEGER
         Contains the dimension of the non-deflated matrix,
         This is the order of the related secular equation. 1 <= K <=N.
[in]C
          C is DOUBLE PRECISION
         C contains garbage if SQRE =0 and the C-value of a Givens
         rotation related to the right null space if SQRE = 1.
[in]S
          S is DOUBLE PRECISION
         S contains garbage if SQRE =0 and the S-value of a Givens
         rotation related to the right null space if SQRE = 1.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension ( K )
[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.
Date
September 2012
Contributors:
Ming Gu and Ren-Cang Li, Computer Science Division, University of California at Berkeley, USA
Osni Marques, LBNL/NERSC, USA

Definition at line 270 of file dlals0.f.

270 *
271 * -- LAPACK computational routine (version 3.4.2) --
272 * -- LAPACK is a software package provided by Univ. of Tennessee, --
273 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274 * September 2012
275 *
276 * .. Scalar Arguments ..
277  INTEGER givptr, icompq, info, k, ldb, ldbx, ldgcol,
278  $ ldgnum, nl, nr, nrhs, sqre
279  DOUBLE PRECISION c, s
280 * ..
281 * .. Array Arguments ..
282  INTEGER givcol( ldgcol, * ), perm( * )
283  DOUBLE PRECISION b( ldb, * ), bx( ldbx, * ), difl( * ),
284  $ difr( ldgnum, * ), givnum( ldgnum, * ),
285  $ poles( ldgnum, * ), work( * ), z( * )
286 * ..
287 *
288 * =====================================================================
289 *
290 * .. Parameters ..
291  DOUBLE PRECISION one, zero, negone
292  parameter( one = 1.0d0, zero = 0.0d0, negone = -1.0d0 )
293 * ..
294 * .. Local Scalars ..
295  INTEGER i, j, m, n, nlp1
296  DOUBLE PRECISION diflj, difrj, dj, dsigj, dsigjp, temp
297 * ..
298 * .. External Subroutines ..
299  EXTERNAL dcopy, dgemv, dlacpy, dlascl, drot, dscal,
300  $ xerbla
301 * ..
302 * .. External Functions ..
303  DOUBLE PRECISION dlamc3, dnrm2
304  EXTERNAL dlamc3, dnrm2
305 * ..
306 * .. Intrinsic Functions ..
307  INTRINSIC max
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( nl.LT.1 ) THEN
318  info = -2
319  ELSE IF( nr.LT.1 ) THEN
320  info = -3
321  ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
322  info = -4
323  END IF
324 *
325  n = nl + nr + 1
326 *
327  IF( nrhs.LT.1 ) THEN
328  info = -5
329  ELSE IF( ldb.LT.n ) THEN
330  info = -7
331  ELSE IF( ldbx.LT.n ) THEN
332  info = -9
333  ELSE IF( givptr.LT.0 ) THEN
334  info = -11
335  ELSE IF( ldgcol.LT.n ) THEN
336  info = -13
337  ELSE IF( ldgnum.LT.n ) THEN
338  info = -15
339  ELSE IF( k.LT.1 ) THEN
340  info = -20
341  END IF
342  IF( info.NE.0 ) THEN
343  CALL xerbla( 'DLALS0', -info )
344  RETURN
345  END IF
346 *
347  m = n + sqre
348  nlp1 = nl + 1
349 *
350  IF( icompq.EQ.0 ) THEN
351 *
352 * Apply back orthogonal transformations from the left.
353 *
354 * Step (1L): apply back the Givens rotations performed.
355 *
356  DO 10 i = 1, givptr
357  CALL drot( nrhs, b( givcol( i, 2 ), 1 ), ldb,
358  $ b( givcol( i, 1 ), 1 ), ldb, givnum( i, 2 ),
359  $ givnum( i, 1 ) )
360  10 CONTINUE
361 *
362 * Step (2L): permute rows of B.
363 *
364  CALL dcopy( nrhs, b( nlp1, 1 ), ldb, bx( 1, 1 ), ldbx )
365  DO 20 i = 2, n
366  CALL dcopy( nrhs, b( perm( i ), 1 ), ldb, bx( i, 1 ), ldbx )
367  20 CONTINUE
368 *
369 * Step (3L): apply the inverse of the left singular vector
370 * matrix to BX.
371 *
372  IF( k.EQ.1 ) THEN
373  CALL dcopy( nrhs, bx, ldbx, b, ldb )
374  IF( z( 1 ).LT.zero ) THEN
375  CALL dscal( nrhs, negone, b, ldb )
376  END IF
377  ELSE
378  DO 50 j = 1, k
379  diflj = difl( j )
380  dj = poles( j, 1 )
381  dsigj = -poles( j, 2 )
382  IF( j.LT.k ) THEN
383  difrj = -difr( j, 1 )
384  dsigjp = -poles( j+1, 2 )
385  END IF
386  IF( ( z( j ).EQ.zero ) .OR. ( poles( j, 2 ).EQ.zero ) )
387  $ THEN
388  work( j ) = zero
389  ELSE
390  work( j ) = -poles( j, 2 )*z( j ) / diflj /
391  $ ( poles( j, 2 )+dj )
392  END IF
393  DO 30 i = 1, j - 1
394  IF( ( z( i ).EQ.zero ) .OR.
395  $ ( poles( i, 2 ).EQ.zero ) ) THEN
396  work( i ) = zero
397  ELSE
398  work( i ) = poles( i, 2 )*z( i ) /
399  $ ( dlamc3( poles( i, 2 ), dsigj )-
400  $ diflj ) / ( poles( i, 2 )+dj )
401  END IF
402  30 CONTINUE
403  DO 40 i = j + 1, k
404  IF( ( z( i ).EQ.zero ) .OR.
405  $ ( poles( i, 2 ).EQ.zero ) ) THEN
406  work( i ) = zero
407  ELSE
408  work( i ) = poles( i, 2 )*z( i ) /
409  $ ( dlamc3( poles( i, 2 ), dsigjp )+
410  $ difrj ) / ( poles( i, 2 )+dj )
411  END IF
412  40 CONTINUE
413  work( 1 ) = negone
414  temp = dnrm2( k, work, 1 )
415  CALL dgemv( 'T', k, nrhs, one, bx, ldbx, work, 1, zero,
416  $ b( j, 1 ), ldb )
417  CALL dlascl( 'G', 0, 0, temp, one, 1, nrhs, b( j, 1 ),
418  $ ldb, info )
419  50 CONTINUE
420  END IF
421 *
422 * Move the deflated rows of BX to B also.
423 *
424  IF( k.LT.max( m, n ) )
425  $ CALL dlacpy( 'A', n-k, nrhs, bx( k+1, 1 ), ldbx,
426  $ b( k+1, 1 ), ldb )
427  ELSE
428 *
429 * Apply back the right orthogonal transformations.
430 *
431 * Step (1R): apply back the new right singular vector matrix
432 * to B.
433 *
434  IF( k.EQ.1 ) THEN
435  CALL dcopy( nrhs, b, ldb, bx, ldbx )
436  ELSE
437  DO 80 j = 1, k
438  dsigj = poles( j, 2 )
439  IF( z( j ).EQ.zero ) THEN
440  work( j ) = zero
441  ELSE
442  work( j ) = -z( j ) / difl( j ) /
443  $ ( dsigj+poles( j, 1 ) ) / difr( j, 2 )
444  END IF
445  DO 60 i = 1, j - 1
446  IF( z( j ).EQ.zero ) THEN
447  work( i ) = zero
448  ELSE
449  work( i ) = z( j ) / ( dlamc3( dsigj, -poles( i+1,
450  $ 2 ) )-difr( i, 1 ) ) /
451  $ ( dsigj+poles( i, 1 ) ) / difr( i, 2 )
452  END IF
453  60 CONTINUE
454  DO 70 i = j + 1, k
455  IF( z( j ).EQ.zero ) THEN
456  work( i ) = zero
457  ELSE
458  work( i ) = z( j ) / ( dlamc3( dsigj, -poles( i,
459  $ 2 ) )-difl( i ) ) /
460  $ ( dsigj+poles( i, 1 ) ) / difr( i, 2 )
461  END IF
462  70 CONTINUE
463  CALL dgemv( 'T', k, nrhs, one, b, ldb, work, 1, zero,
464  $ bx( j, 1 ), ldbx )
465  80 CONTINUE
466  END IF
467 *
468 * Step (2R): if SQRE = 1, apply back the rotation that is
469 * related to the right null space of the subproblem.
470 *
471  IF( sqre.EQ.1 ) THEN
472  CALL dcopy( nrhs, b( m, 1 ), ldb, bx( m, 1 ), ldbx )
473  CALL drot( nrhs, bx( 1, 1 ), ldbx, bx( m, 1 ), ldbx, c, s )
474  END IF
475  IF( k.LT.max( m, n ) )
476  $ CALL dlacpy( 'A', n-k, nrhs, b( k+1, 1 ), ldb, bx( k+1, 1 ),
477  $ ldbx )
478 *
479 * Step (3R): permute rows of B.
480 *
481  CALL dcopy( nrhs, bx( 1, 1 ), ldbx, b( nlp1, 1 ), ldb )
482  IF( sqre.EQ.1 ) THEN
483  CALL dcopy( nrhs, bx( m, 1 ), ldbx, b( m, 1 ), ldb )
484  END IF
485  DO 90 i = 2, n
486  CALL dcopy( nrhs, bx( i, 1 ), ldbx, b( perm( i ), 1 ), ldb )
487  90 CONTINUE
488 *
489 * Step (4R): apply back the Givens rotations performed.
490 *
491  DO 100 i = givptr, 1, -1
492  CALL drot( nrhs, b( givcol( i, 2 ), 1 ), ldb,
493  $ b( givcol( i, 1 ), 1 ), ldb, givnum( i, 2 ),
494  $ -givnum( i, 1 ) )
495  100 CONTINUE
496  END IF
497 *
498  RETURN
499 *
500 * End of DLALS0
501 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
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:141
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
double precision function dlamc3(A, B)
DLAMC3
Definition: dlamch.f:173
subroutine drot(N, DX, INCX, DY, INCY, C, S)
DROT
Definition: drot.f:53
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
double precision function dnrm2(N, X, INCX)
DNRM2
Definition: dnrm2.f:56
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158

Here is the call graph for this function:

Here is the caller graph for this function: