LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slalsa.f
Go to the documentation of this file.
1*> \brief \b SLALSA computes the SVD of the coefficient matrix in compact form. Used by sgelsd.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SLALSA + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slalsa.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slalsa.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slalsa.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SLALSA( ICOMPQ, SMLSIZ, N, NRHS, B, LDB, BX, LDBX, U,
20* LDU, VT, K, DIFL, DIFR, Z, POLES, GIVPTR,
21* GIVCOL, LDGCOL, PERM, GIVNUM, C, S, WORK,
22* IWORK, INFO )
23*
24* .. Scalar Arguments ..
25* INTEGER ICOMPQ, INFO, LDB, LDBX, LDGCOL, LDU, N, NRHS,
26* $ SMLSIZ
27* ..
28* .. Array Arguments ..
29* INTEGER GIVCOL( LDGCOL, * ), GIVPTR( * ), IWORK( * ),
30* $ K( * ), PERM( LDGCOL, * )
31* REAL B( LDB, * ), BX( LDBX, * ), C( * ),
32* $ DIFL( LDU, * ), DIFR( LDU, * ),
33* $ GIVNUM( LDU, * ), POLES( LDU, * ), S( * ),
34* $ U( LDU, * ), VT( LDU, * ), WORK( * ),
35* $ Z( LDU, * )
36* ..
37*
38*
39*> \par Purpose:
40* =============
41*>
42*> \verbatim
43*>
44*> SLALSA is an intermediate step in solving the least squares problem
45*> by computing the SVD of the coefficient matrix in compact form (The
46*> singular vectors are computed as products of simple orthogonal
47*> matrices.).
48*>
49*> If ICOMPQ = 0, SLALSA applies the inverse of the left singular vector
50*> matrix of an upper bidiagonal matrix to the right hand side; and if
51*> ICOMPQ = 1, SLALSA applies the right singular vector matrix to the
52*> right hand side. The singular vector matrices were generated in
53*> compact form by SLALSA.
54*> \endverbatim
55*
56* Arguments:
57* ==========
58*
59*> \param[in] ICOMPQ
60*> \verbatim
61*> ICOMPQ is INTEGER
62*> Specifies whether the left or the right singular vector
63*> matrix is involved.
64*> = 0: Left singular vector matrix
65*> = 1: Right singular vector matrix
66*> \endverbatim
67*>
68*> \param[in] SMLSIZ
69*> \verbatim
70*> SMLSIZ is INTEGER
71*> The maximum size of the subproblems at the bottom of the
72*> computation tree.
73*> \endverbatim
74*>
75*> \param[in] N
76*> \verbatim
77*> N is INTEGER
78*> The row and column dimensions of the upper bidiagonal matrix.
79*> \endverbatim
80*>
81*> \param[in] NRHS
82*> \verbatim
83*> NRHS is INTEGER
84*> The number of columns of B and BX. NRHS must be at least 1.
85*> \endverbatim
86*>
87*> \param[in,out] B
88*> \verbatim
89*> B is REAL array, dimension ( LDB, NRHS )
90*> On input, B contains the right hand sides of the least
91*> squares problem in rows 1 through M.
92*> On output, B contains the solution X in rows 1 through N.
93*> \endverbatim
94*>
95*> \param[in] LDB
96*> \verbatim
97*> LDB is INTEGER
98*> The leading dimension of B in the calling subprogram.
99*> LDB must be at least max(1,MAX( M, N ) ).
100*> \endverbatim
101*>
102*> \param[out] BX
103*> \verbatim
104*> BX is REAL array, dimension ( LDBX, NRHS )
105*> On exit, the result of applying the left or right singular
106*> vector matrix to B.
107*> \endverbatim
108*>
109*> \param[in] LDBX
110*> \verbatim
111*> LDBX is INTEGER
112*> The leading dimension of BX.
113*> \endverbatim
114*>
115*> \param[in] U
116*> \verbatim
117*> U is REAL array, dimension ( LDU, SMLSIZ ).
118*> On entry, U contains the left singular vector matrices of all
119*> subproblems at the bottom level.
120*> \endverbatim
121*>
122*> \param[in] LDU
123*> \verbatim
124*> LDU is INTEGER, LDU = > N.
125*> The leading dimension of arrays U, VT, DIFL, DIFR,
126*> POLES, GIVNUM, and Z.
127*> \endverbatim
128*>
129*> \param[in] VT
130*> \verbatim
131*> VT is REAL array, dimension ( LDU, SMLSIZ+1 ).
132*> On entry, VT**T contains the right singular vector matrices of
133*> all subproblems at the bottom level.
134*> \endverbatim
135*>
136*> \param[in] K
137*> \verbatim
138*> K is INTEGER array, dimension ( N ).
139*> \endverbatim
140*>
141*> \param[in] DIFL
142*> \verbatim
143*> DIFL is REAL array, dimension ( LDU, NLVL ).
144*> where NLVL = INT(log_2 (N/(SMLSIZ+1))) + 1.
145*> \endverbatim
146*>
147*> \param[in] DIFR
148*> \verbatim
149*> DIFR is REAL array, dimension ( LDU, 2 * NLVL ).
150*> On entry, DIFL(*, I) and DIFR(*, 2 * I -1) record
151*> distances between singular values on the I-th level and
152*> singular values on the (I -1)-th level, and DIFR(*, 2 * I)
153*> record the normalizing factors of the right singular vectors
154*> matrices of subproblems on I-th level.
155*> \endverbatim
156*>
157*> \param[in] Z
158*> \verbatim
159*> Z is REAL array, dimension ( LDU, NLVL ).
160*> On entry, Z(1, I) contains the components of the deflation-
161*> adjusted updating row vector for subproblems on the I-th
162*> level.
163*> \endverbatim
164*>
165*> \param[in] POLES
166*> \verbatim
167*> POLES is REAL array, dimension ( LDU, 2 * NLVL ).
168*> On entry, POLES(*, 2 * I -1: 2 * I) contains the new and old
169*> singular values involved in the secular equations on the I-th
170*> level.
171*> \endverbatim
172*>
173*> \param[in] GIVPTR
174*> \verbatim
175*> GIVPTR is INTEGER array, dimension ( N ).
176*> On entry, GIVPTR( I ) records the number of Givens
177*> rotations performed on the I-th problem on the computation
178*> tree.
179*> \endverbatim
180*>
181*> \param[in] GIVCOL
182*> \verbatim
183*> GIVCOL is INTEGER array, dimension ( LDGCOL, 2 * NLVL ).
184*> On entry, for each I, GIVCOL(*, 2 * I - 1: 2 * I) records the
185*> locations of Givens rotations performed on the I-th level on
186*> the computation tree.
187*> \endverbatim
188*>
189*> \param[in] LDGCOL
190*> \verbatim
191*> LDGCOL is INTEGER, LDGCOL = > N.
192*> The leading dimension of arrays GIVCOL and PERM.
193*> \endverbatim
194*>
195*> \param[in] PERM
196*> \verbatim
197*> PERM is INTEGER array, dimension ( LDGCOL, NLVL ).
198*> On entry, PERM(*, I) records permutations done on the I-th
199*> level of the computation tree.
200*> \endverbatim
201*>
202*> \param[in] GIVNUM
203*> \verbatim
204*> GIVNUM is REAL array, dimension ( LDU, 2 * NLVL ).
205*> On entry, GIVNUM(*, 2 *I -1 : 2 * I) records the C- and S-
206*> values of Givens rotations performed on the I-th level on the
207*> computation tree.
208*> \endverbatim
209*>
210*> \param[in] C
211*> \verbatim
212*> C is REAL array, dimension ( N ).
213*> On entry, if the I-th subproblem is not square,
214*> C( I ) contains the C-value of a Givens rotation related to
215*> the right null space of the I-th subproblem.
216*> \endverbatim
217*>
218*> \param[in] S
219*> \verbatim
220*> S is REAL array, dimension ( N ).
221*> On entry, if the I-th subproblem is not square,
222*> S( I ) contains the S-value of a Givens rotation related to
223*> the right null space of the I-th subproblem.
224*> \endverbatim
225*>
226*> \param[out] WORK
227*> \verbatim
228*> WORK is REAL array, dimension (N)
229*> \endverbatim
230*>
231*> \param[out] IWORK
232*> \verbatim
233*> IWORK is INTEGER array, dimension (3*N)
234*> \endverbatim
235*>
236*> \param[out] INFO
237*> \verbatim
238*> INFO is INTEGER
239*> = 0: successful exit.
240*> < 0: if INFO = -i, the i-th argument had an illegal value.
241*> \endverbatim
242*
243* Authors:
244* ========
245*
246*> \author Univ. of Tennessee
247*> \author Univ. of California Berkeley
248*> \author Univ. of Colorado Denver
249*> \author NAG Ltd.
250*
251*> \ingroup lalsa
252*
253*> \par Contributors:
254* ==================
255*>
256*> Ming Gu and Ren-Cang Li, Computer Science Division, University of
257*> California at Berkeley, USA \n
258*> Osni Marques, LBNL/NERSC, USA \n
259*
260* =====================================================================
261 SUBROUTINE slalsa( ICOMPQ, SMLSIZ, N, NRHS, B, LDB, BX, LDBX,
262 $ U,
263 $ LDU, VT, K, DIFL, DIFR, Z, POLES, GIVPTR,
264 $ GIVCOL, LDGCOL, PERM, GIVNUM, C, S, WORK,
265 $ IWORK, INFO )
266*
267* -- LAPACK computational routine --
268* -- LAPACK is a software package provided by Univ. of Tennessee, --
269* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
270*
271* .. Scalar Arguments ..
272 INTEGER ICOMPQ, INFO, LDB, LDBX, LDGCOL, LDU, N, NRHS,
273 $ SMLSIZ
274* ..
275* .. Array Arguments ..
276 INTEGER GIVCOL( LDGCOL, * ), GIVPTR( * ), IWORK( * ),
277 $ K( * ), PERM( LDGCOL, * )
278 REAL B( LDB, * ), BX( LDBX, * ), C( * ),
279 $ DIFL( LDU, * ), DIFR( LDU, * ),
280 $ givnum( ldu, * ), poles( ldu, * ), s( * ),
281 $ u( ldu, * ), vt( ldu, * ), work( * ),
282 $ z( ldu, * )
283* ..
284*
285* =====================================================================
286*
287* .. Parameters ..
288 REAL ZERO, ONE
289 PARAMETER ( ZERO = 0.0e0, one = 1.0e0 )
290* ..
291* .. Local Scalars ..
292 INTEGER I, I1, IC, IM1, INODE, J, LF, LL, LVL, LVL2,
293 $ nd, ndb1, ndiml, ndimr, nl, nlf, nlp1, nlvl,
294 $ nr, nrf, nrp1, sqre
295* ..
296* .. External Subroutines ..
297 EXTERNAL scopy, sgemm, slals0, slasdt,
298 $ 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( 'SLALSA', -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 slasdt( 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 SLASDQ. 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 sgemm( 'T', 'N', nl, nrhs, nl, one, u( nlf, 1 ), ldu,
365 $ b( nlf, 1 ), ldb, zero, bx( nlf, 1 ), ldbx )
366 CALL sgemm( '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 scopy( 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 slals0( icompq, nl, nr, sqre, nrhs, bx( nlf, 1 ),
406 $ ldbx,
407 $ b( nlf, 1 ), ldb, perm( nlf, lvl ),
408 $ givptr( j ), givcol( nlf, lvl2 ), ldgcol,
409 $ givnum( nlf, lvl2 ), ldu, poles( nlf, lvl2 ),
410 $ difl( nlf, lvl ), difr( nlf, lvl2 ),
411 $ z( nlf, lvl ), k( j ), c( j ), s( j ), work,
412 $ info )
413 30 CONTINUE
414 40 CONTINUE
415 GO TO 90
416*
417* ICOMPQ = 1: applying back the right singular vector factors.
418*
419 50 CONTINUE
420*
421* First now go through the right singular vector matrices of all
422* the tree nodes top-down.
423*
424 j = 0
425 DO 70 lvl = 1, nlvl
426 lvl2 = 2*lvl - 1
427*
428* Find the first node LF and last node LL on
429* the current level LVL.
430*
431 IF( lvl.EQ.1 ) THEN
432 lf = 1
433 ll = 1
434 ELSE
435 lf = 2**( lvl-1 )
436 ll = 2*lf - 1
437 END IF
438 DO 60 i = ll, lf, -1
439 im1 = i - 1
440 ic = iwork( inode+im1 )
441 nl = iwork( ndiml+im1 )
442 nr = iwork( ndimr+im1 )
443 nlf = ic - nl
444 nrf = ic + 1
445 IF( i.EQ.ll ) THEN
446 sqre = 0
447 ELSE
448 sqre = 1
449 END IF
450 j = j + 1
451 CALL slals0( icompq, nl, nr, sqre, nrhs, b( nlf, 1 ),
452 $ ldb,
453 $ bx( nlf, 1 ), ldbx, perm( nlf, lvl ),
454 $ givptr( j ), givcol( nlf, lvl2 ), ldgcol,
455 $ givnum( nlf, lvl2 ), ldu, poles( nlf, lvl2 ),
456 $ difl( nlf, lvl ), difr( nlf, lvl2 ),
457 $ z( nlf, lvl ), k( j ), c( j ), s( j ), work,
458 $ info )
459 60 CONTINUE
460 70 CONTINUE
461*
462* The nodes on the bottom level of the tree were solved
463* by SLASDQ. The corresponding right singular vector
464* matrices are in explicit form. Apply them back.
465*
466 ndb1 = ( nd+1 ) / 2
467 DO 80 i = ndb1, nd
468 i1 = i - 1
469 ic = iwork( inode+i1 )
470 nl = iwork( ndiml+i1 )
471 nr = iwork( ndimr+i1 )
472 nlp1 = nl + 1
473 IF( i.EQ.nd ) THEN
474 nrp1 = nr
475 ELSE
476 nrp1 = nr + 1
477 END IF
478 nlf = ic - nl
479 nrf = ic + 1
480 CALL sgemm( 'T', 'N', nlp1, nrhs, nlp1, one, vt( nlf, 1 ),
481 $ ldu,
482 $ b( nlf, 1 ), ldb, zero, bx( nlf, 1 ), ldbx )
483 CALL sgemm( 'T', 'N', nrp1, nrhs, nrp1, one, vt( nrf, 1 ),
484 $ ldu,
485 $ b( nrf, 1 ), ldb, zero, bx( nrf, 1 ), ldbx )
486 80 CONTINUE
487*
488 90 CONTINUE
489*
490 RETURN
491*
492* End of SLALSA
493*
494 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine slals0(icompq, nl, nr, sqre, nrhs, b, ldb, bx, ldbx, perm, givptr, givcol, ldgcol, givnum, ldgnum, poles, difl, difr, z, k, c, s, work, info)
SLALS0 applies back multiplying factors in solving the least squares problem using divide and conquer...
Definition slals0.f:267
subroutine slalsa(icompq, smlsiz, n, nrhs, b, ldb, bx, ldbx, u, ldu, vt, k, difl, difr, z, poles, givptr, givcol, ldgcol, perm, givnum, c, s, work, iwork, info)
SLALSA computes the SVD of the coefficient matrix in compact form. Used by sgelsd.
Definition slalsa.f:266
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:103