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