LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sspt21 ( integer  ITYPE,
character  UPLO,
integer  N,
integer  KBAND,
real, dimension( * )  AP,
real, dimension( * )  D,
real, dimension( * )  E,
real, dimension( ldu, * )  U,
integer  LDU,
real, dimension( * )  VP,
real, dimension( * )  TAU,
real, dimension( * )  WORK,
real, dimension( 2 )  RESULT 
)

SSPT21

Purpose:
 SSPT21  generally checks a decomposition of the form

         A = U S U'

 where ' means transpose, A is symmetric (stored in packed format), U
 is orthogonal, and S is diagonal (if KBAND=0) or symmetric
 tridiagonal (if KBAND=1).  If ITYPE=1, then U is represented as a
 dense matrix, otherwise the U is expressed as a product of
 Householder transformations, whose vectors are stored in the array
 "V" and whose scaling constants are in "TAU"; we shall use the
 letter "V" to refer to the product of Householder transformations
 (which should be equal to U).

 Specifically, if ITYPE=1, then:

         RESULT(1) = | A - U S U' | / ( |A| n ulp ) *andC>         RESULT(2) = | I - UU' | / ( n ulp )

 If ITYPE=2, then:

         RESULT(1) = | A - V S V' | / ( |A| n ulp )

 If ITYPE=3, then:

         RESULT(1) = | I - VU' | / ( n ulp )

 Packed storage means that, for example, if UPLO='U', then the columns
 of the upper triangle of A are stored one after another, so that
 A(1,j+1) immediately follows A(j,j) in the array AP.  Similarly, if
 UPLO='L', then the columns of the lower triangle of A are stored one
 after another in AP, so that A(j+1,j+1) immediately follows A(n,j)
 in the array AP.  This means that A(i,j) is stored in:

    AP( i + j*(j-1)/2 )                 if UPLO='U'

    AP( i + (2*n-j)*(j-1)/2 )           if UPLO='L'

 The array VP bears the same relation to the matrix V that A does to
 AP.

 For ITYPE > 1, the transformation U is expressed as a product
 of Householder transformations:

    If UPLO='U', then  V = H(n-1)...H(1),  where

        H(j) = I  -  tau(j) v(j) v(j)'

    and the first j-1 elements of v(j) are stored in V(1:j-1,j+1),
    (i.e., VP( j*(j+1)/2 + 1 : j*(j+1)/2 + j-1 ) ),
    the j-th element is 1, and the last n-j elements are 0.

    If UPLO='L', then  V = H(1)...H(n-1),  where

        H(j) = I  -  tau(j) v(j) v(j)'

    and the first j elements of v(j) are 0, the (j+1)-st is 1, and the
    (j+2)-nd through n-th elements are stored in V(j+2:n,j) (i.e.,
    in VP( (2*n-j)*(j-1)/2 + j+2 : (2*n-j)*(j-1)/2 + n ) .)
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the type of tests to be performed.
          1: U expressed as a dense orthogonal matrix:
             RESULT(1) = | A - U S U' | / ( |A| n ulp )   *andC>             RESULT(2) = | I - UU' | / ( n ulp )

          2: U expressed as a product V of Housholder transformations:
             RESULT(1) = | A - V S V' | / ( |A| n ulp )

          3: U expressed both as a dense orthogonal matrix and
             as a product of Housholder transformations:
             RESULT(1) = | I - VU' | / ( n ulp )
[in]UPLO
          UPLO is CHARACTER
          If UPLO='U', AP and VP are considered to contain the upper
          triangle of A and V.
          If UPLO='L', AP and VP are considered to contain the lower
          triangle of A and V.
[in]N
          N is INTEGER
          The size of the matrix.  If it is zero, SSPT21 does nothing.
          It must be at least zero.
[in]KBAND
          KBAND is INTEGER
          The bandwidth of the matrix.  It may only be zero or one.
          If zero, then S is diagonal, and E is not referenced.  If
          one, then S is symmetric tri-diagonal.
[in]AP
          AP is REAL array, dimension (N*(N+1)/2)
          The original (unfactored) matrix.  It is assumed to be
          symmetric, and contains the columns of just the upper
          triangle (UPLO='U') or only the lower triangle (UPLO='L'),
          packed one after another.
[in]D
          D is REAL array, dimension (N)
          The diagonal of the (symmetric tri-) diagonal matrix.
[in]E
          E is REAL array, dimension (N-1)
          The off-diagonal of the (symmetric tri-) diagonal matrix.
          E(1) is the (1,2) and (2,1) element, E(2) is the (2,3) and
          (3,2) element, etc.
          Not referenced if KBAND=0.
[in]U
          U is REAL array, dimension (LDU, N)
          If ITYPE=1 or 3, this contains the orthogonal matrix in
          the decomposition, expressed as a dense matrix.  If ITYPE=2,
          then it is not referenced.
[in]LDU
          LDU is INTEGER
          The leading dimension of U.  LDU must be at least N and
          at least 1.
[in]VP
          VP is REAL array, dimension (N*(N+1)/2)
          If ITYPE=2 or 3, the columns of this array contain the
          Householder vectors used to describe the orthogonal matrix
          in the decomposition, as described in purpose.
          *NOTE* If ITYPE=2 or 3, V is modified and restored.  The
          subdiagonal (if UPLO='L') or the superdiagonal (if UPLO='U')
          is set to one, and later reset to its original value, during
          the course of the calculation.
          If ITYPE=1, then it is neither referenced nor modified.
[in]TAU
          TAU is REAL array, dimension (N)
          If ITYPE >= 2, then TAU(j) is the scalar factor of
          v(j) v(j)' in the Householder transformation H(j) of
          the product  U = H(1)...H(n-2)
          If ITYPE < 2, then TAU is not referenced.
[out]WORK
          WORK is REAL array, dimension (N**2+N)
          Workspace.
[out]RESULT
          RESULT is REAL array, dimension (2)
          The values computed by the two tests described above.  The
          values are currently limited to 1/ulp, to avoid overflow.
          RESULT(1) is always modified.  RESULT(2) is modified only
          if ITYPE=1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 221 of file sspt21.f.

221 *
222 * -- LAPACK test routine (version 3.4.0) --
223 * -- LAPACK is a software package provided by Univ. of Tennessee, --
224 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
225 * November 2011
226 *
227 * .. Scalar Arguments ..
228  CHARACTER uplo
229  INTEGER itype, kband, ldu, n
230 * ..
231 * .. Array Arguments ..
232  REAL ap( * ), d( * ), e( * ), result( 2 ), tau( * ),
233  $ u( ldu, * ), vp( * ), work( * )
234 * ..
235 *
236 * =====================================================================
237 *
238 * .. Parameters ..
239  REAL zero, one, ten
240  parameter ( zero = 0.0e0, one = 1.0e0, ten = 10.0e0 )
241  REAL half
242  parameter ( half = 1.0e+0 / 2.0e+0 )
243 * ..
244 * .. Local Scalars ..
245  LOGICAL lower
246  CHARACTER cuplo
247  INTEGER iinfo, j, jp, jp1, jr, lap
248  REAL anorm, temp, ulp, unfl, vsave, wnorm
249 * ..
250 * .. External Functions ..
251  LOGICAL lsame
252  REAL sdot, slamch, slange, slansp
253  EXTERNAL lsame, sdot, slamch, slange, slansp
254 * ..
255 * .. External Subroutines ..
256  EXTERNAL saxpy, scopy, sgemm, slacpy, slaset, sopmtr,
257  $ sspmv, sspr, sspr2
258 * ..
259 * .. Intrinsic Functions ..
260  INTRINSIC max, min, real
261 * ..
262 * .. Executable Statements ..
263 *
264 * 1) Constants
265 *
266  result( 1 ) = zero
267  IF( itype.EQ.1 )
268  $ result( 2 ) = zero
269  IF( n.LE.0 )
270  $ RETURN
271 *
272  lap = ( n*( n+1 ) ) / 2
273 *
274  IF( lsame( uplo, 'U' ) ) THEN
275  lower = .false.
276  cuplo = 'U'
277  ELSE
278  lower = .true.
279  cuplo = 'L'
280  END IF
281 *
282  unfl = slamch( 'Safe minimum' )
283  ulp = slamch( 'Epsilon' )*slamch( 'Base' )
284 *
285 * Some Error Checks
286 *
287  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
288  result( 1 ) = ten / ulp
289  RETURN
290  END IF
291 *
292 * Do Test 1
293 *
294 * Norm of A:
295 *
296  IF( itype.EQ.3 ) THEN
297  anorm = one
298  ELSE
299  anorm = max( slansp( '1', cuplo, n, ap, work ), unfl )
300  END IF
301 *
302 * Compute error matrix:
303 *
304  IF( itype.EQ.1 ) THEN
305 *
306 * ITYPE=1: error = A - U S U'
307 *
308  CALL slaset( 'Full', n, n, zero, zero, work, n )
309  CALL scopy( lap, ap, 1, work, 1 )
310 *
311  DO 10 j = 1, n
312  CALL sspr( cuplo, n, -d( j ), u( 1, j ), 1, work )
313  10 CONTINUE
314 *
315  IF( n.GT.1 .AND. kband.EQ.1 ) THEN
316  DO 20 j = 1, n - 1
317  CALL sspr2( cuplo, n, -e( j ), u( 1, j ), 1, u( 1, j+1 ),
318  $ 1, work )
319  20 CONTINUE
320  END IF
321  wnorm = slansp( '1', cuplo, n, work, work( n**2+1 ) )
322 *
323  ELSE IF( itype.EQ.2 ) THEN
324 *
325 * ITYPE=2: error = V S V' - A
326 *
327  CALL slaset( 'Full', n, n, zero, zero, work, n )
328 *
329  IF( lower ) THEN
330  work( lap ) = d( n )
331  DO 40 j = n - 1, 1, -1
332  jp = ( ( 2*n-j )*( j-1 ) ) / 2
333  jp1 = jp + n - j
334  IF( kband.EQ.1 ) THEN
335  work( jp+j+1 ) = ( one-tau( j ) )*e( j )
336  DO 30 jr = j + 2, n
337  work( jp+jr ) = -tau( j )*e( j )*vp( jp+jr )
338  30 CONTINUE
339  END IF
340 *
341  IF( tau( j ).NE.zero ) THEN
342  vsave = vp( jp+j+1 )
343  vp( jp+j+1 ) = one
344  CALL sspmv( 'L', n-j, one, work( jp1+j+1 ),
345  $ vp( jp+j+1 ), 1, zero, work( lap+1 ), 1 )
346  temp = -half*tau( j )*sdot( n-j, work( lap+1 ), 1,
347  $ vp( jp+j+1 ), 1 )
348  CALL saxpy( n-j, temp, vp( jp+j+1 ), 1, work( lap+1 ),
349  $ 1 )
350  CALL sspr2( 'L', n-j, -tau( j ), vp( jp+j+1 ), 1,
351  $ work( lap+1 ), 1, work( jp1+j+1 ) )
352  vp( jp+j+1 ) = vsave
353  END IF
354  work( jp+j ) = d( j )
355  40 CONTINUE
356  ELSE
357  work( 1 ) = d( 1 )
358  DO 60 j = 1, n - 1
359  jp = ( j*( j-1 ) ) / 2
360  jp1 = jp + j
361  IF( kband.EQ.1 ) THEN
362  work( jp1+j ) = ( one-tau( j ) )*e( j )
363  DO 50 jr = 1, j - 1
364  work( jp1+jr ) = -tau( j )*e( j )*vp( jp1+jr )
365  50 CONTINUE
366  END IF
367 *
368  IF( tau( j ).NE.zero ) THEN
369  vsave = vp( jp1+j )
370  vp( jp1+j ) = one
371  CALL sspmv( 'U', j, one, work, vp( jp1+1 ), 1, zero,
372  $ work( lap+1 ), 1 )
373  temp = -half*tau( j )*sdot( j, work( lap+1 ), 1,
374  $ vp( jp1+1 ), 1 )
375  CALL saxpy( j, temp, vp( jp1+1 ), 1, work( lap+1 ),
376  $ 1 )
377  CALL sspr2( 'U', j, -tau( j ), vp( jp1+1 ), 1,
378  $ work( lap+1 ), 1, work )
379  vp( jp1+j ) = vsave
380  END IF
381  work( jp1+j+1 ) = d( j+1 )
382  60 CONTINUE
383  END IF
384 *
385  DO 70 j = 1, lap
386  work( j ) = work( j ) - ap( j )
387  70 CONTINUE
388  wnorm = slansp( '1', cuplo, n, work, work( lap+1 ) )
389 *
390  ELSE IF( itype.EQ.3 ) THEN
391 *
392 * ITYPE=3: error = U V' - I
393 *
394  IF( n.LT.2 )
395  $ RETURN
396  CALL slacpy( ' ', n, n, u, ldu, work, n )
397  CALL sopmtr( 'R', cuplo, 'T', n, n, vp, tau, work, n,
398  $ work( n**2+1 ), iinfo )
399  IF( iinfo.NE.0 ) THEN
400  result( 1 ) = ten / ulp
401  RETURN
402  END IF
403 *
404  DO 80 j = 1, n
405  work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - one
406  80 CONTINUE
407 *
408  wnorm = slange( '1', n, n, work, n, work( n**2+1 ) )
409  END IF
410 *
411  IF( anorm.GT.wnorm ) THEN
412  result( 1 ) = ( wnorm / anorm ) / ( n*ulp )
413  ELSE
414  IF( anorm.LT.one ) THEN
415  result( 1 ) = ( min( wnorm, n*anorm ) / anorm ) / ( n*ulp )
416  ELSE
417  result( 1 ) = min( wnorm / anorm, REAL( N ) ) / ( n*ulp )
418  END IF
419  END IF
420 *
421 * Do Test 2
422 *
423 * Compute UU' - I
424 *
425  IF( itype.EQ.1 ) THEN
426  CALL sgemm( 'N', 'C', n, n, n, one, u, ldu, u, ldu, zero, work,
427  $ n )
428 *
429  DO 90 j = 1, n
430  work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - one
431  90 CONTINUE
432 *
433  result( 2 ) = min( slange( '1', n, n, work, n,
434  $ work( n**2+1 ) ), REAL( N ) ) / ( n*ulp )
435  END IF
436 *
437  RETURN
438 *
439 * End of SSPT21
440 *
subroutine sspr2(UPLO, N, ALPHA, X, INCX, Y, INCY, AP)
SSPR2
Definition: sspr2.f:144
subroutine sspmv(UPLO, N, ALPHA, AP, X, INCX, BETA, Y, INCY)
SSPMV
Definition: sspmv.f:149
real function sdot(N, SX, INCX, SY, INCY)
SDOT
Definition: sdot.f:53
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:189
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
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:112
real function slange(NORM, M, N, A, LDA, WORK)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: slange.f:116
subroutine sopmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
SOPMTR
Definition: sopmtr.f:152
subroutine saxpy(N, SA, SX, INCX, SY, INCY)
SAXPY
Definition: saxpy.f:54
real function slansp(NORM, UPLO, N, AP, WORK)
SLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: slansp.f:116
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine sspr(UPLO, N, ALPHA, X, INCX, AP)
SSPR
Definition: sspr.f:129
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: