xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision ad1121b0a2ffe41d56b34b046c15202e9c1073cd)
1 #define PETSCMAT_DLL
2 
3 
4 #include "../src/mat/impls/aij/seq/aij.h"
5 #include "petscbt.h"
6 #include "../src/mat/utils/freespace.h"
7 
8 EXTERN_C_BEGIN
9 #undef __FUNCT__
10 #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11 /*
12       Computes an ordering to get most of the large numerical values in the lower triangular part of the matrix
13 */
14 PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
15 {
16   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)mat->data;
17   PetscErrorCode    ierr;
18   PetscInt          i,j,jj,k, kk,n = mat->rmap->n, current = 0, newcurrent = 0,*order;
19   const PetscInt    *ai = a->i, *aj = a->j;
20   const PetscScalar *aa = a->a;
21   PetscTruth        *done;
22   PetscReal         best,past = 0,future;
23 
24   PetscFunctionBegin;
25   /* pick initial row */
26   best = -1;
27   for (i=0; i<n; i++) {
28     future = 0;
29     for (j=ai[i]; j<ai[i+1]; j++) {
30       if (aj[j] != i) future  += PetscAbsScalar(aa[j]); else past = PetscAbsScalar(aa[j]);
31     }
32     if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
33     if (past/future > best) {
34       best = past/future;
35       current = i;
36     }
37   }
38 
39   ierr = PetscMalloc(n*sizeof(PetscTruth),&done);CHKERRQ(ierr);
40   ierr = PetscMalloc(n*sizeof(PetscInt),&order);CHKERRQ(ierr);
41   ierr = PetscMemzero(done,n*sizeof(PetscTruth));CHKERRQ(ierr);
42   order[0] = current;
43   for (i=0; i<n-1; i++) {
44     done[current] = PETSC_TRUE;
45     best          = -1;
46     /* loop over all neighbors of current pivot */
47     for (j=ai[current]; j<ai[current+1]; j++) {
48       jj = aj[j];
49       if (done[jj]) continue;
50       /* loop over columns of potential next row computing weights for below and above diagonal */
51       past = future = 0.0;
52       for (k=ai[jj]; k<ai[jj+1]; k++) {
53         kk = aj[k];
54         if (done[kk]) past += PetscAbsScalar(aa[k]);
55         else if (kk != jj) future  += PetscAbsScalar(aa[k]);
56       }
57       if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
58       if (past/future > best) {
59         best = past/future;
60         newcurrent = jj;
61       }
62     }
63     if (best == -1) { /* no neighbors to select from so select best of all that remain */
64       best = -1;
65       for (k=0; k<n; k++) {
66         if (done[k]) continue;
67         future = 0;
68         past   = 0;
69         for (j=ai[k]; j<ai[k+1]; j++) {
70           kk = aj[j];
71           if (done[kk]) past += PetscAbsScalar(aa[j]);
72           else if (kk != k) future  += PetscAbsScalar(aa[j]);
73         }
74         if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
75         if (past/future > best) {
76           best = past/future;
77           newcurrent = k;
78         }
79       }
80     }
81     if (current == newcurrent) SETERRQ(PETSC_ERR_PLIB,"newcurrent cannot be current");
82     current = newcurrent;
83     order[i+1] = current;
84   }
85   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,order,irow);CHKERRQ(ierr);
86   *icol = *irow;
87   ierr = PetscObjectReference((PetscObject)*irow);CHKERRQ(ierr);
88   ierr = PetscFree(done);CHKERRQ(ierr);
89   ierr = PetscFree(order);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 EXTERN_C_END
93 
94 EXTERN_C_BEGIN
95 #undef __FUNCT__
96 #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
97 PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
98 {
99   PetscFunctionBegin;
100   *flg = PETSC_TRUE;
101   PetscFunctionReturn(0);
102 }
103 EXTERN_C_END
104 
105 EXTERN_C_BEGIN
106 #undef __FUNCT__
107 #define __FUNCT__ "MatGetFactor_seqaij_petsc"
108 PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
109 {
110   PetscInt           n = A->rmap->n;
111   PetscErrorCode     ierr;
112 
113   PetscFunctionBegin;
114   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
115   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
116   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT){
117     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
118     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJ;
119     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ;
120     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqAIJ;
121   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
122     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
123     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
124     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
125     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
126   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
127   (*B)->factor = ftype;
128   PetscFunctionReturn(0);
129 }
130 EXTERN_C_END
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
134 PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
135 {
136   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
137   IS                 isicol;
138   PetscErrorCode     ierr;
139   const PetscInt     *r,*ic;
140   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
141   PetscInt           *bi,*bj,*ajtmp;
142   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
143   PetscReal          f;
144   PetscInt           nlnk,*lnk,k,**bi_ptr;
145   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
146   PetscBT            lnkbt;
147 
148   PetscFunctionBegin;
149   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
150   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
151   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
152   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
153 
154   /* get new row pointers */
155   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
156   bi[0] = 0;
157 
158   /* bdiag is location of diagonal in factor */
159   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
160   bdiag[0] = 0;
161 
162   /* linked list for storing column indices of the active row */
163   nlnk = n + 1;
164   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
165 
166   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
167 
168   /* initial FreeSpace size is f*(ai[n]+1) */
169   f = info->fill;
170   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
171   current_space = free_space;
172 
173   for (i=0; i<n; i++) {
174     /* copy previous fill into linked list */
175     nzi = 0;
176     nnz = ai[r[i]+1] - ai[r[i]];
177     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
178     ajtmp = aj + ai[r[i]];
179     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
180     nzi += nlnk;
181 
182     /* add pivot rows into linked list */
183     row = lnk[n];
184     while (row < i) {
185       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
186       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
187       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
188       nzi += nlnk;
189       row  = lnk[row];
190     }
191     bi[i+1] = bi[i] + nzi;
192     im[i]   = nzi;
193 
194     /* mark bdiag */
195     nzbd = 0;
196     nnz  = nzi;
197     k    = lnk[n];
198     while (nnz-- && k < i){
199       nzbd++;
200       k = lnk[k];
201     }
202     bdiag[i] = bi[i] + nzbd;
203 
204     /* if free space is not available, make more free space */
205     if (current_space->local_remaining<nzi) {
206       nnz = (n - i)*nzi; /* estimated and max additional space needed */
207       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
208       reallocs++;
209     }
210 
211     /* copy data into free space, then initialize lnk */
212     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
213     bi_ptr[i] = current_space->array;
214     current_space->array           += nzi;
215     current_space->local_used      += nzi;
216     current_space->local_remaining -= nzi;
217   }
218 #if defined(PETSC_USE_INFO)
219   if (ai[n] != 0) {
220     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
221     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
222     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
223     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
224     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
225   } else {
226     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
227   }
228 #endif
229 
230   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
231   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
232 
233   /* destroy list of free space and other temporary array(s) */
234   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
235   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
236   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
237   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
238 
239   /* put together the new matrix */
240   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
241   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
242   b    = (Mat_SeqAIJ*)(B)->data;
243   b->free_a       = PETSC_TRUE;
244   b->free_ij      = PETSC_TRUE;
245   b->singlemalloc = PETSC_FALSE;
246   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
247   b->j          = bj;
248   b->i          = bi;
249   b->diag       = bdiag;
250   b->ilen       = 0;
251   b->imax       = 0;
252   b->row        = isrow;
253   b->col        = iscol;
254   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
255   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
256   b->icol       = isicol;
257   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
258 
259   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
260   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
261   b->maxnz = b->nz = bi[n] ;
262 
263   (B)->factor                = MAT_FACTOR_LU;
264   (B)->info.factor_mallocs   = reallocs;
265   (B)->info.fill_ratio_given = f;
266 
267   if (ai[n]) {
268     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
269   } else {
270     (B)->info.fill_ratio_needed = 0.0;
271   }
272   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
273   (B)->ops->solve            = MatSolve_SeqAIJ;
274   (B)->ops->solvetranspose   = MatSolveTranspose_SeqAIJ;
275   /* switch to inodes if appropriate */
276   ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 /*
281     Trouble in factorization, should we dump the original matrix?
282 */
283 #undef __FUNCT__
284 #define __FUNCT__ "MatFactorDumpMatrix"
285 PetscErrorCode MatFactorDumpMatrix(Mat A)
286 {
287   PetscErrorCode ierr;
288   PetscTruth     flg = PETSC_FALSE;
289 
290   PetscFunctionBegin;
291   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_factor_dump_on_error",&flg,PETSC_NULL);CHKERRQ(ierr);
292   if (flg) {
293     PetscViewer viewer;
294     char        filename[PETSC_MAX_PATH_LEN];
295 
296     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
297     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
298     ierr = MatView(A,viewer);CHKERRQ(ierr);
299     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
300   }
301   PetscFunctionReturn(0);
302 }
303 
304 extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec);
305 
306 /* ----------------------------------------------------------- */
307 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat,Vec,Vec);
308 extern PetscErrorCode MatSolve_SeqAIJ_iludt(Mat,Vec,Vec);
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_newdatastruct"
312 PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
313 {
314   Mat            C=B;
315   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
316   IS             isrow = b->row,isicol = b->icol;
317   PetscErrorCode ierr;
318   const PetscInt *r,*ic,*ics;
319   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
320   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
321   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
322   PetscReal      shift=info->shiftinblocks;
323   PetscTruth     row_identity, col_identity;
324 
325   PetscFunctionBegin;
326   /* printf("MatLUFactorNumeric_SeqAIJ_newdatastruct is called ...\n"); */
327   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
328   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
329   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
330   ics  = ic;
331 
332   for (i=0; i<n; i++){
333     /* zero rtmp */
334     /* L part */
335     nz    = bi[i+1] - bi[i];
336     bjtmp = bj + bi[i];
337     for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
338 
339     /* U part */
340     nz = bi[2*n-i+1] - bi[2*n-i];
341     bjtmp = bj + bi[2*n-i];
342     for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
343 
344     /* load in initial (unfactored row) */
345     nz    = ai[r[i]+1] - ai[r[i]];
346     ajtmp = aj + ai[r[i]];
347     v     = aa + ai[r[i]];
348     for (j=0; j<nz; j++) {
349       rtmp[ics[ajtmp[j]]] = v[j];
350     }
351     if (rtmp[ics[r[i]]] == 0.0){
352       rtmp[ics[r[i]]] += shift; /* shift the diagonal of the matrix */
353       /* printf("row %d, shift %g\n",i,shift); */
354     }
355 
356     /* elimination */
357     bjtmp = bj + bi[i];
358     row   = *bjtmp++;
359     nzL   = bi[i+1] - bi[i];
360     k   = 0;
361     while  (k < nzL) {
362       pc = rtmp + row;
363       if (*pc != 0.0) {
364         pv         = b->a + bdiag[row];
365         multiplier = *pc * (*pv);
366         *pc        = multiplier;
367         pj         = b->j + bi[2*n-row]; /* begining of U(row,:) */
368         pv         = b->a + bi[2*n-row];
369         nz         = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries in U(row,:), excluding diag */
370         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
371         ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
372       }
373       row = *bjtmp++; k++;
374     }
375 
376     /* finished row so stick it into b->a */
377     /* L part */
378     pv   = b->a + bi[i] ;
379     pj   = b->j + bi[i] ;
380     nz   = bi[i+1] - bi[i];
381     for (j=0; j<nz; j++) {
382       pv[j] = rtmp[pj[j]];
383     }
384 
385     /* Mark diagonal and invert diagonal for simplier triangular solves */
386     pv  = b->a + bdiag[i];
387     pj  = b->j + bdiag[i];
388     /* if (*pj != i)SETERRQ2(PETSC_ERR_SUP,"row %d != *pj %d",i,*pj) */
389     *pv = 1.0/rtmp[*pj];
390 
391     /* U part */
392     pv = b->a + bi[2*n-i];
393     pj = b->j + bi[2*n-i];
394     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
395     for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]];
396   }
397   ierr = PetscFree(rtmp);CHKERRQ(ierr);
398   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
399   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
400 
401   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
402   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
403   if (row_identity && col_identity) {
404     C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt;
405   } else {
406     C->ops->solve = MatSolve_SeqAIJ_iludt;
407   }
408 
409   C->ops->solveadd           = 0;
410   C->ops->solvetranspose     = 0;
411   C->ops->solvetransposeadd  = 0;
412   C->ops->matsolve           = 0;
413   C->assembled    = PETSC_TRUE;
414   C->preallocated = PETSC_TRUE;
415   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
421 PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
422 {
423   Mat             C=B;
424   Mat_SeqAIJ      *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
425   IS              isrow = b->row,isicol = b->icol;
426   PetscErrorCode  ierr;
427   const PetscInt   *r,*ic,*ics;
428   PetscInt        nz,row,i,j,n=A->rmap->n,diag;
429   const PetscInt  *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
430   const PetscInt  *ajtmp,*bjtmp,*diag_offset = b->diag,*pj;
431   MatScalar       *pv,*rtmp,*pc,multiplier,d;
432   const MatScalar *v,*aa=a->a;
433   PetscReal       rs=0.0;
434   LUShift_Ctx     sctx;
435   PetscInt        newshift,*ddiag;
436 
437   PetscFunctionBegin;
438   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
439   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
440   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
441   ics  = ic;
442 
443   sctx.shift_top      = 0;
444   sctx.nshift_max     = 0;
445   sctx.shift_lo       = 0;
446   sctx.shift_hi       = 0;
447   sctx.shift_fraction = 0;
448 
449   /* if both shift schemes are chosen by user, only use info->shiftpd */
450   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
451     ddiag          = a->diag;
452     sctx.shift_top = info->zeropivot;
453     for (i=0; i<n; i++) {
454       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
455       d  = (aa)[ddiag[i]];
456       rs = -PetscAbsScalar(d) - PetscRealPart(d);
457       v  = aa+ai[i];
458       nz = ai[i+1] - ai[i];
459       for (j=0; j<nz; j++)
460 	rs += PetscAbsScalar(v[j]);
461       if (rs>sctx.shift_top) sctx.shift_top = rs;
462     }
463     sctx.shift_top   *= 1.1;
464     sctx.nshift_max   = 5;
465     sctx.shift_lo     = 0.;
466     sctx.shift_hi     = 1.;
467   }
468 
469   sctx.shift_amount = 0.0;
470   sctx.nshift       = 0;
471   do {
472     sctx.lushift = PETSC_FALSE;
473     for (i=0; i<n; i++){
474       nz    = bi[i+1] - bi[i];
475       bjtmp = bj + bi[i];
476       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
477 
478       /* load in initial (unfactored row) */
479       nz    = ai[r[i]+1] - ai[r[i]];
480       ajtmp = aj + ai[r[i]];
481       v     = aa + ai[r[i]];
482       for (j=0; j<nz; j++) {
483         rtmp[ics[ajtmp[j]]] = v[j];
484       }
485       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
486       /* if (sctx.shift_amount > 0.0) printf("row %d, shift %g\n",i,sctx.shift_amount); */
487 
488       row = *bjtmp++;
489       while  (row < i) {
490         pc = rtmp + row;
491         if (*pc != 0.0) {
492           pv         = b->a + diag_offset[row];
493           pj         = b->j + diag_offset[row] + 1;
494           multiplier = *pc / *pv++;
495           *pc        = multiplier;
496           nz         = bi[row+1] - diag_offset[row] - 1;
497           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
498           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
499         }
500         row = *bjtmp++;
501       }
502       /* finished row so stick it into b->a */
503       pv   = b->a + bi[i] ;
504       pj   = b->j + bi[i] ;
505       nz   = bi[i+1] - bi[i];
506       diag = diag_offset[i] - bi[i];
507       rs   = 0.0;
508       for (j=0; j<nz; j++) {
509         pv[j] = rtmp[pj[j]];
510         rs   += PetscAbsScalar(pv[j]);
511       }
512       rs   -= PetscAbsScalar(pv[diag]);
513 
514       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
515       sctx.rs  = rs;
516       sctx.pv  = pv[diag];
517       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
518       if (newshift == 1) break;
519     }
520 
521     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
522       /*
523        * if no shift in this attempt & shifting & started shifting & can refine,
524        * then try lower shift
525        */
526       sctx.shift_hi       = sctx.shift_fraction;
527       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
528       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
529       sctx.lushift        = PETSC_TRUE;
530       sctx.nshift++;
531     }
532   } while (sctx.lushift);
533 
534   /* invert diagonal entries for simplier triangular solves */
535   for (i=0; i<n; i++) {
536     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
537   }
538   ierr = PetscFree(rtmp);CHKERRQ(ierr);
539   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
540   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
541   if (b->inode.use) {
542     C->ops->solve   = MatSolve_Inode;
543   } else {
544     PetscTruth row_identity, col_identity;
545     ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
546     ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
547     if (row_identity && col_identity) {
548       C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering;
549     } else {
550       C->ops->solve   = MatSolve_SeqAIJ;
551     }
552   }
553   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
554   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
555   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
556   C->ops->matsolve           = MatMatSolve_SeqAIJ;
557   C->assembled    = PETSC_TRUE;
558   C->preallocated = PETSC_TRUE;
559   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
560   if (sctx.nshift){
561      if (info->shiftpd) {
562       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
563     } else if (info->shiftnz) {
564       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
565     }
566   }
567   PetscFunctionReturn(0);
568 }
569 
570 /*
571    This routine implements inplace ILU(0) with row or/and column permutations.
572    Input:
573      A - original matrix
574    Output;
575      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
576          a->j (col index) is permuted by the inverse of colperm, then sorted
577          a->a reordered accordingly with a->j
578          a->diag (ptr to diagonal elements) is updated.
579 */
580 #undef __FUNCT__
581 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
582 PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
583 {
584   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
585   IS             isrow = a->row,isicol = a->icol;
586   PetscErrorCode ierr;
587   const PetscInt *r,*ic,*ics;
588   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
589   PetscInt       *ajtmp,nz,row;
590   PetscInt       *diag = a->diag,nbdiag,*pj;
591   PetscScalar    *rtmp,*pc,multiplier,d;
592   MatScalar      *v,*pv;
593   PetscReal      rs;
594   LUShift_Ctx    sctx;
595   PetscInt       newshift;
596 
597   PetscFunctionBegin;
598   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
599   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
600   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
601   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
602   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
603   ics = ic;
604 
605   sctx.shift_top      = 0;
606   sctx.nshift_max     = 0;
607   sctx.shift_lo       = 0;
608   sctx.shift_hi       = 0;
609   sctx.shift_fraction = 0;
610 
611   /* if both shift schemes are chosen by user, only use info->shiftpd */
612   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
613     sctx.shift_top = 0;
614     for (i=0; i<n; i++) {
615       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
616       d  = (a->a)[diag[i]];
617       rs = -PetscAbsScalar(d) - PetscRealPart(d);
618       v  = a->a+ai[i];
619       nz = ai[i+1] - ai[i];
620       for (j=0; j<nz; j++)
621 	rs += PetscAbsScalar(v[j]);
622       if (rs>sctx.shift_top) sctx.shift_top = rs;
623     }
624     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
625     sctx.shift_top    *= 1.1;
626     sctx.nshift_max   = 5;
627     sctx.shift_lo     = 0.;
628     sctx.shift_hi     = 1.;
629   }
630 
631   sctx.shift_amount = 0;
632   sctx.nshift       = 0;
633   do {
634     sctx.lushift = PETSC_FALSE;
635     for (i=0; i<n; i++){
636       /* load in initial unfactored row */
637       nz    = ai[r[i]+1] - ai[r[i]];
638       ajtmp = aj + ai[r[i]];
639       v     = a->a + ai[r[i]];
640       /* sort permuted ajtmp and values v accordingly */
641       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
642       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
643 
644       diag[r[i]] = ai[r[i]];
645       for (j=0; j<nz; j++) {
646         rtmp[ajtmp[j]] = v[j];
647         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
648       }
649       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
650 
651       row = *ajtmp++;
652       while  (row < i) {
653         pc = rtmp + row;
654         if (*pc != 0.0) {
655           pv         = a->a + diag[r[row]];
656           pj         = aj + diag[r[row]] + 1;
657 
658           multiplier = *pc / *pv++;
659           *pc        = multiplier;
660           nz         = ai[r[row]+1] - diag[r[row]] - 1;
661           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
662           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
663         }
664         row = *ajtmp++;
665       }
666       /* finished row so overwrite it onto a->a */
667       pv   = a->a + ai[r[i]] ;
668       pj   = aj + ai[r[i]] ;
669       nz   = ai[r[i]+1] - ai[r[i]];
670       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
671 
672       rs   = 0.0;
673       for (j=0; j<nz; j++) {
674         pv[j] = rtmp[pj[j]];
675         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
676       }
677 
678       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
679       sctx.rs  = rs;
680       sctx.pv  = pv[nbdiag];
681       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
682       if (newshift == 1) break;
683     }
684 
685     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
686       /*
687        * if no shift in this attempt & shifting & started shifting & can refine,
688        * then try lower shift
689        */
690       sctx.shift_hi        = sctx.shift_fraction;
691       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
692       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
693       sctx.lushift         = PETSC_TRUE;
694       sctx.nshift++;
695     }
696   } while (sctx.lushift);
697 
698   /* invert diagonal entries for simplier triangular solves */
699   for (i=0; i<n; i++) {
700     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
701   }
702 
703   ierr = PetscFree(rtmp);CHKERRQ(ierr);
704   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
705   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
706   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
707   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
708   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
709   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
710   A->assembled = PETSC_TRUE;
711   A->preallocated = PETSC_TRUE;
712   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
713   if (sctx.nshift){
714     if (info->shiftpd) {
715       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
716     } else if (info->shiftnz) {
717       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
718     }
719   }
720   PetscFunctionReturn(0);
721 }
722 
723 /* ----------------------------------------------------------- */
724 #undef __FUNCT__
725 #define __FUNCT__ "MatLUFactor_SeqAIJ"
726 PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
727 {
728   PetscErrorCode ierr;
729   Mat            C;
730 
731   PetscFunctionBegin;
732   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
733   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
734   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
735   A->ops->solve            = C->ops->solve;
736   A->ops->solvetranspose   = C->ops->solvetranspose;
737   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
738   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
739   PetscFunctionReturn(0);
740 }
741 /* ----------------------------------------------------------- */
742 
743 
744 #undef __FUNCT__
745 #define __FUNCT__ "MatSolve_SeqAIJ"
746 PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
747 {
748   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
749   IS                iscol = a->col,isrow = a->row;
750   PetscErrorCode    ierr;
751   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
752   PetscInt          nz;
753   const PetscInt    *rout,*cout,*r,*c;
754   PetscScalar       *x,*tmp,*tmps,sum;
755   const PetscScalar *b;
756   const MatScalar   *aa = a->a,*v;
757 
758   PetscFunctionBegin;
759   if (!n) PetscFunctionReturn(0);
760 
761   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
762   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
763   tmp  = a->solve_work;
764 
765   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
766   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
767 
768   /* forward solve the lower triangular */
769   tmp[0] = b[*r++];
770   tmps   = tmp;
771   for (i=1; i<n; i++) {
772     v   = aa + ai[i] ;
773     vi  = aj + ai[i] ;
774     nz  = a->diag[i] - ai[i];
775     sum = b[*r++];
776     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
777     tmp[i] = sum;
778   }
779 
780   /* backward solve the upper triangular */
781   for (i=n-1; i>=0; i--){
782     v   = aa + a->diag[i] + 1;
783     vi  = aj + a->diag[i] + 1;
784     nz  = ai[i+1] - a->diag[i] - 1;
785     sum = tmp[i];
786     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
787     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
788   }
789 
790   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
791   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
792   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
793   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
794   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
795   PetscFunctionReturn(0);
796 }
797 
798 #undef __FUNCT__
799 #define __FUNCT__ "MatMatSolve_SeqAIJ"
800 PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
801 {
802   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
803   IS              iscol = a->col,isrow = a->row;
804   PetscErrorCode  ierr;
805   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
806   PetscInt        nz,neq;
807   const PetscInt  *rout,*cout,*r,*c;
808   PetscScalar     *x,*b,*tmp,*tmps,sum;
809   const MatScalar *aa = a->a,*v;
810   PetscTruth      bisdense,xisdense;
811 
812   PetscFunctionBegin;
813   if (!n) PetscFunctionReturn(0);
814 
815   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
816   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
817   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
818   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
819 
820   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
821   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
822 
823   tmp  = a->solve_work;
824   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
825   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
826 
827   for (neq=0; neq<B->cmap->n; neq++){
828     /* forward solve the lower triangular */
829     tmp[0] = b[r[0]];
830     tmps   = tmp;
831     for (i=1; i<n; i++) {
832       v   = aa + ai[i] ;
833       vi  = aj + ai[i] ;
834       nz  = a->diag[i] - ai[i];
835       sum = b[r[i]];
836       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
837       tmp[i] = sum;
838     }
839     /* backward solve the upper triangular */
840     for (i=n-1; i>=0; i--){
841       v   = aa + a->diag[i] + 1;
842       vi  = aj + a->diag[i] + 1;
843       nz  = ai[i+1] - a->diag[i] - 1;
844       sum = tmp[i];
845       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
846       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
847     }
848 
849     b += n;
850     x += n;
851   }
852   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
853   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
854   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
855   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
856   ierr = PetscLogFlops(B->cmap->n*(2.0*a->nz - n));CHKERRQ(ierr);
857   PetscFunctionReturn(0);
858 }
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
862 PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
863 {
864   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
865   IS              iscol = a->col,isrow = a->row;
866   PetscErrorCode  ierr;
867   const PetscInt  *r,*c,*rout,*cout;
868   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
869   PetscInt        nz,row;
870   PetscScalar     *x,*b,*tmp,*tmps,sum;
871   const MatScalar *aa = a->a,*v;
872 
873   PetscFunctionBegin;
874   if (!n) PetscFunctionReturn(0);
875 
876   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
877   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
878   tmp  = a->solve_work;
879 
880   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
881   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
882 
883   /* forward solve the lower triangular */
884   tmp[0] = b[*r++];
885   tmps   = tmp;
886   for (row=1; row<n; row++) {
887     i   = rout[row]; /* permuted row */
888     v   = aa + ai[i] ;
889     vi  = aj + ai[i] ;
890     nz  = a->diag[i] - ai[i];
891     sum = b[*r++];
892     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
893     tmp[row] = sum;
894   }
895 
896   /* backward solve the upper triangular */
897   for (row=n-1; row>=0; row--){
898     i   = rout[row]; /* permuted row */
899     v   = aa + a->diag[i] + 1;
900     vi  = aj + a->diag[i] + 1;
901     nz  = ai[i+1] - a->diag[i] - 1;
902     sum = tmp[row];
903     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
904     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
905   }
906 
907   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
908   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
909   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
910   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
911   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
912   PetscFunctionReturn(0);
913 }
914 
915 /* ----------------------------------------------------------- */
916 #include "../src/mat/impls/aij/seq/ftn-kernels/fsolve.h"
917 #undef __FUNCT__
918 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
919 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
920 {
921   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
922   PetscErrorCode    ierr;
923   PetscInt          n = A->rmap->n;
924   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag;
925   PetscScalar       *x;
926   const PetscScalar *b;
927   const MatScalar   *aa = a->a;
928 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
929   PetscInt          adiag_i,i,nz,ai_i;
930   const PetscInt    *vi;
931   const MatScalar   *v;
932   PetscScalar       sum;
933 #endif
934 
935   PetscFunctionBegin;
936   if (!n) PetscFunctionReturn(0);
937 
938   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
939   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
940 
941 #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
942   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
943 #else
944   /* forward solve the lower triangular */
945   x[0] = b[0];
946   for (i=1; i<n; i++) {
947     ai_i = ai[i];
948     v    = aa + ai_i;
949     vi   = aj + ai_i;
950     nz   = adiag[i] - ai_i;
951     sum  = b[i];
952     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
953     x[i] = sum;
954   }
955 
956   /* backward solve the upper triangular */
957   for (i=n-1; i>=0; i--){
958     adiag_i = adiag[i];
959     v       = aa + adiag_i + 1;
960     vi      = aj + adiag_i + 1;
961     nz      = ai[i+1] - adiag_i - 1;
962     sum     = x[i];
963     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
964     x[i]    = sum*aa[adiag_i];
965   }
966 #endif
967   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
968   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
969   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "MatSolveAdd_SeqAIJ"
975 PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
976 {
977   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
978   IS              iscol = a->col,isrow = a->row;
979   PetscErrorCode  ierr;
980   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
981   PetscInt        nz;
982   const PetscInt  *rout,*cout,*r,*c;
983   PetscScalar     *x,*b,*tmp,sum;
984   const MatScalar *aa = a->a,*v;
985 
986   PetscFunctionBegin;
987   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
988 
989   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
990   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
991   tmp  = a->solve_work;
992 
993   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
994   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
995 
996   /* forward solve the lower triangular */
997   tmp[0] = b[*r++];
998   for (i=1; i<n; i++) {
999     v   = aa + ai[i] ;
1000     vi  = aj + ai[i] ;
1001     nz  = a->diag[i] - ai[i];
1002     sum = b[*r++];
1003     while (nz--) sum -= *v++ * tmp[*vi++ ];
1004     tmp[i] = sum;
1005   }
1006 
1007   /* backward solve the upper triangular */
1008   for (i=n-1; i>=0; i--){
1009     v   = aa + a->diag[i] + 1;
1010     vi  = aj + a->diag[i] + 1;
1011     nz  = ai[i+1] - a->diag[i] - 1;
1012     sum = tmp[i];
1013     while (nz--) sum -= *v++ * tmp[*vi++ ];
1014     tmp[i] = sum*aa[a->diag[i]];
1015     x[*c--] += tmp[i];
1016   }
1017 
1018   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1019   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1020   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1021   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1022   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1023 
1024   PetscFunctionReturn(0);
1025 }
1026 /* -------------------------------------------------------------------*/
1027 #undef __FUNCT__
1028 #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1029 PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1030 {
1031   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1032   IS              iscol = a->col,isrow = a->row;
1033   PetscErrorCode  ierr;
1034   const PetscInt  *rout,*cout,*r,*c;
1035   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1036   PetscInt        nz,*diag = a->diag;
1037   PetscScalar     *x,*b,*tmp,s1;
1038   const MatScalar *aa = a->a,*v;
1039 
1040   PetscFunctionBegin;
1041   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1042   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1043   tmp  = a->solve_work;
1044 
1045   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1046   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1047 
1048   /* copy the b into temp work space according to permutation */
1049   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1050 
1051   /* forward solve the U^T */
1052   for (i=0; i<n; i++) {
1053     v   = aa + diag[i] ;
1054     vi  = aj + diag[i] + 1;
1055     nz  = ai[i+1] - diag[i] - 1;
1056     s1  = tmp[i];
1057     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1058     while (nz--) {
1059       tmp[*vi++ ] -= (*v++)*s1;
1060     }
1061     tmp[i] = s1;
1062   }
1063 
1064   /* backward solve the L^T */
1065   for (i=n-1; i>=0; i--){
1066     v   = aa + diag[i] - 1 ;
1067     vi  = aj + diag[i] - 1 ;
1068     nz  = diag[i] - ai[i];
1069     s1  = tmp[i];
1070     while (nz--) {
1071       tmp[*vi-- ] -= (*v--)*s1;
1072     }
1073   }
1074 
1075   /* copy tmp into x according to permutation */
1076   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1077 
1078   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1079   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1080   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1081   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1082 
1083   ierr = PetscLogFlops(2.0*a->nz-A->cmap->n);CHKERRQ(ierr);
1084   PetscFunctionReturn(0);
1085 }
1086 
1087 #undef __FUNCT__
1088 #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1089 PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1090 {
1091   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1092   IS              iscol = a->col,isrow = a->row;
1093   PetscErrorCode  ierr;
1094   const PetscInt  *r,*c,*rout,*cout;
1095   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1096   PetscInt        nz,*diag = a->diag;
1097   PetscScalar     *x,*b,*tmp;
1098   const MatScalar *aa = a->a,*v;
1099 
1100   PetscFunctionBegin;
1101   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
1102 
1103   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1104   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1105   tmp = a->solve_work;
1106 
1107   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1108   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1109 
1110   /* copy the b into temp work space according to permutation */
1111   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1112 
1113   /* forward solve the U^T */
1114   for (i=0; i<n; i++) {
1115     v   = aa + diag[i] ;
1116     vi  = aj + diag[i] + 1;
1117     nz  = ai[i+1] - diag[i] - 1;
1118     tmp[i] *= *v++;
1119     while (nz--) {
1120       tmp[*vi++ ] -= (*v++)*tmp[i];
1121     }
1122   }
1123 
1124   /* backward solve the L^T */
1125   for (i=n-1; i>=0; i--){
1126     v   = aa + diag[i] - 1 ;
1127     vi  = aj + diag[i] - 1 ;
1128     nz  = diag[i] - ai[i];
1129     while (nz--) {
1130       tmp[*vi-- ] -= (*v--)*tmp[i];
1131     }
1132   }
1133 
1134   /* copy tmp into x according to permutation */
1135   for (i=0; i<n; i++) x[r[i]] += tmp[i];
1136 
1137   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1138   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1139   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1140   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1141 
1142   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1143   PetscFunctionReturn(0);
1144 }
1145 /* ----------------------------------------------------------------*/
1146 EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1147 EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscTruth);
1148 
1149 /*
1150    ilu(0) with natural ordering under new data structure.
1151    Factored arrays bj and ba are stored as
1152      L(0,:), L(1,:), ...,L(n-1,:),  U(n-1,:),...,U(i,:),U(i-1,:),...,U(0,:)
1153 
1154    bi=fact->i is an array of size 2n+2, in which
1155    bi+
1156      bi[i]      ->  1st entry of L(i,:),i=0,...,i-1
1157      bi[n]      ->  points to L(n-1,:)+1
1158      bi[n+1]    ->  1st entry of U(n-1,:)
1159      bi[2n-i]   ->  1st entry of U(i,:)
1160      bi[2n-i+1] ->  end of U(i,:)+1, the 1st entry of U(i-1,:)
1161      bi[2n]     ->  1st entry of U(0,:)
1162      bi[2n+1]   ->  points to U(0,:)+1
1163 
1164    U(i,:) contains diag[i] as its last entry, i.e.,
1165     U(i,:) = (u[i,i+1],...,u[i,n-1],diag[i])
1166 */
1167 #undef __FUNCT__
1168 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct"
1169 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1170 {
1171 
1172   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1173   PetscErrorCode     ierr;
1174   PetscInt           n=A->rmap->n,*ai=a->i,*aj,*adiag=a->diag;
1175   PetscInt           i,j,nz,*bi,*bj,*bdiag;
1176 
1177   PetscFunctionBegin;
1178   /* printf("MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct ...\n"); */
1179   ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
1180   b = (Mat_SeqAIJ*)(fact)->data;
1181 
1182   /* allocate matrix arrays for new data structure */
1183   ierr = PetscMalloc3(ai[n]+1,PetscScalar,&b->a,ai[n]+1,PetscInt,&b->j,2*n+2,PetscInt,&b->i);CHKERRQ(ierr);
1184   ierr = PetscLogObjectMemory(fact,ai[n]*(sizeof(PetscScalar)+sizeof(PetscInt))+(2*n+2)*sizeof(PetscInt));CHKERRQ(ierr);
1185   b->singlemalloc = PETSC_TRUE;
1186   if (!b->diag){
1187     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&b->diag);CHKERRQ(ierr);
1188   }
1189   bdiag = b->diag;
1190 
1191   if (n > 0) {
1192     ierr = PetscMemzero(b->a,(ai[n])*sizeof(PetscScalar));CHKERRQ(ierr);
1193   }
1194 
1195   /* set bi and bj with new data structure */
1196   bi = b->i;
1197   bj = b->j;
1198 
1199   /* L part */
1200   bi[0] = 0;
1201   for (i=0; i<n; i++){
1202     nz = adiag[i] - ai[i];
1203     bi[i+1] = bi[i] + nz;
1204     aj = a->j + ai[i];
1205     for (j=0; j<nz; j++){
1206       *bj = aj[j]; bj++;
1207     }
1208   }
1209 
1210   /* U part */
1211   bi[n+1] = bi[n];
1212   for (i=n-1; i>=0; i--){
1213     nz = ai[i+1] - adiag[i] - 1;
1214     bi[2*n-i+1] = bi[2*n-i] + nz + 1;
1215     aj = a->j + adiag[i] + 1;
1216     for (j=0; j<nz; j++){
1217       *bj = aj[j]; bj++;
1218     }
1219     /* diag[i] */
1220     *bj = i; bj++;
1221     bdiag[i] = bi[2*n-i+1]-1;
1222   }
1223   PetscFunctionReturn(0);
1224 }
1225 
1226 #undef __FUNCT__
1227 #define __FUNCT__ "PetscFreeSpaceContiguous_newdatastruct"
1228 PetscErrorCode PetscFreeSpaceContiguous_newdatastruct(PetscFreeSpaceList *head,PetscInt *space,PetscInt n,PetscInt *bi,PetscInt *bdiag)
1229 {
1230   PetscFreeSpaceList a;
1231   PetscErrorCode     ierr;
1232   PetscInt           row,nnz,*bj,*array,total;
1233   PetscInt           nnzL,nnzU;
1234 
1235   PetscFunctionBegin;
1236   bi[2*n+1] = bi[n];
1237   row   = 1;
1238   total = 0;
1239   nnzL  = bdiag[0];
1240   while ((*head)!=NULL) {
1241     total += (*head)->local_used;
1242     array  = (*head)->array_head;
1243 
1244     while (bi[row] <= total && row <=n){
1245       /* copy array entries into bj for row-1 */
1246       nnz  = bi[row] - bi[row-1];
1247       /* set bi[row-1] for new datastruct */
1248       if (row -1 <= 1 ){
1249         bi[row -1] = 0;
1250       } else {
1251         bi[row-1] = bi[row-2] + nnzL; /* nnzL of previous row */
1252       }
1253 
1254       /* L part */
1255       nnzL = bdiag[row-1];
1256       bj   = space+bi[row-1];
1257       ierr = PetscMemcpy(bj,array,nnzL*sizeof(PetscInt));CHKERRQ(ierr);
1258 
1259       /* diagonal entry */
1260       bdiag[row-1]        = bi[2*n-(row-1)+1]-1;
1261       space[bdiag[row-1]] = row-1;
1262 
1263       /* U part */
1264       nnzU = nnz - nnzL;
1265       bi[2*n-(row-1)] = bi[2*n-(row-1)+1] - nnzU;
1266       nnzU --; /* exclude diagonal */
1267       bj  = space + bi[2*n-(row-1)];
1268       ierr = PetscMemcpy(bj,array+nnzL+1,nnzU*sizeof(PetscInt));CHKERRQ(ierr);
1269 
1270       array += nnz;
1271       row++;
1272     }
1273 
1274     a     =  (*head)->more_space;
1275     ierr  =  PetscFree((*head)->array_head);CHKERRQ(ierr);
1276     ierr  =  PetscFree(*head);CHKERRQ(ierr);
1277     *head =  a;
1278   }
1279   bi[n] = bi[n-1] + nnzL;
1280   if (bi[n] != bi[n+1]) SETERRQ2(1,"bi[n] %d != bi[n+1] %d",bi[n],bi[n+1]);
1281   PetscFunctionReturn(0);
1282 }
1283 
1284 #undef __FUNCT__
1285 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_newdatastruct"
1286 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1287 {
1288   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1289   IS                 isicol;
1290   PetscErrorCode     ierr;
1291   const PetscInt     *r,*ic;
1292   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
1293   PetscInt           *bi,*cols,nnz,*cols_lvl;
1294   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1295   PetscInt           i,levels,diagonal_fill;
1296   PetscTruth         col_identity,row_identity;
1297   PetscReal          f;
1298   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1299   PetscBT            lnkbt;
1300   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1301   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1302   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1303   PetscTruth         missing;
1304 
1305   PetscFunctionBegin;
1306   //printf("MatILUFactorSymbolic_SeqAIJ_newdatastruct ...\n");
1307   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
1308   f             = info->fill;
1309   levels        = (PetscInt)info->levels;
1310   diagonal_fill = (PetscInt)info->diagonal_fill;
1311   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1312 
1313   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1314   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
1315 
1316   if (!levels && row_identity && col_identity) {
1317     /* special case: ilu(0) with natural ordering */
1318     ierr = MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1319     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ_newdatastruct;
1320 
1321     fact->factor = MAT_FACTOR_ILU;
1322     (fact)->info.factor_mallocs    = 0;
1323     (fact)->info.fill_ratio_given  = info->fill;
1324     (fact)->info.fill_ratio_needed = 1.0;
1325     b               = (Mat_SeqAIJ*)(fact)->data;
1326     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1327     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1328     b->row              = isrow;
1329     b->col              = iscol;
1330     b->icol             = isicol;
1331     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1332     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1333     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1334     /* ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); */
1335     PetscFunctionReturn(0);
1336   }
1337 
1338   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1339   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1340 
1341   /* get new row pointers */
1342   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1343   bi[0] = 0;
1344   /* bdiag is location of diagonal in factor */
1345   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1346   bdiag[0]  = 0;
1347 
1348   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
1349   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
1350 
1351   /* create a linked list for storing column indices of the active row */
1352   nlnk = n + 1;
1353   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1354 
1355   /* initial FreeSpace size is f*(ai[n]+1) */
1356   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
1357   current_space = free_space;
1358   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
1359   current_space_lvl = free_space_lvl;
1360 
1361   for (i=0; i<n; i++) {
1362     nzi = 0;
1363     /* copy current row into linked list */
1364     nnz  = ai[r[i]+1] - ai[r[i]];
1365     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1366     cols = aj + ai[r[i]];
1367     lnk[i] = -1; /* marker to indicate if diagonal exists */
1368     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1369     nzi += nlnk;
1370 
1371     /* make sure diagonal entry is included */
1372     if (diagonal_fill && lnk[i] == -1) {
1373       fm = n;
1374       while (lnk[fm] < i) fm = lnk[fm];
1375       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1376       lnk[fm]    = i;
1377       lnk_lvl[i] = 0;
1378       nzi++; dcount++;
1379     }
1380 
1381     /* add pivot rows into the active row */
1382     nzbd = 0;
1383     prow = lnk[n];
1384     while (prow < i) {
1385       nnz      = bdiag[prow];
1386       cols     = bj_ptr[prow] + nnz + 1;
1387       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1388       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1389       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
1390       nzi += nlnk;
1391       prow = lnk[prow];
1392       nzbd++;
1393     }
1394     bdiag[i] = nzbd;
1395     bi[i+1]  = bi[i] + nzi;
1396 
1397     /* if free space is not available, make more free space */
1398     if (current_space->local_remaining<nzi) {
1399       nnz = 2*nzi*(n - i); /* estimated and max additional space needed */
1400       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1401       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
1402       reallocs++;
1403     }
1404 
1405     /* copy data into free_space and free_space_lvl, then initialize lnk */
1406     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1407     bj_ptr[i]    = current_space->array;
1408     bjlvl_ptr[i] = current_space_lvl->array;
1409 
1410     /* make sure the active row i has diagonal entry */
1411     if (*(bj_ptr[i]+bdiag[i]) != i) {
1412       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1413     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
1414     }
1415 
1416     current_space->array           += nzi;
1417     current_space->local_used      += nzi;
1418     current_space->local_remaining -= nzi;
1419     current_space_lvl->array           += nzi;
1420     current_space_lvl->local_used      += nzi;
1421     current_space_lvl->local_remaining -= nzi;
1422   }
1423 
1424   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1425   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1426 
1427   /* destroy list of free space and other temporary arrays */
1428   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1429 
1430   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
1431   ierr = PetscFreeSpaceContiguous_newdatastruct(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
1432 
1433   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1434   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1435   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
1436 
1437 #if defined(PETSC_USE_INFO)
1438   {
1439     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1440     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1441     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1442     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1443     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1444     if (diagonal_fill) {
1445       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1446     }
1447   }
1448 #endif
1449 
1450   /* put together the new matrix */
1451   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1452   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1453   b = (Mat_SeqAIJ*)(fact)->data;
1454   b->free_a       = PETSC_TRUE;
1455   b->free_ij      = PETSC_TRUE;
1456   b->singlemalloc = PETSC_FALSE;
1457   ierr = PetscMalloc( (bi[2*n+1] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
1458   b->j          = bj;
1459   b->i          = bi;
1460   b->diag       = bdiag;
1461   b->ilen       = 0;
1462   b->imax       = 0;
1463   b->row        = isrow;
1464   b->col        = iscol;
1465   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1466   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1467   b->icol       = isicol;
1468   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1469   /* In b structure:  Free imax, ilen, old a, old j.
1470      Allocate bdiag, solve_work, new a, new j */
1471   ierr = PetscLogObjectMemory(fact,bi[2*n+1] * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1472   b->maxnz = b->nz = bi[2*n+1] ;
1473   (fact)->info.factor_mallocs    = reallocs;
1474   (fact)->info.fill_ratio_given  = f;
1475   (fact)->info.fill_ratio_needed = ((PetscReal)bi[2*n+1])/((PetscReal)ai[n]);
1476   (fact)->ops->lufactornumeric   = MatLUFactorNumeric_SeqAIJ_newdatastruct;
1477   /* ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); */
1478   PetscFunctionReturn(0);
1479 }
1480 
1481 #undef __FUNCT__
1482 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
1483 PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1484 {
1485   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1486   IS                 isicol;
1487   PetscErrorCode     ierr;
1488   const PetscInt     *r,*ic;
1489   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
1490   PetscInt           *bi,*cols,nnz,*cols_lvl;
1491   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1492   PetscInt           i,levels,diagonal_fill;
1493   PetscTruth         col_identity,row_identity;
1494   PetscReal          f;
1495   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1496   PetscBT            lnkbt;
1497   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1498   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1499   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1500   PetscTruth         missing;
1501   PetscTruth         newdatastruct=PETSC_FALSE;
1502 
1503   PetscFunctionBegin;
1504   ierr = PetscOptionsGetTruth(PETSC_NULL,"-ilu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
1505   if (newdatastruct){
1506     ierr = MatILUFactorSymbolic_SeqAIJ_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1507     PetscFunctionReturn(0);
1508   }
1509 
1510   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
1511   f             = info->fill;
1512   levels        = (PetscInt)info->levels;
1513   diagonal_fill = (PetscInt)info->diagonal_fill;
1514   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1515 
1516   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1517   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
1518   if (!levels && row_identity && col_identity) { /* special case: ilu(0) with natural ordering */
1519     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
1520     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1521 
1522     fact->factor = MAT_FACTOR_ILU;
1523     (fact)->info.factor_mallocs    = 0;
1524     (fact)->info.fill_ratio_given  = info->fill;
1525     (fact)->info.fill_ratio_needed = 1.0;
1526     b               = (Mat_SeqAIJ*)(fact)->data;
1527     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1528     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1529     b->row              = isrow;
1530     b->col              = iscol;
1531     b->icol             = isicol;
1532     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1533     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1534     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1535     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1536     PetscFunctionReturn(0);
1537   }
1538 
1539   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1540   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1541 
1542   /* get new row pointers */
1543   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1544   bi[0] = 0;
1545   /* bdiag is location of diagonal in factor */
1546   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1547   bdiag[0]  = 0;
1548 
1549   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
1550   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
1551 
1552   /* create a linked list for storing column indices of the active row */
1553   nlnk = n + 1;
1554   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1555 
1556   /* initial FreeSpace size is f*(ai[n]+1) */
1557   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
1558   current_space = free_space;
1559   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
1560   current_space_lvl = free_space_lvl;
1561 
1562   for (i=0; i<n; i++) {
1563     nzi = 0;
1564     /* copy current row into linked list */
1565     nnz  = ai[r[i]+1] - ai[r[i]];
1566     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1567     cols = aj + ai[r[i]];
1568     lnk[i] = -1; /* marker to indicate if diagonal exists */
1569     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1570     nzi += nlnk;
1571 
1572     /* make sure diagonal entry is included */
1573     if (diagonal_fill && lnk[i] == -1) {
1574       fm = n;
1575       while (lnk[fm] < i) fm = lnk[fm];
1576       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1577       lnk[fm]    = i;
1578       lnk_lvl[i] = 0;
1579       nzi++; dcount++;
1580     }
1581 
1582     /* add pivot rows into the active row */
1583     nzbd = 0;
1584     prow = lnk[n];
1585     while (prow < i) {
1586       nnz      = bdiag[prow];
1587       cols     = bj_ptr[prow] + nnz + 1;
1588       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1589       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1590       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
1591       nzi += nlnk;
1592       prow = lnk[prow];
1593       nzbd++;
1594     }
1595     bdiag[i] = nzbd;
1596     bi[i+1]  = bi[i] + nzi;
1597 
1598     /* if free space is not available, make more free space */
1599     if (current_space->local_remaining<nzi) {
1600       nnz = nzi*(n - i); /* estimated and max additional space needed */
1601       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1602       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
1603       reallocs++;
1604     }
1605 
1606     /* copy data into free_space and free_space_lvl, then initialize lnk */
1607     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1608     bj_ptr[i]    = current_space->array;
1609     bjlvl_ptr[i] = current_space_lvl->array;
1610 
1611     /* make sure the active row i has diagonal entry */
1612     if (*(bj_ptr[i]+bdiag[i]) != i) {
1613       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1614     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
1615     }
1616 
1617     current_space->array           += nzi;
1618     current_space->local_used      += nzi;
1619     current_space->local_remaining -= nzi;
1620     current_space_lvl->array           += nzi;
1621     current_space_lvl->local_used      += nzi;
1622     current_space_lvl->local_remaining -= nzi;
1623   }
1624 
1625   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1626   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1627 
1628   /* destroy list of free space and other temporary arrays */
1629   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1630   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); /* copy free_space -> bj */
1631   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1632   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1633   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
1634 
1635 #if defined(PETSC_USE_INFO)
1636   {
1637     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1638     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1639     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1640     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1641     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1642     if (diagonal_fill) {
1643       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1644     }
1645   }
1646 #endif
1647 
1648   /* put together the new matrix */
1649   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1650   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1651   b = (Mat_SeqAIJ*)(fact)->data;
1652   b->free_a       = PETSC_TRUE;
1653   b->free_ij      = PETSC_TRUE;
1654   b->singlemalloc = PETSC_FALSE;
1655   ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
1656   b->j          = bj;
1657   b->i          = bi;
1658   for (i=0; i<n; i++) bdiag[i] += bi[i];
1659   b->diag       = bdiag;
1660   b->ilen       = 0;
1661   b->imax       = 0;
1662   b->row        = isrow;
1663   b->col        = iscol;
1664   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1665   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1666   b->icol       = isicol;
1667   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1668   /* In b structure:  Free imax, ilen, old a, old j.
1669      Allocate bdiag, solve_work, new a, new j */
1670   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1671   b->maxnz             = b->nz = bi[n] ;
1672   (fact)->info.factor_mallocs    = reallocs;
1673   (fact)->info.fill_ratio_given  = f;
1674   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1675   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1676   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 #include "../src/mat/impls/sbaij/seq/sbaij.h"
1681 #undef __FUNCT__
1682 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1683 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1684 {
1685   Mat            C = B;
1686   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1687   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
1688   IS             ip=b->row,iip = b->icol;
1689   PetscErrorCode ierr;
1690   const PetscInt *rip,*riip;
1691   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
1692   PetscInt       *ai=a->i,*aj=a->j;
1693   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1694   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1695   PetscReal      zeropivot,rs,shiftnz;
1696   PetscReal      shiftpd;
1697   ChShift_Ctx    sctx;
1698   PetscInt       newshift;
1699   PetscTruth     perm_identity;
1700 
1701   PetscFunctionBegin;
1702 
1703   shiftnz   = info->shiftnz;
1704   shiftpd   = info->shiftpd;
1705   zeropivot = info->zeropivot;
1706 
1707   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1708   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
1709 
1710   /* initialization */
1711   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1712   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1713   jl   = il + mbs;
1714   rtmp = (MatScalar*)(jl + mbs);
1715 
1716   sctx.shift_amount = 0;
1717   sctx.nshift       = 0;
1718   do {
1719     sctx.chshift = PETSC_FALSE;
1720     for (i=0; i<mbs; i++) {
1721       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1722     }
1723 
1724     for (k = 0; k<mbs; k++){
1725       bval = ba + bi[k];
1726       /* initialize k-th row by the perm[k]-th row of A */
1727       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1728       for (j = jmin; j < jmax; j++){
1729         col = riip[aj[j]];
1730         if (col >= k){ /* only take upper triangular entry */
1731           rtmp[col] = aa[j];
1732           *bval++  = 0.0; /* for in-place factorization */
1733         }
1734       }
1735       /* shift the diagonal of the matrix */
1736       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1737 
1738       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1739       dk = rtmp[k];
1740       i = jl[k]; /* first row to be added to k_th row  */
1741 
1742       while (i < k){
1743         nexti = jl[i]; /* next row to be added to k_th row */
1744 
1745         /* compute multiplier, update diag(k) and U(i,k) */
1746         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1747         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1748         dk += uikdi*ba[ili];
1749         ba[ili] = uikdi; /* -U(i,k) */
1750 
1751         /* add multiple of row i to k-th row */
1752         jmin = ili + 1; jmax = bi[i+1];
1753         if (jmin < jmax){
1754           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1755           /* update il and jl for row i */
1756           il[i] = jmin;
1757           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1758         }
1759         i = nexti;
1760       }
1761 
1762       /* shift the diagonals when zero pivot is detected */
1763       /* compute rs=sum of abs(off-diagonal) */
1764       rs   = 0.0;
1765       jmin = bi[k]+1;
1766       nz   = bi[k+1] - jmin;
1767       bcol = bj + jmin;
1768       while (nz--){
1769         rs += PetscAbsScalar(rtmp[*bcol]);
1770         bcol++;
1771       }
1772 
1773       sctx.rs = rs;
1774       sctx.pv = dk;
1775       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1776 
1777       if (newshift == 1) {
1778         if (!sctx.shift_amount) {
1779           sctx.shift_amount = 1e-5;
1780         }
1781         break;
1782       }
1783 
1784       /* copy data into U(k,:) */
1785       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1786       jmin = bi[k]+1; jmax = bi[k+1];
1787       if (jmin < jmax) {
1788         for (j=jmin; j<jmax; j++){
1789           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1790         }
1791         /* add the k-th row into il and jl */
1792         il[k] = jmin;
1793         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1794       }
1795     }
1796   } while (sctx.chshift);
1797   ierr = PetscFree(il);CHKERRQ(ierr);
1798 
1799   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1800   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1801 
1802   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1803   if (perm_identity){
1804     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1805     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1806     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1807     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1808   } else {
1809     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1810     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1811     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1812     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1813   }
1814 
1815   C->assembled    = PETSC_TRUE;
1816   C->preallocated = PETSC_TRUE;
1817   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1818   if (sctx.nshift){
1819     if (shiftnz) {
1820       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1821     } else if (shiftpd) {
1822       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1823     }
1824   }
1825   PetscFunctionReturn(0);
1826 }
1827 
1828 #undef __FUNCT__
1829 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1830 PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1831 {
1832   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1833   Mat_SeqSBAIJ       *b;
1834   PetscErrorCode     ierr;
1835   PetscTruth         perm_identity,missing;
1836   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
1837   const PetscInt     *rip,*riip;
1838   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1839   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
1840   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1841   PetscReal          fill=info->fill,levels=info->levels;
1842   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1843   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1844   PetscBT            lnkbt;
1845   IS                 iperm;
1846 
1847   PetscFunctionBegin;
1848   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
1849   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1850   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1851   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1852   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
1853 
1854   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1855   ui[0] = 0;
1856 
1857   /* ICC(0) without matrix ordering: simply copies fill pattern */
1858   if (!levels && perm_identity) {
1859 
1860     for (i=0; i<am; i++) {
1861       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1862     }
1863     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1864     cols = uj;
1865     for (i=0; i<am; i++) {
1866       aj    = a->j + a->diag[i];
1867       ncols = ui[i+1] - ui[i];
1868       for (j=0; j<ncols; j++) *cols++ = *aj++;
1869     }
1870   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1871     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1872     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1873 
1874     /* initialization */
1875     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
1876 
1877     /* jl: linked list for storing indices of the pivot rows
1878        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1879     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
1880     il         = jl + am;
1881     uj_ptr     = (PetscInt**)(il + am);
1882     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
1883     for (i=0; i<am; i++){
1884       jl[i] = am; il[i] = 0;
1885     }
1886 
1887     /* create and initialize a linked list for storing column indices of the active row k */
1888     nlnk = am + 1;
1889     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1890 
1891     /* initial FreeSpace size is fill*(ai[am]+1) */
1892     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
1893     current_space = free_space;
1894     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
1895     current_space_lvl = free_space_lvl;
1896 
1897     for (k=0; k<am; k++){  /* for each active row k */
1898       /* initialize lnk by the column indices of row rip[k] of A */
1899       nzk   = 0;
1900       ncols = ai[rip[k]+1] - ai[rip[k]];
1901       if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
1902       ncols_upper = 0;
1903       for (j=0; j<ncols; j++){
1904         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1905         if (riip[i] >= k){ /* only take upper triangular entry */
1906           ajtmp[ncols_upper] = i;
1907           ncols_upper++;
1908         }
1909       }
1910       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1911       nzk += nlnk;
1912 
1913       /* update lnk by computing fill-in for each pivot row to be merged in */
1914       prow = jl[k]; /* 1st pivot row */
1915 
1916       while (prow < k){
1917         nextprow = jl[prow];
1918 
1919         /* merge prow into k-th row */
1920         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
1921         jmax = ui[prow+1];
1922         ncols = jmax-jmin;
1923         i     = jmin - ui[prow];
1924         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1925         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
1926         j     = *(uj - 1);
1927         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
1928         nzk += nlnk;
1929 
1930         /* update il and jl for prow */
1931         if (jmin < jmax){
1932           il[prow] = jmin;
1933           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
1934         }
1935         prow = nextprow;
1936       }
1937 
1938       /* if free space is not available, make more free space */
1939       if (current_space->local_remaining<nzk) {
1940         i = am - k + 1; /* num of unfactored rows */
1941         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1942         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1943         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
1944         reallocs++;
1945       }
1946 
1947       /* copy data into free_space and free_space_lvl, then initialize lnk */
1948       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
1949       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1950 
1951       /* add the k-th row into il and jl */
1952       if (nzk > 1){
1953         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
1954         jl[k] = jl[i]; jl[i] = k;
1955         il[k] = ui[k] + 1;
1956       }
1957       uj_ptr[k]     = current_space->array;
1958       uj_lvl_ptr[k] = current_space_lvl->array;
1959 
1960       current_space->array           += nzk;
1961       current_space->local_used      += nzk;
1962       current_space->local_remaining -= nzk;
1963 
1964       current_space_lvl->array           += nzk;
1965       current_space_lvl->local_used      += nzk;
1966       current_space_lvl->local_remaining -= nzk;
1967 
1968       ui[k+1] = ui[k] + nzk;
1969     }
1970 
1971 #if defined(PETSC_USE_INFO)
1972     if (ai[am] != 0) {
1973       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1974       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1975       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1976       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
1977     } else {
1978       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1979     }
1980 #endif
1981 
1982     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1983     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
1984     ierr = PetscFree(jl);CHKERRQ(ierr);
1985     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
1986 
1987     /* destroy list of free space and other temporary array(s) */
1988     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1989     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1990     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1991     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1992 
1993   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
1994 
1995   /* put together the new matrix in MATSEQSBAIJ format */
1996 
1997   b    = (Mat_SeqSBAIJ*)(fact)->data;
1998   b->singlemalloc = PETSC_FALSE;
1999   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
2000   b->j    = uj;
2001   b->i    = ui;
2002   b->diag = 0;
2003   b->ilen = 0;
2004   b->imax = 0;
2005   b->row  = perm;
2006   b->col  = perm;
2007   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2008   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2009   b->icol = iperm;
2010   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2011   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2012   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2013   b->maxnz   = b->nz = ui[am];
2014   b->free_a  = PETSC_TRUE;
2015   b->free_ij = PETSC_TRUE;
2016 
2017   (fact)->info.factor_mallocs    = reallocs;
2018   (fact)->info.fill_ratio_given  = fill;
2019   if (ai[am] != 0) {
2020     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2021   } else {
2022     (fact)->info.fill_ratio_needed = 0.0;
2023   }
2024   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2025   PetscFunctionReturn(0);
2026 }
2027 
2028 #undef __FUNCT__
2029 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
2030 PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2031 {
2032   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2033   Mat_SeqSBAIJ       *b;
2034   PetscErrorCode     ierr;
2035   PetscTruth         perm_identity;
2036   PetscReal          fill = info->fill;
2037   const PetscInt     *rip,*riip;
2038   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
2039   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
2040   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
2041   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2042   PetscBT            lnkbt;
2043   IS                 iperm;
2044 
2045   PetscFunctionBegin;
2046   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
2047   /* check whether perm is the identity mapping */
2048   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
2049   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
2050   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
2051   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
2052 
2053   /* initialization */
2054   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
2055   ui[0] = 0;
2056 
2057   /* jl: linked list for storing indices of the pivot rows
2058      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2059   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
2060   il     = jl + am;
2061   cols   = il + am;
2062   ui_ptr = (PetscInt**)(cols + am);
2063   for (i=0; i<am; i++){
2064     jl[i] = am; il[i] = 0;
2065   }
2066 
2067   /* create and initialize a linked list for storing column indices of the active row k */
2068   nlnk = am + 1;
2069   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2070 
2071   /* initial FreeSpace size is fill*(ai[am]+1) */
2072   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
2073   current_space = free_space;
2074 
2075   for (k=0; k<am; k++){  /* for each active row k */
2076     /* initialize lnk by the column indices of row rip[k] of A */
2077     nzk   = 0;
2078     ncols = ai[rip[k]+1] - ai[rip[k]];
2079     if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2080     ncols_upper = 0;
2081     for (j=0; j<ncols; j++){
2082       i = riip[*(aj + ai[rip[k]] + j)];
2083       if (i >= k){ /* only take upper triangular entry */
2084         cols[ncols_upper] = i;
2085         ncols_upper++;
2086       }
2087     }
2088     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2089     nzk += nlnk;
2090 
2091     /* update lnk by computing fill-in for each pivot row to be merged in */
2092     prow = jl[k]; /* 1st pivot row */
2093 
2094     while (prow < k){
2095       nextprow = jl[prow];
2096       /* merge prow into k-th row */
2097       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
2098       jmax = ui[prow+1];
2099       ncols = jmax-jmin;
2100       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2101       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2102       nzk += nlnk;
2103 
2104       /* update il and jl for prow */
2105       if (jmin < jmax){
2106         il[prow] = jmin;
2107         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
2108       }
2109       prow = nextprow;
2110     }
2111 
2112     /* if free space is not available, make more free space */
2113     if (current_space->local_remaining<nzk) {
2114       i = am - k + 1; /* num of unfactored rows */
2115       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2116       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
2117       reallocs++;
2118     }
2119 
2120     /* copy data into free space, then initialize lnk */
2121     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
2122 
2123     /* add the k-th row into il and jl */
2124     if (nzk-1 > 0){
2125       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2126       jl[k] = jl[i]; jl[i] = k;
2127       il[k] = ui[k] + 1;
2128     }
2129     ui_ptr[k] = current_space->array;
2130     current_space->array           += nzk;
2131     current_space->local_used      += nzk;
2132     current_space->local_remaining -= nzk;
2133 
2134     ui[k+1] = ui[k] + nzk;
2135   }
2136 
2137 #if defined(PETSC_USE_INFO)
2138   if (ai[am] != 0) {
2139     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
2140     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
2141     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
2142     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
2143   } else {
2144      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
2145   }
2146 #endif
2147 
2148   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
2149   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
2150   ierr = PetscFree(jl);CHKERRQ(ierr);
2151 
2152   /* destroy list of free space and other temporary array(s) */
2153   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2154   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
2155   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2156 
2157   /* put together the new matrix in MATSEQSBAIJ format */
2158 
2159   b = (Mat_SeqSBAIJ*)(fact)->data;
2160   b->singlemalloc = PETSC_FALSE;
2161   b->free_a       = PETSC_TRUE;
2162   b->free_ij      = PETSC_TRUE;
2163   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
2164   b->j    = uj;
2165   b->i    = ui;
2166   b->diag = 0;
2167   b->ilen = 0;
2168   b->imax = 0;
2169   b->row  = perm;
2170   b->col  = perm;
2171   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2172   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2173   b->icol = iperm;
2174   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2175   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2176   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2177   b->maxnz = b->nz = ui[am];
2178 
2179   (fact)->info.factor_mallocs    = reallocs;
2180   (fact)->info.fill_ratio_given  = fill;
2181   if (ai[am] != 0) {
2182     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2183   } else {
2184     (fact)->info.fill_ratio_needed = 0.0;
2185   }
2186   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2187   PetscFunctionReturn(0);
2188 }
2189 
2190 #undef __FUNCT__
2191 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_iludt"
2192 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat A,Vec bb,Vec xx)
2193 {
2194   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2195   PetscErrorCode    ierr;
2196   PetscInt          n = A->rmap->n;
2197   const PetscInt    *ai = a->i,*aj = a->j,*vi;
2198   PetscScalar       *x,sum;
2199   const PetscScalar *b;
2200   const MatScalar   *aa = a->a,*v;
2201   PetscInt          i,nz;
2202 
2203   PetscFunctionBegin;
2204   if (!n) PetscFunctionReturn(0);
2205 
2206   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2207   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2208 
2209   /* forward solve the lower triangular */
2210   x[0] = b[0];
2211   v    = aa;
2212   vi   = aj;
2213   for (i=1; i<n; i++) {
2214     nz  = ai[i+1] - ai[i];
2215     sum = b[i];
2216     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
2217     /*    while (nz--) sum -= *v++ * x[*vi++];*/
2218     v  += nz;
2219     vi += nz;
2220     x[i] = sum;
2221   }
2222 
2223   /* backward solve the upper triangular */
2224   v   = aa + ai[n+1];
2225   vi  = aj + ai[n+1];
2226   for (i=n-1; i>=0; i--){
2227     nz = ai[2*n-i +1] - ai[2*n-i]-1;
2228     sum = x[i];
2229     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
2230     /* while (nz--) sum -= *v++ * x[*vi++]; */
2231     v   += nz;
2232     vi  += nz; vi++;
2233     x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */
2234   }
2235 
2236   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
2237   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2238   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2239   PetscFunctionReturn(0);
2240 }
2241 
2242 #undef __FUNCT__
2243 #define __FUNCT__ "MatSolve_SeqAIJ_iludt"
2244 PetscErrorCode MatSolve_SeqAIJ_iludt(Mat A,Vec bb,Vec xx)
2245 {
2246   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2247   IS                iscol = a->col,isrow = a->row;
2248   PetscErrorCode    ierr;
2249   PetscInt          i,n=A->rmap->n,*vi,*ai = a->i,*aj = a->j,*adiag=a->diag;
2250   PetscInt          nz;
2251   const PetscInt    *rout,*cout,*r,*c;
2252   PetscScalar       *x,*tmp,*tmps;
2253   const PetscScalar *b;
2254   const MatScalar   *aa = a->a,*v;
2255 
2256   PetscFunctionBegin;
2257   if (!n) PetscFunctionReturn(0);
2258 
2259   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2260   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2261   tmp  = a->solve_work;
2262 
2263   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
2264   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
2265 
2266   /* forward solve the lower triangular */
2267   tmp[0] = b[*r++];
2268   tmps   = tmp;
2269   v      = aa;
2270   vi     = aj;
2271   for (i=1; i<n; i++) {
2272     nz  = ai[i+1] - ai[i];
2273     tmp[i] = b[*r++];
2274     PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz);
2275     v += nz; vi += nz;
2276   }
2277 
2278   /* backward solve the upper triangular */
2279   v   = aa + adiag[n] + 1;
2280   vi  = aj + adiag[n] + 1;
2281   for (i=n-1; i>=0; i--){
2282     nz  = adiag[i] - adiag[i+1] - 1;
2283     PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz);
2284     x[*c--] = tmp[i] = tmp[i]*aa[adiag[i]];
2285     v += nz+1; vi += nz+1;
2286   }
2287 
2288   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
2289   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
2290   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2291   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2292   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
2293   PetscFunctionReturn(0);
2294 }
2295 
2296 #undef __FUNCT__
2297 #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2298 PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
2299 {
2300   Mat                B = *fact;
2301   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b;
2302   IS                 isicol;
2303   PetscErrorCode     ierr;
2304   const PetscInt     *r,*ic;
2305   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
2306   PetscInt           *bi,*bj,*bdiag,*bdiag_rev;
2307   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
2308   PetscInt           nlnk,*lnk;
2309   PetscBT            lnkbt;
2310   PetscTruth         row_identity,icol_identity,both_identity;
2311   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp;
2312   const PetscInt     *ics;
2313   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
2314   PetscReal          dt=info->dt,shift=info->shiftinblocks;
2315   PetscInt           nnz_max;
2316   PetscTruth         missing;
2317 
2318   PetscFunctionBegin;
2319   /* ------- symbolic factorization, can be reused ---------*/
2320   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2321   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2322   adiag=a->diag;
2323 
2324   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2325 
2326   /* bdiag is location of diagonal in factor */
2327   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2328   bdiag_rev = bdiag + n+1;
2329 
2330   /* allocate row pointers bi */
2331   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2332 
2333   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
2334   dtcount = (PetscInt)info->dtcount;
2335   if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */
2336   nnz_max  = ai[n]+2*n*dtcount+2;
2337 
2338   ierr = PetscMalloc((nnz_max+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
2339   ierr = PetscMalloc((nnz_max+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
2340 
2341   /* put together the new matrix */
2342   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2343   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
2344   b    = (Mat_SeqAIJ*)(B)->data;
2345   b->free_a       = PETSC_TRUE;
2346   b->free_ij      = PETSC_TRUE;
2347   b->singlemalloc = PETSC_FALSE;
2348   b->a          = ba;
2349   b->j          = bj;
2350   b->i          = bi;
2351   b->diag       = bdiag;
2352   b->ilen       = 0;
2353   b->imax       = 0;
2354   b->row        = isrow;
2355   b->col        = iscol;
2356   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
2357   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
2358   b->icol       = isicol;
2359   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2360 
2361   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2362   b->maxnz = nnz_max;
2363 
2364   (B)->factor                = MAT_FACTOR_ILUDT;
2365   (B)->info.factor_mallocs   = 0;
2366   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]);
2367   CHKMEMQ;
2368   /* ------- end of symbolic factorization ---------*/
2369 
2370   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2371   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
2372   ics  = ic;
2373 
2374   /* linked list for storing column indices of the active row */
2375   nlnk = n + 1;
2376   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2377 
2378   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
2379   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
2380   jtmp = im + n;
2381   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
2382   ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2383   ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
2384   vtmp = rtmp + n;
2385 
2386   bi[0]    = 0;
2387   bdiag[0] = nnz_max-1; /* location of diag[0] in factor B */
2388   bdiag_rev[n] = bdiag[0];
2389   bi[2*n+1] = bdiag[0]+1; /* endof bj and ba array */
2390   for (i=0; i<n; i++) {
2391     /* copy initial fill into linked list */
2392     nzi = 0; /* nonzeros for active row i */
2393     nzi = ai[r[i]+1] - ai[r[i]];
2394     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
2395     nzi_al = adiag[r[i]] - ai[r[i]];
2396     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
2397     ajtmp = aj + ai[r[i]];
2398     ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2399 
2400     /* load in initial (unfactored row) */
2401     aatmp = a->a + ai[r[i]];
2402     for (j=0; j<nzi; j++) {
2403       rtmp[ics[*ajtmp++]] = *aatmp++;
2404     }
2405 
2406     /* add pivot rows into linked list */
2407     row = lnk[n];
2408     while (row < i ) {
2409       nzi_bl = bi[row+1] - bi[row] + 1;
2410       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
2411       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
2412       nzi  += nlnk;
2413       row   = lnk[row];
2414     }
2415 
2416     /* copy data from lnk into jtmp, then initialize lnk */
2417     ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
2418 
2419     /* numerical factorization */
2420     bjtmp = jtmp;
2421     row   = *bjtmp++; /* 1st pivot row */
2422     while  ( row < i ) {
2423       pc         = rtmp + row;
2424       pv         = ba + bdiag[row]; /* 1./(diag of the pivot row) */
2425       multiplier = (*pc) * (*pv);
2426       *pc        = multiplier;
2427       if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */
2428         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
2429         pv         = ba + bdiag[row+1] + 1;
2430         /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */
2431         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
2432         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
2433         ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
2434       }
2435       row = *bjtmp++;
2436     }
2437 
2438     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
2439     diag_tmp = rtmp[i];  /* save diagonal value - may not needed?? */
2440     nzi_bl = 0; j = 0;
2441     while (jtmp[j] < i){ /* Note: jtmp is sorted */
2442       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2443       nzi_bl++; j++;
2444     }
2445     nzi_bu = nzi - nzi_bl -1;
2446     while (j < nzi){
2447       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2448       j++;
2449     }
2450 
2451     bjtmp = bj + bi[i];
2452     batmp = ba + bi[i];
2453     /* apply level dropping rule to L part */
2454     ncut = nzi_al + dtcount;
2455     if (ncut < nzi_bl){
2456       ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr);
2457       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
2458     } else {
2459       ncut = nzi_bl;
2460     }
2461     for (j=0; j<ncut; j++){
2462       bjtmp[j] = jtmp[j];
2463       batmp[j] = vtmp[j];
2464       /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */
2465     }
2466     bi[i+1] = bi[i] + ncut;
2467     nzi = ncut + 1;
2468 
2469     /* apply level dropping rule to U part */
2470     ncut = nzi_au + dtcount;
2471     if (ncut < nzi_bu){
2472       ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
2473       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
2474     } else {
2475       ncut = nzi_bu;
2476     }
2477     nzi += ncut;
2478 
2479     /* mark bdiagonal */
2480     bdiag[i+1]       = bdiag[i] - (ncut + 1);
2481     bdiag_rev[n-i-1] = bdiag[i+1];
2482     bi[2*n - i]      = bi[2*n - i +1] - (ncut + 1);
2483     bjtmp = bj + bdiag[i];
2484     batmp = ba + bdiag[i];
2485     *bjtmp = i;
2486     *batmp = diag_tmp; /* rtmp[i]; */
2487     if (*batmp == 0.0) {
2488       *batmp = dt+shift;
2489       /* printf(" row %d add shift %g\n",i,shift); */
2490     }
2491     *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */
2492     /* printf(" (%d,%g),",*bjtmp,*batmp); */
2493 
2494     bjtmp = bj + bdiag[i+1]+1;
2495     batmp = ba + bdiag[i+1]+1;
2496     for (k=0; k<ncut; k++){
2497       bjtmp[k] = jtmp[nzi_bl+1+k];
2498       batmp[k] = vtmp[nzi_bl+1+k];
2499       /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */
2500     }
2501     /* printf("\n"); */
2502 
2503     im[i]   = nzi; /* used by PetscLLAddSortedLU() */
2504     /*
2505     printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]);
2506     printf(" ----------------------------\n");
2507     */
2508   } /* for (i=0; i<n; i++) */
2509   /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */
2510   if (bi[n] >= bdiag[n]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[n],bdiag[n]);
2511 
2512   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2513   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2514 
2515   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2516   ierr = PetscFree(im);CHKERRQ(ierr);
2517   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2518 
2519   ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr);
2520   b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n];
2521 
2522   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
2523   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
2524   both_identity = (PetscTruth) (row_identity && icol_identity);
2525   if (row_identity && icol_identity) {
2526     B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt;
2527   } else {
2528     B->ops->solve = MatSolve_SeqAIJ_iludt;
2529   }
2530 
2531   B->ops->lufactorsymbolic  = MatILUDTFactorSymbolic_SeqAIJ;
2532   B->ops->lufactornumeric   = MatILUDTFactorNumeric_SeqAIJ;
2533   B->ops->solveadd          = 0;
2534   B->ops->solvetranspose    = 0;
2535   B->ops->solvetransposeadd = 0;
2536   B->ops->matsolve          = 0;
2537   B->assembled              = PETSC_TRUE;
2538   B->preallocated           = PETSC_TRUE;
2539   PetscFunctionReturn(0);
2540 }
2541 
2542 /* a wraper of MatILUDTFactor_SeqAIJ() */
2543 #undef __FUNCT__
2544 #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ"
2545 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
2546 {
2547   PetscErrorCode     ierr;
2548 
2549   PetscFunctionBegin;
2550   ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr);
2551 
2552   fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ;
2553   PetscFunctionReturn(0);
2554 }
2555 
2556 /*
2557    same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors
2558    - intend to replace existing MatLUFactorNumeric_SeqAIJ()
2559 */
2560 #undef __FUNCT__
2561 #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ"
2562 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info)
2563 {
2564   Mat            C=fact;
2565   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
2566   IS             isrow = b->row,isicol = b->icol;
2567   PetscErrorCode ierr;
2568   const PetscInt *r,*ic,*ics;
2569   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2570   PetscInt       *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj;
2571   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
2572   PetscReal      dt=info->dt,shift=info->shiftinblocks;
2573   PetscTruth     row_identity, col_identity;
2574 
2575   PetscFunctionBegin;
2576   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2577   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
2578   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2579   ics  = ic;
2580 
2581   for (i=0; i<n; i++){
2582     /* initialize rtmp array */
2583     nzl   = bi[i+1] - bi[i];       /* num of nozeros in L(i,:) */
2584     bjtmp = bj + bi[i];
2585     for  (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0;
2586     rtmp[i] = 0.0;
2587     nzu   = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */
2588     bjtmp = bj + bdiag[i+1] + 1;
2589     for  (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0;
2590 
2591     /* load in initial unfactored row of A */
2592     /* printf("row %d\n",i); */
2593     nz    = ai[r[i]+1] - ai[r[i]];
2594     ajtmp = aj + ai[r[i]];
2595     v     = aa + ai[r[i]];
2596     for (j=0; j<nz; j++) {
2597       rtmp[ics[*ajtmp++]] = v[j];
2598       /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */
2599     }
2600     /* printf("\n"); */
2601 
2602     /* numerical factorization */
2603     bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */
2604     nzl   = bi[i+1] - bi[i]; /* num of entries in L(i,:) */
2605     k = 0;
2606     while (k < nzl){
2607       row   = *bjtmp++;
2608       /* printf("  prow %d\n",row); */
2609       pc         = rtmp + row;
2610       pv         = b->a + bdiag[row]; /* 1./(diag of the pivot row) */
2611       multiplier = (*pc) * (*pv);
2612       *pc        = multiplier;
2613       if (PetscAbsScalar(multiplier) > dt){
2614         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
2615         pv         = b->a + bdiag[row+1] + 1;
2616         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
2617         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
2618         /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */
2619       }
2620       k++;
2621     }
2622 
2623     /* finished row so stick it into b->a */
2624     /* L-part */
2625     pv = b->a + bi[i] ;
2626     pj = bj + bi[i] ;
2627     nzl = bi[i+1] - bi[i];
2628     for (j=0; j<nzl; j++) {
2629       pv[j] = rtmp[pj[j]];
2630       /* printf(" (%d,%g),",pj[j],pv[j]); */
2631     }
2632 
2633     /* diagonal: invert diagonal entries for simplier triangular solves */
2634     if (rtmp[i] == 0.0) rtmp[i] = dt+shift;
2635     b->a[bdiag[i]] = 1.0/rtmp[i];
2636     /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */
2637 
2638     /* U-part */
2639     pv = b->a + bdiag[i+1] + 1;
2640     pj = bj + bdiag[i+1] + 1;
2641     nzu = bdiag[i] - bdiag[i+1] - 1;
2642     for (j=0; j<nzu; j++) {
2643       pv[j] = rtmp[pj[j]];
2644       /* printf(" (%d,%g),",pj[j],pv[j]); */
2645     }
2646     /* printf("\n"); */
2647   }
2648 
2649   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2650   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2651   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2652 
2653   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
2654   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
2655   if (row_identity && col_identity) {
2656     C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering_iludt;
2657   } else {
2658     C->ops->solve   = MatSolve_SeqAIJ_iludt;
2659   }
2660   C->ops->solveadd           = 0;
2661   C->ops->solvetranspose     = 0;
2662   C->ops->solvetransposeadd  = 0;
2663   C->ops->matsolve           = 0;
2664   C->assembled    = PETSC_TRUE;
2665   C->preallocated = PETSC_TRUE;
2666   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
2667   PetscFunctionReturn(0);
2668 }
2669