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

DLASD0 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 DLASD0 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 Using a divide and conquer approach, DLASD0 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, DLASDA, 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (M-1)
         Contains the subdiagonal entries of the bidiagonal matrix.
         On exit, E has been destroyed.
[out]U
          U is DOUBLE PRECISION 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 DOUBLE PRECISION 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 work array.
         Dimension must be at least (8 * N)
[out]WORK
          WORK is DOUBLE PRECISION work array.
         Dimension must be at least (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 154 of file dlasd0.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: