xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 010a20ca055c35d67bd4f840442d6048eb8be78a)
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 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
51b380c88cSHong Zhang int MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
5560ec86bdSBarry Smith   SETERRQ(1,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS           iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth   reorder;
62273d9f13SBarry Smith   int          *c,*r,*ic,ierr,i,n = A->m;
63329f5518SBarry Smith   int          *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
64313ae024SBarry Smith   int          *ordcol,*iwk,*iperm,*jw;
65b19713cbSBarry Smith   int          jmax,lfill,job,*o_i,*o_j;
6687828ca2SBarry Smith   PetscScalar  *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
67329f5518SBarry Smith   PetscReal    permtol,af;
6886bacbd2SBarry Smith 
6986bacbd2SBarry Smith   PetscFunctionBegin;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7286bacbd2SBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax);
7386bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7433259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
756faa4630SBarry Smith   lfill   = (int)(info->dtcount/2.0);
766faa4630SBarry Smith   jmax    = (int)(info->fill*a->nz);
7786bacbd2SBarry Smith   permtol = info->dtcol;
7886bacbd2SBarry Smith 
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith   /* ------------------------------------------------------------
8186bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8286bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8386bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8486bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8586bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8686bacbd2SBarry Smith      ------------------------------------------------------------
8786bacbd2SBarry Smith   */
8886bacbd2SBarry Smith 
8986bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9086bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9186bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9286bacbd2SBarry Smith   reorder = PetscNot(reorder);
9386bacbd2SBarry Smith 
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith   /* storage for ilu factor */
96b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr);
97b0a32e0cSBarry Smith   ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr);
9887828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
99b0a32e0cSBarry Smith   ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr);
10086bacbd2SBarry Smith 
10186bacbd2SBarry Smith   /* ------------------------------------------------------------
10286bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10386bacbd2SBarry Smith      ------------------------------------------------------------
10486bacbd2SBarry Smith   */
105b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
106b19713cbSBarry Smith     old_j[i]++;
10786bacbd2SBarry Smith   }
108b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
109b19713cbSBarry Smith     old_i[i]++;
110b19713cbSBarry Smith   };
111*010a20caSHong Zhang 
11286bacbd2SBarry Smith 
11386bacbd2SBarry Smith   if (reorder) {
114c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);           CHKERRQ(ierr);
115c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);           CHKERRQ(ierr);
116c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
117c0c2c81eSBarry Smith       r[i]  = r[i]+1;
118c0c2c81eSBarry Smith       c[i]  = c[i]+1;
119c0c2c81eSBarry Smith     }
120b0a32e0cSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr);
121b0a32e0cSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr);
12287828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1235bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
124c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
125c0c2c81eSBarry Smith       r[i]  = r[i]-1;
126c0c2c81eSBarry Smith       c[i]  = c[i]-1;
127c0c2c81eSBarry Smith     }
128c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
129c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
130b19713cbSBarry Smith     o_a = old_a2;
131b19713cbSBarry Smith     o_j = old_j2;
132b19713cbSBarry Smith     o_i = old_i2;
133b19713cbSBarry Smith   } else {
134b19713cbSBarry Smith     o_a = old_a;
135b19713cbSBarry Smith     o_j = old_j;
136b19713cbSBarry Smith     o_i = old_i;
13786bacbd2SBarry Smith   }
13886bacbd2SBarry Smith 
13986bacbd2SBarry Smith   /* ------------------------------------------------------------
14086bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14186bacbd2SBarry Smith      ------------------------------------------------------------
14286bacbd2SBarry Smith   */
14386bacbd2SBarry Smith 
144b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr);
145b0a32e0cSBarry Smith   ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr);
14687828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14786bacbd2SBarry Smith 
14887828ca2SBarry 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);
14986bacbd2SBarry Smith   if (ierr) {
15086bacbd2SBarry Smith     switch (ierr) {
15129bbc08cSBarry Smith       case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15229bbc08cSBarry Smith       case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax);
15329bbc08cSBarry Smith       case -5: SETERRQ(1,"ilutp(), zero row encountered");
15429bbc08cSBarry Smith       case -1: SETERRQ(1,"ilutp(), input matrix may be wrong");
15529bbc08cSBarry Smith       case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax);
15629bbc08cSBarry Smith       default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr);
15786bacbd2SBarry Smith     }
15886bacbd2SBarry Smith   }
15986bacbd2SBarry Smith 
16086bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16186bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16286bacbd2SBarry Smith 
16386bacbd2SBarry Smith   /* ------------------------------------------------------------
16486bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16586bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16686bacbd2SBarry Smith      ------------------------------------------------------------
16786bacbd2SBarry Smith   */
16886bacbd2SBarry Smith 
16987828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
170b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr);
17186bacbd2SBarry Smith 
1725bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17386bacbd2SBarry Smith 
17486bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17586bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17686bacbd2SBarry Smith 
17786bacbd2SBarry Smith   if (reorder) {
17886bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18086bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
181313ae024SBarry Smith   } else {
182b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
183f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
184b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
185b19713cbSBarry Smith     }
186b19713cbSBarry Smith   }
187b19713cbSBarry Smith 
188b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
190b19713cbSBarry Smith     old_i[i]--;
191b19713cbSBarry Smith   }
192b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
193b19713cbSBarry Smith     old_j[i]--;
194b19713cbSBarry Smith   }
19586bacbd2SBarry Smith 
196b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
198b19713cbSBarry Smith     new_i[i]--;
199b19713cbSBarry Smith   }
200b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
201b19713cbSBarry Smith     new_j[i]--;
202b19713cbSBarry Smith   }
20386bacbd2SBarry Smith 
20486bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20586bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2064c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2074c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
208c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20986bacbd2SBarry Smith   for(i=0; i<n; i++) {
21086bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21186bacbd2SBarry Smith   };
21286bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21307d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21486bacbd2SBarry Smith 
215c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
216c0c2c81eSBarry Smith 
217329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21807d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21986bacbd2SBarry Smith 
22086bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22186bacbd2SBarry Smith 
22286bacbd2SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
22386bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22486bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22586bacbd2SBarry Smith 
22686bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
22786bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
22886bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
22907d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
230b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23186bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
23286bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
23386bacbd2SBarry Smith   b->a             = new_a;
23486bacbd2SBarry Smith   b->j             = new_j;
23586bacbd2SBarry Smith   b->i             = new_i;
23686bacbd2SBarry Smith   b->ilen          = 0;
23786bacbd2SBarry Smith   b->imax          = 0;
2381f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
239313ae024SBarry Smith   b->row           = isirow;
24086bacbd2SBarry Smith   b->col           = iscolf;
24187828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24286bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24386bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24486bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24586bacbd2SBarry Smith 
24686bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24786bacbd2SBarry Smith 
24886bacbd2SBarry Smith   /* check out for identical nodes. If found, use inode functions */
2493a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
25086bacbd2SBarry Smith 
251431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
252b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
253b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
254b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
255b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
256431e710aSBarry Smith 
25786bacbd2SBarry Smith   PetscFunctionReturn(0);
25860ec86bdSBarry Smith #endif
25986bacbd2SBarry Smith }
26086bacbd2SBarry Smith 
261289bc588SBarry Smith /*
262289bc588SBarry Smith     Factorization code for AIJ format.
263289bc588SBarry Smith */
2644a2ae208SSatish Balay #undef __FUNCT__
265b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
266b380c88cSHong Zhang int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
267289bc588SBarry Smith {
268416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
269289bc588SBarry Smith   IS         isicol;
270273d9f13SBarry Smith   int        *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j;
271*010a20caSHong Zhang   int        *ainew,*ajnew,jmax,*fill,*ajtmp,nz;
27291bf9adeSBarry Smith   int        *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im;
2739e046878SBarry Smith   PetscReal  f;
274289bc588SBarry Smith 
2753a40ed3dSBarry Smith   PetscFunctionBegin;
27629bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2774c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
278ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
279ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
280289bc588SBarry Smith 
281289bc588SBarry Smith   /* get new row pointers */
282b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
283*010a20caSHong Zhang   ainew[0] = 0;
284289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
285d3d32019SBarry Smith   f = info->fill;
286*010a20caSHong Zhang   jmax  = (int)(f*ai[n]+1);
287b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
288289bc588SBarry Smith   /* fill is a linked list of nonzeros in active row */
289b0a32e0cSBarry Smith   ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr);
2902fb3ed76SBarry Smith   im   = fill + n + 1;
291289bc588SBarry Smith   /* idnew is location of diagonal in factor */
292b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr);
293*010a20caSHong Zhang   idnew[0] = 0;
294289bc588SBarry Smith 
295289bc588SBarry Smith   for (i=0; i<n; i++) {
296289bc588SBarry Smith     /* first copy previous fill into linked list */
297289bc588SBarry Smith     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
29829bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
299*010a20caSHong Zhang     ajtmp   = aj + ai[r[i]];
300289bc588SBarry Smith     fill[n] = n;
301289bc588SBarry Smith     while (nz--) {
302289bc588SBarry Smith       fm  = n;
303*010a20caSHong Zhang       idx = ic[*ajtmp++];
304289bc588SBarry Smith       do {
305289bc588SBarry Smith         m  = fm;
306289bc588SBarry Smith         fm = fill[m];
307289bc588SBarry Smith       } while (fm < idx);
308289bc588SBarry Smith       fill[m]   = idx;
309289bc588SBarry Smith       fill[idx] = fm;
310289bc588SBarry Smith     }
311289bc588SBarry Smith     row = fill[n];
312289bc588SBarry Smith     while (row < i) {
313*010a20caSHong Zhang       ajtmp = ajnew + idnew[row] + 1;
3142fb3ed76SBarry Smith       nzbd  = 1 + idnew[row] - ainew[row];
3152fb3ed76SBarry Smith       nz    = im[row] - nzbd;
316289bc588SBarry Smith       fm    = row;
3172fb3ed76SBarry Smith       while (nz-- > 0) {
318*010a20caSHong Zhang         idx = *ajtmp++ ;
3192fb3ed76SBarry Smith         nzbd++;
3202fb3ed76SBarry Smith         if (idx == i) im[row] = nzbd;
321289bc588SBarry Smith         do {
322289bc588SBarry Smith           m  = fm;
323289bc588SBarry Smith           fm = fill[m];
324289bc588SBarry Smith         } while (fm < idx);
325289bc588SBarry Smith         if (fm != idx) {
326289bc588SBarry Smith           fill[m]   = idx;
327289bc588SBarry Smith           fill[idx] = fm;
328289bc588SBarry Smith           fm        = idx;
329289bc588SBarry Smith           nnz++;
330289bc588SBarry Smith         }
331289bc588SBarry Smith       }
332289bc588SBarry Smith       row = fill[row];
333289bc588SBarry Smith     }
334289bc588SBarry Smith     /* copy new filled row into permanent storage */
335289bc588SBarry Smith     ainew[i+1] = ainew[i] + nnz;
336496e697eSBarry Smith     if (ainew[i+1] > jmax) {
337ecf371e4SBarry Smith 
338ecf371e4SBarry Smith       /* estimate how much additional space we will need */
339ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
340ecf371e4SBarry Smith       /* just double the memory each time */
341ecf371e4SBarry Smith       int maxadd = jmax;
342ecf371e4SBarry Smith       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
343bbb6d6a8SBarry Smith       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
3442fb3ed76SBarry Smith       jmax += maxadd;
345ecf371e4SBarry Smith 
346ecf371e4SBarry Smith       /* allocate a longer ajnew */
347b0a32e0cSBarry Smith       ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr);
348*010a20caSHong Zhang       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i])*sizeof(int));CHKERRQ(ierr);
349606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
350289bc588SBarry Smith       ajnew = ajtmp;
3512fb3ed76SBarry Smith       realloc++; /* count how many times we realloc */
352289bc588SBarry Smith     }
353*010a20caSHong Zhang     ajtmp = ajnew + ainew[i];
354289bc588SBarry Smith     fm    = fill[n];
355289bc588SBarry Smith     nzi   = 0;
3562fb3ed76SBarry Smith     im[i] = nnz;
357289bc588SBarry Smith     while (nnz--) {
358289bc588SBarry Smith       if (fm < i) nzi++;
359*010a20caSHong Zhang       *ajtmp++ = fm ;
360289bc588SBarry Smith       fm       = fill[fm];
361289bc588SBarry Smith     }
362289bc588SBarry Smith     idnew[i] = ainew[i] + nzi;
363289bc588SBarry Smith   }
364464e8d44SSatish Balay   if (ai[n] != 0) {
365329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
366b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
367b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
368b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
369b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
37005bf559cSSatish Balay   } else {
371b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
37205bf559cSSatish Balay   }
3732fb3ed76SBarry Smith 
374898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
375898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3761987afe7SBarry Smith 
377606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
378289bc588SBarry Smith 
379289bc588SBarry Smith   /* put together the new matrix */
380b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr);
381b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
382416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*B)->data;
383606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3847c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
385e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
386606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
387606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
388*010a20caSHong Zhang   ierr          = PetscMalloc((ainew[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
389416022c9SBarry Smith   b->j          = ajnew;
390416022c9SBarry Smith   b->i          = ainew;
391416022c9SBarry Smith   b->diag       = idnew;
392416022c9SBarry Smith   b->ilen       = 0;
393416022c9SBarry Smith   b->imax       = 0;
394416022c9SBarry Smith   b->row        = isrow;
395416022c9SBarry Smith   b->col        = iscol;
396085256b3SBarry Smith   b->lu_damping        = info->damping;
39787828ca2SBarry Smith   b->lu_zeropivot      = info->zeropivot;
3986cc28720Svictorle   b->lu_shift          = info->lu_shift;
3996cc28720Svictorle   b->lu_shift_fraction = info->lu_shift_fraction;
400c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
401c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40282bf6240SBarry Smith   b->icol       = isicol;
40387828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
404416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4057fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
406*010a20caSHong Zhang   PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(int)+sizeof(PetscScalar)));
407*010a20caSHong Zhang   b->maxnz = b->nz = ainew[n] ;
4087fd98bd6SLois Curfman McInnes 
409329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
410ae068f58SLois Curfman McInnes   (*B)->info.factor_mallocs    = realloc;
411ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4123a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4137dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
414703deb49SSatish Balay 
415e93ccc4dSBarry Smith   if (ai[n] != 0) {
416329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
417eea2913aSSatish Balay   } else {
418eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
419eea2913aSSatish Balay   }
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
421289bc588SBarry Smith }
4220a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4233a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
42441c01911SSatish Balay 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
427416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B)
428289bc588SBarry Smith {
429416022c9SBarry Smith   Mat          C = *B;
430416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data;
43182bf6240SBarry Smith   IS           isrow = b->row,isicol = b->icol;
432273d9f13SBarry Smith   int          *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j;
433*010a20caSHong Zhang   int          *ajtmpold,*ajtmp,nz,row,*ics;
4340cf777b8SBarry Smith   int          *diag_offset = b->diag,diag,*pj,ndamp = 0, nshift=0;
43587828ca2SBarry Smith   PetscScalar  *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4360cf777b8SBarry Smith   PetscReal    damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs,d;
4370cf777b8SBarry Smith   PetscReal    row_shift,shift_fraction,shift_amount,shift_lo=0., shift_hi=1., shift_top=0.;
4380cf777b8SBarry Smith   PetscTruth   damp,lushift;
439289bc588SBarry Smith 
4403a40ed3dSBarry Smith   PetscFunctionBegin;
44178b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44278b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44387828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44487828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
445*010a20caSHong Zhang   rtmps = rtmp; ics = ic;
446289bc588SBarry Smith 
4476cc28720Svictorle   if (!a->diag) {
4480cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A); CHKERRQ(ierr);
4490cf777b8SBarry Smith   }
4500cf777b8SBarry Smith 
4510cf777b8SBarry Smith   if (b->lu_shift) { /* set max shift */
4520cf777b8SBarry Smith     int *aai = a->i,*ddiag = a->diag;
4530cf777b8SBarry Smith     shift_top = 0;
4546cc28720Svictorle     for (i=0; i<n; i++) {
455*010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4566cc28720Svictorle       /* calculate amt of shift needed for this row */
4573aeef889SHong Zhang       if (d<0) {
4583aeef889SHong Zhang         row_shift = 0;
4593aeef889SHong Zhang       } else {
4603aeef889SHong Zhang         row_shift = -2*d;
4613aeef889SHong Zhang       }
462*010a20caSHong Zhang       v = a->a+aai[i];
4638a2e2d9cSvictorle       for (j=0; j<aai[i+1]-aai[i]; j++)
4646cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4656cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4666cc28720Svictorle     }
4676cc28720Svictorle   }
4686cc28720Svictorle 
4696cc28720Svictorle   shift_fraction = 0; shift_amount = 0;
470085256b3SBarry Smith   do {
4718a5eb818SBarry Smith     damp    = PETSC_FALSE;
4726cc28720Svictorle     lushift = PETSC_FALSE;
473289bc588SBarry Smith     for (i=0; i<n; i++) {
474289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
475*010a20caSHong Zhang       ajtmp = aj + ai[i];
47648e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
477289bc588SBarry Smith 
478289bc588SBarry Smith       /* load in initial (unfactored row) */
479416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
480*010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
481*010a20caSHong Zhang       v        = a->a + a->i[r[i]];
482085256b3SBarry Smith       for (j=0; j<nz; j++) {
483085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
484085256b3SBarry Smith       }
485e2dd4fc4Svictorle       rtmp[ics[r[i]]] += damping + shift_amount; /* damp the diagonal of the matrix */
486289bc588SBarry Smith 
487*010a20caSHong Zhang       row = *ajtmp++ ;
488f2582111SSatish Balay       while  (row < i) {
4898c37ef55SBarry Smith         pc = rtmp + row;
4908c37ef55SBarry Smith         if (*pc != 0.0) {
491*010a20caSHong Zhang           pv         = b->a + diag_offset[row] ;
492*010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
49335aab85fSBarry Smith           multiplier = *pc / *pv++;
4948c37ef55SBarry Smith           *pc        = multiplier;
495f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
496f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
497b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
498289bc588SBarry Smith         }
499*010a20caSHong Zhang         row = *ajtmp++;
500289bc588SBarry Smith       }
501416022c9SBarry Smith       /* finished row so stick it into b->a */
502*010a20caSHong Zhang       pv   = b->a + ai[i] ;
503*010a20caSHong Zhang       pj   = b->j + ai[i] ;
5048c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
50535aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
506d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
507d3d32019SBarry Smith       rs   = 0.0;
508d3d32019SBarry Smith       for (j=0; j<nz; j++) {
509d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
510d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
511d3d32019SBarry Smith       }
5126cc28720Svictorle #define MAX_NSHIFT 5
5138a2e2d9cSvictorle       if (PetscRealPart(pv[diag]) < zeropivot*rs && b->lu_shift) {
5146cc28720Svictorle 	if (nshift>MAX_NSHIFT) {
5150cf777b8SBarry Smith 	  SETERRQ(1,"Unable to determine shift to enforce positive definite preconditioner");
5166cc28720Svictorle 	} else if (nshift==MAX_NSHIFT) {
5176cc28720Svictorle 	  shift_fraction = shift_hi;
5186cc28720Svictorle 	} else {
5196cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5206cc28720Svictorle 	}
5216cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5220cf777b8SBarry Smith 	lushift      = PETSC_TRUE;
5230cf777b8SBarry Smith         nshift++;
5240cf777b8SBarry Smith         break;
5256cc28720Svictorle       }
526d3d32019SBarry Smith       if (PetscAbsScalar(pv[diag]) < zeropivot*rs) {
527d3d32019SBarry Smith         if (damping) {
5288a5eb818SBarry Smith           if (ndamp) damping *= 2.0;
529085256b3SBarry Smith           damp = PETSC_TRUE;
530085256b3SBarry Smith           ndamp++;
531085256b3SBarry Smith           break;
532085256b3SBarry Smith         } else {
53363bb2aa6SBarry Smith           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
534d708749eSBarry Smith         }
53535aab85fSBarry Smith       }
53635aab85fSBarry Smith     }
5376cc28720Svictorle     if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5386cc28720Svictorle       /*
5396cc28720Svictorle        * if not already shifting up & shifting & started shifting & can refine,
5406cc28720Svictorle        * then try lower shift
5416cc28720Svictorle        */
5420cf777b8SBarry Smith       shift_hi       = shift_fraction;
5430cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5446cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5450cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5460cf777b8SBarry Smith       nshift++;
5476cc28720Svictorle     }
5486cc28720Svictorle   } while (damp || lushift);
5490f11f9aeSSatish Balay 
55035aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55135aab85fSBarry Smith   for (i=0; i<n; i++) {
552*010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55335aab85fSBarry Smith   }
55435aab85fSBarry Smith 
555606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55678b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
558416022c9SBarry Smith   C->factor = FACTOR_LU;
5597dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
560c456f294SBarry Smith   C->assembled = PETSC_TRUE;
561b0a32e0cSBarry Smith   PetscLogFlops(C->n);
562d3d32019SBarry Smith   if (ndamp) {
563b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping);
564085256b3SBarry Smith   }
5656cc28720Svictorle   if (nshift) {
5666cc28720Svictorle     b->lu_shift_fraction = shift_fraction;
5676cc28720Svictorle     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction\n",shift_fraction);
5686cc28720Svictorle   }
5693a40ed3dSBarry Smith   PetscFunctionReturn(0);
570289bc588SBarry Smith }
5710a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5724a2ae208SSatish Balay #undef __FUNCT__
5734a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
574b380c88cSHong Zhang int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
575da3a660dSBarry Smith {
576273d9f13SBarry Smith   int ierr;
577416022c9SBarry Smith   Mat C;
578416022c9SBarry Smith 
5793a40ed3dSBarry Smith   PetscFunctionBegin;
5809e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
581f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
582273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
583440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5843a40ed3dSBarry Smith   PetscFunctionReturn(0);
585da3a660dSBarry Smith }
5860a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5874a2ae208SSatish Balay #undef __FUNCT__
5884a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
589416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5908c37ef55SBarry Smith {
591416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
592416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
593273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
594*010a20caSHong Zhang   int          nz,*rout,*cout;
59587828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
5968c37ef55SBarry Smith 
5973a40ed3dSBarry Smith   PetscFunctionBegin;
5983a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
59991bf9adeSBarry Smith 
600e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
601e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
602416022c9SBarry Smith   tmp  = a->solve_work;
6038c37ef55SBarry Smith 
6043b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6053b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6068c37ef55SBarry Smith 
6078c37ef55SBarry Smith   /* forward solve the lower triangular */
6088c37ef55SBarry Smith   tmp[0] = b[*r++];
609*010a20caSHong Zhang   tmps   = tmp;
6108c37ef55SBarry Smith   for (i=1; i<n; i++) {
611*010a20caSHong Zhang     v   = aa + ai[i] ;
612*010a20caSHong Zhang     vi  = aj + ai[i] ;
61353e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
61453e7425aSBarry Smith     sum = b[*r++];
615ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6168c37ef55SBarry Smith     tmp[i] = sum;
6178c37ef55SBarry Smith   }
6188c37ef55SBarry Smith 
6198c37ef55SBarry Smith   /* backward solve the upper triangular */
6208c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
621*010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
622*010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
623416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6248c37ef55SBarry Smith     sum = tmp[i];
625ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
626*010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6278c37ef55SBarry Smith   }
6288c37ef55SBarry Smith 
6293b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6303b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
631cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
632cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
633b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6343a40ed3dSBarry Smith   PetscFunctionReturn(0);
6358c37ef55SBarry Smith }
636026e39d0SSatish Balay 
637930ae53cSBarry Smith /* ----------------------------------------------------------- */
6384a2ae208SSatish Balay #undef __FUNCT__
6394a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
640930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
641930ae53cSBarry Smith {
642930ae53cSBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
643273d9f13SBarry Smith   int          n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr;
644362ced78SSatish Balay   PetscScalar  *x,*b,*aa = a->a;
645aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
646d85afc42SSatish Balay   int          adiag_i,i,*vi,nz,ai_i;
647362ced78SSatish Balay   PetscScalar  *v,sum;
648d85afc42SSatish Balay #endif
649930ae53cSBarry Smith 
6503a40ed3dSBarry Smith   PetscFunctionBegin;
6513a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
652930ae53cSBarry Smith 
653e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
654e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
655930ae53cSBarry Smith 
656aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6571c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6581c4feecaSSatish Balay #else
659930ae53cSBarry Smith   /* forward solve the lower triangular */
660930ae53cSBarry Smith   x[0] = b[0];
661930ae53cSBarry Smith   for (i=1; i<n; i++) {
662930ae53cSBarry Smith     ai_i = ai[i];
663930ae53cSBarry Smith     v    = aa + ai_i;
664930ae53cSBarry Smith     vi   = aj + ai_i;
665930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
666930ae53cSBarry Smith     sum  = b[i];
667930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
668930ae53cSBarry Smith     x[i] = sum;
669930ae53cSBarry Smith   }
670930ae53cSBarry Smith 
671930ae53cSBarry Smith   /* backward solve the upper triangular */
672930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
673930ae53cSBarry Smith     adiag_i = adiag[i];
674d708749eSBarry Smith     v       = aa + adiag_i + 1;
675d708749eSBarry Smith     vi      = aj + adiag_i + 1;
676930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
677930ae53cSBarry Smith     sum     = x[i];
678930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
679930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
680930ae53cSBarry Smith   }
6811c4feecaSSatish Balay #endif
682b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
683cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
684cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6853a40ed3dSBarry Smith   PetscFunctionReturn(0);
686930ae53cSBarry Smith }
687930ae53cSBarry Smith 
6884a2ae208SSatish Balay #undef __FUNCT__
6894a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
690416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
691da3a660dSBarry Smith {
692416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
693416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
694273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
695*010a20caSHong Zhang   int          nz,*rout,*cout;
69687828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,sum,*v;
697da3a660dSBarry Smith 
6983a40ed3dSBarry Smith   PetscFunctionBegin;
69978b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
700da3a660dSBarry Smith 
701e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
702e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
703416022c9SBarry Smith   tmp  = a->solve_work;
704da3a660dSBarry Smith 
7053b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7063b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
707da3a660dSBarry Smith 
708da3a660dSBarry Smith   /* forward solve the lower triangular */
709da3a660dSBarry Smith   tmp[0] = b[*r++];
710da3a660dSBarry Smith   for (i=1; i<n; i++) {
711*010a20caSHong Zhang     v   = aa + ai[i] ;
712*010a20caSHong Zhang     vi  = aj + ai[i] ;
713416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
714da3a660dSBarry Smith     sum = b[*r++];
715*010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
716da3a660dSBarry Smith     tmp[i] = sum;
717da3a660dSBarry Smith   }
718da3a660dSBarry Smith 
719da3a660dSBarry Smith   /* backward solve the upper triangular */
720da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
721*010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
722*010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
723416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
724da3a660dSBarry Smith     sum = tmp[i];
725*010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
726*010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
727da3a660dSBarry Smith     x[*c--] += tmp[i];
728da3a660dSBarry Smith   }
729da3a660dSBarry Smith 
7303b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7313b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
732cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
733cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
734b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
735898183ecSLois Curfman McInnes 
7363a40ed3dSBarry Smith   PetscFunctionReturn(0);
737da3a660dSBarry Smith }
738da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7394a2ae208SSatish Balay #undef __FUNCT__
7404a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
7417c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
742da3a660dSBarry Smith {
743416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
7442235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
745273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
746*010a20caSHong Zhang   int          nz,*rout,*cout,*diag = a->diag;
74787828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v,s1;
748da3a660dSBarry Smith 
7493a40ed3dSBarry Smith   PetscFunctionBegin;
750e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
751e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
752416022c9SBarry Smith   tmp  = a->solve_work;
753da3a660dSBarry Smith 
7542235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7552235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
756da3a660dSBarry Smith 
757da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7582235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
759da3a660dSBarry Smith 
760da3a660dSBarry Smith   /* forward solve the U^T */
761da3a660dSBarry Smith   for (i=0; i<n; i++) {
762*010a20caSHong Zhang     v   = aa + diag[i] ;
763*010a20caSHong Zhang     vi  = aj + diag[i] + 1;
764f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
765c38d4ed2SBarry Smith     s1  = tmp[i];
766ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
767da3a660dSBarry Smith     while (nz--) {
768*010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
769da3a660dSBarry Smith     }
770c38d4ed2SBarry Smith     tmp[i] = s1;
771da3a660dSBarry Smith   }
772da3a660dSBarry Smith 
773da3a660dSBarry Smith   /* backward solve the L^T */
774da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
775*010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
776*010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
777f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
778f1af5d2fSBarry Smith     s1  = tmp[i];
779da3a660dSBarry Smith     while (nz--) {
780*010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
781da3a660dSBarry Smith     }
782da3a660dSBarry Smith   }
783da3a660dSBarry Smith 
784da3a660dSBarry Smith   /* copy tmp into x according to permutation */
785da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
786da3a660dSBarry Smith 
7872235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7882235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
789cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
790cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
791da3a660dSBarry Smith 
792b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
7933a40ed3dSBarry Smith   PetscFunctionReturn(0);
794da3a660dSBarry Smith }
795da3a660dSBarry Smith 
7964a2ae208SSatish Balay #undef __FUNCT__
7974a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
7987c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
799da3a660dSBarry Smith {
800416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
8012235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
802273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
803*010a20caSHong Zhang   int          nz,*rout,*cout,*diag = a->diag;
80487828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v;
8056abc6512SBarry Smith 
8063a40ed3dSBarry Smith   PetscFunctionBegin;
8076abc6512SBarry Smith   if (zz != xx) VecCopy(zz,xx);
8086abc6512SBarry Smith 
809e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
810e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
811416022c9SBarry Smith   tmp = a->solve_work;
8126abc6512SBarry Smith 
8132235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8142235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8156abc6512SBarry Smith 
8166abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8172235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8186abc6512SBarry Smith 
8196abc6512SBarry Smith   /* forward solve the U^T */
8206abc6512SBarry Smith   for (i=0; i<n; i++) {
821*010a20caSHong Zhang     v   = aa + diag[i] ;
822*010a20caSHong Zhang     vi  = aj + diag[i] + 1;
823f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8246abc6512SBarry Smith     tmp[i] *= *v++;
8256abc6512SBarry Smith     while (nz--) {
826*010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8276abc6512SBarry Smith     }
8286abc6512SBarry Smith   }
8296abc6512SBarry Smith 
8306abc6512SBarry Smith   /* backward solve the L^T */
8316abc6512SBarry Smith   for (i=n-1; i>=0; i--){
832*010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
833*010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
834f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8356abc6512SBarry Smith     while (nz--) {
836*010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8376abc6512SBarry Smith     }
8386abc6512SBarry Smith   }
8396abc6512SBarry Smith 
8406abc6512SBarry Smith   /* copy tmp into x according to permutation */
8416abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8426abc6512SBarry Smith 
8432235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8442235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
845cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
846cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8476abc6512SBarry Smith 
848b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8493a40ed3dSBarry Smith   PetscFunctionReturn(0);
850da3a660dSBarry Smith }
8519e25ed09SBarry Smith /* ----------------------------------------------------------------*/
852ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat);
85308480c60SBarry Smith 
8544a2ae208SSatish Balay #undef __FUNCT__
8554a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
856b380c88cSHong Zhang int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8579e25ed09SBarry Smith {
858416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
8599e25ed09SBarry Smith   IS         isicol;
860273d9f13SBarry Smith   int        *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j;
86156beaf04SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
862335d9088SBarry Smith   int        *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0;
863*010a20caSHong Zhang   int        incrlev,nnz,i,levels,diagonal_fill;
86477c4ece6SBarry Smith   PetscTruth col_identity,row_identity;
865329f5518SBarry Smith   PetscReal  f;
8669e25ed09SBarry Smith 
8673a40ed3dSBarry Smith   PetscFunctionBegin;
8686cf997b0SBarry Smith   f             = info->fill;
8690cd17293SBarry Smith   levels        = (int)info->levels;
8700cd17293SBarry Smith   diagonal_fill = (int)info->diagonal_fill;
8714c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
87282bf6240SBarry Smith 
87308480c60SBarry Smith   /* special case that simply copies fill pattern */
874be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
875be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
87686bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8772e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
87808480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
87908480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
88008480c60SBarry Smith     if (!b->diag) {
8817c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
88208480c60SBarry Smith     }
8837c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
88408480c60SBarry Smith     b->row              = isrow;
88508480c60SBarry Smith     b->col              = iscol;
88682bf6240SBarry Smith     b->icol             = isicol;
887a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
88887828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
8896cc28720Svictorle     b->lu_shift         = info->lu_shift;
8906cc28720Svictorle     b->lu_shift_fraction= info->lu_shift_fraction;
89187828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
892f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
893c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
894c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
8953a40ed3dSBarry Smith     PetscFunctionReturn(0);
89608480c60SBarry Smith   }
89708480c60SBarry Smith 
898898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
899898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9009e25ed09SBarry Smith 
9019e25ed09SBarry Smith   /* get new row pointers */
902b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
903*010a20caSHong Zhang   ainew[0] = 0;
9049e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
905*010a20caSHong Zhang   jmax = (int)(f*(ai[n]+1));
906b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
9079e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
908b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr);
9099e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
910b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr);
91156beaf04SBarry Smith   /* im is level for each filled value */
912b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr);
91356beaf04SBarry Smith   /* dloc is location of diagonal in factor */
914b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr);
91556beaf04SBarry Smith   dloc[0]  = 0;
91656beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9176cf997b0SBarry Smith 
9186cf997b0SBarry Smith     /* copy current row into linked list */
91956beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
92029bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
921*010a20caSHong Zhang     xi      = aj + ai[r[prow]] ;
9229e25ed09SBarry Smith     fill[n]    = n;
923435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9249e25ed09SBarry Smith     while (nz--) {
9259e25ed09SBarry Smith       fm  = n;
926*010a20caSHong Zhang       idx = ic[*xi++ ];
9279e25ed09SBarry Smith       do {
9289e25ed09SBarry Smith         m  = fm;
9299e25ed09SBarry Smith         fm = fill[m];
9309e25ed09SBarry Smith       } while (fm < idx);
9319e25ed09SBarry Smith       fill[m]   = idx;
9329e25ed09SBarry Smith       fill[idx] = fm;
93356beaf04SBarry Smith       im[idx]   = 0;
9349e25ed09SBarry Smith     }
9356cf997b0SBarry Smith 
9366cf997b0SBarry Smith     /* make sure diagonal entry is included */
937435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9386cf997b0SBarry Smith       fm = n;
939435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
940435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
941435faa5fSBarry Smith       fill[fm]   = prow;
9426cf997b0SBarry Smith       im[prow]   = 0;
9436cf997b0SBarry Smith       nzf++;
944335d9088SBarry Smith       dcount++;
9456cf997b0SBarry Smith     }
9466cf997b0SBarry Smith 
94756beaf04SBarry Smith     nzi = 0;
9489e25ed09SBarry Smith     row = fill[n];
94956beaf04SBarry Smith     while (row < prow) {
95056beaf04SBarry Smith       incrlev = im[row] + 1;
95156beaf04SBarry Smith       nz      = dloc[row];
952*010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
953*010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
954416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
95556beaf04SBarry Smith       fm      = row;
95656beaf04SBarry Smith       while (nnz-- > 0) {
957*010a20caSHong Zhang         idx = *xi++ ;
95856beaf04SBarry Smith         if (*flev + incrlev > levels) {
95956beaf04SBarry Smith           flev++;
96056beaf04SBarry Smith           continue;
96156beaf04SBarry Smith         }
96256beaf04SBarry Smith         do {
9639e25ed09SBarry Smith           m  = fm;
9649e25ed09SBarry Smith           fm = fill[m];
96556beaf04SBarry Smith         } while (fm < idx);
9669e25ed09SBarry Smith         if (fm != idx) {
96756beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9689e25ed09SBarry Smith           fill[m]   = idx;
9699e25ed09SBarry Smith           fill[idx] = fm;
9709e25ed09SBarry Smith           fm        = idx;
97156beaf04SBarry Smith           nzf++;
972ecf371e4SBarry Smith         } else {
97356beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9749e25ed09SBarry Smith         }
97556beaf04SBarry Smith         flev++;
9769e25ed09SBarry Smith       }
9779e25ed09SBarry Smith       row = fill[row];
97856beaf04SBarry Smith       nzi++;
9799e25ed09SBarry Smith     }
9809e25ed09SBarry Smith     /* copy new filled row into permanent storage */
98156beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
982*010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
983ecf371e4SBarry Smith 
984ecf371e4SBarry Smith       /* estimate how much additional space we will need */
985ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
986ecf371e4SBarry Smith       /* just double the memory each time */
987ecf371e4SBarry Smith       /*  maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */
988ecf371e4SBarry Smith       int maxadd = jmax;
98956beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
990bbb6d6a8SBarry Smith       jmax += maxadd;
991ecf371e4SBarry Smith 
992ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
993b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
994*010a20caSHong Zhang       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(int));CHKERRQ(ierr);
995606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
99656beaf04SBarry Smith       ajnew  = xi;
997b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
998*010a20caSHong Zhang       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(int));CHKERRQ(ierr);
999606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
100056beaf04SBarry Smith       ajfill = xi;
1001ecf371e4SBarry Smith       realloc++; /* count how many times we realloc */
10029e25ed09SBarry Smith     }
1003*010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1004*010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
100556beaf04SBarry Smith     dloc[prow]  = nzi;
10069e25ed09SBarry Smith     fm          = fill[n];
100756beaf04SBarry Smith     while (nzf--) {
1008*010a20caSHong Zhang       *xi++   = fm ;
100956beaf04SBarry Smith       *flev++ = im[fm];
10109e25ed09SBarry Smith       fm      = fill[fm];
10119e25ed09SBarry Smith     }
10126cf997b0SBarry Smith     /* make sure row has diagonal entry */
1013*010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
101429bbc08cSBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\
10156cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10166cf997b0SBarry Smith     }
10179e25ed09SBarry Smith   }
1018606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1019898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1020898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1021606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1022606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10239e25ed09SBarry Smith 
1024f64065bbSBarry Smith   {
1025329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1026b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
1027c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1028c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1029b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1030335d9088SBarry Smith     if (diagonal_fill) {
1031b1bcba4aSBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount);
1032335d9088SBarry Smith     }
1033f64065bbSBarry Smith   }
1034bbb6d6a8SBarry Smith 
10359e25ed09SBarry Smith   /* put together the new matrix */
1036b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
1037b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1038416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1039606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10407c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1041*010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10429e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1043606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1044606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1045b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1046416022c9SBarry Smith   b->j          = ajnew;
1047416022c9SBarry Smith   b->i          = ainew;
104856beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1049416022c9SBarry Smith   b->diag       = dloc;
1050416022c9SBarry Smith   b->ilen       = 0;
1051416022c9SBarry Smith   b->imax       = 0;
1052416022c9SBarry Smith   b->row        = isrow;
1053416022c9SBarry Smith   b->col        = iscol;
1054c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1055c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
105682bf6240SBarry Smith   b->icol       = isicol;
105787828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1058416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
105956beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
1060*010a20caSHong Zhang   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(int)+sizeof(PetscScalar)));
1061*010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
1062a2e34c3dSBarry Smith   b->lu_damping        = info->damping;
10639fc2f943Svictorle   b->lu_shift          = info->lu_shift;
10649fc2f943Svictorle   b->lu_shift_fraction = info->lu_shift_fraction;
106587828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10669e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10673a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10683a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1069ae068f58SLois Curfman McInnes 
1070ae068f58SLois Curfman McInnes   (*fact)->info.factor_mallocs    = realloc;
1071ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1072329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10733a40ed3dSBarry Smith   PetscFunctionReturn(0);
10749e25ed09SBarry Smith }
1075d5d45c9bSBarry Smith 
1076a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1077a6175056SHong Zhang #undef __FUNCT__
1078a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1079a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact)
1080a6175056SHong Zhang {
10810968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1082a6175056SHong Zhang   int                 ierr;
1083d5d45c9bSBarry Smith 
1084a6175056SHong Zhang   PetscFunctionBegin;
10850968510aSHong Zhang   if (!a->sbaijMat){
10860968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1087a6175056SHong Zhang   }
108803aac1b8SHong Zhang 
1089b45a75daSHong Zhang   ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr);
10900968510aSHong Zhang   ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
1091064503c5SHong Zhang   a->sbaijMat = PETSC_NULL;
1092653a6975SHong Zhang 
1093a6175056SHong Zhang   PetscFunctionReturn(0);
1094a6175056SHong Zhang }
1095a6175056SHong Zhang 
1096a6175056SHong Zhang #undef __FUNCT__
1097a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
109815e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1099a6175056SHong Zhang {
11000968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
110103aac1b8SHong Zhang   int                 ierr;
1102653a6975SHong Zhang   PetscTruth          perm_identity;
1103a6175056SHong Zhang 
1104a6175056SHong Zhang   PetscFunctionBegin;
1105653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1106653a6975SHong Zhang   if (!perm_identity){
1107653a6975SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1108653a6975SHong Zhang   }
11090968510aSHong Zhang   if (!a->sbaijMat){
11100968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
11110968510aSHong Zhang   }
1112a6175056SHong Zhang 
1113b00f7748SHong Zhang   ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1114653a6975SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1115653a6975SHong Zhang 
1116a6175056SHong Zhang   PetscFunctionReturn(0);
1117a6175056SHong Zhang }
1118d5d45c9bSBarry Smith 
1119f76d2b81SHong Zhang #undef __FUNCT__
1120f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1121f76d2b81SHong Zhang int MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1122f76d2b81SHong Zhang {
1123f76d2b81SHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
112403aac1b8SHong Zhang   int                 ierr;
1125f76d2b81SHong Zhang   PetscTruth          perm_identity;
1126f76d2b81SHong Zhang 
1127f76d2b81SHong Zhang   PetscFunctionBegin;
1128f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1129f76d2b81SHong Zhang   if (!perm_identity){
1130f76d2b81SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1131f76d2b81SHong Zhang   }
1132f76d2b81SHong Zhang   if (!a->sbaijMat){
1133f76d2b81SHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1134f76d2b81SHong Zhang   }
1135f76d2b81SHong Zhang 
1136f76d2b81SHong Zhang   ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1137f76d2b81SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1138f76d2b81SHong Zhang 
1139f76d2b81SHong Zhang   PetscFunctionReturn(0);
1140f76d2b81SHong Zhang }
1141