SuperLU_DIST  4.0
superlu_dist on CPU and GPU clusters
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Classes | Macros | Functions
superlu_zdefs.h File Reference

Distributed SuperLU data types and function prototypes. More...

#include "superlu_defs.h"
#include "dcomplex.h"

Go to the source code of this file.

Classes

struct  Ucb_indptr_t
 
struct  LocalLU_t
 
struct  LUstruct_t
 
struct  pzgsmv_comm_t
 
struct  pxgstrs_comm_t
 
struct  SOLVEstruct_t
 

Macros

#define MAX_LOOKAHEADS   50
 

Functions

void zCreate_CompCol_Matrix_dist (SuperMatrix *, int_t, int_t, int_t, doublecomplex *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
 
void zCreate_CompRowLoc_Matrix_dist (SuperMatrix *, int_t, int_t, int_t, int_t, int_t, doublecomplex *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
 
void zCompRow_to_CompCol_dist (int_t, int_t, int_t, doublecomplex *, int_t *, int_t *, doublecomplex **, int_t **, int_t **)
 Convert a row compressed storage into a column compressed storage. More...
 
int pzCompRow_loc_to_CompCol_global (int_t, SuperMatrix *, gridinfo_t *, SuperMatrix *)
 Gather A from the distributed compressed row format to global A in compressed column format. More...
 
void zCopy_CompCol_Matrix_dist (SuperMatrix *, SuperMatrix *)
 Copy matrix A into matrix B. More...
 
void zCreate_Dense_Matrix_dist (SuperMatrix *, int_t, int_t, doublecomplex *, int_t, Stype_t, Dtype_t, Mtype_t)
 
void zCreate_SuperNode_Matrix_dist (SuperMatrix *, int_t, int_t, int_t, doublecomplex *, int_t *, int_t *, int_t *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
 
void zCopy_Dense_Matrix_dist (int_t, int_t, doublecomplex *, int_t, doublecomplex *, int_t)
 
void zallocateA_dist (int_t, int_t, doublecomplex **, int_t **, int_t **)
 
void zGenXtrue_dist (int_t, int_t, doublecomplex *, int_t)
 
void zFillRHS_dist (char *, int_t, doublecomplex *, int_t, SuperMatrix *, doublecomplex *, int_t)
 Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's. More...
 
int zcreate_matrix (SuperMatrix *, int, doublecomplex **, int *, doublecomplex **, int *, FILE *, gridinfo_t *)
 
int zcreate_matrix_rb (SuperMatrix *, int, doublecomplex **, int *, doublecomplex **, int *, FILE *, gridinfo_t *)
 
int zcreate_matrix_dat (SuperMatrix *, int, doublecomplex **, int *, doublecomplex **, int *, FILE *, gridinfo_t *)
 
void zgsequ_dist (SuperMatrix *, double *, double *, double *, double *, double *, int_t *)
 
double zlangs_dist (char *, SuperMatrix *)
 
void zlaqgs_dist (SuperMatrix *, double *, double *, double, double, double, char *)
 
void pzgsequ (SuperMatrix *, double *, double *, double *, double *, double *, int_t *, gridinfo_t *)
 
double pzlangs (char *, SuperMatrix *, gridinfo_t *)
 
void pzlaqgs (SuperMatrix *, double *, double *, double, double, double, char *)
 
int pzPermute_Dense_Matrix (int_t, int_t, int_t[], int_t[], doublecomplex[], int, doublecomplex[], int, int, gridinfo_t *)
 Permute the distributed dense matrix: B <= perm(X). perm[i] = j means the i-th row of X is in the j-th row of B. More...
 
int sp_ztrsv_dist (char *, char *, char *, SuperMatrix *, SuperMatrix *, doublecomplex *, int *)
 
int sp_zgemv_dist (char *, doublecomplex, SuperMatrix *, doublecomplex *, int, doublecomplex, doublecomplex *, int)
 
int sp_zgemm_dist (char *, int, doublecomplex, SuperMatrix *, doublecomplex *, int, doublecomplex, doublecomplex *, int)
 
float zdistribute (fact_t, int_t, SuperMatrix *, Glu_freeable_t *, LUstruct_t *, gridinfo_t *)
 
void pzgssvx_ABglobal (superlu_options_t *, SuperMatrix *, ScalePermstruct_t *, doublecomplex *, int, int, gridinfo_t *, LUstruct_t *, double *, SuperLUStat_t *, int *)
 
float pzdistribute (fact_t, int_t, SuperMatrix *, ScalePermstruct_t *, Glu_freeable_t *, LUstruct_t *, gridinfo_t *)
 
void pzgssvx (superlu_options_t *, SuperMatrix *, ScalePermstruct_t *, doublecomplex *, int, int, gridinfo_t *, LUstruct_t *, SOLVEstruct_t *, double *, SuperLUStat_t *, int *)
 
int zSolveInit (superlu_options_t *, SuperMatrix *, int_t[], int_t[], int_t, LUstruct_t *, gridinfo_t *, SOLVEstruct_t *)
 Initialize the data structure for the solution phase. More...
 
int_t pxgstrs_init (int_t, int_t, int_t, int_t, int_t[], int_t[], gridinfo_t *grid, Glu_persist_t *, SOLVEstruct_t *)
 
void pxgstrs_finalize (pxgstrs_comm_t *)
 
void zSolveFinalize (superlu_options_t *, SOLVEstruct_t *)
 Release the resources used for the solution phase. More...
 
void zldperm_dist (int_t, int_t, int_t, int_t[], int_t[], doublecomplex[], int_t *, double[], double[])
 
int static_schedule (superlu_options_t *, int, int, LUstruct_t *, gridinfo_t *, SuperLUStat_t *, int_t *, int_t *, int *)
 
int_t pzgstrf (superlu_options_t *, int, int, double, LUstruct_t *, gridinfo_t *, SuperLUStat_t *, int *)
 
void pzgstrs_Bglobal (int_t, LUstruct_t *, gridinfo_t *, doublecomplex *, int_t, int, SuperLUStat_t *, int *)
 
void pzgstrs (int_t, LUstruct_t *, ScalePermstruct_t *, gridinfo_t *, doublecomplex *, int_t, int_t, int_t, int, SOLVEstruct_t *, SuperLUStat_t *, int *)
 
void zlsum_fmod (doublecomplex *, doublecomplex *, doublecomplex *, doublecomplex *, int, int, int_t, int_t *, int_t, int_t, int_t, int_t *, gridinfo_t *, LocalLU_t *, MPI_Request[], SuperLUStat_t *)
 
void zlsum_bmod (doublecomplex *, doublecomplex *, doublecomplex *, int, int_t, int_t *, int_t *, Ucb_indptr_t **, int_t **, int_t *, gridinfo_t *, LocalLU_t *, MPI_Request[], SuperLUStat_t *)
 
void pzgsrfs (int_t, SuperMatrix *, double, LUstruct_t *, ScalePermstruct_t *, gridinfo_t *, doublecomplex[], int_t, doublecomplex[], int_t, int, SOLVEstruct_t *, double *, SuperLUStat_t *, int *)
 
void pzgsrfs_ABXglobal (int_t, SuperMatrix *, double, LUstruct_t *, gridinfo_t *, doublecomplex *, int_t, doublecomplex *, int_t, int, double *, SuperLUStat_t *, int *)
 
int pzgsmv_AXglobal_setup (SuperMatrix *, Glu_persist_t *, gridinfo_t *, int_t *, int_t *[], doublecomplex *[], int_t *[], int_t[])
 
int pzgsmv_AXglobal (int_t, int_t[], doublecomplex[], int_t[], doublecomplex[], doublecomplex[])
 
int pzgsmv_AXglobal_abs (int_t, int_t[], doublecomplex[], int_t[], doublecomplex[], double[])
 
void pzgsmv_init (SuperMatrix *, int_t *, gridinfo_t *, pzgsmv_comm_t *)
 
void pzgsmv (int_t, SuperMatrix *, gridinfo_t *, pzgsmv_comm_t *, doublecomplex x[], doublecomplex ax[])
 
void pzgsmv_finalize (pzgsmv_comm_t *)
 
doublecomplexdoublecomplexMalloc_dist (int_t)
 
doublecomplexdoublecomplexCalloc_dist (int_t)
 
double * doubleMalloc_dist (int_t)
 
double * doubleCalloc_dist (int_t)
 
void * duser_malloc_dist (int_t, int_t)
 
void duser_free_dist (int_t, int_t)
 
int_t zQuerySpace_dist (int_t, LUstruct_t *, gridinfo_t *, SuperLUStat_t *, mem_usage_t *)
 
void Destroy_LU (int_t, gridinfo_t *, LUstruct_t *)
 Destroy distributed L & U matrices. More...
 
void LUstructInit (const int_t, LUstruct_t *)
 Allocate storage in LUstruct. More...
 
void LUstructFree (LUstruct_t *)
 Deallocate LUstruct. More...
 
void zfill_dist (doublecomplex *, int_t, doublecomplex)
 Fills a doublecomplex precision array with a given value. More...
 
void zinf_norm_error_dist (int_t, int_t, doublecomplex *, int_t, doublecomplex *, int_t, gridinfo_t *)
 Check the inf-norm of the error vector. More...
 
void pzinf_norm_error (int, int_t, int_t, doublecomplex[], int_t, doublecomplex[], int_t, gridinfo_t *)
 Check the inf-norm of the error vector. More...
 
void zreadhb_dist (int, FILE *, int_t *, int_t *, int_t *, doublecomplex **, int_t **, int_t **)
 
void zreadtriple (FILE *, int_t *, int_t *, int_t *, doublecomplex **, int_t **, int_t **)
 
void zreadrb_dist (FILE *, int *, int *, int *, doublecomplex **, int **, int **)
 
float zdist_psymbtonum (fact_t, int_t, SuperMatrix *, ScalePermstruct_t *, Pslu_freeable_t *, LUstruct_t *, gridinfo_t *)
 
void zPrintLblocks (int, int_t, gridinfo_t *, Glu_persist_t *, LocalLU_t *)
 Print the blocks in the factored matrix L. More...
 
void zPrintUblocks (int, int_t, gridinfo_t *, Glu_persist_t *, LocalLU_t *)
 Print the blocks in the factored matrix U. More...
 
void zPrint_CompCol_Matrix_dist (SuperMatrix *)
 
void zPrint_Dense_Matrix_dist (SuperMatrix *)
 
int zPrint_CompRowLoc_Matrix_dist (SuperMatrix *)
 
void PrintDoublecomplex (char *, int_t, doublecomplex *)
 
int file_PrintDoublecomplex (FILE *fp, char *, int_t, doublecomplex *)
 
int zgemm_ (const char *, const char *, const int *, const int *, const int *, const doublecomplex *, const doublecomplex *, const int *, const doublecomplex *, const int *, const doublecomplex *, doublecomplex *, const int *)
 
int ztrsv_ (char *, char *, char *, int *, doublecomplex *, int *, doublecomplex *, int *)
 
int ztrsm_ (char *, char *, char *, char *, int *, int *, doublecomplex *, doublecomplex *, int *, doublecomplex *, int *)
 
int zgemv_ (char *, int *, int *, doublecomplex *, doublecomplex *a, int *, doublecomplex *, int *, doublecomplex *, doublecomplex *, int *)
 
int zgeru_ (int *, int *, doublecomplex *, doublecomplex *, int *, doublecomplex *, int *, doublecomplex *, int *)
 

Detailed Description

Distributed SuperLU data types and function prototypes.

– Distributed SuperLU routine (version 2.5) –
Lawrence Berkeley National Lab, Univ. of California Berkeley.
November 1, 2007

Macro Definition Documentation

#define MAX_LOOKAHEADS   50

Function Documentation

void Destroy_LU ( int_t  ,
gridinfo_t ,
LUstruct_t  
)

Destroy distributed L & U matrices.

double* doubleCalloc_dist ( int_t  )
doublecomplex* doublecomplexCalloc_dist ( int_t  )
doublecomplex* doublecomplexMalloc_dist ( int_t  )
double* doubleMalloc_dist ( int_t  )
void duser_free_dist ( int_t  ,
int_t   
)
void* duser_malloc_dist ( int_t  ,
int_t   
)
int file_PrintDoublecomplex ( FILE *  fp,
char *  ,
int_t  ,
doublecomplex  
)
void LUstructFree ( LUstruct_t )

Deallocate LUstruct.

void LUstructInit ( const int_t  ,
LUstruct_t  
)

Allocate storage in LUstruct.

void PrintDoublecomplex ( char *  ,
int_t  ,
doublecomplex  
)
void pxgstrs_finalize ( pxgstrs_comm_t )
int_t pxgstrs_init ( int_t  n,
int_t  m_loc,
int_t  nrhs,
int_t  fst_row,
int_t  perm_r[],
int_t  perm_c[],
gridinfo_t grid,
Glu_persist_t Glu_persist,
SOLVEstruct_t SOLVEstruct 
)

Purpose

  Set up the communication pattern for the triangular solution.

Arguments

n      (input) int (global)
       The dimension of the linear system.
m_loc  (input) int (local)
       The local row dimension of the distributed input matrix.
nrhs   (input) int (global)
       Number of right-hand sides.
fst_row (input) int (global)
       The row number of matrix B's first row in the global matrix.
perm_r (input) int* (global)
       The row permutation vector.
perm_c (input) int* (global)
       The column permutation vector.
grid   (input) gridinfo_t*
       The 2D process mesh.
int pzCompRow_loc_to_CompCol_global ( int_t  ,
SuperMatrix ,
gridinfo_t ,
SuperMatrix  
)

Gather A from the distributed compressed row format to global A in compressed column format.

float pzdistribute ( fact_t  ,
int_t  ,
SuperMatrix ,
ScalePermstruct_t ,
Glu_freeable_t ,
LUstruct_t ,
gridinfo_t  
)
void pzgsequ ( SuperMatrix A,
double *  r,
double *  c,
double *  rowcnd,
double *  colcnd,
double *  amax,
int_t info,
gridinfo_t grid 
)
    

Purpose

    PZGSEQU computes row and column scalings intended to equilibrate an   
    M-by-N sparse matrix A and reduce its condition number. R returns the row
    scale factors and C the column scale factors, chosen to try to make   
    the largest element in each row and column of the matrix B with   
    elements B(i,j)=R(i)*A(i,j)*C(j) have absolute value 1.
    R(i) and C(j) are restricted to be between SMLNUM = smallest safe   
    number and BIGNUM = largest safe number.  Use of these scaling   
    factors is not guaranteed to reduce the condition number of A but   
    works well in practice.
    See supermatrix.h for the definition of 'SuperMatrix' structure.

Arguments

    A       (input) SuperMatrix*
            The matrix of dimension (A->nrow, A->ncol) whose equilibration
            factors are to be computed. The type of A can be:
            Stype = SLU_NR_loc; Dtype = SLU_Z; Mtype = SLU_GE.
    R       (output) double*, size A->nrow
            If INFO = 0 or INFO > M, R contains the row scale factors   
            for A.
    C       (output) double*, size A->ncol
            If INFO = 0,  C contains the column scale factors for A.
    ROWCND  (output) double*
            If INFO = 0 or INFO > M, ROWCND contains the ratio of the   
            smallest R(i) to the largest R(i).  If ROWCND >= 0.1 and   
            AMAX is neither too large nor too small, it is not worth   
            scaling by R.
    COLCND  (output) double*
            If INFO = 0, COLCND contains the ratio of the smallest   
            C(i) to the largest C(i).  If COLCND >= 0.1, it is not   
            worth scaling by C.
    AMAX    (output) double*
            Absolute value of largest matrix element.  If AMAX is very   
            close to overflow or very close to underflow, the matrix   
            should be scaled.
    INFO    (output) int*
            = 0:  successful exit   
            < 0:  if INFO = -i, the i-th argument had an illegal value   
            > 0:  if INFO = i,  and i is   
                  <= M:  the i-th row of A is exactly zero   
                  >  M:  the (i-M)-th column of A is exactly zero
    GRID    (input) gridinof_t*

The 2D process mesh.

 
void pzgsmv ( int_t  ,
SuperMatrix ,
gridinfo_t ,
pzgsmv_comm_t ,
doublecomplex  x[],
doublecomplex  ax[] 
)
int pzgsmv_AXglobal ( int_t  m,
int_t  update[],
doublecomplex  val[],
int_t  bindx[],
doublecomplex  X[],
doublecomplex  ax[] 
)
Performs sparse matrix-vector multiplication.

  • val/bindx stores the distributed MSR matrix A
  • X is global
  • ax product is distributed the same way as A
int pzgsmv_AXglobal_abs ( int_t  ,
int_t  [],
doublecomplex  [],
int_t  [],
doublecomplex  [],
double  [] 
)
int pzgsmv_AXglobal_setup ( SuperMatrix ,
Glu_persist_t ,
gridinfo_t ,
int_t ,
int_t [],
doublecomplex [],
int_t [],
int_t  [] 
)
void pzgsmv_finalize ( pzgsmv_comm_t )
void pzgsmv_init ( SuperMatrix ,
int_t ,
gridinfo_t ,
pzgsmv_comm_t  
)
void pzgsrfs ( int_t  ,
SuperMatrix ,
double  ,
LUstruct_t ,
ScalePermstruct_t ,
gridinfo_t ,
doublecomplex  [],
int_t  ,
doublecomplex  [],
int_t  ,
int  ,
SOLVEstruct_t ,
double *  ,
SuperLUStat_t ,
int *   
)
void pzgsrfs_ABXglobal ( int_t  n,
SuperMatrix A,
double  anorm,
LUstruct_t LUstruct,
gridinfo_t grid,
doublecomplex B,
int_t  ldb,
doublecomplex X,
int_t  ldx,
int  nrhs,
double *  berr,
SuperLUStat_t stat,
int *  info 
)

Purpose

pzgsrfs_ABXglobal 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 the form Pc*Pr*A*Pc', where Pr and Pc
       are permutation matrices. The type of A can be:
       Stype = SLU_NCP; Dtype = SLU_Z; Mtype = SLU_GE.
       NOTE: Currently, A must reside in all processes when calling
             this routine.
anorm  (input) double
       The norm of the original matrix A, or the scaled A if
       equilibration was done.
LUstruct (input) LUstruct_t*
       The distributed data structures storing L and U factors.
       The L and U factors are obtained from pzgstrf for
       the possibly scaled and permuted matrix A.
       See superlu_ddefs.h for the definition of 'LUstruct_t'.
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_ddefs.h for the definition of 'gridinfo_t'.
B      (input) doublecomplex* (global)
       The N-by-NRHS right-hand side matrix of the possibly equilibrated
       and row permuted system.
       NOTE: Currently, B must reside on all processes when calling
             this routine.
ldb    (input) int (global)
       Leading dimension of matrix B.
X      (input/output) doublecomplex* (global)
       On entry, the solution matrix X, as computed by PZGSTRS.
       On exit, the improved solution matrix X.
       If DiagScale = COL or BOTH, X should be premultiplied by diag(C)
       in order to obtain the solution to the original system.
       NOTE: Currently, X must reside on all processes when calling
             this routine.
ldx    (input) int (global)
       Leading dimension of matrix X.
nrhs   (input) int
       Number of right-hand sides.
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.   
void pzgssvx ( superlu_options_t ,
SuperMatrix ,
ScalePermstruct_t ,
doublecomplex ,
int  ,
int  ,
gridinfo_t ,
LUstruct_t ,
SOLVEstruct_t ,
double *  ,
SuperLUStat_t ,
int *   
)
void pzgssvx_ABglobal ( superlu_options_t ,
SuperMatrix ,
ScalePermstruct_t ,
doublecomplex ,
int  ,
int  ,
gridinfo_t ,
LUstruct_t ,
double *  ,
SuperLUStat_t ,
int *   
)
int_t pzgstrf ( superlu_options_t options,
int  m,
int  n,
double  anorm,
LUstruct_t LUstruct,
gridinfo_t grid,
SuperLUStat_t stat,
int *  info 
)

Purpose

 PZGSTRF performs the LU factorization in parallel.

Arguments

options (input) superlu_options_t*
        The structure defines the input parameters to control
        how the LU decomposition will be performed.
        The following field should be defined:
        o ReplaceTinyPivot (yes_no_t)
          Specifies whether to replace the tiny diagonals by
          sqrt(epsilon)*norm(A) during LU factorization.
m      (input) int
       Number of rows in the matrix.
n      (input) int
       Number of columns in the matrix.
anorm  (input) double
       The norm of the original matrix A, or the scaled A if
       equilibration was done.
LUstruct (input/output) LUstruct_t*
        The data structures to store the distributed L and U factors.
        The following fields should be defined:
        o Glu_persist (input) Glu_persist_t*
          Global data structure (xsup, supno) replicated on all processes,
          describing the supernode partition in the factored matrices
          L and U:
        xsup[s] is the leading column of the s-th supernode,
            supno[i] is the supernode number to which column i belongs.
        o Llu (input/output) LocalLU_t*
          The distributed data structures to store L and U factors.
          See superlu_zdefs.h for the definition of 'LocalLU_t'.
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_zdefs.h for the definition of 'gridinfo_t'.
stat   (output) SuperLUStat_t*
       Record the statistics on runtime and floating-point operation count.
       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
       > 0: if info = i, U(i,i) is exactly zero. The factorization has
            been completed, but the factor U is exactly singular,
            and division by zero will occur if it is used to solve a
            system of equations.
void pzgstrs ( int_t  n,
LUstruct_t LUstruct,
ScalePermstruct_t ScalePermstruct,
gridinfo_t grid,
doublecomplex B,
int_t  m_loc,
int_t  fst_row,
int_t  ldb,
int  nrhs,
SOLVEstruct_t SOLVEstruct,
SuperLUStat_t stat,
int *  info 
)
Purpose
=======
PZGSTRS solves a system of distributed linear equations
A*X = B with a general N-by-N matrix A using the LU factorization
computed by PZGSTRF.
If the equilibration, and row and column permutations were performed,
the LU factorization was performed for A1 where
    A1 = Pc*Pr*diag(R)*A*diag(C)*Pc^T = L*U
and the linear system solved is
    A1 * Y = Pc*Pr*B1, where B was overwritten by B1 = diag(R)*B, and
the permutation to B1 by Pc*Pr is applied internally in this routine.
Arguments
=========
n      (input) int (global)
       The order of the system of linear equations.
LUstruct (input) LUstruct_t*
       The distributed data structures storing L and U factors.
       The L and U factors are obtained from PZGSTRF for
       the possibly scaled and permuted matrix A.
       See superlu_zdefs.h for the definition of 'LUstruct_t'.
       A may be scaled and permuted into A1, so that
       A1 = Pc*Pr*diag(R)*A*diag(C)*Pc^T = L*U
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/output) doublecomplex*
       On entry, the distributed right-hand side matrix of the possibly
       equilibrated system. That is, B may be overwritten by diag(R)*B.
       On exit, the distributed solution matrix Y of the possibly
       equilibrated system if info = 0, where Y = Pc*diag(C)^(-1)*X,
       and X is the solution of the original system.
m_loc  (input) int (local)
       The local row dimension of matrix B.
fst_row (input) int (global)
       The row number of B's first row in the global matrix.
ldb    (input) int (local)
       The leading dimension of matrix B.
nrhs   (input) int (global)
       Number of right-hand sides.
SOLVEstruct (input) SOLVEstruct_t* (global)
       Contains the information for the communication during the
       solution phase.
stat   (output) SuperLUStat_t*
       Record the statistics about the triangular solves.
       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
void pzgstrs_Bglobal ( int_t  n,
LUstruct_t LUstruct,
gridinfo_t grid,
doublecomplex B,
int_t  ldb,
int  nrhs,
SuperLUStat_t stat,
int *  info 
)
Purpose
=======
pzgstrs_Bglobal solves a system of distributed linear equations
A*X = B with a general N-by-N matrix A using the LU factorization
computed by pzgstrf.
Arguments
=========
n      (input) int (global)
       The order of the system of linear equations.
LUstruct (input) LUstruct_t*
       The distributed data structures storing L and U factors.
       The L and U factors are obtained from pzgstrf for
       the possibly scaled and permuted matrix A.
       See superlu_ddefs.h for the definition of 'LUstruct_t'.
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_ddefs.h for the definition of 'gridinfo_t'.
B      (input/output) doublecomplex*
       On entry, the right-hand side matrix of the possibly equilibrated
       and row permuted system.
       On exit, the solution matrix of the possibly equilibrated
       and row permuted system if info = 0;
       NOTE: Currently, the N-by-NRHS  matrix B must reside on all 
             processes when calling this routine.
ldb    (input) int (global)
       Leading dimension of matrix B.
nrhs   (input) int (global)
       Number of right-hand sides.
stat   (output) SuperLUStat_t*
       Record the statistics about the triangular solves.
       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

Purpose

pzgstrs_Bglobal solves a system of distributed linear equations
A*X = B with a general N-by-N matrix A using the LU factorization
computed by pzgstrf.

Arguments

n      (input) int (global)
       The order of the system of linear equations.
LUstruct (input) LUstruct_t*
       The distributed data structures storing L and U factors.
       The L and U factors are obtained from pzgstrf for
       the possibly scaled and permuted matrix A.
       See superlu_ddefs.h for the definition of 'LUstruct_t'.
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_ddefs.h for the definition of 'gridinfo_t'.
B      (input/output) doublecomplex*
       On entry, the right-hand side matrix of the possibly equilibrated
       and row permuted system.
       On exit, the solution matrix of the possibly equilibrated
       and row permuted system if info = 0;
       NOTE: Currently, the N-by-NRHS  matrix B must reside on all 
             processes when calling this routine.
ldb    (input) int (global)
       Leading dimension of matrix B.
nrhs   (input) int (global)
       Number of right-hand sides.
stat   (output) SuperLUStat_t*
       Record the statistics about the triangular solves.
       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
void pzinf_norm_error ( int  ,
int_t  ,
int_t  ,
doublecomplex  [],
int_t  ,
doublecomplex  [],
int_t  ,
gridinfo_t  
)

Check the inf-norm of the error vector.

double pzlangs ( char *  norm,
SuperMatrix A,
gridinfo_t grid 
)
 

Purpose

    PZLANGS returns the value of the one norm, or the Frobenius norm, or 
    the infinity norm, or the element of largest absolute value of a 
    real matrix A.

Description

    PZLANGE returns the value
       PZLANGE = ( max(abs(A(i,j))), NORM = 'M' or 'm'   
                 (   
                 ( norm1(A),         NORM = '1', 'O' or 'o'   
                 (   
                 ( normI(A),         NORM = 'I' or 'i'   
                 (   
                 ( normF(A),         NORM = 'F', 'f', 'E' or 'e'
    where  norm1  denotes the  one norm of a matrix (maximum column sum), 
    normI  denotes the  infinity norm  of a matrix  (maximum row sum) and 
    normF  denotes the  Frobenius norm of a matrix (square root of sum of 
    squares).  Note that  max(abs(A(i,j)))  is not a  matrix norm.

Arguments

    NORM    (input) CHARACTER*1   
            Specifies the value to be returned in DLANGE as described above.   
    A       (input) SuperMatrix*
            The M by N sparse matrix A. 
    GRID    (input) gridinof_t*

The 2D process mesh.

 
void pzlaqgs ( SuperMatrix A,
double *  r,
double *  c,
double  rowcnd,
double  colcnd,
double  amax,
char *  equed 
)

Purpose

    PZLAQGS equilibrates a general sparse M by N matrix A using the row
    and column scaling factors in the vectors R and C.
    See supermatrix.h for the definition of 'SuperMatrix' structure.

Arguments

    A       (input/output) SuperMatrix*
            On exit, the equilibrated matrix.  See EQUED for the form of 
            the equilibrated matrix. The type of A can be:
        Stype = SLU_NR_loc; Dtype = SLU_Z; Mtype = SLU_GE.
    R       (input) double*, dimension (A->nrow)
            The row scale factors for A.
    C       (input) double*, dimension (A->ncol)
            The column scale factors for A.
    ROWCND  (input) double
            Ratio of the smallest R(i) to the largest R(i).
    COLCND  (input) double
            Ratio of the smallest C(i) to the largest C(i).
    AMAX    (input) double
            Absolute value of largest matrix entry.
    EQUED   (output) char*
            Specifies the form of equilibration that was done.   
            = 'N':  No equilibration   
            = 'R':  Row equilibration, i.e., A has been premultiplied by  
                    diag(R).   
            = 'C':  Column equilibration, i.e., A has been postmultiplied  
                    by diag(C).   
            = 'B':  Both row and column equilibration, i.e., A has been
                    replaced by diag(R) * A * diag(C).

Internal Parameters

    THRESH is a threshold value used to decide if row or column scaling   
    should be done based on the ratio of the row or column scaling   
    factors.  If ROWCND < THRESH, row scaling is done, and if   
    COLCND < THRESH, column scaling is done.
    LARGE and SMALL are threshold values used to decide if row scaling   
    should be done based on the absolute size of the largest matrix   
    element.  If AMAX > LARGE or AMAX < SMALL, row scaling is done.   


int pzPermute_Dense_Matrix ( int_t  ,
int_t  ,
int_t  [],
int_t  [],
doublecomplex  [],
int  ,
doublecomplex  [],
int  ,
int  ,
gridinfo_t  
)

Permute the distributed dense matrix: B <= perm(X). perm[i] = j means the i-th row of X is in the j-th row of B.

int sp_zgemm_dist ( char *  transa,
int  n,
doublecomplex  alpha,
SuperMatrix A,
doublecomplex b,
int  ldb,
doublecomplex  beta,
doublecomplex c,
int  ldc 
)
  Purpose   
    =======
    sp_z performs one of the matrix-matrix operations
       C := alpha*op( A )*op( B ) + beta*C,
    where  op( X ) is one of
       op( X ) = X   or   op( X ) = X'   or   op( X ) = conjg( X' ),
    alpha and beta are scalars, and A, B and C are matrices, with op( A ) 
    an m by k matrix,  op( B )  a  k by n matrix and  C an m by n matrix.
    Parameters   
    ==========
    TRANSA - (input) char*
             On entry, TRANSA specifies the form of op( A ) to be used in 
             the matrix multiplication as follows:   
                TRANSA = 'N' or 'n',  op( A ) = A.   
                TRANSA = 'T' or 't',  op( A ) = A'.   
                TRANSA = 'C' or 'c',  op( A ) = conjg( A' ).   
             Unchanged on exit.
    TRANSB - (input) char*
             On entry, TRANSB specifies the form of op( B ) to be used in 
             the matrix multiplication as follows:   
                TRANSB = 'N' or 'n',  op( B ) = B.   
                TRANSB = 'T' or 't',  op( B ) = B'.   
                TRANSB = 'C' or 'c',  op( B ) = conjg( B' ).   
             Unchanged on exit.
    M      - (input) int   
             On entry,  M  specifies  the number of rows of the matrix 
         op( A ) and of the matrix C.  M must be at least zero. 
         Unchanged on exit.
    N      - (input) int
             On entry,  N specifies the number of columns of the matrix 
         op( B ) and the number of columns of the matrix C. N must be 
         at least zero.
         Unchanged on exit.
    K      - (input) int
             On entry, K specifies the number of columns of the matrix 
         op( A ) and the number of rows of the matrix op( B ). K must 
         be at least  zero.   
             Unchanged on exit.
    ALPHA  - (input) doublecomplex
             On entry, ALPHA specifies the scalar alpha.
    A      - (input) SuperMatrix*
             Matrix A with a sparse format, of dimension (A->nrow, A->ncol).
             Currently, the type of A can be:
                 Stype = NC or NCP; Dtype = Z; Mtype = GE. 
             In the future, more general A can be handled.
    B      - DOUBLE COMPLEX PRECISION array of DIMENSION ( LDB, kb ), where kb is 
             n when TRANSB = 'N' or 'n',  and is  k otherwise.   
             Before entry with  TRANSB = 'N' or 'n',  the leading k by n 
             part of the array B must contain the matrix B, otherwise 
             the leading n by k part of the array B must contain the 
             matrix B.   
             Unchanged on exit.
    LDB    - (input) int
             On entry, LDB specifies the first dimension of B as declared 
             in the calling (sub) program. LDB must be at least max( 1, n ).  
             Unchanged on exit.
    BETA   - (input) doublecomplex
             On entry, BETA specifies the scalar beta. When BETA is   
             supplied as zero then C need not be set on input.
    C      - DOUBLE COMPLEX PRECISION array of DIMENSION ( LDC, n ).   
             Before entry, the leading m by n part of the array C must 
             contain the matrix C,  except when beta is zero, in which 
             case C need not be set on entry.   
             On exit, the array C is overwritten by the m by n matrix 
         ( alpha*op( A )*B + beta*C ).
    LDC    - (input) int
             On entry, LDC specifies the first dimension of C as declared 
             in the calling (sub)program. LDC must be at least max(1,m).   
             Unchanged on exit.
    ==== Sparse Level 3 Blas routine.  
int sp_zgemv_dist ( char *  trans,
doublecomplex  alpha,
SuperMatrix A,
doublecomplex x,
int  incx,
doublecomplex  beta,
doublecomplex y,
int  incy 
)

Purpose

    sp_zgemv()  performs one of the matrix-vector operations   
       y := alpha*A*x + beta*y,   or   y := alpha*A'*x + beta*y,   
    where alpha and beta are scalars, x and y are vectors and A is a
    sparse A->nrow by A->ncol matrix.

Parameters

    TRANS  - (input) char*
             On entry, TRANS specifies the operation to be performed as   
             follows:   
                TRANS = 'N' or 'n'   y := alpha*A*x + beta*y.   
                TRANS = 'T' or 't'   y := alpha*A'*x + beta*y.   
                TRANS = 'C' or 'c'   y := alpha*A'*x + beta*y.
    ALPHA  - (input) doublecomplex
             On entry, ALPHA specifies the scalar alpha.
    A      - (input) SuperMatrix*
             Before entry, the leading m by n part of the array A must   
             contain the matrix of coefficients.
    X      - (input) doublecomplex*, array of DIMENSION at least   
             ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'   
             and at least   
             ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.   
             Before entry, the incremented array X must contain the   
             vector x.
    INCX   - (input) int
             On entry, INCX specifies the increment for the elements of   
             X. INCX must not be zero.
    BETA   - (input) doublecomplex
             On entry, BETA specifies the scalar beta. When BETA is   
             supplied as zero then Y need not be set on input.
    Y      - (output) doublecomplex*,  array of DIMENSION at least   
             ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'   
             and at least   
             ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.   
             Before entry with BETA non-zero, the incremented array Y   
             must contain the vector y. On exit, Y is overwritten by the 
             updated vector y.
    INCY   - (input) int
             On entry, INCY specifies the increment for the elements of   
             Y. INCY must not be zero.
    ==== Sparse Level 2 Blas routine.   
int sp_ztrsv_dist ( char *  uplo,
char *  trans,
char *  diag,
SuperMatrix L,
SuperMatrix U,
doublecomplex x,
int *  info 
)

Purpose

  sp_ztrsv() solves one of the systems of equations   
      A*x = b,   or   A'*x = b,
  where b and x are n element vectors and A is a sparse unit , or   
  non-unit, upper or lower triangular matrix.   
  No test for singularity or near-singularity is included in this   
  routine. Such tests must be performed before calling this routine.

Parameters

  uplo   - (input) char*
           On entry, uplo specifies whether the matrix is an upper or   
            lower triangular matrix as follows:   
               uplo = 'U' or 'u'   A is an upper triangular matrix.   
               uplo = 'L' or 'l'   A is a lower triangular matrix.
  trans  - (input) char*
            On entry, trans specifies the equations to be solved as   
            follows:   
               trans = 'N' or 'n'   A*x = b.   
               trans = 'T' or 't'   A'*x = b.   
               trans = 'C' or 'c'   A'*x = b.
  diag   - (input) char*
            On entry, diag specifies whether or not A is unit   
            triangular as follows:   
               diag = 'U' or 'u'   A is assumed to be unit triangular.   
               diag = 'N' or 'n'   A is not assumed to be unit   
                                   triangular.
  L       - (input) SuperMatrix*
        The factor L from the factorization Pr*A*Pc=L*U. Use
            compressed row subscripts storage for supernodes,
            i.e., L has types: Stype = SC, Dtype = Z, Mtype = TRLU.
  U       - (input) SuperMatrix*
         The factor U from the factorization Pr*A*Pc=L*U.
         U has types: Stype = NC, Dtype = Z, Mtype = TRU.
  x       - (input/output) doublecomplex*
            Before entry, the incremented array X must contain the n   
            element right-hand side vector b. On exit, X is overwritten 
            with the solution vector x.
  info    - (output) int*
            If *info = -i, the i-th argument had an illegal value.
int static_schedule ( superlu_options_t ,
int  ,
int  ,
LUstruct_t ,
gridinfo_t ,
SuperLUStat_t ,
int_t ,
int_t ,
int *   
)
void zallocateA_dist ( int_t  ,
int_t  ,
doublecomplex **  ,
int_t **  ,
int_t **   
)
void zCompRow_to_CompCol_dist ( int_t  ,
int_t  ,
int_t  ,
doublecomplex ,
int_t ,
int_t ,
doublecomplex **  ,
int_t **  ,
int_t **   
)

Convert a row compressed storage into a column compressed storage.

void zCopy_CompCol_Matrix_dist ( SuperMatrix ,
SuperMatrix  
)

Copy matrix A into matrix B.

void zCopy_Dense_Matrix_dist ( int_t  ,
int_t  ,
doublecomplex ,
int_t  ,
doublecomplex ,
int_t   
)
 Purpose
 =======
 Copies a two-dimensional matrix X to another matrix Y.
void zCreate_CompCol_Matrix_dist ( SuperMatrix ,
int_t  ,
int_t  ,
int_t  ,
doublecomplex ,
int_t ,
int_t ,
Stype_t  ,
Dtype_t  ,
Mtype_t   
)
void zCreate_CompRowLoc_Matrix_dist ( SuperMatrix ,
int_t  ,
int_t  ,
int_t  ,
int_t  ,
int_t  ,
doublecomplex ,
int_t ,
int_t ,
Stype_t  ,
Dtype_t  ,
Mtype_t   
)
void zCreate_Dense_Matrix_dist ( SuperMatrix ,
int_t  ,
int_t  ,
doublecomplex ,
int_t  ,
Stype_t  ,
Dtype_t  ,
Mtype_t   
)
int zcreate_matrix ( SuperMatrix ,
int  ,
doublecomplex **  ,
int *  ,
doublecomplex **  ,
int *  ,
FILE *  ,
gridinfo_t  
)
int zcreate_matrix_dat ( SuperMatrix ,
int  ,
doublecomplex **  ,
int *  ,
doublecomplex **  ,
int *  ,
FILE *  ,
gridinfo_t  
)
int zcreate_matrix_rb ( SuperMatrix ,
int  ,
doublecomplex **  ,
int *  ,
doublecomplex **  ,
int *  ,
FILE *  ,
gridinfo_t  
)
void zCreate_SuperNode_Matrix_dist ( SuperMatrix ,
int_t  ,
int_t  ,
int_t  ,
doublecomplex ,
int_t ,
int_t ,
int_t ,
int_t ,
int_t ,
Stype_t  ,
Dtype_t  ,
Mtype_t   
)
float zdist_psymbtonum ( fact_t  fact,
int_t  n,
SuperMatrix A,
ScalePermstruct_t ScalePermstruct,
Pslu_freeable_t Pslu_freeable,
LUstruct_t LUstruct,
gridinfo_t grid 
)

Purpose

  Distribute the input matrix onto the 2D process mesh.

Arguments

fact (input) fact_t
       Specifies whether or not the L and U structures will be re-used.
       = SamePattern_SameRowPerm: L and U structures are input, and
                                  unchanged on exit.
         This routine should not be called for this case, an error
         is generated.  Instead, pddistribute routine should be called.
       = DOFACT or SamePattern: L and U structures are computed and output.
n      (Input) int
       Dimension of the matrix.
A      (Input) SuperMatrix*
   The distributed input matrix A of dimension (A->nrow, A->ncol).
       A may be overwritten by diag(R)*A*diag(C)*Pc^T.
       The type of A can be: Stype = NR; Dtype = SLU_D; Mtype = GE.
ScalePermstruct (Input) ScalePermstruct_t*
       The data structure to store the scaling and permutation vectors
       describing the transformations performed to the original matrix A.
Glu_freeable (Input) *Glu_freeable_t
       The global structure describing the graph of L and U.
LUstruct (Input) LUstruct_t*
       Data structures for L and U factors.
grid   (Input) gridinfo_t*
       The 2D process mesh.

Return value

  < 0, number of bytes allocated on return from the dist_symbLU

0, number of bytes allocated for performing the distribution

of the data, when out of memory. (an approximation).

float zdistribute ( fact_t  fact,
int_t  n,
SuperMatrix A,
Glu_freeable_t Glu_freeable,
LUstruct_t LUstruct,
gridinfo_t grid 
)

Purpose

  Distribute the matrix onto the 2D process mesh.

Arguments

fact (input) fact_t
       Specifies whether or not the L and U structures will be re-used.
       = SamePattern_SameRowPerm: L and U structures are input, and
                                  unchanged on exit.
       = DOFACT or SamePattern: L and U structures are computed and output.
n      (input) int
       Dimension of the matrix.
A      (input) SuperMatrix*
   The original matrix A, permuted by columns, of dimension
       (A->nrow, A->ncol). The type of A can be:
       Stype = SLU_NCP; Dtype = SLU_Z; Mtype = SLU_GE.
LUstruct (input) LUstruct_t*
       Data structures for L and U factors.
grid   (input) gridinfo_t*
       The 2D process mesh.

Return value

0, working storage required (in bytes).

– Distributed SuperLU routine (version 1.0) –
Lawrence Berkeley National Lab, Univ. of California Berkeley.
September 1, 1999

Purpose

  Distribute the matrix onto the 2D process mesh.

Arguments

fact (input) fact_t
       Specifies whether or not the L and U structures will be re-used.
       = SamePattern_SameRowPerm: L and U structures are input, and
                                  unchanged on exit.
       = DOFACT or SamePattern: L and U structures are computed and output.
n      (input) int
       Dimension of the matrix.
A      (input) SuperMatrix*
   The original matrix A, permuted by columns, of dimension
       (A->nrow, A->ncol). The type of A can be:
       Stype = NCP; Dtype = Z; Mtype = GE.
LUstruct (input) LUstruct_t*
       Data structures for L and U factors.
grid   (input) gridinfo_t*
       The 2D process mesh.
void zfill_dist ( doublecomplex ,
int_t  ,
doublecomplex   
)

Fills a doublecomplex precision array with a given value.

void zFillRHS_dist ( char *  ,
int_t  ,
doublecomplex ,
int_t  ,
SuperMatrix ,
doublecomplex ,
int_t   
)

Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's.

int zgemm_ ( const char *  ,
const char *  ,
const int *  ,
const int *  ,
const int *  ,
const doublecomplex ,
const doublecomplex ,
const int *  ,
const doublecomplex ,
const int *  ,
const doublecomplex ,
doublecomplex ,
const int *   
)
int zgemv_ ( char *  ,
int *  ,
int *  ,
doublecomplex ,
doublecomplex a,
int *  ,
doublecomplex ,
int *  ,
doublecomplex ,
doublecomplex ,
int *   
)
void zGenXtrue_dist ( int_t  ,
int_t  ,
doublecomplex ,
int_t   
)
int zgeru_ ( int *  ,
int *  ,
doublecomplex ,
doublecomplex ,
int *  ,
doublecomplex ,
int *  ,
doublecomplex ,
int *   
)
void zgsequ_dist ( SuperMatrix A,
double *  r,
double *  c,
double *  rowcnd,
double *  colcnd,
double *  amax,
int_t info 
)
   

Purpose

    ZGSEQU_DIST computes row and column scalings intended to equilibrate an   
    M-by-N sparse matrix A and reduce its condition number. R returns the row
    scale factors and C the column scale factors, chosen to try to make   
    the largest element in each row and column of the matrix B with   
    elements B(i,j)=R(i)*A(i,j)*C(j) have absolute value 1.
    R(i) and C(j) are restricted to be between SMLNUM = smallest safe   
    number and BIGNUM = largest safe number.  Use of these scaling   
    factors is not guaranteed to reduce the condition number of A but   
    works well in practice.
    See supermatrix.h for the definition of 'SuperMatrix' structure.

Arguments

    A       (input) SuperMatrix*
            The matrix of dimension (A->nrow, A->ncol) whose equilibration
            factors are to be computed. The type of A can be:
            Stype = SLU_NC; Dtype = SLU_Z; Mtype = SLU_GE.
    R       (output) double*, size A->nrow
            If INFO = 0 or INFO > M, R contains the row scale factors   
            for A.
    C       (output) double*, size A->ncol
            If INFO = 0,  C contains the column scale factors for A.
    ROWCND  (output) double*
            If INFO = 0 or INFO > M, ROWCND contains the ratio of the   
            smallest R(i) to the largest R(i).  If ROWCND >= 0.1 and   
            AMAX is neither too large nor too small, it is not worth   
            scaling by R.
    COLCND  (output) double*
            If INFO = 0, COLCND contains the ratio of the smallest   
            C(i) to the largest C(i).  If COLCND >= 0.1, it is not   
            worth scaling by C.
    AMAX    (output) double*
            Absolute value of largest matrix element.  If AMAX is very   
            close to overflow or very close to underflow, the matrix   
            should be scaled.
    INFO    (output) int*
            = 0:  successful exit   
            < 0:  if INFO = -i, the i-th argument had an illegal value   
            > 0:  if INFO = i,  and i is   
                  <= M:  the i-th row of A is exactly zero   
                  >  M:  the (i-M)-th column of A is exactly zero   


void zinf_norm_error_dist ( int_t  ,
int_t  ,
doublecomplex ,
int_t  ,
doublecomplex ,
int_t  ,
gridinfo_t  
)

Check the inf-norm of the error vector.

double zlangs_dist ( char *  norm,
SuperMatrix A 
)
 

Purpose

    ZLANGS_DIST returns the value of the one norm, or the Frobenius norm, or 
    the infinity norm, or the element of largest absolute value of a 
    real matrix A.

Description

    ZLANGE returns the value
       ZLANGE = ( max(abs(A(i,j))), NORM = 'M' or 'm'   
                (   
                ( norm1(A),         NORM = '1', 'O' or 'o'   
                (   
                ( normI(A),         NORM = 'I' or 'i'   
                (   
                ( normF(A),         NORM = 'F', 'f', 'E' or 'e'
    where  norm1  denotes the  one norm of a matrix (maximum column sum), 
    normI  denotes the  infinity norm  of a matrix  (maximum row sum) and 
    normF  denotes the  Frobenius norm of a matrix (square root of sum of 
    squares).  Note that  max(abs(A(i,j)))  is not a  matrix norm.

Arguments

    NORM    (input) CHARACTER*1   
            Specifies the value to be returned in ZLANGE as described above.   
    A       (input) SuperMatrix*
            The M by N sparse matrix A. 


void zlaqgs_dist ( SuperMatrix A,
double *  r,
double *  c,
double  rowcnd,
double  colcnd,
double  amax,
char *  equed 
)

Purpose

    ZLAQGS_DIST equilibrates a general sparse M by N matrix A using the row
    and column scaling factors in the vectors R and C.
    See supermatrix.h for the definition of 'SuperMatrix' structure.

Arguments

    A       (input/output) SuperMatrix*
            On exit, the equilibrated matrix.  See EQUED for the form of 
            the equilibrated matrix. The type of A can be:
        Stype = SLU_NC; Dtype = SLU_Z; Mtype = SLU_GE.
    R       (input) double*, dimension (A->nrow)
            The row scale factors for A.
    C       (input) double*, dimension (A->ncol)
            The column scale factors for A.
    ROWCND  (input) double
            Ratio of the smallest R(i) to the largest R(i).
    COLCND  (input) double
            Ratio of the smallest C(i) to the largest C(i).
    AMAX    (input) double
            Absolute value of largest matrix entry.
    EQUED   (output) char*
            Specifies the form of equilibration that was done.   
            = 'N':  No equilibration   
            = 'R':  Row equilibration, i.e., A has been premultiplied by  
                    diag(R).   
            = 'C':  Column equilibration, i.e., A has been postmultiplied  
                    by diag(C).   
            = 'B':  Both row and column equilibration, i.e., A has been
                    replaced by diag(R) * A * diag(C).

Internal Parameters

    THRESH is a threshold value used to decide if row or column scaling   
    should be done based on the ratio of the row or column scaling   
    factors.  If ROWCND < THRESH, row scaling is done, and if   
    COLCND < THRESH, column scaling is done.
    LARGE and SMALL are threshold values used to decide if row scaling   
    should be done based on the absolute size of the largest matrix   
    element.  If AMAX > LARGE or AMAX < SMALL, row scaling is done.   


void zldperm_dist ( int_t  job,
int_t  n,
int_t  nnz,
int_t  colptr[],
int_t  adjncy[],
doublecomplex  nzval[],
int_t perm,
double  u[],
double  v[] 
)

Purpose

  ZLDPERM finds a row permutation so that the matrix has large
  entries on the diagonal.

Arguments

job    (input) int
       Control the action. Possible values for JOB are:
       = 1 : Compute a row permutation of the matrix so that the
             permuted matrix has as many entries on its diagonal as
             possible. The values on the diagonal are of arbitrary size.
             HSL subroutine MC21A/AD is used for this.
       = 2 : Compute a row permutation of the matrix so that the smallest 
             value on the diagonal of the permuted matrix is maximized.
       = 3 : Compute a row permutation of the matrix so that the smallest
             value on the diagonal of the permuted matrix is maximized.
             The algorithm differs from the one used for JOB = 2 and may
             have quite a different performance.
       = 4 : Compute a row permutation of the matrix so that the sum
             of the diagonal entries of the permuted matrix is maximized.
       = 5 : Compute a row permutation of the matrix so that the product
             of the diagonal entries of the permuted matrix is maximized
             and vectors to scale the matrix so that the nonzero diagonal 
             entries of the permuted matrix are one in absolute value and 
             all the off-diagonal entries are less than or equal to one in 
             absolute value.
       Restriction: 1 <= JOB <= 5.
n      (input) int
       The order of the matrix.
nnz    (input) int
       The number of nonzeros in the matrix.
adjncy (input) int*, of size nnz
       The adjacency structure of the matrix, which contains the row
       indices of the nonzeros.
colptr (input) int*, of size n+1
       The pointers to the beginning of each column in ADJNCY.
nzval  (input) doublecomplex*, of size nnz
       The nonzero values of the matrix. nzval[k] is the value of
       the entry corresponding to adjncy[k].
       It is not used if job = 1.
perm   (output) int*, of size n
       The permutation vector. perm[i] = j means row i in the
       original matrix is in row j of the permuted matrix.
u      (output) double*, of size n
       If job = 5, the natural logarithms of the row scaling factors.
v      (output) double*, of size n
       If job = 5, the natural logarithms of the column scaling factors. 
       The scaled matrix B has entries b_ij = a_ij * exp(u_i + v_j).
void zlsum_bmod ( doublecomplex lsum,
doublecomplex x,
doublecomplex xk,
int  nrhs,
int_t  k,
int_t bmod,
int_t Urbs,
Ucb_indptr_t **  Ucb_indptr,
int_t **  Ucb_valptr,
int_t xsup,
gridinfo_t grid,
LocalLU_t Llu,
MPI_Request  send_req[],
SuperLUStat_t stat 
)

Purpose

  Perform local block modifications: lsum[i] -= U_i,k * X[k].
void zlsum_fmod ( doublecomplex lsum,
doublecomplex x,
doublecomplex xk,
doublecomplex rtemp,
int  nrhs,
int  knsupc,
int_t  k,
int_t fmod,
int_t  nlb,
int_t  lptr,
int_t  luptr,
int_t xsup,
gridinfo_t grid,
LocalLU_t Llu,
MPI_Request  send_req[],
SuperLUStat_t stat 
)

Purpose

  Perform local block modifications: lsum[i] -= L_i,k * X[k].
void zPrint_CompCol_Matrix_dist ( SuperMatrix )
int zPrint_CompRowLoc_Matrix_dist ( SuperMatrix )
void zPrint_Dense_Matrix_dist ( SuperMatrix )
void zPrintLblocks ( int  ,
int_t  ,
gridinfo_t ,
Glu_persist_t ,
LocalLU_t  
)

Print the blocks in the factored matrix L.

void zPrintUblocks ( int  ,
int_t  ,
gridinfo_t ,
Glu_persist_t ,
LocalLU_t  
)

Print the blocks in the factored matrix U.

int_t zQuerySpace_dist ( int_t  n,
LUstruct_t LUstruct,
gridinfo_t grid,
SuperLUStat_t stat,
mem_usage_t mem_usage 
)
mem_usage consists of the following fields:

  • for_lu (float) The amount of space used in bytes for the L data structures.
  • total (float) The amount of space needed in bytes to perform factorization.
  • expansions (int) Number of memory expansions during the LU factorization.
void zreadhb_dist ( int  iam,
FILE *  fp,
int_t nrow,
int_t ncol,
int_t nonz,
doublecomplex **  nzval,
int_t **  rowind,
int_t **  colptr 
)

Purpose

Read a DOUBLE COMPLEX PRECISION matrix stored in Harwell-Boeing format 
as described below.
Line 1 (A72,A8) 
    Col. 1 - 72   Title (TITLE) 
 Col. 73 - 80  Key (KEY)
Line 2 (5I14) 
    Col. 1 - 14   Total number of lines excluding header (TOTCRD) 
    Col. 15 - 28  Number of lines for pointers (PTRCRD) 
    Col. 29 - 42  Number of lines for row (or variable) indices (INDCRD) 
    Col. 43 - 56  Number of lines for numerical values (VALCRD) 
 Col. 57 - 70  Number of lines for right-hand sides (RHSCRD) 
                   (including starting guesses and solution vectors 
           if present) 
                  (zero indicates no right-hand side data is present)
Line 3 (A3, 11X, 4I14) 
    Col. 1 - 3    Matrix type (see below) (MXTYPE) 
    Col. 15 - 28  Number of rows (or variables) (NROW) 
    Col. 29 - 42  Number of columns (or elements) (NCOL) 
 Col. 43 - 56  Number of row (or variable) indices (NNZERO) 
               (equal to number of entries for assembled matrices) 
    Col. 57 - 70  Number of elemental matrix entries (NELTVL) 
               (zero in the case of assembled matrices) 
Line 4 (2A16, 2A20) 
    Col. 1 - 16   Format for pointers (PTRFMT) 
 Col. 17 - 32  Format for row (or variable) indices (INDFMT) 
 Col. 33 - 52  Format for numerical values of coefficient matrix (VALFMT) 
    Col. 53 - 72 Format for numerical values of right-hand sides (RHSFMT)
Line 5 (A3, 11X, 2I14) Only present if there are right-hand sides present 
    Col. 1        Right-hand side type: 
              F for full storage or M for same format as matrix 
    Col. 2        G if a starting vector(s) (Guess) is supplied. (RHSTYP) 
    Col. 3        X if an exact solution vector(s) is supplied. 
 Col. 15 - 28  Number of right-hand sides (NRHS) 
 Col. 29 - 42  Number of row indices (NRHSIX) 
                  (ignored in case of unassembled matrices)
The three character type field on line 3 describes the matrix type. 
The following table lists the permitted values for each of the three 
characters. As an example of the type field, RSA denotes that the matrix 
is real, symmetric, and assembled.
First Character: 
 R Real matrix 
 C Complex matrix 
 P Pattern only (no numerical values supplied)
Second Character: 
 S Symmetric 
 U Unsymmetric 
 H Hermitian 
 Z Skew symmetric 
 R Rectangular
Third Character: 
 A Assembled 
 E Elemental matrices (unassembled) 
void zreadrb_dist ( FILE *  ,
int *  ,
int *  ,
int *  ,
doublecomplex **  ,
int **  ,
int **   
)
void zreadtriple ( FILE *  fp,
int_t m,
int_t n,
int_t nonz,
doublecomplex **  nzval,
int_t **  rowind,
int_t **  colptr 
)

brief

Output parameters

  (nzval, rowind, colptr): (*rowind)[*] contains the row subscripts of
     nonzeros in columns of matrix A; (*nzval)[*] the numerical values;
 column i of A is given by (*nzval)[k], k = (*rowind)[i],...,
     (*rowind)[i+1]-1.
void zSolveFinalize ( superlu_options_t ,
SOLVEstruct_t  
)

Release the resources used for the solution phase.

int zSolveInit ( superlu_options_t ,
SuperMatrix ,
int_t  [],
int_t  [],
int_t  ,
LUstruct_t ,
gridinfo_t ,
SOLVEstruct_t  
)

Initialize the data structure for the solution phase.

int ztrsm_ ( char *  ,
char *  ,
char *  ,
char *  ,
int *  ,
int *  ,
doublecomplex ,
doublecomplex ,
int *  ,
doublecomplex ,
int *   
)
int ztrsv_ ( char *  ,
char *  ,
char *  ,
int *  ,
doublecomplex ,
int *  ,
doublecomplex ,
int *   
)