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

◆ zglmts()

subroutine zglmts ( integer n,
integer m,
integer p,
complex*16, dimension( lda, * ) a,
complex*16, dimension( lda, * ) af,
integer lda,
complex*16, dimension( ldb, * ) b,
complex*16, dimension( ldb, * ) bf,
integer ldb,
complex*16, dimension( * ) d,
complex*16, dimension( * ) df,
complex*16, dimension( * ) x,
complex*16, dimension( * ) u,
complex*16, dimension( lwork ) work,
integer lwork,
double precision, dimension( * ) rwork,
double precision result )

ZGLMTS

Purpose:
!>
!> ZGLMTS tests ZGGGLM - a subroutine for solving the generalized
!> linear model problem.
!> 
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,M)
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the arrays A, AF. LDA >= max(M,N).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension (LDB,P)
!>          The N-by-P matrix A.
!> 
[out]BF
!>          BF is COMPLEX*16 array, dimension (LDB,P)
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the arrays B, BF. LDB >= max(P,N).
!> 
[in]D
!>          D is COMPLEX*16 array, dimension( N )
!>          On input, the left hand side of the GLM.
!> 
[out]DF
!>          DF is COMPLEX*16 array, dimension( N )
!> 
[out]X
!>          X is COMPLEX*16 array, dimension( M )
!>          solution vector X in the GLM problem.
!> 
[out]U
!>          U is COMPLEX*16 array, dimension( P )
!>          solution vector U in the GLM problem.
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (LWORK)
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (M)
!> 
[out]RESULT
!>          RESULT is DOUBLE PRECISION
!>          The test ratio:
!>                           norm( d - A*x - B*u )
!>            RESULT = -----------------------------------------
!>                     (norm(A)+norm(B))*(norm(x)+norm(u))*EPS
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 144 of file zglmts.f.

146*
147* -- LAPACK test routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 INTEGER LDA, LDB, LWORK, M, N, P
153 DOUBLE PRECISION RESULT
154* ..
155* .. Array Arguments ..
156*
157* ====================================================================
158*
159 DOUBLE PRECISION RWORK( * )
160 COMPLEX*16 A( LDA, * ), AF( LDA, * ), B( LDB, * ),
161 $ BF( LDB, * ), D( * ), DF( * ), U( * ),
162 $ WORK( LWORK ), X( * )
163* ..
164* .. Parameters ..
165 DOUBLE PRECISION ZERO
166 parameter( zero = 0.0d+0 )
167 COMPLEX*16 CONE
168 parameter( cone = 1.0d+0 )
169* ..
170* .. Local Scalars ..
171 INTEGER INFO
172 DOUBLE PRECISION ANORM, BNORM, DNORM, EPS, UNFL, XNORM, YNORM
173* ..
174* .. External Functions ..
175 DOUBLE PRECISION DLAMCH, DZASUM, ZLANGE
176 EXTERNAL dlamch, dzasum, zlange
177* ..
178* .. External Subroutines ..
179*
180 EXTERNAL zcopy, zgemv, zggglm, zlacpy
181* ..
182* .. Intrinsic Functions ..
183 INTRINSIC max
184* ..
185* .. Executable Statements ..
186*
187 eps = dlamch( 'Epsilon' )
188 unfl = dlamch( 'Safe minimum' )
189 anorm = max( zlange( '1', n, m, a, lda, rwork ), unfl )
190 bnorm = max( zlange( '1', n, p, b, ldb, rwork ), unfl )
191*
192* Copy the matrices A and B to the arrays AF and BF,
193* and the vector D the array DF.
194*
195 CALL zlacpy( 'Full', n, m, a, lda, af, lda )
196 CALL zlacpy( 'Full', n, p, b, ldb, bf, ldb )
197 CALL zcopy( n, d, 1, df, 1 )
198*
199* Solve GLM problem
200*
201 CALL zggglm( n, m, p, af, lda, bf, ldb, df, x, u, work, lwork,
202 $ info )
203*
204* Test the residual for the solution of LSE
205*
206* norm( d - A*x - B*u )
207* RESULT = -----------------------------------------
208* (norm(A)+norm(B))*(norm(x)+norm(u))*EPS
209*
210 CALL zcopy( n, d, 1, df, 1 )
211 CALL zgemv( 'No transpose', n, m, -cone, a, lda, x, 1, cone, df,
212 $ 1 )
213*
214 CALL zgemv( 'No transpose', n, p, -cone, b, ldb, u, 1, cone, df,
215 $ 1 )
216*
217 dnorm = dzasum( n, df, 1 )
218 xnorm = dzasum( m, x, 1 ) + dzasum( p, u, 1 )
219 ynorm = anorm + bnorm
220*
221 IF( xnorm.LE.zero ) THEN
222 result = zero
223 ELSE
224 result = ( ( dnorm / ynorm ) / xnorm ) / eps
225 END IF
226*
227 RETURN
228*
229* End of ZGLMTS
230*
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
subroutine zggglm(n, m, p, a, lda, b, ldb, d, x, y, work, lwork, info)
ZGGGLM
Definition zggglm.f:194
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
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:113
Here is the call graph for this function:
Here is the caller graph for this function: