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

◆ dsgesv()

subroutine dsgesv ( integer n,
integer nrhs,
double precision, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) ipiv,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldx, * ) x,
integer ldx,
double precision, dimension( n, * ) work,
real, dimension( * ) swork,
integer iter,
integer info )

DSGESV computes the solution to system of linear equations A * X = B for GE matrices (mixed precision with iterative refinement)

Download DSGESV + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DSGESV computes the solution to a real system of linear equations
!>    A * X = B,
!> where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
!>
!> DSGESV first attempts to factorize the matrix in SINGLE PRECISION
!> and use this factorization within an iterative refinement procedure
!> to produce a solution with DOUBLE PRECISION normwise backward error
!> quality (see below). If the approach fails the method switches to a
!> DOUBLE PRECISION factorization and solve.
!>
!> The iterative refinement is not going to be a winning strategy if
!> the ratio SINGLE PRECISION performance over DOUBLE PRECISION
!> performance is too small. A reasonable strategy should take the
!> number of right-hand sides and the size of the matrix into account.
!> This might be done with a call to ILAENV in the future. Up to now, we
!> always try iterative refinement.
!>
!> The iterative refinement process is stopped if
!>     ITER > ITERMAX
!> or for all the RHS we have:
!>     RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX
!> where
!>     o ITER is the number of the current iteration in the iterative
!>       refinement process
!>     o RNRM is the infinity-norm of the residual
!>     o XNRM is the infinity-norm of the solution
!>     o ANRM is the infinity-operator-norm of the matrix A
!>     o EPS is the machine epsilon returned by DLAMCH('Epsilon')
!> The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00
!> respectively.
!> 
Parameters
[in]N
!>          N is INTEGER
!>          The number of linear equations, i.e., the order of the
!>          matrix A.  N >= 0.
!> 
[in]NRHS
!>          NRHS is INTEGER
!>          The number of right hand sides, i.e., the number of columns
!>          of the matrix B.  NRHS >= 0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array,
!>          dimension (LDA,N)
!>          On entry, the N-by-N coefficient matrix A.
!>          On exit, if iterative refinement has been successfully used
!>          (INFO = 0 and ITER >= 0, see description below), then A is
!>          unchanged, if double precision factorization has been used
!>          (INFO = 0 and ITER < 0, see description below), then the
!>          array A contains the factors L and U from the factorization
!>          A = P*L*U; the unit diagonal elements of L are not stored.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          The pivot indices that define the permutation matrix P;
!>          row i of the matrix was interchanged with row IPIV(i).
!>          Corresponds either to the single precision factorization
!>          (if INFO = 0 and ITER >= 0) or the double precision
!>          factorization (if INFO = 0 and ITER < 0).
!> 
[in]B
!>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
!>          The N-by-NRHS right hand side matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
!>          If INFO = 0, the N-by-NRHS solution matrix X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (N,NRHS)
!>          This array is used to hold the residual vectors.
!> 
[out]SWORK
!>          SWORK is REAL array, dimension (N*(N+NRHS))
!>          This array is used to use the single precision matrix and the
!>          right-hand sides or solutions in single precision.
!> 
[out]ITER
!>          ITER is INTEGER
!>          < 0: iterative refinement has failed, double precision
!>               factorization has been performed
!>               -1 : the routine fell back to full precision for
!>                    implementation- or machine-specific reasons
!>               -2 : narrowing the precision induced an overflow,
!>                    the routine fell back to full precision
!>               -3 : failure of SGETRF
!>               -31: stop the iterative refinement after the 30th
!>                    iterations
!>          > 0: iterative refinement has been successfully used.
!>               Returns the number of iterations
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  if INFO = i, U(i,i) computed in DOUBLE PRECISION is
!>                exactly zero.  The factorization has been completed,
!>                but the factor U is exactly singular, so the solution
!>                could not be computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 191 of file dsgesv.f.

