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

◆ sgqrts()

subroutine sgqrts ( integer  n,
integer  m,
integer  p,
real, dimension( lda, * )  a,
real, dimension( lda, * )  af,
real, dimension( lda, * )  q,
real, dimension( lda, * )  r,
integer  lda,
real, dimension( * )  taua,
real, dimension( ldb, * )  b,
real, dimension( ldb, * )  bf,
real, dimension( ldb, * )  z,
real, dimension( ldb, * )  t,
real, dimension( ldb, * )  bwk,
integer  ldb,
real, dimension( * )  taub,
real, dimension( lwork )  work,
integer  lwork,
real, dimension( * )  rwork,
real, dimension( 4 )  result 
)

SGQRTS

Purpose:
 SGQRTS tests SGGQRF, which computes the GQR factorization of an
 N-by-M matrix A and a N-by-P matrix B: A = Q*R and B = Q*T*Z.
Parameters
[in]N
          N is INTEGER
          The number of rows of the matrices A and B.  N >= 0.
[in]M
          M is INTEGER
          The number of columns of the matrix A.  M >= 0.
[in]P
          P is INTEGER
          The number of columns of the matrix B.  P >= 0.
[in]A
          A is REAL array, dimension (LDA,M)
          The N-by-M matrix A.
[out]AF
          AF is REAL array, dimension (LDA,N)
          Details of the GQR factorization of A and B, as returned
          by SGGQRF, see SGGQRF for further details.
[out]Q
          Q is REAL array, dimension (LDA,N)
          The M-by-M orthogonal matrix Q.
[out]R
          R is REAL array, dimension (LDA,MAX(M,N))
[in]LDA
          LDA is INTEGER
          The leading dimension of the arrays A, AF, R and Q.
          LDA >= max(M,N).
[out]TAUA
          TAUA is REAL array, dimension (min(M,N))
          The scalar factors of the elementary reflectors, as returned
          by SGGQRF.
[in]B
          B is REAL array, dimension (LDB,P)
          On entry, the N-by-P matrix A.
[out]BF
          BF is REAL array, dimension (LDB,N)
          Details of the GQR factorization of A and B, as returned
          by SGGQRF, see SGGQRF for further details.
[out]Z
          Z is REAL array, dimension (LDB,P)
          The P-by-P orthogonal matrix Z.
[out]T
          T is REAL array, dimension (LDB,max(P,N))
[out]BWK
          BWK is REAL array, dimension (LDB,N)
[in]LDB
          LDB is INTEGER
          The leading dimension of the arrays B, BF, Z and T.
          LDB >= max(P,N).
[out]TAUB
          TAUB is REAL array, dimension (min(P,N))
          The scalar factors of the elementary reflectors, as returned
          by SGGRQF.
[out]WORK
          WORK is REAL array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK, LWORK >= max(N,M,P)**2.
[out]RWORK
          RWORK is REAL array, dimension (max(N,M,P))
[out]RESULT
          RESULT is REAL array, dimension (4)
          The test ratios:
            RESULT(1) = norm( R - Q'*A ) / ( MAX(M,N)*norm(A)*ULP)
            RESULT(2) = norm( T*Z - Q'*B ) / (MAX(P,N)*norm(B)*ULP)
            RESULT(3) = norm( I - Q'*Q ) / ( M*ULP )
            RESULT(4) = norm( I - Z'*Z ) / ( P*ULP )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 174 of file sgqrts.f.

176*
177* -- LAPACK test routine --
178* -- LAPACK is a software package provided by Univ. of Tennessee, --
179* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
180*
181* .. Scalar Arguments ..
182 INTEGER LDA, LDB, LWORK, M, P, N
183* ..
184* .. Array Arguments ..
185 REAL A( LDA, * ), AF( LDA, * ), R( LDA, * ),
186 $ Q( LDA, * ), B( LDB, * ), BF( LDB, * ),
187 $ T( LDB, * ), Z( LDB, * ), BWK( LDB, * ),
188 $ TAUA( * ), TAUB( * ), RESULT( 4 ),
189 $ RWORK( * ), WORK( LWORK )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 REAL ZERO, ONE
196 parameter( zero = 0.0e+0, one = 1.0e+0 )
197 REAL ROGUE
198 parameter( rogue = -1.0e+10 )
199* ..
200* .. Local Scalars ..
201 INTEGER INFO
202 REAL ANORM, BNORM, ULP, UNFL, RESID
203* ..
204* .. External Functions ..
205 REAL SLAMCH, SLANGE, SLANSY
206 EXTERNAL slamch, slange, slansy
207* ..
208* .. External Subroutines ..
209 EXTERNAL sgemm, slacpy, slaset, sorgqr,
210 $ sorgrq, ssyrk
211* ..
212* .. Intrinsic Functions ..
213 INTRINSIC max, min, real
214* ..
215* .. Executable Statements ..
216*
217 ulp = slamch( 'Precision' )
218 unfl = slamch( 'Safe minimum' )
219*
220* Copy the matrix A to the array AF.
221*
222 CALL slacpy( 'Full', n, m, a, lda, af, lda )
223 CALL slacpy( 'Full', n, p, b, ldb, bf, ldb )
224*
225 anorm = max( slange( '1', n, m, a, lda, rwork ), unfl )
226 bnorm = max( slange( '1', n, p, b, ldb, rwork ), unfl )
227*
228* Factorize the matrices A and B in the arrays AF and BF.
229*
230 CALL sggqrf( n, m, p, af, lda, taua, bf, ldb, taub, work,
231 $ lwork, info )
232*
233* Generate the N-by-N matrix Q
234*
235 CALL slaset( 'Full', n, n, rogue, rogue, q, lda )
236 CALL slacpy( 'Lower', n-1, m, af( 2,1 ), lda, q( 2,1 ), lda )
237 CALL sorgqr( n, n, min( n, m ), q, lda, taua, work, lwork, info )
238*
239* Generate the P-by-P matrix Z
240*
241 CALL slaset( 'Full', p, p, rogue, rogue, z, ldb )
242 IF( n.LE.p ) THEN
243 IF( n.GT.0 .AND. n.LT.p )
244 $ CALL slacpy( 'Full', n, p-n, bf, ldb, z( p-n+1, 1 ), ldb )
245 IF( n.GT.1 )
246 $ CALL slacpy( 'Lower', n-1, n-1, bf( 2, p-n+1 ), ldb,
247 $ z( p-n+2, p-n+1 ), ldb )
248 ELSE
249 IF( p.GT.1)
250 $ CALL slacpy( 'Lower', p-1, p-1, bf( n-p+2, 1 ), ldb,
251 $ z( 2, 1 ), ldb )
252 END IF
253 CALL sorgrq( p, p, min( n, p ), z, ldb, taub, work, lwork, info )
254*
255* Copy R
256*
257 CALL slaset( 'Full', n, m, zero, zero, r, lda )
258 CALL slacpy( 'Upper', n, m, af, lda, r, lda )
259*
260* Copy T
261*
262 CALL slaset( 'Full', n, p, zero, zero, t, ldb )
263 IF( n.LE.p ) THEN
264 CALL slacpy( 'Upper', n, n, bf( 1, p-n+1 ), ldb, t( 1, p-n+1 ),
265 $ ldb )
266 ELSE
267 CALL slacpy( 'Full', n-p, p, bf, ldb, t, ldb )
268 CALL slacpy( 'Upper', p, p, bf( n-p+1, 1 ), ldb, t( n-p+1, 1 ),
269 $ ldb )
270 END IF
271*
272* Compute R - Q'*A
273*
274 CALL sgemm( 'Transpose', 'No transpose', n, m, n, -one, q, lda, a,
275 $ lda, one, r, lda )
276*
277* Compute norm( R - Q'*A ) / ( MAX(M,N)*norm(A)*ULP ) .
278*
279 resid = slange( '1', n, m, r, lda, rwork )
280 IF( anorm.GT.zero ) THEN
281 result( 1 ) = ( ( resid / real( max(1,m,n) ) ) / anorm ) / ulp
282 ELSE
283 result( 1 ) = zero
284 END IF
285*
286* Compute T*Z - Q'*B
287*
288 CALL sgemm( 'No Transpose', 'No transpose', n, p, p, one, t, ldb,
289 $ z, ldb, zero, bwk, ldb )
290 CALL sgemm( 'Transpose', 'No transpose', n, p, n, -one, q, lda,
291 $ b, ldb, one, bwk, ldb )
292*
293* Compute norm( T*Z - Q'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
294*
295 resid = slange( '1', n, p, bwk, ldb, rwork )
296 IF( bnorm.GT.zero ) THEN
297 result( 2 ) = ( ( resid / real( max(1,p,n ) ) )/bnorm ) / ulp
298 ELSE
299 result( 2 ) = zero
300 END IF
301*
302* Compute I - Q'*Q
303*
304 CALL slaset( 'Full', n, n, zero, one, r, lda )
305 CALL ssyrk( 'Upper', 'Transpose', n, n, -one, q, lda, one, r,
306 $ lda )
307*
308* Compute norm( I - Q'*Q ) / ( N * ULP ) .
309*
310 resid = slansy( '1', 'Upper', n, r, lda, rwork )
311 result( 3 ) = ( resid / real( max( 1, n ) ) ) / ulp
312*
313* Compute I - Z'*Z
314*
315 CALL slaset( 'Full', p, p, zero, one, t, ldb )
316 CALL ssyrk( 'Upper', 'Transpose', p, p, -one, z, ldb, one, t,
317 $ ldb )
318*
319* Compute norm( I - Z'*Z ) / ( P*ULP ) .
320*
321 resid = slansy( '1', 'Upper', p, t, ldb, rwork )
322 result( 4 ) = ( resid / real( max( 1, p ) ) ) / ulp
323*
324 RETURN
325*
326* End of SGQRTS
327*
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine sggqrf(n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
SGGQRF
Definition sggqrf.f:215
subroutine ssyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
SSYRK
Definition ssyrk.f:169
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
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:114
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,...
Definition slansy.f:122
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
subroutine sorgqr(m, n, k, a, lda, tau, work, lwork, info)
SORGQR
Definition sorgqr.f:128
subroutine sorgrq(m, n, k, a, lda, tau, work, lwork, info)
SORGRQ
Definition sorgrq.f:128
Here is the call graph for this function:
Here is the caller graph for this function: