xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision a1db5e057e92a1600c64e24b1f38d51fd8751642)
1 #define PETSCMAT_DLL
2 
3 
4 #include "../src/mat/impls/aij/seq/aij.h"
5 #include "../src/mat/impls/sbaij/seq/sbaij.h"
6 #include "petscbt.h"
7 #include "../src/mat/utils/freespace.h"
8 
9 EXTERN_C_BEGIN
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 = 0, newcurrent = 0,*order;
20   const PetscInt    *ai = a->i, *aj = a->j;
21   const PetscScalar *aa = a->a;
22   PetscTruth        *done;
23   PetscReal         best,past = 0,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 EXTERN_C_END
94 
95 EXTERN_C_BEGIN
96 #undef __FUNCT__
97 #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
98 PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
99 {
100   PetscFunctionBegin;
101   *flg = PETSC_TRUE;
102   PetscFunctionReturn(0);
103 }
104 EXTERN_C_END
105 
106 EXTERN_C_BEGIN
107 #undef __FUNCT__
108 #define __FUNCT__ "MatGetFactor_seqaij_petsc"
109 PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
110 {
111   PetscInt           n = A->rmap->n;
112   PetscErrorCode     ierr;
113 
114   PetscFunctionBegin;
115   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
116   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
117   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT){
118     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
119     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJ;
120     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ;
121     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqAIJ;
122   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
123     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
124     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
125     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
126     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
127   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
128   (*B)->factor = ftype;
129   PetscFunctionReturn(0);
130 }
131 EXTERN_C_END
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
135 PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
136 {
137   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
138   IS                 isicol;
139   PetscErrorCode     ierr;
140   const PetscInt     *r,*ic;
141   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
142   PetscInt           *bi,*bj,*ajtmp;
143   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
144   PetscReal          f;
145   PetscInt           nlnk,*lnk,k,**bi_ptr;
146   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
147   PetscBT            lnkbt;
148   PetscTruth         newdatastruct=PETSC_FALSE;
149 
150   PetscFunctionBegin;
151   ierr = PetscOptionsGetTruth(PETSC_NULL,"-lu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
152   if(newdatastruct){
153     ierr = MatLUFactorSymbolic_SeqAIJ_newdatastruct(B,A,isrow,iscol,info);CHKERRQ(ierr);
154     PetscFunctionReturn(0);
155   }
156 
157   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
158   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
159   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
160   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
161 
162   /* get new row pointers */
163   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
164   bi[0] = 0;
165 
166   /* bdiag is location of diagonal in factor */
167   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
168   bdiag[0] = 0;
169 
170   /* linked list for storing column indices of the active row */
171   nlnk = n + 1;
172   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
173 
174   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
175 
176   /* initial FreeSpace size is f*(ai[n]+1) */
177   f = info->fill;
178   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
179   current_space = free_space;
180 
181   for (i=0; i<n; i++) {
182     /* copy previous fill into linked list */
183     nzi = 0;
184     nnz = ai[r[i]+1] - ai[r[i]];
185     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
186     ajtmp = aj + ai[r[i]];
187     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
188     nzi += nlnk;
189 
190     /* add pivot rows into linked list */
191     row = lnk[n];
192     while (row < i) {
193       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
194       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
195       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
196       nzi += nlnk;
197       row  = lnk[row];
198     }
199     bi[i+1] = bi[i] + nzi;
200     im[i]   = nzi;
201 
202     /* mark bdiag */
203     nzbd = 0;
204     nnz  = nzi;
205     k    = lnk[n];
206     while (nnz-- && k < i){
207       nzbd++;
208       k = lnk[k];
209     }
210     bdiag[i] = bi[i] + nzbd;
211 
212     /* if free space is not available, make more free space */
213     if (current_space->local_remaining<nzi) {
214       nnz = (n - i)*nzi; /* estimated and max additional space needed */
215       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
216       reallocs++;
217     }
218 
219     /* copy data into free space, then initialize lnk */
220     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
221     bi_ptr[i] = current_space->array;
222     current_space->array           += nzi;
223     current_space->local_used      += nzi;
224     current_space->local_remaining -= nzi;
225   }
226 #if defined(PETSC_USE_INFO)
227   if (ai[n] != 0) {
228     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
229     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
230     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
231     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
232     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
233   } else {
234     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
235   }
236 #endif
237 
238   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
239   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
240 
241   /* destroy list of free space and other temporary array(s) */
242   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
243   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
244   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
245   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
246 
247   /* put together the new matrix */
248   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
249   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
250   b    = (Mat_SeqAIJ*)(B)->data;
251   b->free_a       = PETSC_TRUE;
252   b->free_ij      = PETSC_TRUE;
253   b->singlemalloc = PETSC_FALSE;
254   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
255   b->j          = bj;
256   b->i          = bi;
257   b->diag       = bdiag;
258   b->ilen       = 0;
259   b->imax       = 0;
260   b->row        = isrow;
261   b->col        = iscol;
262   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
263   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
264   b->icol       = isicol;
265   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
266 
267   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
268   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
269   b->maxnz = b->nz = bi[n] ;
270 
271   (B)->factor                = MAT_FACTOR_LU;
272   (B)->info.factor_mallocs   = reallocs;
273   (B)->info.fill_ratio_given = f;
274 
275   if (ai[n]) {
276     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
277   } else {
278     (B)->info.fill_ratio_needed = 0.0;
279   }
280   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
281   (B)->ops->solve            = MatSolve_SeqAIJ;
282   (B)->ops->solvetranspose   = MatSolveTranspose_SeqAIJ;
283   /* switch to inodes if appropriate */
284   ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_newdatastruct"
290 PetscErrorCode MatLUFactorSymbolic_SeqAIJ_newdatastruct(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
291 {
292   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
293   IS                 isicol;
294   PetscErrorCode     ierr;
295   const PetscInt     *r,*ic;
296   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
297   PetscInt           *bi,*bj,*ajtmp;
298   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
299   PetscReal          f;
300   PetscInt           nlnk,*lnk,k,**bi_ptr;
301   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
302   PetscBT            lnkbt;
303 
304   PetscFunctionBegin;
305   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
306   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
307   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
308   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
309 
310   /* get new row pointers */
311   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
312   bi[0] = 0;
313 
314   /* bdiag is location of diagonal in factor */
315   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
316   bdiag[0] = 0;
317 
318   /* linked list for storing column indices of the active row */
319   nlnk = n + 1;
320   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
321 
322   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
323 
324   /* initial FreeSpace size is f*(ai[n]+1) */
325   f = info->fill;
326   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
327   current_space = free_space;
328 
329   for (i=0; i<n; i++) {
330     /* copy previous fill into linked list */
331     nzi = 0;
332     nnz = ai[r[i]+1] - ai[r[i]];
333     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
334     ajtmp = aj + ai[r[i]];
335     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
336     nzi += nlnk;
337 
338     /* add pivot rows into linked list */
339     row = lnk[n];
340     while (row < i){
341       nzbd  = bdiag[row] + 1; /* num of entries in the row with column index <= row */
342       ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
343       ierr  = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
344       nzi  += nlnk;
345       row   = lnk[row];
346     }
347     bi[i+1] = bi[i] + nzi;
348     im[i]   = nzi;
349 
350     /* mark bdiag */
351     nzbd = 0;
352     nnz  = nzi;
353     k    = lnk[n];
354     while (nnz-- && k < i){
355       nzbd++;
356       k = lnk[k];
357     }
358     bdiag[i] = nzbd; /* note: bdiag[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */
359 
360     /* if free space is not available, make more free space */
361     if (current_space->local_remaining<nzi) {
362       nnz = 2*(n - i)*nzi; /* estimated and max additional space needed */
363       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
364       reallocs++;
365     }
366 
367     /* copy data into free space, then initialize lnk */
368     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
369     bi_ptr[i] = current_space->array;
370     current_space->array           += nzi;
371     current_space->local_used      += nzi;
372     current_space->local_remaining -= nzi;
373   }
374 #if defined(PETSC_USE_INFO)
375   if (ai[n] != 0) {
376     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
377     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
378     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
379     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
380     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
381   } else {
382     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
383   }
384 #endif
385 
386   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
387   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
388 
389   /* destroy list of free space and other temporary array(s) */
390   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
391   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
392   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
393   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
394 
395   /* put together the new matrix */
396   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
397   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
398   b    = (Mat_SeqAIJ*)(B)->data;
399   b->free_a       = PETSC_TRUE;
400   b->free_ij      = PETSC_TRUE;
401   b->singlemalloc = PETSC_FALSE;
402   ierr          = PetscMalloc((bi[2*n+1])*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
403   b->j          = bj;
404   b->i          = bi;
405   b->diag       = bdiag;
406   b->ilen       = 0;
407   b->imax       = 0;
408   b->row        = isrow;
409   b->col        = iscol;
410   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
411   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
412   b->icol       = isicol;
413   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
414 
415   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
416   ierr = PetscLogObjectMemory(B,bi[2*n+1]*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
417   b->maxnz = b->nz = bi[2*n+1] ;
418 
419   (B)->factor                = MAT_FACTOR_LU;
420   (B)->info.factor_mallocs   = reallocs;
421   (B)->info.fill_ratio_given = f;
422 
423   if (ai[n]) {
424     (B)->info.fill_ratio_needed = ((PetscReal)bi[2*n+1])/((PetscReal)ai[n]);
425   } else {
426     (B)->info.fill_ratio_needed = 0.0;
427   }
428   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ_newdatastruct;
429   PetscFunctionReturn(0);
430 }
431 
432 /*
433     Trouble in factorization, should we dump the original matrix?
434 */
435 #undef __FUNCT__
436 #define __FUNCT__ "MatFactorDumpMatrix"
437 PetscErrorCode MatFactorDumpMatrix(Mat A)
438 {
439   PetscErrorCode ierr;
440   PetscTruth     flg = PETSC_FALSE;
441 
442   PetscFunctionBegin;
443   ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_factor_dump_on_error",&flg,PETSC_NULL);CHKERRQ(ierr);
444   if (flg) {
445     PetscViewer viewer;
446     char        filename[PETSC_MAX_PATH_LEN];
447 
448     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
449     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
450     ierr = MatView(A,viewer);CHKERRQ(ierr);
451     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
452   }
453   PetscFunctionReturn(0);
454 }
455 
456 extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec);
457 
458 /* ----------------------------------------------------------- */
459 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct(Mat,Vec,Vec);
460 extern PetscErrorCode MatSolve_SeqAIJ_newdatastruct(Mat,Vec,Vec);
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_newdatastruct"
464 PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
465 {
466   Mat            C=B;
467   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
468   IS             isrow = b->row,isicol = b->icol;
469   PetscErrorCode ierr;
470   const PetscInt *r,*ic,*ics;
471   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
472   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
473   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
474   PetscReal      shift=info->shiftinblocks;
475   PetscTruth     row_identity,col_identity;
476 
477   PetscFunctionBegin;
478   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
479   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
480   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
481   ics  = ic;
482 
483   for (i=0; i<n; i++){
484     /* zero rtmp */
485     /* L part */
486     nz    = bi[i+1] - bi[i];
487     bjtmp = bj + bi[i];
488     for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
489 
490     /* U part */
491     nz = bi[2*n-i+1] - bi[2*n-i];
492     bjtmp = bj + bi[2*n-i];
493     for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
494 
495     /* load in initial (unfactored row) */
496     nz    = ai[r[i]+1] - ai[r[i]];
497     ajtmp = aj + ai[r[i]];
498     v     = aa + ai[r[i]];
499     for (j=0; j<nz; j++) {
500       rtmp[ics[ajtmp[j]]] = v[j];
501     }
502     if (rtmp[ics[r[i]]] == 0.0){
503       rtmp[ics[r[i]]] += shift; /* shift the diagonal of the matrix */
504       /* printf("row %d, shift %g\n",i,shift); */
505     }
506 
507     /* elimination */
508     bjtmp = bj + bi[i];
509     row   = *bjtmp++;
510     nzL   = bi[i+1] - bi[i];
511     k   = 0;
512     while  (k < nzL) {
513       pc = rtmp + row;
514       if (*pc != 0.0) {
515         pv         = b->a + bdiag[row];
516         multiplier = *pc * (*pv);
517         *pc        = multiplier;
518         pj         = b->j + bi[2*n-row]; /* begining of U(row,:) */
519         pv         = b->a + bi[2*n-row];
520         nz         = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries in U(row,:), excluding diag */
521         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
522         ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
523       }
524       row = *bjtmp++; k++;
525     }
526 
527     /* finished row so stick it into b->a */
528     /* L part */
529     pv   = b->a + bi[i] ;
530     pj   = b->j + bi[i] ;
531     nz   = bi[i+1] - bi[i];
532     for (j=0; j<nz; j++) {
533       pv[j] = rtmp[pj[j]];
534     }
535 
536     /* Mark diagonal and invert diagonal for simplier triangular solves */
537     pv  = b->a + bdiag[i];
538     pj  = b->j + bdiag[i];
539     /* if (*pj != i)SETERRQ2(PETSC_ERR_SUP,"row %d != *pj %d",i,*pj) */
540     *pv = 1.0/rtmp[*pj];
541 
542     /* U part */
543     pv = b->a + bi[2*n-i];
544     pj = b->j + bi[2*n-i];
545     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
546     for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]];
547   }
548   ierr = PetscFree(rtmp);CHKERRQ(ierr);
549   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
550   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
551 
552   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
553   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
554   if (row_identity && col_identity) {
555     C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_newdatastruct;
556   } else {
557     C->ops->solve = MatSolve_SeqAIJ_newdatastruct;
558   }
559 
560   C->ops->solveadd           = 0;
561   C->ops->solvetranspose     = 0;
562   C->ops->solvetransposeadd  = 0;
563   C->ops->matsolve           = 0;
564   C->assembled    = PETSC_TRUE;
565   C->preallocated = PETSC_TRUE;
566   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
567   PetscFunctionReturn(0);
568 }
569 
570 #undef __FUNCT__
571 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
572 PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
573 {
574   Mat             C=B;
575   Mat_SeqAIJ      *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
576   IS              isrow = b->row,isicol = b->icol;
577   PetscErrorCode  ierr;
578   const PetscInt   *r,*ic,*ics;
579   PetscInt        nz,row,i,j,n=A->rmap->n,diag;
580   const PetscInt  *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
581   const PetscInt  *ajtmp,*bjtmp,*diag_offset = b->diag,*pj;
582   MatScalar       *pv,*rtmp,*pc,multiplier,d;
583   const MatScalar *v,*aa=a->a;
584   PetscReal       rs=0.0;
585   LUShift_Ctx     sctx;
586   PetscInt        newshift,*ddiag;
587 
588   PetscFunctionBegin;
589   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
590   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
591   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
592   ics  = ic;
593 
594   sctx.shift_top      = 0;
595   sctx.nshift_max     = 0;
596   sctx.shift_lo       = 0;
597   sctx.shift_hi       = 0;
598   sctx.shift_fraction = 0;
599 
600   /* if both shift schemes are chosen by user, only use info->shiftpd */
601   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
602     ddiag          = a->diag;
603     sctx.shift_top = info->zeropivot;
604     for (i=0; i<n; i++) {
605       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
606       d  = (aa)[ddiag[i]];
607       rs = -PetscAbsScalar(d) - PetscRealPart(d);
608       v  = aa+ai[i];
609       nz = ai[i+1] - ai[i];
610       for (j=0; j<nz; j++)
611 	rs += PetscAbsScalar(v[j]);
612       if (rs>sctx.shift_top) sctx.shift_top = rs;
613     }
614     sctx.shift_top   *= 1.1;
615     sctx.nshift_max   = 5;
616     sctx.shift_lo     = 0.;
617     sctx.shift_hi     = 1.;
618   }
619 
620   sctx.shift_amount = 0.0;
621   sctx.nshift       = 0;
622   do {
623     sctx.lushift = PETSC_FALSE;
624     for (i=0; i<n; i++){
625       nz    = bi[i+1] - bi[i];
626       bjtmp = bj + bi[i];
627       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
628 
629       /* load in initial (unfactored row) */
630       nz    = ai[r[i]+1] - ai[r[i]];
631       ajtmp = aj + ai[r[i]];
632       v     = aa + ai[r[i]];
633       for (j=0; j<nz; j++) {
634         rtmp[ics[ajtmp[j]]] = v[j];
635       }
636       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
637       /* if (sctx.shift_amount > 0.0) printf("row %d, shift %g\n",i,sctx.shift_amount); */
638 
639       row = *bjtmp++;
640       while  (row < i) {
641         pc = rtmp + row;
642         if (*pc != 0.0) {
643           pv         = b->a + diag_offset[row];
644           pj         = b->j + diag_offset[row] + 1;
645           multiplier = *pc / *pv++;
646           *pc        = multiplier;
647           nz         = bi[row+1] - diag_offset[row] - 1;
648           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
649           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
650         }
651         row = *bjtmp++;
652       }
653       /* finished row so stick it into b->a */
654       pv   = b->a + bi[i] ;
655       pj   = b->j + bi[i] ;
656       nz   = bi[i+1] - bi[i];
657       diag = diag_offset[i] - bi[i];
658       rs   = 0.0;
659       for (j=0; j<nz; j++) {
660         pv[j] = rtmp[pj[j]];
661         rs   += PetscAbsScalar(pv[j]);
662       }
663       rs   -= PetscAbsScalar(pv[diag]);
664 
665       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
666       sctx.rs  = rs;
667       sctx.pv  = pv[diag];
668       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
669       if (newshift == 1) break;
670     }
671 
672     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
673       /*
674        * if no shift in this attempt & shifting & started shifting & can refine,
675        * then try lower shift
676        */
677       sctx.shift_hi       = sctx.shift_fraction;
678       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
679       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
680       sctx.lushift        = PETSC_TRUE;
681       sctx.nshift++;
682     }
683   } while (sctx.lushift);
684 
685   /* invert diagonal entries for simplier triangular solves */
686   for (i=0; i<n; i++) {
687     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
688   }
689   ierr = PetscFree(rtmp);CHKERRQ(ierr);
690   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
691   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
692   if (b->inode.use) {
693     C->ops->solve   = MatSolve_Inode;
694   } else {
695     PetscTruth row_identity, col_identity;
696     ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
697     ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
698     if (row_identity && col_identity) {
699       C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering;
700     } else {
701       C->ops->solve   = MatSolve_SeqAIJ;
702     }
703   }
704   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
705   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
706   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
707   C->ops->matsolve           = MatMatSolve_SeqAIJ;
708   C->assembled    = PETSC_TRUE;
709   C->preallocated = PETSC_TRUE;
710   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
711   if (sctx.nshift){
712      if (info->shiftpd) {
713       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);
714     } else if (info->shiftnz) {
715       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
716     }
717   }
718   PetscFunctionReturn(0);
719 }
720 
721 /*
722    This routine implements inplace ILU(0) with row or/and column permutations.
723    Input:
724      A - original matrix
725    Output;
726      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
727          a->j (col index) is permuted by the inverse of colperm, then sorted
728          a->a reordered accordingly with a->j
729          a->diag (ptr to diagonal elements) is updated.
730 */
731 #undef __FUNCT__
732 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
733 PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
734 {
735   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
736   IS             isrow = a->row,isicol = a->icol;
737   PetscErrorCode ierr;
738   const PetscInt *r,*ic,*ics;
739   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
740   PetscInt       *ajtmp,nz,row;
741   PetscInt       *diag = a->diag,nbdiag,*pj;
742   PetscScalar    *rtmp,*pc,multiplier,d;
743   MatScalar      *v,*pv;
744   PetscReal      rs;
745   LUShift_Ctx    sctx;
746   PetscInt       newshift;
747 
748   PetscFunctionBegin;
749   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
750   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
751   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
752   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
753   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
754   ics = ic;
755 
756   sctx.shift_top      = 0;
757   sctx.nshift_max     = 0;
758   sctx.shift_lo       = 0;
759   sctx.shift_hi       = 0;
760   sctx.shift_fraction = 0;
761 
762   /* if both shift schemes are chosen by user, only use info->shiftpd */
763   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
764     sctx.shift_top = 0;
765     for (i=0; i<n; i++) {
766       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
767       d  = (a->a)[diag[i]];
768       rs = -PetscAbsScalar(d) - PetscRealPart(d);
769       v  = a->a+ai[i];
770       nz = ai[i+1] - ai[i];
771       for (j=0; j<nz; j++)
772 	rs += PetscAbsScalar(v[j]);
773       if (rs>sctx.shift_top) sctx.shift_top = rs;
774     }
775     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
776     sctx.shift_top    *= 1.1;
777     sctx.nshift_max   = 5;
778     sctx.shift_lo     = 0.;
779     sctx.shift_hi     = 1.;
780   }
781 
782   sctx.shift_amount = 0;
783   sctx.nshift       = 0;
784   do {
785     sctx.lushift = PETSC_FALSE;
786     for (i=0; i<n; i++){
787       /* load in initial unfactored row */
788       nz    = ai[r[i]+1] - ai[r[i]];
789       ajtmp = aj + ai[r[i]];
790       v     = a->a + ai[r[i]];
791       /* sort permuted ajtmp and values v accordingly */
792       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
793       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
794 
795       diag[r[i]] = ai[r[i]];
796       for (j=0; j<nz; j++) {
797         rtmp[ajtmp[j]] = v[j];
798         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
799       }
800       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
801 
802       row = *ajtmp++;
803       while  (row < i) {
804         pc = rtmp + row;
805         if (*pc != 0.0) {
806           pv         = a->a + diag[r[row]];
807           pj         = aj + diag[r[row]] + 1;
808 
809           multiplier = *pc / *pv++;
810           *pc        = multiplier;
811           nz         = ai[r[row]+1] - diag[r[row]] - 1;
812           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
813           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
814         }
815         row = *ajtmp++;
816       }
817       /* finished row so overwrite it onto a->a */
818       pv   = a->a + ai[r[i]] ;
819       pj   = aj + ai[r[i]] ;
820       nz   = ai[r[i]+1] - ai[r[i]];
821       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
822 
823       rs   = 0.0;
824       for (j=0; j<nz; j++) {
825         pv[j] = rtmp[pj[j]];
826         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
827       }
828 
829       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
830       sctx.rs  = rs;
831       sctx.pv  = pv[nbdiag];
832       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
833       if (newshift == 1) break;
834     }
835 
836     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
837       /*
838        * if no shift in this attempt & shifting & started shifting & can refine,
839        * then try lower shift
840        */
841       sctx.shift_hi        = sctx.shift_fraction;
842       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
843       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
844       sctx.lushift         = PETSC_TRUE;
845       sctx.nshift++;
846     }
847   } while (sctx.lushift);
848 
849   /* invert diagonal entries for simplier triangular solves */
850   for (i=0; i<n; i++) {
851     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
852   }
853 
854   ierr = PetscFree(rtmp);CHKERRQ(ierr);
855   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
856   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
857   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
858   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
859   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
860   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
861   A->assembled = PETSC_TRUE;
862   A->preallocated = PETSC_TRUE;
863   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
864   if (sctx.nshift){
865     if (info->shiftpd) {
866       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);
867     } else if (info->shiftnz) {
868       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
869     }
870   }
871   PetscFunctionReturn(0);
872 }
873 
874 /* ----------------------------------------------------------- */
875 #undef __FUNCT__
876 #define __FUNCT__ "MatLUFactor_SeqAIJ"
877 PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
878 {
879   PetscErrorCode ierr;
880   Mat            C;
881 
882   PetscFunctionBegin;
883   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
884   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
885   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
886   A->ops->solve            = C->ops->solve;
887   A->ops->solvetranspose   = C->ops->solvetranspose;
888   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
889   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
890   PetscFunctionReturn(0);
891 }
892 /* ----------------------------------------------------------- */
893 
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "MatSolve_SeqAIJ"
897 PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
898 {
899   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
900   IS                iscol = a->col,isrow = a->row;
901   PetscErrorCode    ierr;
902   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
903   PetscInt          nz;
904   const PetscInt    *rout,*cout,*r,*c;
905   PetscScalar       *x,*tmp,*tmps,sum;
906   const PetscScalar *b;
907   const MatScalar   *aa = a->a,*v;
908 
909   PetscFunctionBegin;
910   if (!n) PetscFunctionReturn(0);
911 
912   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
913   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
914   tmp  = a->solve_work;
915 
916   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
917   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
918 
919   /* forward solve the lower triangular */
920   tmp[0] = b[*r++];
921   tmps   = tmp;
922   for (i=1; i<n; i++) {
923     v   = aa + ai[i] ;
924     vi  = aj + ai[i] ;
925     nz  = a->diag[i] - ai[i];
926     sum = b[*r++];
927     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
928     tmp[i] = sum;
929   }
930 
931   /* backward solve the upper triangular */
932   for (i=n-1; i>=0; i--){
933     v   = aa + a->diag[i] + 1;
934     vi  = aj + a->diag[i] + 1;
935     nz  = ai[i+1] - a->diag[i] - 1;
936     sum = tmp[i];
937     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
938     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
939   }
940 
941   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
942   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
943   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
944   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
945   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "MatMatSolve_SeqAIJ"
951 PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
952 {
953   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
954   IS              iscol = a->col,isrow = a->row;
955   PetscErrorCode  ierr;
956   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
957   PetscInt        nz,neq;
958   const PetscInt  *rout,*cout,*r,*c;
959   PetscScalar     *x,*b,*tmp,*tmps,sum;
960   const MatScalar *aa = a->a,*v;
961   PetscTruth      bisdense,xisdense;
962 
963   PetscFunctionBegin;
964   if (!n) PetscFunctionReturn(0);
965 
966   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
967   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
968   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
969   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
970 
971   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
972   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
973 
974   tmp  = a->solve_work;
975   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
976   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
977 
978   for (neq=0; neq<B->cmap->n; neq++){
979     /* forward solve the lower triangular */
980     tmp[0] = b[r[0]];
981     tmps   = tmp;
982     for (i=1; i<n; i++) {
983       v   = aa + ai[i] ;
984       vi  = aj + ai[i] ;
985       nz  = a->diag[i] - ai[i];
986       sum = b[r[i]];
987       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
988       tmp[i] = sum;
989     }
990     /* backward solve the upper triangular */
991     for (i=n-1; i>=0; i--){
992       v   = aa + a->diag[i] + 1;
993       vi  = aj + a->diag[i] + 1;
994       nz  = ai[i+1] - a->diag[i] - 1;
995       sum = tmp[i];
996       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
997       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
998     }
999 
1000     b += n;
1001     x += n;
1002   }
1003   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1004   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1005   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
1006   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
1007   ierr = PetscLogFlops(B->cmap->n*(2.0*a->nz - n));CHKERRQ(ierr);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
1013 PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
1014 {
1015   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1016   IS              iscol = a->col,isrow = a->row;
1017   PetscErrorCode  ierr;
1018   const PetscInt  *r,*c,*rout,*cout;
1019   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1020   PetscInt        nz,row;
1021   PetscScalar     *x,*b,*tmp,*tmps,sum;
1022   const MatScalar *aa = a->a,*v;
1023 
1024   PetscFunctionBegin;
1025   if (!n) PetscFunctionReturn(0);
1026 
1027   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1028   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1029   tmp  = a->solve_work;
1030 
1031   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1032   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1033 
1034   /* forward solve the lower triangular */
1035   tmp[0] = b[*r++];
1036   tmps   = tmp;
1037   for (row=1; row<n; row++) {
1038     i   = rout[row]; /* permuted row */
1039     v   = aa + ai[i] ;
1040     vi  = aj + ai[i] ;
1041     nz  = a->diag[i] - ai[i];
1042     sum = b[*r++];
1043     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1044     tmp[row] = sum;
1045   }
1046 
1047   /* backward solve the upper triangular */
1048   for (row=n-1; row>=0; row--){
1049     i   = rout[row]; /* permuted row */
1050     v   = aa + a->diag[i] + 1;
1051     vi  = aj + a->diag[i] + 1;
1052     nz  = ai[i+1] - a->diag[i] - 1;
1053     sum = tmp[row];
1054     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1055     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
1056   }
1057 
1058   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1059   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1060   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1061   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1062   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 /* ----------------------------------------------------------- */
1067 #include "../src/mat/impls/aij/seq/ftn-kernels/fsolve.h"
1068 #undef __FUNCT__
1069 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
1070 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
1071 {
1072   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1073   PetscErrorCode    ierr;
1074   PetscInt          n = A->rmap->n;
1075   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag;
1076   PetscScalar       *x;
1077   const PetscScalar *b;
1078   const MatScalar   *aa = a->a;
1079 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1080   PetscInt          adiag_i,i,nz,ai_i;
1081   const PetscInt    *vi;
1082   const MatScalar   *v;
1083   PetscScalar       sum;
1084 #endif
1085 
1086   PetscFunctionBegin;
1087   if (!n) PetscFunctionReturn(0);
1088 
1089   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
1090   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1091 
1092 #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1093   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
1094 #else
1095   /* forward solve the lower triangular */
1096   x[0] = b[0];
1097   for (i=1; i<n; i++) {
1098     ai_i = ai[i];
1099     v    = aa + ai_i;
1100     vi   = aj + ai_i;
1101     nz   = adiag[i] - ai_i;
1102     sum  = b[i];
1103     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1104     x[i] = sum;
1105   }
1106 
1107   /* backward solve the upper triangular */
1108   for (i=n-1; i>=0; i--){
1109     adiag_i = adiag[i];
1110     v       = aa + adiag_i + 1;
1111     vi      = aj + adiag_i + 1;
1112     nz      = ai[i+1] - adiag_i - 1;
1113     sum     = x[i];
1114     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1115     x[i]    = sum*aa[adiag_i];
1116   }
1117 #endif
1118   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
1119   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
1120   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1121   PetscFunctionReturn(0);
1122 }
1123 
1124 #undef __FUNCT__
1125 #define __FUNCT__ "MatSolveAdd_SeqAIJ"
1126 PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1127 {
1128   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1129   IS              iscol = a->col,isrow = a->row;
1130   PetscErrorCode  ierr;
1131   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1132   PetscInt        nz;
1133   const PetscInt  *rout,*cout,*r,*c;
1134   PetscScalar     *x,*b,*tmp,sum;
1135   const MatScalar *aa = a->a,*v;
1136 
1137   PetscFunctionBegin;
1138   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
1139 
1140   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1141   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1142   tmp  = a->solve_work;
1143 
1144   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1145   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1146 
1147   /* forward solve the lower triangular */
1148   tmp[0] = b[*r++];
1149   for (i=1; i<n; i++) {
1150     v   = aa + ai[i] ;
1151     vi  = aj + ai[i] ;
1152     nz  = a->diag[i] - ai[i];
1153     sum = b[*r++];
1154     while (nz--) sum -= *v++ * tmp[*vi++ ];
1155     tmp[i] = sum;
1156   }
1157 
1158   /* backward solve the upper triangular */
1159   for (i=n-1; i>=0; i--){
1160     v   = aa + a->diag[i] + 1;
1161     vi  = aj + a->diag[i] + 1;
1162     nz  = ai[i+1] - a->diag[i] - 1;
1163     sum = tmp[i];
1164     while (nz--) sum -= *v++ * tmp[*vi++ ];
1165     tmp[i] = sum*aa[a->diag[i]];
1166     x[*c--] += tmp[i];
1167   }
1168 
1169   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1170   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1171   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1172   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1173   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1174 
1175   PetscFunctionReturn(0);
1176 }
1177 /* -------------------------------------------------------------------*/
1178 #undef __FUNCT__
1179 #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1180 PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1181 {
1182   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1183   IS              iscol = a->col,isrow = a->row;
1184   PetscErrorCode  ierr;
1185   const PetscInt  *rout,*cout,*r,*c;
1186   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1187   PetscInt        nz,*diag = a->diag;
1188   PetscScalar     *x,*b,*tmp,s1;
1189   const MatScalar *aa = a->a,*v;
1190 
1191   PetscFunctionBegin;
1192   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1193   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1194   tmp  = a->solve_work;
1195 
1196   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1197   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1198 
1199   /* copy the b into temp work space according to permutation */
1200   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1201 
1202   /* forward solve the U^T */
1203   for (i=0; i<n; i++) {
1204     v   = aa + diag[i] ;
1205     vi  = aj + diag[i] + 1;
1206     nz  = ai[i+1] - diag[i] - 1;
1207     s1  = tmp[i];
1208     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1209     while (nz--) {
1210       tmp[*vi++ ] -= (*v++)*s1;
1211     }
1212     tmp[i] = s1;
1213   }
1214 
1215   /* backward solve the L^T */
1216   for (i=n-1; i>=0; i--){
1217     v   = aa + diag[i] - 1 ;
1218     vi  = aj + diag[i] - 1 ;
1219     nz  = diag[i] - ai[i];
1220     s1  = tmp[i];
1221     while (nz--) {
1222       tmp[*vi-- ] -= (*v--)*s1;
1223     }
1224   }
1225 
1226   /* copy tmp into x according to permutation */
1227   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1228 
1229   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1230   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1231   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1232   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1233 
1234   ierr = PetscLogFlops(2.0*a->nz-A->cmap->n);CHKERRQ(ierr);
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1240 PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1241 {
1242   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1243   IS              iscol = a->col,isrow = a->row;
1244   PetscErrorCode  ierr;
1245   const PetscInt  *r,*c,*rout,*cout;
1246   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1247   PetscInt        nz,*diag = a->diag;
1248   PetscScalar     *x,*b,*tmp;
1249   const MatScalar *aa = a->a,*v;
1250 
1251   PetscFunctionBegin;
1252   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
1253 
1254   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1255   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1256   tmp = a->solve_work;
1257 
1258   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1259   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1260 
1261   /* copy the b into temp work space according to permutation */
1262   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1263 
1264   /* forward solve the U^T */
1265   for (i=0; i<n; i++) {
1266     v   = aa + diag[i] ;
1267     vi  = aj + diag[i] + 1;
1268     nz  = ai[i+1] - diag[i] - 1;
1269     tmp[i] *= *v++;
1270     while (nz--) {
1271       tmp[*vi++ ] -= (*v++)*tmp[i];
1272     }
1273   }
1274 
1275   /* backward solve the L^T */
1276   for (i=n-1; i>=0; i--){
1277     v   = aa + diag[i] - 1 ;
1278     vi  = aj + diag[i] - 1 ;
1279     nz  = diag[i] - ai[i];
1280     while (nz--) {
1281       tmp[*vi-- ] -= (*v--)*tmp[i];
1282     }
1283   }
1284 
1285   /* copy tmp into x according to permutation */
1286   for (i=0; i<n; i++) x[r[i]] += tmp[i];
1287 
1288   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1289   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1290   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1291   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1292 
1293   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1294   PetscFunctionReturn(0);
1295 }
1296 /* ----------------------------------------------------------------*/
1297 EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1298 EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption,PetscTruth);
1299 
1300 /*
1301    ilu(0) with natural ordering under new data structure.
1302    Factored arrays bj and ba are stored as
1303      L(0,:), L(1,:), ...,L(n-1,:),  U(n-1,:),...,U(i,:),U(i-1,:),...,U(0,:)
1304 
1305    bi=fact->i is an array of size 2n+2, in which
1306    bi+
1307      bi[i]      ->  1st entry of L(i,:),i=0,...,i-1
1308      bi[n]      ->  points to L(n-1,:)+1
1309      bi[n+1]    ->  1st entry of U(n-1,:)
1310      bi[2n-i]   ->  1st entry of U(i,:)
1311      bi[2n-i+1] ->  end of U(i,:)+1, the 1st entry of U(i-1,:)
1312      bi[2n]     ->  1st entry of U(0,:)
1313      bi[2n+1]   ->  points to U(0,:)+1
1314 
1315    U(i,:) contains diag[i] as its last entry, i.e.,
1316     U(i,:) = (u[i,i+1],...,u[i,n-1],diag[i])
1317 */
1318 #undef __FUNCT__
1319 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct"
1320 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1321 {
1322 
1323   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1324   PetscErrorCode     ierr;
1325   PetscInt           n=A->rmap->n,*ai=a->i,*aj,*adiag=a->diag;
1326   PetscInt           i,j,nz,*bi,*bj,*bdiag;
1327 
1328   PetscFunctionBegin;
1329   ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
1330   b    = (Mat_SeqAIJ*)(fact)->data;
1331 
1332   /* allocate matrix arrays for new data structure */
1333   ierr = PetscMalloc3(ai[n]+1,PetscScalar,&b->a,ai[n]+1,PetscInt,&b->j,2*n+2,PetscInt,&b->i);CHKERRQ(ierr);
1334   ierr = PetscLogObjectMemory(fact,ai[n]*(sizeof(PetscScalar)+sizeof(PetscInt))+(2*n+2)*sizeof(PetscInt));CHKERRQ(ierr);
1335   b->singlemalloc = PETSC_TRUE;
1336   if (!b->diag){
1337     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&b->diag);CHKERRQ(ierr);
1338   }
1339   bdiag = b->diag;
1340 
1341   if (n > 0) {
1342     ierr = PetscMemzero(b->a,(ai[n])*sizeof(MatScalar));CHKERRQ(ierr);
1343   }
1344 
1345   /* set bi and bj with new data structure */
1346   bi = b->i;
1347   bj = b->j;
1348 
1349   /* L part */
1350   bi[0] = 0;
1351   for (i=0; i<n; i++){
1352     nz = adiag[i] - ai[i];
1353     bi[i+1] = bi[i] + nz;
1354     aj = a->j + ai[i];
1355     for (j=0; j<nz; j++){
1356       *bj = aj[j]; bj++;
1357     }
1358   }
1359 
1360   /* U part */
1361   bi[n+1] = bi[n];
1362   for (i=n-1; i>=0; i--){
1363     nz = ai[i+1] - adiag[i] - 1;
1364     bi[2*n-i+1] = bi[2*n-i] + nz + 1;
1365     aj = a->j + adiag[i] + 1;
1366     for (j=0; j<nz; j++){
1367       *bj = aj[j]; bj++;
1368     }
1369     /* diag[i] */
1370     *bj = i; bj++;
1371     bdiag[i] = bi[2*n-i+1]-1;
1372   }
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 #undef __FUNCT__
1377 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_newdatastruct"
1378 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1379 {
1380   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1381   IS                 isicol;
1382   PetscErrorCode     ierr;
1383   const PetscInt     *r,*ic;
1384   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
1385   PetscInt           *bi,*cols,nnz,*cols_lvl;
1386   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1387   PetscInt           i,levels,diagonal_fill;
1388   PetscTruth         col_identity,row_identity;
1389   PetscReal          f;
1390   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1391   PetscBT            lnkbt;
1392   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1393   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1394   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1395   PetscTruth         missing;
1396 
1397   PetscFunctionBegin;
1398   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);
1399   f             = info->fill;
1400   levels        = (PetscInt)info->levels;
1401   diagonal_fill = (PetscInt)info->diagonal_fill;
1402   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1403 
1404   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1405   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
1406 
1407   if (!levels && row_identity && col_identity) {
1408     /* special case: ilu(0) with natural ordering */
1409     ierr = MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1410     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ_newdatastruct;
1411 
1412     fact->factor = MAT_FACTOR_ILU;
1413     (fact)->info.factor_mallocs    = 0;
1414     (fact)->info.fill_ratio_given  = info->fill;
1415     (fact)->info.fill_ratio_needed = 1.0;
1416     b               = (Mat_SeqAIJ*)(fact)->data;
1417     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1418     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1419     b->row              = isrow;
1420     b->col              = iscol;
1421     b->icol             = isicol;
1422     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1423     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1424     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1425     /* ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); */
1426     PetscFunctionReturn(0);
1427   }
1428 
1429   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1430   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1431 
1432   /* get new row pointers */
1433   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1434   bi[0] = 0;
1435   /* bdiag is location of diagonal in factor */
1436   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1437   bdiag[0]  = 0;
1438 
1439   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
1440   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
1441 
1442   /* create a linked list for storing column indices of the active row */
1443   nlnk = n + 1;
1444   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1445 
1446   /* initial FreeSpace size is f*(ai[n]+1) */
1447   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
1448   current_space = free_space;
1449   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
1450   current_space_lvl = free_space_lvl;
1451 
1452   for (i=0; i<n; i++) {
1453     nzi = 0;
1454     /* copy current row into linked list */
1455     nnz  = ai[r[i]+1] - ai[r[i]];
1456     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1457     cols = aj + ai[r[i]];
1458     lnk[i] = -1; /* marker to indicate if diagonal exists */
1459     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1460     nzi += nlnk;
1461 
1462     /* make sure diagonal entry is included */
1463     if (diagonal_fill && lnk[i] == -1) {
1464       fm = n;
1465       while (lnk[fm] < i) fm = lnk[fm];
1466       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1467       lnk[fm]    = i;
1468       lnk_lvl[i] = 0;
1469       nzi++; dcount++;
1470     }
1471 
1472     /* add pivot rows into the active row */
1473     nzbd = 0;
1474     prow = lnk[n];
1475     while (prow < i) {
1476       nnz      = bdiag[prow];
1477       cols     = bj_ptr[prow] + nnz + 1;
1478       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1479       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1480       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
1481       nzi += nlnk;
1482       prow = lnk[prow];
1483       nzbd++;
1484     }
1485     bdiag[i] = nzbd;
1486     bi[i+1]  = bi[i] + nzi;
1487 
1488     /* if free space is not available, make more free space */
1489     if (current_space->local_remaining<nzi) {
1490       nnz = 2*nzi*(n - i); /* estimated and max additional space needed */
1491       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1492       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
1493       reallocs++;
1494     }
1495 
1496     /* copy data into free_space and free_space_lvl, then initialize lnk */
1497     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1498     bj_ptr[i]    = current_space->array;
1499     bjlvl_ptr[i] = current_space_lvl->array;
1500 
1501     /* make sure the active row i has diagonal entry */
1502     if (*(bj_ptr[i]+bdiag[i]) != i) {
1503       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1504     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
1505     }
1506 
1507     current_space->array           += nzi;
1508     current_space->local_used      += nzi;
1509     current_space->local_remaining -= nzi;
1510     current_space_lvl->array           += nzi;
1511     current_space_lvl->local_used      += nzi;
1512     current_space_lvl->local_remaining -= nzi;
1513   }
1514 
1515   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1516   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1517 
1518   /* destroy list of free space and other temporary arrays */
1519   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1520 
1521   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
1522   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
1523 
1524   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1525   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1526   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
1527 
1528 #if defined(PETSC_USE_INFO)
1529   {
1530     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1531     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1532     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1533     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1534     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1535     if (diagonal_fill) {
1536       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1537     }
1538   }
1539 #endif
1540 
1541   /* put together the new matrix */
1542   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1543   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1544   b = (Mat_SeqAIJ*)(fact)->data;
1545   b->free_a       = PETSC_TRUE;
1546   b->free_ij      = PETSC_TRUE;
1547   b->singlemalloc = PETSC_FALSE;
1548   ierr = PetscMalloc( (bi[2*n+1] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
1549   b->j          = bj;
1550   b->i          = bi;
1551   b->diag       = bdiag;
1552   b->ilen       = 0;
1553   b->imax       = 0;
1554   b->row        = isrow;
1555   b->col        = iscol;
1556   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1557   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1558   b->icol       = isicol;
1559   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1560   /* In b structure:  Free imax, ilen, old a, old j.
1561      Allocate bdiag, solve_work, new a, new j */
1562   ierr = PetscLogObjectMemory(fact,bi[2*n+1] * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1563   b->maxnz = b->nz = bi[2*n+1] ;
1564   (fact)->info.factor_mallocs    = reallocs;
1565   (fact)->info.fill_ratio_given  = f;
1566   (fact)->info.fill_ratio_needed = ((PetscReal)bi[2*n+1])/((PetscReal)ai[n]);
1567   (fact)->ops->lufactornumeric   = MatLUFactorNumeric_SeqAIJ_newdatastruct;
1568   /* ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); */
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 #undef __FUNCT__
1573 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
1574 PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1575 {
1576   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1577   IS                 isicol;
1578   PetscErrorCode     ierr;
1579   const PetscInt     *r,*ic;
1580   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
1581   PetscInt           *bi,*cols,nnz,*cols_lvl;
1582   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1583   PetscInt           i,levels,diagonal_fill;
1584   PetscTruth         col_identity,row_identity;
1585   PetscReal          f;
1586   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1587   PetscBT            lnkbt;
1588   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1589   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1590   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
1591   PetscTruth         missing;
1592   PetscTruth         newdatastruct=PETSC_FALSE;
1593 
1594   PetscFunctionBegin;
1595   ierr = PetscOptionsGetTruth(PETSC_NULL,"-ilu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
1596   if (newdatastruct){
1597     ierr = MatILUFactorSymbolic_SeqAIJ_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1598     PetscFunctionReturn(0);
1599   }
1600 
1601   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);
1602   f             = info->fill;
1603   levels        = (PetscInt)info->levels;
1604   diagonal_fill = (PetscInt)info->diagonal_fill;
1605   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1606 
1607   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1608   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
1609   if (!levels && row_identity && col_identity) { /* special case: ilu(0) with natural ordering */
1610     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
1611     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1612 
1613     fact->factor = MAT_FACTOR_ILU;
1614     (fact)->info.factor_mallocs    = 0;
1615     (fact)->info.fill_ratio_given  = info->fill;
1616     (fact)->info.fill_ratio_needed = 1.0;
1617     b               = (Mat_SeqAIJ*)(fact)->data;
1618     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
1619     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1620     b->row              = isrow;
1621     b->col              = iscol;
1622     b->icol             = isicol;
1623     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1624     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1625     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1626     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1627     PetscFunctionReturn(0);
1628   }
1629 
1630   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1631   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1632 
1633   /* get new row pointers */
1634   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1635   bi[0] = 0;
1636   /* bdiag is location of diagonal in factor */
1637   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1638   bdiag[0]  = 0;
1639 
1640   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
1641   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
1642 
1643   /* create a linked list for storing column indices of the active row */
1644   nlnk = n + 1;
1645   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1646 
1647   /* initial FreeSpace size is f*(ai[n]+1) */
1648   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
1649   current_space = free_space;
1650   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
1651   current_space_lvl = free_space_lvl;
1652 
1653   for (i=0; i<n; i++) {
1654     nzi = 0;
1655     /* copy current row into linked list */
1656     nnz  = ai[r[i]+1] - ai[r[i]];
1657     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1658     cols = aj + ai[r[i]];
1659     lnk[i] = -1; /* marker to indicate if diagonal exists */
1660     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1661     nzi += nlnk;
1662 
1663     /* make sure diagonal entry is included */
1664     if (diagonal_fill && lnk[i] == -1) {
1665       fm = n;
1666       while (lnk[fm] < i) fm = lnk[fm];
1667       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1668       lnk[fm]    = i;
1669       lnk_lvl[i] = 0;
1670       nzi++; dcount++;
1671     }
1672 
1673     /* add pivot rows into the active row */
1674     nzbd = 0;
1675     prow = lnk[n];
1676     while (prow < i) {
1677       nnz      = bdiag[prow];
1678       cols     = bj_ptr[prow] + nnz + 1;
1679       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1680       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1681       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
1682       nzi += nlnk;
1683       prow = lnk[prow];
1684       nzbd++;
1685     }
1686     bdiag[i] = nzbd;
1687     bi[i+1]  = bi[i] + nzi;
1688 
1689     /* if free space is not available, make more free space */
1690     if (current_space->local_remaining<nzi) {
1691       nnz = nzi*(n - i); /* estimated and max additional space needed */
1692       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1693       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
1694       reallocs++;
1695     }
1696 
1697     /* copy data into free_space and free_space_lvl, then initialize lnk */
1698     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1699     bj_ptr[i]    = current_space->array;
1700     bjlvl_ptr[i] = current_space_lvl->array;
1701 
1702     /* make sure the active row i has diagonal entry */
1703     if (*(bj_ptr[i]+bdiag[i]) != i) {
1704       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1705     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
1706     }
1707 
1708     current_space->array           += nzi;
1709     current_space->local_used      += nzi;
1710     current_space->local_remaining -= nzi;
1711     current_space_lvl->array           += nzi;
1712     current_space_lvl->local_used      += nzi;
1713     current_space_lvl->local_remaining -= nzi;
1714   }
1715 
1716   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1717   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1718 
1719   /* destroy list of free space and other temporary arrays */
1720   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1721   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); /* copy free_space -> bj */
1722   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1723   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1724   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
1725 
1726 #if defined(PETSC_USE_INFO)
1727   {
1728     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1729     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1730     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1731     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1732     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1733     if (diagonal_fill) {
1734       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1735     }
1736   }
1737 #endif
1738 
1739   /* put together the new matrix */
1740   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1741   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1742   b = (Mat_SeqAIJ*)(fact)->data;
1743   b->free_a       = PETSC_TRUE;
1744   b->free_ij      = PETSC_TRUE;
1745   b->singlemalloc = PETSC_FALSE;
1746   ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
1747   b->j          = bj;
1748   b->i          = bi;
1749   for (i=0; i<n; i++) bdiag[i] += bi[i];
1750   b->diag       = bdiag;
1751   b->ilen       = 0;
1752   b->imax       = 0;
1753   b->row        = isrow;
1754   b->col        = iscol;
1755   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1756   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1757   b->icol       = isicol;
1758   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1759   /* In b structure:  Free imax, ilen, old a, old j.
1760      Allocate bdiag, solve_work, new a, new j */
1761   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1762   b->maxnz             = b->nz = bi[n] ;
1763   (fact)->info.factor_mallocs    = reallocs;
1764   (fact)->info.fill_ratio_given  = f;
1765   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1766   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1767   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
1768   PetscFunctionReturn(0);
1769 }
1770 
1771 #undef __FUNCT__
1772 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_newdatastruct"
1773 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
1774 {
1775   Mat            C = B;
1776   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1777   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
1778   IS             ip=b->row,iip = b->icol;
1779   PetscErrorCode ierr;
1780   const PetscInt *rip,*riip;
1781   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1782   PetscInt       *ai=a->i,*aj=a->j;
1783   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1784   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1785   PetscReal      zeropivot,shiftnz;
1786   PetscReal      shiftpd;
1787   PetscTruth     perm_identity;
1788 
1789   PetscFunctionBegin;
1790 
1791   shiftnz   = info->shiftnz;
1792   shiftpd   = info->shiftpd;
1793   zeropivot = info->zeropivot;
1794 
1795   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1796   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
1797 
1798   /* allocate working arrays
1799      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1800      il:  for active k row, il[i] gives the index of the 1st nonzero entry in U[i,k:n-1] in bj and ba arrays
1801   */
1802   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1803   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1804   c2r  = il + mbs;
1805   rtmp = (MatScalar*)(c2r + mbs);
1806 
1807   for (i=0; i<mbs; i++) {
1808     c2r[i] = mbs;
1809   }
1810   il[0] = 0;
1811 
1812   for (k = 0; k<mbs; k++){
1813     /* zero rtmp */
1814     nz = bi[k+1] - bi[k];
1815     bjtmp = bj + bi[k];
1816     for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1817 
1818     /* load in initial unfactored row */
1819     bval = ba + bi[k];
1820     jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1821     for (j = jmin; j < jmax; j++){
1822       col = riip[aj[j]];
1823       if (col >= k){ /* only take upper triangular entry */
1824         rtmp[col] = aa[j];
1825         *bval++   = 0.0; /* for in-place factorization */
1826       }
1827     }
1828     /* shift the diagonal of the matrix */
1829 
1830     /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1831     dk = rtmp[k];
1832     i  = c2r[k]; /* first row to be added to k_th row  */
1833 
1834     while (i < k){
1835       nexti = c2r[i]; /* next row to be added to k_th row */
1836 
1837       /* compute multiplier, update diag(k) and U(i,k) */
1838       ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1839       uikdi = - ba[ili]*ba[bdiag[i]];  /* diagonal(k) */
1840       dk   += uikdi*ba[ili]; /* update diag[k] */
1841       ba[ili] = uikdi; /* -U(i,k) */
1842 
1843       /* add multiple of row i to k-th row */
1844       jmin = ili + 1; jmax = bi[i+1];
1845       if (jmin < jmax){
1846         for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1847         /* update il and c2r for row i */
1848         il[i] = jmin;
1849         j = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1850       }
1851       i = nexti;
1852     }
1853 
1854     /* copy data into U(k,:) */
1855     ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1856     jmin = bi[k]; jmax = bi[k+1]-1;
1857     if (jmin < jmax) {
1858       for (j=jmin; j<jmax; j++){
1859         col = bj[j]; ba[j] = rtmp[col];
1860       }
1861       /* add the k-th row into il and c2r */
1862       il[k] = jmin;
1863       i = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1864     }
1865   }
1866 
1867   ierr = PetscFree(il);CHKERRQ(ierr);
1868   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1869   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1870 
1871   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1872   if (perm_identity){
1873     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct;
1874     (B)->ops->solvetranspose  = 0;
1875     (B)->ops->forwardsolve    = 0;
1876     (B)->ops->backwardsolve   = 0;
1877   } else {
1878     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_newdatastruct;
1879     (B)->ops->solvetranspose  = 0;
1880     (B)->ops->forwardsolve    = 0;
1881     (B)->ops->backwardsolve   = 0;
1882   }
1883 
1884   C->assembled    = PETSC_TRUE;
1885   C->preallocated = PETSC_TRUE;
1886   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1887   PetscFunctionReturn(0);
1888 }
1889 
1890 #undef __FUNCT__
1891 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1892 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1893 {
1894   Mat            C = B;
1895   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1896   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
1897   IS             ip=b->row,iip = b->icol;
1898   PetscErrorCode ierr;
1899   const PetscInt *rip,*riip;
1900   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
1901   PetscInt       *ai=a->i,*aj=a->j;
1902   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1903   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1904   PetscReal      zeropivot,rs,shiftnz;
1905   PetscReal      shiftpd;
1906   ChShift_Ctx    sctx;
1907   PetscInt       newshift;
1908   PetscTruth     perm_identity;
1909 
1910   PetscFunctionBegin;
1911   shiftnz   = info->shiftnz;
1912   shiftpd   = info->shiftpd;
1913   zeropivot = info->zeropivot;
1914 
1915   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1916   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
1917 
1918   /* initialization */
1919   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1920   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1921   jl   = il + mbs;
1922   rtmp = (MatScalar*)(jl + mbs);
1923 
1924   sctx.shift_amount = 0;
1925   sctx.nshift       = 0;
1926   do {
1927     sctx.chshift = PETSC_FALSE;
1928     for (i=0; i<mbs; i++) {
1929       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1930     }
1931 
1932     for (k = 0; k<mbs; k++){
1933       bval = ba + bi[k];
1934       /* initialize k-th row by the perm[k]-th row of A */
1935       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1936       for (j = jmin; j < jmax; j++){
1937         col = riip[aj[j]];
1938         if (col >= k){ /* only take upper triangular entry */
1939           rtmp[col] = aa[j];
1940           *bval++  = 0.0; /* for in-place factorization */
1941         }
1942       }
1943       /* shift the diagonal of the matrix */
1944       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1945 
1946       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1947       dk = rtmp[k];
1948       i = jl[k]; /* first row to be added to k_th row  */
1949 
1950       while (i < k){
1951         nexti = jl[i]; /* next row to be added to k_th row */
1952 
1953         /* compute multiplier, update diag(k) and U(i,k) */
1954         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1955         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1956         dk += uikdi*ba[ili];
1957         ba[ili] = uikdi; /* -U(i,k) */
1958 
1959         /* add multiple of row i to k-th row */
1960         jmin = ili + 1; jmax = bi[i+1];
1961         if (jmin < jmax){
1962           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1963           /* update il and jl for row i */
1964           il[i] = jmin;
1965           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1966         }
1967         i = nexti;
1968       }
1969 
1970       /* shift the diagonals when zero pivot is detected */
1971       /* compute rs=sum of abs(off-diagonal) */
1972       rs   = 0.0;
1973       jmin = bi[k]+1;
1974       nz   = bi[k+1] - jmin;
1975       bcol = bj + jmin;
1976       while (nz--){
1977         rs += PetscAbsScalar(rtmp[*bcol]);
1978         bcol++;
1979       }
1980 
1981       sctx.rs = rs;
1982       sctx.pv = dk;
1983       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1984 
1985       if (newshift == 1) {
1986         if (!sctx.shift_amount) {
1987           sctx.shift_amount = 1e-5;
1988         }
1989         break;
1990       }
1991 
1992       /* copy data into U(k,:) */
1993       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1994       jmin = bi[k]+1; jmax = bi[k+1];
1995       if (jmin < jmax) {
1996         for (j=jmin; j<jmax; j++){
1997           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1998         }
1999         /* add the k-th row into il and jl */
2000         il[k] = jmin;
2001         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
2002       }
2003     }
2004   } while (sctx.chshift);
2005   ierr = PetscFree(il);CHKERRQ(ierr);
2006 
2007   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
2008   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
2009 
2010   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
2011   if (perm_identity){
2012     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
2013     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
2014     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
2015     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
2016   } else {
2017     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
2018     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
2019     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
2020     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
2021   }
2022 
2023   C->assembled    = PETSC_TRUE;
2024   C->preallocated = PETSC_TRUE;
2025   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
2026   if (sctx.nshift){
2027     if (shiftnz) {
2028       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
2029     } else if (shiftpd) {
2030       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
2031     }
2032   }
2033   PetscFunctionReturn(0);
2034 }
2035 
2036 #undef __FUNCT__
2037 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ_newdatastruct"
2038 PetscErrorCode MatICCFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2039 {
2040   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2041   Mat_SeqSBAIJ       *b;
2042   PetscErrorCode     ierr;
2043   PetscTruth         perm_identity,missing;
2044   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag;
2045   const PetscInt     *rip,*riip;
2046   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
2047   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
2048   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
2049   PetscReal          fill=info->fill,levels=info->levels;
2050   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2051   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
2052   PetscBT            lnkbt;
2053   IS                 iperm;
2054 
2055   PetscFunctionBegin;
2056   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);
2057   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
2058   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
2059   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
2060   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
2061 
2062   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
2063   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr);
2064   ui[0] = 0;
2065 
2066   /* ICC(0) without matrix ordering: simply copies fill pattern */
2067   if (!levels && perm_identity) {
2068 
2069     for (i=0; i<am; i++) {
2070       ui[i+1]  = ui[i] + ai[i+1] - a->diag[i];
2071       udiag[i] = ui[i+1] - 1; /* points to the last entry of U(i,:) */
2072     }
2073     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2074     cols = uj;
2075     for (i=0; i<am; i++) {
2076       aj    = a->j + a->diag[i] + 1; /* 1st entry of U(i,:) without diagonal */
2077       ncols = ui[i+1] - ui[i] - 1;
2078       for (j=0; j<ncols; j++) *cols++ = aj[j];
2079       *cols++ = i; /* diagoanl is located as the last entry of U(i,:) */
2080     }
2081   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
2082 
2083     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
2084     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
2085 
2086     /* initialization */
2087     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
2088 
2089     /* jl: linked list for storing indices of the pivot rows
2090        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2091     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
2092     il         = jl + am;
2093     uj_ptr     = (PetscInt**)(il + am);
2094     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
2095     for (i=0; i<am; i++){
2096       jl[i] = am; il[i] = 0;
2097     }
2098 
2099     /* create and initialize a linked list for storing column indices of the active row k */
2100     nlnk = am + 1;
2101     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
2102 
2103     /* initial FreeSpace size is fill*(ai[am]+1) */
2104     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
2105     current_space = free_space;
2106     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
2107     current_space_lvl = free_space_lvl;
2108 
2109     for (k=0; k<am; k++){  /* for each active row k */
2110       /* initialize lnk by the column indices of row rip[k] of A */
2111       nzk   = 0;
2112       ncols = ai[rip[k]+1] - ai[rip[k]];
2113       if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2114       ncols_upper = 0;
2115       for (j=0; j<ncols; j++){
2116         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
2117         if (riip[i] >= k){ /* only take upper triangular entry */
2118           ajtmp[ncols_upper] = i;
2119           ncols_upper++;
2120         }
2121       }
2122       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
2123       nzk += nlnk;
2124 
2125       /* update lnk by computing fill-in for each pivot row to be merged in */
2126       prow = jl[k]; /* 1st pivot row */
2127 
2128       while (prow < k){
2129         nextprow = jl[prow];
2130 
2131         /* merge prow into k-th row */
2132         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
2133         jmax = ui[prow+1];
2134         ncols = jmax-jmin;
2135         i     = jmin - ui[prow];
2136         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2137         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
2138         j     = *(uj - 1);
2139         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
2140         nzk += nlnk;
2141 
2142         /* update il and jl for prow */
2143         if (jmin < jmax){
2144           il[prow] = jmin;
2145           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
2146         }
2147         prow = nextprow;
2148       }
2149 
2150       /* if free space is not available, make more free space */
2151       if (current_space->local_remaining<nzk) {
2152         i = am - k + 1; /* num of unfactored rows */
2153         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2154         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
2155         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
2156         reallocs++;
2157       }
2158 
2159       /* copy data into free_space and free_space_lvl, then initialize lnk */
2160       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
2161       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
2162 
2163       /* add the k-th row into il and jl */
2164       if (nzk > 1){
2165         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2166         jl[k] = jl[i]; jl[i] = k;
2167         il[k] = ui[k] + 1;
2168       }
2169       uj_ptr[k]     = current_space->array;
2170       uj_lvl_ptr[k] = current_space_lvl->array;
2171 
2172       current_space->array           += nzk;
2173       current_space->local_used      += nzk;
2174       current_space->local_remaining -= nzk;
2175 
2176       current_space_lvl->array           += nzk;
2177       current_space_lvl->local_used      += nzk;
2178       current_space_lvl->local_remaining -= nzk;
2179 
2180       ui[k+1] = ui[k] + nzk;
2181     }
2182 
2183 #if defined(PETSC_USE_INFO)
2184     if (ai[am] != 0) {
2185       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
2186       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
2187       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
2188       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
2189     } else {
2190       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
2191     }
2192 #endif
2193 
2194     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
2195     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
2196     ierr = PetscFree(jl);CHKERRQ(ierr);
2197     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
2198 
2199     /* destroy list of free space and other temporary array(s) */
2200     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2201     ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,am,ui,udiag);CHKERRQ(ierr);
2202     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2203     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
2204 
2205   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
2206 
2207   /* put together the new matrix in MATSEQSBAIJ format */
2208 
2209   b    = (Mat_SeqSBAIJ*)(fact)->data;
2210   b->singlemalloc = PETSC_FALSE;
2211   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
2212   b->j    = uj;
2213   b->i    = ui;
2214   b->diag = udiag;
2215   b->ilen = 0;
2216   b->imax = 0;
2217   b->row  = perm;
2218   b->col  = perm;
2219   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2220   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2221   b->icol = iperm;
2222   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2223   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2224   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2225   b->maxnz   = b->nz = ui[am];
2226   b->free_a  = PETSC_TRUE;
2227   b->free_ij = PETSC_TRUE;
2228 
2229   (fact)->info.factor_mallocs    = reallocs;
2230   (fact)->info.fill_ratio_given  = fill;
2231   if (ai[am] != 0) {
2232     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2233   } else {
2234     (fact)->info.fill_ratio_needed = 0.0;
2235   }
2236   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_newdatastruct;
2237   PetscFunctionReturn(0);
2238 }
2239 
2240 #undef __FUNCT__
2241 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
2242 PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2243 {
2244   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2245   Mat_SeqSBAIJ       *b;
2246   PetscErrorCode     ierr;
2247   PetscTruth         perm_identity,missing;
2248   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag;
2249   const PetscInt     *rip,*riip;
2250   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
2251   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
2252   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
2253   PetscReal          fill=info->fill,levels=info->levels;
2254   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2255   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
2256   PetscBT            lnkbt;
2257   IS                 iperm;
2258   PetscTruth         newdatastruct=PETSC_FALSE;
2259 
2260   PetscFunctionBegin;
2261   ierr = PetscOptionsGetTruth(PETSC_NULL,"-icc_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr);
2262   if(newdatastruct){
2263     ierr = MatICCFactorSymbolic_SeqAIJ_newdatastruct(fact,A,perm,info);CHKERRQ(ierr);
2264     PetscFunctionReturn(0);
2265   }
2266 
2267   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);
2268   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
2269   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
2270   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
2271   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
2272 
2273   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
2274   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr);
2275   ui[0] = 0;
2276 
2277   /* ICC(0) without matrix ordering: simply copies fill pattern */
2278   if (!levels && perm_identity) {
2279 
2280     for (i=0; i<am; i++) {
2281       ui[i+1]  = ui[i] + ai[i+1] - a->diag[i];
2282       udiag[i] = ui[i];
2283     }
2284     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2285     cols = uj;
2286     for (i=0; i<am; i++) {
2287       aj    = a->j + a->diag[i];
2288       ncols = ui[i+1] - ui[i];
2289       for (j=0; j<ncols; j++) *cols++ = *aj++;
2290     }
2291   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
2292     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
2293     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
2294 
2295     /* initialization */
2296     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
2297 
2298     /* jl: linked list for storing indices of the pivot rows
2299        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2300     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
2301     il         = jl + am;
2302     uj_ptr     = (PetscInt**)(il + am);
2303     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
2304     for (i=0; i<am; i++){
2305       jl[i] = am; il[i] = 0;
2306     }
2307 
2308     /* create and initialize a linked list for storing column indices of the active row k */
2309     nlnk = am + 1;
2310     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
2311 
2312     /* initial FreeSpace size is fill*(ai[am]+1) */
2313     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
2314     current_space = free_space;
2315     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
2316     current_space_lvl = free_space_lvl;
2317 
2318     for (k=0; k<am; k++){  /* for each active row k */
2319       /* initialize lnk by the column indices of row rip[k] of A */
2320       nzk   = 0;
2321       ncols = ai[rip[k]+1] - ai[rip[k]];
2322       if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2323       ncols_upper = 0;
2324       for (j=0; j<ncols; j++){
2325         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
2326         if (riip[i] >= k){ /* only take upper triangular entry */
2327           ajtmp[ncols_upper] = i;
2328           ncols_upper++;
2329         }
2330       }
2331       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
2332       nzk += nlnk;
2333 
2334       /* update lnk by computing fill-in for each pivot row to be merged in */
2335       prow = jl[k]; /* 1st pivot row */
2336 
2337       while (prow < k){
2338         nextprow = jl[prow];
2339 
2340         /* merge prow into k-th row */
2341         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
2342         jmax = ui[prow+1];
2343         ncols = jmax-jmin;
2344         i     = jmin - ui[prow];
2345         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2346         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
2347         j     = *(uj - 1);
2348         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
2349         nzk += nlnk;
2350 
2351         /* update il and jl for prow */
2352         if (jmin < jmax){
2353           il[prow] = jmin;
2354           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
2355         }
2356         prow = nextprow;
2357       }
2358 
2359       /* if free space is not available, make more free space */
2360       if (current_space->local_remaining<nzk) {
2361         i = am - k + 1; /* num of unfactored rows */
2362         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2363         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
2364         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
2365         reallocs++;
2366       }
2367 
2368       /* copy data into free_space and free_space_lvl, then initialize lnk */
2369       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
2370       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
2371 
2372       /* add the k-th row into il and jl */
2373       if (nzk > 1){
2374         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2375         jl[k] = jl[i]; jl[i] = k;
2376         il[k] = ui[k] + 1;
2377       }
2378       uj_ptr[k]     = current_space->array;
2379       uj_lvl_ptr[k] = current_space_lvl->array;
2380 
2381       current_space->array           += nzk;
2382       current_space->local_used      += nzk;
2383       current_space->local_remaining -= nzk;
2384 
2385       current_space_lvl->array           += nzk;
2386       current_space_lvl->local_used      += nzk;
2387       current_space_lvl->local_remaining -= nzk;
2388 
2389       ui[k+1] = ui[k] + nzk;
2390     }
2391 
2392 #if defined(PETSC_USE_INFO)
2393     if (ai[am] != 0) {
2394       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
2395       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
2396       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
2397       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
2398     } else {
2399       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
2400     }
2401 #endif
2402 
2403     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
2404     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
2405     ierr = PetscFree(jl);CHKERRQ(ierr);
2406     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
2407 
2408     /* destroy list of free space and other temporary array(s) */
2409     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2410     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
2411     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2412     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
2413 
2414   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
2415 
2416   /* put together the new matrix in MATSEQSBAIJ format */
2417 
2418   b    = (Mat_SeqSBAIJ*)(fact)->data;
2419   b->singlemalloc = PETSC_FALSE;
2420   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
2421   b->j    = uj;
2422   b->i    = ui;
2423   b->diag = udiag;
2424   b->ilen = 0;
2425   b->imax = 0;
2426   b->row  = perm;
2427   b->col  = perm;
2428   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2429   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2430   b->icol = iperm;
2431   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2432   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2433   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2434   b->maxnz   = b->nz = ui[am];
2435   b->free_a  = PETSC_TRUE;
2436   b->free_ij = PETSC_TRUE;
2437 
2438   (fact)->info.factor_mallocs    = reallocs;
2439   (fact)->info.fill_ratio_given  = fill;
2440   if (ai[am] != 0) {
2441     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2442   } else {
2443     (fact)->info.fill_ratio_needed = 0.0;
2444   }
2445   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2446   PetscFunctionReturn(0);
2447 }
2448 
2449 #undef __FUNCT__
2450 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
2451 PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2452 {
2453   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2454   Mat_SeqSBAIJ       *b;
2455   PetscErrorCode     ierr;
2456   PetscTruth         perm_identity;
2457   PetscReal          fill = info->fill;
2458   const PetscInt     *rip,*riip;
2459   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
2460   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
2461   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
2462   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2463   PetscBT            lnkbt;
2464   IS                 iperm;
2465 
2466   PetscFunctionBegin;
2467   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);
2468   /* check whether perm is the identity mapping */
2469   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
2470   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
2471   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
2472   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
2473 
2474   /* initialization */
2475   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
2476   ui[0] = 0;
2477 
2478   /* jl: linked list for storing indices of the pivot rows
2479      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2480   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
2481   il     = jl + am;
2482   cols   = il + am;
2483   ui_ptr = (PetscInt**)(cols + am);
2484   for (i=0; i<am; i++){
2485     jl[i] = am; il[i] = 0;
2486   }
2487 
2488   /* create and initialize a linked list for storing column indices of the active row k */
2489   nlnk = am + 1;
2490   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2491 
2492   /* initial FreeSpace size is fill*(ai[am]+1) */
2493   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
2494   current_space = free_space;
2495 
2496   for (k=0; k<am; k++){  /* for each active row k */
2497     /* initialize lnk by the column indices of row rip[k] of A */
2498     nzk   = 0;
2499     ncols = ai[rip[k]+1] - ai[rip[k]];
2500     if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2501     ncols_upper = 0;
2502     for (j=0; j<ncols; j++){
2503       i = riip[*(aj + ai[rip[k]] + j)];
2504       if (i >= k){ /* only take upper triangular entry */
2505         cols[ncols_upper] = i;
2506         ncols_upper++;
2507       }
2508     }
2509     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2510     nzk += nlnk;
2511 
2512     /* update lnk by computing fill-in for each pivot row to be merged in */
2513     prow = jl[k]; /* 1st pivot row */
2514 
2515     while (prow < k){
2516       nextprow = jl[prow];
2517       /* merge prow into k-th row */
2518       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
2519       jmax = ui[prow+1];
2520       ncols = jmax-jmin;
2521       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2522       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2523       nzk += nlnk;
2524 
2525       /* update il and jl for prow */
2526       if (jmin < jmax){
2527         il[prow] = jmin;
2528         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
2529       }
2530       prow = nextprow;
2531     }
2532 
2533     /* if free space is not available, make more free space */
2534     if (current_space->local_remaining<nzk) {
2535       i = am - k + 1; /* num of unfactored rows */
2536       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2537       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
2538       reallocs++;
2539     }
2540 
2541     /* copy data into free space, then initialize lnk */
2542     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
2543 
2544     /* add the k-th row into il and jl */
2545     if (nzk-1 > 0){
2546       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2547       jl[k] = jl[i]; jl[i] = k;
2548       il[k] = ui[k] + 1;
2549     }
2550     ui_ptr[k] = current_space->array;
2551     current_space->array           += nzk;
2552     current_space->local_used      += nzk;
2553     current_space->local_remaining -= nzk;
2554 
2555     ui[k+1] = ui[k] + nzk;
2556   }
2557 
2558 #if defined(PETSC_USE_INFO)
2559   if (ai[am] != 0) {
2560     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
2561     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
2562     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
2563     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
2564   } else {
2565      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
2566   }
2567 #endif
2568 
2569   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
2570   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
2571   ierr = PetscFree(jl);CHKERRQ(ierr);
2572 
2573   /* destroy list of free space and other temporary array(s) */
2574   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
2575   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
2576   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2577 
2578   /* put together the new matrix in MATSEQSBAIJ format */
2579 
2580   b = (Mat_SeqSBAIJ*)(fact)->data;
2581   b->singlemalloc = PETSC_FALSE;
2582   b->free_a       = PETSC_TRUE;
2583   b->free_ij      = PETSC_TRUE;
2584   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
2585   b->j    = uj;
2586   b->i    = ui;
2587   b->diag = 0;
2588   b->ilen = 0;
2589   b->imax = 0;
2590   b->row  = perm;
2591   b->col  = perm;
2592   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2593   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
2594   b->icol = iperm;
2595   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2596   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2597   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2598   b->maxnz = b->nz = ui[am];
2599 
2600   (fact)->info.factor_mallocs    = reallocs;
2601   (fact)->info.fill_ratio_given  = fill;
2602   if (ai[am] != 0) {
2603     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2604   } else {
2605     (fact)->info.fill_ratio_needed = 0.0;
2606   }
2607   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2608   PetscFunctionReturn(0);
2609 }
2610 
2611 #undef __FUNCT__
2612 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_newdatastruct"
2613 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx)
2614 {
2615   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2616   PetscErrorCode    ierr;
2617   PetscInt          n = A->rmap->n;
2618   const PetscInt    *ai = a->i,*aj = a->j,*vi;
2619   PetscScalar       *x,sum;
2620   const PetscScalar *b;
2621   const MatScalar   *aa = a->a,*v;
2622   PetscInt          i,nz;
2623 
2624   PetscFunctionBegin;
2625   if (!n) PetscFunctionReturn(0);
2626 
2627   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2628   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2629 
2630   /* forward solve the lower triangular */
2631   x[0] = b[0];
2632   v    = aa;
2633   vi   = aj;
2634   for (i=1; i<n; i++) {
2635     nz  = ai[i+1] - ai[i];
2636     sum = b[i];
2637     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
2638     /*    while (nz--) sum -= *v++ * x[*vi++];*/
2639     v  += nz;
2640     vi += nz;
2641     x[i] = sum;
2642   }
2643 
2644   /* backward solve the upper triangular */
2645   v   = aa + ai[n+1];
2646   vi  = aj + ai[n+1];
2647   for (i=n-1; i>=0; i--){
2648     nz = ai[2*n-i +1] - ai[2*n-i]-1;
2649     sum = x[i];
2650     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
2651     /* while (nz--) sum -= *v++ * x[*vi++]; */
2652     v   += nz;
2653     vi  += nz; vi++;
2654     x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */
2655   }
2656 
2657   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
2658   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2659   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2660   PetscFunctionReturn(0);
2661 }
2662 
2663 #undef __FUNCT__
2664 #define __FUNCT__ "MatSolve_SeqAIJ_newdatastruct"
2665 PetscErrorCode MatSolve_SeqAIJ_newdatastruct(Mat A,Vec bb,Vec xx)
2666 {
2667   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2668   IS                iscol = a->col,isrow = a->row;
2669   PetscErrorCode    ierr;
2670   PetscInt          i,n=A->rmap->n,*vi,*ai=a->i,*aj=a->j,nz,k;
2671   const PetscInt    *rout,*cout,*r,*c;
2672   PetscScalar       *x,*tmp,*tmps,sum;
2673   const PetscScalar *b;
2674   const MatScalar   *aa = a->a,*v;
2675 
2676   PetscFunctionBegin;
2677   if (!n) PetscFunctionReturn(0);
2678 
2679   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2680   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2681   tmp  = a->solve_work;
2682 
2683   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
2684   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
2685 
2686   /* forward solve the lower triangular */
2687   tmp[0] = b[*r++];
2688   tmps   = tmp;
2689   v      = aa;
2690   vi     = aj;
2691   for (i=1; i<n; i++) {
2692     nz  = ai[i+1] - ai[i];
2693     sum = b[*r++];
2694     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
2695     tmp[i] = sum;
2696     v += nz; vi += nz;
2697   }
2698 
2699   /* backward solve the upper triangular */
2700   k  = n+1;
2701   v  = aa + ai[k]; /* 1st entry of U(n-1,:) */
2702   vi = aj + ai[k];
2703   for (i=n-1; i>=0; i--){
2704     k  = 2*n-i;
2705     nz = ai[k +1] - ai[k] - 1;
2706     sum = tmp[i];
2707     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
2708     x[c[i]] = tmp[i] = sum*v[nz]; /* v[nz] = aa[adiag[i]] */
2709     v += nz+1; vi += nz+1;
2710   }
2711 
2712   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
2713   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
2714   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
2715   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2716   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
2717   PetscFunctionReturn(0);
2718 }
2719 
2720 #undef __FUNCT__
2721 #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2722 PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
2723 {
2724   Mat                B = *fact;
2725   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b;
2726   IS                 isicol;
2727   PetscErrorCode     ierr;
2728   const PetscInt     *r,*ic;
2729   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
2730   PetscInt           *bi,*bj,*bdiag,*bdiag_rev;
2731   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,nzi_al,nzi_au;
2732   PetscInt           nlnk,*lnk;
2733   PetscBT            lnkbt;
2734   PetscTruth         row_identity,icol_identity,both_identity;
2735   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp;
2736   const PetscInt     *ics;
2737   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
2738   PetscReal          dt=info->dt,dtcol=info->dtcol,shift=info->shiftinblocks;
2739   PetscInt           dtcount=(PetscInt)info->dtcount,nnz_max;
2740   PetscTruth         missing;
2741 
2742   PetscFunctionBegin;
2743 
2744   if (dt      == PETSC_DEFAULT) dt      = 0.005;
2745   if (dtcol   == PETSC_DEFAULT) dtcol   = 0.01; /* XXX unused! */
2746   if (dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax);
2747 
2748   /* ------- symbolic factorization, can be reused ---------*/
2749   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2750   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2751   adiag=a->diag;
2752 
2753   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2754 
2755   /* bdiag is location of diagonal in factor */
2756   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2757   bdiag_rev = bdiag + n+1;
2758 
2759   /* allocate row pointers bi */
2760   ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2761 
2762   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
2763   if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */
2764   nnz_max  = ai[n]+2*n*dtcount+2;
2765 
2766   ierr = PetscMalloc((nnz_max+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
2767   ierr = PetscMalloc((nnz_max+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
2768 
2769   /* put together the new matrix */
2770   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2771   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
2772   b    = (Mat_SeqAIJ*)(B)->data;
2773   b->free_a       = PETSC_TRUE;
2774   b->free_ij      = PETSC_TRUE;
2775   b->singlemalloc = PETSC_FALSE;
2776   b->a          = ba;
2777   b->j          = bj;
2778   b->i          = bi;
2779   b->diag       = bdiag;
2780   b->ilen       = 0;
2781   b->imax       = 0;
2782   b->row        = isrow;
2783   b->col        = iscol;
2784   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
2785   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
2786   b->icol       = isicol;
2787   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
2788 
2789   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
2790   b->maxnz = nnz_max;
2791 
2792   (B)->factor                = MAT_FACTOR_ILUDT;
2793   (B)->info.factor_mallocs   = 0;
2794   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]);
2795   CHKMEMQ;
2796   /* ------- end of symbolic factorization ---------*/
2797 
2798   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2799   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
2800   ics  = ic;
2801 
2802   /* linked list for storing column indices of the active row */
2803   nlnk = n + 1;
2804   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2805 
2806   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
2807   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
2808   jtmp = im + n;
2809   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
2810   ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2811   ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
2812   vtmp = rtmp + n;
2813 
2814   bi[0]    = 0;
2815   bdiag[0] = nnz_max-1; /* location of diag[0] in factor B */
2816   bdiag_rev[n] = bdiag[0];
2817   bi[2*n+1] = bdiag[0]+1; /* endof bj and ba array */
2818   for (i=0; i<n; i++) {
2819     /* copy initial fill into linked list */
2820     nzi = 0; /* nonzeros for active row i */
2821     nzi = ai[r[i]+1] - ai[r[i]];
2822     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
2823     nzi_al = adiag[r[i]] - ai[r[i]];
2824     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
2825     ajtmp = aj + ai[r[i]];
2826     ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2827 
2828     /* load in initial (unfactored row) */
2829     aatmp = a->a + ai[r[i]];
2830     for (j=0; j<nzi; j++) {
2831       rtmp[ics[*ajtmp++]] = *aatmp++;
2832     }
2833 
2834     /* add pivot rows into linked list */
2835     row = lnk[n];
2836     while (row < i ) {
2837       nzi_bl = bi[row+1] - bi[row] + 1;
2838       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
2839       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
2840       nzi  += nlnk;
2841       row   = lnk[row];
2842     }
2843 
2844     /* copy data from lnk into jtmp, then initialize lnk */
2845     ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
2846 
2847     /* numerical factorization */
2848     bjtmp = jtmp;
2849     row   = *bjtmp++; /* 1st pivot row */
2850     while  ( row < i ) {
2851       pc         = rtmp + row;
2852       pv         = ba + bdiag[row]; /* 1./(diag of the pivot row) */
2853       multiplier = (*pc) * (*pv);
2854       *pc        = multiplier;
2855       if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */
2856         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
2857         pv         = ba + bdiag[row+1] + 1;
2858         /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */
2859         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
2860         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
2861         ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
2862       }
2863       row = *bjtmp++;
2864     }
2865 
2866     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
2867     diag_tmp = rtmp[i];  /* save diagonal value - may not needed?? */
2868     nzi_bl = 0; j = 0;
2869     while (jtmp[j] < i){ /* Note: jtmp is sorted */
2870       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2871       nzi_bl++; j++;
2872     }
2873     nzi_bu = nzi - nzi_bl -1;
2874     while (j < nzi){
2875       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
2876       j++;
2877     }
2878 
2879     bjtmp = bj + bi[i];
2880     batmp = ba + bi[i];
2881     /* apply level dropping rule to L part */
2882     ncut = nzi_al + dtcount;
2883     if (ncut < nzi_bl){
2884       ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr);
2885       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
2886     } else {
2887       ncut = nzi_bl;
2888     }
2889     for (j=0; j<ncut; j++){
2890       bjtmp[j] = jtmp[j];
2891       batmp[j] = vtmp[j];
2892       /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */
2893     }
2894     bi[i+1] = bi[i] + ncut;
2895     nzi = ncut + 1;
2896 
2897     /* apply level dropping rule to U part */
2898     ncut = nzi_au + dtcount;
2899     if (ncut < nzi_bu){
2900       ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
2901       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
2902     } else {
2903       ncut = nzi_bu;
2904     }
2905     nzi += ncut;
2906 
2907     /* mark bdiagonal */
2908     bdiag[i+1]       = bdiag[i] - (ncut + 1);
2909     bdiag_rev[n-i-1] = bdiag[i+1];
2910     bi[2*n - i]      = bi[2*n - i +1] - (ncut + 1);
2911     bjtmp = bj + bdiag[i];
2912     batmp = ba + bdiag[i];
2913     *bjtmp = i;
2914     *batmp = diag_tmp; /* rtmp[i]; */
2915     if (*batmp == 0.0) {
2916       *batmp = dt+shift;
2917       /* printf(" row %d add shift %g\n",i,shift); */
2918     }
2919     *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */
2920     /* printf(" (%d,%g),",*bjtmp,*batmp); */
2921 
2922     bjtmp = bj + bdiag[i+1]+1;
2923     batmp = ba + bdiag[i+1]+1;
2924     for (k=0; k<ncut; k++){
2925       bjtmp[k] = jtmp[nzi_bl+1+k];
2926       batmp[k] = vtmp[nzi_bl+1+k];
2927       /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */
2928     }
2929     /* printf("\n"); */
2930 
2931     im[i]   = nzi; /* used by PetscLLAddSortedLU() */
2932     /*
2933     printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]);
2934     printf(" ----------------------------\n");
2935     */
2936   } /* for (i=0; i<n; i++) */
2937   /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */
2938   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]);
2939 
2940   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2941   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2942 
2943   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
2944   ierr = PetscFree(im);CHKERRQ(ierr);
2945   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2946 
2947   ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr);
2948   b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n];
2949 
2950   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
2951   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
2952   both_identity = (PetscTruth) (row_identity && icol_identity);
2953   if (row_identity && icol_identity) {
2954     B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_newdatastruct;
2955   } else {
2956     B->ops->solve = MatSolve_SeqAIJ_newdatastruct;
2957   }
2958 
2959   B->ops->lufactorsymbolic  = MatILUDTFactorSymbolic_SeqAIJ;
2960   B->ops->lufactornumeric   = MatILUDTFactorNumeric_SeqAIJ;
2961   B->ops->solveadd          = 0;
2962   B->ops->solvetranspose    = 0;
2963   B->ops->solvetransposeadd = 0;
2964   B->ops->matsolve          = 0;
2965   B->assembled              = PETSC_TRUE;
2966   B->preallocated           = PETSC_TRUE;
2967   PetscFunctionReturn(0);
2968 }
2969 
2970 /* a wraper of MatILUDTFactor_SeqAIJ() */
2971 #undef __FUNCT__
2972 #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ"
2973 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
2974 {
2975   PetscErrorCode     ierr;
2976 
2977   PetscFunctionBegin;
2978   ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr);
2979 
2980   fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ;
2981   PetscFunctionReturn(0);
2982 }
2983 
2984 /*
2985    same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors
2986    - intend to replace existing MatLUFactorNumeric_SeqAIJ()
2987 */
2988 #undef __FUNCT__
2989 #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ"
2990 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info)
2991 {
2992   Mat            C=fact;
2993   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
2994   IS             isrow = b->row,isicol = b->icol;
2995   PetscErrorCode ierr;
2996   const PetscInt *r,*ic,*ics;
2997   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
2998   PetscInt       *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj;
2999   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
3000   PetscReal      dt=info->dt,shift=info->shiftinblocks;
3001   PetscTruth     row_identity, col_identity;
3002 
3003   PetscFunctionBegin;
3004   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
3005   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
3006   ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
3007   ics  = ic;
3008 
3009   for (i=0; i<n; i++){
3010     /* initialize rtmp array */
3011     nzl   = bi[i+1] - bi[i];       /* num of nozeros in L(i,:) */
3012     bjtmp = bj + bi[i];
3013     for  (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0;
3014     rtmp[i] = 0.0;
3015     nzu   = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */
3016     bjtmp = bj + bdiag[i+1] + 1;
3017     for  (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0;
3018 
3019     /* load in initial unfactored row of A */
3020     /* printf("row %d\n",i); */
3021     nz    = ai[r[i]+1] - ai[r[i]];
3022     ajtmp = aj + ai[r[i]];
3023     v     = aa + ai[r[i]];
3024     for (j=0; j<nz; j++) {
3025       rtmp[ics[*ajtmp++]] = v[j];
3026       /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */
3027     }
3028     /* printf("\n"); */
3029 
3030     /* numerical factorization */
3031     bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */
3032     nzl   = bi[i+1] - bi[i]; /* num of entries in L(i,:) */
3033     k = 0;
3034     while (k < nzl){
3035       row   = *bjtmp++;
3036       /* printf("  prow %d\n",row); */
3037       pc         = rtmp + row;
3038       pv         = b->a + bdiag[row]; /* 1./(diag of the pivot row) */
3039       multiplier = (*pc) * (*pv);
3040       *pc        = multiplier;
3041       if (PetscAbsScalar(multiplier) > dt){
3042         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
3043         pv         = b->a + bdiag[row+1] + 1;
3044         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
3045         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
3046         /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */
3047       }
3048       k++;
3049     }
3050 
3051     /* finished row so stick it into b->a */
3052     /* L-part */
3053     pv = b->a + bi[i] ;
3054     pj = bj + bi[i] ;
3055     nzl = bi[i+1] - bi[i];
3056     for (j=0; j<nzl; j++) {
3057       pv[j] = rtmp[pj[j]];
3058       /* printf(" (%d,%g),",pj[j],pv[j]); */
3059     }
3060 
3061     /* diagonal: invert diagonal entries for simplier triangular solves */
3062     if (rtmp[i] == 0.0) rtmp[i] = dt+shift;
3063     b->a[bdiag[i]] = 1.0/rtmp[i];
3064     /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */
3065 
3066     /* U-part */
3067     pv = b->a + bdiag[i+1] + 1;
3068     pj = bj + bdiag[i+1] + 1;
3069     nzu = bdiag[i] - bdiag[i+1] - 1;
3070     for (j=0; j<nzu; j++) {
3071       pv[j] = rtmp[pj[j]];
3072       /* printf(" (%d,%g),",pj[j],pv[j]); */
3073     }
3074     /* printf("\n"); */
3075   }
3076 
3077   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3078   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3079   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
3080 
3081   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
3082   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
3083   if (row_identity && col_identity) {
3084     C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering_newdatastruct;
3085   } else {
3086     C->ops->solve   = MatSolve_SeqAIJ_newdatastruct;
3087   }
3088   C->ops->solveadd           = 0;
3089   C->ops->solvetranspose     = 0;
3090   C->ops->solvetransposeadd  = 0;
3091   C->ops->matsolve           = 0;
3092   C->assembled    = PETSC_TRUE;
3093   C->preallocated = PETSC_TRUE;
3094   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
3095   PetscFunctionReturn(0);
3096 }
3097