xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 63bb2aa60f4717cbfb86ebed81540aa885ea575c)
173f4d377SMatthew Knepley /*$Id: aijfact.c,v 1.167 2001/09/11 16:32:26 bsmith Exp $*/
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
490f02eecSBarry Smith #include "src/vec/vecimpl.h"
51c4feecaSSatish Balay #include "src/inline/dot.h"
63b2fbd54SBarry Smith 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
991e9ee9fSBarry Smith int MatOrdering_Flow_SeqAIJ(Mat mat,MatOrderingType type,IS *irow,IS *icol)
103b2fbd54SBarry Smith {
113a40ed3dSBarry Smith   PetscFunctionBegin;
123a40ed3dSBarry Smith 
1329bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
14aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
15d64ed03dSBarry Smith   PetscFunctionReturn(0);
16d64ed03dSBarry Smith #endif
173b2fbd54SBarry Smith }
183b2fbd54SBarry Smith 
1986bacbd2SBarry Smith 
20ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqAIJ(Mat);
213a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
2286bacbd2SBarry Smith 
2387828ca2SBarry Smith EXTERN int SPARSEKIT2dperm(int*,PetscScalar*,int*,int*,PetscScalar*,int*,int*,int*,int*,int*);
2487828ca2SBarry Smith EXTERN int SPARSEKIT2ilutp(int*,PetscScalar*,int*,int*,int*,PetscReal,PetscReal*,int*,PetscScalar*,int*,int*,int*,PetscScalar*,int*,int*,int*);
2587828ca2SBarry Smith EXTERN int SPARSEKIT2msrcsr(int*,PetscScalar*,int*,PetscScalar*,int*,int*,PetscScalar*,int*);
2686bacbd2SBarry Smith 
274a2ae208SSatish Balay #undef __FUNCT__
284a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2986bacbd2SBarry Smith   /* ------------------------------------------------------------
3086bacbd2SBarry Smith 
3186bacbd2SBarry Smith           This interface was contribed by Tony Caola
3286bacbd2SBarry Smith 
3386bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
345bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
355bd2ddc7SBarry Smith      SPARSEKIT2.
365bd2ddc7SBarry Smith 
375bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
385bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
3986bacbd2SBarry Smith 
4086bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4186bacbd2SBarry Smith      help in getting this routine ironed out.
4286bacbd2SBarry Smith 
435bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
445bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
455bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
465bd2ddc7SBarry Smith      Yousef Saad's code.
4786bacbd2SBarry Smith 
485bd2ddc7SBarry Smith      ishift = 0, for indices start at 1
495bd2ddc7SBarry Smith      ishift = 1, for indices starting at 0
5086bacbd2SBarry Smith      ------------------------------------------------------------
5186bacbd2SBarry Smith */
5286bacbd2SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,MatILUInfo *info,IS isrow,IS iscol,Mat *fact)
5386bacbd2SBarry Smith {
5460ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5560ec86bdSBarry Smith   PetscFunctionBegin;
5660ec86bdSBarry Smith   SETERRQ(1,"This distribution does not include GNU Copyright code\n\
5760ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5860ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5960ec86bdSBarry Smith #else
6086bacbd2SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b;
6107d2056aSBarry Smith   IS           iscolf,isicol,isirow;
6286bacbd2SBarry Smith   PetscTruth   reorder;
63273d9f13SBarry Smith   int          *c,*r,*ic,ierr,i,n = A->m;
64329f5518SBarry Smith   int          *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
65313ae024SBarry Smith   int          *ordcol,*iwk,*iperm,*jw;
665bd2ddc7SBarry Smith   int          ishift = !a->indexshift;
67b19713cbSBarry Smith   int          jmax,lfill,job,*o_i,*o_j;
6887828ca2SBarry Smith   PetscScalar  *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
69329f5518SBarry Smith   PetscReal    permtol,af;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   PetscFunctionBegin;
7286bacbd2SBarry Smith 
7386bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7486bacbd2SBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax);
7586bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7633259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
776faa4630SBarry Smith   lfill   = (int)(info->dtcount/2.0);
786faa4630SBarry Smith   jmax    = (int)(info->fill*a->nz);
7986bacbd2SBarry Smith   permtol = info->dtcol;
8086bacbd2SBarry Smith 
8186bacbd2SBarry Smith 
8286bacbd2SBarry Smith   /* ------------------------------------------------------------
8386bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8486bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8586bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8686bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8786bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8886bacbd2SBarry Smith      ------------------------------------------------------------
8986bacbd2SBarry Smith   */
9086bacbd2SBarry Smith 
9186bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9286bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9386bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9486bacbd2SBarry Smith   reorder = PetscNot(reorder);
9586bacbd2SBarry Smith 
9686bacbd2SBarry Smith 
9786bacbd2SBarry Smith   /* storage for ilu factor */
98b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr);
99b0a32e0cSBarry Smith   ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr);
10087828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
101b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr);
10286bacbd2SBarry Smith 
10386bacbd2SBarry Smith   /* ------------------------------------------------------------
10486bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10586bacbd2SBarry Smith      ------------------------------------------------------------
10686bacbd2SBarry Smith   */
107b19713cbSBarry Smith   if (ishift) {
108b19713cbSBarry Smith     for (i=old_i[0];i<old_i[n];i++) {
109b19713cbSBarry Smith       old_j[i]++;
11086bacbd2SBarry Smith     }
111b19713cbSBarry Smith     for(i=0;i<n+1;i++) {
112b19713cbSBarry Smith       old_i[i]++;
113b19713cbSBarry Smith     };
11486bacbd2SBarry Smith   };
11586bacbd2SBarry Smith 
11686bacbd2SBarry Smith   if (reorder) {
117c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);           CHKERRQ(ierr);
118c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);           CHKERRQ(ierr);
119c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
120c0c2c81eSBarry Smith       r[i]  = r[i]+1;
121c0c2c81eSBarry Smith       c[i]  = c[i]+1;
122c0c2c81eSBarry Smith     }
123b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr);
124b0a32e0cSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr);
12587828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1265bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
127c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
128c0c2c81eSBarry Smith       r[i]  = r[i]-1;
129c0c2c81eSBarry Smith       c[i]  = c[i]-1;
130c0c2c81eSBarry Smith     }
131c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
132c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
133b19713cbSBarry Smith     o_a = old_a2;
134b19713cbSBarry Smith     o_j = old_j2;
135b19713cbSBarry Smith     o_i = old_i2;
136b19713cbSBarry Smith   } else {
137b19713cbSBarry Smith     o_a = old_a;
138b19713cbSBarry Smith     o_j = old_j;
139b19713cbSBarry Smith     o_i = old_i;
14086bacbd2SBarry Smith   }
14186bacbd2SBarry Smith 
14286bacbd2SBarry Smith   /* ------------------------------------------------------------
14386bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14486bacbd2SBarry Smith      ------------------------------------------------------------
14586bacbd2SBarry Smith   */
14686bacbd2SBarry Smith 
147b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr);
148b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr);
14987828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
15086bacbd2SBarry Smith 
15187828ca2SBarry Smith   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);
15286bacbd2SBarry Smith   if (ierr) {
15386bacbd2SBarry Smith     switch (ierr) {
15429bbc08cSBarry Smith       case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15529bbc08cSBarry Smith       case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15629bbc08cSBarry Smith       case -5: SETERRQ(1,"ilutp(), zero row encountered");
15729bbc08cSBarry Smith       case -1: SETERRQ(1,"ilutp(), input matrix may be wrong");
15829bbc08cSBarry Smith       case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax);
15929bbc08cSBarry Smith       default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr);
16086bacbd2SBarry Smith     }
16186bacbd2SBarry Smith   }
16286bacbd2SBarry Smith 
16386bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16486bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16586bacbd2SBarry Smith 
16686bacbd2SBarry Smith   /* ------------------------------------------------------------
16786bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16886bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16986bacbd2SBarry Smith      ------------------------------------------------------------
17086bacbd2SBarry Smith   */
17186bacbd2SBarry Smith 
17287828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
173b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr);
17486bacbd2SBarry Smith 
1755bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17686bacbd2SBarry Smith 
17786bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17886bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17986bacbd2SBarry Smith 
18086bacbd2SBarry Smith   if (reorder) {
18186bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
18286bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18386bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
184313ae024SBarry Smith   } else {
185b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
186f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
187b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
188b19713cbSBarry Smith     }
189b19713cbSBarry Smith   }
190b19713cbSBarry Smith 
191b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
192b19713cbSBarry Smith   if (ishift) {
19386bacbd2SBarry Smith     for (i=0; i<n+1; i++) {
194b19713cbSBarry Smith       old_i[i]--;
195b19713cbSBarry Smith     }
196b19713cbSBarry Smith     for (i=old_i[0];i<old_i[n];i++) {
197b19713cbSBarry Smith       old_j[i]--;
198b19713cbSBarry Smith     }
19986bacbd2SBarry Smith   }
20086bacbd2SBarry Smith 
201b19713cbSBarry Smith   /* Make the factored matrix 0-based */
202b19713cbSBarry Smith   if (ishift) {
20386bacbd2SBarry Smith     for (i=0; i<n+1; i++) {
204b19713cbSBarry Smith       new_i[i]--;
205b19713cbSBarry Smith     }
206b19713cbSBarry Smith     for (i=new_i[0];i<new_i[n];i++) {
207b19713cbSBarry Smith       new_j[i]--;
208b19713cbSBarry Smith     }
20986bacbd2SBarry Smith   }
21086bacbd2SBarry Smith 
21186bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
21286bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2134c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2144c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
215c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21686bacbd2SBarry Smith   for(i=0; i<n; i++) {
21786bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21886bacbd2SBarry Smith   };
21986bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
22007d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
22186bacbd2SBarry Smith 
222c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
223c0c2c81eSBarry Smith 
224329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
22507d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
22686bacbd2SBarry Smith 
22786bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22886bacbd2SBarry Smith 
22986bacbd2SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
23086bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
23186bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
23286bacbd2SBarry Smith 
23386bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
23486bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
23586bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23607d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
237b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23886bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
23986bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
24086bacbd2SBarry Smith   b->a             = new_a;
24186bacbd2SBarry Smith   b->j             = new_j;
24286bacbd2SBarry Smith   b->i             = new_i;
24386bacbd2SBarry Smith   b->ilen          = 0;
24486bacbd2SBarry Smith   b->imax          = 0;
2451f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
246313ae024SBarry Smith   b->row           = isirow;
24786bacbd2SBarry Smith   b->col           = iscolf;
24887828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24986bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
25086bacbd2SBarry Smith   b->indexshift    = a->indexshift;
25186bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
25286bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
25386bacbd2SBarry Smith 
25486bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
25586bacbd2SBarry Smith 
25686bacbd2SBarry Smith   /* check out for identical nodes. If found, use inode functions */
2573a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
25886bacbd2SBarry Smith 
259431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
260b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
261b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
262b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
263b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
264431e710aSBarry Smith 
26586bacbd2SBarry Smith   PetscFunctionReturn(0);
26660ec86bdSBarry Smith #endif
26786bacbd2SBarry Smith }
26886bacbd2SBarry Smith 
269289bc588SBarry Smith /*
270289bc588SBarry Smith     Factorization code for AIJ format.
271289bc588SBarry Smith */
2724a2ae208SSatish Balay #undef __FUNCT__
273b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
2749e046878SBarry Smith int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatLUInfo *info,Mat *B)
275289bc588SBarry Smith {
276416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
277289bc588SBarry Smith   IS         isicol;
278273d9f13SBarry Smith   int        *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j;
279416022c9SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*ajtmp,nz,shift = a->indexshift;
28091bf9adeSBarry Smith   int        *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im;
2819e046878SBarry Smith   PetscReal  f;
282289bc588SBarry Smith 
2833a40ed3dSBarry Smith   PetscFunctionBegin;
28429bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2854c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
286ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
287ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
288289bc588SBarry Smith 
289289bc588SBarry Smith   /* get new row pointers */
290b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
291dbb450caSBarry Smith   ainew[0] = -shift;
292289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
293d3d32019SBarry Smith   f = info->fill;
294dbb450caSBarry Smith   jmax  = (int)(f*ai[n]+(!shift));
295b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
296289bc588SBarry Smith   /* fill is a linked list of nonzeros in active row */
297b0a32e0cSBarry Smith   ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr);
2982fb3ed76SBarry Smith   im   = fill + n + 1;
299289bc588SBarry Smith   /* idnew is location of diagonal in factor */
300b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr);
301dbb450caSBarry Smith   idnew[0] = -shift;
302289bc588SBarry Smith 
303289bc588SBarry Smith   for (i=0; i<n; i++) {
304289bc588SBarry Smith     /* first copy previous fill into linked list */
305289bc588SBarry Smith     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
30629bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
307dbb450caSBarry Smith     ajtmp   = aj + ai[r[i]] + shift;
308289bc588SBarry Smith     fill[n] = n;
309289bc588SBarry Smith     while (nz--) {
310289bc588SBarry Smith       fm  = n;
311dbb450caSBarry Smith       idx = ic[*ajtmp++ + shift];
312289bc588SBarry Smith       do {
313289bc588SBarry Smith         m  = fm;
314289bc588SBarry Smith         fm = fill[m];
315289bc588SBarry Smith       } while (fm < idx);
316289bc588SBarry Smith       fill[m]   = idx;
317289bc588SBarry Smith       fill[idx] = fm;
318289bc588SBarry Smith     }
319289bc588SBarry Smith     row = fill[n];
320289bc588SBarry Smith     while (row < i) {
321dbb450caSBarry Smith       ajtmp = ajnew + idnew[row] + (!shift);
3222fb3ed76SBarry Smith       nzbd  = 1 + idnew[row] - ainew[row];
3232fb3ed76SBarry Smith       nz    = im[row] - nzbd;
324289bc588SBarry Smith       fm    = row;
3252fb3ed76SBarry Smith       while (nz-- > 0) {
326dbb450caSBarry Smith         idx = *ajtmp++ + shift;
3272fb3ed76SBarry Smith         nzbd++;
3282fb3ed76SBarry Smith         if (idx == i) im[row] = nzbd;
329289bc588SBarry Smith         do {
330289bc588SBarry Smith           m  = fm;
331289bc588SBarry Smith           fm = fill[m];
332289bc588SBarry Smith         } while (fm < idx);
333289bc588SBarry Smith         if (fm != idx) {
334289bc588SBarry Smith           fill[m]   = idx;
335289bc588SBarry Smith           fill[idx] = fm;
336289bc588SBarry Smith           fm        = idx;
337289bc588SBarry Smith           nnz++;
338289bc588SBarry Smith         }
339289bc588SBarry Smith       }
340289bc588SBarry Smith       row = fill[row];
341289bc588SBarry Smith     }
342289bc588SBarry Smith     /* copy new filled row into permanent storage */
343289bc588SBarry Smith     ainew[i+1] = ainew[i] + nnz;
344496e697eSBarry Smith     if (ainew[i+1] > jmax) {
345ecf371e4SBarry Smith 
346ecf371e4SBarry Smith       /* estimate how much additional space we will need */
347ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
348ecf371e4SBarry Smith       /* just double the memory each time */
349ecf371e4SBarry Smith       int maxadd = jmax;
350ecf371e4SBarry Smith       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
351bbb6d6a8SBarry Smith       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
3522fb3ed76SBarry Smith       jmax += maxadd;
353ecf371e4SBarry Smith 
354ecf371e4SBarry Smith       /* allocate a longer ajnew */
355b0a32e0cSBarry Smith       ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr);
356549d3d68SSatish Balay       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i]+shift)*sizeof(int));CHKERRQ(ierr);
357606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
358289bc588SBarry Smith       ajnew = ajtmp;
3592fb3ed76SBarry Smith       realloc++; /* count how many times we realloc */
360289bc588SBarry Smith     }
361dbb450caSBarry Smith     ajtmp = ajnew + ainew[i] + shift;
362289bc588SBarry Smith     fm    = fill[n];
363289bc588SBarry Smith     nzi   = 0;
3642fb3ed76SBarry Smith     im[i] = nnz;
365289bc588SBarry Smith     while (nnz--) {
366289bc588SBarry Smith       if (fm < i) nzi++;
367dbb450caSBarry Smith       *ajtmp++ = fm - shift;
368289bc588SBarry Smith       fm       = fill[fm];
369289bc588SBarry Smith     }
370289bc588SBarry Smith     idnew[i] = ainew[i] + nzi;
371289bc588SBarry Smith   }
372464e8d44SSatish Balay   if (ai[n] != 0) {
373329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
374b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
375b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
376b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
377b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
37805bf559cSSatish Balay   } else {
379b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
38005bf559cSSatish Balay   }
3812fb3ed76SBarry Smith 
382898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
383898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3841987afe7SBarry Smith 
385606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
386289bc588SBarry Smith 
387289bc588SBarry Smith   /* put together the new matrix */
388b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr);
389b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
390416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*B)->data;
391606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3927c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
393e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
394606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
395606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
39687828ca2SBarry Smith   ierr          = PetscMalloc((ainew[n]+shift+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
397416022c9SBarry Smith   b->j          = ajnew;
398416022c9SBarry Smith   b->i          = ainew;
399416022c9SBarry Smith   b->diag       = idnew;
400416022c9SBarry Smith   b->ilen       = 0;
401416022c9SBarry Smith   b->imax       = 0;
402416022c9SBarry Smith   b->row        = isrow;
403416022c9SBarry Smith   b->col        = iscol;
404085256b3SBarry Smith   b->lu_damping   = info->damping;
40587828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
406c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
407c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40882bf6240SBarry Smith   b->icol       = isicol;
40987828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
410416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4117fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
41287828ca2SBarry Smith   PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(PetscScalar)));
413416022c9SBarry Smith   b->maxnz = b->nz = ainew[n] + shift;
4147fd98bd6SLois Curfman McInnes 
415329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
416ae068f58SLois Curfman McInnes   (*B)->info.factor_mallocs    = realloc;
417ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4183a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4197dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
420703deb49SSatish Balay 
421e93ccc4dSBarry Smith   if (ai[n] != 0) {
422329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
423eea2913aSSatish Balay   } else {
424eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
425eea2913aSSatish Balay   }
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
427289bc588SBarry Smith }
4280a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4293a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
43041c01911SSatish Balay 
4314a2ae208SSatish Balay #undef __FUNCT__
4324a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
433416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B)
434289bc588SBarry Smith {
435416022c9SBarry Smith   Mat          C = *B;
436416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data;
43782bf6240SBarry Smith   IS           isrow = b->row,isicol = b->icol;
438273d9f13SBarry Smith   int          *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j;
4391987afe7SBarry Smith   int          *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift;
440085256b3SBarry Smith   int          *diag_offset = b->diag,diag;
441085256b3SBarry Smith   int          *pj,ndamp = 0;
44287828ca2SBarry Smith   PetscScalar  *rtmp,*v,*pc,multiplier,*pv,*rtmps;
443d3d32019SBarry Smith   PetscReal    damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs;
4448a5eb818SBarry Smith   PetscTruth   damp;
445289bc588SBarry Smith 
4463a40ed3dSBarry Smith   PetscFunctionBegin;
44778b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44878b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44987828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45087828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
4510cf60383SBarry Smith   rtmps = rtmp + shift; ics = ic + shift;
452289bc588SBarry Smith 
453085256b3SBarry Smith   do {
4548a5eb818SBarry Smith     damp = PETSC_FALSE;
455289bc588SBarry Smith     for (i=0; i<n; i++) {
456289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
457dbb450caSBarry Smith       ajtmp = aj + ai[i] + shift;
45848e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
459289bc588SBarry Smith 
460289bc588SBarry Smith       /* load in initial (unfactored row) */
461416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
462416022c9SBarry Smith       ajtmpold = a->j + a->i[r[i]] + shift;
463416022c9SBarry Smith       v        = a->a + a->i[r[i]] + shift;
464085256b3SBarry Smith       for (j=0; j<nz; j++) {
465085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
4668a5eb818SBarry Smith         if (ndamp && ajtmpold[j] == r[i]) { /* damp the diagonal of the matrix */
467085256b3SBarry Smith           rtmp[ics[ajtmpold[j]]] += damping;
468085256b3SBarry Smith         }
469085256b3SBarry Smith       }
470289bc588SBarry Smith 
471dbb450caSBarry Smith       row = *ajtmp++ + shift;
472f2582111SSatish Balay       while  (row < i) {
4738c37ef55SBarry Smith         pc = rtmp + row;
4748c37ef55SBarry Smith         if (*pc != 0.0) {
4751987afe7SBarry Smith           pv         = b->a + diag_offset[row] + shift;
4761987afe7SBarry Smith           pj         = b->j + diag_offset[row] + (!shift);
47735aab85fSBarry Smith           multiplier = *pc / *pv++;
4788c37ef55SBarry Smith           *pc        = multiplier;
479f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
480f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
481b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
482289bc588SBarry Smith         }
483dbb450caSBarry Smith         row = *ajtmp++ + shift;
484289bc588SBarry Smith       }
485416022c9SBarry Smith       /* finished row so stick it into b->a */
486416022c9SBarry Smith       pv   = b->a + ai[i] + shift;
487416022c9SBarry Smith       pj   = b->j + ai[i] + shift;
4888c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
48935aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
490d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
491d3d32019SBarry Smith       rs   = 0.0;
492d3d32019SBarry Smith       for (j=0; j<nz; j++) {
493d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
494d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
495d3d32019SBarry Smith       }
496d3d32019SBarry Smith       if (PetscAbsScalar(pv[diag]) < zeropivot*rs) {
497d3d32019SBarry Smith         if (damping) {
4988a5eb818SBarry Smith           if (ndamp) damping *= 2.0;
499085256b3SBarry Smith           damp = PETSC_TRUE;
500085256b3SBarry Smith           ndamp++;
501085256b3SBarry Smith           break;
502085256b3SBarry Smith         } else {
503*63bb2aa6SBarry Smith           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
504d708749eSBarry Smith         }
50535aab85fSBarry Smith       }
50635aab85fSBarry Smith     }
507085256b3SBarry Smith   } while (damp);
5080f11f9aeSSatish Balay 
50935aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
51035aab85fSBarry Smith   for (i=0; i<n; i++) {
51135aab85fSBarry Smith     b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift];
51235aab85fSBarry Smith   }
51335aab85fSBarry Smith 
514606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
51578b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
51678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
517416022c9SBarry Smith   C->factor = FACTOR_LU;
5187dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
519c456f294SBarry Smith   C->assembled = PETSC_TRUE;
520b0a32e0cSBarry Smith   PetscLogFlops(C->n);
521d3d32019SBarry Smith   if (ndamp) {
522b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping);
523085256b3SBarry Smith   }
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
525289bc588SBarry Smith }
5260a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5274a2ae208SSatish Balay #undef __FUNCT__
5284a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
5299e046878SBarry Smith int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatLUInfo *info)
530da3a660dSBarry Smith {
531273d9f13SBarry Smith   int ierr;
532416022c9SBarry Smith   Mat C;
533416022c9SBarry Smith 
5343a40ed3dSBarry Smith   PetscFunctionBegin;
5359e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
536f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
537273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
538440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5393a40ed3dSBarry Smith   PetscFunctionReturn(0);
540da3a660dSBarry Smith }
5410a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
544416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5458c37ef55SBarry Smith {
546416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
547416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
548273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
5493b2fbd54SBarry Smith   int          nz,shift = a->indexshift,*rout,*cout;
55087828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
5518c37ef55SBarry Smith 
5523a40ed3dSBarry Smith   PetscFunctionBegin;
5533a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
55491bf9adeSBarry Smith 
555e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
556e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
557416022c9SBarry Smith   tmp  = a->solve_work;
5588c37ef55SBarry Smith 
5593b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
5603b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
5618c37ef55SBarry Smith 
5628c37ef55SBarry Smith   /* forward solve the lower triangular */
5638c37ef55SBarry Smith   tmp[0] = b[*r++];
5640cf60383SBarry Smith   tmps   = tmp + shift;
5658c37ef55SBarry Smith   for (i=1; i<n; i++) {
566dbb450caSBarry Smith     v   = aa + ai[i] + shift;
567dbb450caSBarry Smith     vi  = aj + ai[i] + shift;
568416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
5698c37ef55SBarry Smith     sum = b[*r++];
5700cf60383SBarry Smith     while (nz--) sum -= *v++ * tmps[*vi++];
5718c37ef55SBarry Smith     tmp[i] = sum;
5728c37ef55SBarry Smith   }
5738c37ef55SBarry Smith 
5748c37ef55SBarry Smith   /* backward solve the upper triangular */
5758c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
576416022c9SBarry Smith     v   = aa + a->diag[i] + (!shift);
577416022c9SBarry Smith     vi  = aj + a->diag[i] + (!shift);
578416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
5798c37ef55SBarry Smith     sum = tmp[i];
5800cf60383SBarry Smith     while (nz--) sum -= *v++ * tmps[*vi++];
581416022c9SBarry Smith     x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift];
5828c37ef55SBarry Smith   }
5838c37ef55SBarry Smith 
5843b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
5853b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
586cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
587cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
588b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
5908c37ef55SBarry Smith }
591026e39d0SSatish Balay 
592930ae53cSBarry Smith /* ----------------------------------------------------------- */
5934a2ae208SSatish Balay #undef __FUNCT__
5944a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
595930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
596930ae53cSBarry Smith {
597930ae53cSBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
598273d9f13SBarry Smith   int          n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr;
599362ced78SSatish Balay   PetscScalar  *x,*b,*aa = a->a;
600aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
601d85afc42SSatish Balay   int          adiag_i,i,*vi,nz,ai_i;
602362ced78SSatish Balay   PetscScalar  *v,sum;
603d85afc42SSatish Balay #endif
604930ae53cSBarry Smith 
6053a40ed3dSBarry Smith   PetscFunctionBegin;
6063a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
607930ae53cSBarry Smith   if (a->indexshift) {
6083a40ed3dSBarry Smith      ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr);
6093a40ed3dSBarry Smith      PetscFunctionReturn(0);
610930ae53cSBarry Smith   }
611930ae53cSBarry Smith 
612e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
613e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
614930ae53cSBarry Smith 
615aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6161c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6171c4feecaSSatish Balay #else
618930ae53cSBarry Smith   /* forward solve the lower triangular */
619930ae53cSBarry Smith   x[0] = b[0];
620930ae53cSBarry Smith   for (i=1; i<n; i++) {
621930ae53cSBarry Smith     ai_i = ai[i];
622930ae53cSBarry Smith     v    = aa + ai_i;
623930ae53cSBarry Smith     vi   = aj + ai_i;
624930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
625930ae53cSBarry Smith     sum  = b[i];
626930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
627930ae53cSBarry Smith     x[i] = sum;
628930ae53cSBarry Smith   }
629930ae53cSBarry Smith 
630930ae53cSBarry Smith   /* backward solve the upper triangular */
631930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
632930ae53cSBarry Smith     adiag_i = adiag[i];
633d708749eSBarry Smith     v       = aa + adiag_i + 1;
634d708749eSBarry Smith     vi      = aj + adiag_i + 1;
635930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
636930ae53cSBarry Smith     sum     = x[i];
637930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
638930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
639930ae53cSBarry Smith   }
6401c4feecaSSatish Balay #endif
641b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
642cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
643cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6443a40ed3dSBarry Smith   PetscFunctionReturn(0);
645930ae53cSBarry Smith }
646930ae53cSBarry Smith 
6474a2ae208SSatish Balay #undef __FUNCT__
6484a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
649416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
650da3a660dSBarry Smith {
651416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
652416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
653273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
6543b2fbd54SBarry Smith   int          nz,shift = a->indexshift,*rout,*cout;
65587828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,sum,*v;
656da3a660dSBarry Smith 
6573a40ed3dSBarry Smith   PetscFunctionBegin;
65878b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
659da3a660dSBarry Smith 
660e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
661e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
662416022c9SBarry Smith   tmp  = a->solve_work;
663da3a660dSBarry Smith 
6643b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6653b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
666da3a660dSBarry Smith 
667da3a660dSBarry Smith   /* forward solve the lower triangular */
668da3a660dSBarry Smith   tmp[0] = b[*r++];
669da3a660dSBarry Smith   for (i=1; i<n; i++) {
670dbb450caSBarry Smith     v   = aa + ai[i] + shift;
671dbb450caSBarry Smith     vi  = aj + ai[i] + shift;
672416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
673da3a660dSBarry Smith     sum = b[*r++];
674dbb450caSBarry Smith     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
675da3a660dSBarry Smith     tmp[i] = sum;
676da3a660dSBarry Smith   }
677da3a660dSBarry Smith 
678da3a660dSBarry Smith   /* backward solve the upper triangular */
679da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
680416022c9SBarry Smith     v   = aa + a->diag[i] + (!shift);
681416022c9SBarry Smith     vi  = aj + a->diag[i] + (!shift);
682416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
683da3a660dSBarry Smith     sum = tmp[i];
684dbb450caSBarry Smith     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
685416022c9SBarry Smith     tmp[i] = sum*aa[a->diag[i]+shift];
686da3a660dSBarry Smith     x[*c--] += tmp[i];
687da3a660dSBarry Smith   }
688da3a660dSBarry Smith 
6893b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6903b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
691cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
692cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
693b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
694898183ecSLois Curfman McInnes 
6953a40ed3dSBarry Smith   PetscFunctionReturn(0);
696da3a660dSBarry Smith }
697da3a660dSBarry Smith /* -------------------------------------------------------------------*/
6984a2ae208SSatish Balay #undef __FUNCT__
6994a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
7007c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
701da3a660dSBarry Smith {
702416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
7032235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
704273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
705f1af5d2fSBarry Smith   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
70687828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v,s1;
707da3a660dSBarry Smith 
7083a40ed3dSBarry Smith   PetscFunctionBegin;
709e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
710e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
711416022c9SBarry Smith   tmp  = a->solve_work;
712da3a660dSBarry Smith 
7132235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7142235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
715da3a660dSBarry Smith 
716da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7172235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
718da3a660dSBarry Smith 
719da3a660dSBarry Smith   /* forward solve the U^T */
720da3a660dSBarry Smith   for (i=0; i<n; i++) {
721f1af5d2fSBarry Smith     v   = aa + diag[i] + shift;
722f1af5d2fSBarry Smith     vi  = aj + diag[i] + (!shift);
723f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
724c38d4ed2SBarry Smith     s1  = tmp[i];
725ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
726da3a660dSBarry Smith     while (nz--) {
727f1af5d2fSBarry Smith       tmp[*vi++ + shift] -= (*v++)*s1;
728da3a660dSBarry Smith     }
729c38d4ed2SBarry Smith     tmp[i] = s1;
730da3a660dSBarry Smith   }
731da3a660dSBarry Smith 
732da3a660dSBarry Smith   /* backward solve the L^T */
733da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
734f1af5d2fSBarry Smith     v   = aa + diag[i] - 1 + shift;
735f1af5d2fSBarry Smith     vi  = aj + diag[i] - 1 + shift;
736f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
737f1af5d2fSBarry Smith     s1  = tmp[i];
738da3a660dSBarry Smith     while (nz--) {
739f1af5d2fSBarry Smith       tmp[*vi-- + shift] -= (*v--)*s1;
740da3a660dSBarry Smith     }
741da3a660dSBarry Smith   }
742da3a660dSBarry Smith 
743da3a660dSBarry Smith   /* copy tmp into x according to permutation */
744da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
745da3a660dSBarry Smith 
7462235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7472235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
748cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
749cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
750da3a660dSBarry Smith 
751b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
7523a40ed3dSBarry Smith   PetscFunctionReturn(0);
753da3a660dSBarry Smith }
754da3a660dSBarry Smith 
7554a2ae208SSatish Balay #undef __FUNCT__
7564a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
7577c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
758da3a660dSBarry Smith {
759416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
7602235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
761273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
762f1af5d2fSBarry Smith   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
76387828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v;
7646abc6512SBarry Smith 
7653a40ed3dSBarry Smith   PetscFunctionBegin;
7666abc6512SBarry Smith   if (zz != xx) VecCopy(zz,xx);
7676abc6512SBarry Smith 
768e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
769e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
770416022c9SBarry Smith   tmp = a->solve_work;
7716abc6512SBarry Smith 
7722235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7732235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
7746abc6512SBarry Smith 
7756abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
7762235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
7776abc6512SBarry Smith 
7786abc6512SBarry Smith   /* forward solve the U^T */
7796abc6512SBarry Smith   for (i=0; i<n; i++) {
780f1af5d2fSBarry Smith     v   = aa + diag[i] + shift;
781f1af5d2fSBarry Smith     vi  = aj + diag[i] + (!shift);
782f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
7836abc6512SBarry Smith     tmp[i] *= *v++;
7846abc6512SBarry Smith     while (nz--) {
785dbb450caSBarry Smith       tmp[*vi++ + shift] -= (*v++)*tmp[i];
7866abc6512SBarry Smith     }
7876abc6512SBarry Smith   }
7886abc6512SBarry Smith 
7896abc6512SBarry Smith   /* backward solve the L^T */
7906abc6512SBarry Smith   for (i=n-1; i>=0; i--){
791f1af5d2fSBarry Smith     v   = aa + diag[i] - 1 + shift;
792f1af5d2fSBarry Smith     vi  = aj + diag[i] - 1 + shift;
793f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
7946abc6512SBarry Smith     while (nz--) {
795dbb450caSBarry Smith       tmp[*vi-- + shift] -= (*v--)*tmp[i];
7966abc6512SBarry Smith     }
7976abc6512SBarry Smith   }
7986abc6512SBarry Smith 
7996abc6512SBarry Smith   /* copy tmp into x according to permutation */
8006abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8016abc6512SBarry Smith 
8022235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8032235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
804cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
805cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8066abc6512SBarry Smith 
807b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8083a40ed3dSBarry Smith   PetscFunctionReturn(0);
809da3a660dSBarry Smith }
8109e25ed09SBarry Smith /* ----------------------------------------------------------------*/
811ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat);
81208480c60SBarry Smith 
8134a2ae208SSatish Balay #undef __FUNCT__
8144a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
8156cf997b0SBarry Smith int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatILUInfo *info,Mat *fact)
8169e25ed09SBarry Smith {
817416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
8189e25ed09SBarry Smith   IS         isicol;
819273d9f13SBarry Smith   int        *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j;
82056beaf04SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
821335d9088SBarry Smith   int        *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0;
8226cf997b0SBarry Smith   int        incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill;
82377c4ece6SBarry Smith   PetscTruth col_identity,row_identity;
824329f5518SBarry Smith   PetscReal  f;
8259e25ed09SBarry Smith 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
8276cf997b0SBarry Smith   f             = info->fill;
8280cd17293SBarry Smith   levels        = (int)info->levels;
8290cd17293SBarry Smith   diagonal_fill = (int)info->diagonal_fill;
8304c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
83182bf6240SBarry Smith 
83208480c60SBarry Smith   /* special case that simply copies fill pattern */
833be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
834be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
83586bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8362e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
83708480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
83808480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
83908480c60SBarry Smith     if (!b->diag) {
8407c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
84108480c60SBarry Smith     }
8427c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
84308480c60SBarry Smith     b->row              = isrow;
84408480c60SBarry Smith     b->col              = iscol;
84582bf6240SBarry Smith     b->icol             = isicol;
846a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
84787828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
84887828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
849f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
850c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
851c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
8523a40ed3dSBarry Smith     PetscFunctionReturn(0);
85308480c60SBarry Smith   }
85408480c60SBarry Smith 
855898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
856898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
8579e25ed09SBarry Smith 
8589e25ed09SBarry Smith   /* get new row pointers */
859b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
860dbb450caSBarry Smith   ainew[0] = -shift;
8619e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
862dbb450caSBarry Smith   jmax = (int)(f*(ai[n]+!shift));
863b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
8649e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
865b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr);
8669e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
867b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr);
86856beaf04SBarry Smith   /* im is level for each filled value */
869b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr);
87056beaf04SBarry Smith   /* dloc is location of diagonal in factor */
871b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr);
87256beaf04SBarry Smith   dloc[0]  = 0;
87356beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
8746cf997b0SBarry Smith 
8756cf997b0SBarry Smith     /* copy current row into linked list */
87656beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
87729bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
878dbb450caSBarry Smith     xi      = aj + ai[r[prow]] + shift;
8799e25ed09SBarry Smith     fill[n]    = n;
880435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
8819e25ed09SBarry Smith     while (nz--) {
8829e25ed09SBarry Smith       fm  = n;
883dbb450caSBarry Smith       idx = ic[*xi++ + shift];
8849e25ed09SBarry Smith       do {
8859e25ed09SBarry Smith         m  = fm;
8869e25ed09SBarry Smith         fm = fill[m];
8879e25ed09SBarry Smith       } while (fm < idx);
8889e25ed09SBarry Smith       fill[m]   = idx;
8899e25ed09SBarry Smith       fill[idx] = fm;
89056beaf04SBarry Smith       im[idx]   = 0;
8919e25ed09SBarry Smith     }
8926cf997b0SBarry Smith 
8936cf997b0SBarry Smith     /* make sure diagonal entry is included */
894435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
8956cf997b0SBarry Smith       fm = n;
896435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
897435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
898435faa5fSBarry Smith       fill[fm]   = prow;
8996cf997b0SBarry Smith       im[prow]   = 0;
9006cf997b0SBarry Smith       nzf++;
901335d9088SBarry Smith       dcount++;
9026cf997b0SBarry Smith     }
9036cf997b0SBarry Smith 
90456beaf04SBarry Smith     nzi = 0;
9059e25ed09SBarry Smith     row = fill[n];
90656beaf04SBarry Smith     while (row < prow) {
90756beaf04SBarry Smith       incrlev = im[row] + 1;
90856beaf04SBarry Smith       nz      = dloc[row];
9096cf997b0SBarry Smith       xi      = ajnew  + ainew[row] + shift + nz + 1;
910dbb450caSBarry Smith       flev    = ajfill + ainew[row] + shift + nz + 1;
911416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
91256beaf04SBarry Smith       fm      = row;
91356beaf04SBarry Smith       while (nnz-- > 0) {
914dbb450caSBarry Smith         idx = *xi++ + shift;
91556beaf04SBarry Smith         if (*flev + incrlev > levels) {
91656beaf04SBarry Smith           flev++;
91756beaf04SBarry Smith           continue;
91856beaf04SBarry Smith         }
91956beaf04SBarry Smith         do {
9209e25ed09SBarry Smith           m  = fm;
9219e25ed09SBarry Smith           fm = fill[m];
92256beaf04SBarry Smith         } while (fm < idx);
9239e25ed09SBarry Smith         if (fm != idx) {
92456beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9259e25ed09SBarry Smith           fill[m]   = idx;
9269e25ed09SBarry Smith           fill[idx] = fm;
9279e25ed09SBarry Smith           fm        = idx;
92856beaf04SBarry Smith           nzf++;
929ecf371e4SBarry Smith         } else {
93056beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9319e25ed09SBarry Smith         }
93256beaf04SBarry Smith         flev++;
9339e25ed09SBarry Smith       }
9349e25ed09SBarry Smith       row = fill[row];
93556beaf04SBarry Smith       nzi++;
9369e25ed09SBarry Smith     }
9379e25ed09SBarry Smith     /* copy new filled row into permanent storage */
93856beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
939d7e8b826SBarry Smith     if (ainew[prow+1] > jmax-shift) {
940ecf371e4SBarry Smith 
941ecf371e4SBarry Smith       /* estimate how much additional space we will need */
942ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
943ecf371e4SBarry Smith       /* just double the memory each time */
944ecf371e4SBarry Smith       /*  maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */
945ecf371e4SBarry Smith       int maxadd = jmax;
94656beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
947bbb6d6a8SBarry Smith       jmax += maxadd;
948ecf371e4SBarry Smith 
949ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
950b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
951549d3d68SSatish Balay       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
952606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
95356beaf04SBarry Smith       ajnew  = xi;
954b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
955549d3d68SSatish Balay       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
956606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
95756beaf04SBarry Smith       ajfill = xi;
958ecf371e4SBarry Smith       realloc++; /* count how many times we realloc */
9599e25ed09SBarry Smith     }
960dbb450caSBarry Smith     xi          = ajnew + ainew[prow] + shift;
961dbb450caSBarry Smith     flev        = ajfill + ainew[prow] + shift;
96256beaf04SBarry Smith     dloc[prow]  = nzi;
9639e25ed09SBarry Smith     fm          = fill[n];
96456beaf04SBarry Smith     while (nzf--) {
965dbb450caSBarry Smith       *xi++   = fm - shift;
96656beaf04SBarry Smith       *flev++ = im[fm];
9679e25ed09SBarry Smith       fm      = fill[fm];
9689e25ed09SBarry Smith     }
9696cf997b0SBarry Smith     /* make sure row has diagonal entry */
9706cf997b0SBarry Smith     if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) {
97129bbc08cSBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\
9726cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
9736cf997b0SBarry Smith     }
9749e25ed09SBarry Smith   }
975606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
976898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
977898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
978606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
979606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
9809e25ed09SBarry Smith 
981f64065bbSBarry Smith   {
982329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
983b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
984c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
985c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
986b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
987335d9088SBarry Smith     if (diagonal_fill) {
988b1bcba4aSBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount);
989335d9088SBarry Smith     }
990f64065bbSBarry Smith   }
991bbb6d6a8SBarry Smith 
9929e25ed09SBarry Smith   /* put together the new matrix */
993b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
994b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
995416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
996606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
9977c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
99887828ca2SBarry Smith   len = (ainew[n] + shift)*sizeof(PetscScalar);
9999e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1000606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1001606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1002b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1003416022c9SBarry Smith   b->j          = ajnew;
1004416022c9SBarry Smith   b->i          = ainew;
100556beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1006416022c9SBarry Smith   b->diag       = dloc;
1007416022c9SBarry Smith   b->ilen       = 0;
1008416022c9SBarry Smith   b->imax       = 0;
1009416022c9SBarry Smith   b->row        = isrow;
1010416022c9SBarry Smith   b->col        = iscol;
1011c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1012c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
101382bf6240SBarry Smith   b->icol       = isicol;
101487828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1015416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
101656beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
101787828ca2SBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(PetscScalar)));
1018416022c9SBarry Smith   b->maxnz          = b->nz = ainew[n] + shift;
1019a2e34c3dSBarry Smith   b->lu_damping   = info->damping;
102087828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10219e25ed09SBarry Smith   (*fact)->factor   = FACTOR_LU;
10223a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10233a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1024ae068f58SLois Curfman McInnes 
1025ae068f58SLois Curfman McInnes   (*fact)->info.factor_mallocs    = realloc;
1026ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1027329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
1028329f5518SBarry Smith   (*fact)->factor                 =  FACTOR_LU;
10293a40ed3dSBarry Smith   PetscFunctionReturn(0);
10309e25ed09SBarry Smith }
1031d5d45c9bSBarry Smith 
1032d5d45c9bSBarry Smith 
1033d5d45c9bSBarry Smith 
1034d5d45c9bSBarry Smith 
1035