LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlasd1.f
Go to the documentation of this file.
1*> \brief \b DLASD1 computes the SVD of an upper bidiagonal matrix B of the specified size. Used by sbdsdc.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DLASD1 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasd1.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasd1.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasd1.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DLASD1( NL, NR, SQRE, D, ALPHA, BETA, U, LDU, VT, LDVT,
20* IDXQ, IWORK, WORK, INFO )
21*
22* .. Scalar Arguments ..
23* INTEGER INFO, LDU, LDVT, NL, NR, SQRE
24* DOUBLE PRECISION ALPHA, BETA
25* ..
26* .. Array Arguments ..
27* INTEGER IDXQ( * ), IWORK( * )
28* DOUBLE PRECISION D( * ), U( LDU, * ), VT( LDVT, * ), WORK( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DLASD1 computes the SVD of an upper bidiagonal N-by-M matrix B,
38*> where N = NL + NR + 1 and M = N + SQRE. DLASD1 is called from DLASD0.
39*>
40*> A related subroutine DLASD7 handles the case in which the singular
41*> values (and the singular vectors in factored form) are desired.
42*>
43*> DLASD1 computes the SVD as follows:
44*>
45*> ( D1(in) 0 0 0 )
46*> B = U(in) * ( Z1**T a Z2**T b ) * VT(in)
47*> ( 0 0 D2(in) 0 )
48*>
49*> = U(out) * ( D(out) 0) * VT(out)
50*>
51*> where Z**T = (Z1**T a Z2**T b) = u**T VT**T, and u is a vector of dimension M
52*> with ALPHA and BETA in the NL+1 and NL+2 th entries and zeros
53*> elsewhere; and the entry b is empty if SQRE = 0.
54*>
55*> The left singular vectors of the original matrix are stored in U, and
56*> the transpose of the right singular vectors are stored in VT, and the
57*> singular values are in D. The algorithm consists of three stages:
58*>
59*> The first stage consists of deflating the size of the problem
60*> when there are multiple singular values or when there are zeros in
61*> the Z vector. For each such occurrence the dimension of the
62*> secular equation problem is reduced by one. This stage is
63*> performed by the routine DLASD2.
64*>
65*> The second stage consists of calculating the updated
66*> singular values. This is done by finding the square roots of the
67*> roots of the secular equation via the routine DLASD4 (as called
68*> by DLASD3). This routine also calculates the singular vectors of
69*> the current problem.
70*>
71*> The final stage consists of computing the updated singular vectors
72*> directly using the updated singular values. The singular vectors
73*> for the current problem are multiplied with the singular vectors
74*> from the overall problem.
75*> \endverbatim
76*
77* Arguments:
78* ==========
79*
80*> \param[in] NL
81*> \verbatim
82*> NL is INTEGER
83*> The row dimension of the upper block. NL >= 1.
84*> \endverbatim
85*>
86*> \param[in] NR
87*> \verbatim
88*> NR is INTEGER
89*> The row dimension of the lower block. NR >= 1.
90*> \endverbatim
91*>
92*> \param[in] SQRE
93*> \verbatim
94*> SQRE is INTEGER
95*> = 0: the lower block is an NR-by-NR square matrix.
96*> = 1: the lower block is an NR-by-(NR+1) rectangular matrix.
97*>
98*> The bidiagonal matrix has row dimension N = NL + NR + 1,
99*> and column dimension M = N + SQRE.
100*> \endverbatim
101*>
102*> \param[in,out] D
103*> \verbatim
104*> D is DOUBLE PRECISION array,
105*> dimension (N = NL+NR+1).
106*> On entry D(1:NL,1:NL) contains the singular values of the
107*> upper block; and D(NL+2:N) contains the singular values of
108*> the lower block. On exit D(1:N) contains the singular values
109*> of the modified matrix.
110*> \endverbatim
111*>
112*> \param[in,out] ALPHA
113*> \verbatim
114*> ALPHA is DOUBLE PRECISION
115*> Contains the diagonal element associated with the added row.
116*> \endverbatim
117*>
118*> \param[in,out] BETA
119*> \verbatim
120*> BETA is DOUBLE PRECISION
121*> Contains the off-diagonal element associated with the added
122*> row.
123*> \endverbatim
124*>
125*> \param[in,out] U
126*> \verbatim
127*> U is DOUBLE PRECISION array, dimension(LDU,N)
128*> On entry U(1:NL, 1:NL) contains the left singular vectors of
129*> the upper block; U(NL+2:N, NL+2:N) contains the left singular
130*> vectors of the lower block. On exit U contains the left
131*> singular vectors of the bidiagonal matrix.
132*> \endverbatim
133*>
134*> \param[in] LDU
135*> \verbatim
136*> LDU is INTEGER
137*> The leading dimension of the array U. LDU >= max( 1, N ).
138*> \endverbatim
139*>
140*> \param[in,out] VT
141*> \verbatim
142*> VT is DOUBLE PRECISION array, dimension(LDVT,M)
143*> where M = N + SQRE.
144*> On entry VT(1:NL+1, 1:NL+1)**T contains the right singular
145*> vectors of the upper block; VT(NL+2:M, NL+2:M)**T contains
146*> the right singular vectors of the lower block. On exit
147*> VT**T contains the right singular vectors of the
148*> bidiagonal matrix.
149*> \endverbatim
150*>
151*> \param[in] LDVT
152*> \verbatim
153*> LDVT is INTEGER
154*> The leading dimension of the array VT. LDVT >= max( 1, M ).
155*> \endverbatim
156*>
157*> \param[in,out] IDXQ
158*> \verbatim
159*> IDXQ is INTEGER array, dimension(N)
160*> This contains the permutation which will reintegrate the
161*> subproblem just solved back into sorted order, i.e.
162*> D( IDXQ( I = 1, N ) ) will be in ascending order.
163*> \endverbatim
164*>
165*> \param[out] IWORK
166*> \verbatim
167*> IWORK is INTEGER array, dimension( 4 * N )
168*> \endverbatim
169*>
170*> \param[out] WORK
171*> \verbatim
172*> WORK is DOUBLE PRECISION array, dimension( 3*M**2 + 2*M )
173*> \endverbatim
174*>
175*> \param[out] INFO
176*> \verbatim
177*> INFO is INTEGER
178*> = 0: successful exit.
179*> < 0: if INFO = -i, the i-th argument had an illegal value.
180*> > 0: if INFO = 1, a singular value did not converge
181*> \endverbatim
182*
183* Authors:
184* ========
185*
186*> \author Univ. of Tennessee
187*> \author Univ. of California Berkeley
188*> \author Univ. of Colorado Denver
189*> \author NAG Ltd.
190*
191*> \ingroup lasd1
192*
193*> \par Contributors:
194* ==================
195*>
196*> Ming Gu and Huan Ren, Computer Science Division, University of
197*> California at Berkeley, USA
198*>
199* =====================================================================
200 SUBROUTINE dlasd1( NL, NR, SQRE, D, ALPHA, BETA, U, LDU, VT,
201 $ LDVT,
202 $ IDXQ, IWORK, WORK, INFO )
203*
204* -- LAPACK auxiliary routine --
205* -- LAPACK is a software package provided by Univ. of Tennessee, --
206* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
207*
208* .. Scalar Arguments ..
209 INTEGER INFO, LDU, LDVT, NL, NR, SQRE
210 DOUBLE PRECISION ALPHA, BETA
211* ..
212* .. Array Arguments ..
213 INTEGER IDXQ( * ), IWORK( * )
214 DOUBLE PRECISION D( * ), U( LDU, * ), VT( LDVT, * ), WORK( * )
215* ..
216*
217* =====================================================================
218*
219* .. Parameters ..
220*
221 DOUBLE PRECISION ONE, ZERO
222 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
223* ..
224* .. Local Scalars ..
225 INTEGER COLTYP, I, IDX, IDXC, IDXP, IQ, ISIGMA, IU2,
226 $ IVT2, IZ, K, LDQ, LDU2, LDVT2, M, N, N1, N2
227 DOUBLE PRECISION ORGNRM
228* ..
229* .. External Subroutines ..
230 EXTERNAL dlamrg, dlascl, dlasd2, dlasd3,
231 $ xerbla
232* ..
233* .. Intrinsic Functions ..
234 INTRINSIC abs, max
235* ..
236* .. Executable Statements ..
237*
238* Test the input parameters.
239*
240 info = 0
241*
242 IF( nl.LT.1 ) THEN
243 info = -1
244 ELSE IF( nr.LT.1 ) THEN
245 info = -2
246 ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
247 info = -3
248 END IF
249 IF( info.NE.0 ) THEN
250 CALL xerbla( 'DLASD1', -info )
251 RETURN
252 END IF
253*
254 n = nl + nr + 1
255 m = n + sqre
256*
257* The following values are for bookkeeping purposes only. They are
258* integer pointers which indicate the portion of the workspace
259* used by a particular array in DLASD2 and DLASD3.
260*
261 ldu2 = n
262 ldvt2 = m
263*
264 iz = 1
265 isigma = iz + m
266 iu2 = isigma + n
267 ivt2 = iu2 + ldu2*n
268 iq = ivt2 + ldvt2*m
269*
270 idx = 1
271 idxc = idx + n
272 coltyp = idxc + n
273 idxp = coltyp + n
274*
275* Scale.
276*
277 orgnrm = max( abs( alpha ), abs( beta ) )
278 d( nl+1 ) = zero
279 DO 10 i = 1, n
280 IF( abs( d( i ) ).GT.orgnrm ) THEN
281 orgnrm = abs( d( i ) )
282 END IF
283 10 CONTINUE
284 CALL dlascl( 'G', 0, 0, orgnrm, one, n, 1, d, n, info )
285 alpha = alpha / orgnrm
286 beta = beta / orgnrm
287*
288* Deflate singular values.
289*
290 CALL dlasd2( nl, nr, sqre, k, d, work( iz ), alpha, beta, u,
291 $ ldu,
292 $ vt, ldvt, work( isigma ), work( iu2 ), ldu2,
293 $ work( ivt2 ), ldvt2, iwork( idxp ), iwork( idx ),
294 $ iwork( idxc ), idxq, iwork( coltyp ), info )
295*
296* Solve Secular Equation and update singular vectors.
297*
298 ldq = k
299 CALL dlasd3( nl, nr, sqre, k, d, work( iq ), ldq,
300 $ work( isigma ),
301 $ u, ldu, work( iu2 ), ldu2, vt, ldvt, work( ivt2 ),
302 $ ldvt2, iwork( idxc ), iwork( coltyp ), work( iz ),
303 $ info )
304*
305* Report the convergence failure.
306*
307 IF( info.NE.0 ) THEN
308 RETURN
309 END IF
310*
311* Unscale.
312*
313 CALL dlascl( 'G', 0, 0, one, orgnrm, n, 1, d, n, info )
314*
315* Prepare the IDXQ sorting permutation.
316*
317 n1 = k
318 n2 = n - k
319 CALL dlamrg( n1, n2, d, 1, -1, idxq )
320*
321 RETURN
322*
323* End of DLASD1
324*
325 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlamrg(n1, n2, a, dtrd1, dtrd2, index)
DLAMRG creates a permutation list to merge the entries of two independently sorted sets into a single...
Definition dlamrg.f:97
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:142
subroutine dlasd1(nl, nr, sqre, d, alpha, beta, u, ldu, vt, ldvt, idxq, iwork, work, info)
DLASD1 computes the SVD of an upper bidiagonal matrix B of the specified size. Used by sbdsdc.
Definition dlasd1.f:203
subroutine dlasd2(nl, nr, sqre, k, d, z, alpha, beta, u, ldu, vt, ldvt, dsigma, u2, ldu2, vt2, ldvt2, idxp, idx, idxc, idxq, coltyp, info)
DLASD2 merges the two sets of singular values together into a single sorted set. Used by sbdsdc.
Definition dlasd2.f:268
subroutine dlasd3(nl, nr, sqre, k, d, q, ldq, dsigma, u, ldu, u2, ldu2, vt, ldvt, vt2, ldvt2, idxc, ctot, z, info)
DLASD3 finds all square roots of the roots of the secular equation, as defined by the values in D and...
Definition dlasd3.f:216