LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ zgebrd()

subroutine zgebrd ( integer  M,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
complex*16, dimension( * )  TAUQ,
complex*16, dimension( * )  TAUP,
complex*16, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

ZGEBRD

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

Purpose:
 ZGEBRD 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*16 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 DOUBLE PRECISION array, dimension (min(M,N))
          The diagonal elements of the bidiagonal matrix B:
          D(i) = A(i,i).
[out]E
          E is DOUBLE PRECISION 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*16 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*16 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*16 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 2017
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 207 of file zgebrd.f.

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