LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cgebrd ( integer  M,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  D,
real, dimension( * )  E,
complex, dimension( * )  TAUQ,
complex, dimension( * )  TAUP,
complex, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

CGEBRD

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

Purpose:
 CGEBRD reduces a general complex M-by-N matrix A to upper or lower
 bidiagonal form B by a unitary transformation: Q**H * A * P = B.

 If m >= n, B is upper bidiagonal; if m < n, B is lower bidiagonal.
Parameters
[in]M
          M is INTEGER
          The number of rows in the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns in the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the M-by-N general matrix to be reduced.
          On exit,
          if m >= n, the diagonal and the first superdiagonal are
            overwritten with the upper bidiagonal matrix B; the
            elements below the diagonal, with the array TAUQ, represent
            the unitary matrix Q as a product of elementary
            reflectors, and the elements above the first superdiagonal,
            with the array TAUP, represent the unitary matrix P as
            a product of elementary reflectors;
          if m < n, the diagonal and the first subdiagonal are
            overwritten with the lower bidiagonal matrix B; the
            elements below the first subdiagonal, with the array TAUQ,
            represent the unitary matrix Q as a product of
            elementary reflectors, and the elements above the diagonal,
            with the array TAUP, represent the unitary matrix P as
            a product of elementary reflectors.
          See Further Details.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]D
          D is REAL array, dimension (min(M,N))
          The diagonal elements of the bidiagonal matrix B:
          D(i) = A(i,i).
[out]E
          E is REAL array, dimension (min(M,N)-1)
          The off-diagonal elements of the bidiagonal matrix B:
          if m >= n, E(i) = A(i,i+1) for i = 1,2,...,n-1;
          if m < n, E(i) = A(i+1,i) for i = 1,2,...,m-1.
[out]TAUQ
          TAUQ is COMPLEX array dimension (min(M,N))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix Q. See Further Details.
[out]TAUP
          TAUP is COMPLEX array, dimension (min(M,N))
          The scalar factors of the elementary reflectors which
          represent the unitary matrix P. See Further Details.
