LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sqrt05 ( integer  M,
integer  N,
integer  L,
integer  NB,
real, dimension(6)  RESULT 
)

SQRT05

Purpose:
 SQRT05 tests STPQRT and STPMQRT.
Parameters
[in]M
          M is INTEGER
          Number of rows in lower part of the test matrix.
[in]N
          N is INTEGER
          Number of columns in test matrix.
[in]L
          L is INTEGER
          The number of rows of the upper trapezoidal part the
          lower test matrix.  0 <= L <= M.
[in]NB
          NB is INTEGER
          Block size of test matrix.  NB <= N.
[out]RESULT
          RESULT is REAL array, dimension (6)
          Results of each of the six tests below.

          RESULT(1) = | A - Q R |
          RESULT(2) = | I - Q^H Q |
          RESULT(3) = | Q C - Q C |
          RESULT(4) = | Q^H C - Q^H C |
          RESULT(5) = | C Q - C Q | 
          RESULT(6) = | C Q^H - C Q^H |
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
April 2012

Definition at line 82 of file sqrt05.f.

82  IMPLICIT NONE
83 *
84 * -- LAPACK test routine (version 3.4.1) --
85 * -- LAPACK is a software package provided by Univ. of Tennessee, --
86 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
87 * April 2012
88 *
89 * .. Scalar Arguments ..
90  INTEGER lwork, m, n, l, nb, ldt
91 * .. Return values ..
92  REAL result(6)
93 *
94 * =====================================================================
95 *
96 * ..
97 * .. Local allocatable arrays
98  REAL, ALLOCATABLE :: af(:,:), q(:,:),
99  $ r(:,:), rwork(:), work( : ), t(:,:),
100  $ cf(:,:), df(:,:), a(:,:), c(:,:), d(:,:)
101 *
102 * .. Parameters ..
103  REAL zero, one
104  parameter( zero = 0.0, one = 1.0 )
105 * ..
106 * .. Local Scalars ..
107  INTEGER info, j, k, m2, np1
108  REAL anorm, eps, resid, cnorm, dnorm
109 * ..
110 * .. Local Arrays ..
111  INTEGER iseed( 4 )
112 * ..
113 * .. External Functions ..
114  REAL slamch
115  REAL slange, slansy
116  LOGICAL lsame
117  EXTERNAL slamch, slange, slansy, lsame
118 * ..
119 * .. Data statements ..
120  DATA iseed / 1988, 1989, 1990, 1991 /
121 *
122  eps = slamch( 'Epsilon' )
123  k = n
124  m2 = m+n
125  IF( m.GT.0 ) THEN
126  np1 = n+1
127  ELSE
128  np1 = 1
129  END IF
130  lwork = m2*m2*nb
131 *
132 * Dynamically allocate all arrays
133 *
134  ALLOCATE(a(m2,n),af(m2,n),q(m2,m2),r(m2,m2),rwork(m2),
135  $ work(lwork),t(nb,n),c(m2,n),cf(m2,n),
136  $ d(n,m2),df(n,m2) )
137 *
138 * Put random stuff into A
139 *
140  ldt=nb
141  CALL slaset( 'Full', m2, n, zero, zero, a, m2 )
142  CALL slaset( 'Full', nb, n, zero, zero, t, nb )
143  DO j=1,n
144  CALL slarnv( 2, iseed, j, a( 1, j ) )
145  END DO
146  IF( m.GT.0 ) THEN
147  DO j=1,n
148  CALL slarnv( 2, iseed, m-l, a( n+1, j ) )
149  END DO
150  END IF
151  IF( l.GT.0 ) THEN
152  DO j=1,n
153  CALL slarnv( 2, iseed, min(j,l), a( n+m-l+1, j ) )
154  END DO
155  END IF
156 *
157 * Copy the matrix A to the array AF.
158 *
159  CALL slacpy( 'Full', m2, n, a, m2, af, m2 )
160 *
161 * Factor the matrix A in the array AF.
162 *
163  CALL stpqrt( m,n,l,nb,af,m2,af(np1,1),m2,t,ldt,work,info)
164 *
165 * Generate the (M+N)-by-(M+N) matrix Q by applying H to I
166 *
167  CALL slaset( 'Full', m2, m2, zero, one, q, m2 )
168  CALL sgemqrt( 'R', 'N', m2, m2, k, nb, af, m2, t, ldt, q, m2,
169  $ work, info )
170 *
171 * Copy R
172 *
173  CALL slaset( 'Full', m2, n, zero, zero, r, m2 )
174  CALL slacpy( 'Upper', m2, n, af, m2, r, m2 )
175 *
176 * Compute |R - Q'*A| / |A| and store in RESULT(1)
177 *
178  CALL sgemm( 'T', 'N', m2, n, m2, -one, q, m2, a, m2, one, r, m2 )
179  anorm = slange( '1', m2, n, a, m2, rwork )
180  resid = slange( '1', m2, n, r, m2, rwork )
181  IF( anorm.GT.zero ) THEN
182  result( 1 ) = resid / (eps*anorm*max(1,m2))
183  ELSE
184  result( 1 ) = zero
185  END IF
186 *
187 * Compute |I - Q'*Q| and store in RESULT(2)
188 *
189  CALL slaset( 'Full', m2, m2, zero, one, r, m2 )
190  CALL ssyrk( 'U', 'C', m2, m2, -one, q, m2, one,
191  $ r, m2 )
192  resid = slansy( '1', 'Upper', m2, r, m2, rwork )
193  result( 2 ) = resid / (eps*max(1,m2))
194 *
195 * Generate random m-by-n matrix C and a copy CF
196 *
197  DO j=1,n
198  CALL slarnv( 2, iseed, m2, c( 1, j ) )
199  END DO
200  cnorm = slange( '1', m2, n, c, m2, rwork)
201  CALL slacpy( 'Full', m2, n, c, m2, cf, m2 )
202 *
203 * Apply Q to C as Q*C
204 *
205  CALL stpmqrt( 'L','N', m,n,k,l,nb,af(np1,1),m2,t,ldt,cf,
206  $ m2,cf(np1,1),m2,work,info)
207 *
208 * Compute |Q*C - Q*C| / |C|
209 *
210  CALL sgemm( 'N', 'N', m2, n, m2, -one, q,m2,c,m2,one,cf,m2)
211  resid = slange( '1', m2, n, cf, m2, rwork )
212  IF( cnorm.GT.zero ) THEN
213  result( 3 ) = resid / (eps*max(1,m2)*cnorm)
214  ELSE
215  result( 3 ) = zero
216  END IF
217 *
218 * Copy C into CF again
219 *
220  CALL slacpy( 'Full', m2, n, c, m2, cf, m2 )
221 *
222 * Apply Q to C as QT*C
223 *
224  CALL stpmqrt('L','T',m,n,k,l,nb,af(np1,1),m2,t,ldt,cf,m2,
225  $ cf(np1,1),m2,work,info)
226 *
227 * Compute |QT*C - QT*C| / |C|
228 *
229  CALL sgemm('T','N',m2,n,m2,-one,q,m2,c,m2,one,cf,m2)
230  resid = slange( '1', m2, n, cf, m2, rwork )
231  IF( cnorm.GT.zero ) THEN
232  result( 4 ) = resid / (eps*max(1,m2)*cnorm)
233  ELSE
234  result( 4 ) = zero
235  END IF
236 *
237 * Generate random n-by-m matrix D and a copy DF
238 *
239  DO j=1,m2
240  CALL slarnv( 2, iseed, n, d( 1, j ) )
241  END DO
242  dnorm = slange( '1', n, m2, d, n, rwork)
243  CALL slacpy( 'Full', n, m2, d, n, df, n )
244 *
245 * Apply Q to D as D*Q
246 *
247  CALL stpmqrt('R','N',n,m,n,l,nb,af(np1,1),m2,t,ldt,df,n,
248  $ df(1,np1),n,work,info)
249 *
250 * Compute |D*Q - D*Q| / |D|
251 *
252  CALL sgemm('N','N',n,m2,m2,-one,d,n,q,m2,one,df,n)
253  resid = slange('1',n, m2,df,n,rwork )
254  IF( cnorm.GT.zero ) THEN
255  result( 5 ) = resid / (eps*max(1,m2)*dnorm)
256  ELSE
257  result( 5 ) = zero
258  END IF
259 *
260 * Copy D into DF again
261 *
262  CALL slacpy('Full',n,m2,d,n,df,n )
263 *
264 * Apply Q to D as D*QT
265 *
266  CALL stpmqrt('R','T',n,m,n,l,nb,af(np1,1),m2,t,ldt,df,n,
267  $ df(1,np1),n,work,info)
268 
269 *
270 * Compute |D*QT - D*QT| / |D|
271 *
272  CALL sgemm( 'N', 'T', n, m2, m2, -one, d, n, q, m2, one, df, n )
273  resid = slange( '1', n, m2, df, n, rwork )
274  IF( cnorm.GT.zero ) THEN
275  result( 6 ) = resid / (eps*max(1,m2)*dnorm)
276  ELSE
277  result( 6 ) = zero
278  END IF
279 *
280 * Deallocate all arrays
281 *
282  DEALLOCATE ( a, af, q, r, rwork, work, t, c, d, cf, df)
283  RETURN
subroutine ssyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
SSYRK
Definition: ssyrk.f:171
subroutine slarnv(IDIST, ISEED, N, X)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: slarnv.f:99
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:189
subroutine stpmqrt(SIDE, TRANS, M, N, K, L, NB, V, LDV, T, LDT, A, LDA, B, LDB, WORK, INFO)
STPMQRT
Definition: stpmqrt.f:218
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 sgemqrt(SIDE, TRANS, M, N, K, NB, V, LDV, T, LDT, C, LDC, WORK, INFO)
SGEMQRT
Definition: sgemqrt.f:170
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 stpqrt(M, N, L, NB, A, LDA, B, LDB, T, LDT, WORK, INFO)
STPQRT
Definition: stpqrt.f:191
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
real function slansy(NORM, UPLO, N, A, LDA, WORK)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: slansy.f:124

Here is the call graph for this function:

Here is the caller graph for this function: