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

◆ sgghrd()

subroutine sgghrd ( character  compq,
character  compz,
integer  n,
integer  ilo,
integer  ihi,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldb, * )  b,
integer  ldb,
real, dimension( ldq, * )  q,
integer  ldq,
real, dimension( ldz, * )  z,
integer  ldz,
integer  info 
)

SGGHRD

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

Purpose:
 SGGHRD reduces a pair of real matrices (A,B) to generalized upper
 Hessenberg form using orthogonal transformations, where A is a
 general matrix and B is upper triangular.  The form of the
 generalized eigenvalue problem is
    A*x = lambda*B*x,
 and B is typically made upper triangular by computing its QR
 factorization and moving the orthogonal matrix Q to the left side
 of the equation.

 This subroutine simultaneously reduces A to a Hessenberg matrix H:
    Q**T*A*Z = H
 and transforms B to another upper triangular matrix T:
    Q**T*B*Z = T
 in order to reduce the problem to its standard form
    H*y = lambda*T*y
 where y = Z**T*x.

 The orthogonal matrices Q and Z are determined as products of Givens
 rotations.  They may either be formed explicitly, or they may be
 postmultiplied into input matrices Q1 and Z1, so that

      Q1 * A * Z1**T = (Q1*Q) * H * (Z1*Z)**T

      Q1 * B * Z1**T = (Q1*Q) * T * (Z1*Z)**T

 If Q1 is the orthogonal matrix from the QR factorization of B in the
 original equation A*x = lambda*B*x, then SGGHRD reduces the original
 problem to generalized Hessenberg form.
Parameters
[in]COMPQ
          COMPQ is CHARACTER*1
          = 'N': do not compute Q;
          = 'I': Q is initialized to the unit matrix, and the
                 orthogonal matrix Q is returned;
          = 'V': Q must contain an orthogonal matrix Q1 on entry,
                 and the product Q1*Q is returned.
[in]COMPZ
          COMPZ is CHARACTER*1
          = 'N': do not compute Z;
          = 'I': Z is initialized to the unit matrix, and the
                 orthogonal matrix Z is returned;
          = 'V': Z must contain an orthogonal matrix Z1 on entry,
                 and the product Z1*Z is returned.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in]ILO
          ILO is INTEGER
[in]IHI
          IHI is INTEGER

          ILO and IHI mark the rows and columns of A which are to be
          reduced.  It is assumed that A is already upper triangular
          in rows and columns 1:ILO-1 and IHI+1:N.  ILO and IHI are
          normally set by a previous call to SGGBAL; otherwise they
          should be set to 1 and N respectively.
          1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.
[in,out]A
          A is REAL array, dimension (LDA, N)
          On entry, the N-by-N general matrix to be reduced.
          On exit, the upper triangle and the first subdiagonal of A
          are overwritten with the upper Hessenberg matrix H, and the
          rest is set to zero.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension (LDB, N)
          On entry, the N-by-N upper triangular matrix B.
          On exit, the upper triangular matrix T = Q**T B Z.  The
          elements below the diagonal are set to zero.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in,out]Q
          Q is REAL array, dimension (LDQ, N)
          On entry, if COMPQ = 'V', the orthogonal matrix Q1,
          typically from the QR factorization of B.
          On exit, if COMPQ='I', the orthogonal matrix Q, and if
          COMPQ = 'V', the product Q1*Q.
          Not referenced if COMPQ='N'.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.
          LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.
[in,out]Z
          Z is REAL array, dimension (LDZ, N)
          On entry, if COMPZ = 'V', the orthogonal matrix Z1.
          On exit, if COMPZ='I', the orthogonal matrix Z, and if
          COMPZ = 'V', the product Z1*Z.
          Not referenced if COMPZ='N'.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.
          LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.
[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.
Further Details:
  This routine reduces A to Hessenberg and B to triangular form by
  an unblocked reduction, as described in _Matrix_Computations_,
  by Golub and Van Loan (Johns Hopkins Press.)

Definition at line 205 of file sgghrd.f.

207*
208* -- LAPACK computational routine --
209* -- LAPACK is a software package provided by Univ. of Tennessee, --
210* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
211*
212* .. Scalar Arguments ..
213 CHARACTER COMPQ, COMPZ
214 INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, N
215* ..
216* .. Array Arguments ..
217 REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
218 $ Z( LDZ, * )
219* ..
220*
221* =====================================================================
222*
223* .. Parameters ..
224 REAL ONE, ZERO
225 parameter( one = 1.0e+0, zero = 0.0e+0 )
226* ..
227* .. Local Scalars ..
228 LOGICAL ILQ, ILZ
229 INTEGER ICOMPQ, ICOMPZ, JCOL, JROW
230 REAL C, S, TEMP
231* ..
232* .. External Functions ..
233 LOGICAL LSAME
234 EXTERNAL lsame
235* ..
236* .. External Subroutines ..
237 EXTERNAL slartg, slaset, srot, xerbla
238* ..
239* .. Intrinsic Functions ..
240 INTRINSIC max
241* ..
242* .. Executable Statements ..
243*
244* Decode COMPQ
245*
246 IF( lsame( compq, 'N' ) ) THEN
247 ilq = .false.
248 icompq = 1
249 ELSE IF( lsame( compq, 'V' ) ) THEN
250 ilq = .true.
251 icompq = 2
252 ELSE IF( lsame( compq, 'I' ) ) THEN
253 ilq = .true.
254 icompq = 3
255 ELSE
256 icompq = 0
257 END IF
258*
259* Decode COMPZ
260*
261 IF( lsame( compz, 'N' ) ) THEN
262 ilz = .false.
263 icompz = 1
264 ELSE IF( lsame( compz, 'V' ) ) THEN
265 ilz = .true.
266 icompz = 2
267 ELSE IF( lsame( compz, 'I' ) ) THEN
268 ilz = .true.
269 icompz = 3
270 ELSE
271 icompz = 0
272 END IF
273*
274* Test the input parameters.
275*
276 info = 0
277 IF( icompq.LE.0 ) THEN
278 info = -1
279 ELSE IF( icompz.LE.0 ) THEN
280 info = -2
281 ELSE IF( n.LT.0 ) THEN
282 info = -3
283 ELSE IF( ilo.LT.1 ) THEN
284 info = -4
285 ELSE IF( ihi.GT.n .OR. ihi.LT.ilo-1 ) THEN
286 info = -5
287 ELSE IF( lda.LT.max( 1, n ) ) THEN
288 info = -7
289 ELSE IF( ldb.LT.max( 1, n ) ) THEN
290 info = -9
291 ELSE IF( ( ilq .AND. ldq.LT.n ) .OR. ldq.LT.1 ) THEN
292 info = -11
293 ELSE IF( ( ilz .AND. ldz.LT.n ) .OR. ldz.LT.1 ) THEN
294 info = -13
295 END IF
296 IF( info.NE.0 ) THEN
297 CALL xerbla( 'SGGHRD', -info )
298 RETURN
299 END IF
300*
301* Initialize Q and Z if desired.
302*
303 IF( icompq.EQ.3 )
304 $ CALL slaset( 'Full', n, n, zero, one, q, ldq )
305 IF( icompz.EQ.3 )
306 $ CALL slaset( 'Full', n, n, zero, one, z, ldz )
307*
308* Quick return if possible
309*
310 IF( n.LE.1 )
311 $ RETURN
312*
313* Zero out lower triangle of B
314*
315 DO 20 jcol = 1, n - 1
316 DO 10 jrow = jcol + 1, n
317 b( jrow, jcol ) = zero
318 10 CONTINUE
319 20 CONTINUE
320*
321* Reduce A and B
322*
323 DO 40 jcol = ilo, ihi - 2
324*
325 DO 30 jrow = ihi, jcol + 2, -1
326*
327* Step 1: rotate rows JROW-1, JROW to kill A(JROW,JCOL)
328*
329 temp = a( jrow-1, jcol )
330 CALL slartg( temp, a( jrow, jcol ), c, s,
331 $ a( jrow-1, jcol ) )
332 a( jrow, jcol ) = zero
333 CALL srot( n-jcol, a( jrow-1, jcol+1 ), lda,
334 $ a( jrow, jcol+1 ), lda, c, s )
335 CALL srot( n+2-jrow, b( jrow-1, jrow-1 ), ldb,
336 $ b( jrow, jrow-1 ), ldb, c, s )
337 IF( ilq )
338 $ CALL srot( n, q( 1, jrow-1 ), 1, q( 1, jrow ), 1, c, s )
339*
340* Step 2: rotate columns JROW, JROW-1 to kill B(JROW,JROW-1)
341*
342 temp = b( jrow, jrow )
343 CALL slartg( temp, b( jrow, jrow-1 ), c, s,
344 $ b( jrow, jrow ) )
345 b( jrow, jrow-1 ) = zero
346 CALL srot( ihi, a( 1, jrow ), 1, a( 1, jrow-1 ), 1, c, s )
347 CALL srot( jrow-1, b( 1, jrow ), 1, b( 1, jrow-1 ), 1, c,
348 $ s )
349 IF( ilz )
350 $ CALL srot( n, z( 1, jrow ), 1, z( 1, jrow-1 ), 1, c, s )
351 30 CONTINUE
352 40 CONTINUE
353*
354 RETURN
355*
356* End of SGGHRD
357*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine slartg(f, g, c, s, r)
SLARTG generates a plane rotation with real cosine and real sine.
Definition slartg.f90:111
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:110
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: