SuperLU Distributed 8.2.1
Distributed memory sparse direct solver
pzgsrfs.c File Reference

Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates. More...

#include <math.h>
#include "superlu_zdefs.h"
Include dependency graph for pzgsrfs.c:

Macros

#define ITMAX   20
 

Functions

void pzgsrfs (superlu_dist_options_t *options, int_t n, SuperMatrix *A, double anorm, zLUstruct_t *LUstruct, zScalePermstruct_t *ScalePermstruct, gridinfo_t *grid, doublecomplex *B, int_t ldb, doublecomplex *X, int_t ldx, int nrhs, zSOLVEstruct_t *SOLVEstruct, double *berr, SuperLUStat_t *stat, int *info)
 

Detailed Description

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

Macro Definition Documentation

◆ ITMAX

#define ITMAX   20

Function Documentation

◆ pzgsrfs()

void pzgsrfs ( superlu_dist_options_t options,
int_t  n,
SuperMatrix A,
double  anorm,
zLUstruct_t LUstruct,
zScalePermstruct_t ScalePermstruct,
gridinfo_t grid,
doublecomplex B,
int_t  ldb,
doublecomplex X,
int_t  ldx,
int  nrhs,
zSOLVEstruct_t SOLVEstruct,
double *  berr,
SuperLUStat_t stat,
int *  info 
)
Purpose
=======

PZGSRFS improves the computed solution to a system of linear
equations and provides error bounds and backward error estimates
for the solution.

Arguments
=========

options (input) superlu_dist_options_t* (global)
        The structure defines the input parameters to control
        how the LU decomposition and triangular solve are performed.

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_Z; Mtype = SLU_GE.

anorm  (input) double
       The norm of the original matrix A, or the scaled A if
       equilibration was done.

LUstruct (input) zLUstruct_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_zdefs.h for the definition of 'zLUstruct_t'.

ScalePermstruct (input) zScalePermstruct_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) doublecomplex* (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) doublecomplex* (local)
       On entry, the solution matrix Y, as computed by PZGSTRS, 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) zSOLVEstruct_t* (global)
       Contains the information for the communication during the
       solution phase.

berr   (output) double*, 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.
Here is the call graph for this function: