SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
|
Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates. More...
Macros | |
#define | ITMAX 20 |
Functions | |
void | psgsrfs (int_t n, SuperMatrix *A, float anorm, sLUstruct_t *LUstruct, sScalePermstruct_t *ScalePermstruct, gridinfo_t *grid, float *B, int_t ldb, float *X, int_t ldx, int nrhs, sSOLVEstruct_t *SOLVEstruct, float *berr, SuperLUStat_t *stat, int *info) |
Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates.
Copyright (c) 2003, The Regents of the University of California, through Lawrence Berkeley National Laboratory (subject to receipt of any required approvals from U.S. Dept. of Energy)
All rights reserved.
The source code is distributed under BSD license, see the file License.txt at the top-level directory.
-- Distributed SuperLU routine (version 4.3) -- Lawrence Berkeley National Lab, Univ. of California Berkeley. March 15, 2003 Last modified: December 31, 2015
#define ITMAX 20 |
void psgsrfs | ( | int_t | n, |
SuperMatrix * | A, | ||
float | anorm, | ||
sLUstruct_t * | LUstruct, | ||
sScalePermstruct_t * | ScalePermstruct, | ||
gridinfo_t * | grid, | ||
float * | B, | ||
int_t | ldb, | ||
float * | X, | ||
int_t | ldx, | ||
int | nrhs, | ||
sSOLVEstruct_t * | SOLVEstruct, | ||
float * | berr, | ||
SuperLUStat_t * | stat, | ||
int * | info | ||
) |
Purpose ======= PSGSRFS improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution. Arguments ========= n (input) int (global) The order of the system of linear equations. A (input) SuperMatrix* The original matrix A, or the scaled A if equilibration was done. A is also permuted into diag(R)*A*diag(C)*Pc'. The type of A can be: Stype = SLU_NR_loc; Dtype = SLU_S; Mtype = SLU_GE. anorm (input) double The norm of the original matrix A, or the scaled A if equilibration was done. LUstruct (input) sLUstruct_t* The distributed data structures storing L and U factors. The L and U factors are obtained from pdgstrf for the possibly scaled and permuted matrix A. See superlu_sdefs.h for the definition of 'sLUstruct_t'. ScalePermstruct (input) sScalePermstruct_t* (global) The data structure to store the scaling and permutation vectors describing the transformations performed to the matrix A. grid (input) gridinfo_t* The 2D process mesh. It contains the MPI communicator, the number of process rows (NPROW), the number of process columns (NPCOL), and my process rank. It is an input argument to all the parallel routines. Grid can be initialized by subroutine SUPERLU_GRIDINIT. See superlu_defs.h for the definition of 'gridinfo_t'. B (input) float* (local) The m_loc-by-NRHS right-hand side matrix of the possibly equilibrated system. That is, B may be overwritten by diag(R)*B. ldb (input) int (local) Leading dimension of matrix B. X (input/output) float* (local) On entry, the solution matrix Y, as computed by PDGSTRS, of the transformed system A1*Y = Pc*Pr*B. where A1 = Pc*Pr*diag(R)*A*diag(C)*Pc' and Y = Pc*diag(C)^(-1)*X. On exit, the improved solution matrix Y. In order to obtain the solution X to the original system, Y should be permutated by Pc^T, and premultiplied by diag(C) if DiagScale = COL or BOTH. This must be done after this routine is called. ldx (input) int (local) Leading dimension of matrix X. nrhs (input) int Number of right-hand sides. SOLVEstruct (output) sSOLVEstruct_t* (global) Contains the information for the communication during the solution phase. berr (output) float*, dimension (nrhs) The componentwise relative backward error of each solution vector X(j) (i.e., the smallest relative change in any element of A or B that makes X(j) an exact solution). stat (output) SuperLUStat_t* Record the statistics about the refinement steps. See util.h for the definition of SuperLUStat_t. info (output) int* = 0: successful exit < 0: if info = -i, the i-th argument had an illegal value Internal Parameters =================== ITMAX is the maximum number of steps of iterative refinement.