LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgbbrd()

subroutine sgbbrd ( character vect,
integer m,
integer n,
integer ncc,
integer kl,
integer ku,
real, dimension( ldab, * ) ab,
integer ldab,
real, dimension( * ) d,
real, dimension( * ) e,
real, dimension( ldq, * ) q,
integer ldq,
real, dimension( ldpt, * ) pt,
integer ldpt,
real, dimension( ldc, * ) c,
integer ldc,
real, dimension( * ) work,
integer info )

SGBBRD

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

Purpose:
!>
!> SGBBRD reduces a real general m-by-n band matrix A to upper
!> bidiagonal form B by an orthogonal transformation: Q**T * A * P = B.
!>
!> The routine computes B, and optionally forms Q or P**T, or computes
!> Q**T*C for a given matrix C.
!> 
Parameters
[in]VECT
!>          VECT is CHARACTER*1
!>          Specifies whether or not the matrices Q and P**T are to be
!>          formed.
!>          = 'N': do not form Q or P**T;
!>          = 'Q': form Q only;
!>          = 'P': form P**T only;
!>          = 'B': form both.
!> 
[in]M
!>          M is INTEGER
!>          The number of rows of the matrix A.  M >= 0.
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrix A.  N >= 0.
!> 
[in]NCC
!>          NCC is INTEGER
!>          The number of columns of the matrix C.  NCC >= 0.
!> 
[in]KL
!>          KL is INTEGER
!>          The number of subdiagonals of the matrix A. KL >= 0.
!> 
[in]KU
!>          KU is INTEGER
!>          The number of superdiagonals of the matrix A. KU >= 0.
!> 
[in,out]AB
!>          AB is REAL array, dimension (LDAB,N)
!>          On entry, the m-by-n band matrix A, stored in rows 1 to
!>          KL+KU+1. The j-th column of A is stored in the j-th column of
!>          the array AB as follows:
!>          AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku)<=i<=min(m,j+kl).
!>          On exit, A is overwritten by values generated during the
!>          reduction.
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array A. LDAB >= KL+KU+1.
!> 
[out]D
!>          D is REAL array, dimension (min(M,N))
!>          The diagonal elements of the bidiagonal matrix B.
!> 
[out]E
!>          E is REAL array, dimension (min(M,N)-1)
!>          The superdiagonal elements of the bidiagonal matrix B.
!> 
[out]Q
!>          Q is REAL array, dimension (LDQ,M)
!>          If VECT = 'Q' or 'B', the m-by-m orthogonal matrix Q.
!>          If VECT = 'N' or 'P', the array Q is not referenced.
!> 
[in]LDQ
!>          LDQ is INTEGER
!>          The leading dimension of the array Q.
!>          LDQ >= max(1,M) if VECT = 'Q' or 'B'; LDQ >= 1 otherwise.
!> 
[out]PT
!>          PT is REAL array, dimension (LDPT,N)
!>          If VECT = 'P' or 'B', the n-by-n orthogonal matrix P'.
!>          If VECT = 'N' or 'Q', the array PT is not referenced.
!> 
[in]LDPT
!>          LDPT is INTEGER
!>          The leading dimension of the array PT.
!>          LDPT >= max(1,N) if VECT = 'P' or 'B'; LDPT >= 1 otherwise.
!> 
[in,out]C
!>          C is REAL array, dimension (LDC,NCC)
!>          On entry, an m-by-ncc matrix C.
!>          On exit, C is overwritten by Q**T*C.
!>          C is not referenced if NCC = 0.
!> 
[in]LDC
!>          LDC is INTEGER
!>          The leading dimension of the array C.
!>          LDC >= max(1,M) if NCC > 0; LDC >= 1 if NCC = 0.
!> 
[out]WORK
!>          WORK is REAL array, dimension (2*max(M,N))
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit.
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 183 of file sgbbrd.f.

