xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision fd62cfe2599e3e567f75a22cddcef37895525805)
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"
51393bc97SHong Zhang #include "petscbt.h"
61393bc97SHong Zhang #include "src/mat/utils/freespace.h"
73b2fbd54SBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
10dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
113b2fbd54SBarry Smith {
123a40ed3dSBarry Smith   PetscFunctionBegin;
133a40ed3dSBarry Smith 
1429bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
15aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
16d64ed03dSBarry Smith   PetscFunctionReturn(0);
17d64ed03dSBarry Smith #endif
183b2fbd54SBarry Smith }
193b2fbd54SBarry Smith 
2086bacbd2SBarry Smith 
21dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2286bacbd2SBarry Smith 
2338baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
249cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
2686bacbd2SBarry Smith 
274a2ae208SSatish Balay #undef __FUNCT__
284a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2986bacbd2SBarry Smith   /* ------------------------------------------------------------
3086bacbd2SBarry Smith 
3186bacbd2SBarry Smith           This interface was contribed by Tony Caola
3286bacbd2SBarry Smith 
3386bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
345bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
355bd2ddc7SBarry Smith      SPARSEKIT2.
365bd2ddc7SBarry Smith 
375bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
385bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
3986bacbd2SBarry Smith 
4086bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4186bacbd2SBarry Smith      help in getting this routine ironed out.
4286bacbd2SBarry Smith 
435bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
445bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
455bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
465bd2ddc7SBarry Smith      Yousef Saad's code.
4786bacbd2SBarry Smith 
4886bacbd2SBarry Smith      ------------------------------------------------------------
4986bacbd2SBarry Smith */
50dfbe8321SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5186bacbd2SBarry Smith {
5260ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5360ec86bdSBarry Smith   PetscFunctionBegin;
54e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5560ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5660ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5760ec86bdSBarry Smith #else
5886bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
5907d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6086bacbd2SBarry Smith   PetscTruth     reorder;
619cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
629cc1f4e3SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->m;
6338baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6438baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6538baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6687828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6754a8161fSBarry Smith   PetscReal      af;
6886bacbd2SBarry Smith 
6986bacbd2SBarry Smith   PetscFunctionBegin;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7238baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7386bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7433259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7538baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7638baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7786bacbd2SBarry Smith 
7886bacbd2SBarry Smith 
7986bacbd2SBarry Smith   /* ------------------------------------------------------------
8086bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8186bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8286bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8386bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8486bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8586bacbd2SBarry Smith      ------------------------------------------------------------
8686bacbd2SBarry Smith   */
8786bacbd2SBarry Smith 
8886bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
8986bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9086bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9186bacbd2SBarry Smith   reorder = PetscNot(reorder);
9286bacbd2SBarry Smith 
9386bacbd2SBarry Smith 
9486bacbd2SBarry Smith   /* storage for ilu factor */
9538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9638baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
9787828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
9838baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
9986bacbd2SBarry Smith 
10086bacbd2SBarry Smith   /* ------------------------------------------------------------
10186bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10286bacbd2SBarry Smith      ------------------------------------------------------------
10386bacbd2SBarry Smith   */
104b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
105b19713cbSBarry Smith     old_j[i]++;
10686bacbd2SBarry Smith   }
107b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
108b19713cbSBarry Smith     old_i[i]++;
109b19713cbSBarry Smith   };
110010a20caSHong Zhang 
11186bacbd2SBarry Smith 
11286bacbd2SBarry Smith   if (reorder) {
113c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
114c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
115c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
116c0c2c81eSBarry Smith       r[i]  = r[i]+1;
117c0c2c81eSBarry Smith       c[i]  = c[i]+1;
118c0c2c81eSBarry Smith     }
11938baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12038baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12187828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1225bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
123c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
124c0c2c81eSBarry Smith       r[i]  = r[i]-1;
125c0c2c81eSBarry Smith       c[i]  = c[i]-1;
126c0c2c81eSBarry Smith     }
127c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
128c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
129b19713cbSBarry Smith     o_a = old_a2;
130b19713cbSBarry Smith     o_j = old_j2;
131b19713cbSBarry Smith     o_i = old_i2;
132b19713cbSBarry Smith   } else {
133b19713cbSBarry Smith     o_a = old_a;
134b19713cbSBarry Smith     o_j = old_j;
135b19713cbSBarry Smith     o_i = old_i;
13686bacbd2SBarry Smith   }
13786bacbd2SBarry Smith 
13886bacbd2SBarry Smith   /* ------------------------------------------------------------
13986bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14086bacbd2SBarry Smith      ------------------------------------------------------------
14186bacbd2SBarry Smith   */
14286bacbd2SBarry Smith 
14338baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14587828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14686bacbd2SBarry Smith 
14754a8161fSBarry 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);
1486849ba73SBarry Smith   if (sierr) {
1496849ba73SBarry Smith     switch (sierr) {
15077431f27SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
15177431f27SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
152e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
153e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15477431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15577431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15686bacbd2SBarry Smith     }
15786bacbd2SBarry Smith   }
15886bacbd2SBarry Smith 
15986bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16086bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   /* ------------------------------------------------------------
16386bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16486bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16586bacbd2SBarry Smith      ------------------------------------------------------------
16686bacbd2SBarry Smith   */
16786bacbd2SBarry Smith 
16887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16938baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17086bacbd2SBarry Smith 
1715bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17286bacbd2SBarry Smith 
17386bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17486bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   if (reorder) {
17786bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
180313ae024SBarry Smith   } else {
181b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
182f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
183b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
184b19713cbSBarry Smith     }
185b19713cbSBarry Smith   }
186b19713cbSBarry Smith 
187b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18886bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
189b19713cbSBarry Smith     old_i[i]--;
190b19713cbSBarry Smith   }
191b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
192b19713cbSBarry Smith     old_j[i]--;
193b19713cbSBarry Smith   }
19486bacbd2SBarry Smith 
195b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19686bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
197b19713cbSBarry Smith     new_i[i]--;
198b19713cbSBarry Smith   }
199b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
200b19713cbSBarry Smith     new_j[i]--;
201b19713cbSBarry Smith   }
20286bacbd2SBarry Smith 
20386bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20486bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2054c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2064c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
207c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20886bacbd2SBarry Smith   for(i=0; i<n; i++) {
20986bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21086bacbd2SBarry Smith   };
21186bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21207d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21386bacbd2SBarry Smith 
214c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
215c0c2c81eSBarry Smith 
216329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21707d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21886bacbd2SBarry Smith 
21986bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22086bacbd2SBarry Smith 
221f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
222f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
223f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
22486bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22586bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22686bacbd2SBarry Smith 
22786bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
22886bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
22986bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23007d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
231b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23286bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
23386bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
23486bacbd2SBarry Smith   b->a             = new_a;
23586bacbd2SBarry Smith   b->j             = new_j;
23686bacbd2SBarry Smith   b->i             = new_i;
23786bacbd2SBarry Smith   b->ilen          = 0;
23886bacbd2SBarry Smith   b->imax          = 0;
2391f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
240313ae024SBarry Smith   b->row           = isirow;
24186bacbd2SBarry Smith   b->col           = iscolf;
24287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24386bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24486bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24586bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24686bacbd2SBarry Smith 
24786bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24886bacbd2SBarry Smith 
249431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
250b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
251b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
252b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
253b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
254431e710aSBarry Smith 
2554846f1f5SKris Buschelman   ierr = MatILUDTFactor_Inode(A,info,isrow,iscol,fact);CHKERRQ(ierr);
25671c2f376SKris Buschelman 
25786bacbd2SBarry Smith   PetscFunctionReturn(0);
25860ec86bdSBarry Smith #endif
25986bacbd2SBarry Smith }
26086bacbd2SBarry Smith 
2614a2ae208SSatish Balay #undef __FUNCT__
262b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
263dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
264289bc588SBarry Smith {
265416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
266289bc588SBarry Smith   IS             isicol;
2676849ba73SBarry Smith   PetscErrorCode ierr;
26838baddfdSBarry Smith   PetscInt       *r,*ic,i,n=A->m,*ai=a->i,*aj=a->j;
2691393bc97SHong Zhang   PetscInt       *bi,*bj,*ajtmp;
2701393bc97SHong Zhang   PetscInt       *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2719e046878SBarry Smith   PetscReal      f;
2721393bc97SHong Zhang   PetscInt       nlnk,*lnk,k,*cols,**bi_ptr;
2731393bc97SHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
2741393bc97SHong Zhang   PetscBT        lnkbt;
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 */
2831393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2841393bc97SHong Zhang   bi[0] = 0;
2851393bc97SHong Zhang 
2861393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2871393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2881393bc97SHong Zhang   bdiag[0] = 0;
2891393bc97SHong Zhang 
2901393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2911393bc97SHong Zhang   nlnk = n + 1;
2921393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2931393bc97SHong Zhang 
2941393bc97SHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&cols);CHKERRQ(ierr);
2951393bc97SHong Zhang   im     = cols + n;
2961393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2971393bc97SHong Zhang 
2981393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
299d3d32019SBarry Smith   f = info->fill;
3001393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3011393bc97SHong Zhang   current_space = free_space;
302289bc588SBarry Smith 
303289bc588SBarry Smith   for (i=0; i<n; i++) {
3041393bc97SHong Zhang     /* copy previous fill into linked list */
305289bc588SBarry Smith     nzi = 0;
3061393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3071393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3081393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
3091393bc97SHong Zhang     for (k=0; k<nnz; k++) cols[k] = ic[*(ajtmp+k)]; /* note: cols is not sorted when iscol!=indentity */
3101393bc97SHong Zhang     ierr = PetscLLAdd(nnz,cols,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3111393bc97SHong Zhang     nzi += nlnk;
3121393bc97SHong Zhang 
3131393bc97SHong Zhang     /* add pivot rows into linked list */
3141393bc97SHong Zhang     row = lnk[n];
3151393bc97SHong Zhang     while (row < i) {
3161393bc97SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1;
3171393bc97SHong Zhang       ajtmp   = bi_ptr[row] + nzbd;
3181393bc97SHong Zhang       nnz     = im[row] - nzbd; /* num of columns with row<indices<=i */
3191393bc97SHong Zhang       im[row] = nzbd;
3201393bc97SHong Zhang       ierr = PetscLLAddSortedLU(nnz,ajtmp,row,nlnk,lnk,lnkbt,i,nzbd);CHKERRQ(ierr);
3211393bc97SHong Zhang       nzi     += nlnk;
3221393bc97SHong Zhang       im[row] += nzbd;  /* update im[row]: num of cols with index<=i */
3231393bc97SHong Zhang 
3241393bc97SHong Zhang       row = lnk[row];
325289bc588SBarry Smith     }
3261393bc97SHong Zhang 
3271393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3281393bc97SHong Zhang     im[i]   = nzi;
3291393bc97SHong Zhang 
3301393bc97SHong Zhang     /* mark bdiag */
3311393bc97SHong Zhang     nzbd = 0;
3321393bc97SHong Zhang     nnz  = nzi;
3331393bc97SHong Zhang     k    = lnk[n];
3341393bc97SHong Zhang     while (nnz-- && k < i){
3351393bc97SHong Zhang       nzbd++;
3361393bc97SHong Zhang       k = lnk[k];
3371393bc97SHong Zhang     }
3381393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3391393bc97SHong Zhang 
3401393bc97SHong Zhang     /* if free space is not available, make more free space */
3411393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3421393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
3431393bc97SHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
3441393bc97SHong Zhang       reallocs++;
3451393bc97SHong Zhang     }
3461393bc97SHong Zhang 
3471393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3481393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3491393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3501393bc97SHong Zhang     current_space->array           += nzi;
3511393bc97SHong Zhang     current_space->local_used      += nzi;
3521393bc97SHong Zhang     current_space->local_remaining -= nzi;
353289bc588SBarry Smith   }
354464e8d44SSatish Balay   if (ai[n] != 0) {
3551393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
356418422e8SSatish Balay     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
357b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
358b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
359b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
36005bf559cSSatish Balay   } else {
361b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
36205bf559cSSatish Balay   }
3632fb3ed76SBarry Smith 
364898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
365898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3661987afe7SBarry Smith 
3671393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3681393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
3691393bc97SHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3701393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3711393bc97SHong Zhang   ierr = PetscFree(cols);CHKERRQ(ierr);
372289bc588SBarry Smith 
373289bc588SBarry Smith   /* put together the new matrix */
374f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr);
375f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
376f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr);
37752e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
378416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
379606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3807c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
381e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
382606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
383606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
3841393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3851393bc97SHong Zhang   b->j          = bj;
3861393bc97SHong Zhang   b->i          = bi;
3871393bc97SHong Zhang   b->diag       = bdiag;
388416022c9SBarry Smith   b->ilen       = 0;
389416022c9SBarry Smith   b->imax       = 0;
390416022c9SBarry Smith   b->row        = isrow;
391416022c9SBarry Smith   b->col        = iscol;
392c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
393c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
39482bf6240SBarry Smith   b->icol       = isicol;
39587828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3961393bc97SHong Zhang 
3971393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3981393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3991393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4007fd98bd6SLois Curfman McInnes 
401329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
402418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
403ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
404703deb49SSatish Balay 
405e93ccc4dSBarry Smith   if (ai[n] != 0) {
4061393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
407eea2913aSSatish Balay   } else {
408eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
409eea2913aSSatish Balay   }
4104846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
41171c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
413289bc588SBarry Smith }
4141393bc97SHong Zhang 
4150a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4164a2ae208SSatish Balay #undef __FUNCT__
4174a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
418af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
419289bc588SBarry Smith {
420416022c9SBarry Smith   Mat            C=*B;
421416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
42282bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4236849ba73SBarry Smith   PetscErrorCode ierr;
424b3bf805bSHong Zhang   PetscInt       *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j;
425b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
426d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
42787828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
428540703acSHong Zhang   PetscReal      zeropivot,rs,d,shift_nz;
429d2276718SHong Zhang   PetscReal      row_shift,shift_top=0.;
430d2276718SHong Zhang   PetscTruth     shift_pd;
431b3bf805bSHong Zhang   LUShift_Ctx    sctx;
432d2276718SHong Zhang   PetscInt       newshift;
433289bc588SBarry Smith 
4343a40ed3dSBarry Smith   PetscFunctionBegin;
435540703acSHong Zhang   shift_nz  = info->shiftnz;
4360a29876aSHong Zhang   shift_pd  = info->shiftpd;
4370a29876aSHong Zhang   zeropivot = info->zeropivot;
4380a29876aSHong Zhang 
43978b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44078b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44187828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44287828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
443010a20caSHong Zhang   rtmps = rtmp; ics = ic;
444289bc588SBarry Smith 
4456cc28720Svictorle   if (!a->diag) {
4460cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4470cf777b8SBarry Smith   }
448e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
449540703acSHong Zhang   if (shift_pd && shift_nz) shift_nz = 0.0;
450e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45138baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4520cf777b8SBarry Smith     shift_top = 0;
4536cc28720Svictorle     for (i=0; i<n; i++) {
454010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4556cc28720Svictorle       /* calculate amt of shift needed for this row */
4566c037d1bSvictorle       if (d<=0) {
4573aeef889SHong Zhang         row_shift = 0;
4583aeef889SHong Zhang       } else {
4593aeef889SHong Zhang         row_shift = -2*d;
4603aeef889SHong Zhang       }
461010a20caSHong Zhang       v  = a->a+aai[i];
462e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
463e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4646cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4656cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4666cc28720Svictorle     }
4676cc28720Svictorle   }
4686cc28720Svictorle 
469d2276718SHong Zhang   sctx.shift_amount = 0;
470d2276718SHong Zhang   sctx.shift_top    = shift_top;
471d2276718SHong Zhang   sctx.nshift       = 0;
472d2276718SHong Zhang   sctx.nshift_max   = 5;
473d2276718SHong Zhang   sctx.shift_lo     = 0.;
474d2276718SHong Zhang   sctx.shift_hi     = 1.;
475085256b3SBarry Smith   do {
476d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
477289bc588SBarry Smith     for (i=0; i<n; i++){
478b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
479b3bf805bSHong Zhang       bjtmp = bj + bi[i];
480b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
481289bc588SBarry Smith 
482289bc588SBarry Smith       /* load in initial (unfactored row) */
483416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
484b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
485010a20caSHong Zhang       v     = a->a + a->i[r[i]];
486085256b3SBarry Smith       for (j=0; j<nz; j++) {
487b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
488085256b3SBarry Smith       }
489d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
490289bc588SBarry Smith 
491b3bf805bSHong Zhang       row = *bjtmp++;
492f2582111SSatish Balay       while  (row < i) {
4938c37ef55SBarry Smith         pc = rtmp + row;
4948c37ef55SBarry Smith         if (*pc != 0.0) {
495010a20caSHong Zhang           pv         = b->a + diag_offset[row];
496010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
49735aab85fSBarry Smith           multiplier = *pc / *pv++;
4988c37ef55SBarry Smith           *pc        = multiplier;
499b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
500f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
501efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
502289bc588SBarry Smith         }
503b3bf805bSHong Zhang         row = *bjtmp++;
504289bc588SBarry Smith       }
505416022c9SBarry Smith       /* finished row so stick it into b->a */
506b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
507b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
508b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
509b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
510d3d32019SBarry Smith       rs   = 0.0;
511d3d32019SBarry Smith       for (j=0; j<nz; j++) {
512d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
513d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
514d3d32019SBarry Smith       }
515d2276718SHong Zhang 
516b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
517d2276718SHong Zhang       sctx.rs  = rs;
518d2276718SHong Zhang       sctx.pv  = pv[diag];
5193e8c821fSHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
520d2276718SHong Zhang       if (newshift == 1){
521b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
522d2276718SHong Zhang       } else if (newshift == -1){
523d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
524d708749eSBarry Smith       }
52535aab85fSBarry Smith     }
526d2276718SHong Zhang 
527d2276718SHong Zhang     if (shift_pd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5286cc28720Svictorle       /*
5296c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5306cc28720Svictorle        * then try lower shift
5316cc28720Svictorle        */
532d2276718SHong Zhang       sctx.shift_hi       = info->shift_fraction;
533d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
534d2276718SHong Zhang       sctx.shift_amount   = info->shift_fraction * sctx.shift_top;
535d2276718SHong Zhang       sctx.lushift        = PETSC_TRUE;
536d2276718SHong Zhang       sctx.nshift++;
5376cc28720Svictorle     }
538d2276718SHong Zhang   } while (sctx.lushift);
5390f11f9aeSSatish Balay 
54035aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
54135aab85fSBarry Smith   for (i=0; i<n; i++) {
542010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
54335aab85fSBarry Smith   }
54435aab85fSBarry Smith 
545606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
54678b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
54778b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
548416022c9SBarry Smith   C->factor = FACTOR_LU;
5497dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
550c456f294SBarry Smith   C->assembled = PETSC_TRUE;
551efee365bSSatish Balay   ierr = PetscLogFlops(C->n);CHKERRQ(ierr);
552d2276718SHong Zhang   if (sctx.nshift){
553540703acSHong Zhang     if (shift_nz) {
554540703acSHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
5550697b6caSHong Zhang     } else if (shift_pd) {
556d2276718SHong Zhang       b->lu_shift_fraction = info->shift_fraction;
557d2276718SHong 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);
5586cc28720Svictorle     }
5590697b6caSHong Zhang   }
5603a40ed3dSBarry Smith   PetscFunctionReturn(0);
561289bc588SBarry Smith }
5620889b5dcSvictorle 
5630889b5dcSvictorle #undef __FUNCT__
5640889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
565dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5660889b5dcSvictorle {
5670889b5dcSvictorle   PetscFunctionBegin;
5680889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5690889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5700889b5dcSvictorle   PetscFunctionReturn(0);
5710889b5dcSvictorle }
5720889b5dcSvictorle 
5730889b5dcSvictorle 
5740a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5754a2ae208SSatish Balay #undef __FUNCT__
5764a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
577dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
578da3a660dSBarry Smith {
579dfbe8321SBarry Smith   PetscErrorCode ierr;
580416022c9SBarry Smith   Mat            C;
581416022c9SBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
5839e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
584af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
585273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
58652e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5873a40ed3dSBarry Smith   PetscFunctionReturn(0);
588da3a660dSBarry Smith }
5890a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5904a2ae208SSatish Balay #undef __FUNCT__
5914a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
592dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5938c37ef55SBarry Smith {
594416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
595416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
5966849ba73SBarry Smith   PetscErrorCode ierr;
59738baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
59838baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
59987828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6008c37ef55SBarry Smith 
6013a40ed3dSBarry Smith   PetscFunctionBegin;
6023a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
60391bf9adeSBarry Smith 
6041ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6051ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
606416022c9SBarry Smith   tmp  = a->solve_work;
6078c37ef55SBarry Smith 
6083b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6093b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6108c37ef55SBarry Smith 
6118c37ef55SBarry Smith   /* forward solve the lower triangular */
6128c37ef55SBarry Smith   tmp[0] = b[*r++];
613010a20caSHong Zhang   tmps   = tmp;
6148c37ef55SBarry Smith   for (i=1; i<n; i++) {
615010a20caSHong Zhang     v   = aa + ai[i] ;
616010a20caSHong Zhang     vi  = aj + ai[i] ;
61753e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
61853e7425aSBarry Smith     sum = b[*r++];
619ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6208c37ef55SBarry Smith     tmp[i] = sum;
6218c37ef55SBarry Smith   }
6228c37ef55SBarry Smith 
6238c37ef55SBarry Smith   /* backward solve the upper triangular */
6248c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
625010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
626010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
627416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6288c37ef55SBarry Smith     sum = tmp[i];
629ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
630010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6318c37ef55SBarry Smith   }
6328c37ef55SBarry Smith 
6333b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6343b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6351ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6361ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
637efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
6383a40ed3dSBarry Smith   PetscFunctionReturn(0);
6398c37ef55SBarry Smith }
640026e39d0SSatish Balay 
641930ae53cSBarry Smith /* ----------------------------------------------------------- */
6424a2ae208SSatish Balay #undef __FUNCT__
6434a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
644dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
645930ae53cSBarry Smith {
646930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6476849ba73SBarry Smith   PetscErrorCode ierr;
64838baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
649362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
650aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
65138baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
652362ced78SSatish Balay   PetscScalar    *v,sum;
653d85afc42SSatish Balay #endif
654930ae53cSBarry Smith 
6553a40ed3dSBarry Smith   PetscFunctionBegin;
6563a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
657930ae53cSBarry Smith 
6581ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6591ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
660930ae53cSBarry Smith 
661aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6621c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6631c4feecaSSatish Balay #else
664930ae53cSBarry Smith   /* forward solve the lower triangular */
665930ae53cSBarry Smith   x[0] = b[0];
666930ae53cSBarry Smith   for (i=1; i<n; i++) {
667930ae53cSBarry Smith     ai_i = ai[i];
668930ae53cSBarry Smith     v    = aa + ai_i;
669930ae53cSBarry Smith     vi   = aj + ai_i;
670930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
671930ae53cSBarry Smith     sum  = b[i];
672930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
673930ae53cSBarry Smith     x[i] = sum;
674930ae53cSBarry Smith   }
675930ae53cSBarry Smith 
676930ae53cSBarry Smith   /* backward solve the upper triangular */
677930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
678930ae53cSBarry Smith     adiag_i = adiag[i];
679d708749eSBarry Smith     v       = aa + adiag_i + 1;
680d708749eSBarry Smith     vi      = aj + adiag_i + 1;
681930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
682930ae53cSBarry Smith     sum     = x[i];
683930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
684930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
685930ae53cSBarry Smith   }
6861c4feecaSSatish Balay #endif
687efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
6881ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6903a40ed3dSBarry Smith   PetscFunctionReturn(0);
691930ae53cSBarry Smith }
692930ae53cSBarry Smith 
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
695dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
696da3a660dSBarry Smith {
697416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
698416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6996849ba73SBarry Smith   PetscErrorCode ierr;
70038baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
70138baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
70287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
703da3a660dSBarry Smith 
7043a40ed3dSBarry Smith   PetscFunctionBegin;
70578b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
706da3a660dSBarry Smith 
7071ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7081ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
709416022c9SBarry Smith   tmp  = a->solve_work;
710da3a660dSBarry Smith 
7113b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7123b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
713da3a660dSBarry Smith 
714da3a660dSBarry Smith   /* forward solve the lower triangular */
715da3a660dSBarry Smith   tmp[0] = b[*r++];
716da3a660dSBarry Smith   for (i=1; i<n; i++) {
717010a20caSHong Zhang     v   = aa + ai[i] ;
718010a20caSHong Zhang     vi  = aj + ai[i] ;
719416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
720da3a660dSBarry Smith     sum = b[*r++];
721010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
722da3a660dSBarry Smith     tmp[i] = sum;
723da3a660dSBarry Smith   }
724da3a660dSBarry Smith 
725da3a660dSBarry Smith   /* backward solve the upper triangular */
726da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
727010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
728010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
729416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
730da3a660dSBarry Smith     sum = tmp[i];
731010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
732010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
733da3a660dSBarry Smith     x[*c--] += tmp[i];
734da3a660dSBarry Smith   }
735da3a660dSBarry Smith 
7363b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7373b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7381ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7391ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
740efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
741898183ecSLois Curfman McInnes 
7423a40ed3dSBarry Smith   PetscFunctionReturn(0);
743da3a660dSBarry Smith }
744da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7454a2ae208SSatish Balay #undef __FUNCT__
7464a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
747dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
748da3a660dSBarry Smith {
749416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7502235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7516849ba73SBarry Smith   PetscErrorCode ierr;
75238baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
75338baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
75487828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
755da3a660dSBarry Smith 
7563a40ed3dSBarry Smith   PetscFunctionBegin;
7571ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
759416022c9SBarry Smith   tmp  = a->solve_work;
760da3a660dSBarry Smith 
7612235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7622235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
763da3a660dSBarry Smith 
764da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7652235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
766da3a660dSBarry Smith 
767da3a660dSBarry Smith   /* forward solve the U^T */
768da3a660dSBarry Smith   for (i=0; i<n; i++) {
769010a20caSHong Zhang     v   = aa + diag[i] ;
770010a20caSHong Zhang     vi  = aj + diag[i] + 1;
771f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
772c38d4ed2SBarry Smith     s1  = tmp[i];
773ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
774da3a660dSBarry Smith     while (nz--) {
775010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
776da3a660dSBarry Smith     }
777c38d4ed2SBarry Smith     tmp[i] = s1;
778da3a660dSBarry Smith   }
779da3a660dSBarry Smith 
780da3a660dSBarry Smith   /* backward solve the L^T */
781da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
782010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
783010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
784f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
785f1af5d2fSBarry Smith     s1  = tmp[i];
786da3a660dSBarry Smith     while (nz--) {
787010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
788da3a660dSBarry Smith     }
789da3a660dSBarry Smith   }
790da3a660dSBarry Smith 
791da3a660dSBarry Smith   /* copy tmp into x according to permutation */
792da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
793da3a660dSBarry Smith 
7942235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7952235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7961ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7971ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
798da3a660dSBarry Smith 
799efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr);
8003a40ed3dSBarry Smith   PetscFunctionReturn(0);
801da3a660dSBarry Smith }
802da3a660dSBarry Smith 
8034a2ae208SSatish Balay #undef __FUNCT__
8044a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
805dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
806da3a660dSBarry Smith {
807416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8082235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8096849ba73SBarry Smith   PetscErrorCode ierr;
81038baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
81138baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
81287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8136abc6512SBarry Smith 
8143a40ed3dSBarry Smith   PetscFunctionBegin;
81523bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8166abc6512SBarry Smith 
8171ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
819416022c9SBarry Smith   tmp = a->solve_work;
8206abc6512SBarry Smith 
8212235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8222235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8236abc6512SBarry Smith 
8246abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8252235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8266abc6512SBarry Smith 
8276abc6512SBarry Smith   /* forward solve the U^T */
8286abc6512SBarry Smith   for (i=0; i<n; i++) {
829010a20caSHong Zhang     v   = aa + diag[i] ;
830010a20caSHong Zhang     vi  = aj + diag[i] + 1;
831f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8326abc6512SBarry Smith     tmp[i] *= *v++;
8336abc6512SBarry Smith     while (nz--) {
834010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8356abc6512SBarry Smith     }
8366abc6512SBarry Smith   }
8376abc6512SBarry Smith 
8386abc6512SBarry Smith   /* backward solve the L^T */
8396abc6512SBarry Smith   for (i=n-1; i>=0; i--){
840010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
841010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
842f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8436abc6512SBarry Smith     while (nz--) {
844010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8456abc6512SBarry Smith     }
8466abc6512SBarry Smith   }
8476abc6512SBarry Smith 
8486abc6512SBarry Smith   /* copy tmp into x according to permutation */
8496abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8506abc6512SBarry Smith 
8512235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8522235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8531ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8541ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8556abc6512SBarry Smith 
856efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
858da3a660dSBarry Smith }
8599e25ed09SBarry Smith /* ----------------------------------------------------------------*/
860dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
86108480c60SBarry Smith 
8624a2ae208SSatish Balay #undef __FUNCT__
8634a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
864dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8659e25ed09SBarry Smith {
866416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8679e25ed09SBarry Smith   IS             isicol;
8686849ba73SBarry Smith   PetscErrorCode ierr;
8695a8e39fbSHong Zhang   PetscInt       *r,*ic,n=A->m,*ai=a->i,*aj=a->j;
8705a8e39fbSHong Zhang   PetscInt       *bi,*cols,nnz,*cols_lvl;
8715a8e39fbSHong Zhang   PetscInt       *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
8725a8e39fbSHong Zhang   PetscInt       i,levels,diagonal_fill;
87377c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
874329f5518SBarry Smith   PetscReal      f;
8755a8e39fbSHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
8765a8e39fbSHong Zhang   PetscBT        lnkbt;
8775a8e39fbSHong Zhang   PetscInt       nzi,*bj,**bj_ptr,**bjlvl_ptr;
8785a8e39fbSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
8795a8e39fbSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
8809e25ed09SBarry Smith 
8813a40ed3dSBarry Smith   PetscFunctionBegin;
8826cf997b0SBarry Smith   f             = info->fill;
88338baddfdSBarry Smith   levels        = (PetscInt)info->levels;
88438baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8854c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
88682bf6240SBarry Smith 
88708480c60SBarry Smith   /* special case that simply copies fill pattern */
888be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
889be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
89086bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8912e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
89208480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
89308480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
89408480c60SBarry Smith     if (!b->diag) {
8957c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89608480c60SBarry Smith     }
8977c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
89808480c60SBarry Smith     b->row              = isrow;
89908480c60SBarry Smith     b->col              = iscol;
90082bf6240SBarry Smith     b->icol             = isicol;
90187828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
902f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
903c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
904c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9053a40ed3dSBarry Smith     PetscFunctionReturn(0);
90608480c60SBarry Smith   }
90708480c60SBarry Smith 
908898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
909898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9109e25ed09SBarry Smith 
9119e25ed09SBarry Smith   /* get new row pointers */
9125a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9135a8e39fbSHong Zhang   bi[0] = 0;
9145a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9155a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9165a8e39fbSHong Zhang   bdiag[0]  = 0;
9176cf997b0SBarry Smith 
9185a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9195a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9205a8e39fbSHong Zhang 
9215a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9225a8e39fbSHong Zhang   nlnk = n + 1;
9235a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9245a8e39fbSHong Zhang 
9255a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
9265a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9275a8e39fbSHong Zhang   current_space = free_space;
9285a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9295a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
9305a8e39fbSHong Zhang 
9315a8e39fbSHong Zhang   for (i=0; i<n; i++) {
9325a8e39fbSHong Zhang     nzi = 0;
9336cf997b0SBarry Smith     /* copy current row into linked list */
9345a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
9355a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
9365a8e39fbSHong Zhang     cols = aj + ai[r[i]];
9375a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
9385a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9395a8e39fbSHong Zhang     nzi += nlnk;
9406cf997b0SBarry Smith 
9416cf997b0SBarry Smith     /* make sure diagonal entry is included */
9425a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
9436cf997b0SBarry Smith       fm = n;
9445a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
9455a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
9465a8e39fbSHong Zhang       lnk[fm]    = i;
9475a8e39fbSHong Zhang       lnk_lvl[i] = 0;
9485a8e39fbSHong Zhang       nzi++; dcount++;
9496cf997b0SBarry Smith     }
9506cf997b0SBarry Smith 
9515a8e39fbSHong Zhang     /* add pivot rows into the active row */
9525a8e39fbSHong Zhang     nzbd = 0;
9535a8e39fbSHong Zhang     prow = lnk[n];
9545a8e39fbSHong Zhang     while (prow < i) {
9555a8e39fbSHong Zhang       nnz      = bdiag[prow];
9565a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
9575a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
9585a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
9595a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
9605a8e39fbSHong Zhang       nzi += nlnk;
9615a8e39fbSHong Zhang       prow = lnk[prow];
9625a8e39fbSHong Zhang       nzbd++;
96356beaf04SBarry Smith     }
9645a8e39fbSHong Zhang     bdiag[i] = nzbd;
9655a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
966ecf371e4SBarry Smith 
9675a8e39fbSHong Zhang     /* if free space is not available, make more free space */
9685a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
9695a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
9705a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
9715a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space_lvl);CHKERRQ(ierr);
9725a8e39fbSHong Zhang       reallocs++;
9735a8e39fbSHong Zhang     }
974ecf371e4SBarry Smith 
9755a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
9765a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
9775a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
9785a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
9795a8e39fbSHong Zhang 
9805a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
9815a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
98277431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
9835a8e39fbSHong Zhang     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i);
9846cf997b0SBarry Smith     }
9855a8e39fbSHong Zhang 
9865a8e39fbSHong Zhang     current_space->array           += nzi;
9875a8e39fbSHong Zhang     current_space->local_used      += nzi;
9885a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
9895a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
9905a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
9915a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
9929e25ed09SBarry Smith   }
9935a8e39fbSHong Zhang 
994898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
995898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
9965a8e39fbSHong Zhang 
9975a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
9985a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
9995a8e39fbSHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
10005a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
10015a8e39fbSHong Zhang   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
10025a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10039e25ed09SBarry Smith 
1004f64065bbSBarry Smith   {
10055a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1006418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1007c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1008c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1009b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1010335d9088SBarry Smith     if (diagonal_fill) {
101177431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1012335d9088SBarry Smith     }
1013f64065bbSBarry Smith   }
1014bbb6d6a8SBarry Smith 
10159e25ed09SBarry Smith   /* put together the new matrix */
1016f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1017f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1018f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
101952e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1020416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1021606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10227c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10235a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
10249e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1025606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1026606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1027b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10285a8e39fbSHong Zhang   b->j          = bj;
10295a8e39fbSHong Zhang   b->i          = bi;
10305a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
10315a8e39fbSHong Zhang   b->diag       = bdiag;
1032416022c9SBarry Smith   b->ilen       = 0;
1033416022c9SBarry Smith   b->imax       = 0;
1034416022c9SBarry Smith   b->row        = isrow;
1035416022c9SBarry Smith   b->col        = iscol;
1036c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1037c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
103882bf6240SBarry Smith   b->icol       = isicol;
103987828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1040416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
10415a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
10425a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
10435a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
10449e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1045418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1046ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
10475a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
104871c2f376SKris Buschelman 
1049*fd62cfe2SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,info,isrow,iscol,fact);CHKERRQ(ierr);
105071c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
10514846f1f5SKris Buschelman 
10523a40ed3dSBarry Smith   PetscFunctionReturn(0);
10539e25ed09SBarry Smith }
1054d5d45c9bSBarry Smith 
10553c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1056a6175056SHong Zhang #undef __FUNCT__
1057a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1058af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1059a6175056SHong Zhang {
10602ed133dbSHong Zhang   Mat            C = *B;
10610968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10622ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
10632ed133dbSHong Zhang   IS             ip=b->row;
10642ed133dbSHong Zhang   PetscErrorCode ierr;
10652ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
10662ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
10672ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1068622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1069540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1070540703acSHong Zhang   PetscTruth     shiftpd;
1071540703acSHong Zhang   ChShift_Ctx    sctx;
1072540703acSHong Zhang   PetscInt       newshift;
1073d5d45c9bSBarry Smith 
1074a6175056SHong Zhang   PetscFunctionBegin;
1075540703acSHong Zhang   shiftnz   = info->shiftnz;
1076540703acSHong Zhang   shiftpd   = info->shiftpd;
1077ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1078ee45ca4aSHong Zhang 
10792ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
10802ed133dbSHong Zhang 
10812ed133dbSHong Zhang   /* initialization */
10822ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
10832ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
10842ed133dbSHong Zhang   jl   = il + mbs;
10852ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
10862ed133dbSHong Zhang 
1087540703acSHong Zhang   sctx.shift_amount = 0;
1088540703acSHong Zhang   sctx.nshift       = 0;
10892ed133dbSHong Zhang   do {
1090540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
10912ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
10922ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1093a6175056SHong Zhang     }
10942ed133dbSHong Zhang 
10952ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1096c05c3958SHong Zhang       bval = ba + bi[k];
10972ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
10982ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
10992ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11002ed133dbSHong Zhang         col = rip[aj[j]];
11012ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11022ed133dbSHong Zhang           rtmp[col] = aa[j];
1103c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11042ed133dbSHong Zhang         }
11052ed133dbSHong Zhang       }
110639e3d630SHong Zhang       /* shift the diagonal of the matrix */
1107540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11082ed133dbSHong Zhang 
11092ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11102ed133dbSHong Zhang       dk = rtmp[k];
11112ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11122ed133dbSHong Zhang 
11132ed133dbSHong Zhang       while (i < k){
11142ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11152ed133dbSHong Zhang 
11162ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11172ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11182ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11192ed133dbSHong Zhang         dk += uikdi*ba[ili];
11202ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11212ed133dbSHong Zhang 
11222ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11232ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11242ed133dbSHong Zhang         if (jmin < jmax){
11252ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11262ed133dbSHong Zhang           /* update il and jl for row i */
11272ed133dbSHong Zhang           il[i] = jmin;
11282ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11292ed133dbSHong Zhang         }
11302ed133dbSHong Zhang         i = nexti;
11312ed133dbSHong Zhang       }
11322ed133dbSHong Zhang 
11330697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11340697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11350697b6caSHong Zhang       rs   = 0.0;
11363c9e8598SHong Zhang       jmin = bi[k]+1;
11373c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11383c9e8598SHong Zhang       if (nz){
11393c9e8598SHong Zhang         bcol = bj + jmin;
11403c9e8598SHong Zhang         while (nz--){
114189f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
114289f845c8SHong Zhang           bcol++;
11433c9e8598SHong Zhang         }
11443c9e8598SHong Zhang       }
1145540703acSHong Zhang 
1146540703acSHong Zhang       sctx.rs = rs;
1147540703acSHong Zhang       sctx.pv = dk;
11483e8c821fSHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1149540703acSHong Zhang       if (newshift == 1){
1150540703acSHong Zhang         break;    /* sctx.shift_amount is updated */
1151540703acSHong Zhang       } else if (newshift == -1){
11520697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11533c9e8598SHong Zhang       }
11543c9e8598SHong Zhang 
11553c9e8598SHong Zhang       /* copy data into U(k,:) */
115639e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
115739e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
115839e3d630SHong Zhang       if (jmin < jmax) {
115939e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
116039e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
11613c9e8598SHong Zhang         }
116239e3d630SHong Zhang         /* add the k-th row into il and jl */
11633c9e8598SHong Zhang         il[k] = jmin;
11643c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
11653c9e8598SHong Zhang       }
11663c9e8598SHong Zhang     }
1167540703acSHong Zhang   } while (sctx.chshift);
11683c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
11693c9e8598SHong Zhang 
117039e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
11713c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
11723c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
11733c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1174efee365bSSatish Balay   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
1175540703acSHong Zhang   if (sctx.nshift){
1176540703acSHong Zhang     if (shiftnz) {
1177540703acSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
1178540703acSHong Zhang     } else if (shiftpd) {
1179540703acSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
11800697b6caSHong Zhang     }
11813c9e8598SHong Zhang   }
11823c9e8598SHong Zhang   PetscFunctionReturn(0);
11833c9e8598SHong Zhang }
1184a6175056SHong Zhang 
1185a6175056SHong Zhang #undef __FUNCT__
1186a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1187dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1188a6175056SHong Zhang {
11890968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1190ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1191ed59401aSHong Zhang   Mat            B;
1192dfbe8321SBarry Smith   PetscErrorCode ierr;
1193653a6975SHong Zhang   PetscTruth     perm_identity;
1194622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1195ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1196391eac55SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
11975a8e39fbSHong Zhang   PetscInt       ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1198ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
11997a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
12007a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12017a48dd6fSHong Zhang   PetscBT        lnkbt;
1202a6175056SHong Zhang 
1203a6175056SHong Zhang   PetscFunctionBegin;
1204653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12057a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12067a48dd6fSHong Zhang 
120739e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
120839e3d630SHong Zhang   ui[0] = 0;
120939e3d630SHong Zhang 
1210622e2028SHong Zhang   /* special case that simply copies fill pattern */
1211622e2028SHong Zhang   if (!levels && perm_identity) {
1212622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1213ed59401aSHong Zhang     for (i=0; i<am; i++) {
121439e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1215ed59401aSHong Zhang     }
121639e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
121739e3d630SHong Zhang     cols = uj;
1218ed59401aSHong Zhang     for (i=0; i<am; i++) {
1219ed59401aSHong Zhang       aj    = a->j + a->diag[i];
122039e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
122139e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1222ed59401aSHong Zhang     }
122339e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12247a48dd6fSHong Zhang     /* initialization */
12255a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12267a48dd6fSHong Zhang 
12277a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12287a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12291393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12307a48dd6fSHong Zhang     il         = jl + am;
12317a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12327a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12337a48dd6fSHong Zhang     for (i=0; i<am; i++){
12347a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12357a48dd6fSHong Zhang     }
12367a48dd6fSHong Zhang 
12377a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12387a48dd6fSHong Zhang     nlnk = am + 1;
12392ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12407a48dd6fSHong Zhang 
12417a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12427a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12437a48dd6fSHong Zhang     current_space = free_space;
12447a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12457a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12467a48dd6fSHong Zhang 
12477a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12487a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12497a48dd6fSHong Zhang       nzk   = 0;
1250622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1251622e2028SHong Zhang       ncols_upper = 0;
1252622e2028SHong Zhang       for (j=0; j<ncols; j++){
12535a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
12545a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
12555a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1256622e2028SHong Zhang           ncols_upper++;
1257622e2028SHong Zhang         }
1258622e2028SHong Zhang       }
12595a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12607a48dd6fSHong Zhang       nzk += nlnk;
12617a48dd6fSHong Zhang 
12627a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
12637a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
12647a48dd6fSHong Zhang 
12657a48dd6fSHong Zhang       while (prow < k){
12667a48dd6fSHong Zhang         nextprow = jl[prow];
12677a48dd6fSHong Zhang 
12687a48dd6fSHong Zhang         /* merge prow into k-th row */
12697a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
12707a48dd6fSHong Zhang         jmax = ui[prow+1];
12717a48dd6fSHong Zhang         ncols = jmax-jmin;
1272ed59401aSHong Zhang         i     = jmin - ui[prow];
1273ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
12745a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
12755a8e39fbSHong Zhang         j     = *(uj - 1);
12765a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
12777a48dd6fSHong Zhang         nzk += nlnk;
12787a48dd6fSHong Zhang 
12797a48dd6fSHong Zhang         /* update il and jl for prow */
12807a48dd6fSHong Zhang         if (jmin < jmax){
12817a48dd6fSHong Zhang           il[prow] = jmin;
12827a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
12837a48dd6fSHong Zhang         }
12847a48dd6fSHong Zhang         prow = nextprow;
12857a48dd6fSHong Zhang       }
12867a48dd6fSHong Zhang 
12877a48dd6fSHong Zhang       /* if free space is not available, make more free space */
12887a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
12897a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
12907a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
12917a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
12927a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
12937a48dd6fSHong Zhang         reallocs++;
12947a48dd6fSHong Zhang       }
12957a48dd6fSHong Zhang 
12967a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
12972ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
12987a48dd6fSHong Zhang 
12997a48dd6fSHong Zhang       /* add the k-th row into il and jl */
130039e3d630SHong Zhang       if (nzk > 1){
13017a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13027a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13037a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13047a48dd6fSHong Zhang       }
13057a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13067a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13077a48dd6fSHong Zhang 
13087a48dd6fSHong Zhang       current_space->array           += nzk;
13097a48dd6fSHong Zhang       current_space->local_used      += nzk;
13107a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13117a48dd6fSHong Zhang 
13127a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13137a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13147a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13157a48dd6fSHong Zhang 
13167a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13177a48dd6fSHong Zhang     }
13187a48dd6fSHong Zhang 
13197a48dd6fSHong Zhang     if (ai[am] != 0) {
132039e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
13217a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
13227a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
13237a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
13247a48dd6fSHong Zhang     } else {
1325ed59401aSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
13267a48dd6fSHong Zhang     }
13277a48dd6fSHong Zhang 
13287a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13297a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13305a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13317a48dd6fSHong Zhang 
13327a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13337a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13347a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13352ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13367a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13377a48dd6fSHong Zhang 
133839e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
133939e3d630SHong Zhang 
13407a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13417a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1342ed59401aSHong Zhang   B = *fact;
1343ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1344ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
13457a48dd6fSHong Zhang 
1346ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13477a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
13487a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13497a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
13507a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
13517a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
13527a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13537a48dd6fSHong Zhang   b->j    = uj;
13547a48dd6fSHong Zhang   b->i    = ui;
13557a48dd6fSHong Zhang   b->diag = 0;
13567a48dd6fSHong Zhang   b->ilen = 0;
13577a48dd6fSHong Zhang   b->imax = 0;
13587a48dd6fSHong Zhang   b->row  = perm;
13597a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
13607a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13617a48dd6fSHong Zhang   b->icol = perm;
13627a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13637a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
136452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13657a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
13667a48dd6fSHong Zhang 
1367ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1368ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1369ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
13707a48dd6fSHong Zhang   if (ai[am] != 0) {
1371ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
13727a48dd6fSHong Zhang   } else {
1373ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
13747a48dd6fSHong Zhang   }
137539e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1376622e2028SHong Zhang   if (perm_identity){
1377ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1378ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1379622e2028SHong Zhang   }
1380a6175056SHong Zhang   PetscFunctionReturn(0);
1381a6175056SHong Zhang }
1382d5d45c9bSBarry Smith 
1383f76d2b81SHong Zhang #undef __FUNCT__
1384f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1385dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1386f76d2b81SHong Zhang {
1387f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
138810c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
13892ed133dbSHong Zhang   Mat            B;
1390dfbe8321SBarry Smith   PetscErrorCode ierr;
1391f76d2b81SHong Zhang   PetscTruth     perm_identity;
139210c27e3dSHong Zhang   PetscReal      fill = info->fill;
13931393bc97SHong Zhang   PetscInt       *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
139410c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
13952ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
139610c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
139710c27e3dSHong Zhang   PetscBT        lnkbt;
13982ed133dbSHong Zhang   IS             iperm;
1399f76d2b81SHong Zhang 
1400f76d2b81SHong Zhang   PetscFunctionBegin;
14012ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1402f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14032ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14042ed133dbSHong Zhang 
1405f76d2b81SHong Zhang   if (!perm_identity){
14062ed133dbSHong Zhang     /* check if perm is symmetric! */
14072ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14082ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14091393bc97SHong Zhang     for (i=0; i<am; i++) {
14102ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14112ed133dbSHong Zhang     }
14122ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14132ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1414f76d2b81SHong Zhang   }
141510c27e3dSHong Zhang 
141610c27e3dSHong Zhang   /* initialization */
14171393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
141810c27e3dSHong Zhang   ui[0] = 0;
141910c27e3dSHong Zhang 
142010c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14211393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14221393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14231393bc97SHong Zhang   il     = jl + am;
14241393bc97SHong Zhang   cols   = il + am;
14251393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14261393bc97SHong Zhang   for (i=0; i<am; i++){
14271393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1428f76d2b81SHong Zhang   }
142910c27e3dSHong Zhang 
143010c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14311393bc97SHong Zhang   nlnk = am + 1;
14321393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
143310c27e3dSHong Zhang 
14341393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
14351393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
143610c27e3dSHong Zhang   current_space = free_space;
143710c27e3dSHong Zhang 
14381393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
143910c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
144010c27e3dSHong Zhang     nzk   = 0;
14412ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14422ed133dbSHong Zhang     ncols_upper = 0;
14432ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1444622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14452ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14462ed133dbSHong Zhang         cols[ncols_upper] = i;
14472ed133dbSHong Zhang         ncols_upper++;
14482ed133dbSHong Zhang       }
14492ed133dbSHong Zhang     }
14501393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
145110c27e3dSHong Zhang     nzk += nlnk;
145210c27e3dSHong Zhang 
145310c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
145410c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
145510c27e3dSHong Zhang 
145610c27e3dSHong Zhang     while (prow < k){
145710c27e3dSHong Zhang       nextprow = jl[prow];
145810c27e3dSHong Zhang       /* merge prow into k-th row */
14591393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
146010c27e3dSHong Zhang       jmax = ui[prow+1];
146110c27e3dSHong Zhang       ncols = jmax-jmin;
14621393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
14635a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
146410c27e3dSHong Zhang       nzk += nlnk;
146510c27e3dSHong Zhang 
146610c27e3dSHong Zhang       /* update il and jl for prow */
146710c27e3dSHong Zhang       if (jmin < jmax){
146810c27e3dSHong Zhang         il[prow] = jmin;
14692ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
147010c27e3dSHong Zhang       }
147110c27e3dSHong Zhang       prow = nextprow;
147210c27e3dSHong Zhang     }
147310c27e3dSHong Zhang 
147410c27e3dSHong Zhang     /* if free space is not available, make more free space */
147510c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
14761393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
147710c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
147810c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
147910c27e3dSHong Zhang       reallocs++;
148010c27e3dSHong Zhang     }
148110c27e3dSHong Zhang 
148210c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
14831393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
148410c27e3dSHong Zhang 
148510c27e3dSHong Zhang     /* add the k-th row into il and jl */
148610c27e3dSHong Zhang     if (nzk-1 > 0){
14871393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
148810c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
148910c27e3dSHong Zhang       il[k] = ui[k] + 1;
149010c27e3dSHong Zhang     }
14912ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
149210c27e3dSHong Zhang     current_space->array           += nzk;
149310c27e3dSHong Zhang     current_space->local_used      += nzk;
149410c27e3dSHong Zhang     current_space->local_remaining -= nzk;
149510c27e3dSHong Zhang 
149610c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
149710c27e3dSHong Zhang   }
149810c27e3dSHong Zhang 
14991393bc97SHong Zhang   if (ai[am] != 0) {
15001393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
150178910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
150278910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
150378910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
150410c27e3dSHong Zhang   } else {
150578910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
150610c27e3dSHong Zhang   }
150710c27e3dSHong Zhang 
150810c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
150910c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
151010c27e3dSHong Zhang 
151110c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15121393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
151310c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
151410c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
151510c27e3dSHong Zhang 
151610c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15171393bc97SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
15182ed133dbSHong Zhang   B    = *fact;
15192ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
15202ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
152110c27e3dSHong Zhang 
15222ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
152310c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
152410c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
152510c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
152610c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
152710c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15281393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
152910c27e3dSHong Zhang   b->j    = uj;
153010c27e3dSHong Zhang   b->i    = ui;
153110c27e3dSHong Zhang   b->diag = 0;
153210c27e3dSHong Zhang   b->ilen = 0;
153310c27e3dSHong Zhang   b->imax = 0;
153410c27e3dSHong Zhang   b->row  = perm;
153510c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
153610c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
153710c27e3dSHong Zhang   b->icol = perm;
153810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15391393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15401393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
15411393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
154210c27e3dSHong Zhang 
15432ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15442ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15452ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15461393bc97SHong Zhang   if (ai[am] != 0) {
15471393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
154810c27e3dSHong Zhang   } else {
15492ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
155010c27e3dSHong Zhang   }
155139e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15522ed133dbSHong Zhang   if (perm_identity){
155310c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
155410c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15552ed133dbSHong Zhang   }
1556f76d2b81SHong Zhang   PetscFunctionReturn(0);
1557f76d2b81SHong Zhang }
1558