193*
194* -- LAPACK driver routine --
195* -- LAPACK is a software package provided by Univ. of Tennessee, --
196* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
197*
198* .. Scalar Arguments ..
199 INTEGER INFO, ITER, LDA, LDB, LDX, N, NRHS
200* ..
201* .. Array Arguments ..
202 INTEGER IPIV( * )
203 REAL SWORK( * )
204 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( N, * ),
205 $ X( LDX, * )
206* ..
207*
208* =====================================================================
209*
210* .. Parameters ..
211 LOGICAL DOITREF
212 parameter( doitref = .true. )
213*
214 INTEGER ITERMAX
215 parameter( itermax = 30 )
216*
217 DOUBLE PRECISION BWDMAX
218 parameter( bwdmax = 1.0e+00 )
219*
220 DOUBLE PRECISION NEGONE, ONE
221 parameter( negone = -1.0d+0, one = 1.0d+0 )
222*
223* .. Local Scalars ..
224 INTEGER I, IITER, PTSA, PTSX
225 DOUBLE PRECISION ANRM, CTE, EPS, RNRM, XNRM
226*
227* .. External Subroutines ..
228 EXTERNAL daxpy, dgemm, dlacpy, dlag2s, dgetrf,
229 $ dgetrs,
231* ..
232* .. External Functions ..
233 INTEGER IDAMAX
234 DOUBLE PRECISION DLAMCH, DLANGE
235 EXTERNAL idamax, dlamch, dlange
236* ..
237* .. Intrinsic Functions ..
238 INTRINSIC abs, dble, max, sqrt
239* ..
240* .. Executable Statements ..
241*
242 info = 0
243 iter = 0
244*
245* Test the input parameters.
246*
247 IF( n.LT.0 ) THEN
248 info = -1
249 ELSE IF( nrhs.LT.0 ) THEN
250 info = -2
251 ELSE IF( lda.LT.max( 1, n ) ) THEN
252 info = -4
253 ELSE IF( ldb.LT.max( 1, n ) ) THEN
254 info = -7
255 ELSE IF( ldx.LT.max( 1, n ) ) THEN
256 info = -9
257 END IF
258 IF( info.NE.0 ) THEN
259 CALL xerbla( 'DSGESV', -info )
260 RETURN
261 END IF
262*
263* Quick return if (N.EQ.0).
264*
265 IF( n.EQ.0 )
266 $ RETURN
267*
268* Skip single precision iterative refinement if a priori slower
269* than double precision factorization.
270*
271 IF( .NOT.doitref ) THEN
272 iter = -1
273 GO TO 40
274 END IF
275*
276* Compute some constants.
277*
278 anrm = dlange( 'I', n, n, a, lda, work )
279 eps = dlamch( 'Epsilon' )
280 cte = anrm*eps*sqrt( dble( n ) )*bwdmax
281*
282* Set the indices PTSA, PTSX for referencing SA and SX in SWORK.
283*
284 ptsa = 1
285 ptsx = ptsa + n*n
286*
287* Convert B from double precision to single precision and store the
288* result in SX.
289*
290 CALL dlag2s( n, nrhs, b, ldb, swork( ptsx ), n, info )
291*
292 IF( info.NE.0 ) THEN
293 iter = -2
294 GO TO 40
295 END IF
296*
297* Convert A from double precision to single precision and store the
298* result in SA.
299*
300 CALL dlag2s( n, n, a, lda, swork( ptsa ), n, info )
301*
302 IF( info.NE.0 ) THEN
303 iter = -2
304 GO TO 40
305 END IF
306*
307* Compute the LU factorization of SA.
308*
309 CALL sgetrf( n, n, swork( ptsa ), n, ipiv, info )
310*
311 IF( info.NE.0 ) THEN
312 iter = -3
313 GO TO 40
314 END IF
315*
316* Solve the system SA*SX = SB.
317*
318 CALL sgetrs( 'No transpose', n, nrhs, swork( ptsa ), n, ipiv,
319 $ swork( ptsx ), n, info )
320*
321* Convert SX back to double precision
322*
323 CALL slag2d( n, nrhs, swork( ptsx ), n, x, ldx, info )
324*
325* Compute R = B - AX (R is WORK).
326*
327 CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
328*
329 CALL dgemm( 'No Transpose', 'No Transpose', n, nrhs, n, negone,
330 $ a,
331 $ lda, x, ldx, one, work, n )
332*
333* Check whether the NRHS normwise backward errors satisfy the
334* stopping criterion. If yes, set ITER=0 and return.
335*
336 DO i = 1, nrhs
337 xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
338 rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
339 IF( rnrm.GT.xnrm*cte )
340 $ GO TO 10
341 END DO
342*
343* If we are here, the NRHS normwise backward errors satisfy the
344* stopping criterion. We are good to exit.
345*
346 iter = 0
347 RETURN
348*
349 10 CONTINUE
350*
351 DO 30 iiter = 1, itermax
352*
353* Convert R (in WORK) from double precision to single precision
354* and store the result in SX.
355*
356 CALL dlag2s( n, nrhs, work, n, swork( ptsx ), n, info )
357*
358 IF( info.NE.0 ) THEN
359 iter = -2
360 GO TO 40
361 END IF
362*
363* Solve the system SA*SX = SR.
364*
365 CALL sgetrs( 'No transpose', n, nrhs, swork( ptsa ), n,
366 $ ipiv,
367 $ swork( ptsx ), n, info )
368*
369* Convert SX back to double precision and update the current
370* iterate.
371*
372 CALL slag2d( n, nrhs, swork( ptsx ), n, work, n, info )
373*
374 DO i = 1, nrhs
375 CALL daxpy( n, one, work( 1, i ), 1, x( 1, i ), 1 )
376 END DO
377*
378* Compute R = B - AX (R is WORK).
379*
380 CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
381*
382 CALL dgemm( 'No Transpose', 'No Transpose', n, nrhs, n,
383 $ negone,
384 $ a, lda, x, ldx, one, work, n )
385*
386* Check whether the NRHS normwise backward errors satisfy the
387* stopping criterion. If yes, set ITER=IITER>0 and return.
388*
389 DO i = 1, nrhs
390 xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
391 rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
392 IF( rnrm.GT.xnrm*cte )
393 $ GO TO 20
394 END DO
395*
396* If we are here, the NRHS normwise backward errors satisfy the
397* stopping criterion, we are good to exit.
398*
399 iter = iiter
400*
401 RETURN
402*
403 20 CONTINUE
404*
405 30 CONTINUE
406*
407* If we are at this place of the code, this is because we have
408* performed ITER=ITERMAX iterations and never satisfied the
409* stopping criterion, set up the ITER flag accordingly and follow up
410* on double precision routine.
411*
412 iter = -itermax - 1
413*
414 40 CONTINUE
415*
416* Single-precision iterative refinement failed to converge to a
417* satisfactory solution, so we resort to double precision.
418*
419 CALL dgetrf( n, n, a, lda, ipiv, info )
420*
421 IF( info.NE.0 )
422 $ RETURN
423*
424 CALL dlacpy( 'All', n, nrhs, b, ldb, x, ldx )
425 CALL dgetrs( 'No transpose', n, nrhs, a, lda, ipiv, x, ldx,
426 $ info )
427*
428 RETURN
429*
430* End of DSGESV
431*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlag2s(m, n, a, lda, sa, ldsa, info)
DLAG2S converts a double precision matrix to a single precision matrix.
Definition dlag2s.f:106
subroutine slag2d(m, n, sa, ldsa, a, lda, info)
SLAG2D converts a single precision matrix to a double precision matrix.
Definition slag2d.f:102
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine sgetrf(m, n, a, lda, ipiv, info)
SGETRF
Definition sgetrf.f:106
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
Definition dgetrf.f:106
subroutine sgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
SGETRS
Definition sgetrs.f:119
subroutine dgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
DGETRS
Definition dgetrs.f:119
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
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:112
Here is the call graph for this function:
Here is the caller graph for this function: