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