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