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

◆ dlahilb()

subroutine dlahilb ( integer n,
integer nrhs,
double precision, dimension(lda, n) a,
integer lda,
double precision, dimension(ldx, nrhs) x,
integer ldx,
double precision, dimension(ldb, nrhs) b,
integer ldb,
double precision, dimension(n) work,
integer info )

DLAHILB

Purpose:
!>
!> DLAHILB generates an N by N scaled Hilbert matrix in A along with
!> NRHS right-hand sides in B and solutions in X such that A*X=B.
!>
!> The Hilbert matrix is scaled by M = LCM(1, 2, ..., 2*N-1) so that all
!> entries are integers.  The right-hand sides are the first NRHS
!> columns of M * the identity matrix, and the solutions are the
!> first NRHS columns of the inverse Hilbert matrix.
!>
!> The condition number of the Hilbert matrix grows exponentially with
!> its size, roughly as O(e ** (3.5*N)).  Additionally, the inverse
!> Hilbert matrices beyond a relatively small dimension cannot be
!> generated exactly without extra precision.  Precision is exhausted
!> when the largest entry in the inverse Hilbert matrix is greater than
!> 2 to the power of the number of bits in the fraction of the data type
!> used plus one, which is 24 for single precision.
!>
!> In single, the generated solution is exact for N <= 6 and has
!> small componentwise error for 7 <= N <= 11.
!> 
Parameters
[in]N
!>          N is INTEGER
!>          The dimension of the matrix A.
!> 
[in]NRHS
!>          NRHS is INTEGER
!>          The requested number of right-hand sides.
!> 
[out]A
!>          A is DOUBLE PRECISION array, dimension (LDA, N)
!>          The generated scaled Hilbert matrix.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= N.
!> 
[out]X
!>          X is DOUBLE PRECISION array, dimension (LDX, NRHS)
!>          The generated exact solutions.  Currently, the first NRHS
!>          columns of the inverse Hilbert matrix.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= N.
!> 
[out]B
!>          B is DOUBLE PRECISION array, dimension (LDB, NRHS)
!>          The generated right-hand sides.  Currently, the first NRHS
!>          columns of LCM(1, 2, ..., 2*N-1) * the identity matrix.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= N.
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          = 1: N is too large; the data is still generated but may not
!>               be not exact.
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 123 of file dlahilb.f.

124*
125* -- LAPACK test routine --
126* -- LAPACK is a software package provided by Univ. of Tennessee, --
127* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
128*
129* .. Scalar Arguments ..
130 INTEGER N, NRHS, LDA, LDX, LDB, INFO
131* .. Array Arguments ..
132 DOUBLE PRECISION A(LDA, N), X(LDX, NRHS), B(LDB, NRHS), WORK(N)
133* ..
134*
135* =====================================================================
136* .. Local Scalars ..
137 INTEGER TM, TI, R
138 INTEGER M
139 INTEGER I, J
140* ..
141* .. Parameters ..
142* NMAX_EXACT the largest dimension where the generated data is
143* exact.
144* NMAX_APPROX the largest dimension where the generated data has
145* a small componentwise relative error.
146 INTEGER NMAX_EXACT, NMAX_APPROX
147 parameter(nmax_exact = 6, nmax_approx = 11)
148
149* ..
150* .. External Functions
151 EXTERNAL dlaset
152 INTRINSIC dble
153* ..
154* .. Executable Statements ..
155*
156* Test the input arguments
157*
158 info = 0
159 IF (n .LT. 0 .OR. n .GT. nmax_approx) THEN
160 info = -1
161 ELSE IF (nrhs .LT. 0) THEN
162 info = -2
163 ELSE IF (lda .LT. n) THEN
164 info = -4
165 ELSE IF (ldx .LT. n) THEN
166 info = -6
167 ELSE IF (ldb .LT. n) THEN
168 info = -8
169 END IF
170 IF (info .LT. 0) THEN
171 CALL xerbla('DLAHILB', -info)
172 RETURN
173 END IF
174 IF (n .GT. nmax_exact) THEN
175 info = 1
176 END IF
177*
178* Compute M = the LCM of the integers [1, 2*N-1]. The largest
179* reasonable N is small enough that integers suffice (up to N = 11).
180 m = 1
181 DO i = 2, (2*n-1)
182 tm = m
183 ti = i
184 r = mod(tm, ti)
185 DO WHILE (r .NE. 0)
186 tm = ti
187 ti = r
188 r = mod(tm, ti)
189 END DO
190 m = (m / ti) * i
191 END DO
192*
193* Generate the scaled Hilbert matrix in A
194 DO j = 1, n
195 DO i = 1, n
196 a(i, j) = dble(m) / (i + j - 1)
197 END DO
198 END DO
199*
200* Generate matrix B as simply the first NRHS columns of M * the
201* identity.
202 CALL dlaset('Full', n, nrhs, 0.0d+0, dble(m), b, ldb)
203
204* Generate the true solutions in X. Because B = the first NRHS
205* columns of M*I, the true solutions are just the first NRHS columns
206* of the inverse Hilbert matrix.
207 work(1) = n
208 DO j = 2, n
209 work(j) = ( ( (work(j-1)/(j-1)) * (j-1 - n) ) /(j-1) )
210 $ * (n +j -1)
211 END DO
212*
213 DO j = 1, nrhs
214 DO i = 1, n
215 x(i, j) = (work(i)*work(j)) / (i + j - 1)
216 END DO
217 END DO
218*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:108
Here is the call graph for this function:
Here is the caller graph for this function: