xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision fb10cecf1927e102c85c17a39b04f58da85ca3ce)
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   };
111010a20caSHong 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;
271010a20caSHong 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);
283010a20caSHong Zhang   ainew[0] = 0;
284289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
285d3d32019SBarry Smith   f = info->fill;
286010a20caSHong 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);
293010a20caSHong 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");
299010a20caSHong Zhang     ajtmp   = aj + ai[r[i]];
300289bc588SBarry Smith     fill[n] = n;
301289bc588SBarry Smith     while (nz--) {
302289bc588SBarry Smith       fm  = n;
303010a20caSHong 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) {
313010a20caSHong 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) {
318010a20caSHong 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);
348010a20caSHong 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     }
353010a20caSHong 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++;
359010a20caSHong 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);
388010a20caSHong 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;
3982cea7109SBarry Smith   b->lu_shift          = info->shift;
3992cea7109SBarry Smith   b->lu_shift_fraction = info->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 */
406010a20caSHong Zhang   PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(int)+sizeof(PetscScalar)));
407010a20caSHong 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;
433010a20caSHong 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);
445010a20caSHong 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++) {
455010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4566cc28720Svictorle       /* calculate amt of shift needed for this row */
4576c037d1bSvictorle       if (d<=0) {
4583aeef889SHong Zhang         row_shift = 0;
4593aeef889SHong Zhang       } else {
4603aeef889SHong Zhang         row_shift = -2*d;
4613aeef889SHong Zhang       }
462010a20caSHong 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];
475010a20caSHong 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]];
480010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
481010a20caSHong 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 
487010a20caSHong Zhang       row = *ajtmp++ ;
488f2582111SSatish Balay       while  (row < i) {
4898c37ef55SBarry Smith         pc = rtmp + row;
4908c37ef55SBarry Smith         if (*pc != 0.0) {
491010a20caSHong Zhang           pv         = b->a + diag_offset[row] ;
492010a20caSHong 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         }
499010a20caSHong Zhang         row = *ajtmp++;
500289bc588SBarry Smith       }
501416022c9SBarry Smith       /* finished row so stick it into b->a */
502010a20caSHong Zhang       pv   = b->a + ai[i] ;
503010a20caSHong 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;
5186c037d1bSvictorle 	  lushift      = PETSC_FALSE;
5196cc28720Svictorle 	} else {
5206cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5216c037d1bSvictorle 	  lushift      = PETSC_TRUE;
5226cc28720Svictorle 	}
5236cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5240cf777b8SBarry Smith         nshift++;
5250cf777b8SBarry Smith         break;
5266cc28720Svictorle       }
527d3d32019SBarry Smith       if (PetscAbsScalar(pv[diag]) < zeropivot*rs) {
528d3d32019SBarry Smith         if (damping) {
5298a5eb818SBarry Smith           if (ndamp) damping *= 2.0;
530085256b3SBarry Smith           damp = PETSC_TRUE;
531085256b3SBarry Smith           ndamp++;
532085256b3SBarry Smith           break;
533085256b3SBarry Smith         } else {
53463bb2aa6SBarry Smith           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
535d708749eSBarry Smith         }
53635aab85fSBarry Smith       }
53735aab85fSBarry Smith     }
5386cc28720Svictorle     if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5396cc28720Svictorle       /*
5406c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5416cc28720Svictorle        * then try lower shift
5426cc28720Svictorle        */
5430cf777b8SBarry Smith       shift_hi       = shift_fraction;
5440cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5456cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5460cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5470cf777b8SBarry Smith       nshift++;
5486cc28720Svictorle     }
5496cc28720Svictorle   } while (damp || lushift);
5500f11f9aeSSatish Balay 
55135aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55235aab85fSBarry Smith   for (i=0; i<n; i++) {
553010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55435aab85fSBarry Smith   }
55535aab85fSBarry Smith 
556606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55778b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
559416022c9SBarry Smith   C->factor = FACTOR_LU;
5607dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
561c456f294SBarry Smith   C->assembled = PETSC_TRUE;
562b0a32e0cSBarry Smith   PetscLogFlops(C->n);
563d3d32019SBarry Smith   if (ndamp) {
564b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping);
565085256b3SBarry Smith   }
5666cc28720Svictorle   if (nshift) {
5676cc28720Svictorle     b->lu_shift_fraction = shift_fraction;
568*fb10cecfSBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %d\n",shift_fraction,shift_top,nshift);
5696cc28720Svictorle   }
5703a40ed3dSBarry Smith   PetscFunctionReturn(0);
571289bc588SBarry Smith }
5720889b5dcSvictorle 
5730889b5dcSvictorle #undef __FUNCT__
5740889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
5750889b5dcSvictorle int MatUsePETSc_SeqAIJ(Mat A)
5760889b5dcSvictorle {
5770889b5dcSvictorle   PetscFunctionBegin;
5780889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5790889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5800889b5dcSvictorle   PetscFunctionReturn(0);
5810889b5dcSvictorle }
5820889b5dcSvictorle 
5830889b5dcSvictorle 
5840a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5854a2ae208SSatish Balay #undef __FUNCT__
5864a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
587b380c88cSHong Zhang int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
588da3a660dSBarry Smith {
589273d9f13SBarry Smith   int ierr;
590416022c9SBarry Smith   Mat C;
591416022c9SBarry Smith 
5923a40ed3dSBarry Smith   PetscFunctionBegin;
5939e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
594f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
595273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
596440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5973a40ed3dSBarry Smith   PetscFunctionReturn(0);
598da3a660dSBarry Smith }
5990a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
6004a2ae208SSatish Balay #undef __FUNCT__
6014a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
602416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6038c37ef55SBarry Smith {
604416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
605416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
606273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
607010a20caSHong Zhang   int          nz,*rout,*cout;
60887828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6098c37ef55SBarry Smith 
6103a40ed3dSBarry Smith   PetscFunctionBegin;
6113a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
61291bf9adeSBarry Smith 
613e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
614e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
615416022c9SBarry Smith   tmp  = a->solve_work;
6168c37ef55SBarry Smith 
6173b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6183b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6198c37ef55SBarry Smith 
6208c37ef55SBarry Smith   /* forward solve the lower triangular */
6218c37ef55SBarry Smith   tmp[0] = b[*r++];
622010a20caSHong Zhang   tmps   = tmp;
6238c37ef55SBarry Smith   for (i=1; i<n; i++) {
624010a20caSHong Zhang     v   = aa + ai[i] ;
625010a20caSHong Zhang     vi  = aj + ai[i] ;
62653e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
62753e7425aSBarry Smith     sum = b[*r++];
628ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6298c37ef55SBarry Smith     tmp[i] = sum;
6308c37ef55SBarry Smith   }
6318c37ef55SBarry Smith 
6328c37ef55SBarry Smith   /* backward solve the upper triangular */
6338c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
634010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
635010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
636416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6378c37ef55SBarry Smith     sum = tmp[i];
638ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
639010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6408c37ef55SBarry Smith   }
6418c37ef55SBarry Smith 
6423b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6433b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
644cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
645cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
646b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6473a40ed3dSBarry Smith   PetscFunctionReturn(0);
6488c37ef55SBarry Smith }
649026e39d0SSatish Balay 
650930ae53cSBarry Smith /* ----------------------------------------------------------- */
6514a2ae208SSatish Balay #undef __FUNCT__
6524a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
653930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
654930ae53cSBarry Smith {
655930ae53cSBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
656273d9f13SBarry Smith   int          n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr;
657362ced78SSatish Balay   PetscScalar  *x,*b,*aa = a->a;
658aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
659d85afc42SSatish Balay   int          adiag_i,i,*vi,nz,ai_i;
660362ced78SSatish Balay   PetscScalar  *v,sum;
661d85afc42SSatish Balay #endif
662930ae53cSBarry Smith 
6633a40ed3dSBarry Smith   PetscFunctionBegin;
6643a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
665930ae53cSBarry Smith 
666e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
667e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
668930ae53cSBarry Smith 
669aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6701c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6711c4feecaSSatish Balay #else
672930ae53cSBarry Smith   /* forward solve the lower triangular */
673930ae53cSBarry Smith   x[0] = b[0];
674930ae53cSBarry Smith   for (i=1; i<n; i++) {
675930ae53cSBarry Smith     ai_i = ai[i];
676930ae53cSBarry Smith     v    = aa + ai_i;
677930ae53cSBarry Smith     vi   = aj + ai_i;
678930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
679930ae53cSBarry Smith     sum  = b[i];
680930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
681930ae53cSBarry Smith     x[i] = sum;
682930ae53cSBarry Smith   }
683930ae53cSBarry Smith 
684930ae53cSBarry Smith   /* backward solve the upper triangular */
685930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
686930ae53cSBarry Smith     adiag_i = adiag[i];
687d708749eSBarry Smith     v       = aa + adiag_i + 1;
688d708749eSBarry Smith     vi      = aj + adiag_i + 1;
689930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
690930ae53cSBarry Smith     sum     = x[i];
691930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
692930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
693930ae53cSBarry Smith   }
6941c4feecaSSatish Balay #endif
695b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
696cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
697cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
699930ae53cSBarry Smith }
700930ae53cSBarry Smith 
7014a2ae208SSatish Balay #undef __FUNCT__
7024a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
703416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
704da3a660dSBarry Smith {
705416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
706416022c9SBarry Smith   IS           iscol = a->col,isrow = a->row;
707273d9f13SBarry Smith   int          *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
708010a20caSHong Zhang   int          nz,*rout,*cout;
70987828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,sum,*v;
710da3a660dSBarry Smith 
7113a40ed3dSBarry Smith   PetscFunctionBegin;
71278b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
713da3a660dSBarry Smith 
714e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
715e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
716416022c9SBarry Smith   tmp  = a->solve_work;
717da3a660dSBarry Smith 
7183b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7193b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
720da3a660dSBarry Smith 
721da3a660dSBarry Smith   /* forward solve the lower triangular */
722da3a660dSBarry Smith   tmp[0] = b[*r++];
723da3a660dSBarry Smith   for (i=1; i<n; i++) {
724010a20caSHong Zhang     v   = aa + ai[i] ;
725010a20caSHong Zhang     vi  = aj + ai[i] ;
726416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
727da3a660dSBarry Smith     sum = b[*r++];
728010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
729da3a660dSBarry Smith     tmp[i] = sum;
730da3a660dSBarry Smith   }
731da3a660dSBarry Smith 
732da3a660dSBarry Smith   /* backward solve the upper triangular */
733da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
734010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
735010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
736416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
737da3a660dSBarry Smith     sum = tmp[i];
738010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
739010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
740da3a660dSBarry Smith     x[*c--] += tmp[i];
741da3a660dSBarry Smith   }
742da3a660dSBarry Smith 
7433b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7443b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
745cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
746cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
747b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
748898183ecSLois Curfman McInnes 
7493a40ed3dSBarry Smith   PetscFunctionReturn(0);
750da3a660dSBarry Smith }
751da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7524a2ae208SSatish Balay #undef __FUNCT__
7534a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
7547c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
755da3a660dSBarry Smith {
756416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
7572235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
758273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
759010a20caSHong Zhang   int          nz,*rout,*cout,*diag = a->diag;
76087828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v,s1;
761da3a660dSBarry Smith 
7623a40ed3dSBarry Smith   PetscFunctionBegin;
763e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
764e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
765416022c9SBarry Smith   tmp  = a->solve_work;
766da3a660dSBarry Smith 
7672235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7682235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
769da3a660dSBarry Smith 
770da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7712235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
772da3a660dSBarry Smith 
773da3a660dSBarry Smith   /* forward solve the U^T */
774da3a660dSBarry Smith   for (i=0; i<n; i++) {
775010a20caSHong Zhang     v   = aa + diag[i] ;
776010a20caSHong Zhang     vi  = aj + diag[i] + 1;
777f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
778c38d4ed2SBarry Smith     s1  = tmp[i];
779ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
780da3a660dSBarry Smith     while (nz--) {
781010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
782da3a660dSBarry Smith     }
783c38d4ed2SBarry Smith     tmp[i] = s1;
784da3a660dSBarry Smith   }
785da3a660dSBarry Smith 
786da3a660dSBarry Smith   /* backward solve the L^T */
787da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
788010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
789010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
790f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
791f1af5d2fSBarry Smith     s1  = tmp[i];
792da3a660dSBarry Smith     while (nz--) {
793010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
794da3a660dSBarry Smith     }
795da3a660dSBarry Smith   }
796da3a660dSBarry Smith 
797da3a660dSBarry Smith   /* copy tmp into x according to permutation */
798da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
799da3a660dSBarry Smith 
8002235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8012235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
802cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
803cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
804da3a660dSBarry Smith 
805b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8063a40ed3dSBarry Smith   PetscFunctionReturn(0);
807da3a660dSBarry Smith }
808da3a660dSBarry Smith 
8094a2ae208SSatish Balay #undef __FUNCT__
8104a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
8117c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
812da3a660dSBarry Smith {
813416022c9SBarry Smith   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
8142235e667SBarry Smith   IS           iscol = a->col,isrow = a->row;
815273d9f13SBarry Smith   int          *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
816010a20caSHong Zhang   int          nz,*rout,*cout,*diag = a->diag;
81787828ca2SBarry Smith   PetscScalar  *x,*b,*tmp,*aa = a->a,*v;
8186abc6512SBarry Smith 
8193a40ed3dSBarry Smith   PetscFunctionBegin;
8206abc6512SBarry Smith   if (zz != xx) VecCopy(zz,xx);
8216abc6512SBarry Smith 
822e1311b90SBarry Smith   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
823e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
824416022c9SBarry Smith   tmp = a->solve_work;
8256abc6512SBarry Smith 
8262235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8272235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8286abc6512SBarry Smith 
8296abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8302235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8316abc6512SBarry Smith 
8326abc6512SBarry Smith   /* forward solve the U^T */
8336abc6512SBarry Smith   for (i=0; i<n; i++) {
834010a20caSHong Zhang     v   = aa + diag[i] ;
835010a20caSHong Zhang     vi  = aj + diag[i] + 1;
836f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8376abc6512SBarry Smith     tmp[i] *= *v++;
8386abc6512SBarry Smith     while (nz--) {
839010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8406abc6512SBarry Smith     }
8416abc6512SBarry Smith   }
8426abc6512SBarry Smith 
8436abc6512SBarry Smith   /* backward solve the L^T */
8446abc6512SBarry Smith   for (i=n-1; i>=0; i--){
845010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
846010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
847f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8486abc6512SBarry Smith     while (nz--) {
849010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8506abc6512SBarry Smith     }
8516abc6512SBarry Smith   }
8526abc6512SBarry Smith 
8536abc6512SBarry Smith   /* copy tmp into x according to permutation */
8546abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8556abc6512SBarry Smith 
8562235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8572235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
858cb5b572fSBarry Smith   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
859cb5b572fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8606abc6512SBarry Smith 
861b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8623a40ed3dSBarry Smith   PetscFunctionReturn(0);
863da3a660dSBarry Smith }
8649e25ed09SBarry Smith /* ----------------------------------------------------------------*/
865ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat);
86608480c60SBarry Smith 
8674a2ae208SSatish Balay #undef __FUNCT__
8684a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
869b380c88cSHong Zhang int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8709e25ed09SBarry Smith {
871416022c9SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
8729e25ed09SBarry Smith   IS         isicol;
873273d9f13SBarry Smith   int        *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j;
87456beaf04SBarry Smith   int        *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
875335d9088SBarry Smith   int        *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0;
876010a20caSHong Zhang   int        incrlev,nnz,i,levels,diagonal_fill;
87777c4ece6SBarry Smith   PetscTruth col_identity,row_identity;
878329f5518SBarry Smith   PetscReal  f;
8799e25ed09SBarry Smith 
8803a40ed3dSBarry Smith   PetscFunctionBegin;
8816cf997b0SBarry Smith   f             = info->fill;
8820cd17293SBarry Smith   levels        = (int)info->levels;
8830cd17293SBarry Smith   diagonal_fill = (int)info->diagonal_fill;
8844c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
88582bf6240SBarry Smith 
88608480c60SBarry Smith   /* special case that simply copies fill pattern */
887be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
888be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
88986bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8902e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
89108480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
89208480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
89308480c60SBarry Smith     if (!b->diag) {
8947c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89508480c60SBarry Smith     }
8967c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89708480c60SBarry Smith     b->row              = isrow;
89808480c60SBarry Smith     b->col              = iscol;
89982bf6240SBarry Smith     b->icol             = isicol;
900a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
90187828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
9022cea7109SBarry Smith     b->lu_shift         = info->shift;
9032cea7109SBarry Smith     b->lu_shift_fraction= info->shift_fraction;
90487828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
905f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
906c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
907c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9083a40ed3dSBarry Smith     PetscFunctionReturn(0);
90908480c60SBarry Smith   }
91008480c60SBarry Smith 
911898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
912898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9139e25ed09SBarry Smith 
9149e25ed09SBarry Smith   /* get new row pointers */
915b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr);
916010a20caSHong Zhang   ainew[0] = 0;
9179e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
918010a20caSHong Zhang   jmax = (int)(f*(ai[n]+1));
919b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr);
9209e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
921b0a32e0cSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr);
9229e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
923b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr);
92456beaf04SBarry Smith   /* im is level for each filled value */
925b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr);
92656beaf04SBarry Smith   /* dloc is location of diagonal in factor */
927b0a32e0cSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr);
92856beaf04SBarry Smith   dloc[0]  = 0;
92956beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9306cf997b0SBarry Smith 
9316cf997b0SBarry Smith     /* copy current row into linked list */
93256beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
93329bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
934010a20caSHong Zhang     xi      = aj + ai[r[prow]] ;
9359e25ed09SBarry Smith     fill[n]    = n;
936435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9379e25ed09SBarry Smith     while (nz--) {
9389e25ed09SBarry Smith       fm  = n;
939010a20caSHong Zhang       idx = ic[*xi++ ];
9409e25ed09SBarry Smith       do {
9419e25ed09SBarry Smith         m  = fm;
9429e25ed09SBarry Smith         fm = fill[m];
9439e25ed09SBarry Smith       } while (fm < idx);
9449e25ed09SBarry Smith       fill[m]   = idx;
9459e25ed09SBarry Smith       fill[idx] = fm;
94656beaf04SBarry Smith       im[idx]   = 0;
9479e25ed09SBarry Smith     }
9486cf997b0SBarry Smith 
9496cf997b0SBarry Smith     /* make sure diagonal entry is included */
950435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9516cf997b0SBarry Smith       fm = n;
952435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
953435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
954435faa5fSBarry Smith       fill[fm]   = prow;
9556cf997b0SBarry Smith       im[prow]   = 0;
9566cf997b0SBarry Smith       nzf++;
957335d9088SBarry Smith       dcount++;
9586cf997b0SBarry Smith     }
9596cf997b0SBarry Smith 
96056beaf04SBarry Smith     nzi = 0;
9619e25ed09SBarry Smith     row = fill[n];
96256beaf04SBarry Smith     while (row < prow) {
96356beaf04SBarry Smith       incrlev = im[row] + 1;
96456beaf04SBarry Smith       nz      = dloc[row];
965010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
966010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
967416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
96856beaf04SBarry Smith       fm      = row;
96956beaf04SBarry Smith       while (nnz-- > 0) {
970010a20caSHong Zhang         idx = *xi++ ;
97156beaf04SBarry Smith         if (*flev + incrlev > levels) {
97256beaf04SBarry Smith           flev++;
97356beaf04SBarry Smith           continue;
97456beaf04SBarry Smith         }
97556beaf04SBarry Smith         do {
9769e25ed09SBarry Smith           m  = fm;
9779e25ed09SBarry Smith           fm = fill[m];
97856beaf04SBarry Smith         } while (fm < idx);
9799e25ed09SBarry Smith         if (fm != idx) {
98056beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9819e25ed09SBarry Smith           fill[m]   = idx;
9829e25ed09SBarry Smith           fill[idx] = fm;
9839e25ed09SBarry Smith           fm        = idx;
98456beaf04SBarry Smith           nzf++;
985ecf371e4SBarry Smith         } else {
98656beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9879e25ed09SBarry Smith         }
98856beaf04SBarry Smith         flev++;
9899e25ed09SBarry Smith       }
9909e25ed09SBarry Smith       row = fill[row];
99156beaf04SBarry Smith       nzi++;
9929e25ed09SBarry Smith     }
9939e25ed09SBarry Smith     /* copy new filled row into permanent storage */
99456beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
995010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
996ecf371e4SBarry Smith 
997ecf371e4SBarry Smith       /* estimate how much additional space we will need */
998ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
999ecf371e4SBarry Smith       /* just double the memory each time */
1000ecf371e4SBarry Smith       /*  maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */
1001ecf371e4SBarry Smith       int maxadd = jmax;
100256beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1003bbb6d6a8SBarry Smith       jmax += maxadd;
1004ecf371e4SBarry Smith 
1005ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
1006b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
1007010a20caSHong Zhang       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(int));CHKERRQ(ierr);
1008606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
100956beaf04SBarry Smith       ajnew  = xi;
1010b0a32e0cSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr);
1011010a20caSHong Zhang       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(int));CHKERRQ(ierr);
1012606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
101356beaf04SBarry Smith       ajfill = xi;
1014ecf371e4SBarry Smith       realloc++; /* count how many times we realloc */
10159e25ed09SBarry Smith     }
1016010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1017010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
101856beaf04SBarry Smith     dloc[prow]  = nzi;
10199e25ed09SBarry Smith     fm          = fill[n];
102056beaf04SBarry Smith     while (nzf--) {
1021010a20caSHong Zhang       *xi++   = fm ;
102256beaf04SBarry Smith       *flev++ = im[fm];
10239e25ed09SBarry Smith       fm      = fill[fm];
10249e25ed09SBarry Smith     }
10256cf997b0SBarry Smith     /* make sure row has diagonal entry */
1026010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
102729bbc08cSBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\
10286cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10296cf997b0SBarry Smith     }
10309e25ed09SBarry Smith   }
1031606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1032898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1033898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1034606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1035606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10369e25ed09SBarry Smith 
1037f64065bbSBarry Smith   {
1038329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1039b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af);
1040c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1041c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1042b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1043335d9088SBarry Smith     if (diagonal_fill) {
1044b1bcba4aSBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount);
1045335d9088SBarry Smith     }
1046f64065bbSBarry Smith   }
1047bbb6d6a8SBarry Smith 
10489e25ed09SBarry Smith   /* put together the new matrix */
1049b4fd4287SBarry Smith   ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr);
1050b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1051416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1052606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10537c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1054010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10559e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1056606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1057606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1058b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1059416022c9SBarry Smith   b->j          = ajnew;
1060416022c9SBarry Smith   b->i          = ainew;
106156beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1062416022c9SBarry Smith   b->diag       = dloc;
1063416022c9SBarry Smith   b->ilen       = 0;
1064416022c9SBarry Smith   b->imax       = 0;
1065416022c9SBarry Smith   b->row        = isrow;
1066416022c9SBarry Smith   b->col        = iscol;
1067c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1068c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
106982bf6240SBarry Smith   b->icol       = isicol;
107087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1071416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
107256beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
1073010a20caSHong Zhang   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(int)+sizeof(PetscScalar)));
1074010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
1075a2e34c3dSBarry Smith   b->lu_damping        = info->damping;
10762cea7109SBarry Smith   b->lu_shift          = info->shift;
10772cea7109SBarry Smith   b->lu_shift_fraction = info->shift_fraction;
107887828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10799e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10803a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10813a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1082ae068f58SLois Curfman McInnes 
1083ae068f58SLois Curfman McInnes   (*fact)->info.factor_mallocs    = realloc;
1084ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1085329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10863a40ed3dSBarry Smith   PetscFunctionReturn(0);
10879e25ed09SBarry Smith }
1088d5d45c9bSBarry Smith 
1089a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1090a6175056SHong Zhang #undef __FUNCT__
1091a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1092a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact)
1093a6175056SHong Zhang {
10940968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1095a6175056SHong Zhang   int                 ierr;
1096d5d45c9bSBarry Smith 
1097a6175056SHong Zhang   PetscFunctionBegin;
10980968510aSHong Zhang   if (!a->sbaijMat){
10990968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1100a6175056SHong Zhang   }
110103aac1b8SHong Zhang 
1102b45a75daSHong Zhang   ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr);
11030968510aSHong Zhang   ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
1104064503c5SHong Zhang   a->sbaijMat = PETSC_NULL;
1105653a6975SHong Zhang 
1106a6175056SHong Zhang   PetscFunctionReturn(0);
1107a6175056SHong Zhang }
1108a6175056SHong Zhang 
1109a6175056SHong Zhang #undef __FUNCT__
1110a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
111115e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1112a6175056SHong Zhang {
11130968510aSHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
111403aac1b8SHong Zhang   int                 ierr;
1115653a6975SHong Zhang   PetscTruth          perm_identity;
1116a6175056SHong Zhang 
1117a6175056SHong Zhang   PetscFunctionBegin;
1118653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1119653a6975SHong Zhang   if (!perm_identity){
1120653a6975SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1121653a6975SHong Zhang   }
11220968510aSHong Zhang   if (!a->sbaijMat){
11230968510aSHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
11240968510aSHong Zhang   }
1125a6175056SHong Zhang 
1126b00f7748SHong Zhang   ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1127653a6975SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1128653a6975SHong Zhang 
1129a6175056SHong Zhang   PetscFunctionReturn(0);
1130a6175056SHong Zhang }
1131d5d45c9bSBarry Smith 
1132f76d2b81SHong Zhang #undef __FUNCT__
1133f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1134f76d2b81SHong Zhang int MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1135f76d2b81SHong Zhang {
1136f76d2b81SHong Zhang   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
113703aac1b8SHong Zhang   int                 ierr;
1138f76d2b81SHong Zhang   PetscTruth          perm_identity;
1139f76d2b81SHong Zhang 
1140f76d2b81SHong Zhang   PetscFunctionBegin;
1141f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1142f76d2b81SHong Zhang   if (!perm_identity){
1143f76d2b81SHong Zhang     SETERRQ(1,"Non-identity permutation is not supported yet");
1144f76d2b81SHong Zhang   }
1145f76d2b81SHong Zhang   if (!a->sbaijMat){
1146f76d2b81SHong Zhang     ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr);
1147f76d2b81SHong Zhang   }
1148f76d2b81SHong Zhang 
1149f76d2b81SHong Zhang   ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
1150f76d2b81SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1151f76d2b81SHong Zhang 
1152f76d2b81SHong Zhang   PetscFunctionReturn(0);
1153f76d2b81SHong Zhang }
1154