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

◆ slahilb()

subroutine slahilb ( integer  n,
integer  nrhs,
real, dimension(lda, n)  a,
integer  lda,
real, dimension(ldx, nrhs)  x,
integer  ldx,
real, dimension(ldb, nrhs)  b,
integer  ldb,
real, dimension(n)  work,
integer  info 
)

SLAHILB

Purpose:
 SLAHILB 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 REAL 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 REAL 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 REAL 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 REAL 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 slahilb.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 REAL 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* .. External Subroutines ..
150 EXTERNAL xerbla
151* ..
152* .. External Functions
153 EXTERNAL slaset
154 INTRINSIC real
155* ..
156* .. Executable Statements ..
157*
158* Test the input arguments
159*
160 info = 0
161 IF (n .LT. 0 .OR. n .GT. nmax_approx) THEN
162 info = -1
163 ELSE IF (nrhs .LT. 0) THEN
164 info = -2
165 ELSE IF (lda .LT. n) THEN
166 info = -4
167 ELSE IF (ldx .LT. n) THEN
168 info = -6
169 ELSE IF (ldb .LT. n) THEN
170 info = -8
171 END IF
172 IF (info .LT. 0) THEN
173 CALL xerbla('SLAHILB', -info)
174 RETURN
175 END IF
176 IF (n .GT. nmax_exact) THEN
177 info = 1
178 END IF
179*
180* Compute M = the LCM of the integers [1, 2*N-1]. The largest
181* reasonable N is small enough that integers suffice (up to N = 11).
182 m = 1
183 DO i = 2, (2*n-1)
184 tm = m
185 ti = i
186 r = mod(tm, ti)
187 DO WHILE (r .NE. 0)
188 tm = ti
189 ti = r
190 r = mod(tm, ti)
191 END DO
192 m = (m / ti) * i
193 END DO
194*
195* Generate the scaled Hilbert matrix in A
196 DO j = 1, n
197 DO i = 1, n
198 a(i, j) = real(m) / (i + j - 1)
199 END DO
200 END DO
201*
202* Generate matrix B as simply the first NRHS columns of M * the
203* identity.
204 CALL slaset('Full', n, nrhs, 0.0, real(m), b, ldb)
205*
206* Generate the true solutions in X. Because B = the first NRHS
207* columns of M*I, the true solutions are just the first NRHS columns
208* of the inverse Hilbert matrix.
209 work(1) = n
210 DO j = 2, n
211 work(j) = ( ( (work(j-1)/(j-1)) * (j-1 - n) ) /(j-1) )
212 $ * (n +j -1)
213 END DO
214*
215 DO j = 1, nrhs
216 DO i = 1, n
217 x(i, j) = (work(i)*work(j)) / (i + j - 1)
218 END DO
219 END DO
220*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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
Here is the call graph for this function: