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

◆ dglmts()

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

DGLMTS

Purpose:
 DGLMTS tests DGGGLM - 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 DOUBLE PRECISION array, dimension (LDA,M)
          The N-by-M matrix A.
[out]AF
          AF is DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (LDB,P)
          The N-by-P matrix A.
[out]BF
          BF is DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension( N )
          On input, the left hand side of the GLM.
[out]DF
          DF is DOUBLE PRECISION array, dimension( N )
[out]X
          X is DOUBLE PRECISION array, dimension( M )
          solution vector X in the GLM problem.
[out]U
          U is DOUBLE PRECISION array, dimension( P )
          solution vector U in the GLM problem.
[out]WORK
          WORK is DOUBLE PRECISION 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 dglmts.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 A( LDA, * ), AF( LDA, * ), B( LDB, * ),
160 $ BF( LDB, * ), D( * ), DF( * ), RWORK( * ),
161 $ U( * ), WORK( LWORK ), X( * )
162* ..
163* .. Parameters ..
164 DOUBLE PRECISION ZERO, ONE
165 parameter( zero = 0.0d+0, one = 1.0d+0 )
166* ..
167* .. Local Scalars ..
168 INTEGER INFO
169 DOUBLE PRECISION ANORM, BNORM, DNORM, EPS, UNFL, XNORM, YNORM
170* ..
171* .. External Functions ..
172 DOUBLE PRECISION DASUM, DLAMCH, DLANGE
173 EXTERNAL dasum, dlamch, dlange
174* ..
175* .. External Subroutines ..
176*
177 EXTERNAL dcopy, dgemv, dggglm, dlacpy
178* ..
179* .. Intrinsic Functions ..
180 INTRINSIC max
181* ..
182* .. Executable Statements ..
183*
184 eps = dlamch( 'Epsilon' )
185 unfl = dlamch( 'Safe minimum' )
186 anorm = max( dlange( '1', n, m, a, lda, rwork ), unfl )
187 bnorm = max( dlange( '1', n, p, b, ldb, rwork ), unfl )
188*
189* Copy the matrices A and B to the arrays AF and BF,
190* and the vector D the array DF.
191*
192 CALL dlacpy( 'Full', n, m, a, lda, af, lda )
193 CALL dlacpy( 'Full', n, p, b, ldb, bf, ldb )
194 CALL dcopy( n, d, 1, df, 1 )
195*
196* Solve GLM problem
197*
198 CALL dggglm( n, m, p, af, lda, bf, ldb, df, x, u, work, lwork,
199 $ info )
200*
201* Test the residual for the solution of LSE
202*
203* norm( d - A*x - B*u )
204* RESULT = -----------------------------------------
205* (norm(A)+norm(B))*(norm(x)+norm(u))*EPS
206*
207 CALL dcopy( n, d, 1, df, 1 )
208 CALL dgemv( 'No transpose', n, m, -one, a, lda, x, 1, one, df, 1 )
209*
210 CALL dgemv( 'No transpose', n, p, -one, b, ldb, u, 1, one, df, 1 )
211*
212 dnorm = dasum( n, df, 1 )
213 xnorm = dasum( m, x, 1 ) + dasum( p, u, 1 )
214 ynorm = anorm + bnorm
215*
216 IF( xnorm.LE.zero ) THEN
217 result = zero
218 ELSE
219 result = ( ( dnorm / ynorm ) / xnorm ) / eps
220 END IF
221*
222 RETURN
223*
224* End of DGLMTS
225*
double precision function dasum(n, dx, incx)
DASUM
Definition dasum.f:71
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine dggglm(n, m, p, a, lda, b, ldb, d, x, y, work, lwork, info)
DGGGLM
Definition dggglm.f:185
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
Here is the call graph for this function:
Here is the caller graph for this function: