LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine slasd0 ( integer  N,
integer  SQRE,
real, dimension( * )  D,
real, dimension( * )  E,
real, dimension( ldu, * )  U,
integer  LDU,
real, dimension( ldvt, * )  VT,
integer  LDVT,
integer  SMLSIZ,
integer, dimension( * )  IWORK,
real, dimension( * )  WORK,
integer  INFO 
)

SLASD0 computes the singular values of a real upper bidiagonal n-by-m matrix B with diagonal d and off-diagonal e. Used by sbdsdc.

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

Purpose:
 Using a divide and conquer approach, SLASD0 computes the singular
 value decomposition (SVD) of a real upper bidiagonal N-by-M
 matrix B with diagonal D and offdiagonal E, where M = N + SQRE.
 The algorithm computes orthogonal matrices U and VT such that
 B = U * S * VT. The singular values S are overwritten on D.

 A related subroutine, SLASDA, computes only the singular values,
 and optionally, the singular vectors in compact form.
Parameters
[in]N
          N is INTEGER
         On entry, the row dimension of the upper bidiagonal matrix.
         This is also the dimension of the main diagonal array D.
[in]SQRE
          SQRE is INTEGER
         Specifies the column dimension of the bidiagonal matrix.
         = 0: The bidiagonal matrix has column dimension M = N;
         = 1: The bidiagonal matrix has column dimension M = N+1;
[in,out]D
          D is REAL array, dimension (N)
         On entry D contains the main diagonal of the bidiagonal
         matrix.
         On exit D, if INFO = 0, contains its singular values.
[in]E
          E is REAL array, dimension (M-1)
         Contains the subdiagonal entries of the bidiagonal matrix.
         On exit, E has been destroyed.
[out]U
          U is REAL array, dimension at least (LDQ, N)
         On exit, U contains the left singular vectors.
[in]LDU
          LDU is INTEGER
         On entry, leading dimension of U.
[out]VT
          VT is REAL array, dimension at least (LDVT, M)
         On exit, VT**T contains the right singular vectors.
[in]LDVT
          LDVT is INTEGER
         On entry, leading dimension of VT.
[in]SMLSIZ
          SMLSIZ is INTEGER
         On entry, maximum size of the subproblems at the
         bottom of the computation tree.
[out]IWORK
          IWORK is INTEGER array, dimension (8*N)
[out]WORK
          WORK is REAL array, dimension (3*M**2+2*M)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = 1, a singular value did not converge
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015
Contributors:
Ming Gu and Huan Ren, Computer Science Division, University of California at Berkeley, USA

Definition at line 152 of file slasd0.f.

