xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision cee9d6f295ef30f06a70cc94141c40a93dc7ea43)
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__ "MatILUFactorSymbolic_SeqAIJ"
1228 PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1229 {
1230   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1231   IS                 isicol;
1232   PetscErrorCode     ierr;
1233   const PetscInt     *r,*ic;
1234   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
1235   PetscInt           *bi,*cols,nnz,*cols_lvl;
1236   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1237   PetscInt           i,levels,diagonal_fill;
1238   PetscTruth         col_identity,row_identity;
1239   PetscReal          f;
1240   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1241   PetscBT            lnkbt;
1242   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1243   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1244   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1245   PetscTruth         missing;
1246 
1247   PetscFunctionBegin;
1248   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);
1249   f             = info->fill;
1250   levels        = (PetscInt)info->levels;
1251   diagonal_fill = (PetscInt)info->diagonal_fill;
1252   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1253 
1254   /* special case that simply copies fill pattern */
1255   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1256   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
1257   if (!levels && row_identity && col_identity) { /* ilu(0) with natural ordering */
1258 
1259     PetscTruth newdatastruct=PETSC_FALSE;
1260     ierr = PetscOptionsGetTruth(PETSC_NULL,"-ilu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
1261     if (newdatastruct){
1262       ierr = MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1263       (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ_newdatastruct;
1264     } else {
1265       ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
1266       (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1267     }
1268 
1269     fact->factor = MAT_FACTOR_ILU;
1270     (fact)->info.factor_mallocs    = 0;
1271     (fact)->info.fill_ratio_given  = info->fill;
1272     (fact)->info.fill_ratio_needed = 1.0;
1273     b               = (Mat_SeqAIJ*)(fact)->data;
1274     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1275     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1276     b->row              = isrow;
1277     b->col              = iscol;
1278     b->icol             = isicol;
1279     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1280     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1281     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1282     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1283     PetscFunctionReturn(0);
1284   }
1285 
1286   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1287   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1288 
1289   /* get new row pointers */
1290   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1291   bi[0] = 0;
1292   /* bdiag is location of diagonal in factor */
1293   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1294   bdiag[0]  = 0;
1295 
1296   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
1297   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
1298 
1299   /* create a linked list for storing column indices of the active row */
1300   nlnk = n + 1;
1301   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1302 
1303   /* initial FreeSpace size is f*(ai[n]+1) */
1304   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
1305   current_space = free_space;
1306   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
1307   current_space_lvl = free_space_lvl;
1308 
1309   for (i=0; i<n; i++) {
1310     nzi = 0;
1311     /* copy current row into linked list */
1312     nnz  = ai[r[i]+1] - ai[r[i]];
1313     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1314     cols = aj + ai[r[i]];
1315     lnk[i] = -1; /* marker to indicate if diagonal exists */
1316     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1317     nzi += nlnk;
1318 
1319     /* make sure diagonal entry is included */
1320     if (diagonal_fill && lnk[i] == -1) {
1321       fm = n;
1322       while (lnk[fm] < i) fm = lnk[fm];
1323       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1324       lnk[fm]    = i;
1325       lnk_lvl[i] = 0;
1326       nzi++; dcount++;
1327     }
1328 
1329     /* add pivot rows into the active row */
1330     nzbd = 0;
1331     prow = lnk[n];
1332     while (prow < i) {
1333       nnz      = bdiag[prow];
1334       cols     = bj_ptr[prow] + nnz + 1;
1335       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1336       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1337       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
1338       nzi += nlnk;
1339       prow = lnk[prow];
1340       nzbd++;
1341     }
1342     bdiag[i] = nzbd;
1343     bi[i+1]  = bi[i] + nzi;
1344 
1345     /* if free space is not available, make more free space */
1346     if (current_space->local_remaining<nzi) {
1347       nnz = nzi*(n - i); /* estimated and max additional space needed */
1348       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1349       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
1350       reallocs++;
1351     }
1352 
1353     /* copy data into free_space and free_space_lvl, then initialize lnk */
1354     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1355     bj_ptr[i]    = current_space->array;
1356     bjlvl_ptr[i] = current_space_lvl->array;
1357 
1358     /* make sure the active row i has diagonal entry */
1359     if (*(bj_ptr[i]+bdiag[i]) != i) {
1360       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1361     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
1362     }
1363 
1364     current_space->array           += nzi;
1365     current_space->local_used      += nzi;
1366     current_space->local_remaining -= nzi;
1367     current_space_lvl->array           += nzi;
1368     current_space_lvl->local_used      += nzi;
1369     current_space_lvl->local_remaining -= nzi;
1370   }
1371 
1372   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1373   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1374 
1375   /* destroy list of free space and other temporary arrays */
1376   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1377   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1378   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1379   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1380   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
1381 
1382 #if defined(PETSC_USE_INFO)
1383   {
1384     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1385     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1386     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1387     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1388     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1389     if (diagonal_fill) {
1390       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1391     }
1392   }
1393 #endif
1394 
1395   /* put together the new matrix */
1396   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1397   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1398   b = (Mat_SeqAIJ*)(fact)->data;
1399   b->free_a       = PETSC_TRUE;
1400   b->free_ij      = PETSC_TRUE;
1401   b->singlemalloc = PETSC_FALSE;
1402   ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
1403   b->j          = bj;
1404   b->i          = bi;
1405   for (i=0; i<n; i++) bdiag[i] += bi[i];
1406   b->diag       = bdiag;
1407   b->ilen       = 0;
1408   b->imax       = 0;
1409   b->row        = isrow;
1410   b->col        = iscol;
1411   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1412   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1413   b->icol       = isicol;
1414   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1415   /* In b structure:  Free imax, ilen, old a, old j.
1416      Allocate bdiag, solve_work, new a, new j */
1417   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1418   b->maxnz             = b->nz = bi[n] ;
1419   (fact)->info.factor_mallocs    = reallocs;
1420   (fact)->info.fill_ratio_given  = f;
1421   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1422   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1423   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1424   PetscFunctionReturn(0);
1425 }
1426 
1427 #include "../src/mat/impls/sbaij/seq/sbaij.h"
1428 #undef __FUNCT__
1429 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1430 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1431 {
1432   Mat            C = B;
1433   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1434   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
1435   IS             ip=b->row,iip = b->icol;
1436   PetscErrorCode ierr;
1437   const PetscInt *rip,*riip;
1438   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
1439   PetscInt       *ai=a->i,*aj=a->j;
1440   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1441   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1442   PetscReal      zeropivot,rs,shiftnz;
1443   PetscReal      shiftpd;
1444   ChShift_Ctx    sctx;
1445   PetscInt       newshift;
1446   PetscTruth     perm_identity;
1447 
1448   PetscFunctionBegin;
1449 
1450   shiftnz   = info->shiftnz;
1451   shiftpd   = info->shiftpd;
1452   zeropivot = info->zeropivot;
1453 
1454   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1455   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
1456 
1457   /* initialization */
1458   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1459   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1460   jl   = il + mbs;
1461   rtmp = (MatScalar*)(jl + mbs);
1462 
1463   sctx.shift_amount = 0;
1464   sctx.nshift       = 0;
1465   do {
1466     sctx.chshift = PETSC_FALSE;
1467     for (i=0; i<mbs; i++) {
1468       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1469     }
1470 
1471     for (k = 0; k<mbs; k++){
1472       bval = ba + bi[k];
1473       /* initialize k-th row by the perm[k]-th row of A */
1474       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1475       for (j = jmin; j < jmax; j++){
1476         col = riip[aj[j]];
1477         if (col >= k){ /* only take upper triangular entry */
1478           rtmp[col] = aa[j];
1479           *bval++  = 0.0; /* for in-place factorization */
1480         }
1481       }
1482       /* shift the diagonal of the matrix */
1483       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1484 
1485       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1486       dk = rtmp[k];
1487       i = jl[k]; /* first row to be added to k_th row  */
1488 
1489       while (i < k){
1490         nexti = jl[i]; /* next row to be added to k_th row */
1491 
1492         /* compute multiplier, update diag(k) and U(i,k) */
1493         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1494         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1495         dk += uikdi*ba[ili];
1496         ba[ili] = uikdi; /* -U(i,k) */
1497 
1498         /* add multiple of row i to k-th row */
1499         jmin = ili + 1; jmax = bi[i+1];
1500         if (jmin < jmax){
1501           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1502           /* update il and jl for row i */
1503           il[i] = jmin;
1504           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1505         }
1506         i = nexti;
1507       }
1508 
1509       /* shift the diagonals when zero pivot is detected */
1510       /* compute rs=sum of abs(off-diagonal) */
1511       rs   = 0.0;
1512       jmin = bi[k]+1;
1513       nz   = bi[k+1] - jmin;
1514       bcol = bj + jmin;
1515       while (nz--){
1516         rs += PetscAbsScalar(rtmp[*bcol]);
1517         bcol++;
1518       }
1519 
1520       sctx.rs = rs;
1521       sctx.pv = dk;
1522       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1523 
1524       if (newshift == 1) {
1525         if (!sctx.shift_amount) {
1526           sctx.shift_amount = 1e-5;
1527         }
1528         break;
1529       }
1530 
1531       /* copy data into U(k,:) */
1532       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1533       jmin = bi[k]+1; jmax = bi[k+1];
1534       if (jmin < jmax) {
1535         for (j=jmin; j<jmax; j++){
1536           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1537         }
1538         /* add the k-th row into il and jl */
1539         il[k] = jmin;
1540         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1541       }
1542     }
1543   } while (sctx.chshift);
1544   ierr = PetscFree(il);CHKERRQ(ierr);
1545 
1546   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1547   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1548 
1549   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1550   if (perm_identity){
1551     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1552     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1553     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1554     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1555   } else {
1556     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1557     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1558     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1559     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1560   }
1561 
1562   C->assembled    = PETSC_TRUE;
1563   C->preallocated = PETSC_TRUE;
1564   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1565   if (sctx.nshift){
1566     if (shiftnz) {
1567       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1568     } else if (shiftpd) {
1569       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1570     }
1571   }
1572   PetscFunctionReturn(0);
1573 }
1574 
1575 #undef __FUNCT__
1576 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1577 PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1578 {
1579   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1580   Mat_SeqSBAIJ       *b;
1581   PetscErrorCode     ierr;
1582   PetscTruth         perm_identity,missing;
1583   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
1584   const PetscInt     *rip,*riip;
1585   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1586   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
1587   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1588   PetscReal          fill=info->fill,levels=info->levels;
1589   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1590   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1591   PetscBT            lnkbt;
1592   IS                 iperm;
1593 
1594   PetscFunctionBegin;
1595   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);
1596   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1597   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1598   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1599   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
1600 
1601   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1602   ui[0] = 0;
1603 
1604   /* ICC(0) without matrix ordering: simply copies fill pattern */
1605   if (!levels && perm_identity) {
1606 
1607     for (i=0; i<am; i++) {
1608       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1609     }
1610     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1611     cols = uj;
1612     for (i=0; i<am; i++) {
1613       aj    = a->j + a->diag[i];
1614       ncols = ui[i+1] - ui[i];
1615       for (j=0; j<ncols; j++) *cols++ = *aj++;
1616     }
1617   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1618     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1619     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1620 
1621     /* initialization */
1622     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
1623 
1624     /* jl: linked list for storing indices of the pivot rows
1625        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1626     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
1627     il         = jl + am;
1628     uj_ptr     = (PetscInt**)(il + am);
1629     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
1630     for (i=0; i<am; i++){
1631       jl[i] = am; il[i] = 0;
1632     }
1633 
1634     /* create and initialize a linked list for storing column indices of the active row k */
1635     nlnk = am + 1;
1636     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1637 
1638     /* initial FreeSpace size is fill*(ai[am]+1) */
1639     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
1640     current_space = free_space;
1641     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
1642     current_space_lvl = free_space_lvl;
1643 
1644     for (k=0; k<am; k++){  /* for each active row k */
1645       /* initialize lnk by the column indices of row rip[k] of A */
1646       nzk   = 0;
1647       ncols = ai[rip[k]+1] - ai[rip[k]];
1648       if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
1649       ncols_upper = 0;
1650       for (j=0; j<ncols; j++){
1651         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1652         if (riip[i] >= k){ /* only take upper triangular entry */
1653           ajtmp[ncols_upper] = i;
1654           ncols_upper++;
1655         }
1656       }
1657       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1658       nzk += nlnk;
1659 
1660       /* update lnk by computing fill-in for each pivot row to be merged in */
1661       prow = jl[k]; /* 1st pivot row */
1662 
1663       while (prow < k){
1664         nextprow = jl[prow];
1665 
1666         /* merge prow into k-th row */
1667         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
1668         jmax = ui[prow+1];
1669         ncols = jmax-jmin;
1670         i     = jmin - ui[prow];
1671         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1672         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
1673         j     = *(uj - 1);
1674         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
1675         nzk += nlnk;
1676 
1677         /* update il and jl for prow */
1678         if (jmin < jmax){
1679           il[prow] = jmin;
1680           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
1681         }
1682         prow = nextprow;
1683       }
1684 
1685       /* if free space is not available, make more free space */
1686       if (current_space->local_remaining<nzk) {
1687         i = am - k + 1; /* num of unfactored rows */
1688         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1689         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1690         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
1691         reallocs++;
1692       }
1693 
1694       /* copy data into free_space and free_space_lvl, then initialize lnk */
1695       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
1696       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1697 
1698       /* add the k-th row into il and jl */
1699       if (nzk > 1){
1700         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
1701         jl[k] = jl[i]; jl[i] = k;
1702         il[k] = ui[k] + 1;
1703       }
1704       uj_ptr[k]     = current_space->array;
1705       uj_lvl_ptr[k] = current_space_lvl->array;
1706 
1707       current_space->array           += nzk;
1708       current_space->local_used      += nzk;
1709       current_space->local_remaining -= nzk;
1710 
1711       current_space_lvl->array           += nzk;
1712       current_space_lvl->local_used      += nzk;
1713       current_space_lvl->local_remaining -= nzk;
1714 
1715       ui[k+1] = ui[k] + nzk;
1716     }
1717 
1718 #if defined(PETSC_USE_INFO)
1719     if (ai[am] != 0) {
1720       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1721       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1722       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1723       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
1724     } else {
1725       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1726     }
1727 #endif
1728 
1729     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1730     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
1731     ierr = PetscFree(jl);CHKERRQ(ierr);
1732     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
1733 
1734     /* destroy list of free space and other temporary array(s) */
1735     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1736     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1737     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1738     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1739 
1740   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
1741 
1742   /* put together the new matrix in MATSEQSBAIJ format */
1743 
1744   b    = (Mat_SeqSBAIJ*)(fact)->data;
1745   b->singlemalloc = PETSC_FALSE;
1746   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
1747   b->j    = uj;
1748   b->i    = ui;
1749   b->diag = 0;
1750   b->ilen = 0;
1751   b->imax = 0;
1752   b->row  = perm;
1753   b->col  = perm;
1754   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1755   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1756   b->icol = iperm;
1757   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1758   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1759   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1760   b->maxnz   = b->nz = ui[am];
1761   b->free_a  = PETSC_TRUE;
1762   b->free_ij = PETSC_TRUE;
1763 
1764   (fact)->info.factor_mallocs    = reallocs;
1765   (fact)->info.fill_ratio_given  = fill;
1766   if (ai[am] != 0) {
1767     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
1768   } else {
1769     (fact)->info.fill_ratio_needed = 0.0;
1770   }
1771   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1772   PetscFunctionReturn(0);
1773 }
1774 
1775 #undef __FUNCT__
1776 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1777 PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1778 {
1779   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1780   Mat_SeqSBAIJ       *b;
1781   PetscErrorCode     ierr;
1782   PetscTruth         perm_identity;
1783   PetscReal          fill = info->fill;
1784   const PetscInt     *rip,*riip;
1785   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
1786   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1787   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1788   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1789   PetscBT            lnkbt;
1790   IS                 iperm;
1791 
1792   PetscFunctionBegin;
1793   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);
1794   /* check whether perm is the identity mapping */
1795   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1796   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
1797   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1798   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1799 
1800   /* initialization */
1801   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1802   ui[0] = 0;
1803 
1804   /* jl: linked list for storing indices of the pivot rows
1805      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1806   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
1807   il     = jl + am;
1808   cols   = il + am;
1809   ui_ptr = (PetscInt**)(cols + am);
1810   for (i=0; i<am; i++){
1811     jl[i] = am; il[i] = 0;
1812   }
1813 
1814   /* create and initialize a linked list for storing column indices of the active row k */
1815   nlnk = am + 1;
1816   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1817 
1818   /* initial FreeSpace size is fill*(ai[am]+1) */
1819   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
1820   current_space = free_space;
1821 
1822   for (k=0; k<am; k++){  /* for each active row k */
1823     /* initialize lnk by the column indices of row rip[k] of A */
1824     nzk   = 0;
1825     ncols = ai[rip[k]+1] - ai[rip[k]];
1826     if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
1827     ncols_upper = 0;
1828     for (j=0; j<ncols; j++){
1829       i = riip[*(aj + ai[rip[k]] + j)];
1830       if (i >= k){ /* only take upper triangular entry */
1831         cols[ncols_upper] = i;
1832         ncols_upper++;
1833       }
1834     }
1835     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1836     nzk += nlnk;
1837 
1838     /* update lnk by computing fill-in for each pivot row to be merged in */
1839     prow = jl[k]; /* 1st pivot row */
1840 
1841     while (prow < k){
1842       nextprow = jl[prow];
1843       /* merge prow into k-th row */
1844       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
1845       jmax = ui[prow+1];
1846       ncols = jmax-jmin;
1847       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1848       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1849       nzk += nlnk;
1850 
1851       /* update il and jl for prow */
1852       if (jmin < jmax){
1853         il[prow] = jmin;
1854         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
1855       }
1856       prow = nextprow;
1857     }
1858 
1859     /* if free space is not available, make more free space */
1860     if (current_space->local_remaining<nzk) {
1861       i = am - k + 1; /* num of unfactored rows */
1862       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1863       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1864       reallocs++;
1865     }
1866 
1867     /* copy data into free space, then initialize lnk */
1868     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1869 
1870     /* add the k-th row into il and jl */
1871     if (nzk-1 > 0){
1872       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
1873       jl[k] = jl[i]; jl[i] = k;
1874       il[k] = ui[k] + 1;
1875     }
1876     ui_ptr[k] = current_space->array;
1877     current_space->array           += nzk;
1878     current_space->local_used      += nzk;
1879     current_space->local_remaining -= nzk;
1880 
1881     ui[k+1] = ui[k] + nzk;
1882   }
1883 
1884 #if defined(PETSC_USE_INFO)
1885   if (ai[am] != 0) {
1886     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1887     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1888     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1889     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
1890   } else {
1891      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1892   }
1893 #endif
1894 
1895   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1896   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
1897   ierr = PetscFree(jl);CHKERRQ(ierr);
1898 
1899   /* destroy list of free space and other temporary array(s) */
1900   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1901   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1902   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1903 
1904   /* put together the new matrix in MATSEQSBAIJ format */
1905 
1906   b = (Mat_SeqSBAIJ*)(fact)->data;
1907   b->singlemalloc = PETSC_FALSE;
1908   b->free_a       = PETSC_TRUE;
1909   b->free_ij      = PETSC_TRUE;
1910   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
1911   b->j    = uj;
1912   b->i    = ui;
1913   b->diag = 0;
1914   b->ilen = 0;
1915   b->imax = 0;
1916   b->row  = perm;
1917   b->col  = perm;
1918   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1919   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1920   b->icol = iperm;
1921   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1922   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1923   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1924   b->maxnz = b->nz = ui[am];
1925 
1926   (fact)->info.factor_mallocs    = reallocs;
1927   (fact)->info.fill_ratio_given  = fill;
1928   if (ai[am] != 0) {
1929     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
1930   } else {
1931     (fact)->info.fill_ratio_needed = 0.0;
1932   }
1933   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1934   PetscFunctionReturn(0);
1935 }
1936 
1937 #undef __FUNCT__
1938 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_iludt"
1939 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat A,Vec bb,Vec xx)
1940 {
1941   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1942   PetscErrorCode    ierr;
1943   PetscInt          n = A->rmap->n;
1944   const PetscInt    *ai = a->i,*aj = a->j,*vi;
1945   PetscScalar       *x,sum;
1946   const PetscScalar *b;
1947   const MatScalar   *aa = a->a,*v;
1948   PetscInt          i,nz;
1949 
1950   PetscFunctionBegin;
1951   if (!n) PetscFunctionReturn(0);
1952 
1953   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
1954   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1955 
1956   /* forward solve the lower triangular */
1957   x[0] = b[0];
1958   v    = aa;
1959   vi   = aj;
1960   for (i=1; i<n; i++) {
1961     nz  = ai[i+1] - ai[i];
1962     sum = b[i];
1963     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1964     /*    while (nz--) sum -= *v++ * x[*vi++];*/
1965     v  += nz;
1966     vi += nz;
1967     x[i] = sum;
1968   }
1969 
1970   /* backward solve the upper triangular */
1971   v   = aa + ai[n+1];
1972   vi  = aj + ai[n+1];
1973   for (i=n-1; i>=0; i--){
1974     nz = ai[2*n-i +1] - ai[2*n-i]-1;
1975     sum = x[i];
1976     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1977     /* while (nz--) sum -= *v++ * x[*vi++]; */
1978     v   += nz;
1979     vi  += nz; vi++;
1980     x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */
1981   }
1982 
1983   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
1984   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
1985   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1986   PetscFunctionReturn(0);
1987 }
1988 
1989 #undef __FUNCT__
1990 #define __FUNCT__ "MatSolve_SeqAIJ_iludt"
1991 PetscErrorCode MatSolve_SeqAIJ_iludt(Mat A,Vec bb,Vec xx)
1992 {
1993   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1994   IS                iscol = a->col,isrow = a->row;
1995   PetscErrorCode    ierr;
1996   PetscInt          i,n=A->rmap->n,*vi,*ai = a->i,*aj = a->j,*adiag=a->diag;
1997   PetscInt          nz;
1998   const PetscInt    *rout,*cout,*r,*c;
1999   PetscScalar       *x,*tmp,*tmps;
2000   const PetscScalar *b;
2001   const MatScalar   *aa = a->a,*v;
2002 
2003   PetscFunctionBegin;
2004   if (!n) PetscFunctionReturn(0);
2005 
2006   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2007   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2008   tmp  = a->solve_work;
2009 
2010   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
2011   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
2012 
2013   /* forward solve the lower triangular */
2014   tmp[0] = b[*r++];
2015   tmps   = tmp;
2016   v      = aa;
2017   vi     = aj;
2018   for (i=1; i<n; i++) {
2019     nz  = ai[i+1] - ai[i];
2020     tmp[i] = b[*r++];
2021     PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz);
2022     v += nz; vi += nz;
2023   }
2024 
2025   /* backward solve the upper triangular */
2026   v   = aa + adiag[n] + 1;
2027   vi  = aj + adiag[n] + 1;
2028   for (i=n-1; i>=0; i--){
2029     nz  = adiag[i] - adiag[i+1] - 1;
2030     PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz);
2031     x[*c--] = tmp[i] = tmp[i]*aa[adiag[i]];
2032     v += nz+1; vi += nz+1;
2033   }
2034 
2035   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
2036   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
2037   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2038   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2039   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
2040   PetscFunctionReturn(0);
2041 }
2042 
2043 #undef __FUNCT__
2044 #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2045 PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
2046 {
2047   Mat                B = *fact;
2048   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b;
2049   IS                 isicol;
2050   PetscErrorCode     ierr;
2051   const PetscInt     *r,*ic;
2052   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
2053   PetscInt           *bi,*bj,*bdiag,*bdiag_rev;
2054   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
2055   PetscInt           nlnk,*lnk;
2056   PetscBT            lnkbt;
2057   PetscTruth         row_identity,icol_identity,both_identity;
2058   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp;
2059   const PetscInt     *ics;
2060   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
2061   PetscReal          dt=info->dt,shift=info->shiftinblocks;
2062   PetscInt           nnz_max;
2063   PetscTruth         missing;
2064 
2065   PetscFunctionBegin;
2066   /* ------- symbolic factorization, can be reused ---------*/
2067   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2068   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2069   adiag=a->diag;
2070 
2071   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2072 
2073   /* bdiag is location of diagonal in factor */
2074   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2075   bdiag_rev = bdiag + n+1;
2076 
2077   /* allocate row pointers bi */
2078   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2079 
2080   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
2081   dtcount = (PetscInt)info->dtcount;
2082   if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */
2083   nnz_max  = ai[n]+2*n*dtcount+2;
2084 
2085   ierr = PetscMalloc((nnz_max+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
2086   ierr = PetscMalloc((nnz_max+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
2087 
2088   /* put together the new matrix */
2089   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2090   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
2091   b    = (Mat_SeqAIJ*)(B)->data;
2092   b->free_a       = PETSC_TRUE;
2093   b->free_ij      = PETSC_TRUE;
2094   b->singlemalloc = PETSC_FALSE;
2095   b->a          = ba;
2096   b->j          = bj;
2097   b->i          = bi;
2098   b->diag       = bdiag;
2099   b->ilen       = 0;
2100   b->imax       = 0;
2101   b->row        = isrow;
2102   b->col        = iscol;
2103   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
2104   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
2105   b->icol       = isicol;
2106   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2107 
2108   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2109   b->maxnz = nnz_max;
2110 
2111   (B)->factor                = MAT_FACTOR_ILUDT;
2112   (B)->info.factor_mallocs   = 0;
2113   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]);
2114   CHKMEMQ;
2115   /* ------- end of symbolic factorization ---------*/
2116 
2117   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2118   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
2119   ics  = ic;
2120 
2121   /* linked list for storing column indices of the active row */
2122   nlnk = n + 1;
2123   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2124 
2125   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
2126   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
2127   jtmp = im + n;
2128   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
2129   ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2130   ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
2131   vtmp = rtmp + n;
2132 
2133   bi[0]    = 0;
2134   bdiag[0] = nnz_max-1; /* location of diag[0] in factor B */
2135   bdiag_rev[n] = bdiag[0];
2136   bi[2*n+1] = bdiag[0]+1; /* endof bj and ba array */
2137   for (i=0; i<n; i++) {
2138     /* copy initial fill into linked list */
2139     nzi = 0; /* nonzeros for active row i */
2140     nzi = ai[r[i]+1] - ai[r[i]];
2141     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
2142     nzi_al = adiag[r[i]] - ai[r[i]];
2143     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
2144     ajtmp = aj + ai[r[i]];
2145     ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2146 
2147     /* load in initial (unfactored row) */
2148     aatmp = a->a + ai[r[i]];
2149     for (j=0; j<nzi; j++) {
2150       rtmp[ics[*ajtmp++]] = *aatmp++;
2151     }
2152 
2153     /* add pivot rows into linked list */
2154     row = lnk[n];
2155     while (row < i ) {
2156       nzi_bl = bi[row+1] - bi[row] + 1;
2157       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
2158       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
2159       nzi  += nlnk;
2160       row   = lnk[row];
2161     }
2162 
2163     /* copy data from lnk into jtmp, then initialize lnk */
2164     ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
2165 
2166     /* numerical factorization */
2167     bjtmp = jtmp;
2168     row   = *bjtmp++; /* 1st pivot row */
2169     while  ( row < i ) {
2170       pc         = rtmp + row;
2171       pv         = ba + bdiag[row]; /* 1./(diag of the pivot row) */
2172       multiplier = (*pc) * (*pv);
2173       *pc        = multiplier;
2174       if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */
2175         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
2176         pv         = ba + bdiag[row+1] + 1;
2177         /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */
2178         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
2179         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
2180         ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
2181       }
2182       row = *bjtmp++;
2183     }
2184 
2185     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
2186     diag_tmp = rtmp[i];  /* save diagonal value - may not needed?? */
2187     nzi_bl = 0; j = 0;
2188     while (jtmp[j] < i){ /* Note: jtmp is sorted */
2189       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2190       nzi_bl++; j++;
2191     }
2192     nzi_bu = nzi - nzi_bl -1;
2193     while (j < nzi){
2194       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2195       j++;
2196     }
2197 
2198     bjtmp = bj + bi[i];
2199     batmp = ba + bi[i];
2200     /* apply level dropping rule to L part */
2201     ncut = nzi_al + dtcount;
2202     if (ncut < nzi_bl){
2203       ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr);
2204       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
2205     } else {
2206       ncut = nzi_bl;
2207     }
2208     for (j=0; j<ncut; j++){
2209       bjtmp[j] = jtmp[j];
2210       batmp[j] = vtmp[j];
2211       /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */
2212     }
2213     bi[i+1] = bi[i] + ncut;
2214     nzi = ncut + 1;
2215 
2216     /* apply level dropping rule to U part */
2217     ncut = nzi_au + dtcount;
2218     if (ncut < nzi_bu){
2219       ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
2220       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
2221     } else {
2222       ncut = nzi_bu;
2223     }
2224     nzi += ncut;
2225 
2226     /* mark bdiagonal */
2227     bdiag[i+1]       = bdiag[i] - (ncut + 1);
2228     bdiag_rev[n-i-1] = bdiag[i+1];
2229     bi[2*n - i]      = bi[2*n - i +1] - (ncut + 1);
2230     bjtmp = bj + bdiag[i];
2231     batmp = ba + bdiag[i];
2232     *bjtmp = i;
2233     *batmp = diag_tmp; /* rtmp[i]; */
2234     if (*batmp == 0.0) {
2235       *batmp = dt+shift;
2236       /* printf(" row %d add shift %g\n",i,shift); */
2237     }
2238     *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */
2239     /* printf(" (%d,%g),",*bjtmp,*batmp); */
2240 
2241     bjtmp = bj + bdiag[i+1]+1;
2242     batmp = ba + bdiag[i+1]+1;
2243     for (k=0; k<ncut; k++){
2244       bjtmp[k] = jtmp[nzi_bl+1+k];
2245       batmp[k] = vtmp[nzi_bl+1+k];
2246       /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */
2247     }
2248     /* printf("\n"); */
2249 
2250     im[i]   = nzi; /* used by PetscLLAddSortedLU() */
2251     /*
2252     printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]);
2253     printf(" ----------------------------\n");
2254     */
2255   } /* for (i=0; i<n; i++) */
2256   /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */
2257   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]);
2258 
2259   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2260   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2261 
2262   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2263   ierr = PetscFree(im);CHKERRQ(ierr);
2264   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2265 
2266   ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr);
2267   b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n];
2268 
2269   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
2270   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
2271   both_identity = (PetscTruth) (row_identity && icol_identity);
2272   if (row_identity && icol_identity) {
2273     B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt;
2274   } else {
2275     B->ops->solve = MatSolve_SeqAIJ_iludt;
2276   }
2277 
2278   B->ops->lufactorsymbolic  = MatILUDTFactorSymbolic_SeqAIJ;
2279   B->ops->lufactornumeric   = MatILUDTFactorNumeric_SeqAIJ;
2280   B->ops->solveadd          = 0;
2281   B->ops->solvetranspose    = 0;
2282   B->ops->solvetransposeadd = 0;
2283   B->ops->matsolve          = 0;
2284   B->assembled              = PETSC_TRUE;
2285   B->preallocated           = PETSC_TRUE;
2286   PetscFunctionReturn(0);
2287 }
2288 
2289 /* a wraper of MatILUDTFactor_SeqAIJ() */
2290 #undef __FUNCT__
2291 #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ"
2292 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
2293 {
2294   PetscErrorCode     ierr;
2295 
2296   PetscFunctionBegin;
2297   ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr);
2298 
2299   fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ;
2300   PetscFunctionReturn(0);
2301 }
2302 
2303 /*
2304    same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors
2305    - intend to replace existing MatLUFactorNumeric_SeqAIJ()
2306 */
2307 #undef __FUNCT__
2308 #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ"
2309 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info)
2310 {
2311   Mat            C=fact;
2312   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
2313   IS             isrow = b->row,isicol = b->icol;
2314   PetscErrorCode ierr;
2315   const PetscInt *r,*ic,*ics;
2316   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2317   PetscInt       *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj;
2318   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
2319   PetscReal      dt=info->dt,shift=info->shiftinblocks;
2320   PetscTruth     row_identity, col_identity;
2321 
2322   PetscFunctionBegin;
2323   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2324   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
2325   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2326   ics  = ic;
2327 
2328   for (i=0; i<n; i++){
2329     /* initialize rtmp array */
2330     nzl   = bi[i+1] - bi[i];       /* num of nozeros in L(i,:) */
2331     bjtmp = bj + bi[i];
2332     for  (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0;
2333     rtmp[i] = 0.0;
2334     nzu   = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */
2335     bjtmp = bj + bdiag[i+1] + 1;
2336     for  (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0;
2337 
2338     /* load in initial unfactored row of A */
2339     /* printf("row %d\n",i); */
2340     nz    = ai[r[i]+1] - ai[r[i]];
2341     ajtmp = aj + ai[r[i]];
2342     v     = aa + ai[r[i]];
2343     for (j=0; j<nz; j++) {
2344       rtmp[ics[*ajtmp++]] = v[j];
2345       /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */
2346     }
2347     /* printf("\n"); */
2348 
2349     /* numerical factorization */
2350     bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */
2351     nzl   = bi[i+1] - bi[i]; /* num of entries in L(i,:) */
2352     k = 0;
2353     while (k < nzl){
2354       row   = *bjtmp++;
2355       /* printf("  prow %d\n",row); */
2356       pc         = rtmp + row;
2357       pv         = b->a + bdiag[row]; /* 1./(diag of the pivot row) */
2358       multiplier = (*pc) * (*pv);
2359       *pc        = multiplier;
2360       if (PetscAbsScalar(multiplier) > dt){
2361         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
2362         pv         = b->a + bdiag[row+1] + 1;
2363         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
2364         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
2365         /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */
2366       }
2367       k++;
2368     }
2369 
2370     /* finished row so stick it into b->a */
2371     /* L-part */
2372     pv = b->a + bi[i] ;
2373     pj = bj + bi[i] ;
2374     nzl = bi[i+1] - bi[i];
2375     for (j=0; j<nzl; j++) {
2376       pv[j] = rtmp[pj[j]];
2377       /* printf(" (%d,%g),",pj[j],pv[j]); */
2378     }
2379 
2380     /* diagonal: invert diagonal entries for simplier triangular solves */
2381     if (rtmp[i] == 0.0) rtmp[i] = dt+shift;
2382     b->a[bdiag[i]] = 1.0/rtmp[i];
2383     /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */
2384 
2385     /* U-part */
2386     pv = b->a + bdiag[i+1] + 1;
2387     pj = bj + bdiag[i+1] + 1;
2388     nzu = bdiag[i] - bdiag[i+1] - 1;
2389     for (j=0; j<nzu; j++) {
2390       pv[j] = rtmp[pj[j]];
2391       /* printf(" (%d,%g),",pj[j],pv[j]); */
2392     }
2393     /* printf("\n"); */
2394   }
2395 
2396   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2397   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2398   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2399 
2400   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
2401   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
2402   if (row_identity && col_identity) {
2403     C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering_iludt;
2404   } else {
2405     C->ops->solve   = MatSolve_SeqAIJ_iludt;
2406   }
2407   C->ops->solveadd           = 0;
2408   C->ops->solvetranspose     = 0;
2409   C->ops->solvetransposeadd  = 0;
2410   C->ops->matsolve           = 0;
2411   C->assembled    = PETSC_TRUE;
2412   C->preallocated = PETSC_TRUE;
2413   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
2414   PetscFunctionReturn(0);
2415 }
2416