IFPACK  Development
 All Classes Files Functions Variables Enumerations Friends
Ifpack_IKLU_Utils.cpp
00001 /*@HEADER
00002 // ***********************************************************************
00003 //
00004 //       Ifpack: Object-Oriented Algebraic Preconditioner Package
00005 //                 Copyright (2002) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ***********************************************************************
00040 //@HEADER
00041 */
00042 
00043 #include "Ifpack_ConfigDefs.h"
00044 #include "Ifpack_IKLU_Utils.h"
00045 
00046 //-----------------------------------------------------------------
00047 // Allocation and Destruction
00048 //-----------------------------------------------------------------
00049 
00050 /* allocate a sparse matrix (triplet form or compressed-column form) */
00051 csr *csr_spalloc (int m, int n, int nzmax, int values, int triplet)
00052 {
00053   csr *A = (csr*)calloc (1, sizeof (csr)) ;    /* allocate the csr struct */
00054   if (!A) return (NULL) ;                 /* out of memory */
00055   A->m = m ;                              /* define dimensions and nzmax */
00056   A->n = n ;
00057   A->nzmax = nzmax = CS_MAX (nzmax, 1) ;
00058   A->nz = triplet ? 0 : -1 ;              /* allocate triplet or comp.row */
00059   A->p = (int*)malloc (triplet ? CS_MAX(nzmax,1) * sizeof (int) : CS_MAX(m+1,1) * sizeof (int)) ;
00060   A->j = (int*)malloc (CS_MAX(nzmax,1) * sizeof (int)) ;
00061   A->x = values ? (double*)malloc (CS_MAX(nzmax,1) * sizeof (double)) : NULL ;
00062   return ((!A->p || !A->j || (values && !A->x)) ? csr_spfree (A) : A) ;
00063 }
00064 
00065 /* change the max # of entries sparse matrix */
00066 int csr_sprealloc (csr *A, int nzmax)
00067 {
00068     int ok, oki, okj = 1, okx = 1 ;
00069     if (!A) return (0) ;
00070     nzmax = (nzmax <= 0) ? (A->p [A->m]) : nzmax ;
00071     A->j = (int*)csr_realloc (A->j, nzmax, sizeof (int), &oki) ;
00072     if (CS_TRIPLET (A)) A->p = (int*)csr_realloc (A->p, nzmax, sizeof (int), &okj) ;
00073     if (A->x) A->x = (double*)csr_realloc (A->x, nzmax, sizeof (double), &okx) ;
00074     ok = (oki && okj && okx) ;
00075     if (ok) A->nzmax = nzmax ;
00076     return (ok) ;
00077 }
00078 
00079 /* wrapper for realloc */
00080 void *csr_realloc (void *p, int n, size_t size, int *ok)
00081 {
00082     void *pnew ;
00083     pnew = realloc (p, CS_MAX (n,1) * size) ; /* realloc the block */
00084     *ok = (pnew != NULL) ;                  /* realloc fails if pnew is NULL */
00085     return ((*ok) ? pnew : p) ;             /* return original p if failure */
00086 }
00087 
00088 /* free a sparse matrix */
00089 csr *csr_spfree (csr *A)
00090 {
00091   if (!A) return (NULL);     /* do nothing if A already NULL */
00092   if (A->p) free(A->p);
00093   if (A->j) free(A->j);
00094   if (A->x) free(A->x);
00095   if (A) free(A);
00096   return (NULL) ;      /* free the csr struct and return NULL */
00097 }
00098 
00099 /* free a symbolic factorization */
00100 css *csr_sfree (css *S)
00101 {
00102   if (!S) return (NULL) ;     /* do nothing if S already NULL */
00103   if (S->pinv) free(S->pinv);
00104   if (S->q) free(S->q);
00105   if (S->parent) free(S->parent);
00106   if (S->cp) free(S->cp);
00107   if (S->leftmost) free(S->leftmost);
00108   if (S) free(S);
00109   return (NULL) ;      /* free the css struct and return NULL */
00110 }
00111 
00112 csrn *csr_nfree (csrn *N)
00113 {
00114   if (!N) return (NULL) ;     /* do nothing if N already NULL */
00115   csr_spfree (N->L) ;
00116   csr_spfree (N->U) ;
00117   if (N->pinv) free(N->pinv);
00118   if (N->perm) free(N->perm);
00119   if (N->B) free(N->B);
00120   if (N) free(N);
00121   return (NULL) ;      /* free the csn struct and return NULL */
00122 }
00123 
00124 /* free workspace and return a sparse matrix result */
00125 csr *csr_done (csr *C, void *w, void *x, int ok)
00126 {
00127   if (w) free(w);                       /* free workspace */
00128   if (x) free(x);
00129   return (ok ? C : csr_spfree (C)) ;  /* return result if OK, else free it */
00130 }
00131 
00132 /* free workspace and return a numeric factorization (Cholesky, LU, or QR) */
00133 csrn *csr_ndone (csrn *N, csr *C, void *w, void *x, int ok)
00134 {
00135   csr_spfree (C) ;                    /* free temporary matrix */
00136   if (w) free(w);                     /* free workspace */
00137   if (x) free(x);
00138   return (ok ? N : csr_nfree (N)) ;   /* return result if OK, else free it */
00139 }
00140 
00141 /* free workspace and return int array result */
00142 int *csr_idone (int *p, csr *C, void *w, int ok)
00143 {
00144   csr_spfree (C) ;                    /* free temporary matrix */
00145   if (w) free(w);                     /* free workspace */
00146   /* return result if OK, else free it */
00147   if (ok)
00148     return p;
00149   else {
00150     if (p) free(p);
00151     return NULL;
00152   }
00153 }
00154 
00155 //-----------------------------------------------------------------
00156 // Basic CSR math routines
00157 //-----------------------------------------------------------------
00158 
00159 /* p [0..n] = cumulative sum of c [0..n-1], and then copy p [0..n-1] into c */
00160 double csr_cumsum (int *p, int *c, int n)
00161 {
00162   int i, nz = 0 ;
00163   double nz2 = 0 ;
00164   if (!p || !c) return (-1) ;     /* check inputs */
00165   for (i = 0 ; i < n ; i++)
00166     {
00167       p [i] = nz ;
00168       nz += c [i] ;
00169       nz2 += c [i] ;              /* also in double to avoid int overflow */
00170       c [i] = p [i] ;             /* also copy p[0..n-1] back into c[0..n-1]*/
00171     }
00172   p [n] = nz ;
00173   return (nz2) ;                  /* return sum (c [0..n-1]) */
00174 }
00175 
00176 /* x = x + alpha * B(i,:), where x is a dense vector and B(i,:) is sparse */
00177 int csr_scatter (const csr *B, int i, double alpha, int *w, double *x, int mark,
00178          csr *C, int nz)
00179 {
00180   int j, p, *Bp, *Bj, *Cj ;
00181   double *Bx ;
00182   if (!CS_CSC (B) || !w || !CS_CSC (C)) return (-1) ;     /* check inputs */
00183   Bp = B->p ; Bj = B->j ; Bx = B->x ; Cj = C->j ;
00184   for (p = Bp [i] ; p < Bp [i+1] ; p++)
00185     {
00186       j = Bj [p] ;                            /* B(i,j) is nonzero */
00187       if (w [j] < mark)
00188         {
00189       w [j] = mark ;                      /* j is new entry in row i */
00190       Cj [nz++] = j ;                     /* add j to pattern of C(i,:) */
00191       if (x) x [j] = alpha * Bx [p] ;     /* x(j) = alpha*B(i,j) */
00192         }
00193       else if (x) x [j] += alpha * Bx [p] ;   /* j exists in C(i,:) already */
00194     }
00195   return (nz) ;
00196 }
00197 
00198 /* C = alpha*A + beta*B */
00199 csr *csr_add (const csr *A, const csr *B, double alpha, double beta)
00200 {
00201   int p, j, nz = 0, anz, *Cp, *Cj, *Bp, m, n, bnz, *w, values ;
00202   double *x, *Bx, *Cx ;
00203   csr *C ;
00204   if (!CS_CSC (A) || !CS_CSC (B)) return (NULL) ;         /* check inputs */
00205   if ( (A->m != B->m) || (A->n != B->n) ) return (NULL);
00206   m = A->m ; anz = A->p [m] ;
00207   n = B->n ; Bp = B->p ; Bx = B->x ; bnz = Bp [m] ;
00208   w = (int*)calloc (CS_MAX(n,1), sizeof (int)) ;                       /* get workspace */
00209   values = (A->x != NULL) && (Bx != NULL) ;
00210   x = values ? (double*)malloc (n * sizeof (double)) : NULL ;    /* get workspace */
00211   C = csr_spalloc (m, n, anz + bnz, values, 0) ;          /* allocate result*/
00212   if (!C || !w || (values && !x)) return (csr_done (C, w, x, 0)) ;
00213   Cp = C->p ; Cj = C->j ; Cx = C->x ;
00214   for (j = 0 ; j < n ; j++)
00215     {
00216       Cp [j] = nz ;                   /* row j of C starts here */
00217       nz = csr_scatter (A, j, alpha, w, x, j+1, C, nz) ;   /* alpha*A(j,:)*/
00218       nz = csr_scatter (B, j, beta, w, x, j+1, C, nz) ;    /* beta*B(j,:) */
00219       if (values) for (p = Cp [j] ; p < nz ; p++) Cx [p] = x [Cj [p]] ;
00220     }
00221   Cp [m] = nz ;                       /* finalize the last row of C */
00222   csr_sprealloc (C, 0) ;              /* remove extra space from C */
00223   return (csr_done (C, w, x, 1)) ;    /* success; free workspace, return C */
00224 }
00225 
00226 /* C = A' */
00227 csr *csr_transpose (const csr *A, int values)
00228 {
00229   int p, q, i, *Cp, *Cj, n, m, *Ap, *Aj, *w ;
00230   double *Cx, *Ax ;
00231   csr *C ;
00232   if (!CS_CSC (A)) return (NULL) ;    /* check inputs */
00233   m = A->m ; n = A->n ; Ap = A->p ; Aj = A->j ; Ax = A->x ;
00234   C = csr_spalloc (n, m, Ap [m], values && Ax, 0) ;      /* allocate result */
00235   w = (int*)calloc (CS_MAX(n,1), sizeof (int)) ;         /* get workspace */
00236   if (!C || !w) return (csr_done (C, w, NULL, 0)) ;      /* out of memory */
00237   Cp = C->p ; Cj = C->j ; Cx = C->x ;
00238   for (p = 0 ; p < Ap [m] ; p++) w [Aj [p]]++ ;          /* col counts */
00239   csr_cumsum (Cp, w, n) ;                                 /* col pointers */
00240   for (i = 0 ; i < m ; i++)
00241     {
00242       for (p = Ap [i] ; p < Ap [i+1] ; p++)
00243         {
00244       Cj [q = w [Aj [p]]++] = i ; /* place A(i,j) as entry C(j,i) */
00245       if (Cx) Cx [q] = Ax [p] ;
00246         }
00247     }
00248   return (csr_done (C, w, NULL, 1)) ; /* success; free w and return C */
00249 }
00250 
00251 /* C = A*B */
00252 csr *csr_multiply (const csr *A, const csr *B)
00253 {
00254   int p, j, nz = 0, anz, *Cp, *Cj, *Ap, m, k, n, bnz, *w, values, *Aj;
00255   double *x, *Ax, *Cx ;
00256   csr *C ;
00257   if (!CS_CSC (A) || !CS_CSC (B)) return (NULL) ;      /* check inputs */
00258   k = A->n;
00259   if (k != B->m ) return(NULL);
00260   m = A->m ; anz = A->p [A->m] ;
00261   n = B->n ; Ap = A->p ; Aj = A->j ; Ax = A->x ; bnz = B->p [k] ;
00262   w = (int*)calloc (CS_MAX(n,1), sizeof (int)) ;            /* get workspace */
00263   values = (Ax != NULL) && (B->x != NULL) ;
00264   x = values ? (double*)malloc (n * sizeof (double)) : NULL ; /* get workspace */
00265   C = csr_spalloc (m, n, anz + bnz, values, 0) ;       /* allocate result */
00266   if (!C || !w || (values && !x)) return (csr_done (C, w, x, 0)) ;
00267   Cp = C->p ;
00268   for (j = 0 ; j < m ; j++)
00269     {
00270       if (nz + n > C->nzmax && !csr_sprealloc (C, 2*(C->nzmax)+m))
00271         {
00272       return (csr_done (C, w, x, 0)) ;            /* out of memory */
00273         }
00274       Cj = C->j ; Cx = C->x ;         /* C->j and C->x may be reallocated */
00275       Cp [j] = nz ;                   /* row j of C starts here */
00276       for (p = Ap [j] ; p < Ap [j+1] ; p++)
00277         {
00278       nz = csr_scatter (B, Aj [p], Ax ? Ax [p] : 1, w, x, j+1, C, nz) ;
00279         }
00280       if (values) for (p = Cp [j] ; p < nz ; p++) Cx [p] = x [Cj [p]] ;
00281     }
00282   Cp [m] = nz ;                       /* finalize the last column of C */
00283   csr_sprealloc (C, 0) ;              /* remove extra space from C */
00284   return (csr_done (C, w, x, 1)) ;    /* success; free workspace, return C */
00285 }
00286 
00287 
00288 //-----------------------------------------------------------------
00289 // Symbolic analysis
00290 //-----------------------------------------------------------------
00291 
00292 /* symbolic ordering and analysis for LU */
00293 css *csr_sqr (int order, const csr *A )
00294 {
00295     int n, ok = 1;
00296     css *S ;
00297     if (!CS_CSC (A)) return (NULL) ;        /* check inputs */
00298     n = A->n ;
00299     S = (css*)calloc(1, sizeof (css)) ;       /* allocate result S */
00300     if (!S) return (NULL) ;                 /* out of memory */
00301     S->q = csr_amd (order, A) ;             /* fill-reducing ordering */
00302     if (!S->q)
00303     {
00304         printf(" csr_sqr error no permutation\n");
00305     }
00306     if (order && !S->q) return (csr_sfree (S)) ;
00307 
00308     /* LU factorization */
00309     S->unz = (double) CS_MIN(4*(A->p [n]) + n, n * n );
00310     S->lnz = S->unz ;                       /* guess nnz(L) and nnz(U) */
00311     return (ok ? S : csr_sfree (S)) ;        /* return result S */
00312 }
00313 
00314 
00315 /* xi [top...n-1] = nodes reachable from graph of G*P' via nodes in B(:,k).
00316  * xi [n...2n-1] used as workspace */
00317 int csr_reach (csr *G, const csr *B, int k, int *xi, const int *pinv)
00318 {
00319     int p, m, top, *Bp, *Bj, *Gp ;
00320     if (!CS_CSC (G) || !CS_CSC (B) || !xi) return (-1) ;    /* check inputs */
00321     m = G->m ; Bp = B->p ; Bj = B->j ; Gp = G->p ;
00322     top = m ;
00323     for (p = Bp [k] ; p < Bp [k+1] ; p++)
00324     {
00325         if (!CS_MARKED (Gp, Bj [p]))    /* start a dfs at unmarked node i */
00326         {
00327             top = csr_dfs (Bj [p], G, top, xi, xi+m, pinv) ;
00328         }
00329     }
00330     for (p = top ; p < m ; p++) CS_MARK (Gp, xi [p]) ;  /* restore G */
00331     return (top) ;
00332 }
00333 
00334 /* depth-first-search of the graph of a csr matrix, starting at node j */
00335 int csr_dfs (int j, csr *G, int top, int *xi, int *pstack, const int *pinv)
00336 {
00337     int i, p, p2, done, jnew, head = 0, *Gp, *Gj ;
00338     if (!CS_CSC (G) || !xi || !pstack) return (-1) ;    /* check inputs */
00339     Gp = G->p ; Gj = G->j ;
00340     xi [0] = j ;                /* initialize the recursion stack */
00341     while (head >= 0)
00342     {
00343         j = xi [head] ;         /* get j from the top of the recursion stack */
00344         jnew = pinv ? (pinv [j]) : j ;
00345         if (!CS_MARKED (Gp, j))
00346         {
00347             CS_MARK (Gp, j) ;       /* mark node j as visited */
00348             pstack [head] = (jnew < 0) ? 0 : CS_UNFLIP (Gp [jnew]) ;
00349         }
00350         done = 1 ;                  /* node j done if no unvisited neighbors */
00351         p2 = (jnew < 0) ? 0 : CS_UNFLIP (Gp [jnew+1]) ;
00352         for (p = pstack [head] ; p < p2 ; p++)  /* examine all neighbors of j */
00353         {
00354             i = Gj [p] ;            /* consider neighbor node i */
00355             if (CS_MARKED (Gp, i)) continue ;   /* skip visited node i */
00356             pstack [head] = p ;     /* pause depth-first search of node j */
00357             xi [++head] = i ;       /* start dfs at node i */
00358             done = 0 ;              /* node j is not done */
00359             break ;                 /* break, to start dfs (i) */
00360         }
00361         if (done)               /* depth-first search at node j is done */
00362         {
00363             head-- ;            /* remove j from the recursion stack */
00364             xi [--top] = j ;    /* and place in the output stack */
00365         }
00366     }
00367     return (top) ;
00368 }
00369 
00370 /* depth-first search and postorder of a tree rooted at node j */
00371 int csr_tdfs (int j, int k, int *head, const int *next, int *post, int *stack)
00372 {
00373   int i, p, top = 0 ;
00374   if (!head || !next || !post || !stack) return (-1) ;    /* check inputs */
00375   stack [0] = j ;                 /* place j on the stack */
00376   while (top >= 0)                /* while (stack is not empty) */
00377     {
00378       p = stack [top] ;           /* p = top of stack */
00379       i = head [p] ;              /* i = youngest child of p */
00380       if (i == -1)
00381         {
00382       top-- ;                 /* p has no unordered children left */
00383       post [k++] = p ;        /* node p is the kth postordered node */
00384         }
00385       else
00386         {
00387       head [p] = next [i] ;   /* remove i from children of p */
00388       stack [++top] = i ;     /* start dfs on child node i */
00389         }
00390     }
00391   return (k) ;
00392 }
00393 
00394 //-----------------------------------------------------------------
00395 // LU factorization
00396 //-----------------------------------------------------------------
00397 
00398 /*
00399  * Given sparse A,  
00400  * [L,U,pinv]=lu(A, [q lnz unz]). lnz and unz can be guesses 
00401  * Hypotheses: m=n, Lj[Lp[i+1]-1]=i and Uj[Lp[i]]=i
00402  */
00403 csrn *csr_lu (const csr *A, const css *S, double tol)
00404 {
00405     csr *L, *U ;
00406     csrn *N ;
00407     double pivot, *Lx, *Ux, *x,  a, t;
00408     int *Lp, *Lj, *Up, *Uj, *pinv, *xi, *q, n, ipiv, k, top, p, i, row, lnz, unz;
00409     int debug = 0;
00410 
00411     if (!CS_CSC (A) ) printf(" error csrlu: A not csc\n");
00412     if (!CS_CSC (A) || !S) return (NULL) ;      /* check inputs */
00413     n = A->n ;
00414     if (n != A->m) return (NULL) ;      /* check inputs */
00415     q = S->q ; lnz = (int)S->lnz ; unz = (int)S->unz ;
00416     x = (double*)malloc(n * sizeof (double)) ;      /* get double workspace */
00417     xi = (int*)malloc (2*n * sizeof (int)) ;        /* get int workspace */
00418     N = (csrn*)calloc (1, sizeof (csrn)) ;          /* allocate result */
00419     if (!(A) ) printf(" error csrlu: allocation of N failed\n");
00420     if (!x || !xi || !N) return (csr_ndone (N, NULL, xi, x, 0)) ;
00421 
00422     N->L = L = csr_spalloc (n, n, lnz, 1, 0) ;      /* allocate result L */
00423     N->U = U = csr_spalloc (n, n, unz, 1, 0) ;      /* allocate result U */
00424     N->pinv = pinv = (int*)malloc (n * sizeof (int)) ;  /* allocate result pinv */
00425     N->perm = (int*)malloc (n * sizeof (int)) ;         /* allocate result perm */
00426     if (!L || !U || !pinv) return (csr_ndone (N, NULL, xi, x, 0)) ;
00427     Lp = L->p ; Up = U->p ;
00428     for (i = 0 ; i < n ; i++) x [i] = 0 ;       /* clear workspace */
00429     for (i = 0 ; i < n ; i++) pinv [i] = -1 ;       /* no rows pivotal yet */
00430     for (k = 1 ; k <= n ; k++) Up [k] = 0 ;     /* no rows of U yet */
00431     for (k = 1 ; k <= n ; k++) Lp [k] = 0 ;     /* no rows of L yet either */
00432     lnz = unz = 0 ;
00433     if( debug )
00434     {
00435        printf ("A:\n") ; csr_print (A, 0) ;  
00436     }
00437     for (k = 0 ; k < n ; k++)       /* compute L(:,k) and U(:,k) */
00438     {
00439     /* --- Triangular solve --------------------------------------------- */
00440     Lp [k] = lnz ;          /* L(:,k) starts here */
00441     Up [k] = unz ;          /* U(:,k) starts here */
00442     if ((lnz + n > L->nzmax && !csr_sprealloc (L, 2*L->nzmax + n)) ||
00443         (unz + n > U->nzmax && !csr_sprealloc (U, 2*U->nzmax + n)))
00444     {
00445         return (csr_ndone (N, NULL, xi, x, 0)) ;
00446     }
00447         Lj = L->j ; Lx = L->x ; Uj = U->j ; Ux = U->x ;
00448         row = q ? (q [k]) : k ;
00449         if( debug > 1 )
00450         {
00451             printf("--------------------------------\n");
00452             printf(" %d spsolve row=%d \n",k, row);
00453             printf(" pinv = %d %d %d %d \n", pinv[0], pinv[1], pinv[2], pinv[3]);
00454         }
00455         top = csr_spsolve (U, A, row, xi, x, pinv, 1) ; /* x = U\A(row,:) */
00456         if( debug > 1 ) printf("top=%d  x = %g %g %g %g \n", top,x[0],x[1],x[2],x[3]);
00457         /* --- Find pivot --------------------------------------------------- */
00458         ipiv = -1 ;
00459         a = -1 ;
00460         for (p = top ; p < n ; p++)
00461         {
00462             i = xi [p] ;            /* x(i) is nonzero */
00463             if (pinv [i] < 0)       /* row i is not yet pivotal */
00464             {
00465                 if ((t = fabs (x [i])) > a)
00466                 {
00467                     a = t ;         /* largest pivot candidate so far */
00468                     ipiv = i ;
00469                 }
00470             }
00471             else                    /* x(i) is the entry L(pinv[i],k) */
00472             {
00473                 Lj [lnz] = pinv [i] ;
00474                 Lx [lnz++] = x [i] ;
00475             }
00476         }
00477         if (ipiv == -1 || a <= 0) return (csr_ndone (N, NULL, xi, x, 0)) ;
00478         if (pinv [row] < 0 && fabs (x [row]) >= a*tol) ipiv = row;
00479         pivot = x [ipiv] ;          /* the chosen pivot */
00480         Lj [lnz] = k ;              /* last entry in L(:,k) is L(k,k) */
00481         Lx [lnz++] = pivot ;
00482         if( debug > 1 ) { printf ("L:") ; csr_print (L, 0) ;  }
00483 
00484         /* --- Divide by pivot ---------------------------------------------- */
00485         pinv [ipiv] = k ;           /* ipiv is the kth pivot row */
00486         Uj [unz] = ipiv ;           /* first entry in U(:,k) is U(k,k) = 1 */
00487         Ux [unz++] = 1 ;
00488         for (p = top ; p < n ; p++) /* U(k+1:n,k) = x / pivot */
00489         {
00490             i = xi [p] ;
00491             if (pinv [i] < 0)       /* x(i) is an entry in U(:,k) */
00492             {
00493                 Uj [unz] = i ;      /* save unpermuted row in U */
00494                 Ux [unz++] = x [i] / pivot ;    /* scale pivot row */
00495             }
00496             x [i] = 0 ;             /* x [0..n-1] = 0 for next k */
00497         }
00498         if( debug > 1 )
00499         {
00500             printf ("U:") ; csr_print (U, 0) ;  
00501             printf("------------------------------------\n");
00502         }
00503     }
00504     /* --- Finalize U and L ------------------------------------------------- */
00505     Lp [n] = lnz ;
00506     if( debug ) { printf ("L:") ; csr_print (L, 0) ;  }
00507     Up [n] = unz ;
00508     Uj = U->j ;                     /* fix column indices of U for final pinv */
00509     for (p = 0 ; p < unz ; p++) Uj [p] = pinv [Uj [p]] ;
00510 
00511     csr_sprealloc (L, 0) ;      /* remove extra space from L and U */
00512     csr_sprealloc (U, 0) ;
00513     if( debug ) { printf ("U:") ; csr_print (U, 0) ;  }
00514     return (csr_ndone (N, NULL, xi, x, 1)) ;    /* success */
00515 }
00516 
00517 //-----------------------------------------------------------------
00518 // Triangular solves
00519 //-----------------------------------------------------------------
00520 
00521 /* solve xG=b(k,:), where G is either upper (up=1) or lower (up=0) triangular */
00522 int csr_spsolve (csr *G, const csr *B, int k, int *xi, double *x, const int *pinv, int up)
00523 {
00524   int i, I, p, q, px, top, n, *Gp, *Gj, *Bp, *Bj ;
00525   int debug = 0;
00526   double *Gx, *Bx ;
00527   if (!CS_CSC (G) || !CS_CSC (B) || !xi || !x) return (-1) ;
00528   Gp = G->p ; Gj = G->j ; Gx = G->x ; n = G->n ;
00529   Bp = B->p ; Bj = B->j ; Bx = B->x ;
00530   top = csr_reach (G, B, k, xi, pinv) ;       /* xi[top..n-1]=Reach(B(:,k)) */
00531   for (p = top ; p < n ; p++) x [xi [p]] = 0 ;    /* clear x */
00532   for (p = Bp [k] ; p < Bp [k+1] ; p++) x [Bj [p]] = Bx [p] ; /* scatter B */
00533   if( debug )
00534     printf("solve k=%d   x= %g %g %g %g \n", k, x[0], x[1], x[2], x[3]);
00535   if( debug )
00536     printf("top=%d   xi= %d %d %d %d \n", top , xi[0], xi[1], xi[2], xi[3]);
00537   for (px = top ; px < n ; px++)
00538     {
00539       i = xi [px] ;                               /* x(i) is nonzero */
00540       I = pinv ? (pinv [i]) : i ;                 /* i maps to col I of G */
00541       if (I < 0) continue ;                       /* row I is empty */
00542       /* dmd */
00543       x [i] /= Gx [up ? (Gp [I]) : (Gp [I+1]-1)] ;/* x(i) /= G(i,i) */
00544       p = up ? (Gp [I]+1) : (Gp [I]) ;            /* up: L(i,i) 1st entry */
00545       q = up ? (Gp [I+1]) : (Gp [I+1]-1) ;        /* up: U(i,i) last entry */
00546       for ( ; p < q ; p++)
00547         {
00548       if( debug )
00549         printf("%d %d solve %d %g %g \n", px, i ,p,  Gx [p] , x [Gj [p]]  );
00550       
00551       x [Gj[p]] -= Gx [p] * x [i] ;           /* x(i) -= G(i,j) * x(j) */
00552         }
00553       if( debug )
00554     printf(" x= %g %g %g %g \n", x[0], x[1], x[2], x[3]);
00555     }
00556   return (top) ;                                  /* return top of stack */
00557 }
00558 
00559 //-----------------------------------------------------------------
00560 // AMD routine
00561 //-----------------------------------------------------------------
00562 
00563 /*
00564  * amd for csr matrices , and  order =1
00565  * Hypothesis:  m = n
00566  */
00567 /* clear w */
00568 static int csr_wclear (int mark, int lemax, int *w, int n)
00569 {
00570     int k ;
00571     if (mark < 2 || (mark + lemax < 0))
00572     {
00573         for (k = 0 ; k < n ; k++) if (w [k] != 0) w [k] = 1 ;
00574         mark = 2 ;
00575     }
00576     return (mark) ;     /* at this point, w [0..n-1] < mark holds */
00577 }
00578 
00579 /* keep off-diagonal entries; drop diagonal entries */
00580 static int csr_diag (int i, int j, double aij, void *other) { return (i != j) ; }
00581 
00582 /* p = amd(A+A') if symmetric is true, or amd(A'A) otherwise */
00583 int *csr_amd (int order, const csr *A)  /* order 0:natural, 1:Chol, 2:LU, 3:QR */
00584 {
00585   csr *C, *A2, *AT ;
00586   int *Cp, *Cj, *last, *W, *len, *nv, *next, *P, *head, *elen, *degree, *w,
00587     *hhead, *ATp, *ATj, d, dk, dext, lemax = 0, e, elenk, eln, i, j, k, k1,
00588     k2, k3, jlast, ln, dense, nzmax, mindeg = 0, nvi, nvj, nvk, mark, wnvi,
00589     ok, cnz, nel = 0, p, p1, p2, p3, p4, pj, pk, pk1, pk2, pn, q, n, m, t ;
00590   unsigned int h ;
00591   /* --- Construct matrix C ----------------------------------------------- */
00592   if (!CS_CSC (A) || order <= 0 || order > 3) return (NULL) ; /* check */
00593   AT = csr_transpose (A, 0) ;             /* compute A' */
00594   if (!AT) return (NULL) ;
00595   m = A->m ; n = A->n ;
00596   if ( n != m) return(NULL); /* For rectangular matrices, use csr_amd */
00597   dense = (int)CS_MAX (16, 10 * sqrt ((double) n)) ;   /* find dense threshold */
00598   dense = CS_MIN (n-2, dense) ;
00599   if (order == 1 && n == m)
00600     {
00601       C = csr_add (A, AT, 0, 0) ;         /* C = A+A' */
00602     }
00603   else if (order == 2)
00604     {
00605       ATp = AT->p ;                       /* drop dense columns from AT */
00606       ATj = AT->j ;
00607       for (p2 = 0, j = 0 ; j < m ; j++)
00608         {
00609       p = ATp [j] ;                   /* column j of AT starts here */
00610       ATp [j] = p2 ;                  /* new column j starts here */
00611       if (ATp [j+1] - p > dense) continue ;   /* skip dense col j */
00612       for ( ; p < ATp [j+1] ; p++) ATj [p2++] = ATj [p] ;
00613         }
00614       ATp [m] = p2 ;                      /* finalize AT */
00615       A2 = csr_transpose (AT, 0) ;        /* A2 = AT' */
00616       C = A2 ? csr_multiply (AT, A2) : NULL ; /* C=A'*A with no dense rows */
00617       csr_spfree (A2) ;
00618     }
00619   else
00620     {
00621       C = csr_multiply (AT, A) ;          /* C=A'*A */
00622     }
00623   csr_spfree (AT) ;
00624   if (!C) return (NULL) ;
00625   csr_fkeep (C, &csr_diag, NULL) ;         /* drop diagonal entries */
00626   Cp = C->p ;
00627   cnz = Cp [n] ;
00628   P = (int*)malloc (CS_MAX(n+1,1) * sizeof (int)) ;     /* allocate result */
00629   W = (int*)malloc (CS_MAX(8*(n+1),1) * sizeof (int)) ; /* get workspace */
00630   t = cnz + cnz/5 + 2*n ;                 /* add elbow room to C */
00631   if (!P || !W || !csr_sprealloc (C, t)) return (csr_idone (P, C, W, 0)) ;
00632   len  = W           ; nv     = W +   (n+1) ; next   = W + 2*(n+1) ;
00633   head = W + 3*(n+1) ; elen   = W + 4*(n+1) ; degree = W + 5*(n+1) ;
00634   w    = W + 6*(n+1) ; hhead  = W + 7*(n+1) ;
00635   last = P ;                              /* use P as workspace for last */
00636   /* --- Initialize quotient graph ---------------------------------------- */
00637   for (k = 0 ; k < n ; k++) len [k] = Cp [k+1] - Cp [k] ;
00638   len [n] = 0 ;
00639   nzmax = C->nzmax ;
00640   Cj = C->j ;
00641   for (i = 0 ; i <= n ; i++)
00642     {
00643       head [i] = -1 ;                     /* degree list i is empty */
00644       last [i] = -1 ;
00645       next [i] = -1 ;
00646       hhead [i] = -1 ;                    /* hash list i is empty */
00647       nv [i] = 1 ;                        /* node i is just one node */
00648       w [i] = 1 ;                         /* node i is alive */
00649       elen [i] = 0 ;                      /* Ek of node i is empty */
00650       degree [i] = len [i] ;              /* degree of node i */
00651     }
00652   mark = csr_wclear (0, 0, w, n) ;         /* clear w */
00653   elen [n] = -2 ;                         /* n is a dead element */
00654   Cp [n] = -1 ;                           /* n is a root of assembly tree */
00655   w [n] = 0 ;                             /* n is a dead element */
00656   /* --- Initialize degree lists ------------------------------------------ */
00657   for (i = 0 ; i < n ; i++)
00658     {
00659       d = degree [i] ;
00660       if (d == 0)                         /* node i is empty */
00661         {
00662       elen [i] = -2 ;                 /* element i is dead */
00663       nel++ ;
00664       Cp [i] = -1 ;                   /* i is a root of assemby tree */
00665       w [i] = 0 ;
00666         }
00667       else if (d > dense)                 /* node i is dense */
00668         {
00669       nv [i] = 0 ;                    /* absorb i into element n */
00670       elen [i] = -1 ;                 /* node i is dead */
00671       nel++ ;
00672       Cp [i] = CS_FLIP (n) ;
00673       nv [n]++ ;
00674         }
00675       else
00676         {
00677       if (head [d] != -1) last [head [d]] = i ;
00678       next [i] = head [d] ;           /* put node i in degree list d */
00679       head [d] = i ;
00680         }
00681     }
00682   while (nel < n)                         /* while (selecting pivots) do */
00683     {
00684       /* --- Select node of minimum approximate degree -------------------- */
00685       for (k = -1 ; mindeg < n && (k = head [mindeg]) == -1 ; mindeg++) ;
00686       if (next [k] != -1) last [next [k]] = -1 ;
00687       head [mindeg] = next [k] ;          /* remove k from degree list */
00688       elenk = elen [k] ;                  /* elenk = |Ek| */
00689       nvk = nv [k] ;                      /* # of nodes k represents */
00690       nel += nvk ;                        /* nv[k] nodes of A eliminated */
00691       /* --- Garbage collection ------------------------------------------- */
00692       if (elenk > 0 && cnz + mindeg >= nzmax)
00693         {
00694       for (j = 0 ; j < n ; j++)
00695             {
00696           if ((p = Cp [j]) >= 0)      /* j is a live node or element */
00697                 {
00698           Cp [j] = Cj [p] ;       /* save first entry of object */
00699           Cj [p] = CS_FLIP (j) ;  /* first entry is now CS_FLIP(j) */
00700                 }
00701             }
00702       for (q = 0, p = 0 ; p < cnz ; ) /* scan all of memory */
00703             {
00704           if ((j = CS_FLIP (Cj [p++])) >= 0)  /* found object j */
00705                 {
00706           Cj [q] = Cp [j] ;       /* restore first entry of object */
00707           Cp [j] = q++ ;          /* new pointer to object j */
00708           for (k3 = 0 ; k3 < len [j]-1 ; k3++) Cj [q++] = Cj [p++] ;
00709                 }
00710             }
00711       cnz = q ;                       /* Cj [cnz...nzmax-1] now free */
00712         }
00713       /* --- Construct new element ---------------------------------------- */
00714       dk = 0 ;
00715       nv [k] = -nvk ;                     /* flag k as in Lk */
00716       p = Cp [k] ;
00717       pk1 = (elenk == 0) ? p : cnz ;      /* do in place if elen[k] == 0 */
00718       pk2 = pk1 ;
00719       for (k1 = 1 ; k1 <= elenk + 1 ; k1++)
00720         {
00721       if (k1 > elenk)
00722             {
00723           e = k ;                     /* search the nodes in k */
00724           pj = p ;                    /* list of nodes starts at Cj[pj]*/
00725           ln = len [k] - elenk ;      /* length of list of nodes in k */
00726             }
00727       else
00728             {
00729           e = Cj [p++] ;              /* search the nodes in e */
00730           pj = Cp [e] ;
00731           ln = len [e] ;              /* length of list of nodes in e */
00732             }
00733       for (k2 = 1 ; k2 <= ln ; k2++)
00734             {
00735           i = Cj [pj++] ;
00736           if ((nvi = nv [i]) <= 0) continue ; /* node i dead, or seen */
00737           dk += nvi ;                 /* degree[Lk] += size of node i */
00738           nv [i] = -nvi ;             /* negate nv[i] to denote i in Lk*/
00739           Cj [pk2++] = i ;            /* place i in Lk */
00740           if (next [i] != -1) last [next [i]] = last [i] ;
00741           if (last [i] != -1)         /* remove i from degree list */
00742                 {
00743           next [last [i]] = next [i] ;
00744                 }
00745           else
00746                 {
00747           head [degree [i]] = next [i] ;
00748                 }
00749             }
00750       if (e != k)
00751             {
00752           Cp [e] = CS_FLIP (k) ;      /* absorb e into k */
00753           w [e] = 0 ;                 /* e is now a dead element */
00754             }
00755         }
00756       if (elenk != 0) cnz = pk2 ;         /* Cj [cnz...nzmax] is free */
00757       degree [k] = dk ;                   /* external degree of k - |Lk\i| */
00758       Cp [k] = pk1 ;                      /* element k is in Cj[pk1..pk2-1] */
00759       len [k] = pk2 - pk1 ;
00760       elen [k] = -2 ;                     /* k is now an element */
00761       /* --- Find set differences ----------------------------------------- */
00762       mark = csr_wclear (mark, lemax, w, n) ;  /* clear w if necessary */
00763       for (pk = pk1 ; pk < pk2 ; pk++)    /* scan 1: find |Le\Lk| */
00764         {
00765       i = Cj [pk] ;
00766       if ((eln = elen [i]) <= 0) continue ;/* skip if elen[i] empty */
00767       nvi = -nv [i] ;                      /* nv [i] was negated */
00768       wnvi = mark - nvi ;
00769       for (p = Cp [i] ; p <= Cp [i] + eln - 1 ; p++)  /* scan Ei */
00770             {
00771           e = Cj [p] ;
00772           if (w [e] >= mark)
00773                 {
00774           w [e] -= nvi ;          /* decrement |Le\Lk| */
00775                 }
00776           else if (w [e] != 0)        /* ensure e is a live element */
00777                 {
00778           w [e] = degree [e] + wnvi ; /* 1st time e seen in scan 1 */
00779                 }
00780             }
00781         }
00782       /* --- Degree update ------------------------------------------------ */
00783       for (pk = pk1 ; pk < pk2 ; pk++)    /* scan2: degree update */
00784         {
00785       i = Cj [pk] ;                   /* consider node i in Lk */
00786       p1 = Cp [i] ;
00787       p2 = p1 + elen [i] - 1 ;
00788       pn = p1 ;
00789       for (h = 0, d = 0, p = p1 ; p <= p2 ; p++)    /* scan Ei */
00790             {
00791           e = Cj [p] ;
00792           if (w [e] != 0)             /* e is an unabsorbed element */
00793                 {
00794           dext = w [e] - mark ;   /* dext = |Le\Lk| */
00795           if (dext > 0)
00796                     {
00797               d += dext ;         /* sum up the set differences */
00798               Cj [pn++] = e ;     /* keep e in Ei */
00799               h += e ;            /* compute the hash of node i */
00800                     }
00801           else
00802                     {
00803               Cp [e] = CS_FLIP (k) ;  /* aggressive absorb. e->k */
00804               w [e] = 0 ;             /* e is a dead element */
00805                     }
00806                 }
00807             }
00808       elen [i] = pn - p1 + 1 ;        /* elen[i] = |Ei| */
00809       p3 = pn ;
00810       p4 = p1 + len [i] ;
00811       for (p = p2 + 1 ; p < p4 ; p++) /* prune edges in Ai */
00812             {
00813           j = Cj [p] ;
00814           if ((nvj = nv [j]) <= 0) continue ; /* node j dead or in Lk */
00815           d += nvj ;                  /* degree(i) += |j| */
00816           Cj [pn++] = j ;             /* place j in node list of i */
00817           h += j ;                    /* compute hash for node i */
00818             }
00819       if (d == 0)                     /* check for mass elimination */
00820             {
00821           Cp [i] = CS_FLIP (k) ;      /* absorb i into k */
00822           nvi = -nv [i] ;
00823           dk -= nvi ;                 /* |Lk| -= |i| */
00824           nvk += nvi ;                /* |k| += nv[i] */
00825           nel += nvi ;
00826           nv [i] = 0 ;
00827           elen [i] = -1 ;             /* node i is dead */
00828             }
00829       else
00830             {
00831           degree [i] = CS_MIN (degree [i], d) ;   /* update degree(i) */
00832           Cj [pn] = Cj [p3] ;         /* move first node to end */
00833           Cj [p3] = Cj [p1] ;         /* move 1st el. to end of Ei */
00834           Cj [p1] = k ;               /* add k as 1st element in of Ei */
00835           
00836           len [i] = pn - p1 + 1 ;     /* new len of adj. list of node i */
00837           h %= n ;                    /* finalize hash of i */
00838           next [i] = hhead [h] ;      /* place i in hash bucket */
00839           hhead [h] = i ;
00840           last [i] = h ;              /* save hash of i in last[i] */
00841             }
00842         }                                   /* scan2 is done */
00843       degree [k] = dk ;                   /* finalize |Lk| */
00844       lemax = CS_MAX (lemax, dk) ;
00845       mark = csr_wclear (mark+lemax, lemax, w, n) ;    /* clear w */
00846       /* --- Supernode detection ------------------------------------------ */
00847       for (pk = pk1 ; pk < pk2 ; pk++)
00848         {
00849       i = Cj [pk] ;
00850       if (nv [i] >= 0) continue ;         /* skip if i is dead */
00851       h = last [i] ;                      /* scan hash bucket of node i */
00852       i = hhead [h] ;
00853       hhead [h] = -1 ;                    /* hash bucket will be empty */
00854       for ( ; i != -1 && next [i] != -1 ; i = next [i], mark++)
00855             {
00856           ln = len [i] ;
00857           eln = elen [i] ;
00858           for (p = Cp [i]+1 ; p <= Cp [i] + ln-1 ; p++) w [Cj [p]] = mark;
00859           jlast = i ;
00860           for (j = next [i] ; j != -1 ; ) /* compare i with all j */
00861                 {
00862           ok = (len [j] == ln) && (elen [j] == eln) ;
00863           for (p = Cp [j] + 1 ; ok && p <= Cp [j] + ln - 1 ; p++)
00864             {
00865               if (w [Cj [p]] != mark) ok = 0 ;    /* compare i and j*/
00866                     }
00867           if (ok)                     /* i and j are identical */
00868                     {
00869               Cp [j] = CS_FLIP (i) ;  /* absorb j into i */
00870               nv [i] += nv [j] ;
00871               nv [j] = 0 ;
00872               elen [j] = -1 ;         /* node j is dead */
00873               j = next [j] ;          /* delete j from hash bucket */
00874               next [jlast] = j ;
00875                     }
00876           else
00877                     {
00878               jlast = j ;             /* j and i are different */
00879               j = next [j] ;
00880                     }
00881                 }
00882             }
00883         }
00884       /* --- Finalize new element------------------------------------------ */
00885       for (p = pk1, pk = pk1 ; pk < pk2 ; pk++)   /* finalize Lk */
00886         {
00887       i = Cj [pk] ;
00888       if ((nvi = -nv [i]) <= 0) continue ;/* skip if i is dead */
00889       nv [i] = nvi ;                      /* restore nv[i] */
00890       d = degree [i] + dk - nvi ;         /* compute external degree(i) */
00891       d = CS_MIN (d, n - nel - nvi) ;
00892       if (head [d] != -1) last [head [d]] = i ;
00893       next [i] = head [d] ;               /* put i back in degree list */
00894       last [i] = -1 ;
00895       head [d] = i ;
00896       mindeg = CS_MIN (mindeg, d) ;       /* find new minimum degree */
00897       degree [i] = d ;
00898       Cj [p++] = i ;                      /* place i in Lk */
00899         }
00900       nv [k] = nvk ;                      /* # nodes absorbed into k */
00901       if ((len [k] = p-pk1) == 0)         /* length of adj list of element k*/
00902         {
00903       Cp [k] = -1 ;                   /* k is a root of the tree */
00904       w [k] = 0 ;                     /* k is now a dead element */
00905         }
00906       if (elenk != 0) cnz = p ;           /* free unused space in Lk */
00907     }
00908   /* --- Postordering ----------------------------------------------------- */
00909   for (i = 0 ; i < n ; i++) Cp [i] = CS_FLIP (Cp [i]) ;/* fix assembly tree */
00910   for (j = 0 ; j <= n ; j++) head [j] = -1 ;
00911   for (j = n ; j >= 0 ; j--)              /* place unordered nodes in lists */
00912     {
00913       if (nv [j] > 0) continue ;          /* skip if j is an element */
00914       next [j] = head [Cp [j]] ;          /* place j in list of its parent */
00915       head [Cp [j]] = j ;
00916     }
00917   for (e = n ; e >= 0 ; e--)              /* place elements in lists */
00918     {
00919       if (nv [e] <= 0) continue ;         /* skip unless e is an element */
00920       if (Cp [e] != -1)
00921         {
00922       next [e] = head [Cp [e]] ;      /* place e in list of its parent */
00923       head [Cp [e]] = e ;
00924         }
00925     }
00926   for (k = 0, i = 0 ; i <= n ; i++)       /* postorder the assembly tree */
00927     {
00928       if (Cp [i] == -1) k = csr_tdfs (i, k, head, next, P, w) ;
00929     }
00930   return (csr_idone (P, C, W, 1)) ;
00931 }
00932 
00933 
00934 //-----------------------------------------------------------------
00935 // Misc utilities
00936 //-----------------------------------------------------------------
00937 
00938 /* print a sparse matrix */
00939 int csr_print (const csr *A, int brief)
00940 {
00941   int p, j, m, n, nzmax, nz, nnz, *Ap, *Aj ;
00942   double *Ax ;
00943   if (!A) { printf ("(null)\n") ; return (0) ; }
00944   m = A->m ; n = A->n ; Ap = A->p ; Aj = A->j ; Ax = A->x ;
00945   nzmax = A->nzmax ; nz = A->nz ; nnz = Ap [m];
00946   if (nz < 0)
00947     {
00948       if( nnz == 0)
00949         {  /* print intermeidate matrices from csr_lu */
00950       while ((Ap[m] == 0) && (m > 0))
00951         {
00952               --m;
00953         }
00954       nnz = Ap [m];
00955         }
00956       if( nnz > 0)
00957         {
00958       printf ("%d-by-%d, nzmax: %d nnz: %d, mxnorm: %g\n", m, n, nzmax,
00959           Ap [m], csr_norm (A)) ;
00960       for (j = 0 ; j < m ; j++)
00961             {
00962           printf ("    row %d : locations %d to %d\n", j, Ap [j], Ap [j+1]-1);
00963           for (p = Ap [j] ; p < Ap [j+1] ; p++)
00964                 {
00965           printf ("      %d : %g\n", Aj [p], Ax ? Ax [p] : 1) ;
00966           if (brief && p > 20) { printf ("  ...\n") ; return (1) ; }
00967                 }
00968             }
00969         }else{
00970     printf ("%d-by-%d, zero matrix with nzmax: %d\n", m, n, nzmax);
00971       }
00972     }
00973   else
00974     {
00975       printf ("triplet: %d-by-%d, nzmax: %d nnz: %d\n", m, n, nzmax, nz) ;
00976       for (p = 0 ; p < nz ; p++)
00977         {
00978       printf ("    %d %d : %g\n", Aj [p], Ap [p], Ax ? Ax [p] : 1) ;
00979       if (brief && p > 20) { printf ("  ...\n") ; return (1) ; }
00980         }
00981     }
00982   return (1) ;
00983 }
00984 
00985 /* infinity-norm of a sparse matrix = norm(A,inf), max row sum */
00986 double csr_norm (const csr *A)
00987 {
00988     int p, j, m, *Ap ;
00989     double *Ax,  norm = 0, s ;
00990     if (!CS_CSC (A) || !A->x) return (-1) ;        /* check inputs */
00991     m = A->m ; Ap = A->p ; Ax = A->x ;
00992     for (j = 0 ; j < m ; j++)
00993     {
00994         for (s = 0, p = Ap [j] ; p < Ap [j+1] ; p++) s += fabs (Ax [p]);
00995         norm = CS_MAX (norm, s) ;
00996     }
00997     return (norm) ;
00998 }
00999 
01000 /* drop entries for which fkeep(A(i,j)) is false; return nz if OK, else -1 */
01001 int csr_fkeep (csr *A, int (*fkeep) (int, int, double, void *), void *other)
01002 {
01003     int j, p, nz = 0, m, *Ap, *Aj ;
01004     double *Ax ;
01005     if (!CS_CSC (A) || !fkeep) return (-1) ;    /* check inputs */
01006     m = A->m ; Ap = A->p ; Aj = A->j ; Ax = A->x ;
01007     for (j = 0 ; j < m ; j++)
01008     {
01009         p = Ap [j] ;                        /* get current location of col j */
01010         Ap [j] = nz ;                       /* record new location of col j */
01011         for ( ; p < Ap [j+1] ; p++)
01012         {
01013             if (fkeep (Aj [p], j, Ax ? Ax [p] : 1, other))
01014             {
01015                 if (Ax) Ax [nz] = Ax [p] ;  /* keep A(i,j) */
01016                 Aj [nz++] = Aj [p] ;
01017             }
01018         }
01019     }
01020     Ap [m] = nz ;                           /* finalize A */
01021     csr_sprealloc (A, 0) ;                  /* remove extra space from A */
01022     return (nz) ;
01023 }
 All Classes Files Functions Variables Enumerations Friends