152 *
153 * -- LAPACK auxiliary routine (version 3.6.0) --
154 * -- LAPACK is a software package provided by Univ. of Tennessee, --
155 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
156 * November 2015
157 *
158 * .. Scalar Arguments ..
159  INTEGER info, ldu, ldvt, n, smlsiz, sqre
160 * ..
161 * .. Array Arguments ..
162  INTEGER iwork( * )
163  REAL d( * ), e( * ), u( ldu, * ), vt( ldvt, * ),
164  $ work( * )
165 * ..
166 *
167 * =====================================================================
168 *
169 * .. Local Scalars ..
170  INTEGER i, i1, ic, idxq, idxqc, im1, inode, itemp, iwk,
171  $ j, lf, ll, lvl, m, ncc, nd, ndb1, ndiml, ndimr,
172  $ nl, nlf, nlp1, nlvl, nr, nrf, nrp1, sqrei
173  REAL alpha, beta
174 * ..
175 * .. External Subroutines ..
176  EXTERNAL slasd1, slasdq, slasdt, xerbla
177 * ..
178 * .. Executable Statements ..
179 *
180 * Test the input parameters.
181 *
182  info = 0
183 *
184  IF( n.LT.0 ) THEN
185  info = -1
186  ELSE IF( ( sqre.LT.0 ) .OR. ( sqre.GT.1 ) ) THEN
187  info = -2
188  END IF
189 *
190  m = n + sqre
191 *
192  IF( ldu.LT.n ) THEN
193  info = -6
194  ELSE IF( ldvt.LT.m ) THEN
195  info = -8
196  ELSE IF( smlsiz.LT.3 ) THEN
197  info = -9
198  END IF
199  IF( info.NE.0 ) THEN
200  CALL xerbla( 'SLASD0', -info )
201  RETURN
202  END IF
203 *
204 * If the input matrix is too small, call SLASDQ to find the SVD.
205 *
206  IF( n.LE.smlsiz ) THEN
207  CALL slasdq( 'U', sqre, n, m, n, 0, d, e, vt, ldvt, u, ldu, u,
208  $ ldu, work, info )
209  RETURN
210  END IF
211 *
212 * Set up the computation tree.
213 *
214  inode = 1
215  ndiml = inode + n
216  ndimr = ndiml + n
217  idxq = ndimr + n
218  iwk = idxq + n
219  CALL slasdt( n, nlvl, nd, iwork( inode ), iwork( ndiml ),
220  $ iwork( ndimr ), smlsiz )
221 *
222 * For the nodes on bottom level of the tree, solve
223 * their subproblems by SLASDQ.
224 *
225  ndb1 = ( nd+1 ) / 2
226  ncc = 0
227  DO 30 i = ndb1, nd
228 *
229 * IC : center row of each node
230 * NL : number of rows of left subproblem
231 * NR : number of rows of right subproblem
232 * NLF: starting row of the left subproblem
233 * NRF: starting row of the right subproblem
234 *
235  i1 = i - 1
236  ic = iwork( inode+i1 )
237  nl = iwork( ndiml+i1 )
238  nlp1 = nl + 1
239  nr = iwork( ndimr+i1 )
240  nrp1 = nr + 1
241  nlf = ic - nl
242  nrf = ic + 1
243  sqrei = 1
244  CALL slasdq( 'U', sqrei, nl, nlp1, nl, ncc, d( nlf ), e( nlf ),
245  $ vt( nlf, nlf ), ldvt, u( nlf, nlf ), ldu,
246  $ u( nlf, nlf ), ldu, work, info )
247  IF( info.NE.0 ) THEN
248  RETURN
249  END IF
250  itemp = idxq + nlf - 2
251  DO 10 j = 1, nl
252  iwork( itemp+j ) = j
253  10 CONTINUE
254  IF( i.EQ.nd ) THEN
255  sqrei = sqre
256  ELSE
257  sqrei = 1
258  END IF
259  nrp1 = nr + sqrei
260  CALL slasdq( 'U', sqrei, nr, nrp1, nr, ncc, d( nrf ), e( nrf ),
261  $ vt( nrf, nrf ), ldvt, u( nrf, nrf ), ldu,
262  $ u( nrf, nrf ), ldu, work, info )
263  IF( info.NE.0 ) THEN
264  RETURN
265  END IF
266  itemp = idxq + ic
267  DO 20 j = 1, nr
268  iwork( itemp+j-1 ) = j
269  20 CONTINUE
270  30 CONTINUE
271 *
272 * Now conquer each subproblem bottom-up.
273 *
274  DO 50 lvl = nlvl, 1, -1
275 *
276 * Find the first node LF and last node LL on the
277 * current level LVL.
278 *
279  IF( lvl.EQ.1 ) THEN
280  lf = 1
281  ll = 1
282  ELSE
283  lf = 2**( lvl-1 )
284  ll = 2*lf - 1
285  END IF
286  DO 40 i = lf, ll
287  im1 = i - 1
288  ic = iwork( inode+im1 )
289  nl = iwork( ndiml+im1 )
290  nr = iwork( ndimr+im1 )
291  nlf = ic - nl
292  IF( ( sqre.EQ.0 ) .AND. ( i.EQ.ll ) ) THEN
293  sqrei = sqre
294  ELSE
295  sqrei = 1
296  END IF
297  idxqc = idxq + nlf - 1
298  alpha = d( ic )
299  beta = e( ic )
300  CALL slasd1( nl, nr, sqrei, d( nlf ), alpha, beta,
301  $ u( nlf, nlf ), ldu, vt( nlf, nlf ), ldvt,
302  $ iwork( idxqc ), iwork( iwk ), work, info )
303 *
304 * Report the possible convergence failure.
305 *
306  IF( info.NE.0 ) THEN
307  RETURN
308  END IF
309  40 CONTINUE
310  50 CONTINUE
311 *
312  RETURN
313 *
314 * End of SLASD0
315 *
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:107
subroutine slasd1(NL, NR, SQRE, D, ALPHA, BETA, U, LDU, VT, LDVT, IDXQ, IWORK, WORK, INFO)
SLASD1 computes the SVD of an upper bidiagonal matrix B of the specified size. Used by sbdsdc...
Definition: slasd1.f:206
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
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:213

Here is the call graph for this function:

Here is the caller graph for this function: