LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlals0.f
Go to the documentation of this file.
1*> \brief \b DLALS0 applies back multiplying factors in solving the least squares problem using divide and conquer SVD approach. Used by sgelsd.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DLALS0 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlals0.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlals0.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlals0.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DLALS0( ICOMPQ, NL, NR, SQRE, NRHS, B, LDB, BX, LDBX,
22* PERM, GIVPTR, GIVCOL, LDGCOL, GIVNUM, LDGNUM,
23* POLES, DIFL, DIFR, Z, K, C, S, WORK, INFO )
24*
25* .. Scalar Arguments ..
26* INTEGER GIVPTR, ICOMPQ, INFO, K, LDB, LDBX, LDGCOL,
27* $ LDGNUM, NL, NR, NRHS, SQRE
28* DOUBLE PRECISION C, S
29* ..
30* .. Array Arguments ..
31* INTEGER GIVCOL( LDGCOL, * ), PERM( * )
32* DOUBLE PRECISION B( LDB, * ), BX( LDBX, * ), DIFL( * ),
33* $ DIFR( LDGNUM, * ), GIVNUM( LDGNUM, * ),
34* $ POLES( LDGNUM, * ), WORK( * ), Z( * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> DLALS0 applies back the multiplying factors of either the left or the
44*> right singular vector matrix of a diagonal matrix appended by a row
45*> to the right hand side matrix B in solving the least squares problem
46*> using the divide-and-conquer SVD approach.
47*>
48*> For the left singular vector matrix, three types of orthogonal
49*> matrices are involved:
50*>
51*> (1L) Givens rotations: the number of such rotations is GIVPTR; the
52*> pairs of columns/rows they were applied to are stored in GIVCOL;
53*> and the C- and S-values of these rotations are stored in GIVNUM.
54*>
55*> (2L) Permutation. The (NL+1)-st row of B is to be moved to the first
56*> row, and for J=2:N, PERM(J)-th row of B is to be moved to the
57*> J-th row.
58*>
59*> (3L) The left singular vector matrix of the remaining matrix.
60*>
61*> For the right singular vector matrix, four types of orthogonal
62*> matrices are involved:
63*>
64*> (1R) The right singular vector matrix of the remaining matrix.
65*>
66*> (2R) If SQRE = 1, one extra Givens rotation to generate the right
67*> null space.
68*>
69*> (3R) The inverse transformation of (2L).
70*>
71*> (4R) The inverse transformation of (1L).
72*> \endverbatim
73*
74* Arguments:
75* ==========
76*
77*> \param[in] ICOMPQ
78*> \verbatim
79*> ICOMPQ is INTEGER
80*> Specifies whether singular vectors are to be computed in
81*> factored form:
82*> = 0: Left singular vector matrix.
83*> = 1: Right singular vector matrix.
84*> \endverbatim
85*>
86*> \param[in] NL
87*> \verbatim
88*> NL is INTEGER
89*> The row dimension of the upper block. NL >= 1.
90*> \endverbatim
91*>
92*> \param[in] NR
93*> \verbatim
94*> NR is INTEGER
95*> The row dimension of the lower block. NR >= 1.
96*> \endverbatim
97*>
98*> \param[in] SQRE
99*> \verbatim
100*> SQRE is INTEGER
101*> = 0: the lower block is an NR-by-NR square matrix.
102*> = 1: the lower block is an NR-by-(NR+1) rectangular matrix.
103*>
104*> The bidiagonal matrix has row dimension N = NL + NR + 1,
105*> and column dimension M = N + SQRE.
106*> \endverbatim
107*>
108*> \param[in] NRHS
109*> \verbatim
110*> NRHS is INTEGER
111*> The number of columns of B and BX. NRHS must be at least 1.
112*> \endverbatim
113*>
114*> \param[in,out] B
115*> \verbatim
116*> B is DOUBLE PRECISION array, dimension ( LDB, NRHS )
117*> On input, B contains the right hand sides of the least
118*> squares problem in rows 1 through M. On output, B contains
119*> the solution X in rows 1 through N.
120*> \endverbatim
121*>
122*> \param[in] LDB
123*> \verbatim
124*> LDB is INTEGER
125*> The leading dimension of B. LDB must be at least
126*> max(1,MAX( M, N ) ).
127*> \endverbatim
128*>
129*> \param[out] BX
130*> \verbatim
131*> BX is DOUBLE PRECISION array, dimension ( LDBX, NRHS )
132*> \endverbatim
133*>
134*> \param[in] LDBX
135*> \verbatim
136*> LDBX is INTEGER
137*> The leading dimension of BX.
138*> \endverbatim
139*>
140*> \param[in] PERM
141*> \verbatim
142*> PERM is INTEGER array, dimension ( N )
143*> The permutations (from deflation and sorting) applied
144*> to the two blocks.
145*> \endverbatim
146*>
147*> \param[in] GIVPTR
148*> \verbatim
149*> GIVPTR is INTEGER
150*> The number of Givens rotations which took place in this
151*> subproblem.
152*> \endverbatim
153*>
154*> \param[in] GIVCOL
155*> \verbatim
156*> GIVCOL is INTEGER array, dimension ( LDGCOL, 2 )
157*> Each pair of numbers indicates a pair of rows/columns
158*> involved in a Givens rotation.
159*> \endverbatim
160*>
161*> \param[in] LDGCOL
162*> \verbatim
163*> LDGCOL is INTEGER
164*> The leading dimension of GIVCOL, must be at least N.
165*> \endverbatim
166*>
167*> \param[in] GIVNUM
168*> \verbatim
169*> GIVNUM is DOUBLE PRECISION array, dimension ( LDGNUM, 2 )
170*> Each number indicates the C or S value used in the
171*> corresponding Givens rotation.
172*> \endverbatim
173*>
174*> \param[in] LDGNUM
175*> \verbatim
176*> LDGNUM is INTEGER
177*> The leading dimension of arrays DIFR, POLES and
178*> GIVNUM, must be at least K.
179*> \endverbatim
180*>
181*> \param[in] POLES
182*> \verbatim
183*> POLES is DOUBLE PRECISION array, dimension ( LDGNUM, 2 )
184*> On entry, POLES(1:K, 1) contains the new singular
185*> values obtained from solving the secular equation, and
186*> POLES(1:K, 2) is an array containing the poles in the secular
187*> equation.
188*> \endverbatim
189*>
190*> \param[in] DIFL
191*> \verbatim
192*> DIFL is DOUBLE PRECISION array, dimension ( K ).
193*> On entry, DIFL(I) is the distance between I-th updated
194*> (undeflated) singular value and the I-th (undeflated) old
195*> singular value.
196*> \endverbatim
197*>
198*> \param[in] DIFR
199*> \verbatim
200*> DIFR is DOUBLE PRECISION array, dimension ( LDGNUM, 2 ).
201*> On entry, DIFR(I, 1) contains the distances between I-th
202*> updated (undeflated) singular value and the I+1-th
203*> (undeflated) old singular value. And DIFR(I, 2) is the
204*> normalizing factor for the I-th right singular vector.
205*> \endverbatim
206*>
207*> \param[in] Z
208*> \verbatim
209*> Z is DOUBLE PRECISION array, dimension ( K )
210*> Contain the components of the deflation-adjusted updating row
211*> vector.
212*> \endverbatim
213*>
214*> \param[in] K
215*> \verbatim
216*> K is INTEGER
217*> Contains the dimension of the non-deflated matrix,
218*> This is the order of the related secular equation. 1 <= K <=N.
219*> \endverbatim
220*>
221*> \param[in] C
222*> \verbatim
223*> C is DOUBLE PRECISION
224*> C contains garbage if SQRE =0 and the C-value of a Givens
225*> rotation related to the right null space if SQRE = 1.
226*> \endverbatim
227*>
228*> \param[in] S
229*> \verbatim
230*> S is DOUBLE PRECISION
231*> S contains garbage if SQRE =0 and the S-value of a Givens
232*> rotation related to the right null space if SQRE = 1.
233*> \endverbatim
234*>
235*> \param[out] WORK
236*> \verbatim
237*> WORK is DOUBLE PRECISION array, dimension ( K )
238*> \endverbatim
239*>
240*> \param[out] INFO
241*> \verbatim
242*> INFO is INTEGER
243*> = 0: successful exit.
244*> < 0: if INFO = -i, the i-th argument had an illegal value.
245*> \endverbatim
246*
247* Authors:
248* ========
249*
250*> \author Univ. of Tennessee
251*> \author Univ. of California Berkeley
252*> \author Univ. of Colorado Denver
253*> \author NAG Ltd.
254*
255*> \ingroup lals0
256*
257*> \par Contributors:
258* ==================
259*>
260*> Ming Gu and Ren-Cang Li, Computer Science Division, University of
261*> California at Berkeley, USA \n
262*> Osni Marques, LBNL/NERSC, USA \n
263*
264* =====================================================================
265 SUBROUTINE dlals0( ICOMPQ, NL, NR, SQRE, NRHS, B, LDB, BX, LDBX,
266 $ PERM, GIVPTR, GIVCOL, LDGCOL, GIVNUM, LDGNUM,
267 $ POLES, DIFL, DIFR, Z, K, C, S, WORK, INFO )
268*
269* -- LAPACK computational routine --
270* -- LAPACK is a software package provided by Univ. of Tennessee, --
271* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
272*
273* .. Scalar Arguments ..
274 INTEGER GIVPTR, ICOMPQ, INFO, K, LDB, LDBX, LDGCOL,
275 $ LDGNUM, NL, NR, NRHS, SQRE
276 DOUBLE PRECISION C, S
277* ..
278* .. Array Arguments ..
279 INTEGER GIVCOL( LDGCOL, * ), PERM( * )
280 DOUBLE PRECISION B( LDB, * ), BX( LDBX, * ), DIFL( * ),
281 $ difr( ldgnum, * ), givnum( ldgnum, * ),
282 $ poles( ldgnum, * ), work( * ), z( * )
283* ..
284*
285* =====================================================================
286*
287* .. Parameters ..
288 DOUBLE PRECISION ONE, ZERO, NEGONE
289 PARAMETER ( ONE = 1.0d0, zero = 0.0d0, negone = -1.0d0 )
290* ..
291* .. Local Scalars ..
292 INTEGER I, J, M, N, NLP1
293 DOUBLE PRECISION DIFLJ, DIFRJ, DJ, DSIGJ, DSIGJP, TEMP
294* ..
295* .. External Subroutines ..
296 EXTERNAL dcopy, dgemv, dlacpy, dlascl, drot, dscal,
297 $ xerbla
298* ..
299* .. External Functions ..
300 DOUBLE PRECISION DLAMC3, DNRM2
301 EXTERNAL DLAMC3, DNRM2
302* ..
303* .. Intrinsic Functions ..
304 INTRINSIC max
305* ..
306* .. Executable Statements ..
307*
308* Test the input parameters.
309*
310 info = 0
311 n = nl + nr + 1
312*
313 IF( ( icompq.LT.0 ) .OR. ( icompq.GT.1 ) ) THEN
314 info = -1
315 ELSE IF( nl.LT.1 ) THEN
316 info = -2
317 ELSE IF( nr.LT.1 ) THEN
318 info = -3
319 ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
320 info = -4
321 ELSE IF( nrhs.LT.1 ) THEN
322 info = -5
323 ELSE IF( ldb.LT.n ) THEN
324 info = -7
325 ELSE IF( ldbx.LT.n ) THEN
326 info = -9
327 ELSE IF( givptr.LT.0 ) THEN
328 info = -11
329 ELSE IF( ldgcol.LT.n ) THEN
330 info = -13
331 ELSE IF( ldgnum.LT.n ) THEN
332 info = -15
333 ELSE IF( k.LT.1 ) THEN
334 info = -20
335 END IF
336 IF( info.NE.0 ) THEN
337 CALL xerbla( 'DLALS0', -info )
338 RETURN
339 END IF
340*
341 m = n + sqre
342 nlp1 = nl + 1
343*
344 IF( icompq.EQ.0 ) THEN
345*
346* Apply back orthogonal transformations from the left.
347*
348* Step (1L): apply back the Givens rotations performed.
349*
350 DO 10 i = 1, givptr
351 CALL drot( nrhs, b( givcol( i, 2 ), 1 ), ldb,
352 $ b( givcol( i, 1 ), 1 ), ldb, givnum( i, 2 ),
353 $ givnum( i, 1 ) )
354 10 CONTINUE
355*
356* Step (2L): permute rows of B.
357*
358 CALL dcopy( nrhs, b( nlp1, 1 ), ldb, bx( 1, 1 ), ldbx )
359 DO 20 i = 2, n
360 CALL dcopy( nrhs, b( perm( i ), 1 ), ldb, bx( i, 1 ), ldbx )
361 20 CONTINUE
362*
363* Step (3L): apply the inverse of the left singular vector
364* matrix to BX.
365*
366 IF( k.EQ.1 ) THEN
367 CALL dcopy( nrhs, bx, ldbx, b, ldb )
368 IF( z( 1 ).LT.zero ) THEN
369 CALL dscal( nrhs, negone, b, ldb )
370 END IF
371 ELSE
372 DO 50 j = 1, k
373 diflj = difl( j )
374 dj = poles( j, 1 )
375 dsigj = -poles( j, 2 )
376 IF( j.LT.k ) THEN
377 difrj = -difr( j, 1 )
378 dsigjp = -poles( j+1, 2 )
379 END IF
380 IF( ( z( j ).EQ.zero ) .OR. ( poles( j, 2 ).EQ.zero ) )
381 $ THEN
382 work( j ) = zero
383 ELSE
384 work( j ) = -poles( j, 2 )*z( j ) / diflj /
385 $ ( poles( j, 2 )+dj )
386 END IF
387 DO 30 i = 1, j - 1
388 IF( ( z( i ).EQ.zero ) .OR.
389 $ ( poles( i, 2 ).EQ.zero ) ) THEN
390 work( i ) = zero
391 ELSE
392*
393* Use calls to the subroutine DLAMC3 to enforce the
394* parentheses (x+y)+z. The goal is to prevent
395* optimizing compilers from doing x+(y+z).
396*
397 work( i ) = poles( i, 2 )*z( i ) /
398 $ ( dlamc3( poles( i, 2 ), dsigj )-
399 $ diflj ) / ( poles( i, 2 )+dj )
400 END IF
401 30 CONTINUE
402 DO 40 i = j + 1, k
403 IF( ( z( i ).EQ.zero ) .OR.
404 $ ( poles( i, 2 ).EQ.zero ) ) THEN
405 work( i ) = zero
406 ELSE
407 work( i ) = poles( i, 2 )*z( i ) /
408 $ ( dlamc3( poles( i, 2 ), dsigjp )+
409 $ difrj ) / ( poles( i, 2 )+dj )
410 END IF
411 40 CONTINUE
412 work( 1 ) = negone
413 temp = dnrm2( k, work, 1 )
414 CALL dgemv( 'T', k, nrhs, one, bx, ldbx, work, 1, zero,
415 $ b( j, 1 ), ldb )
416 CALL dlascl( 'G', 0, 0, temp, one, 1, nrhs, b( j, 1 ),
417 $ ldb, info )
418 50 CONTINUE
419 END IF
420*
421* Move the deflated rows of BX to B also.
422*
423 IF( k.LT.max( m, n ) )
424 $ CALL dlacpy( 'A', n-k, nrhs, bx( k+1, 1 ), ldbx,
425 $ b( k+1, 1 ), ldb )
426 ELSE
427*
428* Apply back the right orthogonal transformations.
429*
430* Step (1R): apply back the new right singular vector matrix
431* to B.
432*
433 IF( k.EQ.1 ) THEN
434 CALL dcopy( nrhs, b, ldb, bx, ldbx )
435 ELSE
436 DO 80 j = 1, k
437 dsigj = poles( j, 2 )
438 IF( z( j ).EQ.zero ) THEN
439 work( j ) = zero
440 ELSE
441 work( j ) = -z( j ) / difl( j ) /
442 $ ( dsigj+poles( j, 1 ) ) / difr( j, 2 )
443 END IF
444 DO 60 i = 1, j - 1
445 IF( z( j ).EQ.zero ) THEN
446 work( i ) = zero
447 ELSE
448*
449* Use calls to the subroutine DLAMC3 to enforce the
450* parentheses (x+y)+z. The goal is to prevent
451* optimizing compilers from doing x+(y+z).
452*
453 work( i ) = z( j ) / ( dlamc3( dsigj, -poles( i+1,
454 $ 2 ) )-difr( i, 1 ) ) /
455 $ ( dsigj+poles( i, 1 ) ) / difr( i, 2 )
456 END IF
457 60 CONTINUE
458 DO 70 i = j + 1, k
459 IF( z( j ).EQ.zero ) THEN
460 work( i ) = zero
461 ELSE
462 work( i ) = z( j ) / ( dlamc3( dsigj, -poles( i,
463 $ 2 ) )-difl( i ) ) /
464 $ ( dsigj+poles( i, 1 ) ) / difr( i, 2 )
465 END IF
466 70 CONTINUE
467 CALL dgemv( 'T', k, nrhs, one, b, ldb, work, 1, zero,
468 $ bx( j, 1 ), ldbx )
469 80 CONTINUE
470 END IF
471*
472* Step (2R): if SQRE = 1, apply back the rotation that is
473* related to the right null space of the subproblem.
474*
475 IF( sqre.EQ.1 ) THEN
476 CALL dcopy( nrhs, b( m, 1 ), ldb, bx( m, 1 ), ldbx )
477 CALL drot( nrhs, bx( 1, 1 ), ldbx, bx( m, 1 ), ldbx, c, s )
478 END IF
479 IF( k.LT.max( m, n ) )
480 $ CALL dlacpy( 'A', n-k, nrhs, b( k+1, 1 ), ldb, bx( k+1, 1 ),
481 $ ldbx )
482*
483* Step (3R): permute rows of B.
484*
485 CALL dcopy( nrhs, bx( 1, 1 ), ldbx, b( nlp1, 1 ), ldb )
486 IF( sqre.EQ.1 ) THEN
487 CALL dcopy( nrhs, bx( m, 1 ), ldbx, b( m, 1 ), ldb )
488 END IF
489 DO 90 i = 2, n
490 CALL dcopy( nrhs, bx( i, 1 ), ldbx, b( perm( i ), 1 ), ldb )
491 90 CONTINUE
492*
493* Step (4R): apply back the Givens rotations performed.
494*
495 DO 100 i = givptr, 1, -1
496 CALL drot( nrhs, b( givcol( i, 2 ), 1 ), ldb,
497 $ b( givcol( i, 1 ), 1 ), ldb, givnum( i, 2 ),
498 $ -givnum( i, 1 ) )
499 100 CONTINUE
500 END IF
501*
502 RETURN
503*
504* End of DLALS0
505*
506 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
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
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 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
subroutine drot(n, dx, incx, dy, incy, c, s)
DROT
Definition drot.f:92
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79