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