[out]WORK
          WORK is COMPLEX array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.  LWORK >= max(1,M,N).
          For optimum performance LWORK >= (M+N)*NB, where NB
          is the optimal blocksize.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[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.
Date
November 2011
Further Details:
  The matrices Q and P are represented as products of elementary
  reflectors:

  If m >= n,

     Q = H(1) H(2) . . . H(n)  and  P = G(1) G(2) . . . G(n-1)

  Each H(i) and G(i) has the form:

     H(i) = I - tauq * v * v**H  and G(i) = I - taup * u * u**H

  where tauq and taup are complex scalars, and v and u are complex
  vectors; v(1:i-1) = 0, v(i) = 1, and v(i+1:m) is stored on exit in
  A(i+1:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+2:n) is stored on exit in
  A(i,i+2:n); tauq is stored in TAUQ(i) and taup in TAUP(i).

  If m < n,

     Q = H(1) H(2) . . . H(m-1)  and  P = G(1) G(2) . . . G(m)

  Each H(i) and G(i) has the form:

     H(i) = I - tauq * v * v**H  and G(i) = I - taup * u * u**H

  where tauq and taup are complex scalars, and v and u are complex
  vectors; v(1:i) = 0, v(i+1) = 1, and v(i+2:m) is stored on exit in
  A(i+2:m,i); u(1:i-1) = 0, u(i) = 1, and u(i+1:n) is stored on exit in
  A(i,i+1:n); tauq is stored in TAUQ(i) and taup in TAUP(i).

  The contents of A on exit are illustrated by the following examples:

  m = 6 and n = 5 (m > n):          m = 5 and n = 6 (m < n):

    (  d   e   u1  u1  u1 )           (  d   u1  u1  u1  u1  u1 )
    (  v1  d   e   u2  u2 )           (  e   d   u2  u2  u2  u2 )
    (  v1  v2  d   e   u3 )           (  v1  e   d   u3  u3  u3 )
    (  v1  v2  v3  d   e  )           (  v1  v2  e   d   u4  u4 )
    (  v1  v2  v3  v4  d  )           (  v1  v2  v3  e   d   u5 )
    (  v1  v2  v3  v4  v5 )

  where d and e denote diagonal and off-diagonal elements of B, vi
  denotes an element of the vector defining H(i), and ui an element of
  the vector defining G(i).

Definition at line 208 of file cgebrd.f.

208 *
209 * -- LAPACK computational routine (version 3.4.0) --
210 * -- LAPACK is a software package provided by Univ. of Tennessee, --
211 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
212 * November 2011
213 *
214 * .. Scalar Arguments ..
215  INTEGER info, lda, lwork, m, n
216 * ..
217 * .. Array Arguments ..
218  REAL d( * ), e( * )
219  COMPLEX a( lda, * ), taup( * ), tauq( * ),
220  $ work( * )
221 * ..
222 *
223 * =====================================================================
224 *
225 * .. Parameters ..
226  COMPLEX one
227  parameter ( one = ( 1.0e+0, 0.0e+0 ) )
228 * ..
229 * .. Local Scalars ..
230  LOGICAL lquery
231  INTEGER i, iinfo, j, ldwrkx, ldwrky, lwkopt, minmn, nb,
232  $ nbmin, nx
233  REAL ws
234 * ..
235 * .. External Subroutines ..
236  EXTERNAL cgebd2, cgemm, clabrd, xerbla
237 * ..
238 * .. Intrinsic Functions ..
239  INTRINSIC max, min, real
240 * ..
241 * .. External Functions ..
242  INTEGER ilaenv
243  EXTERNAL ilaenv
244 * ..
245 * .. Executable Statements ..
246 *
247 * Test the input parameters
248 *
249  info = 0
250  nb = max( 1, ilaenv( 1, 'CGEBRD', ' ', m, n, -1, -1 ) )
251  lwkopt = ( m+n )*nb
252  work( 1 ) = REAL( lwkopt )
253  lquery = ( lwork.EQ.-1 )
254  IF( m.LT.0 ) THEN
255  info = -1
256  ELSE IF( n.LT.0 ) THEN
257  info = -2
258  ELSE IF( lda.LT.max( 1, m ) ) THEN
259  info = -4
260  ELSE IF( lwork.LT.max( 1, m, n ) .AND. .NOT.lquery ) THEN
261  info = -10
262  END IF
263  IF( info.LT.0 ) THEN
264  CALL xerbla( 'CGEBRD', -info )
265  RETURN
266  ELSE IF( lquery ) THEN
267  RETURN
268  END IF
269 *
270 * Quick return if possible
271 *
272  minmn = min( m, n )
273  IF( minmn.EQ.0 ) THEN
274  work( 1 ) = 1
275  RETURN
276  END IF
277 *
278  ws = max( m, n )
279  ldwrkx = m
280  ldwrky = n
281 *
282  IF( nb.GT.1 .AND. nb.LT.minmn ) THEN
283 *
284 * Set the crossover point NX.
285 *
286  nx = max( nb, ilaenv( 3, 'CGEBRD', ' ', m, n, -1, -1 ) )
287 *
288 * Determine when to switch from blocked to unblocked code.
289 *
290  IF( nx.LT.minmn ) THEN
291  ws = ( m+n )*nb
292  IF( lwork.LT.ws ) THEN
293 *
294 * Not enough work space for the optimal NB, consider using
295 * a smaller block size.
296 *
297  nbmin = ilaenv( 2, 'CGEBRD', ' ', m, n, -1, -1 )
298  IF( lwork.GE.( m+n )*nbmin ) THEN
299  nb = lwork / ( m+n )
300  ELSE
301  nb = 1
302  nx = minmn
303  END IF
304  END IF
305  END IF
306  ELSE
307  nx = minmn
308  END IF
309 *
310  DO 30 i = 1, minmn - nx, nb
311 *
312 * Reduce rows and columns i:i+ib-1 to bidiagonal form and return
313 * the matrices X and Y which are needed to update the unreduced
314 * part of the matrix
315 *
316  CALL clabrd( m-i+1, n-i+1, nb, a( i, i ), lda, d( i ), e( i ),
317  $ tauq( i ), taup( i ), work, ldwrkx,
318  $ work( ldwrkx*nb+1 ), ldwrky )
319 *
320 * Update the trailing submatrix A(i+ib:m,i+ib:n), using
321 * an update of the form A := A - V*Y**H - X*U**H
322 *
323  CALL cgemm( 'No transpose', 'Conjugate transpose', m-i-nb+1,
324  $ n-i-nb+1, nb, -one, a( i+nb, i ), lda,
325  $ work( ldwrkx*nb+nb+1 ), ldwrky, one,
326  $ a( i+nb, i+nb ), lda )
327  CALL cgemm( 'No transpose', 'No transpose', m-i-nb+1, n-i-nb+1,
328  $ nb, -one, work( nb+1 ), ldwrkx, a( i, i+nb ), lda,
329  $ one, a( i+nb, i+nb ), lda )
330 *
331 * Copy diagonal and off-diagonal elements of B back into A
332 *
333  IF( m.GE.n ) THEN
334  DO 10 j = i, i + nb - 1
335  a( j, j ) = d( j )
336  a( j, j+1 ) = e( j )
337  10 CONTINUE
338  ELSE
339  DO 20 j = i, i + nb - 1
340  a( j, j ) = d( j )
341  a( j+1, j ) = e( j )
342  20 CONTINUE
343  END IF
344  30 CONTINUE
345 *
346 * Use unblocked code to reduce the remainder of the matrix
347 *
348  CALL cgebd2( m-i+1, n-i+1, a( i, i ), lda, d( i ), e( i ),
349  $ tauq( i ), taup( i ), work, iinfo )
350  work( 1 ) = ws
351  RETURN
352 *
353 * End of CGEBRD
354 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine clabrd(M, N, NB, A, LDA, D, E, TAUQ, TAUP, X, LDX, Y, LDY)
CLABRD reduces the first nb rows and columns of a general matrix to a bidiagonal form.
Definition: clabrd.f:214
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
subroutine cgebd2(M, N, A, LDA, D, E, TAUQ, TAUP, WORK, INFO)
CGEBD2 reduces a general matrix to bidiagonal form using an unblocked algorithm.
Definition: cgebd2.f:192

Here is the call graph for this function:

Here is the caller graph for this function: