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

◆ zgqrts()

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

ZGQRTS

Purpose:
 ZGQRTS tests ZGGQRF, 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 COMPLEX*16 array, dimension (LDA,M)
          The N-by-M matrix A.
[out]AF
          AF is COMPLEX*16 array, dimension (LDA,N)
          Details of the GQR factorization of A and B, as returned
          by ZGGQRF, see CGGQRF for further details.
[out]Q
          Q is COMPLEX*16 array, dimension (LDA,N)
          The M-by-M unitary matrix Q.
[out]R
          R is COMPLEX*16 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 COMPLEX*16 array, dimension (min(M,N))
          The scalar factors of the elementary reflectors, as returned
          by ZGGQRF.
[in]B
          B is COMPLEX*16 array, dimension (LDB,P)
          On entry, the N-by-P matrix A.
[out]BF
          BF is COMPLEX*16 array, dimension (LDB,N)
          Details of the GQR factorization of A and B, as returned
          by ZGGQRF, see CGGQRF for further details.
[out]Z
          Z is COMPLEX*16 array, dimension (LDB,P)
          The P-by-P unitary matrix Z.
[out]T
          T is COMPLEX*16 array, dimension (LDB,max(P,N))
[out]BWK
          BWK is COMPLEX*16 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 COMPLEX*16 array, dimension (min(P,N))
          The scalar factors of the elementary reflectors, as returned
          by DGGRQF.
[out]WORK
          WORK is COMPLEX*16 array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK, LWORK >= max(N,M,P)**2.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(N,M,P))
[out]RESULT
          RESULT is DOUBLE PRECISION 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 zgqrts.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, N, P
183* ..
184* .. Array Arguments ..
185 DOUBLE PRECISION RESULT( 4 ), RWORK( * )
186 COMPLEX*16 A( LDA, * ), AF( LDA, * ), B( LDB, * ),
187 $ BF( LDB, * ), BWK( LDB, * ), Q( LDA, * ),
188 $ R( LDA, * ), T( LDB, * ), TAUA( * ), TAUB( * ),
189 $ WORK( LWORK ), Z( LDB, * )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 DOUBLE PRECISION ZERO, ONE
196 parameter( zero = 0.0d+0, one = 1.0d+0 )
197 COMPLEX*16 CZERO, CONE
198 parameter( czero = ( 0.0d+0, 0.0d+0 ),
199 $ cone = ( 1.0d+0, 0.0d+0 ) )
200 COMPLEX*16 CROGUE
201 parameter( crogue = ( -1.0d+10, 0.0d+0 ) )
202* ..
203* .. Local Scalars ..
204 INTEGER INFO
205 DOUBLE PRECISION ANORM, BNORM, RESID, ULP, UNFL
206* ..
207* .. External Functions ..
208 DOUBLE PRECISION DLAMCH, ZLANGE, ZLANHE
209 EXTERNAL dlamch, zlange, zlanhe
210* ..
211* .. External Subroutines ..
212 EXTERNAL zgemm, zggqrf, zherk, zlacpy, zlaset, zungqr,
213 $ zungrq
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC dble, max, min
217* ..
218* .. Executable Statements ..
219*
220 ulp = dlamch( 'Precision' )
221 unfl = dlamch( 'Safe minimum' )
222*
223* Copy the matrix A to the array AF.
224*
225 CALL zlacpy( 'Full', n, m, a, lda, af, lda )
226 CALL zlacpy( 'Full', n, p, b, ldb, bf, ldb )
227*
228 anorm = max( zlange( '1', n, m, a, lda, rwork ), unfl )
229 bnorm = max( zlange( '1', n, p, b, ldb, rwork ), unfl )
230*
231* Factorize the matrices A and B in the arrays AF and BF.
232*
233 CALL zggqrf( n, m, p, af, lda, taua, bf, ldb, taub, work, lwork,
234 $ info )
235*
236* Generate the N-by-N matrix Q
237*
238 CALL zlaset( 'Full', n, n, crogue, crogue, q, lda )
239 CALL zlacpy( 'Lower', n-1, m, af( 2, 1 ), lda, q( 2, 1 ), lda )
240 CALL zungqr( n, n, min( n, m ), q, lda, taua, work, lwork, info )
241*
242* Generate the P-by-P matrix Z
243*
244 CALL zlaset( 'Full', p, p, crogue, crogue, z, ldb )
245 IF( n.LE.p ) THEN
246 IF( n.GT.0 .AND. n.LT.p )
247 $ CALL zlacpy( 'Full', n, p-n, bf, ldb, z( p-n+1, 1 ), ldb )
248 IF( n.GT.1 )
249 $ CALL zlacpy( 'Lower', n-1, n-1, bf( 2, p-n+1 ), ldb,
250 $ z( p-n+2, p-n+1 ), ldb )
251 ELSE
252 IF( p.GT.1 )
253 $ CALL zlacpy( 'Lower', p-1, p-1, bf( n-p+2, 1 ), ldb,
254 $ z( 2, 1 ), ldb )
255 END IF
256 CALL zungrq( p, p, min( n, p ), z, ldb, taub, work, lwork, info )
257*
258* Copy R
259*
260 CALL zlaset( 'Full', n, m, czero, czero, r, lda )
261 CALL zlacpy( 'Upper', n, m, af, lda, r, lda )
262*
263* Copy T
264*
265 CALL zlaset( 'Full', n, p, czero, czero, t, ldb )
266 IF( n.LE.p ) THEN
267 CALL zlacpy( 'Upper', n, n, bf( 1, p-n+1 ), ldb, t( 1, p-n+1 ),
268 $ ldb )
269 ELSE
270 CALL zlacpy( 'Full', n-p, p, bf, ldb, t, ldb )
271 CALL zlacpy( 'Upper', p, p, bf( n-p+1, 1 ), ldb, t( n-p+1, 1 ),
272 $ ldb )
273 END IF
274*
275* Compute R - Q'*A
276*
277 CALL zgemm( 'Conjugate transpose', 'No transpose', n, m, n, -cone,
278 $ q, lda, a, lda, cone, r, lda )
279*
280* Compute norm( R - Q'*A ) / ( MAX(M,N)*norm(A)*ULP ) .
281*
282 resid = zlange( '1', n, m, r, lda, rwork )
283 IF( anorm.GT.zero ) THEN
284 result( 1 ) = ( ( resid / dble( max( 1, m, n ) ) ) / anorm ) /
285 $ ulp
286 ELSE
287 result( 1 ) = zero
288 END IF
289*
290* Compute T*Z - Q'*B
291*
292 CALL zgemm( 'No Transpose', 'No transpose', n, p, p, cone, t, ldb,
293 $ z, ldb, czero, bwk, ldb )
294 CALL zgemm( 'Conjugate transpose', 'No transpose', n, p, n, -cone,
295 $ q, lda, b, ldb, cone, bwk, ldb )
296*
297* Compute norm( T*Z - Q'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
298*
299 resid = zlange( '1', n, p, bwk, ldb, rwork )
300 IF( bnorm.GT.zero ) THEN
301 result( 2 ) = ( ( resid / dble( max( 1, p, n ) ) ) / bnorm ) /
302 $ ulp
303 ELSE
304 result( 2 ) = zero
305 END IF
306*
307* Compute I - Q'*Q
308*
309 CALL zlaset( 'Full', n, n, czero, cone, r, lda )
310 CALL zherk( 'Upper', 'Conjugate transpose', n, n, -one, q, lda,
311 $ one, r, lda )
312*
313* Compute norm( I - Q'*Q ) / ( N * ULP ) .
314*
315 resid = zlanhe( '1', 'Upper', n, r, lda, rwork )
316 result( 3 ) = ( resid / dble( max( 1, n ) ) ) / ulp
317*
318* Compute I - Z'*Z
319*
320 CALL zlaset( 'Full', p, p, czero, cone, t, ldb )
321 CALL zherk( 'Upper', 'Conjugate transpose', p, p, -one, z, ldb,
322 $ one, t, ldb )
323*
324* Compute norm( I - Z'*Z ) / ( P*ULP ) .
325*
326 resid = zlanhe( '1', 'Upper', p, t, ldb, rwork )
327 result( 4 ) = ( resid / dble( max( 1, p ) ) ) / ulp
328*
329 RETURN
330*
331* End of ZGQRTS
332*
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zggqrf(n, m, p, a, lda, taua, b, ldb, taub, work, lwork, info)
ZGGQRF
Definition zggqrf.f:215
subroutine zherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
ZHERK
Definition zherk.f:173
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
double precision function zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:124
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:128
subroutine zungrq(m, n, k, a, lda, tau, work, lwork, info)
ZUNGRQ
Definition zungrq.f:128
Here is the call graph for this function:
Here is the caller graph for this function: