xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision e2dd4fc4bb028d23435829a5bbdb20e479bb1d04)
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"
6ed480e8bSBarry Smith #include "src/inline/spops.h"
73b2fbd54SBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
1091e9ee9fSBarry Smith int MatOrdering_Flow_SeqAIJ(Mat mat,MatOrderingType type,IS *irow,IS *icol)
113b2fbd54SBarry Smith {
123a40ed3dSBarry Smith   PetscFunctionBegin;
133a40ed3dSBarry Smith 
1429bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
15aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
16d64ed03dSBarry Smith   PetscFunctionReturn(0);
17d64ed03dSBarry Smith #endif
183b2fbd54SBarry Smith }
193b2fbd54SBarry Smith 
2086bacbd2SBarry Smith 
21ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqAIJ(Mat);
223a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
2386bacbd2SBarry Smith 
2487828ca2SBarry Smith EXTERN int SPARSEKIT2dperm(int*,PetscScalar*,int*,int*,PetscScalar*,int*,int*,int*,int*,int*);
2587828ca2SBarry Smith EXTERN int SPARSEKIT2ilutp(int*,PetscScalar*,int*,int*,int*,PetscReal,PetscReal*,int*,PetscScalar*,int*,int*,int*,PetscScalar*,int*,int*,int*);
2687828ca2SBarry Smith EXTERN int SPARSEKIT2msrcsr(int*,PetscScalar*,int*,PetscScalar*,int*,int*,PetscScalar*,int*);
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
495bd2ddc7SBarry Smith      ishift = 0, for indices start at 1
505bd2ddc7SBarry Smith      ishift = 1, for indices starting at 0
5186bacbd2SBarry Smith      ------------------------------------------------------------
5286bacbd2SBarry Smith */
53b380c88cSHong Zhang int MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5486bacbd2SBarry Smith {
5560ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5660ec86bdSBarry Smith   PetscFunctionBegin;
5760ec86bdSBarry Smith   SETERRQ(1,"This distribution does not include GNU Copyright code\n\
5860ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5960ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
6060ec86bdSBarry Smith #else
6186bacbd2SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b;
6207d2056aSBarry Smith   IS           iscolf,isicol,isirow;
6386bacbd2SBarry Smith   PetscTruth   reorder;
64273d9f13SBarry Smith   int          *c,*r,*ic,ierr,i,n = A->m;
65329f5518SBarry Smith   int          *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
66313ae024SBarry Smith   int          *ordcol,*iwk,*iperm,*jw;
675bd2ddc7SBarry Smith   int          ishift = !a->indexshift;
68b19713cbSBarry Smith   int          jmax,lfill,job,*o_i,*o_j;
6987828ca2SBarry Smith   PetscScalar  *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
70329f5518SBarry Smith   PetscReal    permtol,af;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   PetscFunctionBegin;
7386bacbd2SBarry Smith 
7486bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7586bacbd2SBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax);
7686bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7733259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
786faa4630SBarry Smith   lfill   = (int)(info->dtcount/2.0);
796faa4630SBarry Smith   jmax    = (int)(info->fill*a->nz);
8086bacbd2SBarry Smith   permtol = info->dtcol;
8186bacbd2SBarry Smith 
8286bacbd2SBarry Smith 
8386bacbd2SBarry Smith   /* ------------------------------------------------------------
8486bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8586bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8686bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8786bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8886bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8986bacbd2SBarry Smith      ------------------------------------------------------------
9086bacbd2SBarry Smith   */
9186bacbd2SBarry Smith 
9286bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9386bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9486bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9586bacbd2SBarry Smith   reorder = PetscNot(reorder);
9686bacbd2SBarry Smith 
9786bacbd2SBarry Smith 
9886bacbd2SBarry Smith   /* storage for ilu factor */
99b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr);
100b0a32e0cSBarry Smith   ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr);
10187828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
102b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr);
10386bacbd2SBarry Smith 
10486bacbd2SBarry Smith   /* ------------------------------------------------------------
10586bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10686bacbd2SBarry Smith      ------------------------------------------------------------
10786bacbd2SBarry Smith   */
108b19713cbSBarry Smith   if (ishift) {
109b19713cbSBarry Smith     for (i=old_i[0];i<old_i[n];i++) {
110b19713cbSBarry Smith       old_j[i]++;
11186bacbd2SBarry Smith     }
112b19713cbSBarry Smith     for(i=0;i<n+1;i++) {
113b19713cbSBarry Smith       old_i[i]++;
114b19713cbSBarry Smith     };
11586bacbd2SBarry Smith   };
11686bacbd2SBarry Smith 
11786bacbd2SBarry Smith   if (reorder) {
118c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);           CHKERRQ(ierr);
119c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);           CHKERRQ(ierr);
120c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
121c0c2c81eSBarry Smith       r[i]  = r[i]+1;
122c0c2c81eSBarry Smith       c[i]  = c[i]+1;
123c0c2c81eSBarry Smith     }
124b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr);
125b0a32e0cSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr);
12687828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1275bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
128c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
129c0c2c81eSBarry Smith       r[i]  = r[i]-1;
130c0c2c81eSBarry Smith       c[i]  = c[i]-1;
131c0c2c81eSBarry Smith     }
132c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
133c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
134b19713cbSBarry Smith     o_a = old_a2;
135b19713cbSBarry Smith     o_j = old_j2;
136b19713cbSBarry Smith     o_i = old_i2;
137b19713cbSBarry Smith   } else {
138b19713cbSBarry Smith     o_a = old_a;
139b19713cbSBarry Smith     o_j = old_j;
140b19713cbSBarry Smith     o_i = old_i;
14186bacbd2SBarry Smith   }
14286bacbd2SBarry Smith 
14386bacbd2SBarry Smith   /* ------------------------------------------------------------
14486bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14586bacbd2SBarry Smith      ------------------------------------------------------------
14686bacbd2SBarry Smith   */
14786bacbd2SBarry Smith 
148b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr);
149b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr);
15087828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
15186bacbd2SBarry Smith 
15287828ca2SBarry 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);
15386bacbd2SBarry Smith   if (ierr) {
15486bacbd2SBarry Smith     switch (ierr) {
15529bbc08cSBarry Smith       case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15629bbc08cSBarry Smith       case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15729bbc08cSBarry Smith       case -5: SETERRQ(1,"ilutp(), zero row encountered");
15829bbc08cSBarry Smith       case -1: SETERRQ(1,"ilutp(), input matrix may be wrong");
15929bbc08cSBarry Smith       case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax);
16029bbc08cSBarry Smith       default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr);
16186bacbd2SBarry Smith     }
16286bacbd2SBarry Smith   }
16386bacbd2SBarry Smith 
16486bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16586bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16686bacbd2SBarry Smith 
16786bacbd2SBarry Smith   /* ------------------------------------------------------------
16886bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16986bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
17086bacbd2SBarry Smith      ------------------------------------------------------------
17186bacbd2SBarry Smith   */
17286bacbd2SBarry Smith 
17387828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
174b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr);
17586bacbd2SBarry Smith 
1765bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17786bacbd2SBarry Smith 
17886bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17986bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
18086bacbd2SBarry Smith 
18186bacbd2SBarry Smith   if (reorder) {
18286bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
18386bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18486bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
185313ae024SBarry Smith   } else {
186b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
187f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
188b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
189b19713cbSBarry Smith     }
190b19713cbSBarry Smith   }
191b19713cbSBarry Smith 
192b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
193b19713cbSBarry Smith   if (ishift) {
19486bacbd2SBarry Smith     for (i=0; i<n+1; i++) {
195b19713cbSBarry Smith       old_i[i]--;
196b19713cbSBarry Smith     }
197b19713cbSBarry Smith     for (i=old_i[0];i<old_i[n];i++) {
198b19713cbSBarry Smith       old_j[i]--;
199b19713cbSBarry Smith     }
20086bacbd2SBarry Smith   }
20186bacbd2SBarry Smith 
202b19713cbSBarry Smith   /* Make the factored matrix 0-based */
203b19713cbSBarry Smith   if (ishift) {
20486bacbd2SBarry Smith     for (i=0; i<n+1; i++) {
205b19713cbSBarry Smith       new_i[i]--;
206b19713cbSBarry Smith     }
207b19713cbSBarry Smith     for (i=new_i[0];i<new_i[n];i++) {
208b19713cbSBarry Smith       new_j[i]--;
209b19713cbSBarry Smith     }
21086bacbd2SBarry Smith   }
21186bacbd2SBarry Smith 
21286bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
21386bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2144c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2154c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
216c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21786bacbd2SBarry Smith   for(i=0; i<n; i++) {
21886bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21986bacbd2SBarry Smith   };
22086bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
22107d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
22286bacbd2SBarry Smith 
223c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
224c0c2c81eSBarry Smith 
225329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
22607d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
22786bacbd2SBarry Smith 
22886bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22986bacbd2SBarry Smith 
23086bacbd2SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
23186bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
23286bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
23386bacbd2SBarry Smith 
23486bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
23586bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
23686bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23707d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
238b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23986bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
24086bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
24186bacbd2SBarry Smith   b->a             = new_a;
24286bacbd2SBarry Smith   b->j             = new_j;
24386bacbd2SBarry Smith   b->i             = new_i;
24486bacbd2SBarry Smith   b->ilen          = 0;
24586bacbd2SBarry Smith   b->imax          = 0;
2461f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
247313ae024SBarry Smith   b->row           = isirow;
24886bacbd2SBarry Smith   b->col           = iscolf;
24987828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
25086bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
25186bacbd2SBarry Smith   b->indexshift    = a->indexshift;
25286bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
25386bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
25486bacbd2SBarry Smith 
25586bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
25686bacbd2SBarry Smith 
25786bacbd2SBarry Smith   /* check out for identical nodes. If found, use inode functions */
2583a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
25986bacbd2SBarry Smith 
260431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
261b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
262b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
263b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
264b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
265431e710aSBarry Smith 
26686bacbd2SBarry Smith   PetscFunctionReturn(0);
26760ec86bdSBarry Smith #endif
26886bacbd2SBarry Smith }
26986bacbd2SBarry Smith 
270289bc588SBarry Smith /*
271289bc588SBarry Smith     Factorization code for AIJ format.
272289bc588SBarry Smith */
2734a2ae208SSatish Balay #undef __FUNCT__
274b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
275b380c88cSHong Zhang int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
276289bc588SBarry Smith {
277416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
278289bc588SBarry Smith   IS         isicol;
279273d9f13SBarry Smith   int        *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j;
280416022c9SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*ajtmp,nz,shift = a->indexshift;
28191bf9adeSBarry Smith   int        *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im;
2829e046878SBarry Smith   PetscReal  f;
283289bc588SBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
28529bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2864c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
287ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
288ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
289289bc588SBarry Smith 
290289bc588SBarry Smith   /* get new row pointers */
291b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
292dbb450caSBarry Smith   ainew[0] = -shift;
293289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
294d3d32019SBarry Smith   f = info->fill;
295dbb450caSBarry Smith   jmax  = (int)(f*ai[n]+(!shift));
296b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
297289bc588SBarry Smith   /* fill is a linked list of nonzeros in active row */
298b0a32e0cSBarry Smith   ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr);
2992fb3ed76SBarry Smith   im   = fill + n + 1;
300289bc588SBarry Smith   /* idnew is location of diagonal in factor */
301b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr);
302dbb450caSBarry Smith   idnew[0] = -shift;
303289bc588SBarry Smith 
304289bc588SBarry Smith   for (i=0; i<n; i++) {
305289bc588SBarry Smith     /* first copy previous fill into linked list */
306289bc588SBarry Smith     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
30729bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
308dbb450caSBarry Smith     ajtmp   = aj + ai[r[i]] + shift;
309289bc588SBarry Smith     fill[n] = n;
310289bc588SBarry Smith     while (nz--) {
311289bc588SBarry Smith       fm  = n;
312dbb450caSBarry Smith       idx = ic[*ajtmp++ + shift];
313289bc588SBarry Smith       do {
314289bc588SBarry Smith         m  = fm;
315289bc588SBarry Smith         fm = fill[m];
316289bc588SBarry Smith       } while (fm < idx);
317289bc588SBarry Smith       fill[m]   = idx;
318289bc588SBarry Smith       fill[idx] = fm;
319289bc588SBarry Smith     }
320289bc588SBarry Smith     row = fill[n];
321289bc588SBarry Smith     while (row < i) {
322dbb450caSBarry Smith       ajtmp = ajnew + idnew[row] + (!shift);
3232fb3ed76SBarry Smith       nzbd  = 1 + idnew[row] - ainew[row];
3242fb3ed76SBarry Smith       nz    = im[row] - nzbd;
325289bc588SBarry Smith       fm    = row;
3262fb3ed76SBarry Smith       while (nz-- > 0) {
327dbb450caSBarry Smith         idx = *ajtmp++ + shift;
3282fb3ed76SBarry Smith         nzbd++;
3292fb3ed76SBarry Smith         if (idx == i) im[row] = nzbd;
330289bc588SBarry Smith         do {
331289bc588SBarry Smith           m  = fm;
332289bc588SBarry Smith           fm = fill[m];
333289bc588SBarry Smith         } while (fm < idx);
334289bc588SBarry Smith         if (fm != idx) {
335289bc588SBarry Smith           fill[m]   = idx;
336289bc588SBarry Smith           fill[idx] = fm;
337289bc588SBarry Smith           fm        = idx;
338289bc588SBarry Smith           nnz++;
339289bc588SBarry Smith         }
340289bc588SBarry Smith       }
341289bc588SBarry Smith       row = fill[row];
342289bc588SBarry Smith     }
343289bc588SBarry Smith     /* copy new filled row into permanent storage */
344289bc588SBarry Smith     ainew[i+1] = ainew[i] + nnz;
345496e697eSBarry Smith     if (ainew[i+1] > jmax) {
346ecf371e4SBarry Smith 
347ecf371e4SBarry Smith       /* estimate how much additional space we will need */
348ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
349ecf371e4SBarry Smith       /* just double the memory each time */
350ecf371e4SBarry Smith       int maxadd = jmax;
351ecf371e4SBarry Smith       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
352bbb6d6a8SBarry Smith       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
3532fb3ed76SBarry Smith       jmax += maxadd;
354ecf371e4SBarry Smith 
355ecf371e4SBarry Smith       /* allocate a longer ajnew */
356b0a32e0cSBarry Smith       ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr);
357549d3d68SSatish Balay       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i]+shift)*sizeof(int));CHKERRQ(ierr);
358606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
359289bc588SBarry Smith       ajnew = ajtmp;
3602fb3ed76SBarry Smith       realloc++; /* count how many times we realloc */
361289bc588SBarry Smith     }
362dbb450caSBarry Smith     ajtmp = ajnew + ainew[i] + shift;
363289bc588SBarry Smith     fm    = fill[n];
364289bc588SBarry Smith     nzi   = 0;
3652fb3ed76SBarry Smith     im[i] = nnz;
366289bc588SBarry Smith     while (nnz--) {
367289bc588SBarry Smith       if (fm < i) nzi++;
368dbb450caSBarry Smith       *ajtmp++ = fm - shift;
369289bc588SBarry Smith       fm       = fill[fm];
370289bc588SBarry Smith     }
371289bc588SBarry Smith     idnew[i] = ainew[i] + nzi;
372289bc588SBarry Smith   }
373464e8d44SSatish Balay   if (ai[n] != 0) {
374329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
375b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
376b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
377b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
378b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
37905bf559cSSatish Balay   } else {
380b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
38105bf559cSSatish Balay   }
3822fb3ed76SBarry Smith 
383898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
384898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3851987afe7SBarry Smith 
386606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
387289bc588SBarry Smith 
388289bc588SBarry Smith   /* put together the new matrix */
389b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr);
390b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
391416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*B)->data;
392606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3937c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
394e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
395606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
396606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
39787828ca2SBarry Smith   ierr          = PetscMalloc((ainew[n]+shift+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
398416022c9SBarry Smith   b->j          = ajnew;
399416022c9SBarry Smith   b->i          = ainew;
400416022c9SBarry Smith   b->diag       = idnew;
401416022c9SBarry Smith   b->ilen       = 0;
402416022c9SBarry Smith   b->imax       = 0;
403416022c9SBarry Smith   b->row        = isrow;
404416022c9SBarry Smith   b->col        = iscol;
405085256b3SBarry Smith   b->lu_damping        = info->damping;
40687828ca2SBarry Smith   b->lu_zeropivot      = info->zeropivot;
4076cc28720Svictorle   b->lu_shift          = info->lu_shift;
4086cc28720Svictorle   b->lu_shift_fraction = info->lu_shift_fraction;
409c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
410c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
41182bf6240SBarry Smith   b->icol       = isicol;
41287828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
413416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4147fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
41587828ca2SBarry Smith   PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(PetscScalar)));
416416022c9SBarry Smith   b->maxnz = b->nz = ainew[n] + shift;
4177fd98bd6SLois Curfman McInnes 
418329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
419ae068f58SLois Curfman McInnes   (*B)->info.factor_mallocs    = realloc;
420ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4213a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4227dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
423703deb49SSatish Balay 
424e93ccc4dSBarry Smith   if (ai[n] != 0) {
425329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
426eea2913aSSatish Balay   } else {
427eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
428eea2913aSSatish Balay   }
4293a40ed3dSBarry Smith   PetscFunctionReturn(0);
430289bc588SBarry Smith }
4310a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4323a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
43341c01911SSatish Balay 
4344a2ae208SSatish Balay #undef __FUNCT__
4354a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
436416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B)
437289bc588SBarry Smith {
438416022c9SBarry Smith   Mat          C = *B;
439416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data;
44082bf6240SBarry Smith   IS           isrow = b->row,isicol = b->icol;
441273d9f13SBarry Smith   int          *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j;
4421987afe7SBarry Smith   int          *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift;
4430cf777b8SBarry Smith   int          *diag_offset = b->diag,diag,*pj,ndamp = 0, nshift=0;
44487828ca2SBarry Smith   PetscScalar  *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4450cf777b8SBarry Smith   PetscReal    damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs,d;
4460cf777b8SBarry Smith   PetscReal    row_shift,shift_fraction,shift_amount,shift_lo=0., shift_hi=1., shift_top=0.;
4470cf777b8SBarry Smith   PetscTruth   damp,lushift;
448289bc588SBarry Smith 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
45078b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
45178b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
45287828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45387828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
4540cf60383SBarry Smith   rtmps = rtmp + shift; ics = ic + shift;
455289bc588SBarry Smith 
4566cc28720Svictorle   if (!a->diag) {
4570cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A); CHKERRQ(ierr);
4580cf777b8SBarry Smith   }
4590cf777b8SBarry Smith 
4600cf777b8SBarry Smith   if (b->lu_shift) { /* set max shift */
4610cf777b8SBarry Smith     int *aai = a->i,*ddiag = a->diag;
4620cf777b8SBarry Smith     shift_top = 0;
4636cc28720Svictorle     for (i=0; i<n; i++) {
4640cf777b8SBarry Smith       d = PetscAbsScalar((a->a)[ddiag[i]+shift]);
4656cc28720Svictorle       /* calculate amt of shift needed for this row */
46628bb5da7SSatish Balay       if (d<0) row_shift = 0; else row_shift = -2*d;
4678a2e2d9cSvictorle       v = a->a+aai[i]+shift;
4688a2e2d9cSvictorle       for (j=0; j<aai[i+1]-aai[i]; j++)
4696cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4706cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4716cc28720Svictorle     }
4726cc28720Svictorle   }
4736cc28720Svictorle 
4746cc28720Svictorle   shift_fraction = 0; shift_amount = 0;
475085256b3SBarry Smith   do {
4768a5eb818SBarry Smith     damp    = PETSC_FALSE;
4776cc28720Svictorle     lushift = PETSC_FALSE;
478289bc588SBarry Smith     for (i=0; i<n; i++) {
479289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
480dbb450caSBarry Smith       ajtmp = aj + ai[i] + shift;
48148e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
482289bc588SBarry Smith 
483289bc588SBarry Smith       /* load in initial (unfactored row) */
484416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
485416022c9SBarry Smith       ajtmpold = a->j + a->i[r[i]] + shift;
486416022c9SBarry Smith       v        = a->a + a->i[r[i]] + shift;
487085256b3SBarry Smith       for (j=0; j<nz; j++) {
488085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
489085256b3SBarry Smith       }
490*e2dd4fc4Svictorle       rtmp[ics[r[i]]] += damping + shift_amount; /* damp the diagonal of the matrix */
491289bc588SBarry Smith 
492dbb450caSBarry Smith       row = *ajtmp++ + shift;
493f2582111SSatish Balay       while  (row < i) {
4948c37ef55SBarry Smith         pc = rtmp + row;
4958c37ef55SBarry Smith         if (*pc != 0.0) {
4961987afe7SBarry Smith           pv         = b->a + diag_offset[row] + shift;
4971987afe7SBarry Smith           pj         = b->j + diag_offset[row] + (!shift);
49835aab85fSBarry Smith           multiplier = *pc / *pv++;
4998c37ef55SBarry Smith           *pc        = multiplier;
500f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
501f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
502b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
503289bc588SBarry Smith         }
504dbb450caSBarry Smith         row = *ajtmp++ + shift;
505289bc588SBarry Smith       }
506416022c9SBarry Smith       /* finished row so stick it into b->a */
507416022c9SBarry Smith       pv   = b->a + ai[i] + shift;
508416022c9SBarry Smith       pj   = b->j + ai[i] + shift;
5098c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
51035aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
511d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
512d3d32019SBarry Smith       rs   = 0.0;
513d3d32019SBarry Smith       for (j=0; j<nz; j++) {
514d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
515d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
516d3d32019SBarry Smith       }
5176cc28720Svictorle #define MAX_NSHIFT 5
5188a2e2d9cSvictorle       if (PetscRealPart(pv[diag]) < zeropivot*rs && b->lu_shift) {
5196cc28720Svictorle 	if (nshift>MAX_NSHIFT) {
5200cf777b8SBarry Smith 	  SETERRQ(1,"Unable to determine shift to enforce positive definite preconditioner");
5216cc28720Svictorle 	} else if (nshift==MAX_NSHIFT) {
5226cc28720Svictorle 	  shift_fraction = shift_hi;
5236cc28720Svictorle 	} else {
5246cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5256cc28720Svictorle 	}
5266cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5270cf777b8SBarry Smith 	lushift      = PETSC_TRUE;
5280cf777b8SBarry Smith         nshift++;
5290cf777b8SBarry Smith         break;
5306cc28720Svictorle       }
531d3d32019SBarry Smith       if (PetscAbsScalar(pv[diag]) < zeropivot*rs) {
532d3d32019SBarry Smith         if (damping) {
5338a5eb818SBarry Smith           if (ndamp) damping *= 2.0;
534085256b3SBarry Smith           damp = PETSC_TRUE;
535085256b3SBarry Smith           ndamp++;
536085256b3SBarry Smith           break;
537085256b3SBarry Smith         } else {
53863bb2aa6SBarry Smith           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
539d708749eSBarry Smith         }
54035aab85fSBarry Smith       }
54135aab85fSBarry Smith     }
5426cc28720Svictorle     if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5436cc28720Svictorle       /*
5446cc28720Svictorle        * if not already shifting up & shifting & started shifting & can refine,
5456cc28720Svictorle        * then try lower shift
5466cc28720Svictorle        */
5470cf777b8SBarry Smith       shift_hi       = shift_fraction;
5480cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5496cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5500cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5510cf777b8SBarry Smith       nshift++;
5526cc28720Svictorle     }
5536cc28720Svictorle   } while (damp || lushift);
5540f11f9aeSSatish Balay 
55535aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55635aab85fSBarry Smith   for (i=0; i<n; i++) {
55735aab85fSBarry Smith     b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift];
55835aab85fSBarry Smith   }
55935aab85fSBarry Smith 
560606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
56178b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
56278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
563416022c9SBarry Smith   C->factor = FACTOR_LU;
5647dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
565c456f294SBarry Smith   C->assembled = PETSC_TRUE;
566b0a32e0cSBarry Smith   PetscLogFlops(C->n);
567d3d32019SBarry Smith   if (ndamp) {
568b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping);
569085256b3SBarry Smith   }
5706cc28720Svictorle   if (nshift) {
5716cc28720Svictorle     b->lu_shift_fraction = shift_fraction;
5726cc28720Svictorle     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction\n",shift_fraction);
5736cc28720Svictorle   }
5743a40ed3dSBarry Smith   PetscFunctionReturn(0);
575289bc588SBarry Smith }
5760a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5774a2ae208SSatish Balay #undef __FUNCT__
5784a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
579b380c88cSHong Zhang int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
580da3a660dSBarry Smith {
581273d9f13SBarry Smith   int ierr;
582416022c9SBarry Smith   Mat C;
583416022c9SBarry Smith 
5843a40ed3dSBarry Smith   PetscFunctionBegin;
5859e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
586f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
587273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
588440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
590da3a660dSBarry Smith }
5910a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5924a2ae208SSatish Balay #undef __FUNCT__
5934a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
594416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5958c37ef55SBarry Smith {
596416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
597416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
598273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
5993b2fbd54SBarry Smith   int          nz,shift = a->indexshift,*rout,*cout;
60087828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6018c37ef55SBarry Smith 
6023a40ed3dSBarry Smith   PetscFunctionBegin;
6033a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
60491bf9adeSBarry Smith 
605e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
606e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
607416022c9SBarry Smith   tmp  = a->solve_work;
6088c37ef55SBarry Smith 
6093b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6103b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6118c37ef55SBarry Smith 
6128c37ef55SBarry Smith   /* forward solve the lower triangular */
6138c37ef55SBarry Smith   tmp[0] = b[*r++];
6140cf60383SBarry Smith   tmps   = tmp + shift;
6158c37ef55SBarry Smith   for (i=1; i<n; i++) {
616dbb450caSBarry Smith     v   = aa + ai[i] + shift;
617dbb450caSBarry Smith     vi  = aj + ai[i] + shift;
61853e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
61953e7425aSBarry Smith     sum = b[*r++];
620ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6218c37ef55SBarry Smith     tmp[i] = sum;
6228c37ef55SBarry Smith   }
6238c37ef55SBarry Smith 
6248c37ef55SBarry Smith   /* backward solve the upper triangular */
6258c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
626416022c9SBarry Smith     v   = aa + a->diag[i] + (!shift);
627416022c9SBarry Smith     vi  = aj + a->diag[i] + (!shift);
628416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6298c37ef55SBarry Smith     sum = tmp[i];
630ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
631416022c9SBarry Smith     x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift];
6328c37ef55SBarry Smith   }
6338c37ef55SBarry Smith 
6343b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6353b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
636cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
637cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
638b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6393a40ed3dSBarry Smith   PetscFunctionReturn(0);
6408c37ef55SBarry Smith }
641026e39d0SSatish Balay 
642930ae53cSBarry Smith /* ----------------------------------------------------------- */
6434a2ae208SSatish Balay #undef __FUNCT__
6444a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
645930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
646930ae53cSBarry Smith {
647930ae53cSBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
648273d9f13SBarry Smith   int          n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr;
649362ced78SSatish Balay   PetscScalar  *x,*b,*aa = a->a;
650aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
651d85afc42SSatish Balay   int          adiag_i,i,*vi,nz,ai_i;
652362ced78SSatish Balay   PetscScalar  *v,sum;
653d85afc42SSatish Balay #endif
654930ae53cSBarry Smith 
6553a40ed3dSBarry Smith   PetscFunctionBegin;
6563a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
657930ae53cSBarry Smith   if (a->indexshift) {
6583a40ed3dSBarry Smith      ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr);
6593a40ed3dSBarry Smith      PetscFunctionReturn(0);
660930ae53cSBarry Smith   }
661930ae53cSBarry Smith 
662e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
663e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
664930ae53cSBarry Smith 
665aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6661c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6671c4feecaSSatish Balay #else
668930ae53cSBarry Smith   /* forward solve the lower triangular */
669930ae53cSBarry Smith   x[0] = b[0];
670930ae53cSBarry Smith   for (i=1; i<n; i++) {
671930ae53cSBarry Smith     ai_i = ai[i];
672930ae53cSBarry Smith     v    = aa + ai_i;
673930ae53cSBarry Smith     vi   = aj + ai_i;
674930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
675930ae53cSBarry Smith     sum  = b[i];
676930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
677930ae53cSBarry Smith     x[i] = sum;
678930ae53cSBarry Smith   }
679930ae53cSBarry Smith 
680930ae53cSBarry Smith   /* backward solve the upper triangular */
681930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
682930ae53cSBarry Smith     adiag_i = adiag[i];
683d708749eSBarry Smith     v       = aa + adiag_i + 1;
684d708749eSBarry Smith     vi      = aj + adiag_i + 1;
685930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
686930ae53cSBarry Smith     sum     = x[i];
687930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
688930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
689930ae53cSBarry Smith   }
6901c4feecaSSatish Balay #endif
691b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
692cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
693cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
695930ae53cSBarry Smith }
696930ae53cSBarry Smith 
6974a2ae208SSatish Balay #undef __FUNCT__
6984a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
699416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
700da3a660dSBarry Smith {
701416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
702416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
703273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
7043b2fbd54SBarry Smith   int          nz,shift = a->indexshift,*rout,*cout;
70587828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,sum,*v;
706da3a660dSBarry Smith 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
70878b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
709da3a660dSBarry Smith 
710e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
711e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
712416022c9SBarry Smith   tmp  = a->solve_work;
713da3a660dSBarry Smith 
7143b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7153b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
716da3a660dSBarry Smith 
717da3a660dSBarry Smith   /* forward solve the lower triangular */
718da3a660dSBarry Smith   tmp[0] = b[*r++];
719da3a660dSBarry Smith   for (i=1; i<n; i++) {
720dbb450caSBarry Smith     v   = aa + ai[i] + shift;
721dbb450caSBarry Smith     vi  = aj + ai[i] + shift;
722416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
723da3a660dSBarry Smith     sum = b[*r++];
724dbb450caSBarry Smith     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
725da3a660dSBarry Smith     tmp[i] = sum;
726da3a660dSBarry Smith   }
727da3a660dSBarry Smith 
728da3a660dSBarry Smith   /* backward solve the upper triangular */
729da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
730416022c9SBarry Smith     v   = aa + a->diag[i] + (!shift);
731416022c9SBarry Smith     vi  = aj + a->diag[i] + (!shift);
732416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
733da3a660dSBarry Smith     sum = tmp[i];
734dbb450caSBarry Smith     while (nz--) sum -= *v++ * tmp[*vi++ + shift];
735416022c9SBarry Smith     tmp[i] = sum*aa[a->diag[i]+shift];
736da3a660dSBarry Smith     x[*c--] += tmp[i];
737da3a660dSBarry Smith   }
738da3a660dSBarry Smith 
7393b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7403b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
741cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
742cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
743b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
744898183ecSLois Curfman McInnes 
7453a40ed3dSBarry Smith   PetscFunctionReturn(0);
746da3a660dSBarry Smith }
747da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7484a2ae208SSatish Balay #undef __FUNCT__
7494a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
7507c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
751da3a660dSBarry Smith {
752416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
7532235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
754273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
755f1af5d2fSBarry Smith   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
75687828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v,s1;
757da3a660dSBarry Smith 
7583a40ed3dSBarry Smith   PetscFunctionBegin;
759e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
760e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
761416022c9SBarry Smith   tmp  = a->solve_work;
762da3a660dSBarry Smith 
7632235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7642235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
765da3a660dSBarry Smith 
766da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7672235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
768da3a660dSBarry Smith 
769da3a660dSBarry Smith   /* forward solve the U^T */
770da3a660dSBarry Smith   for (i=0; i<n; i++) {
771f1af5d2fSBarry Smith     v   = aa + diag[i] + shift;
772f1af5d2fSBarry Smith     vi  = aj + diag[i] + (!shift);
773f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
774c38d4ed2SBarry Smith     s1  = tmp[i];
775ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
776da3a660dSBarry Smith     while (nz--) {
777f1af5d2fSBarry Smith       tmp[*vi++ + shift] -= (*v++)*s1;
778da3a660dSBarry Smith     }
779c38d4ed2SBarry Smith     tmp[i] = s1;
780da3a660dSBarry Smith   }
781da3a660dSBarry Smith 
782da3a660dSBarry Smith   /* backward solve the L^T */
783da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
784f1af5d2fSBarry Smith     v   = aa + diag[i] - 1 + shift;
785f1af5d2fSBarry Smith     vi  = aj + diag[i] - 1 + shift;
786f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
787f1af5d2fSBarry Smith     s1  = tmp[i];
788da3a660dSBarry Smith     while (nz--) {
789f1af5d2fSBarry Smith       tmp[*vi-- + shift] -= (*v--)*s1;
790da3a660dSBarry Smith     }
791da3a660dSBarry Smith   }
792da3a660dSBarry Smith 
793da3a660dSBarry Smith   /* copy tmp into x according to permutation */
794da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
795da3a660dSBarry Smith 
7962235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7972235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
798cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
799cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
800da3a660dSBarry Smith 
801b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8023a40ed3dSBarry Smith   PetscFunctionReturn(0);
803da3a660dSBarry Smith }
804da3a660dSBarry Smith 
8054a2ae208SSatish Balay #undef __FUNCT__
8064a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
8077c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
808da3a660dSBarry Smith {
809416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
8102235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
811273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
812f1af5d2fSBarry Smith   int          nz,shift = a->indexshift,*rout,*cout,*diag = a->diag;
81387828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v;
8146abc6512SBarry Smith 
8153a40ed3dSBarry Smith   PetscFunctionBegin;
8166abc6512SBarry Smith   if (zz != xx) VecCopy(zz,xx);
8176abc6512SBarry Smith 
818e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
819e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
820416022c9SBarry Smith   tmp = a->solve_work;
8216abc6512SBarry Smith 
8222235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8232235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8246abc6512SBarry Smith 
8256abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8262235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8276abc6512SBarry Smith 
8286abc6512SBarry Smith   /* forward solve the U^T */
8296abc6512SBarry Smith   for (i=0; i<n; i++) {
830f1af5d2fSBarry Smith     v   = aa + diag[i] + shift;
831f1af5d2fSBarry Smith     vi  = aj + diag[i] + (!shift);
832f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8336abc6512SBarry Smith     tmp[i] *= *v++;
8346abc6512SBarry Smith     while (nz--) {
835dbb450caSBarry Smith       tmp[*vi++ + shift] -= (*v++)*tmp[i];
8366abc6512SBarry Smith     }
8376abc6512SBarry Smith   }
8386abc6512SBarry Smith 
8396abc6512SBarry Smith   /* backward solve the L^T */
8406abc6512SBarry Smith   for (i=n-1; i>=0; i--){
841f1af5d2fSBarry Smith     v   = aa + diag[i] - 1 + shift;
842f1af5d2fSBarry Smith     vi  = aj + diag[i] - 1 + shift;
843f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8446abc6512SBarry Smith     while (nz--) {
845dbb450caSBarry Smith       tmp[*vi-- + shift] -= (*v--)*tmp[i];
8466abc6512SBarry Smith     }
8476abc6512SBarry Smith   }
8486abc6512SBarry Smith 
8496abc6512SBarry Smith   /* copy tmp into x according to permutation */
8506abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8516abc6512SBarry Smith 
8522235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8532235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
854cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
855cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8566abc6512SBarry Smith 
857b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8583a40ed3dSBarry Smith   PetscFunctionReturn(0);
859da3a660dSBarry Smith }
8609e25ed09SBarry Smith /* ----------------------------------------------------------------*/
861ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat);
86208480c60SBarry Smith 
8634a2ae208SSatish Balay #undef __FUNCT__
8644a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
865b380c88cSHong Zhang int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8669e25ed09SBarry Smith {
867416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
8689e25ed09SBarry Smith   IS         isicol;
869273d9f13SBarry Smith   int        *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j;
87056beaf04SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
871335d9088SBarry Smith   int        *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0;
8726cf997b0SBarry Smith   int        incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill;
87377c4ece6SBarry Smith   PetscTruth col_identity,row_identity;
874329f5518SBarry Smith   PetscReal  f;
8759e25ed09SBarry Smith 
8763a40ed3dSBarry Smith   PetscFunctionBegin;
8776cf997b0SBarry Smith   f             = info->fill;
8780cd17293SBarry Smith   levels        = (int)info->levels;
8790cd17293SBarry Smith   diagonal_fill = (int)info->diagonal_fill;
8804c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
88182bf6240SBarry Smith 
88208480c60SBarry Smith   /* special case that simply copies fill pattern */
883be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
884be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
88586bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8862e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
88708480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
88808480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
88908480c60SBarry Smith     if (!b->diag) {
8907c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89108480c60SBarry Smith     }
8927c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89308480c60SBarry Smith     b->row              = isrow;
89408480c60SBarry Smith     b->col              = iscol;
89582bf6240SBarry Smith     b->icol             = isicol;
896a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
89787828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
8986cc28720Svictorle     b->lu_shift         = info->lu_shift;
8996cc28720Svictorle     b->lu_shift_fraction= info->lu_shift_fraction;
90087828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
901f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
902c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
903c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9043a40ed3dSBarry Smith     PetscFunctionReturn(0);
90508480c60SBarry Smith   }
90608480c60SBarry Smith 
907898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
908898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9099e25ed09SBarry Smith 
9109e25ed09SBarry Smith   /* get new row pointers */
911b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
912dbb450caSBarry Smith   ainew[0] = -shift;
9139e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
914dbb450caSBarry Smith   jmax = (int)(f*(ai[n]+!shift));
915b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
9169e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
917b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr);
9189e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
919b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr);
92056beaf04SBarry Smith   /* im is level for each filled value */
921b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr);
92256beaf04SBarry Smith   /* dloc is location of diagonal in factor */
923b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr);
92456beaf04SBarry Smith   dloc[0]  = 0;
92556beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9266cf997b0SBarry Smith 
9276cf997b0SBarry Smith     /* copy current row into linked list */
92856beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
92929bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
930dbb450caSBarry Smith     xi      = aj + ai[r[prow]] + shift;
9319e25ed09SBarry Smith     fill[n]    = n;
932435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9339e25ed09SBarry Smith     while (nz--) {
9349e25ed09SBarry Smith       fm  = n;
935dbb450caSBarry Smith       idx = ic[*xi++ + shift];
9369e25ed09SBarry Smith       do {
9379e25ed09SBarry Smith         m  = fm;
9389e25ed09SBarry Smith         fm = fill[m];
9399e25ed09SBarry Smith       } while (fm < idx);
9409e25ed09SBarry Smith       fill[m]   = idx;
9419e25ed09SBarry Smith       fill[idx] = fm;
94256beaf04SBarry Smith       im[idx]   = 0;
9439e25ed09SBarry Smith     }
9446cf997b0SBarry Smith 
9456cf997b0SBarry Smith     /* make sure diagonal entry is included */
946435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9476cf997b0SBarry Smith       fm = n;
948435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
949435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
950435faa5fSBarry Smith       fill[fm]   = prow;
9516cf997b0SBarry Smith       im[prow]   = 0;
9526cf997b0SBarry Smith       nzf++;
953335d9088SBarry Smith       dcount++;
9546cf997b0SBarry Smith     }
9556cf997b0SBarry Smith 
95656beaf04SBarry Smith     nzi = 0;
9579e25ed09SBarry Smith     row = fill[n];
95856beaf04SBarry Smith     while (row < prow) {
95956beaf04SBarry Smith       incrlev = im[row] + 1;
96056beaf04SBarry Smith       nz      = dloc[row];
9616cf997b0SBarry Smith       xi      = ajnew  + ainew[row] + shift + nz + 1;
962dbb450caSBarry Smith       flev    = ajfill + ainew[row] + shift + nz + 1;
963416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
96456beaf04SBarry Smith       fm      = row;
96556beaf04SBarry Smith       while (nnz-- > 0) {
966dbb450caSBarry Smith         idx = *xi++ + shift;
96756beaf04SBarry Smith         if (*flev + incrlev > levels) {
96856beaf04SBarry Smith           flev++;
96956beaf04SBarry Smith           continue;
97056beaf04SBarry Smith         }
97156beaf04SBarry Smith         do {
9729e25ed09SBarry Smith           m  = fm;
9739e25ed09SBarry Smith           fm = fill[m];
97456beaf04SBarry Smith         } while (fm < idx);
9759e25ed09SBarry Smith         if (fm != idx) {
97656beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9779e25ed09SBarry Smith           fill[m]   = idx;
9789e25ed09SBarry Smith           fill[idx] = fm;
9799e25ed09SBarry Smith           fm        = idx;
98056beaf04SBarry Smith           nzf++;
981ecf371e4SBarry Smith         } else {
98256beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9839e25ed09SBarry Smith         }
98456beaf04SBarry Smith         flev++;
9859e25ed09SBarry Smith       }
9869e25ed09SBarry Smith       row = fill[row];
98756beaf04SBarry Smith       nzi++;
9889e25ed09SBarry Smith     }
9899e25ed09SBarry Smith     /* copy new filled row into permanent storage */
99056beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
991d7e8b826SBarry Smith     if (ainew[prow+1] > jmax-shift) {
992ecf371e4SBarry Smith 
993ecf371e4SBarry Smith       /* estimate how much additional space we will need */
994ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
995ecf371e4SBarry Smith       /* just double the memory each time */
996ecf371e4SBarry Smith       /*  maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */
997ecf371e4SBarry Smith       int maxadd = jmax;
99856beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
999bbb6d6a8SBarry Smith       jmax += maxadd;
1000ecf371e4SBarry Smith 
1001ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
1002b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
1003549d3d68SSatish Balay       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
1004606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
100556beaf04SBarry Smith       ajnew  = xi;
1006b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
1007549d3d68SSatish Balay       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr);
1008606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
100956beaf04SBarry Smith       ajfill = xi;
1010ecf371e4SBarry Smith       realloc++; /* count how many times we realloc */
10119e25ed09SBarry Smith     }
1012dbb450caSBarry Smith     xi          = ajnew + ainew[prow] + shift;
1013dbb450caSBarry Smith     flev        = ajfill + ainew[prow] + shift;
101456beaf04SBarry Smith     dloc[prow]  = nzi;
10159e25ed09SBarry Smith     fm          = fill[n];
101656beaf04SBarry Smith     while (nzf--) {
1017dbb450caSBarry Smith       *xi++   = fm - shift;
101856beaf04SBarry Smith       *flev++ = im[fm];
10199e25ed09SBarry Smith       fm      = fill[fm];
10209e25ed09SBarry Smith     }
10216cf997b0SBarry Smith     /* make sure row has diagonal entry */
10226cf997b0SBarry Smith     if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) {
102329bbc08cSBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\
10246cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10256cf997b0SBarry Smith     }
10269e25ed09SBarry Smith   }
1027606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1028898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1029898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1030606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1031606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10329e25ed09SBarry Smith 
1033f64065bbSBarry Smith   {
1034329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1035b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
1036c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1037c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1038b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1039335d9088SBarry Smith     if (diagonal_fill) {
1040b1bcba4aSBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount);
1041335d9088SBarry Smith     }
1042f64065bbSBarry Smith   }
1043bbb6d6a8SBarry Smith 
10449e25ed09SBarry Smith   /* put together the new matrix */
1045b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
1046b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1047416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1048606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10497c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
105087828ca2SBarry Smith   len = (ainew[n] + shift)*sizeof(PetscScalar);
10519e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1052606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1053606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1054b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1055416022c9SBarry Smith   b->j          = ajnew;
1056416022c9SBarry Smith   b->i          = ainew;
105756beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1058416022c9SBarry Smith   b->diag       = dloc;
1059416022c9SBarry Smith   b->ilen       = 0;
1060416022c9SBarry Smith   b->imax       = 0;
1061416022c9SBarry Smith   b->row        = isrow;
1062416022c9SBarry Smith   b->col        = iscol;
1063c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1064c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
106582bf6240SBarry Smith   b->icol       = isicol;
106687828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1067416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
106856beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
106987828ca2SBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(PetscScalar)));
1070416022c9SBarry Smith   b->maxnz        = b->nz = ainew[n] + shift;
1071a2e34c3dSBarry Smith   b->lu_damping   = info->damping;
107287828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10739e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10743a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10753a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1076ae068f58SLois Curfman McInnes 
1077ae068f58SLois Curfman McInnes   (*fact)->info.factor_mallocs    = realloc;
1078ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1079329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10803a40ed3dSBarry Smith   PetscFunctionReturn(0);
10819e25ed09SBarry Smith }
1082d5d45c9bSBarry Smith 
1083a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1084a6175056SHong Zhang #undef __FUNCT__
1085a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1086a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact)
1087a6175056SHong Zhang {
10880968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1089a6175056SHong Zhang   int                 ierr;
1090d5d45c9bSBarry Smith 
1091a6175056SHong Zhang   PetscFunctionBegin;
10920968510aSHong Zhang   if (!a->sbaijMat){
10930968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1094a6175056SHong Zhang   }
109503aac1b8SHong Zhang 
1096b45a75daSHong Zhang   ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr);
10970968510aSHong Zhang   ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
1098064503c5SHong Zhang   a->sbaijMat = PETSC_NULL;
1099653a6975SHong Zhang 
1100a6175056SHong Zhang   PetscFunctionReturn(0);
1101a6175056SHong Zhang }
1102a6175056SHong Zhang 
1103a6175056SHong Zhang #undef __FUNCT__
1104a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
110515e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1106a6175056SHong Zhang {
11070968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
110803aac1b8SHong Zhang   int                 ierr;
1109653a6975SHong Zhang   PetscTruth          perm_identity;
1110a6175056SHong Zhang 
1111a6175056SHong Zhang   PetscFunctionBegin;
1112653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1113653a6975SHong Zhang   if (!perm_identity){
1114653a6975SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1115653a6975SHong Zhang   }
11160968510aSHong Zhang   if (!a->sbaijMat){
11170968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
11180968510aSHong Zhang   }
1119a6175056SHong Zhang 
1120b00f7748SHong Zhang   ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1121653a6975SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1122653a6975SHong Zhang 
1123a6175056SHong Zhang   PetscFunctionReturn(0);
1124a6175056SHong Zhang }
1125d5d45c9bSBarry Smith 
1126f76d2b81SHong Zhang #undef __FUNCT__
1127f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1128f76d2b81SHong Zhang int MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1129f76d2b81SHong Zhang {
1130f76d2b81SHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
113103aac1b8SHong Zhang   int                 ierr;
1132f76d2b81SHong Zhang   PetscTruth          perm_identity;
1133f76d2b81SHong Zhang 
1134f76d2b81SHong Zhang   PetscFunctionBegin;
1135f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1136f76d2b81SHong Zhang   if (!perm_identity){
1137f76d2b81SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1138f76d2b81SHong Zhang   }
1139f76d2b81SHong Zhang   if (!a->sbaijMat){
1140f76d2b81SHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1141f76d2b81SHong Zhang   }
1142f76d2b81SHong Zhang 
1143f76d2b81SHong Zhang   ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1144f76d2b81SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1145f76d2b81SHong Zhang 
1146f76d2b81SHong Zhang   PetscFunctionReturn(0);
1147f76d2b81SHong Zhang }
1148