xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision d4e9f3b6443ec535fe5e409f785fa3cf019ceb6c)
1 /*$Id: aijfact.c,v 1.167 2001/09/11 16:32:26 bsmith Exp $*/
2 
3 #include "src/mat/impls/aij/seq/aij.h"
4 #include "src/vec/vecimpl.h"
5 #include "src/inline/dot.h"
6 #include "src/inline/spops.h"
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
10 int MatOrdering_Flow_SeqAIJ(Mat mat,MatOrderingType type,IS *irow,IS *icol)
11 {
12   PetscFunctionBegin;
13 
14   SETERRQ(PETSC_ERR_SUP,"Code not written");
15 #if !defined(PETSC_USE_DEBUG)
16   PetscFunctionReturn(0);
17 #endif
18 }
19 
20 
21 EXTERN int MatMarkDiagonal_SeqAIJ(Mat);
22 EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
23 
24 EXTERN int SPARSEKIT2dperm(int*,PetscScalar*,int*,int*,PetscScalar*,int*,int*,int*,int*,int*);
25 EXTERN int SPARSEKIT2ilutp(int*,PetscScalar*,int*,int*,int*,PetscReal,PetscReal*,int*,PetscScalar*,int*,int*,int*,PetscScalar*,int*,int*,int*);
26 EXTERN int SPARSEKIT2msrcsr(int*,PetscScalar*,int*,PetscScalar*,int*,int*,PetscScalar*,int*);
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
30   /* ------------------------------------------------------------
31 
32           This interface was contribed by Tony Caola
33 
34      This routine is an interface to the pivoting drop-tolerance
35      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
36      SPARSEKIT2.
37 
38      The SPARSEKIT2 routines used here are covered by the GNU
39      copyright; see the file gnu in this directory.
40 
41      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
42      help in getting this routine ironed out.
43 
44      The major drawback to this routine is that if info->fill is
45      not large enough it fails rather than allocating more space;
46      this can be fixed by hacking/improving the f2c version of
47      Yousef Saad's code.
48 
49      ishift = 0, for indices start at 1
50      ishift = 1, for indices starting at 0
51      ------------------------------------------------------------
52 */
53 int MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
54 {
55 #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
56   PetscFunctionBegin;
57   SETERRQ(1,"This distribution does not include GNU Copyright code\n\
58   You can obtain the drop tolerance routines by installing PETSc from\n\
59   www.mcs.anl.gov/petsc\n");
60 #else
61   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b;
62   IS           iscolf,isicol,isirow;
63   PetscTruth   reorder;
64   int          *c,*r,*ic,ierr,i,n = A->m;
65   int          *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
66   int          *ordcol,*iwk,*iperm,*jw;
67   int          ishift = !a->indexshift;
68   int          jmax,lfill,job,*o_i,*o_j;
69   PetscScalar  *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
70   PetscReal    permtol,af;
71 
72   PetscFunctionBegin;
73 
74   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
75   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax);
76   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
77   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
78   lfill   = (int)(info->dtcount/2.0);
79   jmax    = (int)(info->fill*a->nz);
80   permtol = info->dtcol;
81 
82 
83   /* ------------------------------------------------------------
84      If reorder=.TRUE., then the original matrix has to be
85      reordered to reflect the user selected ordering scheme, and
86      then de-reordered so it is in it's original format.
87      Because Saad's dperm() is NOT in place, we have to copy
88      the original matrix and allocate more storage. . .
89      ------------------------------------------------------------
90   */
91 
92   /* set reorder to true if either isrow or iscol is not identity */
93   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
94   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
95   reorder = PetscNot(reorder);
96 
97 
98   /* storage for ilu factor */
99   ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr);
100   ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr);
101   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
102   ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr);
103 
104   /* ------------------------------------------------------------
105      Make sure that everything is Fortran formatted (1-Based)
106      ------------------------------------------------------------
107   */
108   if (ishift) {
109     for (i=old_i[0];i<old_i[n];i++) {
110       old_j[i]++;
111     }
112     for(i=0;i<n+1;i++) {
113       old_i[i]++;
114     };
115   };
116 
117   if (reorder) {
118     ierr = ISGetIndices(iscol,&c);           CHKERRQ(ierr);
119     ierr = ISGetIndices(isrow,&r);           CHKERRQ(ierr);
120     for(i=0;i<n;i++) {
121       r[i]  = r[i]+1;
122       c[i]  = c[i]+1;
123     }
124     ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr);
125     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr);
126     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
127     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
128     for (i=0;i<n;i++) {
129       r[i]  = r[i]-1;
130       c[i]  = c[i]-1;
131     }
132     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
133     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
134     o_a = old_a2;
135     o_j = old_j2;
136     o_i = old_i2;
137   } else {
138     o_a = old_a;
139     o_j = old_j;
140     o_i = old_i;
141   }
142 
143   /* ------------------------------------------------------------
144      Call Saad's ilutp() routine to generate the factorization
145      ------------------------------------------------------------
146   */
147 
148   ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr);
149   ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr);
150   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
151 
152   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&permtol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&ierr);
153   if (ierr) {
154     switch (ierr) {
155       case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
156       case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
157       case -5: SETERRQ(1,"ilutp(), zero row encountered");
158       case -1: SETERRQ(1,"ilutp(), input matrix may be wrong");
159       case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax);
160       default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr);
161     }
162   }
163 
164   ierr = PetscFree(w);CHKERRQ(ierr);
165   ierr = PetscFree(jw);CHKERRQ(ierr);
166 
167   /* ------------------------------------------------------------
168      Saad's routine gives the result in Modified Sparse Row (msr)
169      Convert to Compressed Sparse Row format (csr)
170      ------------------------------------------------------------
171   */
172 
173   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
174   ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr);
175 
176   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
177 
178   ierr = PetscFree(iwk);CHKERRQ(ierr);
179   ierr = PetscFree(wk);CHKERRQ(ierr);
180 
181   if (reorder) {
182     ierr = PetscFree(old_a2);CHKERRQ(ierr);
183     ierr = PetscFree(old_j2);CHKERRQ(ierr);
184     ierr = PetscFree(old_i2);CHKERRQ(ierr);
185   } else {
186     /* fix permutation of old_j that the factorization introduced */
187     for (i=old_i[0]; i<old_i[n]; i++) {
188       old_j[i-1] = iperm[old_j[i-1]-1];
189     }
190   }
191 
192   /* get rid of the shift to indices starting at 1 */
193   if (ishift) {
194     for (i=0; i<n+1; i++) {
195       old_i[i]--;
196     }
197     for (i=old_i[0];i<old_i[n];i++) {
198       old_j[i]--;
199     }
200   }
201 
202   /* Make the factored matrix 0-based */
203   if (ishift) {
204     for (i=0; i<n+1; i++) {
205       new_i[i]--;
206     }
207     for (i=new_i[0];i<new_i[n];i++) {
208       new_j[i]--;
209     }
210   }
211 
212   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
213   /*-- permute the right-hand-side and solution vectors           --*/
214   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
215   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
216   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
217   for(i=0; i<n; i++) {
218     ordcol[i] = ic[iperm[i]-1];
219   };
220   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
221   ierr = ISDestroy(isicol);CHKERRQ(ierr);
222 
223   ierr = PetscFree(iperm);CHKERRQ(ierr);
224 
225   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
226   ierr = PetscFree(ordcol);CHKERRQ(ierr);
227 
228   /*----- put together the new matrix -----*/
229 
230   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
231   (*fact)->factor    = FACTOR_LU;
232   (*fact)->assembled = PETSC_TRUE;
233 
234   b = (Mat_SeqAIJ*)(*fact)->data;
235   ierr = PetscFree(b->imax);CHKERRQ(ierr);
236   b->sorted        = PETSC_FALSE;
237   b->singlemalloc  = PETSC_FALSE;
238   /* the next line frees the default space generated by the MatCreate() */
239   ierr             = PetscFree(b->a);CHKERRQ(ierr);
240   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
241   b->a             = new_a;
242   b->j             = new_j;
243   b->i             = new_i;
244   b->ilen          = 0;
245   b->imax          = 0;
246   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
247   b->row           = isirow;
248   b->col           = iscolf;
249   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
250   b->maxnz = b->nz = new_i[n];
251   b->indexshift    = a->indexshift;
252   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
253   (*fact)->info.factor_mallocs = 0;
254 
255   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
256 
257   /* check out for identical nodes. If found, use inode functions */
258   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
259 
260   af = ((double)b->nz)/((double)a->nz) + .001;
261   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
262   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
263   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
264   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
265 
266   PetscFunctionReturn(0);
267 #endif
268 }
269 
270 /*
271     Factorization code for AIJ format.
272 */
273 #undef __FUNCT__
274 #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
275 int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
276 {
277   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
278   IS         isicol;
279   int        *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j;
280   int        *ainew,*ajnew,jmax,*fill,*ajtmp,nz,shift = a->indexshift;
281   int        *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im;
282   PetscReal  f;
283 
284   PetscFunctionBegin;
285   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
286   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
287   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
288   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
289 
290   /* get new row pointers */
291   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
292   ainew[0] = -shift;
293   /* don't know how many column pointers are needed so estimate */
294   f = info->fill;
295   jmax  = (int)(f*ai[n]+(!shift));
296   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
297   /* fill is a linked list of nonzeros in active row */
298   ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr);
299   im   = fill + n + 1;
300   /* idnew is location of diagonal in factor */
301   ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr);
302   idnew[0] = -shift;
303 
304   for (i=0; i<n; i++) {
305     /* first copy previous fill into linked list */
306     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
307     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
308     ajtmp   = aj + ai[r[i]] + shift;
309     fill[n] = n;
310     while (nz--) {
311       fm  = n;
312       idx = ic[*ajtmp++ + shift];
313       do {
314         m  = fm;
315         fm = fill[m];
316       } while (fm < idx);
317       fill[m]   = idx;
318       fill[idx] = fm;
319     }
320     row = fill[n];
321     while (row < i) {
322       ajtmp = ajnew + idnew[row] + (!shift);
323       nzbd  = 1 + idnew[row] - ainew[row];
324       nz    = im[row] - nzbd;
325       fm    = row;
326       while (nz-- > 0) {
327         idx = *ajtmp++ + shift;
328         nzbd++;
329         if (idx == i) im[row] = nzbd;
330         do {
331           m  = fm;
332           fm = fill[m];
333         } while (fm < idx);
334         if (fm != idx) {
335           fill[m]   = idx;
336           fill[idx] = fm;
337           fm        = idx;
338           nnz++;
339         }
340       }
341       row = fill[row];
342     }
343     /* copy new filled row into permanent storage */
344     ainew[i+1] = ainew[i] + nnz;
345     if (ainew[i+1] > jmax) {
346 
347       /* estimate how much additional space we will need */
348       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
349       /* just double the memory each time */
350       int maxadd = jmax;
351       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
352       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
353       jmax += maxadd;
354 
355       /* allocate a longer ajnew */
356       ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr);
357       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i]+shift)*sizeof(int));CHKERRQ(ierr);
358       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
359       ajnew = ajtmp;
360       realloc++; /* count how many times we realloc */
361     }
362     ajtmp = ajnew + ainew[i] + shift;
363     fm    = fill[n];
364     nzi   = 0;
365     im[i] = nnz;
366     while (nnz--) {
367       if (fm < i) nzi++;
368       *ajtmp++ = fm - shift;
369       fm       = fill[fm];
370     }
371     idnew[i] = ainew[i] + nzi;
372   }
373   if (ai[n] != 0) {
374     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
375     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
376     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
377     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
378     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
379   } else {
380     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
381   }
382 
383   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
384   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
385 
386   ierr = PetscFree(fill);CHKERRQ(ierr);
387 
388   /* put together the new matrix */
389   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr);
390   PetscLogObjectParent(*B,isicol);
391   b = (Mat_SeqAIJ*)(*B)->data;
392   ierr = PetscFree(b->imax);CHKERRQ(ierr);
393   b->singlemalloc = PETSC_FALSE;
394   /* the next line frees the default space generated by the Create() */
395   ierr = PetscFree(b->a);CHKERRQ(ierr);
396   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
397   ierr          = PetscMalloc((ainew[n]+shift+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
398   b->j          = ajnew;
399   b->i          = ainew;
400   b->diag       = idnew;
401   b->ilen       = 0;
402   b->imax       = 0;
403   b->row        = isrow;
404   b->col        = iscol;
405   b->lu_damping   = info->damping;
406   b->lu_zeropivot = info->zeropivot;
407   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
408   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
409   b->icol       = isicol;
410   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
411   /* In b structure:  Free imax, ilen, old a, old j.
412      Allocate idnew, solve_work, new a, new j */
413   PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(PetscScalar)));
414   b->maxnz = b->nz = ainew[n] + shift;
415 
416   (*B)->factor                 =  FACTOR_LU;
417   (*B)->info.factor_mallocs    = realloc;
418   (*B)->info.fill_ratio_given  = f;
419   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
420   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
421 
422   if (ai[n] != 0) {
423     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
424   } else {
425     (*B)->info.fill_ratio_needed = 0.0;
426   }
427   PetscFunctionReturn(0);
428 }
429 /* ----------------------------------------------------------- */
430 EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
431 
432 #undef __FUNCT__
433 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
434 int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B)
435 {
436   Mat          C = *B;
437   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data;
438   IS           isrow = b->row,isicol = b->icol;
439   int          *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j;
440   int          *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift;
441   int          *diag_offset = b->diag,diag;
442   int          *pj,ndamp = 0;
443   PetscScalar  *rtmp,*v,*pc,multiplier,*pv,*rtmps;
444   PetscReal    damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs;
445   PetscTruth   damp;
446 
447   PetscFunctionBegin;
448   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
449   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
450   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
451   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
452   rtmps = rtmp + shift; ics = ic + shift;
453 
454   do {
455     damp = PETSC_FALSE;
456     for (i=0; i<n; i++) {
457       nz    = ai[i+1] - ai[i];
458       ajtmp = aj + ai[i] + shift;
459       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
460 
461       /* load in initial (unfactored row) */
462       nz       = a->i[r[i]+1] - a->i[r[i]];
463       ajtmpold = a->j + a->i[r[i]] + shift;
464       v        = a->a + a->i[r[i]] + shift;
465       for (j=0; j<nz; j++) {
466         rtmp[ics[ajtmpold[j]]] = v[j];
467         if (ndamp && ajtmpold[j] == r[i]) { /* damp the diagonal of the matrix */
468           rtmp[ics[ajtmpold[j]]] += damping;
469         }
470       }
471 
472       row = *ajtmp++ + shift;
473       while  (row < i) {
474         pc = rtmp + row;
475         if (*pc != 0.0) {
476           pv         = b->a + diag_offset[row] + shift;
477           pj         = b->j + diag_offset[row] + (!shift);
478           multiplier = *pc / *pv++;
479           *pc        = multiplier;
480           nz         = ai[row+1] - diag_offset[row] - 1;
481           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
482           PetscLogFlops(2*nz);
483         }
484         row = *ajtmp++ + shift;
485       }
486       /* finished row so stick it into b->a */
487       pv   = b->a + ai[i] + shift;
488       pj   = b->j + ai[i] + shift;
489       nz   = ai[i+1] - ai[i];
490       diag = diag_offset[i] - ai[i];
491       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
492       rs   = 0.0;
493       for (j=0; j<nz; j++) {
494         pv[j] = rtmps[pj[j]];
495         if (j != diag) rs += PetscAbsScalar(pv[j]);
496       }
497       if (PetscAbsScalar(pv[diag]) < zeropivot*rs) {
498         if (damping) {
499           if (ndamp) damping *= 2.0;
500           damp = PETSC_TRUE;
501           ndamp++;
502           break;
503         } else {
504           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
505         }
506       }
507     }
508   } while (damp);
509 
510   /* invert diagonal entries for simplier triangular solves */
511   for (i=0; i<n; i++) {
512     b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift];
513   }
514 
515   ierr = PetscFree(rtmp);CHKERRQ(ierr);
516   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
517   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
518   C->factor = FACTOR_LU;
519   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
520   C->assembled = PETSC_TRUE;
521   PetscLogFlops(C->n);
522   if (ndamp) {
523     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping);
524   }
525   PetscFunctionReturn(0);
526 }
527 /* ----------------------------------------------------------- */
528 #undef __FUNCT__
529 #define __FUNCT__ "MatLUFactor_SeqAIJ"
530 int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
531 {
532   int ierr;
533   Mat C;
534 
535   PetscFunctionBegin;
536   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
537   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
538   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
539   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
540   PetscFunctionReturn(0);
541 }
542 /* ----------------------------------------------------------- */
543 #undef __FUNCT__
544 #define __FUNCT__ "MatSolve_SeqAIJ"
545 int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
546 {
547   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
548   IS           iscol = a->col,isrow = a->row;
549   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
550   int          nz,shift = a->indexshift,*rout,*cout;
551   PetscScalar  *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
552 
553   PetscFunctionBegin;
554   if (!n) PetscFunctionReturn(0);
555 
556   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
557   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
558   tmp  = a->solve_work;
559 
560   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
561   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
562 
563   /* forward solve the lower triangular */
564   tmp[0] = b[*r++];
565   tmps   = tmp + shift;
566   for (i=1; i<n; i++) {
567     v   = aa + ai[i] + shift;
568     vi  = aj + ai[i] + shift;
569     nz  = a->diag[i] - ai[i];
570     sum = b[*r++];
571     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
572     tmp[i] = sum;
573   }
574 
575   /* backward solve the upper triangular */
576   for (i=n-1; i>=0; i--){
577     v   = aa + a->diag[i] + (!shift);
578     vi  = aj + a->diag[i] + (!shift);
579     nz  = ai[i+1] - a->diag[i] - 1;
580     sum = tmp[i];
581     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
582     x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift];
583   }
584 
585   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
586   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
587   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
588   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
589   PetscLogFlops(2*a->nz - A->n);
590   PetscFunctionReturn(0);
591 }
592 
593 /* ----------------------------------------------------------- */
594 #undef __FUNCT__
595 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
596 int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
597 {
598   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
599   int          n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr;
600   PetscScalar  *x,*b,*aa = a->a;
601 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
602   int          adiag_i,i,*vi,nz,ai_i;
603   PetscScalar  *v,sum;
604 #endif
605 
606   PetscFunctionBegin;
607   if (!n) PetscFunctionReturn(0);
608   if (a->indexshift) {
609      ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr);
610      PetscFunctionReturn(0);
611   }
612 
613   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
614   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
615 
616 #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
617   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
618 #else
619   /* forward solve the lower triangular */
620   x[0] = b[0];
621   for (i=1; i<n; i++) {
622     ai_i = ai[i];
623     v    = aa + ai_i;
624     vi   = aj + ai_i;
625     nz   = adiag[i] - ai_i;
626     sum  = b[i];
627     while (nz--) sum -= *v++ * x[*vi++];
628     x[i] = sum;
629   }
630 
631   /* backward solve the upper triangular */
632   for (i=n-1; i>=0; i--){
633     adiag_i = adiag[i];
634     v       = aa + adiag_i + 1;
635     vi      = aj + adiag_i + 1;
636     nz      = ai[i+1] - adiag_i - 1;
637     sum     = x[i];
638     while (nz--) sum -= *v++ * x[*vi++];
639     x[i]    = sum*aa[adiag_i];
640   }
641 #endif
642   PetscLogFlops(2*a->nz - A->n);
643   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
644   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
645   PetscFunctionReturn(0);
646 }
647 
648 #undef __FUNCT__
649 #define __FUNCT__ "MatSolveAdd_SeqAIJ"
650 int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
651 {
652   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
653   IS           iscol = a->col,isrow = a->row;
654   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
655   int          nz,shift = a->indexshift,*rout,*cout;
656   PetscScalar  *x,*b,*tmp,*aa = a->a,sum,*v;
657 
658   PetscFunctionBegin;
659   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
660 
661   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
662   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
663   tmp  = a->solve_work;
664 
665   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
666   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
667 
668   /* forward solve the lower triangular */
669   tmp[0] = b[*r++];
670   for (i=1; i<n; i++) {
671     v   = aa + ai[i] + shift;
672     vi  = aj + ai[i] + shift;
673     nz  = a->diag[i] - ai[i];
674     sum = b[*r++];
675     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
676     tmp[i] = sum;
677   }
678 
679   /* backward solve the upper triangular */
680   for (i=n-1; i>=0; i--){
681     v   = aa + a->diag[i] + (!shift);
682     vi  = aj + a->diag[i] + (!shift);
683     nz  = ai[i+1] - a->diag[i] - 1;
684     sum = tmp[i];
685     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
686     tmp[i] = sum*aa[a->diag[i]+shift];
687     x[*c--] += tmp[i];
688   }
689 
690   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
691   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
692   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
693   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
694   PetscLogFlops(2*a->nz);
695 
696   PetscFunctionReturn(0);
697 }
698 /* -------------------------------------------------------------------*/
699 #undef __FUNCT__
700 #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
701 int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
702 {
703   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
704   IS           iscol = a->col,isrow = a->row;
705   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
706   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
707   PetscScalar  *x,*b,*tmp,*aa = a->a,*v,s1;
708 
709   PetscFunctionBegin;
710   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
711   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
712   tmp  = a->solve_work;
713 
714   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
715   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
716 
717   /* copy the b into temp work space according to permutation */
718   for (i=0; i<n; i++) tmp[i] = b[c[i]];
719 
720   /* forward solve the U^T */
721   for (i=0; i<n; i++) {
722     v   = aa + diag[i] + shift;
723     vi  = aj + diag[i] + (!shift);
724     nz  = ai[i+1] - diag[i] - 1;
725     s1  = tmp[i];
726     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
727     while (nz--) {
728       tmp[*vi++ + shift] -= (*v++)*s1;
729     }
730     tmp[i] = s1;
731   }
732 
733   /* backward solve the L^T */
734   for (i=n-1; i>=0; i--){
735     v   = aa + diag[i] - 1 + shift;
736     vi  = aj + diag[i] - 1 + shift;
737     nz  = diag[i] - ai[i];
738     s1  = tmp[i];
739     while (nz--) {
740       tmp[*vi-- + shift] -= (*v--)*s1;
741     }
742   }
743 
744   /* copy tmp into x according to permutation */
745   for (i=0; i<n; i++) x[r[i]] = tmp[i];
746 
747   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
748   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
749   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
750   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
751 
752   PetscLogFlops(2*a->nz-A->n);
753   PetscFunctionReturn(0);
754 }
755 
756 #undef __FUNCT__
757 #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
758 int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
759 {
760   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
761   IS           iscol = a->col,isrow = a->row;
762   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
763   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
764   PetscScalar  *x,*b,*tmp,*aa = a->a,*v;
765 
766   PetscFunctionBegin;
767   if (zz != xx) VecCopy(zz,xx);
768 
769   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
770   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
771   tmp = a->solve_work;
772 
773   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
774   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
775 
776   /* copy the b into temp work space according to permutation */
777   for (i=0; i<n; i++) tmp[i] = b[c[i]];
778 
779   /* forward solve the U^T */
780   for (i=0; i<n; i++) {
781     v   = aa + diag[i] + shift;
782     vi  = aj + diag[i] + (!shift);
783     nz  = ai[i+1] - diag[i] - 1;
784     tmp[i] *= *v++;
785     while (nz--) {
786       tmp[*vi++ + shift] -= (*v++)*tmp[i];
787     }
788   }
789 
790   /* backward solve the L^T */
791   for (i=n-1; i>=0; i--){
792     v   = aa + diag[i] - 1 + shift;
793     vi  = aj + diag[i] - 1 + shift;
794     nz  = diag[i] - ai[i];
795     while (nz--) {
796       tmp[*vi-- + shift] -= (*v--)*tmp[i];
797     }
798   }
799 
800   /* copy tmp into x according to permutation */
801   for (i=0; i<n; i++) x[r[i]] += tmp[i];
802 
803   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
804   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
805   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
806   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
807 
808   PetscLogFlops(2*a->nz);
809   PetscFunctionReturn(0);
810 }
811 /* ----------------------------------------------------------------*/
812 EXTERN int MatMissingDiagonal_SeqAIJ(Mat);
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
816 int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
817 {
818   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
819   IS         isicol;
820   int        *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j;
821   int        *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
822   int        *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0;
823   int        incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill;
824   PetscTruth col_identity,row_identity;
825   PetscReal  f;
826 
827   PetscFunctionBegin;
828   f             = info->fill;
829   levels        = (int)info->levels;
830   diagonal_fill = (int)info->diagonal_fill;
831   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
832 
833   /* special case that simply copies fill pattern */
834   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
835   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
836   if (!levels && row_identity && col_identity) {
837     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
838     (*fact)->factor = FACTOR_LU;
839     b               = (Mat_SeqAIJ*)(*fact)->data;
840     if (!b->diag) {
841       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
842     }
843     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
844     b->row              = isrow;
845     b->col              = iscol;
846     b->icol             = isicol;
847     b->lu_damping       = info->damping;
848     b->lu_zeropivot     = info->zeropivot;
849     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
850     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
851     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
852     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
853     PetscFunctionReturn(0);
854   }
855 
856   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
857   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
858 
859   /* get new row pointers */
860   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
861   ainew[0] = -shift;
862   /* don't know how many column pointers are needed so estimate */
863   jmax = (int)(f*(ai[n]+!shift));
864   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
865   /* ajfill is level of fill for each fill entry */
866   ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr);
867   /* fill is a linked list of nonzeros in active row */
868   ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr);
869   /* im is level for each filled value */
870   ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr);
871   /* dloc is location of diagonal in factor */
872   ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr);
873   dloc[0]  = 0;
874   for (prow=0; prow<n; prow++) {
875 
876     /* copy current row into linked list */
877     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
878     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
879     xi      = aj + ai[r[prow]] + shift;
880     fill[n]    = n;
881     fill[prow] = -1; /* marker to indicate if diagonal exists */
882     while (nz--) {
883       fm  = n;
884       idx = ic[*xi++ + shift];
885       do {
886         m  = fm;
887         fm = fill[m];
888       } while (fm < idx);
889       fill[m]   = idx;
890       fill[idx] = fm;
891       im[idx]   = 0;
892     }
893 
894     /* make sure diagonal entry is included */
895     if (diagonal_fill && fill[prow] == -1) {
896       fm = n;
897       while (fill[fm] < prow) fm = fill[fm];
898       fill[prow] = fill[fm]; /* insert diagonal into linked list */
899       fill[fm]   = prow;
900       im[prow]   = 0;
901       nzf++;
902       dcount++;
903     }
904 
905     nzi = 0;
906     row = fill[n];
907     while (row < prow) {
908       incrlev = im[row] + 1;
909       nz      = dloc[row];
910       xi      = ajnew  + ainew[row] + shift + nz + 1;
911       flev    = ajfill + ainew[row] + shift + nz + 1;
912       nnz     = ainew[row+1] - ainew[row] - nz - 1;
913       fm      = row;
914       while (nnz-- > 0) {
915         idx = *xi++ + shift;
916         if (*flev + incrlev > levels) {
917           flev++;
918           continue;
919         }
920         do {
921           m  = fm;
922           fm = fill[m];
923         } while (fm < idx);
924         if (fm != idx) {
925           im[idx]   = *flev + incrlev;
926           fill[m]   = idx;
927           fill[idx] = fm;
928           fm        = idx;
929           nzf++;
930         } else {
931           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
932         }
933         flev++;
934       }
935       row = fill[row];
936       nzi++;
937     }
938     /* copy new filled row into permanent storage */
939     ainew[prow+1] = ainew[prow] + nzf;
940     if (ainew[prow+1] > jmax-shift) {
941 
942       /* estimate how much additional space we will need */
943       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
944       /* just double the memory each time */
945       /*  maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */
946       int maxadd = jmax;
947       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
948       jmax += maxadd;
949 
950       /* allocate a longer ajnew and ajfill */
951       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
952       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
953       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
954       ajnew  = xi;
955       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
956       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
957       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
958       ajfill = xi;
959       realloc++; /* count how many times we realloc */
960     }
961     xi          = ajnew + ainew[prow] + shift;
962     flev        = ajfill + ainew[prow] + shift;
963     dloc[prow]  = nzi;
964     fm          = fill[n];
965     while (nzf--) {
966       *xi++   = fm - shift;
967       *flev++ = im[fm];
968       fm      = fill[fm];
969     }
970     /* make sure row has diagonal entry */
971     if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) {
972       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\
973     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
974     }
975   }
976   ierr = PetscFree(ajfill);CHKERRQ(ierr);
977   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
978   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
979   ierr = PetscFree(fill);CHKERRQ(ierr);
980   ierr = PetscFree(im);CHKERRQ(ierr);
981 
982   {
983     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
984     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
985     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
986     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
987     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
988     if (diagonal_fill) {
989       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount);
990     }
991   }
992 
993   /* put together the new matrix */
994   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
995   PetscLogObjectParent(*fact,isicol);
996   b = (Mat_SeqAIJ*)(*fact)->data;
997   ierr = PetscFree(b->imax);CHKERRQ(ierr);
998   b->singlemalloc = PETSC_FALSE;
999   len = (ainew[n] + shift)*sizeof(PetscScalar);
1000   /* the next line frees the default space generated by the Create() */
1001   ierr = PetscFree(b->a);CHKERRQ(ierr);
1002   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1003   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1004   b->j          = ajnew;
1005   b->i          = ainew;
1006   for (i=0; i<n; i++) dloc[i] += ainew[i];
1007   b->diag       = dloc;
1008   b->ilen       = 0;
1009   b->imax       = 0;
1010   b->row        = isrow;
1011   b->col        = iscol;
1012   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1013   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1014   b->icol       = isicol;
1015   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1016   /* In b structure:  Free imax, ilen, old a, old j.
1017      Allocate dloc, solve_work, new a, new j */
1018   PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(PetscScalar)));
1019   b->maxnz        = b->nz = ainew[n] + shift;
1020   b->lu_damping   = info->damping;
1021   b->lu_zeropivot = info->zeropivot;
1022   (*fact)->factor = FACTOR_LU;
1023   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
1024   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1025 
1026   (*fact)->info.factor_mallocs    = realloc;
1027   (*fact)->info.fill_ratio_given  = f;
1028   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 #include "src/mat/impls/sbaij/seq/sbaij.h"
1033 #undef __FUNCT__
1034 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1035 int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact)
1036 {
1037   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1038   Mat_SeqSBAIJ        *b = (Mat_SeqSBAIJ*)(*fact)->data;
1039   int                 ierr;
1040 
1041   PetscFunctionBegin;
1042   if (!a->sbaijMat){
1043     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1044   }
1045   ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr);
1046   if (b->factor_levels > 0){
1047     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
1048   }
1049   a->sbaijMat = PETSC_NULL;
1050 
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 #undef __FUNCT__
1055 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1056 int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1057 {
1058   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1059   Mat_SeqSBAIJ        *b;
1060   int                 ierr,levels = info->levels;
1061   PetscTruth          perm_identity;
1062 
1063   PetscFunctionBegin;
1064   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1065   if (!perm_identity){
1066     SETERRQ(1,"Non-identity permutation is not supported yet");
1067   }
1068   if (!a->sbaijMat){
1069     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1070   }
1071 
1072   if (levels > 0){
1073     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1074     b    = (Mat_SeqSBAIJ*)(*fact)->data;
1075   } else { /* in-place icc(0) */
1076     (*fact)         = a->sbaijMat;
1077     (*fact)->factor = FACTOR_CHOLESKY;
1078     b               = (Mat_SeqSBAIJ*)(*fact)->data;
1079     b->row          = perm;
1080     b->icol         = perm;
1081     ierr            = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1082     ierr            = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1083     ierr            = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1084   }
1085   b->factor_levels = levels;
1086 
1087   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1088   (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1089 
1090   PetscFunctionReturn(0);
1091 }
1092 
1093