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