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