185*
186* -- LAPACK computational routine --
187* -- LAPACK is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 CHARACTER VECT
192 INTEGER INFO, KL, KU, LDAB, LDC, LDPT, LDQ, M, N, NCC
193* ..
194* .. Array Arguments ..
195 REAL AB( LDAB, * ), C( LDC, * ), D( * ), E( * ),
196 $ PT( LDPT, * ), Q( LDQ, * ), WORK( * )
197* ..
198*
199* =====================================================================
200*
201* .. Parameters ..
202 REAL ZERO, ONE
203 parameter( zero = 0.0e+0, one = 1.0e+0 )
204* ..
205* .. Local Scalars ..
206 LOGICAL WANTB, WANTC, WANTPT, WANTQ
207 INTEGER I, INCA, J, J1, J2, KB, KB1, KK, KLM, KLU1,
208 $ KUN, L, MINMN, ML, ML0, MN, MU, MU0, NR, NRT
209 REAL RA, RB, RC, RS
210* ..
211* .. External Subroutines ..
212 EXTERNAL slargv, slartg, slartv, slaset, srot,
213 $ xerbla
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max, min
217* ..
218* .. External Functions ..
219 LOGICAL LSAME
220 EXTERNAL lsame
221* ..
222* .. Executable Statements ..
223*
224* Test the input parameters
225*
226 wantb = lsame( vect, 'B' )
227 wantq = lsame( vect, 'Q' ) .OR. wantb
228 wantpt = lsame( vect, 'P' ) .OR. wantb
229 wantc = ncc.GT.0
230 klu1 = kl + ku + 1
231 info = 0
232 IF( .NOT.wantq .AND. .NOT.wantpt .AND.
233 $ .NOT.lsame( vect, 'N' ) )
234 $ THEN
235 info = -1
236 ELSE IF( m.LT.0 ) THEN
237 info = -2
238 ELSE IF( n.LT.0 ) THEN
239 info = -3
240 ELSE IF( ncc.LT.0 ) THEN
241 info = -4
242 ELSE IF( kl.LT.0 ) THEN
243 info = -5
244 ELSE IF( ku.LT.0 ) THEN
245 info = -6
246 ELSE IF( ldab.LT.klu1 ) THEN
247 info = -8
248 ELSE IF( ldq.LT.1 .OR. wantq .AND. ldq.LT.max( 1, m ) ) THEN
249 info = -12
250 ELSE IF( ldpt.LT.1 .OR. wantpt .AND. ldpt.LT.max( 1, n ) ) THEN
251 info = -14
252 ELSE IF( ldc.LT.1 .OR. wantc .AND. ldc.LT.max( 1, m ) ) THEN
253 info = -16
254 END IF
255 IF( info.NE.0 ) THEN
256 CALL xerbla( 'SGBBRD', -info )
257 RETURN
258 END IF
259*
260* Initialize Q and P**T to the unit matrix, if needed
261*
262 IF( wantq )
263 $ CALL slaset( 'Full', m, m, zero, one, q, ldq )
264 IF( wantpt )
265 $ CALL slaset( 'Full', n, n, zero, one, pt, ldpt )
266*
267* Quick return if possible.
268*
269 IF( m.EQ.0 .OR. n.EQ.0 )
270 $ RETURN
271*
272 minmn = min( m, n )
273*
274 IF( kl+ku.GT.1 ) THEN
275*
276* Reduce to upper bidiagonal form if KU > 0; if KU = 0, reduce
277* first to lower bidiagonal form and then transform to upper
278* bidiagonal
279*
280 IF( ku.GT.0 ) THEN
281 ml0 = 1
282 mu0 = 2
283 ELSE
284 ml0 = 2
285 mu0 = 1
286 END IF
287*
288* Wherever possible, plane rotations are generated and applied in
289* vector operations of length NR over the index set J1:J2:KLU1.
290*
291* The sines of the plane rotations are stored in WORK(1:max(m,n))
292* and the cosines in WORK(max(m,n)+1:2*max(m,n)).
293*
294 mn = max( m, n )
295 klm = min( m-1, kl )
296 kun = min( n-1, ku )
297 kb = klm + kun
298 kb1 = kb + 1
299 inca = kb1*ldab
300 nr = 0
301 j1 = klm + 2
302 j2 = 1 - kun
303*
304 DO 90 i = 1, minmn
305*
306* Reduce i-th column and i-th row of matrix to bidiagonal form
307*
308 ml = klm + 1
309 mu = kun + 1
310 DO 80 kk = 1, kb
311 j1 = j1 + kb
312 j2 = j2 + kb
313*
314* generate plane rotations to annihilate nonzero elements
315* which have been created below the band
316*
317 IF( nr.GT.0 )
318 $ CALL slargv( nr, ab( klu1, j1-klm-1 ), inca,
319 $ work( j1 ), kb1, work( mn+j1 ), kb1 )
320*
321* apply plane rotations from the left
322*
323 DO 10 l = 1, kb
324 IF( j2-klm+l-1.GT.n ) THEN
325 nrt = nr - 1
326 ELSE
327 nrt = nr
328 END IF
329 IF( nrt.GT.0 )
330 $ CALL slartv( nrt, ab( klu1-l, j1-klm+l-1 ),
331 $ inca,
332 $ ab( klu1-l+1, j1-klm+l-1 ), inca,
333 $ work( mn+j1 ), work( j1 ), kb1 )
334 10 CONTINUE
335*
336 IF( ml.GT.ml0 ) THEN
337 IF( ml.LE.m-i+1 ) THEN
338*
339* generate plane rotation to annihilate a(i+ml-1,i)
340* within the band, and apply rotation from the left
341*
342 CALL slartg( ab( ku+ml-1, i ), ab( ku+ml, i ),
343 $ work( mn+i+ml-1 ), work( i+ml-1 ),
344 $ ra )
345 ab( ku+ml-1, i ) = ra
346 IF( i.LT.n )
347 $ CALL srot( min( ku+ml-2, n-i ),
348 $ ab( ku+ml-2, i+1 ), ldab-1,
349 $ ab( ku+ml-1, i+1 ), ldab-1,
350 $ work( mn+i+ml-1 ), work( i+ml-1 ) )
351 END IF
352 nr = nr + 1
353 j1 = j1 - kb1
354 END IF
355*
356 IF( wantq ) THEN
357*
358* accumulate product of plane rotations in Q
359*
360 DO 20 j = j1, j2, kb1
361 CALL srot( m, q( 1, j-1 ), 1, q( 1, j ), 1,
362 $ work( mn+j ), work( j ) )
363 20 CONTINUE
364 END IF
365*
366 IF( wantc ) THEN
367*
368* apply plane rotations to C
369*
370 DO 30 j = j1, j2, kb1
371 CALL srot( ncc, c( j-1, 1 ), ldc, c( j, 1 ),
372 $ ldc,
373 $ work( mn+j ), work( j ) )
374 30 CONTINUE
375 END IF
376*
377 IF( j2+kun.GT.n ) THEN
378*
379* adjust J2 to keep within the bounds of the matrix
380*
381 nr = nr - 1
382 j2 = j2 - kb1
383 END IF
384*
385 DO 40 j = j1, j2, kb1
386*
387* create nonzero element a(j-1,j+ku) above the band
388* and store it in WORK(n+1:2*n)
389*
390 work( j+kun ) = work( j )*ab( 1, j+kun )
391 ab( 1, j+kun ) = work( mn+j )*ab( 1, j+kun )
392 40 CONTINUE
393*
394* generate plane rotations to annihilate nonzero elements
395* which have been generated above the band
396*
397 IF( nr.GT.0 )
398 $ CALL slargv( nr, ab( 1, j1+kun-1 ), inca,
399 $ work( j1+kun ), kb1, work( mn+j1+kun ),
400 $ kb1 )
401*
402* apply plane rotations from the right
403*
404 DO 50 l = 1, kb
405 IF( j2+l-1.GT.m ) THEN
406 nrt = nr - 1
407 ELSE
408 nrt = nr
409 END IF
410 IF( nrt.GT.0 )
411 $ CALL slartv( nrt, ab( l+1, j1+kun-1 ), inca,
412 $ ab( l, j1+kun ), inca,
413 $ work( mn+j1+kun ), work( j1+kun ),
414 $ kb1 )
415 50 CONTINUE
416*
417 IF( ml.EQ.ml0 .AND. mu.GT.mu0 ) THEN
418 IF( mu.LE.n-i+1 ) THEN
419*
420* generate plane rotation to annihilate a(i,i+mu-1)
421* within the band, and apply rotation from the right
422*
423 CALL slartg( ab( ku-mu+3, i+mu-2 ),
424 $ ab( ku-mu+2, i+mu-1 ),
425 $ work( mn+i+mu-1 ), work( i+mu-1 ),
426 $ ra )
427 ab( ku-mu+3, i+mu-2 ) = ra
428 CALL srot( min( kl+mu-2, m-i ),
429 $ ab( ku-mu+4, i+mu-2 ), 1,
430 $ ab( ku-mu+3, i+mu-1 ), 1,
431 $ work( mn+i+mu-1 ), work( i+mu-1 ) )
432 END IF
433 nr = nr + 1
434 j1 = j1 - kb1
435 END IF
436*
437 IF( wantpt ) THEN
438*
439* accumulate product of plane rotations in P**T
440*
441 DO 60 j = j1, j2, kb1
442 CALL srot( n, pt( j+kun-1, 1 ), ldpt,
443 $ pt( j+kun, 1 ), ldpt, work( mn+j+kun ),
444 $ work( j+kun ) )
445 60 CONTINUE
446 END IF
447*
448 IF( j2+kb.GT.m ) THEN
449*
450* adjust J2 to keep within the bounds of the matrix
451*
452 nr = nr - 1
453 j2 = j2 - kb1
454 END IF
455*
456 DO 70 j = j1, j2, kb1
457*
458* create nonzero element a(j+kl+ku,j+ku-1) below the
459* band and store it in WORK(1:n)
460*
461 work( j+kb ) = work( j+kun )*ab( klu1, j+kun )
462 ab( klu1, j+kun ) = work( mn+j+kun )*ab( klu1, j+kun )
463 70 CONTINUE
464*
465 IF( ml.GT.ml0 ) THEN
466 ml = ml - 1
467 ELSE
468 mu = mu - 1
469 END IF
470 80 CONTINUE
471 90 CONTINUE
472 END IF
473*
474 IF( ku.EQ.0 .AND. kl.GT.0 ) THEN
475*
476* A has been reduced to lower bidiagonal form
477*
478* Transform lower bidiagonal form to upper bidiagonal by applying
479* plane rotations from the left, storing diagonal elements in D
480* and off-diagonal elements in E
481*
482 DO 100 i = 1, min( m-1, n )
483 CALL slartg( ab( 1, i ), ab( 2, i ), rc, rs, ra )
484 d( i ) = ra
485 IF( i.LT.n ) THEN
486 e( i ) = rs*ab( 1, i+1 )
487 ab( 1, i+1 ) = rc*ab( 1, i+1 )
488 END IF
489 IF( wantq )
490 $ CALL srot( m, q( 1, i ), 1, q( 1, i+1 ), 1, rc, rs )
491 IF( wantc )
492 $ CALL srot( ncc, c( i, 1 ), ldc, c( i+1, 1 ), ldc, rc,
493 $ rs )
494 100 CONTINUE
495 IF( m.LE.n )
496 $ d( m ) = ab( 1, m )
497 ELSE IF( ku.GT.0 ) THEN
498*
499* A has been reduced to upper bidiagonal form
500*
501 IF( m.LT.n ) THEN
502*
503* Annihilate a(m,m+1) by applying plane rotations from the
504* right, storing diagonal elements in D and off-diagonal
505* elements in E
506*
507 rb = ab( ku, m+1 )
508 DO 110 i = m, 1, -1
509 CALL slartg( ab( ku+1, i ), rb, rc, rs, ra )
510 d( i ) = ra
511 IF( i.GT.1 ) THEN
512 rb = -rs*ab( ku, i )
513 e( i-1 ) = rc*ab( ku, i )
514 END IF
515 IF( wantpt )
516 $ CALL srot( n, pt( i, 1 ), ldpt, pt( m+1, 1 ), ldpt,
517 $ rc, rs )
518 110 CONTINUE
519 ELSE
520*
521* Copy off-diagonal elements to E and diagonal elements to D
522*
523 DO 120 i = 1, minmn - 1
524 e( i ) = ab( ku, i+1 )
525 120 CONTINUE
526 DO 130 i = 1, minmn
527 d( i ) = ab( ku+1, i )
528 130 CONTINUE
529 END IF
530 ELSE
531*
532* A is diagonal. Set elements of E to zero and copy diagonal
533* elements to D.
534*
535 DO 140 i = 1, minmn - 1
536 e( i ) = zero
537 140 CONTINUE
538 DO 150 i = 1, minmn
539 d( i ) = ab( 1, i )
540 150 CONTINUE
541 END IF
542 RETURN
543*
544* End of SGBBRD
545*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine slargv(n, x, incx, y, incy, c, incc)
SLARGV generates a vector of plane rotations with real cosines and real sines.
Definition slargv.f:102
subroutine slartg(f, g, c, s, r)
SLARTG generates a plane rotation with real cosine and real sine.
Definition slartg.f90:111
subroutine slartv(n, x, incx, y, incy, c, s, incc)
SLARTV applies a vector of plane rotations with real cosines and real sines to the elements of a pair...
Definition slartv.f:106
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:108
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine srot(n, sx, incx, sy, incy, c, s)
SROT
Definition srot.f:92
Here is the call graph for this function:
Here is the caller graph for this function: