xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 391eac555c6c020ce3eb8b85384fe86947e3b5af)
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 
423d2276718SHong Zhang #include "src/ksp/pc/impls/factor/factor.h"
4244a2ae208SSatish Balay #undef __FUNCT__
4254a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
426af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
427289bc588SBarry Smith {
428416022c9SBarry Smith   Mat            C=*B;
429416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
43082bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4316849ba73SBarry Smith   PetscErrorCode ierr;
43238baddfdSBarry Smith   PetscInt       *r,*ic,i,j,n=A->m,*ai=b->i,*aj=b->j;
43338baddfdSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ics;
434d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
43587828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
436ee45ca4aSHong Zhang   PetscReal      zeropivot,rs,d,shift_nonzero;
437d2276718SHong Zhang   PetscReal      row_shift,shift_top=0.;
438d2276718SHong Zhang   PetscTruth     shift_pd;
439d2276718SHong Zhang   Shift_Ctx      sctx;
440d2276718SHong Zhang   PetscInt       newshift;
441289bc588SBarry Smith 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
4430a29876aSHong Zhang   shift_nonzero  = info->shiftnz;
4440a29876aSHong Zhang   shift_pd       = info->shiftpd;
4450a29876aSHong Zhang   zeropivot      = info->zeropivot;
4460a29876aSHong Zhang 
44778b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44878b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44987828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45087828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
451010a20caSHong Zhang   rtmps = rtmp; ics = ic;
452289bc588SBarry Smith 
4536cc28720Svictorle   if (!a->diag) {
4540cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4550cf777b8SBarry Smith   }
456e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
457e24cfe64SHong Zhang   if (shift_pd && shift_nonzero) shift_nonzero = 0.0;
458e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45938baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4600cf777b8SBarry Smith     shift_top = 0;
4616cc28720Svictorle     for (i=0; i<n; i++) {
462010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4636cc28720Svictorle       /* calculate amt of shift needed for this row */
4646c037d1bSvictorle       if (d<=0) {
4653aeef889SHong Zhang         row_shift = 0;
4663aeef889SHong Zhang       } else {
4673aeef889SHong Zhang         row_shift = -2*d;
4683aeef889SHong Zhang       }
469010a20caSHong Zhang       v  = a->a+aai[i];
470e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
471e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4726cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4736cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4746cc28720Svictorle     }
4756cc28720Svictorle   }
4766cc28720Svictorle 
477d2276718SHong Zhang   sctx.shift_amount = 0;
478d2276718SHong Zhang   sctx.shift_top    = shift_top;
479d2276718SHong Zhang   sctx.nshift       = 0;
480d2276718SHong Zhang   sctx.nshift_max   = 5;
481d2276718SHong Zhang   sctx.shift_lo     = 0.;
482d2276718SHong Zhang   sctx.shift_hi     = 1.;
483085256b3SBarry Smith   do {
484d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
485289bc588SBarry Smith     for (i=0; i<n; i++){
486289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
487010a20caSHong Zhang       ajtmp = aj + ai[i];
48848e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
489289bc588SBarry Smith 
490289bc588SBarry Smith       /* load in initial (unfactored row) */
491416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
492010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
493010a20caSHong Zhang       v        = a->a + a->i[r[i]];
494085256b3SBarry Smith       for (j=0; j<nz; j++) {
495085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
496085256b3SBarry Smith       }
497d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
498289bc588SBarry Smith 
499010a20caSHong Zhang       row = *ajtmp++;
500f2582111SSatish Balay       while  (row < i) {
5018c37ef55SBarry Smith         pc = rtmp + row;
5028c37ef55SBarry Smith         if (*pc != 0.0) {
503010a20caSHong Zhang           pv         = b->a + diag_offset[row];
504010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50535aab85fSBarry Smith           multiplier = *pc / *pv++;
5068c37ef55SBarry Smith           *pc        = multiplier;
507f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
508f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
509b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
510289bc588SBarry Smith         }
511010a20caSHong Zhang         row = *ajtmp++;
512289bc588SBarry Smith       }
513416022c9SBarry Smith       /* finished row so stick it into b->a */
514010a20caSHong Zhang       pv   = b->a + ai[i] ;
515010a20caSHong Zhang       pj   = b->j + ai[i] ;
5168c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
51735aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
518d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
519d3d32019SBarry Smith       rs   = 0.0;
520d3d32019SBarry Smith       for (j=0; j<nz; j++) {
521d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
522d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
523d3d32019SBarry Smith       }
524d2276718SHong Zhang 
525d2276718SHong Zhang       sctx.rs  = rs;
526d2276718SHong Zhang       sctx.pv  = pv[diag];
527d2276718SHong Zhang       newshift = 0;
528d2276718SHong Zhang       ierr = PCLUFactorCheckShift(A,info,B,&sctx,&newshift);CHKERRQ(ierr);
529d2276718SHong Zhang       if (newshift == 1){
530e24cfe64SHong Zhang         break;
531d2276718SHong Zhang       } else if (newshift == -1){
532d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
533d708749eSBarry Smith       }
53435aab85fSBarry Smith     }
535d2276718SHong Zhang 
536d2276718SHong Zhang     if (shift_pd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5376cc28720Svictorle       /*
5386c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5396cc28720Svictorle        * then try lower shift
5406cc28720Svictorle        */
541d2276718SHong Zhang       sctx.shift_hi       = info->shift_fraction;
542d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
543d2276718SHong Zhang       sctx.shift_amount   = info->shift_fraction * sctx.shift_top;
544d2276718SHong Zhang       sctx.lushift        = PETSC_TRUE;
545d2276718SHong Zhang       sctx.nshift++;
5466cc28720Svictorle     }
547d2276718SHong Zhang   } while (sctx.lushift);
5480f11f9aeSSatish Balay 
54935aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55035aab85fSBarry Smith   for (i=0; i<n; i++) {
551010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55235aab85fSBarry Smith   }
55335aab85fSBarry Smith 
554606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55578b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
557416022c9SBarry Smith   C->factor = FACTOR_LU;
5587dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
559c456f294SBarry Smith   C->assembled = PETSC_TRUE;
560b0a32e0cSBarry Smith   PetscLogFlops(C->n);
561d2276718SHong Zhang   if (sctx.nshift){
5620697b6caSHong Zhang     if (shift_nonzero) {
563d2276718SHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
5640697b6caSHong Zhang     } else if (shift_pd) {
565d2276718SHong Zhang       b->lu_shift_fraction = info->shift_fraction;
566d2276718SHong 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);
5676cc28720Svictorle     }
5680697b6caSHong Zhang   }
5693a40ed3dSBarry Smith   PetscFunctionReturn(0);
570289bc588SBarry Smith }
5710889b5dcSvictorle 
5720889b5dcSvictorle #undef __FUNCT__
5730889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
574dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5750889b5dcSvictorle {
5760889b5dcSvictorle   PetscFunctionBegin;
5770889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5780889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5790889b5dcSvictorle   PetscFunctionReturn(0);
5800889b5dcSvictorle }
5810889b5dcSvictorle 
5820889b5dcSvictorle 
5830a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5844a2ae208SSatish Balay #undef __FUNCT__
5854a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
586dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
587da3a660dSBarry Smith {
588dfbe8321SBarry Smith   PetscErrorCode ierr;
589416022c9SBarry Smith   Mat            C;
590416022c9SBarry Smith 
5913a40ed3dSBarry Smith   PetscFunctionBegin;
5929e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
593af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
594273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
595440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
5963a40ed3dSBarry Smith   PetscFunctionReturn(0);
597da3a660dSBarry Smith }
5980a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5994a2ae208SSatish Balay #undef __FUNCT__
6004a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
601dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6028c37ef55SBarry Smith {
603416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
604416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6056849ba73SBarry Smith   PetscErrorCode ierr;
60638baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
60738baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
60887828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6098c37ef55SBarry Smith 
6103a40ed3dSBarry Smith   PetscFunctionBegin;
6113a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
61291bf9adeSBarry Smith 
6131ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6141ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
615416022c9SBarry Smith   tmp  = a->solve_work;
6168c37ef55SBarry Smith 
6173b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6183b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6198c37ef55SBarry Smith 
6208c37ef55SBarry Smith   /* forward solve the lower triangular */
6218c37ef55SBarry Smith   tmp[0] = b[*r++];
622010a20caSHong Zhang   tmps   = tmp;
6238c37ef55SBarry Smith   for (i=1; i<n; i++) {
624010a20caSHong Zhang     v   = aa + ai[i] ;
625010a20caSHong Zhang     vi  = aj + ai[i] ;
62653e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
62753e7425aSBarry Smith     sum = b[*r++];
628ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6298c37ef55SBarry Smith     tmp[i] = sum;
6308c37ef55SBarry Smith   }
6318c37ef55SBarry Smith 
6328c37ef55SBarry Smith   /* backward solve the upper triangular */
6338c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
634010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
635010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
636416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6378c37ef55SBarry Smith     sum = tmp[i];
638ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
639010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6408c37ef55SBarry Smith   }
6418c37ef55SBarry Smith 
6423b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6433b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6441ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6451ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
646b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6473a40ed3dSBarry Smith   PetscFunctionReturn(0);
6488c37ef55SBarry Smith }
649026e39d0SSatish Balay 
650930ae53cSBarry Smith /* ----------------------------------------------------------- */
6514a2ae208SSatish Balay #undef __FUNCT__
6524a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
653dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
654930ae53cSBarry Smith {
655930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6566849ba73SBarry Smith   PetscErrorCode ierr;
65738baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
658362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
659aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
66038baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
661362ced78SSatish Balay   PetscScalar    *v,sum;
662d85afc42SSatish Balay #endif
663930ae53cSBarry Smith 
6643a40ed3dSBarry Smith   PetscFunctionBegin;
6653a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
666930ae53cSBarry Smith 
6671ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6681ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
669930ae53cSBarry Smith 
670aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6711c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6721c4feecaSSatish Balay #else
673930ae53cSBarry Smith   /* forward solve the lower triangular */
674930ae53cSBarry Smith   x[0] = b[0];
675930ae53cSBarry Smith   for (i=1; i<n; i++) {
676930ae53cSBarry Smith     ai_i = ai[i];
677930ae53cSBarry Smith     v    = aa + ai_i;
678930ae53cSBarry Smith     vi   = aj + ai_i;
679930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
680930ae53cSBarry Smith     sum  = b[i];
681930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
682930ae53cSBarry Smith     x[i] = sum;
683930ae53cSBarry Smith   }
684930ae53cSBarry Smith 
685930ae53cSBarry Smith   /* backward solve the upper triangular */
686930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
687930ae53cSBarry Smith     adiag_i = adiag[i];
688d708749eSBarry Smith     v       = aa + adiag_i + 1;
689d708749eSBarry Smith     vi      = aj + adiag_i + 1;
690930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
691930ae53cSBarry Smith     sum     = x[i];
692930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
693930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
694930ae53cSBarry Smith   }
6951c4feecaSSatish Balay #endif
696b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6971ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6981ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6993a40ed3dSBarry Smith   PetscFunctionReturn(0);
700930ae53cSBarry Smith }
701930ae53cSBarry Smith 
7024a2ae208SSatish Balay #undef __FUNCT__
7034a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
704dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
705da3a660dSBarry Smith {
706416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
707416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7086849ba73SBarry Smith   PetscErrorCode ierr;
70938baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
71038baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
71187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
712da3a660dSBarry Smith 
7133a40ed3dSBarry Smith   PetscFunctionBegin;
71478b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
715da3a660dSBarry Smith 
7161ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7171ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
718416022c9SBarry Smith   tmp  = a->solve_work;
719da3a660dSBarry Smith 
7203b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7213b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
722da3a660dSBarry Smith 
723da3a660dSBarry Smith   /* forward solve the lower triangular */
724da3a660dSBarry Smith   tmp[0] = b[*r++];
725da3a660dSBarry Smith   for (i=1; i<n; i++) {
726010a20caSHong Zhang     v   = aa + ai[i] ;
727010a20caSHong Zhang     vi  = aj + ai[i] ;
728416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
729da3a660dSBarry Smith     sum = b[*r++];
730010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
731da3a660dSBarry Smith     tmp[i] = sum;
732da3a660dSBarry Smith   }
733da3a660dSBarry Smith 
734da3a660dSBarry Smith   /* backward solve the upper triangular */
735da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
736010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
737010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
738416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
739da3a660dSBarry Smith     sum = tmp[i];
740010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
741010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
742da3a660dSBarry Smith     x[*c--] += tmp[i];
743da3a660dSBarry Smith   }
744da3a660dSBarry Smith 
7453b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7463b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7471ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7481ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
749b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
750898183ecSLois Curfman McInnes 
7513a40ed3dSBarry Smith   PetscFunctionReturn(0);
752da3a660dSBarry Smith }
753da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7544a2ae208SSatish Balay #undef __FUNCT__
7554a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
756dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
757da3a660dSBarry Smith {
758416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7592235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7606849ba73SBarry Smith   PetscErrorCode ierr;
76138baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
76238baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
76387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
764da3a660dSBarry Smith 
7653a40ed3dSBarry Smith   PetscFunctionBegin;
7661ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7671ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
768416022c9SBarry Smith   tmp  = a->solve_work;
769da3a660dSBarry Smith 
7702235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7712235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
772da3a660dSBarry Smith 
773da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7742235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
775da3a660dSBarry Smith 
776da3a660dSBarry Smith   /* forward solve the U^T */
777da3a660dSBarry Smith   for (i=0; i<n; i++) {
778010a20caSHong Zhang     v   = aa + diag[i] ;
779010a20caSHong Zhang     vi  = aj + diag[i] + 1;
780f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
781c38d4ed2SBarry Smith     s1  = tmp[i];
782ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
783da3a660dSBarry Smith     while (nz--) {
784010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
785da3a660dSBarry Smith     }
786c38d4ed2SBarry Smith     tmp[i] = s1;
787da3a660dSBarry Smith   }
788da3a660dSBarry Smith 
789da3a660dSBarry Smith   /* backward solve the L^T */
790da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
791010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
792010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
793f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
794f1af5d2fSBarry Smith     s1  = tmp[i];
795da3a660dSBarry Smith     while (nz--) {
796010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
797da3a660dSBarry Smith     }
798da3a660dSBarry Smith   }
799da3a660dSBarry Smith 
800da3a660dSBarry Smith   /* copy tmp into x according to permutation */
801da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
802da3a660dSBarry Smith 
8032235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8042235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8051ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
807da3a660dSBarry Smith 
808b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8093a40ed3dSBarry Smith   PetscFunctionReturn(0);
810da3a660dSBarry Smith }
811da3a660dSBarry Smith 
8124a2ae208SSatish Balay #undef __FUNCT__
8134a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
814dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
815da3a660dSBarry Smith {
816416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8172235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8186849ba73SBarry Smith   PetscErrorCode ierr;
81938baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
82038baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
82187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8226abc6512SBarry Smith 
8233a40ed3dSBarry Smith   PetscFunctionBegin;
82423bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8256abc6512SBarry Smith 
8261ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8271ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
828416022c9SBarry Smith   tmp = a->solve_work;
8296abc6512SBarry Smith 
8302235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8312235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8326abc6512SBarry Smith 
8336abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8342235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8356abc6512SBarry Smith 
8366abc6512SBarry Smith   /* forward solve the U^T */
8376abc6512SBarry Smith   for (i=0; i<n; i++) {
838010a20caSHong Zhang     v   = aa + diag[i] ;
839010a20caSHong Zhang     vi  = aj + diag[i] + 1;
840f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8416abc6512SBarry Smith     tmp[i] *= *v++;
8426abc6512SBarry Smith     while (nz--) {
843010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8446abc6512SBarry Smith     }
8456abc6512SBarry Smith   }
8466abc6512SBarry Smith 
8476abc6512SBarry Smith   /* backward solve the L^T */
8486abc6512SBarry Smith   for (i=n-1; i>=0; i--){
849010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
850010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
851f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8526abc6512SBarry Smith     while (nz--) {
853010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8546abc6512SBarry Smith     }
8556abc6512SBarry Smith   }
8566abc6512SBarry Smith 
8576abc6512SBarry Smith   /* copy tmp into x according to permutation */
8586abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8596abc6512SBarry Smith 
8602235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8612235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8621ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8631ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8646abc6512SBarry Smith 
865b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8663a40ed3dSBarry Smith   PetscFunctionReturn(0);
867da3a660dSBarry Smith }
8689e25ed09SBarry Smith /* ----------------------------------------------------------------*/
869dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
87008480c60SBarry Smith 
8714a2ae208SSatish Balay #undef __FUNCT__
8724a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
873dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8749e25ed09SBarry Smith {
875416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8769e25ed09SBarry Smith   IS             isicol;
8776849ba73SBarry Smith   PetscErrorCode ierr;
87838baddfdSBarry Smith   PetscInt       *r,*ic,prow,n = A->m,*ai = a->i,*aj = a->j;
87938baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
880418422e8SSatish Balay   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,len, reallocs = 0,dcount = 0;
88138baddfdSBarry Smith   PetscInt       incrlev,nnz,i,levels,diagonal_fill;
88277c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
883329f5518SBarry Smith   PetscReal      f;
8849e25ed09SBarry Smith 
8853a40ed3dSBarry Smith   PetscFunctionBegin;
8866cf997b0SBarry Smith   f             = info->fill;
88738baddfdSBarry Smith   levels        = (PetscInt)info->levels;
88838baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8894c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
89082bf6240SBarry Smith 
89108480c60SBarry Smith   /* special case that simply copies fill pattern */
892be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
893be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
89486bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8952e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
89608480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
89708480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
89808480c60SBarry Smith     if (!b->diag) {
8997c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90008480c60SBarry Smith     }
9017c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90208480c60SBarry Smith     b->row              = isrow;
90308480c60SBarry Smith     b->col              = iscol;
90482bf6240SBarry Smith     b->icol             = isicol;
90587828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
906f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
907c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
908c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9093a40ed3dSBarry Smith     PetscFunctionReturn(0);
91008480c60SBarry Smith   }
91108480c60SBarry Smith 
912898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
913898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9149e25ed09SBarry Smith 
9159e25ed09SBarry Smith   /* get new row pointers */
91638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
917010a20caSHong Zhang   ainew[0] = 0;
9189e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
91938baddfdSBarry Smith   jmax = (PetscInt)(f*(ai[n]+1));
92038baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
9219e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
92238baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajfill);CHKERRQ(ierr);
9239e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
92438baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
92556beaf04SBarry Smith   /* im is level for each filled value */
92638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
92756beaf04SBarry Smith   /* dloc is location of diagonal in factor */
92838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&dloc);CHKERRQ(ierr);
92956beaf04SBarry Smith   dloc[0]  = 0;
93056beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9316cf997b0SBarry Smith 
9326cf997b0SBarry Smith     /* copy current row into linked list */
93356beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
93429bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
935010a20caSHong Zhang     xi      = aj + ai[r[prow]];
9369e25ed09SBarry Smith     fill[n] = n;
937435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9389e25ed09SBarry Smith     while (nz--) {
9399e25ed09SBarry Smith       fm  = n;
940010a20caSHong Zhang       idx = ic[*xi++ ];
9419e25ed09SBarry Smith       do {
9429e25ed09SBarry Smith         m  = fm;
9439e25ed09SBarry Smith         fm = fill[m];
9449e25ed09SBarry Smith       } while (fm < idx);
9459e25ed09SBarry Smith       fill[m]   = idx;
9469e25ed09SBarry Smith       fill[idx] = fm;
94756beaf04SBarry Smith       im[idx]   = 0;
9489e25ed09SBarry Smith     }
9496cf997b0SBarry Smith 
9506cf997b0SBarry Smith     /* make sure diagonal entry is included */
951435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9526cf997b0SBarry Smith       fm = n;
953435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
954435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
955435faa5fSBarry Smith       fill[fm]   = prow;
9566cf997b0SBarry Smith       im[prow]   = 0;
9576cf997b0SBarry Smith       nzf++;
958335d9088SBarry Smith       dcount++;
9596cf997b0SBarry Smith     }
9606cf997b0SBarry Smith 
96156beaf04SBarry Smith     nzi = 0;
9629e25ed09SBarry Smith     row = fill[n];
96356beaf04SBarry Smith     while (row < prow) {
96456beaf04SBarry Smith       incrlev = im[row] + 1;
96556beaf04SBarry Smith       nz      = dloc[row];
966010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
967010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
968416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
96956beaf04SBarry Smith       fm      = row;
97056beaf04SBarry Smith       while (nnz-- > 0) {
971010a20caSHong Zhang         idx = *xi++ ;
97256beaf04SBarry Smith         if (*flev + incrlev > levels) {
97356beaf04SBarry Smith           flev++;
97456beaf04SBarry Smith           continue;
97556beaf04SBarry Smith         }
97656beaf04SBarry Smith         do {
9779e25ed09SBarry Smith           m  = fm;
9789e25ed09SBarry Smith           fm = fill[m];
97956beaf04SBarry Smith         } while (fm < idx);
9809e25ed09SBarry Smith         if (fm != idx) {
98156beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9829e25ed09SBarry Smith           fill[m]   = idx;
9839e25ed09SBarry Smith           fill[idx] = fm;
9849e25ed09SBarry Smith           fm        = idx;
98556beaf04SBarry Smith           nzf++;
986ecf371e4SBarry Smith         } else {
98756beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9889e25ed09SBarry Smith         }
98956beaf04SBarry Smith         flev++;
9909e25ed09SBarry Smith       }
9919e25ed09SBarry Smith       row = fill[row];
99256beaf04SBarry Smith       nzi++;
9939e25ed09SBarry Smith     }
9949e25ed09SBarry Smith     /* copy new filled row into permanent storage */
99556beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
996010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
997ecf371e4SBarry Smith 
998ecf371e4SBarry Smith       /* estimate how much additional space we will need */
999ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
1000ecf371e4SBarry Smith       /* just double the memory each time */
100138baddfdSBarry Smith       /*  maxadd = (PetscInt)((f*(ai[n]+!shift)*(n-prow+5))/n); */
100238baddfdSBarry Smith       PetscInt maxadd = jmax;
100356beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1004bbb6d6a8SBarry Smith       jmax += maxadd;
1005ecf371e4SBarry Smith 
1006ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
100738baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
100838baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1009606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
101056beaf04SBarry Smith       ajnew  = xi;
101138baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
101238baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1013606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
101456beaf04SBarry Smith       ajfill = xi;
1015418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
10169e25ed09SBarry Smith     }
1017010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1018010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
101956beaf04SBarry Smith     dloc[prow]  = nzi;
10209e25ed09SBarry Smith     fm          = fill[n];
102156beaf04SBarry Smith     while (nzf--) {
1022010a20caSHong Zhang       *xi++   = fm ;
102356beaf04SBarry Smith       *flev++ = im[fm];
10249e25ed09SBarry Smith       fm      = fill[fm];
10259e25ed09SBarry Smith     }
10266cf997b0SBarry Smith     /* make sure row has diagonal entry */
1027010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
102877431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10296cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10306cf997b0SBarry Smith     }
10319e25ed09SBarry Smith   }
1032606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1033898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1034898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1035606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1036606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10379e25ed09SBarry Smith 
1038f64065bbSBarry Smith   {
1039329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1040418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1041c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1042c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1043b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1044335d9088SBarry Smith     if (diagonal_fill) {
104577431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1046335d9088SBarry Smith     }
1047f64065bbSBarry Smith   }
1048bbb6d6a8SBarry Smith 
10499e25ed09SBarry Smith   /* put together the new matrix */
1050f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1051f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1052f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
1053b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1054416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1055606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10567c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1057010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10589e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1059606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1060606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1061b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1062416022c9SBarry Smith   b->j          = ajnew;
1063416022c9SBarry Smith   b->i          = ainew;
106456beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1065416022c9SBarry Smith   b->diag       = dloc;
1066416022c9SBarry Smith   b->ilen       = 0;
1067416022c9SBarry Smith   b->imax       = 0;
1068416022c9SBarry Smith   b->row        = isrow;
1069416022c9SBarry Smith   b->col        = iscol;
1070c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1071c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
107282bf6240SBarry Smith   b->icol       = isicol;
107387828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1074416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
107556beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
107638baddfdSBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
1077010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
10789e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10793a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10803a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1081418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1082ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1083329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10843a40ed3dSBarry Smith   PetscFunctionReturn(0);
10859e25ed09SBarry Smith }
1086d5d45c9bSBarry Smith 
10873c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1088a6175056SHong Zhang #undef __FUNCT__
1089a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1090af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1091a6175056SHong Zhang {
10922ed133dbSHong Zhang   Mat            C = *B;
10930968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10942ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
10952ed133dbSHong Zhang   IS             ip=b->row;
10962ed133dbSHong Zhang   PetscErrorCode ierr;
10972ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
10982ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
10992ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1100622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1101ee45ca4aSHong Zhang   PetscReal      zeropivot,shift_amount,rs,shift_nonzero;
1102ee45ca4aSHong Zhang   PetscTruth     chshift,shift_pd;
110339e3d630SHong Zhang   PetscInt       nshift=0;
1104d5d45c9bSBarry Smith 
1105a6175056SHong Zhang   PetscFunctionBegin;
1106ee45ca4aSHong Zhang   shift_nonzero  = info->shiftnz;
1107ee45ca4aSHong Zhang   shift_pd       = info->shiftpd;
1108ee45ca4aSHong Zhang   zeropivot      = info->zeropivot;
1109ee45ca4aSHong Zhang 
11102ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11112ed133dbSHong Zhang 
11122ed133dbSHong Zhang   /* initialization */
11132ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11142ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11152ed133dbSHong Zhang   jl   = il + mbs;
11162ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11172ed133dbSHong Zhang 
11182ed133dbSHong Zhang   shift_amount = 0;
11192ed133dbSHong Zhang   do {
11202ed133dbSHong Zhang     chshift = PETSC_FALSE;
11212ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11222ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1123a6175056SHong Zhang     }
11242ed133dbSHong Zhang 
11252ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1126c05c3958SHong Zhang       bval = ba + bi[k];
11272ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11282ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11292ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11302ed133dbSHong Zhang         col = rip[aj[j]];
11312ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11322ed133dbSHong Zhang           rtmp[col] = aa[j];
1133c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11342ed133dbSHong Zhang         }
11352ed133dbSHong Zhang       }
113639e3d630SHong Zhang       /* shift the diagonal of the matrix */
113739e3d630SHong Zhang       if (nshift) rtmp[k] += shift_amount;
11382ed133dbSHong Zhang 
11392ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11402ed133dbSHong Zhang       dk = rtmp[k];
11412ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11422ed133dbSHong Zhang 
11432ed133dbSHong Zhang       while (i < k){
11442ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11452ed133dbSHong Zhang 
11462ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11472ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11482ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11492ed133dbSHong Zhang         dk += uikdi*ba[ili];
11502ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11512ed133dbSHong Zhang 
11522ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11532ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11542ed133dbSHong Zhang         if (jmin < jmax){
11552ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11562ed133dbSHong Zhang           /* update il and jl for row i */
11572ed133dbSHong Zhang           il[i] = jmin;
11582ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11592ed133dbSHong Zhang         }
11602ed133dbSHong Zhang         i = nexti;
11612ed133dbSHong Zhang       }
11622ed133dbSHong Zhang 
11630697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11640697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11650697b6caSHong Zhang       rs   = 0.0;
11663c9e8598SHong Zhang       jmin = bi[k]+1;
11673c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11683c9e8598SHong Zhang       if (nz){
11693c9e8598SHong Zhang         bcol = bj + jmin;
11703c9e8598SHong Zhang         while (nz--){
11713c9e8598SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol++]);
11723c9e8598SHong Zhang         }
11733c9e8598SHong Zhang       }
11740697b6caSHong Zhang       if (PetscAbsScalar(dk) <= zeropivot*rs && shift_nonzero){
11750697b6caSHong Zhang         /* force |diag(*B)| > zeropivot*rs */
11760697b6caSHong Zhang         if (!nshift){
11770697b6caSHong Zhang           shift_amount = shift_nonzero;
11780697b6caSHong Zhang         } else {
11790697b6caSHong Zhang           shift_amount *= 2.0;
11800697b6caSHong Zhang         }
11810697b6caSHong Zhang         chshift = PETSC_TRUE;
11820697b6caSHong Zhang         nshift++;
11830697b6caSHong Zhang         break;
11840697b6caSHong Zhang       } else if (PetscRealPart(dk) <= zeropivot*rs && shift_pd){
11850697b6caSHong Zhang         /* calculate a shift that would make this row diagonally dominant */
11860697b6caSHong Zhang 	shift_amount = PetscMax(rs+PetscAbs(PetscRealPart(dk)),1.1*shift_amount);
11873c9e8598SHong Zhang 	chshift      = PETSC_TRUE;
11883c9e8598SHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
11893c9e8598SHong Zhang 	   we increase the shift until it converges. There is no guarantee that
11903c9e8598SHong Zhang 	   this algorithm converges faster or slower, or is better or worse
11913c9e8598SHong Zhang 	   than the ILU algorithm. */
11923c9e8598SHong Zhang 	nshift++;
11933c9e8598SHong Zhang 	break;
11940697b6caSHong Zhang       } else if (PetscAbsScalar(dk) <= zeropivot*rs){
11950697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11963c9e8598SHong Zhang       }
11973c9e8598SHong Zhang 
11983c9e8598SHong Zhang       /* copy data into U(k,:) */
119939e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
120039e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
120139e3d630SHong Zhang       if (jmin < jmax) {
120239e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
120339e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12043c9e8598SHong Zhang         }
120539e3d630SHong Zhang         /* add the k-th row into il and jl */
12063c9e8598SHong Zhang         il[k] = jmin;
12073c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12083c9e8598SHong Zhang       }
12093c9e8598SHong Zhang     }
12100697b6caSHong Zhang   } while (chshift);
12113c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12123c9e8598SHong Zhang 
121339e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12143c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12153c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12163c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
12173c9e8598SHong Zhang   PetscLogFlops(C->m);
12183c9e8598SHong Zhang   if (nshift){
12190697b6caSHong Zhang     if (shift_nonzero) {
122039e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",nshift,shift_amount);
12210697b6caSHong Zhang     } else if (shift_pd) {
122239e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_pd tries %D, shift_amount %g\n",nshift,shift_amount);
12230697b6caSHong Zhang     }
12243c9e8598SHong Zhang   }
12253c9e8598SHong Zhang   PetscFunctionReturn(0);
12263c9e8598SHong Zhang }
1227a6175056SHong Zhang 
12287a48dd6fSHong Zhang #include "petscbt.h"
12297a48dd6fSHong Zhang #include "src/mat/utils/freespace.h"
1230a6175056SHong Zhang #undef __FUNCT__
1231a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1232dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1233a6175056SHong Zhang {
12340968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1235ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1236ed59401aSHong Zhang   Mat            B;
1237dfbe8321SBarry Smith   PetscErrorCode ierr;
1238653a6975SHong Zhang   PetscTruth     perm_identity;
1239622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1240ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1241*391eac55SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1242*391eac55SHong Zhang   PetscInt       ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
1243ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
12447a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
12457a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12467a48dd6fSHong Zhang   PetscBT        lnkbt;
1247a6175056SHong Zhang 
1248a6175056SHong Zhang   PetscFunctionBegin;
1249653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12507a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12517a48dd6fSHong Zhang 
125239e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
125339e3d630SHong Zhang   ui[0] = 0;
125439e3d630SHong Zhang 
1255622e2028SHong Zhang   /* special case that simply copies fill pattern */
1256622e2028SHong Zhang   if (!levels && perm_identity) {
1257622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1258ed59401aSHong Zhang     for (i=0; i<am; i++) {
125939e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1260ed59401aSHong Zhang     }
126139e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
126239e3d630SHong Zhang     cols = uj;
1263ed59401aSHong Zhang     for (i=0; i<am; i++) {
1264ed59401aSHong Zhang       aj    = a->j + a->diag[i];
126539e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
126639e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1267ed59401aSHong Zhang     }
126839e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12697a48dd6fSHong Zhang     /* initialization */
1270622e2028SHong Zhang     ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
12717a48dd6fSHong Zhang 
12727a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12737a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12747a48dd6fSHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
12757a48dd6fSHong Zhang     il         = jl + am;
12767a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12777a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12787a48dd6fSHong Zhang     for (i=0; i<am; i++){
12797a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12807a48dd6fSHong Zhang     }
12817a48dd6fSHong Zhang 
12827a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12837a48dd6fSHong Zhang     nlnk = am + 1;
12842ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12857a48dd6fSHong Zhang 
12867a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12877a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12887a48dd6fSHong Zhang     current_space = free_space;
12897a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12907a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12917a48dd6fSHong Zhang 
12927a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12937a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12947a48dd6fSHong Zhang       nzk   = 0;
1295622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1296622e2028SHong Zhang       ncols_upper = 0;
1297622e2028SHong Zhang       cols        = cols_lvl + am;
1298622e2028SHong Zhang       for (j=0; j<ncols; j++){
1299622e2028SHong Zhang         i = rip[*(aj + ai[rip[k]] + j)];
1300622e2028SHong Zhang         if (i >= k){ /* only take upper triangular entry */
1301622e2028SHong Zhang           cols[ncols_upper] = i;
1302622e2028SHong Zhang           cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
1303622e2028SHong Zhang           ncols_upper++;
1304622e2028SHong Zhang         }
1305622e2028SHong Zhang       }
1306622e2028SHong Zhang       ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13077a48dd6fSHong Zhang       nzk += nlnk;
13087a48dd6fSHong Zhang 
13097a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13107a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13117a48dd6fSHong Zhang 
13127a48dd6fSHong Zhang       while (prow < k){
13137a48dd6fSHong Zhang         nextprow = jl[prow];
13147a48dd6fSHong Zhang 
13157a48dd6fSHong Zhang         /* merge prow into k-th row */
13167a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13177a48dd6fSHong Zhang         jmax = ui[prow+1];
13187a48dd6fSHong Zhang         ncols = jmax-jmin;
1319ed59401aSHong Zhang         i     = jmin - ui[prow];
1320ed59401aSHong Zhang         cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1321ed59401aSHong Zhang         for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
13222ed133dbSHong Zhang         ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13237a48dd6fSHong Zhang         nzk += nlnk;
13247a48dd6fSHong Zhang 
13257a48dd6fSHong Zhang         /* update il and jl for prow */
13267a48dd6fSHong Zhang         if (jmin < jmax){
13277a48dd6fSHong Zhang           il[prow] = jmin;
13287a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13297a48dd6fSHong Zhang         }
13307a48dd6fSHong Zhang         prow = nextprow;
13317a48dd6fSHong Zhang       }
13327a48dd6fSHong Zhang 
13337a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13347a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13357a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13367a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
13377a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
13387a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
13397a48dd6fSHong Zhang         reallocs++;
13407a48dd6fSHong Zhang       }
13417a48dd6fSHong Zhang 
13427a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13432ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13447a48dd6fSHong Zhang 
13457a48dd6fSHong Zhang       /* add the k-th row into il and jl */
134639e3d630SHong Zhang       if (nzk > 1){
13477a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13487a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13497a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13507a48dd6fSHong Zhang       }
13517a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13527a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13537a48dd6fSHong Zhang 
13547a48dd6fSHong Zhang       current_space->array           += nzk;
13557a48dd6fSHong Zhang       current_space->local_used      += nzk;
13567a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13577a48dd6fSHong Zhang 
13587a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13597a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13607a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13617a48dd6fSHong Zhang 
13627a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13637a48dd6fSHong Zhang     }
13647a48dd6fSHong Zhang 
13657a48dd6fSHong Zhang     if (ai[am] != 0) {
136639e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
13677a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
13687a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
13697a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
13707a48dd6fSHong Zhang     } else {
1371ed59401aSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
13727a48dd6fSHong Zhang     }
13737a48dd6fSHong Zhang 
13747a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13757a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13767a48dd6fSHong Zhang     ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
13777a48dd6fSHong Zhang 
13787a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13797a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13807a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13812ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13827a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13837a48dd6fSHong Zhang 
138439e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
138539e3d630SHong Zhang 
13867a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13877a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1388ed59401aSHong Zhang   B = *fact;
1389ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1390ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
13917a48dd6fSHong Zhang 
1392ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13937a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
13947a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13957a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
13967a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
13977a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
13987a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13997a48dd6fSHong Zhang   b->j    = uj;
14007a48dd6fSHong Zhang   b->i    = ui;
14017a48dd6fSHong Zhang   b->diag = 0;
14027a48dd6fSHong Zhang   b->ilen = 0;
14037a48dd6fSHong Zhang   b->imax = 0;
14047a48dd6fSHong Zhang   b->row  = perm;
14057a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
14067a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14077a48dd6fSHong Zhang   b->icol = perm;
14087a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14097a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1410ed59401aSHong Zhang   PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
14117a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
14127a48dd6fSHong Zhang 
1413ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1414ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1415ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14167a48dd6fSHong Zhang   if (ai[am] != 0) {
1417ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14187a48dd6fSHong Zhang   } else {
1419ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14207a48dd6fSHong Zhang   }
142139e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1422622e2028SHong Zhang   if (perm_identity){
1423ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1424ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1425622e2028SHong Zhang   }
1426a6175056SHong Zhang   PetscFunctionReturn(0);
1427a6175056SHong Zhang }
1428d5d45c9bSBarry Smith 
1429f76d2b81SHong Zhang #undef __FUNCT__
1430f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1431dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1432f76d2b81SHong Zhang {
1433f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
143410c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
14352ed133dbSHong Zhang   Mat            B;
1436dfbe8321SBarry Smith   PetscErrorCode ierr;
1437f76d2b81SHong Zhang   PetscTruth     perm_identity;
143810c27e3dSHong Zhang   PetscReal      fill = info->fill;
14392ed133dbSHong Zhang   PetscInt       *rip,*riip,i,mbs=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
144010c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14412ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
144210c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
144310c27e3dSHong Zhang   PetscBT        lnkbt;
14442ed133dbSHong Zhang   IS             iperm;
1445f76d2b81SHong Zhang 
1446f76d2b81SHong Zhang   PetscFunctionBegin;
14472ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1448f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14492ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14502ed133dbSHong Zhang 
1451f76d2b81SHong Zhang   if (!perm_identity){
14522ed133dbSHong Zhang     /* check if perm is symmetric! */
14532ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14542ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14552ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14562ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14572ed133dbSHong Zhang     }
14582ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14592ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1460f76d2b81SHong Zhang   }
146110c27e3dSHong Zhang 
146210c27e3dSHong Zhang   /* initialization */
14632ed133dbSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
146410c27e3dSHong Zhang   ui[0] = 0;
146510c27e3dSHong Zhang 
146610c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14672ed133dbSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
14682ed133dbSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
14692ed133dbSHong Zhang   il     = jl + mbs;
14702ed133dbSHong Zhang   cols   = il + mbs;
14712ed133dbSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
14722ed133dbSHong Zhang   for (i=0; i<mbs; i++){
14732ed133dbSHong Zhang     jl[i] = mbs; il[i] = 0;
1474f76d2b81SHong Zhang   }
147510c27e3dSHong Zhang 
147610c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14772ed133dbSHong Zhang   nlnk = mbs + 1;
14782ed133dbSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
147910c27e3dSHong Zhang 
14802ed133dbSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
14812ed133dbSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
148210c27e3dSHong Zhang   current_space = free_space;
148310c27e3dSHong Zhang 
14842ed133dbSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
148510c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
148610c27e3dSHong Zhang     nzk   = 0;
14872ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14882ed133dbSHong Zhang     ncols_upper = 0;
14892ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1490622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14912ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14922ed133dbSHong Zhang         cols[ncols_upper] = i;
14932ed133dbSHong Zhang         ncols_upper++;
14942ed133dbSHong Zhang       }
14952ed133dbSHong Zhang     }
14962ed133dbSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
149710c27e3dSHong Zhang     nzk += nlnk;
149810c27e3dSHong Zhang 
149910c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
150010c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
150110c27e3dSHong Zhang 
150210c27e3dSHong Zhang     while (prow < k){
150310c27e3dSHong Zhang       nextprow = jl[prow];
150410c27e3dSHong Zhang       /* merge prow into k-th row */
15052ed133dbSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
150610c27e3dSHong Zhang       jmax = ui[prow+1];
150710c27e3dSHong Zhang       ncols = jmax-jmin;
15082ed133dbSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
15092ed133dbSHong Zhang       ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
151010c27e3dSHong Zhang       nzk += nlnk;
151110c27e3dSHong Zhang 
151210c27e3dSHong Zhang       /* update il and jl for prow */
151310c27e3dSHong Zhang       if (jmin < jmax){
151410c27e3dSHong Zhang         il[prow] = jmin;
15152ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
151610c27e3dSHong Zhang       }
151710c27e3dSHong Zhang       prow = nextprow;
151810c27e3dSHong Zhang     }
151910c27e3dSHong Zhang 
152010c27e3dSHong Zhang     /* if free space is not available, make more free space */
152110c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15222ed133dbSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
152310c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
152410c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
152510c27e3dSHong Zhang       reallocs++;
152610c27e3dSHong Zhang     }
152710c27e3dSHong Zhang 
152810c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15292ed133dbSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
153010c27e3dSHong Zhang 
153110c27e3dSHong Zhang     /* add the k-th row into il and jl */
153210c27e3dSHong Zhang     if (nzk-1 > 0){
15332ed133dbSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
153410c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
153510c27e3dSHong Zhang       il[k] = ui[k] + 1;
153610c27e3dSHong Zhang     }
15372ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
153810c27e3dSHong Zhang     current_space->array           += nzk;
153910c27e3dSHong Zhang     current_space->local_used      += nzk;
154010c27e3dSHong Zhang     current_space->local_remaining -= nzk;
154110c27e3dSHong Zhang 
154210c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
154310c27e3dSHong Zhang   }
154410c27e3dSHong Zhang 
15452ed133dbSHong Zhang   if (ai[mbs] != 0) {
154678910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[mbs]-mbs))/((PetscReal)ai[mbs]);
154778910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
154878910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
154978910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
155010c27e3dSHong Zhang   } else {
155178910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
155210c27e3dSHong Zhang   }
155310c27e3dSHong Zhang 
155410c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
155510c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
155610c27e3dSHong Zhang 
155710c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15582ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
155910c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
156010c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
156110c27e3dSHong Zhang 
156210c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15632ed133dbSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
15642ed133dbSHong Zhang   B    = *fact;
15652ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
15662ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
156710c27e3dSHong Zhang 
15682ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
156910c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
157010c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
157110c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
157210c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
157310c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15742ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
157510c27e3dSHong Zhang   b->j    = uj;
157610c27e3dSHong Zhang   b->i    = ui;
157710c27e3dSHong Zhang   b->diag = 0;
157810c27e3dSHong Zhang   b->ilen = 0;
157910c27e3dSHong Zhang   b->imax = 0;
158010c27e3dSHong Zhang   b->row  = perm;
158110c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
158210c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
158310c27e3dSHong Zhang   b->icol = perm;
158410c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15852ed133dbSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15862ed133dbSHong Zhang   PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
15872ed133dbSHong Zhang   b->maxnz = b->nz = ui[mbs];
158810c27e3dSHong Zhang 
15892ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15902ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15912ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15922ed133dbSHong Zhang   if (ai[mbs] != 0) {
15932ed133dbSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
159410c27e3dSHong Zhang   } else {
15952ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
159610c27e3dSHong Zhang   }
159739e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15982ed133dbSHong Zhang   if (perm_identity){
159910c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
160010c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
16012ed133dbSHong Zhang   }
1602f76d2b81SHong Zhang   PetscFunctionReturn(0);
1603f76d2b81SHong Zhang }
1604