kusano 7d535a
kusano 7d535a
<meta content="text/html;charset=UTF-8" http-equiv="Content-Type">
kusano 7d535a
<title>SuperLU: SRC/cgsitrf.c File Reference</title>
kusano 7d535a
<link href="doxygen.css" rel="stylesheet" type="text/css">
kusano 7d535a
<link href="tabs.css" rel="stylesheet" type="text/css">
kusano 7d535a
kusano 7d535a
kusano 7d535a
kusano 7d535a
  
kusano 7d535a
    
    kusano 7d535a
          
  • Main Page
  • kusano 7d535a
          
  • Data Structures
  • kusano 7d535a
          
  • Files
  • kusano 7d535a
        
    kusano 7d535a
      
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a

    SRC/cgsitrf.c File Reference

    #include "slu_cdefs.h"
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a

    Functions

    kusano 7d535a
    void cgsitrf (superlu_options_t *options, SuperMatrix *A, int relax, int panel_size, int *etree, void *work, int lwork, int *perm_c, int *perm_r, SuperMatrix *L, SuperMatrix *U, SuperLUStat_t *stat, int *info)
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a

    Function Documentation

    kusano 7d535a
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a
          
    kusano 7d535a
            
    kusano 7d535a
              void cgsitrf           
    kusano 7d535a
              (
    kusano 7d535a
              superlu_options_t
    kusano 7d535a
               options, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              SuperMatrix
    kusano 7d535a
               A, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int 
    kusano 7d535a
               relax, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int 
    kusano 7d535a
               panel_size, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int * 
    kusano 7d535a
               etree, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              void * 
    kusano 7d535a
               work, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int 
    kusano 7d535a
               lwork, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int * 
    kusano 7d535a
               perm_c, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int * 
    kusano 7d535a
               perm_r, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              SuperMatrix
    kusano 7d535a
               L, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              SuperMatrix
    kusano 7d535a
               U, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              SuperLUStat_t
    kusano 7d535a
               stat, 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              
    kusano 7d535a
              int * 
    kusano 7d535a
               info 
    kusano 7d535a
            
    kusano 7d535a
            
    kusano 7d535a
              
    kusano 7d535a
              )
    kusano 7d535a
              
    kusano 7d535a
            
    kusano 7d535a
          
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a
    kusano 7d535a

    kusano 7d535a
    kusano 7d535a
     Purpose
    kusano 7d535a
     =======

    kusano 7d535a
     CGSITRF computes an ILU factorization of a general sparse m-by-n
    kusano 7d535a
     matrix A using partial pivoting with row interchanges.
    kusano 7d535a
     The factorization has the form
    kusano 7d535a
         Pr * A = L * U
    kusano 7d535a
     where Pr is a row permutation matrix, L is lower triangular with unit
    kusano 7d535a
     diagonal elements (lower trapezoidal if A->nrow > A->ncol), and U is upper
    kusano 7d535a
     triangular (upper trapezoidal if A->nrow < A->ncol).

    kusano 7d535a
     See supermatrix.h for the definition of 'SuperMatrix' structure.

    kusano 7d535a
     Arguments
    kusano 7d535a
     =========

    kusano 7d535a
     options (input) superlu_options_t*
    kusano 7d535a
    	   The structure defines the input parameters to control
    kusano 7d535a
    	   how the ILU decomposition will be performed.

    kusano 7d535a
     A	    (input) SuperMatrix*
    kusano 7d535a
    	    Original matrix A, permuted by columns, of dimension
    kusano 7d535a
    	    (A->nrow, A->ncol). The type of A can be:
    kusano 7d535a
    	    Stype = SLU_NCP; Dtype = SLU_C; Mtype = SLU_GE.

    kusano 7d535a
     relax    (input) int
    kusano 7d535a
    	    To control degree of relaxing supernodes. If the number
    kusano 7d535a
    	    of nodes (columns) in a subtree of the elimination tree is less
    kusano 7d535a
    	    than relax, this subtree is considered as one supernode,
    kusano 7d535a
    	    regardless of the row structures of those columns.

    kusano 7d535a
     panel_size (input) int
    kusano 7d535a
    	    A panel consists of at most panel_size consecutive columns.

    kusano 7d535a
     etree    (input) int*, dimension (A->ncol)
    kusano 7d535a
    	    Elimination tree of A'*A.
    kusano 7d535a
    	    Note: etree is a vector of parent pointers for a forest whose
    kusano 7d535a
    	    vertices are the integers 0 to A->ncol-1; etree[root]==A->ncol.
    kusano 7d535a
    	    On input, the columns of A should be permuted so that the
    kusano 7d535a
    	    etree is in a certain postorder.

    kusano 7d535a
     work     (input/output) void*, size (lwork) (in bytes)
    kusano 7d535a
    	    User-supplied work space and space for the output data structures.
    kusano 7d535a
    	    Not referenced if lwork = 0;

    kusano 7d535a
     lwork   (input) int
    kusano 7d535a
    	   Specifies the size of work array in bytes.
    kusano 7d535a
    	   = 0:  allocate space internally by system malloc;
    kusano 7d535a
    	   > 0:  use user-supplied work array of length lwork in bytes,
    kusano 7d535a
    		 returns error if space runs out.
    kusano 7d535a
    	   = -1: the routine guesses the amount of space needed without
    kusano 7d535a
    		 performing the factorization, and returns it in
    kusano 7d535a
    		 *info; no other side effects.

    kusano 7d535a
     perm_c   (input) int*, dimension (A->ncol)
    kusano 7d535a
    	    Column permutation vector, which defines the
    kusano 7d535a
    	    permutation matrix Pc; perm_c[i] = j means column i of A is
    kusano 7d535a
    	    in position j in A*Pc.
    kusano 7d535a
    	    When searching for diagonal, perm_c[*] is applied to the
    kusano 7d535a
    	    row subscripts of A, so that diagonal threshold pivoting
    kusano 7d535a
    	    can find the diagonal of A, rather than that of A*Pc.

    kusano 7d535a
     perm_r   (input/output) int*, dimension (A->nrow)
    kusano 7d535a
    	    Row permutation vector which defines the permutation matrix Pr,
    kusano 7d535a
    	    perm_r[i] = j means row i of A is in position j in Pr*A.
    kusano 7d535a
    	    If options->Fact = SamePattern_SameRowPerm, the pivoting routine
    kusano 7d535a
    	       will try to use the input perm_r, unless a certain threshold
    kusano 7d535a
    	       criterion is violated. In that case, perm_r is overwritten by
    kusano 7d535a
    	       a new permutation determined by partial pivoting or diagonal
    kusano 7d535a
    	       threshold pivoting.
    kusano 7d535a
    	    Otherwise, perm_r is output argument;

    kusano 7d535a
     L	    (output) SuperMatrix*
    kusano 7d535a
    	    The factor L from the factorization Pr*A=L*U; use compressed row
    kusano 7d535a
    	    subscripts storage for supernodes, i.e., L has type:
    kusano 7d535a
    	    Stype = SLU_SC, Dtype = SLU_C, Mtype = SLU_TRLU.

    kusano 7d535a
     U	    (output) SuperMatrix*
    kusano 7d535a
    	    The factor U from the factorization Pr*A*Pc=L*U. Use column-wise
    kusano 7d535a
    	    storage scheme, i.e., U has types: Stype = SLU_NC,
    kusano 7d535a
    	    Dtype = SLU_C, Mtype = SLU_TRU.

    kusano 7d535a
     stat     (output) SuperLUStat_t*
    kusano 7d535a
    	    Record the statistics on runtime and floating-point operation count.
    kusano 7d535a
    	    See slu_util.h for the definition of 'SuperLUStat_t'.

    kusano 7d535a
     info     (output) int*
    kusano 7d535a
    	    = 0: successful exit
    kusano 7d535a
    	    < 0: if info = -i, the i-th argument had an illegal value
    kusano 7d535a
    	    > 0: if info = i, and i is
    kusano 7d535a
    	       <= A->ncol: number of zero pivots. They are replaced by small
    kusano 7d535a
    		  entries according to options->ILU_FillTol.
    kusano 7d535a
    	       > A->ncol: number of bytes allocated when memory allocation
    kusano 7d535a
    		  failure occurred, plus A->ncol. If lwork = -1, it is
    kusano 7d535a
    		  the estimated amount of space needed, plus A->ncol.

    kusano 7d535a
     ======================================================================

    kusano 7d535a
     Local Working Arrays:
    kusano 7d535a
     ======================
    kusano 7d535a
       m = number of rows in the matrix
    kusano 7d535a
       n = number of columns in the matrix

    kusano 7d535a
       marker[0:3*m-1]: marker[i] = j means that node i has been
    kusano 7d535a
    	reached when working on column j.
    kusano 7d535a
    	Storage: relative to original row subscripts
    kusano 7d535a
    	NOTE: There are 4 of them:
    kusano 7d535a
    	      marker/marker1 are used for panel dfs, see (ilu_)dpanel_dfs.c;
    kusano 7d535a
    	      marker2 is used for inner-factorization, see (ilu)_dcolumn_dfs.c;
    kusano 7d535a
    	      marker_relax(has its own space) is used for relaxed supernodes.

    kusano 7d535a
       parent[0:m-1]: parent vector used during dfs
    kusano 7d535a
    	Storage: relative to new row subscripts

    kusano 7d535a
       xplore[0:m-1]: xplore[i] gives the location of the next (dfs)
    kusano 7d535a
    	unexplored neighbor of i in lsub[*]

    kusano 7d535a
       segrep[0:nseg-1]: contains the list of supernodal representatives
    kusano 7d535a
    	in topological order of the dfs. A supernode representative is the
    kusano 7d535a
    	last column of a supernode.
    kusano 7d535a
    	The maximum size of segrep[] is n.

    kusano 7d535a
       repfnz[0:W*m-1]: for a nonzero segment U[*,j] that ends at a
    kusano 7d535a
    	supernodal representative r, repfnz[r] is the location of the first
    kusano 7d535a
    	nonzero in this segment.  It is also used during the dfs: repfnz[r]>0
    kusano 7d535a
    	indicates the supernode r has been explored.
    kusano 7d535a
    	NOTE: There are W of them, each used for one column of a panel.

    kusano 7d535a
       panel_lsub[0:W*m-1]: temporary for the nonzeros row indices below
    kusano 7d535a
    	the panel diagonal. These are filled in during dpanel_dfs(), and are
    kusano 7d535a
    	used later in the inner LU factorization within the panel.
    kusano 7d535a
    	panel_lsub[]/dense[] pair forms the SPA data structure.
    kusano 7d535a
    	NOTE: There are W of them.

    kusano 7d535a
       dense[0:W*m-1]: sparse accumulating (SPA) vector for intermediate values;
    kusano 7d535a
    		   NOTE: there are W of them.

    kusano 7d535a
       tempv[0:*]: real temporary used for dense numeric kernels;
    kusano 7d535a
    	The size of this array is defined by NUM_TEMPV() in slu_util.h.
    kusano 7d535a
    	It is also used by the dropping routine ilu_ddrop_row().
    kusano 7d535a
      
    kusano 7d535a
    kusano 7d535a

    kusano 7d535a
    kusano 7d535a

    <address style="text-align: right;"><small>Generated on Mon Nov 22 10:23:47 2010 for SuperLU by </small></address>
    kusano 7d535a
    kusano 7d535a
    doxygen 1.5.5 
    kusano 7d535a
    kusano 7d535a