LAPACK  3.9.0
LAPACK: Linear Algebra PACKage
sqrt05.f
Go to the documentation of this file.
1 *> \brief \b SQRT05
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE SQRT05(M,N,L,NB,RESULT)
12 *
13 * .. Scalar Arguments ..
14 * INTEGER LWORK, M, N, L, NB, LDT
15 * .. Return values ..
16 * REAL RESULT(6)
17 *
18 *
19 *> \par Purpose:
20 * =============
21 *>
22 *> \verbatim
23 *>
24 *> SQRT05 tests STPQRT and STPMQRT.
25 *> \endverbatim
26 *
27 * Arguments:
28 * ==========
29 *
30 *> \param[in] M
31 *> \verbatim
32 *> M is INTEGER
33 *> Number of rows in lower part of the test matrix.
34 *> \endverbatim
35 *>
36 *> \param[in] N
37 *> \verbatim
38 *> N is INTEGER
39 *> Number of columns in test matrix.
40 *> \endverbatim
41 *>
42 *> \param[in] L
43 *> \verbatim
44 *> L is INTEGER
45 *> The number of rows of the upper trapezoidal part the
46 *> lower test matrix. 0 <= L <= M.
47 *> \endverbatim
48 *>
49 *> \param[in] NB
50 *> \verbatim
51 *> NB is INTEGER
52 *> Block size of test matrix. NB <= N.
53 *> \endverbatim
54 *>
55 *> \param[out] RESULT
56 *> \verbatim
57 *> RESULT is REAL array, dimension (6)
58 *> Results of each of the six tests below.
59 *>
60 *> RESULT(1) = | A - Q R |
61 *> RESULT(2) = | I - Q^H Q |
62 *> RESULT(3) = | Q C - Q C |
63 *> RESULT(4) = | Q^H C - Q^H C |
64 *> RESULT(5) = | C Q - C Q |
65 *> RESULT(6) = | C Q^H - C Q^H |
66 *> \endverbatim
67 *
68 * Authors:
69 * ========
70 *
71 *> \author Univ. of Tennessee
72 *> \author Univ. of California Berkeley
73 *> \author Univ. of Colorado Denver
74 *> \author NAG Ltd.
75 *
76 *> \date April 2012
77 *
78 *> \ingroup single_lin
79 *
80 * =====================================================================
81  SUBROUTINE sqrt05(M,N,L,NB,RESULT)
82  IMPLICIT NONE
83 *
84 * -- LAPACK test routine (version 3.8.0) --
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 Subroutine ..
114  EXTERNAL sgemm, slarnv, stpmqrt, stpqrt, sgemqrt, ssyrk, slacpy,
115  $ slaset
116 * ..
117 * .. External Functions ..
118  REAL SLAMCH
119  REAL SLANGE, SLANSY
120  LOGICAL LSAME
121  EXTERNAL slamch, slange, slansy, lsame
122 * ..
123 * .. Data statements ..
124  DATA iseed / 1988, 1989, 1990, 1991 /
125 *
126  eps = slamch( 'Epsilon' )
127  k = n
128  m2 = m+n
129  IF( m.GT.0 ) THEN
130  np1 = n+1
131  ELSE
132  np1 = 1
133  END IF
134  lwork = m2*m2*nb
135 *
136 * Dynamically allocate all arrays
137 *
138  ALLOCATE(a(m2,n),af(m2,n),q(m2,m2),r(m2,m2),rwork(m2),
139  $ work(lwork),t(nb,n),c(m2,n),cf(m2,n),
140  $ d(n,m2),df(n,m2) )
141 *
142 * Put random stuff into A
143 *
144  ldt=nb
145  CALL slaset( 'Full', m2, n, zero, zero, a, m2 )
146  CALL slaset( 'Full', nb, n, zero, zero, t, nb )
147  DO j=1,n
148  CALL slarnv( 2, iseed, j, a( 1, j ) )
149  END DO
150  IF( m.GT.0 ) THEN
151  DO j=1,n
152  CALL slarnv( 2, iseed, m-l, a( n+1, j ) )
153  END DO
154  END IF
155  IF( l.GT.0 ) THEN
156  DO j=1,n
157  CALL slarnv( 2, iseed, min(j,l), a( n+m-l+1, j ) )
158  END DO
159  END IF
160 *
161 * Copy the matrix A to the array AF.
162 *
163  CALL slacpy( 'Full', m2, n, a, m2, af, m2 )
164 *
165 * Factor the matrix A in the array AF.
166 *
167  CALL stpqrt( m,n,l,nb,af,m2,af(np1,1),m2,t,ldt,work,info)
168 *
169 * Generate the (M+N)-by-(M+N) matrix Q by applying H to I
170 *
171  CALL slaset( 'Full', m2, m2, zero, one, q, m2 )
172  CALL sgemqrt( 'R', 'N', m2, m2, k, nb, af, m2, t, ldt, q, m2,
173  $ work, info )
174 *
175 * Copy R
176 *
177  CALL slaset( 'Full', m2, n, zero, zero, r, m2 )
178  CALL slacpy( 'Upper', m2, n, af, m2, r, m2 )
179 *
180 * Compute |R - Q'*A| / |A| and store in RESULT(1)
181 *
182  CALL sgemm( 'T', 'N', m2, n, m2, -one, q, m2, a, m2, one, r, m2 )
183  anorm = slange( '1', m2, n, a, m2, rwork )
184  resid = slange( '1', m2, n, r, m2, rwork )
185  IF( anorm.GT.zero ) THEN
186  result( 1 ) = resid / (eps*anorm*max(1,m2))
187  ELSE
188  result( 1 ) = zero
189  END IF
190 *
191 * Compute |I - Q'*Q| and store in RESULT(2)
192 *
193  CALL slaset( 'Full', m2, m2, zero, one, r, m2 )
194  CALL ssyrk( 'U', 'C', m2, m2, -one, q, m2, one,
195  $ r, m2 )
196  resid = slansy( '1', 'Upper', m2, r, m2, rwork )
197  result( 2 ) = resid / (eps*max(1,m2))
198 *
199 * Generate random m-by-n matrix C and a copy CF
200 *
201  DO j=1,n
202  CALL slarnv( 2, iseed, m2, c( 1, j ) )
203  END DO
204  cnorm = slange( '1', m2, n, c, m2, rwork)
205  CALL slacpy( 'Full', m2, n, c, m2, cf, m2 )
206 *
207 * Apply Q to C as Q*C
208 *
209  CALL stpmqrt( 'L','N', m,n,k,l,nb,af(np1,1),m2,t,ldt,cf,
210  $ m2,cf(np1,1),m2,work,info)
211 *
212 * Compute |Q*C - Q*C| / |C|
213 *
214  CALL sgemm( 'N', 'N', m2, n, m2, -one, q,m2,c,m2,one,cf,m2)
215  resid = slange( '1', m2, n, cf, m2, rwork )
216  IF( cnorm.GT.zero ) THEN
217  result( 3 ) = resid / (eps*max(1,m2)*cnorm)
218  ELSE
219  result( 3 ) = zero
220  END IF
221 *
222 * Copy C into CF again
223 *
224  CALL slacpy( 'Full', m2, n, c, m2, cf, m2 )
225 *
226 * Apply Q to C as QT*C
227 *
228  CALL stpmqrt('L','T',m,n,k,l,nb,af(np1,1),m2,t,ldt,cf,m2,
229  $ cf(np1,1),m2,work,info)
230 *
231 * Compute |QT*C - QT*C| / |C|
232 *
233  CALL sgemm('T','N',m2,n,m2,-one,q,m2,c,m2,one,cf,m2)
234  resid = slange( '1', m2, n, cf, m2, rwork )
235  IF( cnorm.GT.zero ) THEN
236  result( 4 ) = resid / (eps*max(1,m2)*cnorm)
237  ELSE
238  result( 4 ) = zero
239  END IF
240 *
241 * Generate random n-by-m matrix D and a copy DF
242 *
243  DO j=1,m2
244  CALL slarnv( 2, iseed, n, d( 1, j ) )
245  END DO
246  dnorm = slange( '1', n, m2, d, n, rwork)
247  CALL slacpy( 'Full', n, m2, d, n, df, n )
248 *
249 * Apply Q to D as D*Q
250 *
251  CALL stpmqrt('R','N',n,m,n,l,nb,af(np1,1),m2,t,ldt,df,n,
252  $ df(1,np1),n,work,info)
253 *
254 * Compute |D*Q - D*Q| / |D|
255 *
256  CALL sgemm('N','N',n,m2,m2,-one,d,n,q,m2,one,df,n)
257  resid = slange('1',n, m2,df,n,rwork )
258  IF( cnorm.GT.zero ) THEN
259  result( 5 ) = resid / (eps*max(1,m2)*dnorm)
260  ELSE
261  result( 5 ) = zero
262  END IF
263 *
264 * Copy D into DF again
265 *
266  CALL slacpy('Full',n,m2,d,n,df,n )
267 *
268 * Apply Q to D as D*QT
269 *
270  CALL stpmqrt('R','T',n,m,n,l,nb,af(np1,1),m2,t,ldt,df,n,
271  $ df(1,np1),n,work,info)
272 
273 *
274 * Compute |D*QT - D*QT| / |D|
275 *
276  CALL sgemm( 'N', 'T', n, m2, m2, -one, d, n, q, m2, one, df, n )
277  resid = slange( '1', n, m2, df, n, rwork )
278  IF( cnorm.GT.zero ) THEN
279  result( 6 ) = resid / (eps*max(1,m2)*dnorm)
280  ELSE
281  result( 6 ) = zero
282  END IF
283 *
284 * Deallocate all arrays
285 *
286  DEALLOCATE ( a, af, q, r, rwork, work, t, c, d, cf, df)
287  RETURN
288  END
289 
slarnv
subroutine slarnv(IDIST, ISEED, N, X)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: slarnv.f:99
sgemm
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:189
slacpy
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
stpmqrt
subroutine stpmqrt(SIDE, TRANS, M, N, K, L, NB, V, LDV, T, LDT, A, LDA, B, LDB, WORK, INFO)
STPMQRT
Definition: stpmqrt.f:218
slaset
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
sgemqrt
subroutine sgemqrt(SIDE, TRANS, M, N, K, NB, V, LDV, T, LDT, C, LDC, WORK, INFO)
SGEMQRT
Definition: sgemqrt.f:170
stpqrt
subroutine stpqrt(M, N, L, NB, A, LDA, B, LDB, T, LDT, WORK, INFO)
STPQRT
Definition: stpqrt.f:191
ssyrk
subroutine ssyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
SSYRK
Definition: ssyrk.f:171
sqrt05
subroutine sqrt05(M, N, L, NB, RESULT)
SQRT05
Definition: sqrt05.f:82