xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision b3bf805bc5fa6ada2d57fa222e3d33ae005e5770)
1289bc588SBarry Smith 
270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
31c4feecaSSatish Balay #include "src/inline/dot.h"
4ed480e8bSBarry Smith #include "src/inline/spops.h"
53b2fbd54SBarry Smith 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
8dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
93b2fbd54SBarry Smith {
103a40ed3dSBarry Smith   PetscFunctionBegin;
113a40ed3dSBarry Smith 
1229bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
13aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
14d64ed03dSBarry Smith   PetscFunctionReturn(0);
15d64ed03dSBarry Smith #endif
163b2fbd54SBarry Smith }
173b2fbd54SBarry Smith 
1886bacbd2SBarry Smith 
19dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
20dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
2186bacbd2SBarry Smith 
2238baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
239cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2438baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
2586bacbd2SBarry Smith 
264a2ae208SSatish Balay #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2886bacbd2SBarry Smith   /* ------------------------------------------------------------
2986bacbd2SBarry Smith 
3086bacbd2SBarry Smith           This interface was contribed by Tony Caola
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
335bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
345bd2ddc7SBarry Smith      SPARSEKIT2.
355bd2ddc7SBarry Smith 
365bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
375bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
3886bacbd2SBarry Smith 
3986bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4086bacbd2SBarry Smith      help in getting this routine ironed out.
4186bacbd2SBarry Smith 
425bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
435bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
445bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
455bd2ddc7SBarry Smith      Yousef Saad's code.
4686bacbd2SBarry Smith 
4786bacbd2SBarry Smith      ------------------------------------------------------------
4886bacbd2SBarry Smith */
49dfbe8321SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5086bacbd2SBarry Smith {
5160ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5260ec86bdSBarry Smith   PetscFunctionBegin;
53e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5460ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5560ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5660ec86bdSBarry Smith #else
5786bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
5807d2056aSBarry Smith   IS             iscolf,isicol,isirow;
5986bacbd2SBarry Smith   PetscTruth     reorder;
609cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
619cc1f4e3SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->m;
6238baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6338baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6438baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6587828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6654a8161fSBarry Smith   PetscReal      af;
6786bacbd2SBarry Smith 
6886bacbd2SBarry Smith   PetscFunctionBegin;
6986bacbd2SBarry Smith 
7086bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7138baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7286bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7333259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7438baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7538baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7686bacbd2SBarry Smith 
7786bacbd2SBarry Smith 
7886bacbd2SBarry Smith   /* ------------------------------------------------------------
7986bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8086bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8186bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8286bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8386bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8486bacbd2SBarry Smith      ------------------------------------------------------------
8586bacbd2SBarry Smith   */
8686bacbd2SBarry Smith 
8786bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
8886bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
8986bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9086bacbd2SBarry Smith   reorder = PetscNot(reorder);
9186bacbd2SBarry Smith 
9286bacbd2SBarry Smith 
9386bacbd2SBarry Smith   /* storage for ilu factor */
9438baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9538baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
9687828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
9738baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
9886bacbd2SBarry Smith 
9986bacbd2SBarry Smith   /* ------------------------------------------------------------
10086bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10186bacbd2SBarry Smith      ------------------------------------------------------------
10286bacbd2SBarry Smith   */
103b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
104b19713cbSBarry Smith     old_j[i]++;
10586bacbd2SBarry Smith   }
106b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
107b19713cbSBarry Smith     old_i[i]++;
108b19713cbSBarry Smith   };
109010a20caSHong Zhang 
11086bacbd2SBarry Smith 
11186bacbd2SBarry Smith   if (reorder) {
112c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
113c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
114c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
115c0c2c81eSBarry Smith       r[i]  = r[i]+1;
116c0c2c81eSBarry Smith       c[i]  = c[i]+1;
117c0c2c81eSBarry Smith     }
11838baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
11938baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12087828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1215bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
122c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
123c0c2c81eSBarry Smith       r[i]  = r[i]-1;
124c0c2c81eSBarry Smith       c[i]  = c[i]-1;
125c0c2c81eSBarry Smith     }
126c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
127c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
128b19713cbSBarry Smith     o_a = old_a2;
129b19713cbSBarry Smith     o_j = old_j2;
130b19713cbSBarry Smith     o_i = old_i2;
131b19713cbSBarry Smith   } else {
132b19713cbSBarry Smith     o_a = old_a;
133b19713cbSBarry Smith     o_j = old_j;
134b19713cbSBarry Smith     o_i = old_i;
13586bacbd2SBarry Smith   }
13686bacbd2SBarry Smith 
13786bacbd2SBarry Smith   /* ------------------------------------------------------------
13886bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
13986bacbd2SBarry Smith      ------------------------------------------------------------
14086bacbd2SBarry Smith   */
14186bacbd2SBarry Smith 
14238baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14338baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14487828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14586bacbd2SBarry Smith 
14654a8161fSBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1476849ba73SBarry Smith   if (sierr) {
1486849ba73SBarry Smith     switch (sierr) {
14977431f27SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
15077431f27SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
151e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
152e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15377431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15477431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15586bacbd2SBarry Smith     }
15686bacbd2SBarry Smith   }
15786bacbd2SBarry Smith 
15886bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
15986bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16086bacbd2SBarry Smith 
16186bacbd2SBarry Smith   /* ------------------------------------------------------------
16286bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16386bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16486bacbd2SBarry Smith      ------------------------------------------------------------
16586bacbd2SBarry Smith   */
16686bacbd2SBarry Smith 
16787828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
16986bacbd2SBarry Smith 
1705bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17186bacbd2SBarry Smith 
17286bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17386bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17486bacbd2SBarry Smith 
17586bacbd2SBarry Smith   if (reorder) {
17686bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17786bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
179313ae024SBarry Smith   } else {
180b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
181f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
182b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
183b19713cbSBarry Smith     }
184b19713cbSBarry Smith   }
185b19713cbSBarry Smith 
186b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
188b19713cbSBarry Smith     old_i[i]--;
189b19713cbSBarry Smith   }
190b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
191b19713cbSBarry Smith     old_j[i]--;
192b19713cbSBarry Smith   }
19386bacbd2SBarry Smith 
194b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19586bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
196b19713cbSBarry Smith     new_i[i]--;
197b19713cbSBarry Smith   }
198b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
199b19713cbSBarry Smith     new_j[i]--;
200b19713cbSBarry Smith   }
20186bacbd2SBarry Smith 
20286bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20386bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2044c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2054c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
206c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20786bacbd2SBarry Smith   for(i=0; i<n; i++) {
20886bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
20986bacbd2SBarry Smith   };
21086bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21107d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21286bacbd2SBarry Smith 
213c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
214c0c2c81eSBarry Smith 
215329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21607d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21786bacbd2SBarry Smith 
21886bacbd2SBarry Smith   /*----- put together the new matrix -----*/
21986bacbd2SBarry Smith 
220f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
221f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
222f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);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"
266dfbe8321SBarry Smith PetscErrorCode 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;
2706849ba73SBarry Smith   PetscErrorCode ierr;
27138baddfdSBarry Smith   PetscInt       *r,*ic,i,n = A->m,*ai = a->i,*aj = a->j;
27238baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*ajtmp,nz;
273418422e8SSatish Balay   PetscInt       *idnew,idx,row,m,fm,nnz,nzi,reallocs = 0,nzbd,*im;
2749e046878SBarry Smith   PetscReal      f;
275289bc588SBarry Smith 
2763a40ed3dSBarry Smith   PetscFunctionBegin;
27729bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2784c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
279ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
280ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
281289bc588SBarry Smith 
282289bc588SBarry Smith   /* get new row pointers */
28338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
284010a20caSHong Zhang   ainew[0] = 0;
285289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
286d3d32019SBarry Smith   f = info->fill;
28738baddfdSBarry Smith   jmax  = (PetscInt)(f*ai[n]+1);
28838baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
289289bc588SBarry Smith   /* fill is a linked list of nonzeros in active row */
29038baddfdSBarry Smith   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
2912fb3ed76SBarry Smith   im   = fill + n + 1;
292289bc588SBarry Smith   /* idnew is location of diagonal in factor */
29338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&idnew);CHKERRQ(ierr);
294010a20caSHong Zhang   idnew[0] = 0;
295289bc588SBarry Smith 
296289bc588SBarry Smith   for (i=0; i<n; i++) {
297289bc588SBarry Smith     /* first copy previous fill into linked list */
298289bc588SBarry Smith     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
29929bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
300010a20caSHong Zhang     ajtmp   = aj + ai[r[i]];
301289bc588SBarry Smith     fill[n] = n;
302289bc588SBarry Smith     while (nz--) {
303289bc588SBarry Smith       fm  = n;
304010a20caSHong Zhang       idx = ic[*ajtmp++];
305289bc588SBarry Smith       do {
306289bc588SBarry Smith         m  = fm;
307289bc588SBarry Smith         fm = fill[m];
308289bc588SBarry Smith       } while (fm < idx);
309289bc588SBarry Smith       fill[m]   = idx;
310289bc588SBarry Smith       fill[idx] = fm;
311289bc588SBarry Smith     }
312289bc588SBarry Smith     row = fill[n];
313289bc588SBarry Smith     while (row < i) {
314010a20caSHong Zhang       ajtmp = ajnew + idnew[row] + 1;
3152fb3ed76SBarry Smith       nzbd  = 1 + idnew[row] - ainew[row];
3162fb3ed76SBarry Smith       nz    = im[row] - nzbd;
317289bc588SBarry Smith       fm    = row;
3182fb3ed76SBarry Smith       while (nz-- > 0) {
319010a20caSHong Zhang         idx = *ajtmp++ ;
3202fb3ed76SBarry Smith         nzbd++;
3212fb3ed76SBarry Smith         if (idx == i) im[row] = nzbd;
322289bc588SBarry Smith         do {
323289bc588SBarry Smith           m  = fm;
324289bc588SBarry Smith           fm = fill[m];
325289bc588SBarry Smith         } while (fm < idx);
326289bc588SBarry Smith         if (fm != idx) {
327289bc588SBarry Smith           fill[m]   = idx;
328289bc588SBarry Smith           fill[idx] = fm;
329289bc588SBarry Smith           fm        = idx;
330289bc588SBarry Smith           nnz++;
331289bc588SBarry Smith         }
332289bc588SBarry Smith       }
333289bc588SBarry Smith       row = fill[row];
334289bc588SBarry Smith     }
335289bc588SBarry Smith     /* copy new filled row into permanent storage */
336289bc588SBarry Smith     ainew[i+1] = ainew[i] + nnz;
337496e697eSBarry Smith     if (ainew[i+1] > jmax) {
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 */
34138baddfdSBarry Smith       PetscInt 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 */
34738baddfdSBarry Smith       ierr  = PetscMalloc(jmax*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
34838baddfdSBarry Smith       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i])*sizeof(PetscInt));CHKERRQ(ierr);
349606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
350289bc588SBarry Smith       ajnew = ajtmp;
351418422e8SSatish Balay       reallocs++; /* 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]);
366418422e8SSatish Balay     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,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 */
380f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr);
381f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
382f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr);
383b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
384416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
385606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3867c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
387e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
388606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
389606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
390010a20caSHong Zhang   ierr          = PetscMalloc((ainew[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
391416022c9SBarry Smith   b->j          = ajnew;
392416022c9SBarry Smith   b->i          = ainew;
393416022c9SBarry Smith   b->diag       = idnew;
394416022c9SBarry Smith   b->ilen       = 0;
395416022c9SBarry Smith   b->imax       = 0;
396416022c9SBarry Smith   b->row        = isrow;
397416022c9SBarry Smith   b->col        = iscol;
398c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
399c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40082bf6240SBarry Smith   b->icol       = isicol;
40187828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
402416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4037fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
40438baddfdSBarry Smith   PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));
405010a20caSHong Zhang   b->maxnz = b->nz = ainew[n] ;
4067fd98bd6SLois Curfman McInnes 
407329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
408418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
409ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4103a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4117dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
412703deb49SSatish Balay 
413e93ccc4dSBarry Smith   if (ai[n] != 0) {
414329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
415eea2913aSSatish Balay   } else {
416eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
417eea2913aSSatish Balay   }
4183a40ed3dSBarry Smith   PetscFunctionReturn(0);
419289bc588SBarry Smith }
4200a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
421dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
42241c01911SSatish Balay 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
425af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
426289bc588SBarry Smith {
427416022c9SBarry Smith   Mat            C=*B;
428416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
42982bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4306849ba73SBarry Smith   PetscErrorCode ierr;
431*b3bf805bSHong Zhang   PetscInt       *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j;
432*b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
433d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
43487828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
435ee45ca4aSHong Zhang   PetscReal      zeropivot,rs,d,shift_nonzero;
436d2276718SHong Zhang   PetscReal      row_shift,shift_top=0.;
437d2276718SHong Zhang   PetscTruth     shift_pd;
438*b3bf805bSHong Zhang   LUShift_Ctx    sctx;
439d2276718SHong Zhang   PetscInt       newshift;
440289bc588SBarry Smith 
4413a40ed3dSBarry Smith   PetscFunctionBegin;
4420a29876aSHong Zhang   shift_nonzero  = info->shiftnz;
4430a29876aSHong Zhang   shift_pd       = info->shiftpd;
4440a29876aSHong Zhang   zeropivot      = info->zeropivot;
4450a29876aSHong Zhang 
44678b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44778b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44887828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44987828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
450010a20caSHong Zhang   rtmps = rtmp; ics = ic;
451289bc588SBarry Smith 
4526cc28720Svictorle   if (!a->diag) {
4530cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4540cf777b8SBarry Smith   }
455e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
456e24cfe64SHong Zhang   if (shift_pd && shift_nonzero) shift_nonzero = 0.0;
457e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45838baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4590cf777b8SBarry Smith     shift_top = 0;
4606cc28720Svictorle     for (i=0; i<n; i++) {
461010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4626cc28720Svictorle       /* calculate amt of shift needed for this row */
4636c037d1bSvictorle       if (d<=0) {
4643aeef889SHong Zhang         row_shift = 0;
4653aeef889SHong Zhang       } else {
4663aeef889SHong Zhang         row_shift = -2*d;
4673aeef889SHong Zhang       }
468010a20caSHong Zhang       v  = a->a+aai[i];
469e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
470e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4716cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4726cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4736cc28720Svictorle     }
4746cc28720Svictorle   }
4756cc28720Svictorle 
476d2276718SHong Zhang   sctx.shift_amount = 0;
477d2276718SHong Zhang   sctx.shift_top    = shift_top;
478d2276718SHong Zhang   sctx.nshift       = 0;
479d2276718SHong Zhang   sctx.nshift_max   = 5;
480d2276718SHong Zhang   sctx.shift_lo     = 0.;
481d2276718SHong Zhang   sctx.shift_hi     = 1.;
482085256b3SBarry Smith   do {
483d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
484289bc588SBarry Smith     for (i=0; i<n; i++){
485*b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
486*b3bf805bSHong Zhang       bjtmp = bj + bi[i];
487*b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
488289bc588SBarry Smith 
489289bc588SBarry Smith       /* load in initial (unfactored row) */
490416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
491*b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
492010a20caSHong Zhang       v     = a->a + a->i[r[i]];
493085256b3SBarry Smith       for (j=0; j<nz; j++) {
494*b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
495085256b3SBarry Smith       }
496d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
497289bc588SBarry Smith 
498*b3bf805bSHong Zhang       row = *bjtmp++;
499f2582111SSatish Balay       while  (row < i) {
5008c37ef55SBarry Smith         pc = rtmp + row;
5018c37ef55SBarry Smith         if (*pc != 0.0) {
502010a20caSHong Zhang           pv         = b->a + diag_offset[row];
503010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50435aab85fSBarry Smith           multiplier = *pc / *pv++;
5058c37ef55SBarry Smith           *pc        = multiplier;
506*b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
507f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
508b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
509289bc588SBarry Smith         }
510*b3bf805bSHong Zhang         row = *bjtmp++;
511289bc588SBarry Smith       }
512416022c9SBarry Smith       /* finished row so stick it into b->a */
513*b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
514*b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
515*b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
516*b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
517d3d32019SBarry Smith       rs   = 0.0;
518d3d32019SBarry Smith       for (j=0; j<nz; j++) {
519d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
520d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
521d3d32019SBarry Smith       }
522d2276718SHong Zhang 
523*b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
524d2276718SHong Zhang       sctx.rs  = rs;
525d2276718SHong Zhang       sctx.pv  = pv[diag];
526d2276718SHong Zhang       newshift = 0;
527*b3bf805bSHong Zhang       ierr = Mat_LUFactorCheckShift(info,&sctx,&newshift);CHKERRQ(ierr);
528d2276718SHong Zhang       if (newshift == 1){
529*b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
530d2276718SHong Zhang       } else if (newshift == -1){
531d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
532d708749eSBarry Smith       }
53335aab85fSBarry Smith     }
534d2276718SHong Zhang 
535d2276718SHong Zhang     if (shift_pd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5366cc28720Svictorle       /*
5376c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5386cc28720Svictorle        * then try lower shift
5396cc28720Svictorle        */
540d2276718SHong Zhang       sctx.shift_hi       = info->shift_fraction;
541d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
542d2276718SHong Zhang       sctx.shift_amount   = info->shift_fraction * sctx.shift_top;
543d2276718SHong Zhang       sctx.lushift        = PETSC_TRUE;
544d2276718SHong Zhang       sctx.nshift++;
5456cc28720Svictorle     }
546d2276718SHong Zhang   } while (sctx.lushift);
5470f11f9aeSSatish Balay 
54835aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
54935aab85fSBarry Smith   for (i=0; i<n; i++) {
550010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55135aab85fSBarry Smith   }
55235aab85fSBarry Smith 
553606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55478b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55578b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
556416022c9SBarry Smith   C->factor = FACTOR_LU;
5577dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
558c456f294SBarry Smith   C->assembled = PETSC_TRUE;
559b0a32e0cSBarry Smith   PetscLogFlops(C->n);
560d2276718SHong Zhang   if (sctx.nshift){
5610697b6caSHong Zhang     if (shift_nonzero) {
562d2276718SHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
5630697b6caSHong Zhang     } else if (shift_pd) {
564d2276718SHong Zhang       b->lu_shift_fraction = info->shift_fraction;
565d2276718SHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %D\n",info->shift_fraction,shift_top,sctx.nshift);
5666cc28720Svictorle     }
5670697b6caSHong Zhang   }
5683a40ed3dSBarry Smith   PetscFunctionReturn(0);
569289bc588SBarry Smith }
5700889b5dcSvictorle 
5710889b5dcSvictorle #undef __FUNCT__
5720889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
573dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5740889b5dcSvictorle {
5750889b5dcSvictorle   PetscFunctionBegin;
5760889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5770889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5780889b5dcSvictorle   PetscFunctionReturn(0);
5790889b5dcSvictorle }
5800889b5dcSvictorle 
5810889b5dcSvictorle 
5820a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5834a2ae208SSatish Balay #undef __FUNCT__
5844a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
585dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
586da3a660dSBarry Smith {
587dfbe8321SBarry Smith   PetscErrorCode ierr;
588416022c9SBarry Smith   Mat            C;
589416022c9SBarry Smith 
5903a40ed3dSBarry Smith   PetscFunctionBegin;
5919e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
592af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
593273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
594440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5953a40ed3dSBarry Smith   PetscFunctionReturn(0);
596da3a660dSBarry Smith }
5970a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5984a2ae208SSatish Balay #undef __FUNCT__
5994a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
600dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6018c37ef55SBarry Smith {
602416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
603416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6046849ba73SBarry Smith   PetscErrorCode ierr;
60538baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
60638baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
60787828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6088c37ef55SBarry Smith 
6093a40ed3dSBarry Smith   PetscFunctionBegin;
6103a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
61191bf9adeSBarry Smith 
6121ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6131ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
614416022c9SBarry Smith   tmp  = a->solve_work;
6158c37ef55SBarry Smith 
6163b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6173b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6188c37ef55SBarry Smith 
6198c37ef55SBarry Smith   /* forward solve the lower triangular */
6208c37ef55SBarry Smith   tmp[0] = b[*r++];
621010a20caSHong Zhang   tmps   = tmp;
6228c37ef55SBarry Smith   for (i=1; i<n; i++) {
623010a20caSHong Zhang     v   = aa + ai[i] ;
624010a20caSHong Zhang     vi  = aj + ai[i] ;
62553e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
62653e7425aSBarry Smith     sum = b[*r++];
627ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6288c37ef55SBarry Smith     tmp[i] = sum;
6298c37ef55SBarry Smith   }
6308c37ef55SBarry Smith 
6318c37ef55SBarry Smith   /* backward solve the upper triangular */
6328c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
633010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
634010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
635416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6368c37ef55SBarry Smith     sum = tmp[i];
637ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
638010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6398c37ef55SBarry Smith   }
6408c37ef55SBarry Smith 
6413b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6423b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6431ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6441ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
645b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6463a40ed3dSBarry Smith   PetscFunctionReturn(0);
6478c37ef55SBarry Smith }
648026e39d0SSatish Balay 
649930ae53cSBarry Smith /* ----------------------------------------------------------- */
6504a2ae208SSatish Balay #undef __FUNCT__
6514a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
652dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
653930ae53cSBarry Smith {
654930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6556849ba73SBarry Smith   PetscErrorCode ierr;
65638baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
657362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
658aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
65938baddfdSBarry Smith   PetscInt       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 
6661ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6671ebc52fbSHong Zhang   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);
6961ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6971ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6983a40ed3dSBarry Smith   PetscFunctionReturn(0);
699930ae53cSBarry Smith }
700930ae53cSBarry Smith 
7014a2ae208SSatish Balay #undef __FUNCT__
7024a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
703dfbe8321SBarry Smith PetscErrorCode 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;
7076849ba73SBarry Smith   PetscErrorCode ierr;
70838baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
70938baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
71087828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
711da3a660dSBarry Smith 
7123a40ed3dSBarry Smith   PetscFunctionBegin;
71378b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
714da3a660dSBarry Smith 
7151ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7161ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
717416022c9SBarry Smith   tmp  = a->solve_work;
718da3a660dSBarry Smith 
7193b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7203b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
721da3a660dSBarry Smith 
722da3a660dSBarry Smith   /* forward solve the lower triangular */
723da3a660dSBarry Smith   tmp[0] = b[*r++];
724da3a660dSBarry Smith   for (i=1; i<n; i++) {
725010a20caSHong Zhang     v   = aa + ai[i] ;
726010a20caSHong Zhang     vi  = aj + ai[i] ;
727416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
728da3a660dSBarry Smith     sum = b[*r++];
729010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
730da3a660dSBarry Smith     tmp[i] = sum;
731da3a660dSBarry Smith   }
732da3a660dSBarry Smith 
733da3a660dSBarry Smith   /* backward solve the upper triangular */
734da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
735010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
736010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
737416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
738da3a660dSBarry Smith     sum = tmp[i];
739010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
740010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
741da3a660dSBarry Smith     x[*c--] += tmp[i];
742da3a660dSBarry Smith   }
743da3a660dSBarry Smith 
7443b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7453b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7461ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7471ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
748b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
749898183ecSLois Curfman McInnes 
7503a40ed3dSBarry Smith   PetscFunctionReturn(0);
751da3a660dSBarry Smith }
752da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7534a2ae208SSatish Balay #undef __FUNCT__
7544a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
755dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
756da3a660dSBarry Smith {
757416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7582235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7596849ba73SBarry Smith   PetscErrorCode ierr;
76038baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
76138baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
76287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
763da3a660dSBarry Smith 
7643a40ed3dSBarry Smith   PetscFunctionBegin;
7651ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
767416022c9SBarry Smith   tmp  = a->solve_work;
768da3a660dSBarry Smith 
7692235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7702235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
771da3a660dSBarry Smith 
772da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7732235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
774da3a660dSBarry Smith 
775da3a660dSBarry Smith   /* forward solve the U^T */
776da3a660dSBarry Smith   for (i=0; i<n; i++) {
777010a20caSHong Zhang     v   = aa + diag[i] ;
778010a20caSHong Zhang     vi  = aj + diag[i] + 1;
779f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
780c38d4ed2SBarry Smith     s1  = tmp[i];
781ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
782da3a660dSBarry Smith     while (nz--) {
783010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
784da3a660dSBarry Smith     }
785c38d4ed2SBarry Smith     tmp[i] = s1;
786da3a660dSBarry Smith   }
787da3a660dSBarry Smith 
788da3a660dSBarry Smith   /* backward solve the L^T */
789da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
790010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
791010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
792f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
793f1af5d2fSBarry Smith     s1  = tmp[i];
794da3a660dSBarry Smith     while (nz--) {
795010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
796da3a660dSBarry Smith     }
797da3a660dSBarry Smith   }
798da3a660dSBarry Smith 
799da3a660dSBarry Smith   /* copy tmp into x according to permutation */
800da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
801da3a660dSBarry Smith 
8022235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8032235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8041ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8051ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
806da3a660dSBarry Smith 
807b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8083a40ed3dSBarry Smith   PetscFunctionReturn(0);
809da3a660dSBarry Smith }
810da3a660dSBarry Smith 
8114a2ae208SSatish Balay #undef __FUNCT__
8124a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
813dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
814da3a660dSBarry Smith {
815416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8162235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8176849ba73SBarry Smith   PetscErrorCode ierr;
81838baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
81938baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
82087828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8216abc6512SBarry Smith 
8223a40ed3dSBarry Smith   PetscFunctionBegin;
82323bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8246abc6512SBarry Smith 
8251ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
827416022c9SBarry Smith   tmp = a->solve_work;
8286abc6512SBarry Smith 
8292235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8302235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8316abc6512SBarry Smith 
8326abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8332235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8346abc6512SBarry Smith 
8356abc6512SBarry Smith   /* forward solve the U^T */
8366abc6512SBarry Smith   for (i=0; i<n; i++) {
837010a20caSHong Zhang     v   = aa + diag[i] ;
838010a20caSHong Zhang     vi  = aj + diag[i] + 1;
839f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8406abc6512SBarry Smith     tmp[i] *= *v++;
8416abc6512SBarry Smith     while (nz--) {
842010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8436abc6512SBarry Smith     }
8446abc6512SBarry Smith   }
8456abc6512SBarry Smith 
8466abc6512SBarry Smith   /* backward solve the L^T */
8476abc6512SBarry Smith   for (i=n-1; i>=0; i--){
848010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
849010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
850f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8516abc6512SBarry Smith     while (nz--) {
852010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8536abc6512SBarry Smith     }
8546abc6512SBarry Smith   }
8556abc6512SBarry Smith 
8566abc6512SBarry Smith   /* copy tmp into x according to permutation */
8576abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8586abc6512SBarry Smith 
8592235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8602235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8611ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8636abc6512SBarry Smith 
864b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8653a40ed3dSBarry Smith   PetscFunctionReturn(0);
866da3a660dSBarry Smith }
8679e25ed09SBarry Smith /* ----------------------------------------------------------------*/
868dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
86908480c60SBarry Smith 
8704a2ae208SSatish Balay #undef __FUNCT__
8714a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
872dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8739e25ed09SBarry Smith {
874416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8759e25ed09SBarry Smith   IS             isicol;
8766849ba73SBarry Smith   PetscErrorCode ierr;
87738baddfdSBarry Smith   PetscInt       *r,*ic,prow,n = A->m,*ai = a->i,*aj = a->j;
87838baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
879418422e8SSatish Balay   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,len, reallocs = 0,dcount = 0;
88038baddfdSBarry Smith   PetscInt       incrlev,nnz,i,levels,diagonal_fill;
88177c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
882329f5518SBarry Smith   PetscReal      f;
8839e25ed09SBarry Smith 
8843a40ed3dSBarry Smith   PetscFunctionBegin;
8856cf997b0SBarry Smith   f             = info->fill;
88638baddfdSBarry Smith   levels        = (PetscInt)info->levels;
88738baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8884c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
88982bf6240SBarry Smith 
89008480c60SBarry Smith   /* special case that simply copies fill pattern */
891be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
892be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
89386bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8942e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
89508480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
89608480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
89708480c60SBarry Smith     if (!b->diag) {
8987c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89908480c60SBarry Smith     }
9007c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90108480c60SBarry Smith     b->row              = isrow;
90208480c60SBarry Smith     b->col              = iscol;
90382bf6240SBarry Smith     b->icol             = isicol;
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 */
91538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
916010a20caSHong Zhang   ainew[0] = 0;
9179e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
91838baddfdSBarry Smith   jmax = (PetscInt)(f*(ai[n]+1));
91938baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
9209e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
92138baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajfill);CHKERRQ(ierr);
9229e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
92338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
92456beaf04SBarry Smith   /* im is level for each filled value */
92538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
92656beaf04SBarry Smith   /* dloc is location of diagonal in factor */
92738baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&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 */
100038baddfdSBarry Smith       /*  maxadd = (PetscInt)((f*(ai[n]+!shift)*(n-prow+5))/n); */
100138baddfdSBarry Smith       PetscInt 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 */
100638baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
100738baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1008606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
100956beaf04SBarry Smith       ajnew  = xi;
101038baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
101138baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1012606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
101356beaf04SBarry Smith       ajfill = xi;
1014418422e8SSatish Balay       reallocs++; /* 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) {
102777431f27SBarry 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]);
1039418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,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) {
104477431f27SBarry 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 */
1049f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1050f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1051f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
1052b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1053416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1054606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10557c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1056010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10579e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1058606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1059606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1060b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1061416022c9SBarry Smith   b->j          = ajnew;
1062416022c9SBarry Smith   b->i          = ainew;
106356beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1064416022c9SBarry Smith   b->diag       = dloc;
1065416022c9SBarry Smith   b->ilen       = 0;
1066416022c9SBarry Smith   b->imax       = 0;
1067416022c9SBarry Smith   b->row        = isrow;
1068416022c9SBarry Smith   b->col        = iscol;
1069c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1070c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
107182bf6240SBarry Smith   b->icol       = isicol;
107287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1073416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
107456beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
107538baddfdSBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
1076010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
10779e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10783a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10793a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1080418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1081ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1082329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10833a40ed3dSBarry Smith   PetscFunctionReturn(0);
10849e25ed09SBarry Smith }
1085d5d45c9bSBarry Smith 
10863c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1087a6175056SHong Zhang #undef __FUNCT__
1088a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1089af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1090a6175056SHong Zhang {
10912ed133dbSHong Zhang   Mat            C = *B;
10920968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10932ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
10942ed133dbSHong Zhang   IS             ip=b->row;
10952ed133dbSHong Zhang   PetscErrorCode ierr;
10962ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
10972ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
10982ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1099622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1100ee45ca4aSHong Zhang   PetscReal      zeropivot,shift_amount,rs,shift_nonzero;
1101ee45ca4aSHong Zhang   PetscTruth     chshift,shift_pd;
110239e3d630SHong Zhang   PetscInt       nshift=0;
1103d5d45c9bSBarry Smith 
1104a6175056SHong Zhang   PetscFunctionBegin;
1105ee45ca4aSHong Zhang   shift_nonzero  = info->shiftnz;
1106ee45ca4aSHong Zhang   shift_pd       = info->shiftpd;
1107ee45ca4aSHong Zhang   zeropivot      = info->zeropivot;
1108ee45ca4aSHong Zhang 
11092ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11102ed133dbSHong Zhang 
11112ed133dbSHong Zhang   /* initialization */
11122ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11132ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11142ed133dbSHong Zhang   jl   = il + mbs;
11152ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11162ed133dbSHong Zhang 
11172ed133dbSHong Zhang   shift_amount = 0;
11182ed133dbSHong Zhang   do {
11192ed133dbSHong Zhang     chshift = PETSC_FALSE;
11202ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11212ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1122a6175056SHong Zhang     }
11232ed133dbSHong Zhang 
11242ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1125c05c3958SHong Zhang       bval = ba + bi[k];
11262ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11272ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11282ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11292ed133dbSHong Zhang         col = rip[aj[j]];
11302ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11312ed133dbSHong Zhang           rtmp[col] = aa[j];
1132c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11332ed133dbSHong Zhang         }
11342ed133dbSHong Zhang       }
113539e3d630SHong Zhang       /* shift the diagonal of the matrix */
113639e3d630SHong Zhang       if (nshift) rtmp[k] += shift_amount;
11372ed133dbSHong Zhang 
11382ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11392ed133dbSHong Zhang       dk = rtmp[k];
11402ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11412ed133dbSHong Zhang 
11422ed133dbSHong Zhang       while (i < k){
11432ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11442ed133dbSHong Zhang 
11452ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11462ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11472ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11482ed133dbSHong Zhang         dk += uikdi*ba[ili];
11492ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11502ed133dbSHong Zhang 
11512ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11522ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11532ed133dbSHong Zhang         if (jmin < jmax){
11542ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11552ed133dbSHong Zhang           /* update il and jl for row i */
11562ed133dbSHong Zhang           il[i] = jmin;
11572ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11582ed133dbSHong Zhang         }
11592ed133dbSHong Zhang         i = nexti;
11602ed133dbSHong Zhang       }
11612ed133dbSHong Zhang 
11620697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11630697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11640697b6caSHong Zhang       rs   = 0.0;
11653c9e8598SHong Zhang       jmin = bi[k]+1;
11663c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11673c9e8598SHong Zhang       if (nz){
11683c9e8598SHong Zhang         bcol = bj + jmin;
11693c9e8598SHong Zhang         while (nz--){
11703c9e8598SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol++]);
11713c9e8598SHong Zhang         }
11723c9e8598SHong Zhang       }
11730697b6caSHong Zhang       if (PetscAbsScalar(dk) <= zeropivot*rs && shift_nonzero){
11740697b6caSHong Zhang         /* force |diag(*B)| > zeropivot*rs */
11750697b6caSHong Zhang         if (!nshift){
11760697b6caSHong Zhang           shift_amount = shift_nonzero;
11770697b6caSHong Zhang         } else {
11780697b6caSHong Zhang           shift_amount *= 2.0;
11790697b6caSHong Zhang         }
11800697b6caSHong Zhang         chshift = PETSC_TRUE;
11810697b6caSHong Zhang         nshift++;
11820697b6caSHong Zhang         break;
11830697b6caSHong Zhang       } else if (PetscRealPart(dk) <= zeropivot*rs && shift_pd){
11840697b6caSHong Zhang         /* calculate a shift that would make this row diagonally dominant */
11850697b6caSHong Zhang 	shift_amount = PetscMax(rs+PetscAbs(PetscRealPart(dk)),1.1*shift_amount);
11863c9e8598SHong Zhang 	chshift      = PETSC_TRUE;
11873c9e8598SHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
11883c9e8598SHong Zhang 	   we increase the shift until it converges. There is no guarantee that
11893c9e8598SHong Zhang 	   this algorithm converges faster or slower, or is better or worse
11903c9e8598SHong Zhang 	   than the ILU algorithm. */
11913c9e8598SHong Zhang 	nshift++;
11923c9e8598SHong Zhang 	break;
11930697b6caSHong Zhang       } else if (PetscAbsScalar(dk) <= zeropivot*rs){
11940697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11953c9e8598SHong Zhang       }
11963c9e8598SHong Zhang 
11973c9e8598SHong Zhang       /* copy data into U(k,:) */
119839e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
119939e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
120039e3d630SHong Zhang       if (jmin < jmax) {
120139e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
120239e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12033c9e8598SHong Zhang         }
120439e3d630SHong Zhang         /* add the k-th row into il and jl */
12053c9e8598SHong Zhang         il[k] = jmin;
12063c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12073c9e8598SHong Zhang       }
12083c9e8598SHong Zhang     }
12090697b6caSHong Zhang   } while (chshift);
12103c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12113c9e8598SHong Zhang 
121239e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12133c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12143c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12153c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
12163c9e8598SHong Zhang   PetscLogFlops(C->m);
12173c9e8598SHong Zhang   if (nshift){
12180697b6caSHong Zhang     if (shift_nonzero) {
121939e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",nshift,shift_amount);
12200697b6caSHong Zhang     } else if (shift_pd) {
122139e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_pd tries %D, shift_amount %g\n",nshift,shift_amount);
12220697b6caSHong Zhang     }
12233c9e8598SHong Zhang   }
12243c9e8598SHong Zhang   PetscFunctionReturn(0);
12253c9e8598SHong Zhang }
1226a6175056SHong Zhang 
12277a48dd6fSHong Zhang #include "petscbt.h"
12287a48dd6fSHong Zhang #include "src/mat/utils/freespace.h"
1229a6175056SHong Zhang #undef __FUNCT__
1230a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1231dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1232a6175056SHong Zhang {
12330968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1234ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1235ed59401aSHong Zhang   Mat            B;
1236dfbe8321SBarry Smith   PetscErrorCode ierr;
1237653a6975SHong Zhang   PetscTruth     perm_identity;
1238622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1239ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1240391eac55SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1241391eac55SHong Zhang   PetscInt       ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
1242ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
12437a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
12447a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12457a48dd6fSHong Zhang   PetscBT        lnkbt;
1246a6175056SHong Zhang 
1247a6175056SHong Zhang   PetscFunctionBegin;
1248653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12497a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12507a48dd6fSHong Zhang 
125139e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
125239e3d630SHong Zhang   ui[0] = 0;
125339e3d630SHong Zhang 
1254622e2028SHong Zhang   /* special case that simply copies fill pattern */
1255622e2028SHong Zhang   if (!levels && perm_identity) {
1256622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1257ed59401aSHong Zhang     for (i=0; i<am; i++) {
125839e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1259ed59401aSHong Zhang     }
126039e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
126139e3d630SHong Zhang     cols = uj;
1262ed59401aSHong Zhang     for (i=0; i<am; i++) {
1263ed59401aSHong Zhang       aj    = a->j + a->diag[i];
126439e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
126539e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1266ed59401aSHong Zhang     }
126739e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12687a48dd6fSHong Zhang     /* initialization */
1269622e2028SHong Zhang     ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
12707a48dd6fSHong Zhang 
12717a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12727a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12737a48dd6fSHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
12747a48dd6fSHong Zhang     il         = jl + am;
12757a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12767a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12777a48dd6fSHong Zhang     for (i=0; i<am; i++){
12787a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12797a48dd6fSHong Zhang     }
12807a48dd6fSHong Zhang 
12817a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12827a48dd6fSHong Zhang     nlnk = am + 1;
12832ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12847a48dd6fSHong Zhang 
12857a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12867a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12877a48dd6fSHong Zhang     current_space = free_space;
12887a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12897a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12907a48dd6fSHong Zhang 
12917a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12927a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12937a48dd6fSHong Zhang       nzk   = 0;
1294622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1295622e2028SHong Zhang       ncols_upper = 0;
1296622e2028SHong Zhang       cols        = cols_lvl + am;
1297622e2028SHong Zhang       for (j=0; j<ncols; j++){
1298622e2028SHong Zhang         i = rip[*(aj + ai[rip[k]] + j)];
1299622e2028SHong Zhang         if (i >= k){ /* only take upper triangular entry */
1300622e2028SHong Zhang           cols[ncols_upper] = i;
1301622e2028SHong Zhang           cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
1302622e2028SHong Zhang           ncols_upper++;
1303622e2028SHong Zhang         }
1304622e2028SHong Zhang       }
1305622e2028SHong Zhang       ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13067a48dd6fSHong Zhang       nzk += nlnk;
13077a48dd6fSHong Zhang 
13087a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13097a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13107a48dd6fSHong Zhang 
13117a48dd6fSHong Zhang       while (prow < k){
13127a48dd6fSHong Zhang         nextprow = jl[prow];
13137a48dd6fSHong Zhang 
13147a48dd6fSHong Zhang         /* merge prow into k-th row */
13157a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13167a48dd6fSHong Zhang         jmax = ui[prow+1];
13177a48dd6fSHong Zhang         ncols = jmax-jmin;
1318ed59401aSHong Zhang         i     = jmin - ui[prow];
1319ed59401aSHong Zhang         cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1320ed59401aSHong Zhang         for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
13212ed133dbSHong Zhang         ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13227a48dd6fSHong Zhang         nzk += nlnk;
13237a48dd6fSHong Zhang 
13247a48dd6fSHong Zhang         /* update il and jl for prow */
13257a48dd6fSHong Zhang         if (jmin < jmax){
13267a48dd6fSHong Zhang           il[prow] = jmin;
13277a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13287a48dd6fSHong Zhang         }
13297a48dd6fSHong Zhang         prow = nextprow;
13307a48dd6fSHong Zhang       }
13317a48dd6fSHong Zhang 
13327a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13337a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13347a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13357a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
13367a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
13377a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
13387a48dd6fSHong Zhang         reallocs++;
13397a48dd6fSHong Zhang       }
13407a48dd6fSHong Zhang 
13417a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13422ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13437a48dd6fSHong Zhang 
13447a48dd6fSHong Zhang       /* add the k-th row into il and jl */
134539e3d630SHong Zhang       if (nzk > 1){
13467a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13477a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13487a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13497a48dd6fSHong Zhang       }
13507a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13517a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13527a48dd6fSHong Zhang 
13537a48dd6fSHong Zhang       current_space->array           += nzk;
13547a48dd6fSHong Zhang       current_space->local_used      += nzk;
13557a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13567a48dd6fSHong Zhang 
13577a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13587a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13597a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13607a48dd6fSHong Zhang 
13617a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13627a48dd6fSHong Zhang     }
13637a48dd6fSHong Zhang 
13647a48dd6fSHong Zhang     if (ai[am] != 0) {
136539e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
13667a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
13677a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
13687a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
13697a48dd6fSHong Zhang     } else {
1370ed59401aSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
13717a48dd6fSHong Zhang     }
13727a48dd6fSHong Zhang 
13737a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13747a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13757a48dd6fSHong Zhang     ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
13767a48dd6fSHong Zhang 
13777a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13787a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13797a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13802ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13817a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13827a48dd6fSHong Zhang 
138339e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
138439e3d630SHong Zhang 
13857a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13867a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1387ed59401aSHong Zhang   B = *fact;
1388ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1389ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
13907a48dd6fSHong Zhang 
1391ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13927a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
13937a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13947a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
13957a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
13967a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
13977a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13987a48dd6fSHong Zhang   b->j    = uj;
13997a48dd6fSHong Zhang   b->i    = ui;
14007a48dd6fSHong Zhang   b->diag = 0;
14017a48dd6fSHong Zhang   b->ilen = 0;
14027a48dd6fSHong Zhang   b->imax = 0;
14037a48dd6fSHong Zhang   b->row  = perm;
14047a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
14057a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14067a48dd6fSHong Zhang   b->icol = perm;
14077a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14087a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1409ed59401aSHong Zhang   PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
14107a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
14117a48dd6fSHong Zhang 
1412ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1413ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1414ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14157a48dd6fSHong Zhang   if (ai[am] != 0) {
1416ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14177a48dd6fSHong Zhang   } else {
1418ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14197a48dd6fSHong Zhang   }
142039e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1421622e2028SHong Zhang   if (perm_identity){
1422ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1423ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1424622e2028SHong Zhang   }
1425a6175056SHong Zhang   PetscFunctionReturn(0);
1426a6175056SHong Zhang }
1427d5d45c9bSBarry Smith 
1428f76d2b81SHong Zhang #undef __FUNCT__
1429f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1430dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1431f76d2b81SHong Zhang {
1432f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
143310c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
14342ed133dbSHong Zhang   Mat            B;
1435dfbe8321SBarry Smith   PetscErrorCode ierr;
1436f76d2b81SHong Zhang   PetscTruth     perm_identity;
143710c27e3dSHong Zhang   PetscReal      fill = info->fill;
14382ed133dbSHong Zhang   PetscInt       *rip,*riip,i,mbs=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
143910c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14402ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
144110c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
144210c27e3dSHong Zhang   PetscBT        lnkbt;
14432ed133dbSHong Zhang   IS             iperm;
1444f76d2b81SHong Zhang 
1445f76d2b81SHong Zhang   PetscFunctionBegin;
14462ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1447f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14482ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14492ed133dbSHong Zhang 
1450f76d2b81SHong Zhang   if (!perm_identity){
14512ed133dbSHong Zhang     /* check if perm is symmetric! */
14522ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14532ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14542ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14552ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14562ed133dbSHong Zhang     }
14572ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14582ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1459f76d2b81SHong Zhang   }
146010c27e3dSHong Zhang 
146110c27e3dSHong Zhang   /* initialization */
14622ed133dbSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
146310c27e3dSHong Zhang   ui[0] = 0;
146410c27e3dSHong Zhang 
146510c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14662ed133dbSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
14672ed133dbSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
14682ed133dbSHong Zhang   il     = jl + mbs;
14692ed133dbSHong Zhang   cols   = il + mbs;
14702ed133dbSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
14712ed133dbSHong Zhang   for (i=0; i<mbs; i++){
14722ed133dbSHong Zhang     jl[i] = mbs; il[i] = 0;
1473f76d2b81SHong Zhang   }
147410c27e3dSHong Zhang 
147510c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14762ed133dbSHong Zhang   nlnk = mbs + 1;
14772ed133dbSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
147810c27e3dSHong Zhang 
14792ed133dbSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
14802ed133dbSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
148110c27e3dSHong Zhang   current_space = free_space;
148210c27e3dSHong Zhang 
14832ed133dbSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
148410c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
148510c27e3dSHong Zhang     nzk   = 0;
14862ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14872ed133dbSHong Zhang     ncols_upper = 0;
14882ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1489622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14902ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14912ed133dbSHong Zhang         cols[ncols_upper] = i;
14922ed133dbSHong Zhang         ncols_upper++;
14932ed133dbSHong Zhang       }
14942ed133dbSHong Zhang     }
14952ed133dbSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
149610c27e3dSHong Zhang     nzk += nlnk;
149710c27e3dSHong Zhang 
149810c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
149910c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
150010c27e3dSHong Zhang 
150110c27e3dSHong Zhang     while (prow < k){
150210c27e3dSHong Zhang       nextprow = jl[prow];
150310c27e3dSHong Zhang       /* merge prow into k-th row */
15042ed133dbSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
150510c27e3dSHong Zhang       jmax = ui[prow+1];
150610c27e3dSHong Zhang       ncols = jmax-jmin;
15072ed133dbSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
15082ed133dbSHong Zhang       ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150910c27e3dSHong Zhang       nzk += nlnk;
151010c27e3dSHong Zhang 
151110c27e3dSHong Zhang       /* update il and jl for prow */
151210c27e3dSHong Zhang       if (jmin < jmax){
151310c27e3dSHong Zhang         il[prow] = jmin;
15142ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
151510c27e3dSHong Zhang       }
151610c27e3dSHong Zhang       prow = nextprow;
151710c27e3dSHong Zhang     }
151810c27e3dSHong Zhang 
151910c27e3dSHong Zhang     /* if free space is not available, make more free space */
152010c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15212ed133dbSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
152210c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
152310c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
152410c27e3dSHong Zhang       reallocs++;
152510c27e3dSHong Zhang     }
152610c27e3dSHong Zhang 
152710c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15282ed133dbSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
152910c27e3dSHong Zhang 
153010c27e3dSHong Zhang     /* add the k-th row into il and jl */
153110c27e3dSHong Zhang     if (nzk-1 > 0){
15322ed133dbSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
153310c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
153410c27e3dSHong Zhang       il[k] = ui[k] + 1;
153510c27e3dSHong Zhang     }
15362ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
153710c27e3dSHong Zhang     current_space->array           += nzk;
153810c27e3dSHong Zhang     current_space->local_used      += nzk;
153910c27e3dSHong Zhang     current_space->local_remaining -= nzk;
154010c27e3dSHong Zhang 
154110c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
154210c27e3dSHong Zhang   }
154310c27e3dSHong Zhang 
15442ed133dbSHong Zhang   if (ai[mbs] != 0) {
154578910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[mbs]-mbs))/((PetscReal)ai[mbs]);
154678910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
154778910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
154878910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
154910c27e3dSHong Zhang   } else {
155078910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
155110c27e3dSHong Zhang   }
155210c27e3dSHong Zhang 
155310c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
155410c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
155510c27e3dSHong Zhang 
155610c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15572ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
155810c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
155910c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
156010c27e3dSHong Zhang 
156110c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15622ed133dbSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
15632ed133dbSHong Zhang   B    = *fact;
15642ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
15652ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
156610c27e3dSHong Zhang 
15672ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
156810c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
156910c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
157010c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
157110c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
157210c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15732ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
157410c27e3dSHong Zhang   b->j    = uj;
157510c27e3dSHong Zhang   b->i    = ui;
157610c27e3dSHong Zhang   b->diag = 0;
157710c27e3dSHong Zhang   b->ilen = 0;
157810c27e3dSHong Zhang   b->imax = 0;
157910c27e3dSHong Zhang   b->row  = perm;
158010c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
158110c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
158210c27e3dSHong Zhang   b->icol = perm;
158310c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15842ed133dbSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15852ed133dbSHong Zhang   PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
15862ed133dbSHong Zhang   b->maxnz = b->nz = ui[mbs];
158710c27e3dSHong Zhang 
15882ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15892ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15902ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15912ed133dbSHong Zhang   if (ai[mbs] != 0) {
15922ed133dbSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
159310c27e3dSHong Zhang   } else {
15942ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
159510c27e3dSHong Zhang   }
159639e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15972ed133dbSHong Zhang   if (perm_identity){
159810c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
159910c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
16002ed133dbSHong Zhang   }
1601f76d2b81SHong Zhang   PetscFunctionReturn(0);
1602f76d2b81SHong Zhang }
1603