LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slasdq.f
Go to the documentation of this file.
1*> \brief \b SLASDQ computes the SVD of a real bidiagonal matrix with diagonal d and off-diagonal e. Used by sbdsdc.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SLASDQ + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slasdq.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slasdq.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slasdq.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SLASDQ( UPLO, SQRE, N, NCVT, NRU, NCC, D, E, VT, LDVT,
20* U, LDU, C, LDC, WORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER UPLO
24* INTEGER INFO, LDC, LDU, LDVT, N, NCC, NCVT, NRU, SQRE
25* ..
26* .. Array Arguments ..
27* REAL C( LDC, * ), D( * ), E( * ), U( LDU, * ),
28* $ VT( LDVT, * ), WORK( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> SLASDQ computes the singular value decomposition (SVD) of a real
38*> (upper or lower) bidiagonal matrix with diagonal D and offdiagonal
39*> E, accumulating the transformations if desired. Letting B denote
40*> the input bidiagonal matrix, the algorithm computes orthogonal
41*> matrices Q and P such that B = Q * S * P**T (P**T denotes the transpose
42*> of P). The singular values S are overwritten on D.
43*>
44*> The input matrix U is changed to U * Q if desired.
45*> The input matrix VT is changed to P**T * VT if desired.
46*> The input matrix C is changed to Q**T * C if desired.
47*>
48*> See "Computing Small Singular Values of Bidiagonal Matrices With
49*> Guaranteed High Relative Accuracy," by J. Demmel and W. Kahan,
50*> LAPACK Working Note #3, for a detailed description of the algorithm.
51*> \endverbatim
52*
53* Arguments:
54* ==========
55*
56*> \param[in] UPLO
57*> \verbatim
58*> UPLO is CHARACTER*1
59*> On entry, UPLO specifies whether the input bidiagonal matrix
60*> is upper or lower bidiagonal, and whether it is square are
61*> not.
62*> UPLO = 'U' or 'u' B is upper bidiagonal.
63*> UPLO = 'L' or 'l' B is lower bidiagonal.
64*> \endverbatim
65*>
66*> \param[in] SQRE
67*> \verbatim
68*> SQRE is INTEGER
69*> = 0: then the input matrix is N-by-N.
70*> = 1: then the input matrix is N-by-(N+1) if UPLU = 'U' and
71*> (N+1)-by-N if UPLU = 'L'.
72*>
73*> The bidiagonal matrix has
74*> N = NL + NR + 1 rows and
75*> M = N + SQRE >= N columns.
76*> \endverbatim
77*>
78*> \param[in] N
79*> \verbatim
80*> N is INTEGER
81*> On entry, N specifies the number of rows and columns
82*> in the matrix. N must be at least 0.
83*> \endverbatim
84*>
85*> \param[in] NCVT
86*> \verbatim
87*> NCVT is INTEGER
88*> On entry, NCVT specifies the number of columns of
89*> the matrix VT. NCVT must be at least 0.
90*> \endverbatim
91*>
92*> \param[in] NRU
93*> \verbatim
94*> NRU is INTEGER
95*> On entry, NRU specifies the number of rows of
96*> the matrix U. NRU must be at least 0.
97*> \endverbatim
98*>
99*> \param[in] NCC
100*> \verbatim
101*> NCC is INTEGER
102*> On entry, NCC specifies the number of columns of
103*> the matrix C. NCC must be at least 0.
104*> \endverbatim
105*>
106*> \param[in,out] D
107*> \verbatim
108*> D is REAL array, dimension (N)
109*> On entry, D contains the diagonal entries of the
110*> bidiagonal matrix whose SVD is desired. On normal exit,
111*> D contains the singular values in ascending order.
112*> \endverbatim
113*>
114*> \param[in,out] E
115*> \verbatim
116*> E is REAL array.
117*> dimension is (N-1) if SQRE = 0 and N if SQRE = 1.
118*> On entry, the entries of E contain the offdiagonal entries
119*> of the bidiagonal matrix whose SVD is desired. On normal
120*> exit, E will contain 0. If the algorithm does not converge,
121*> D and E will contain the diagonal and superdiagonal entries
122*> of a bidiagonal matrix orthogonally equivalent to the one
123*> given as input.
124*> \endverbatim
125*>
126*> \param[in,out] VT
127*> \verbatim
128*> VT is REAL array, dimension (LDVT, NCVT)
129*> On entry, contains a matrix which on exit has been
130*> premultiplied by P**T, dimension N-by-NCVT if SQRE = 0
131*> and (N+1)-by-NCVT if SQRE = 1 (not referenced if NCVT=0).
132*> \endverbatim
133*>
134*> \param[in] LDVT
135*> \verbatim
136*> LDVT is INTEGER
137*> On entry, LDVT specifies the leading dimension of VT as
138*> declared in the calling (sub) program. LDVT must be at
139*> least 1. If NCVT is nonzero LDVT must also be at least N.
140*> \endverbatim
141*>
142*> \param[in,out] U
143*> \verbatim
144*> U is REAL array, dimension (LDU, N)
145*> On entry, contains a matrix which on exit has been
146*> postmultiplied by Q, dimension NRU-by-N if SQRE = 0
147*> and NRU-by-(N+1) if SQRE = 1 (not referenced if NRU=0).
148*> \endverbatim
149*>
150*> \param[in] LDU
151*> \verbatim
152*> LDU is INTEGER
153*> On entry, LDU specifies the leading dimension of U as
154*> declared in the calling (sub) program. LDU must be at
155*> least max( 1, NRU ) .
156*> \endverbatim
157*>
158*> \param[in,out] C
159*> \verbatim
160*> C is REAL array, dimension (LDC, NCC)
161*> On entry, contains an N-by-NCC matrix which on exit
162*> has been premultiplied by Q**T dimension N-by-NCC if SQRE = 0
163*> and (N+1)-by-NCC if SQRE = 1 (not referenced if NCC=0).
164*> \endverbatim
165*>
166*> \param[in] LDC
167*> \verbatim
168*> LDC is INTEGER
169*> On entry, LDC specifies the leading dimension of C as
170*> declared in the calling (sub) program. LDC must be at
171*> least 1. If NCC is nonzero, LDC must also be at least N.
172*> \endverbatim
173*>
174*> \param[out] WORK
175*> \verbatim
176*> WORK is REAL array, dimension (4*N)
177*> Workspace. Only referenced if one of NCVT, NRU, or NCC is
178*> nonzero, and if N is at least 2.
179*> \endverbatim
180*>
181*> \param[out] INFO
182*> \verbatim
183*> INFO is INTEGER
184*> On exit, a value of 0 indicates a successful exit.
185*> If INFO < 0, argument number -INFO is illegal.
186*> If INFO > 0, the algorithm did not converge, and INFO
187*> specifies how many superdiagonals did not converge.
188*> \endverbatim
189*
190* Authors:
191* ========
192*
193*> \author Univ. of Tennessee
194*> \author Univ. of California Berkeley
195*> \author Univ. of Colorado Denver
196*> \author NAG Ltd.
197*
198*> \ingroup lasdq
199*
200*> \par Contributors:
201* ==================
202*>
203*> Ming Gu and Huan Ren, Computer Science Division, University of
204*> California at Berkeley, USA
205*>
206* =====================================================================
207 SUBROUTINE slasdq( UPLO, SQRE, N, NCVT, NRU, NCC, D, E, VT,
208 $ LDVT,
209 $ U, LDU, C, LDC, WORK, INFO )
210*
211* -- LAPACK auxiliary routine --
212* -- LAPACK is a software package provided by Univ. of Tennessee, --
213* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
214*
215* .. Scalar Arguments ..
216 CHARACTER UPLO
217 INTEGER INFO, LDC, LDU, LDVT, N, NCC, NCVT, NRU, SQRE
218* ..
219* .. Array Arguments ..
220 REAL C( LDC, * ), D( * ), E( * ), U( LDU, * ),
221 $ VT( LDVT, * ), WORK( * )
222* ..
223*
224* =====================================================================
225*
226* .. Parameters ..
227 REAL ZERO
228 PARAMETER ( ZERO = 0.0e+0 )
229* ..
230* .. Local Scalars ..
231 LOGICAL ROTATE
232 INTEGER I, ISUB, IUPLO, J, NP1, SQRE1
233 REAL CS, R, SMIN, SN
234* ..
235* .. External Subroutines ..
236 EXTERNAL sbdsqr, slartg, slasr, sswap,
237 $ xerbla
238* ..
239* .. External Functions ..
240 LOGICAL LSAME
241 EXTERNAL LSAME
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC max
245* ..
246* .. Executable Statements ..
247*
248* Test the input parameters.
249*
250 info = 0
251 iuplo = 0
252 IF( lsame( uplo, 'U' ) )
253 $ iuplo = 1
254 IF( lsame( uplo, 'L' ) )
255 $ iuplo = 2
256 IF( iuplo.EQ.0 ) THEN
257 info = -1
258 ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
259 info = -2
260 ELSE IF( n.LT.0 ) THEN
261 info = -3
262 ELSE IF( ncvt.LT.0 ) THEN
263 info = -4
264 ELSE IF( nru.LT.0 ) THEN
265 info = -5
266 ELSE IF( ncc.LT.0 ) THEN
267 info = -6
268 ELSE IF( ( ncvt.EQ.0 .AND. ldvt.LT.1 ) .OR.
269 $ ( ncvt.GT.0 .AND. ldvt.LT.max( 1, n ) ) ) THEN
270 info = -10
271 ELSE IF( ldu.LT.max( 1, nru ) ) THEN
272 info = -12
273 ELSE IF( ( ncc.EQ.0 .AND. ldc.LT.1 ) .OR.
274 $ ( ncc.GT.0 .AND. ldc.LT.max( 1, n ) ) ) THEN
275 info = -14
276 END IF
277 IF( info.NE.0 ) THEN
278 CALL xerbla( 'SLASDQ', -info )
279 RETURN
280 END IF
281 IF( n.EQ.0 )
282 $ RETURN
283*
284* ROTATE is true if any singular vectors desired, false otherwise
285*
286 rotate = ( ncvt.GT.0 ) .OR. ( nru.GT.0 ) .OR. ( ncc.GT.0 )
287 np1 = n + 1
288 sqre1 = sqre
289*
290* If matrix non-square upper bidiagonal, rotate to be lower
291* bidiagonal. The rotations are on the right.
292*
293 IF( ( iuplo.EQ.1 ) .AND. ( sqre1.EQ.1 ) ) THEN
294 DO 10 i = 1, n - 1
295 CALL slartg( d( i ), e( i ), cs, sn, r )
296 d( i ) = r
297 e( i ) = sn*d( i+1 )
298 d( i+1 ) = cs*d( i+1 )
299 IF( rotate ) THEN
300 work( i ) = cs
301 work( n+i ) = sn
302 END IF
303 10 CONTINUE
304 CALL slartg( d( n ), e( n ), cs, sn, r )
305 d( n ) = r
306 e( n ) = zero
307 IF( rotate ) THEN
308 work( n ) = cs
309 work( n+n ) = sn
310 END IF
311 iuplo = 2
312 sqre1 = 0
313*
314* Update singular vectors if desired.
315*
316 IF( ncvt.GT.0 )
317 $ CALL slasr( 'L', 'V', 'F', np1, ncvt, work( 1 ),
318 $ work( np1 ), vt, ldvt )
319 END IF
320*
321* If matrix lower bidiagonal, rotate to be upper bidiagonal
322* by applying Givens rotations on the left.
323*
324 IF( iuplo.EQ.2 ) THEN
325 DO 20 i = 1, n - 1
326 CALL slartg( d( i ), e( i ), cs, sn, r )
327 d( i ) = r
328 e( i ) = sn*d( i+1 )
329 d( i+1 ) = cs*d( i+1 )
330 IF( rotate ) THEN
331 work( i ) = cs
332 work( n+i ) = sn
333 END IF
334 20 CONTINUE
335*
336* If matrix (N+1)-by-N lower bidiagonal, one additional
337* rotation is needed.
338*
339 IF( sqre1.EQ.1 ) THEN
340 CALL slartg( d( n ), e( n ), cs, sn, r )
341 d( n ) = r
342 IF( rotate ) THEN
343 work( n ) = cs
344 work( n+n ) = sn
345 END IF
346 END IF
347*
348* Update singular vectors if desired.
349*
350 IF( nru.GT.0 ) THEN
351 IF( sqre1.EQ.0 ) THEN
352 CALL slasr( 'R', 'V', 'F', nru, n, work( 1 ),
353 $ work( np1 ), u, ldu )
354 ELSE
355 CALL slasr( 'R', 'V', 'F', nru, np1, work( 1 ),
356 $ work( np1 ), u, ldu )
357 END IF
358 END IF
359 IF( ncc.GT.0 ) THEN
360 IF( sqre1.EQ.0 ) THEN
361 CALL slasr( 'L', 'V', 'F', n, ncc, work( 1 ),
362 $ work( np1 ), c, ldc )
363 ELSE
364 CALL slasr( 'L', 'V', 'F', np1, ncc, work( 1 ),
365 $ work( np1 ), c, ldc )
366 END IF
367 END IF
368 END IF
369*
370* Call SBDSQR to compute the SVD of the reduced real
371* N-by-N upper bidiagonal matrix.
372*
373 CALL sbdsqr( 'U', n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c,
374 $ ldc, work, info )
375*
376* Sort the singular values into ascending order (insertion sort on
377* singular values, but only one transposition per singular vector)
378*
379 DO 40 i = 1, n
380*
381* Scan for smallest D(I).
382*
383 isub = i
384 smin = d( i )
385 DO 30 j = i + 1, n
386 IF( d( j ).LT.smin ) THEN
387 isub = j
388 smin = d( j )
389 END IF
390 30 CONTINUE
391 IF( isub.NE.i ) THEN
392*
393* Swap singular values and vectors.
394*
395 d( isub ) = d( i )
396 d( i ) = smin
397 IF( ncvt.GT.0 )
398 $ CALL sswap( ncvt, vt( isub, 1 ), ldvt, vt( i, 1 ),
399 $ ldvt )
400 IF( nru.GT.0 )
401 $ CALL sswap( nru, u( 1, isub ), 1, u( 1, i ), 1 )
402 IF( ncc.GT.0 )
403 $ CALL sswap( ncc, c( isub, 1 ), ldc, c( i, 1 ), ldc )
404 END IF
405 40 CONTINUE
406*
407 RETURN
408*
409* End of SLASDQ
410*
411 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, work, info)
SBDSQR
Definition sbdsqr.f:240
subroutine slartg(f, g, c, s, r)
SLARTG generates a plane rotation with real cosine and real sine.
Definition slartg.f90:111
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:210
subroutine slasr(side, pivot, direct, m, n, c, s, a, lda)
SLASR applies a sequence of plane rotations to a general rectangular matrix.
Definition slasr.f:197
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82