xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 0c451bc42d7bade440162d987bcef1c8f26e50ee)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2386bacbd2SBarry Smith 
2438baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
259cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2638baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
639cc1f4e3SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->m;
6438baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6538baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6638baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6787828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6854a8161fSBarry Smith   PetscReal      af;
6986bacbd2SBarry Smith 
7086bacbd2SBarry Smith   PetscFunctionBegin;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7338baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7486bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7533259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7638baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7738baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7886bacbd2SBarry Smith 
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith   /* ------------------------------------------------------------
8186bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8286bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8386bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8486bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8586bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8686bacbd2SBarry Smith      ------------------------------------------------------------
8786bacbd2SBarry Smith   */
8886bacbd2SBarry Smith 
8986bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9086bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9186bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9286bacbd2SBarry Smith   reorder = PetscNot(reorder);
9386bacbd2SBarry Smith 
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith   /* storage for ilu factor */
9638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9738baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
9887828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
9938baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10086bacbd2SBarry Smith 
10186bacbd2SBarry Smith   /* ------------------------------------------------------------
10286bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10386bacbd2SBarry Smith      ------------------------------------------------------------
10486bacbd2SBarry Smith   */
105b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
106b19713cbSBarry Smith     old_j[i]++;
10786bacbd2SBarry Smith   }
108b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
109b19713cbSBarry Smith     old_i[i]++;
110b19713cbSBarry Smith   };
111010a20caSHong Zhang 
11286bacbd2SBarry Smith 
11386bacbd2SBarry Smith   if (reorder) {
114c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
115c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
116c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
117c0c2c81eSBarry Smith       r[i]  = r[i]+1;
118c0c2c81eSBarry Smith       c[i]  = c[i]+1;
119c0c2c81eSBarry Smith     }
12038baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12138baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12287828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1235bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
124c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
125c0c2c81eSBarry Smith       r[i]  = r[i]-1;
126c0c2c81eSBarry Smith       c[i]  = c[i]-1;
127c0c2c81eSBarry Smith     }
128c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
129c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
130b19713cbSBarry Smith     o_a = old_a2;
131b19713cbSBarry Smith     o_j = old_j2;
132b19713cbSBarry Smith     o_i = old_i2;
133b19713cbSBarry Smith   } else {
134b19713cbSBarry Smith     o_a = old_a;
135b19713cbSBarry Smith     o_j = old_j;
136b19713cbSBarry Smith     o_i = old_i;
13786bacbd2SBarry Smith   }
13886bacbd2SBarry Smith 
13986bacbd2SBarry Smith   /* ------------------------------------------------------------
14086bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14186bacbd2SBarry Smith      ------------------------------------------------------------
14286bacbd2SBarry Smith   */
14386bacbd2SBarry Smith 
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14538baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14687828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14786bacbd2SBarry Smith 
14854a8161fSBarry 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);
1496849ba73SBarry Smith   if (sierr) {
1506849ba73SBarry Smith     switch (sierr) {
15177431f27SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
15277431f27SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
153e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
154e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15577431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15677431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15786bacbd2SBarry Smith     }
15886bacbd2SBarry Smith   }
15986bacbd2SBarry Smith 
16086bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16186bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16286bacbd2SBarry Smith 
16386bacbd2SBarry Smith   /* ------------------------------------------------------------
16486bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16586bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16686bacbd2SBarry Smith      ------------------------------------------------------------
16786bacbd2SBarry Smith   */
16886bacbd2SBarry Smith 
16987828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
17038baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17186bacbd2SBarry Smith 
1725bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17386bacbd2SBarry Smith 
17486bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17586bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17686bacbd2SBarry Smith 
17786bacbd2SBarry Smith   if (reorder) {
17886bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18086bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
181313ae024SBarry Smith   } else {
182b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
183f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
184b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
185b19713cbSBarry Smith     }
186b19713cbSBarry Smith   }
187b19713cbSBarry Smith 
188b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
190b19713cbSBarry Smith     old_i[i]--;
191b19713cbSBarry Smith   }
192b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
193b19713cbSBarry Smith     old_j[i]--;
194b19713cbSBarry Smith   }
19586bacbd2SBarry Smith 
196b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
198b19713cbSBarry Smith     new_i[i]--;
199b19713cbSBarry Smith   }
200b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
201b19713cbSBarry Smith     new_j[i]--;
202b19713cbSBarry Smith   }
20386bacbd2SBarry Smith 
20486bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20586bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2064c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2074c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
208c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20986bacbd2SBarry Smith   for(i=0; i<n; i++) {
21086bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21186bacbd2SBarry Smith   };
21286bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21307d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21486bacbd2SBarry Smith 
215c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
216c0c2c81eSBarry Smith 
217329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21807d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21986bacbd2SBarry Smith 
22086bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22186bacbd2SBarry Smith 
222f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
223f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
224ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
22586bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22686bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22786bacbd2SBarry Smith 
22886bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
229a96a251dSBarry Smith   b->freedata      = PETSC_TRUE;
23086bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23107d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23286bacbd2SBarry Smith   b->a             = new_a;
23386bacbd2SBarry Smith   b->j             = new_j;
23486bacbd2SBarry Smith   b->i             = new_i;
23586bacbd2SBarry Smith   b->ilen          = 0;
23686bacbd2SBarry Smith   b->imax          = 0;
2371f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
238313ae024SBarry Smith   b->row           = isirow;
23986bacbd2SBarry Smith   b->col           = iscolf;
24087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24186bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24286bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24386bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24486bacbd2SBarry Smith 
24586bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24686bacbd2SBarry Smith 
247431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
24863ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af));CHKERRQ(ierr);
24963ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af));CHKERRQ(ierr);
25063ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af));CHKERRQ(ierr);
25163ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:for best performance.\n"));CHKERRQ(ierr);
252431e710aSBarry Smith 
2537529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25471c2f376SKris Buschelman 
25586bacbd2SBarry Smith   PetscFunctionReturn(0);
25660ec86bdSBarry Smith #endif
25786bacbd2SBarry Smith }
25886bacbd2SBarry Smith 
2594a2ae208SSatish Balay #undef __FUNCT__
260b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
261dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
262289bc588SBarry Smith {
263416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
264289bc588SBarry Smith   IS             isicol;
2656849ba73SBarry Smith   PetscErrorCode ierr;
26638baddfdSBarry Smith   PetscInt       *r,*ic,i,n=A->m,*ai=a->i,*aj=a->j;
2671393bc97SHong Zhang   PetscInt       *bi,*bj,*ajtmp;
2681393bc97SHong Zhang   PetscInt       *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2699e046878SBarry Smith   PetscReal      f;
2701393bc97SHong Zhang   PetscInt       nlnk,*lnk,k,*cols,**bi_ptr;
2711393bc97SHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
2721393bc97SHong Zhang   PetscBT        lnkbt;
273289bc588SBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
27529bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2764c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
277ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
278ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
279289bc588SBarry Smith 
280289bc588SBarry Smith   /* get new row pointers */
2811393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2821393bc97SHong Zhang   bi[0] = 0;
2831393bc97SHong Zhang 
2841393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2851393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2861393bc97SHong Zhang   bdiag[0] = 0;
2871393bc97SHong Zhang 
2881393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2891393bc97SHong Zhang   nlnk = n + 1;
2901393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2911393bc97SHong Zhang 
2921393bc97SHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&cols);CHKERRQ(ierr);
2931393bc97SHong Zhang   im     = cols + n;
2941393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2951393bc97SHong Zhang 
2961393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
297d3d32019SBarry Smith   f = info->fill;
2981393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
2991393bc97SHong Zhang   current_space = free_space;
300289bc588SBarry Smith 
301289bc588SBarry Smith   for (i=0; i<n; i++) {
3021393bc97SHong Zhang     /* copy previous fill into linked list */
303289bc588SBarry Smith     nzi = 0;
3041393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3051393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3061393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
3071393bc97SHong Zhang     for (k=0; k<nnz; k++) cols[k] = ic[*(ajtmp+k)]; /* note: cols is not sorted when iscol!=indentity */
3081393bc97SHong Zhang     ierr = PetscLLAdd(nnz,cols,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3091393bc97SHong Zhang     nzi += nlnk;
3101393bc97SHong Zhang 
3111393bc97SHong Zhang     /* add pivot rows into linked list */
3121393bc97SHong Zhang     row = lnk[n];
3131393bc97SHong Zhang     while (row < i) {
3141393bc97SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1;
3151393bc97SHong Zhang       ajtmp   = bi_ptr[row] + nzbd;
3161393bc97SHong Zhang       nnz     = im[row] - nzbd; /* num of columns with row<indices<=i */
3171393bc97SHong Zhang       im[row] = nzbd;
3181393bc97SHong Zhang       ierr = PetscLLAddSortedLU(nnz,ajtmp,row,nlnk,lnk,lnkbt,i,nzbd);CHKERRQ(ierr);
3191393bc97SHong Zhang       nzi     += nlnk;
3201393bc97SHong Zhang       im[row] += nzbd;  /* update im[row]: num of cols with index<=i */
3211393bc97SHong Zhang 
3221393bc97SHong Zhang       row = lnk[row];
323289bc588SBarry Smith     }
3241393bc97SHong Zhang 
3251393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3261393bc97SHong Zhang     im[i]   = nzi;
3271393bc97SHong Zhang 
3281393bc97SHong Zhang     /* mark bdiag */
3291393bc97SHong Zhang     nzbd = 0;
3301393bc97SHong Zhang     nnz  = nzi;
3311393bc97SHong Zhang     k    = lnk[n];
3321393bc97SHong Zhang     while (nnz-- && k < i){
3331393bc97SHong Zhang       nzbd++;
3341393bc97SHong Zhang       k = lnk[k];
3351393bc97SHong Zhang     }
3361393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3371393bc97SHong Zhang 
3381393bc97SHong Zhang     /* if free space is not available, make more free space */
3391393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3401393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
3411393bc97SHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
3421393bc97SHong Zhang       reallocs++;
3431393bc97SHong Zhang     }
3441393bc97SHong Zhang 
3451393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3461393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3471393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3481393bc97SHong Zhang     current_space->array           += nzi;
3491393bc97SHong Zhang     current_space->local_used      += nzi;
3501393bc97SHong Zhang     current_space->local_remaining -= nzi;
351289bc588SBarry Smith   }
35263ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
353464e8d44SSatish Balay   if (ai[n] != 0) {
3541393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
35563ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr);
356*0c451bc4SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %G or use \n",af));CHKERRQ(ierr);
35763ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af));CHKERRQ(ierr);
35863ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr);
35905bf559cSSatish Balay   } else {
36063ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"));CHKERRQ(ierr);
36105bf559cSSatish Balay   }
36263ba0a88SBarry Smith #endif
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);
376ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
37752e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
378416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
379a96a251dSBarry Smith   b->freedata     = PETSC_TRUE;
3807c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3811393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3821393bc97SHong Zhang   b->j          = bj;
3831393bc97SHong Zhang   b->i          = bi;
3841393bc97SHong Zhang   b->diag       = bdiag;
385416022c9SBarry Smith   b->ilen       = 0;
386416022c9SBarry Smith   b->imax       = 0;
387416022c9SBarry Smith   b->row        = isrow;
388416022c9SBarry Smith   b->col        = iscol;
389c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
390c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
39182bf6240SBarry Smith   b->icol       = isicol;
39287828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3931393bc97SHong Zhang 
3941393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3951393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3961393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
3977fd98bd6SLois Curfman McInnes 
398329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
399418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
400ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
401703deb49SSatish Balay 
402e93ccc4dSBarry Smith   if (ai[n] != 0) {
4031393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
404eea2913aSSatish Balay   } else {
405eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
406eea2913aSSatish Balay   }
4074846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
40871c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4093a40ed3dSBarry Smith   PetscFunctionReturn(0);
410289bc588SBarry Smith }
4111393bc97SHong Zhang 
4120a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4134a2ae208SSatish Balay #undef __FUNCT__
4144a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
415af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
416289bc588SBarry Smith {
417416022c9SBarry Smith   Mat            C=*B;
418416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
41982bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4206849ba73SBarry Smith   PetscErrorCode ierr;
421b3bf805bSHong Zhang   PetscInt       *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j;
422b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
423d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
42487828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4256a7c8fc2SHong Zhang   PetscScalar    d;
4266a7c8fc2SHong Zhang   PetscReal      rs;
427b3bf805bSHong Zhang   LUShift_Ctx    sctx;
428d2276718SHong Zhang   PetscInt       newshift;
429289bc588SBarry Smith 
4303a40ed3dSBarry Smith   PetscFunctionBegin;
43178b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
43278b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
43387828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
43487828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
435010a20caSHong Zhang   rtmps = rtmp; ics = ic;
436289bc588SBarry Smith 
4376cc28720Svictorle   if (!a->diag) {
4380cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4390cf777b8SBarry Smith   }
440118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
441118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4421a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
44338baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
444118f5789SBarry Smith     sctx.shift_top = 0;
4456cc28720Svictorle     for (i=0; i<n; i++) {
4469f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4479f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4481a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
449010a20caSHong Zhang       v  = a->a+aai[i];
450e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
451e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4521a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4531a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4546cc28720Svictorle     }
455118f5789SBarry Smith     if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12;
456118f5789SBarry Smith     sctx.shift_top    *= 1.1;
457d2276718SHong Zhang     sctx.nshift_max   = 5;
458d2276718SHong Zhang     sctx.shift_lo     = 0.;
459d2276718SHong Zhang     sctx.shift_hi     = 1.;
460a0a356efSHong Zhang   }
461a0a356efSHong Zhang 
462a0a356efSHong Zhang   sctx.shift_amount = 0;
463a0a356efSHong Zhang   sctx.nshift       = 0;
464085256b3SBarry Smith   do {
465d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
466289bc588SBarry Smith     for (i=0; i<n; i++){
467b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
468b3bf805bSHong Zhang       bjtmp = bj + bi[i];
469b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
470289bc588SBarry Smith 
471289bc588SBarry Smith       /* load in initial (unfactored row) */
472416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
473b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
474010a20caSHong Zhang       v     = a->a + a->i[r[i]];
475085256b3SBarry Smith       for (j=0; j<nz; j++) {
476b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
477085256b3SBarry Smith       }
478d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
479289bc588SBarry Smith 
480b3bf805bSHong Zhang       row = *bjtmp++;
481f2582111SSatish Balay       while  (row < i) {
4828c37ef55SBarry Smith         pc = rtmp + row;
4838c37ef55SBarry Smith         if (*pc != 0.0) {
484010a20caSHong Zhang           pv         = b->a + diag_offset[row];
485010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
48635aab85fSBarry Smith           multiplier = *pc / *pv++;
4878c37ef55SBarry Smith           *pc        = multiplier;
488b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
489f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
490efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
491289bc588SBarry Smith         }
492b3bf805bSHong Zhang         row = *bjtmp++;
493289bc588SBarry Smith       }
494416022c9SBarry Smith       /* finished row so stick it into b->a */
495b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
496b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
497b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
498b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
499d3d32019SBarry Smith       rs   = 0.0;
500d3d32019SBarry Smith       for (j=0; j<nz; j++) {
501d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
502d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
503d3d32019SBarry Smith       }
504d2276718SHong Zhang 
505b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
506d2276718SHong Zhang       sctx.rs  = rs;
507d2276718SHong Zhang       sctx.pv  = pv[diag];
5083e8c821fSHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
509d2276718SHong Zhang       if (newshift == 1){
510b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
511d2276718SHong Zhang       } else if (newshift == -1){
512d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
513d708749eSBarry Smith       }
51435aab85fSBarry Smith     }
515d2276718SHong Zhang 
516118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5176cc28720Svictorle       /*
5186c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5196cc28720Svictorle        * then try lower shift
5206cc28720Svictorle        */
521d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
522d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
523d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
524d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
525d2276718SHong Zhang       sctx.nshift++;
5266cc28720Svictorle     }
527d2276718SHong Zhang   } while (sctx.lushift);
5280f11f9aeSSatish Balay 
52935aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
53035aab85fSBarry Smith   for (i=0; i<n; i++) {
531010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
53235aab85fSBarry Smith   }
53335aab85fSBarry Smith 
534606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
53578b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
53678b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
537416022c9SBarry Smith   C->factor = FACTOR_LU;
5387dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
539c456f294SBarry Smith   C->assembled = PETSC_TRUE;
540efee365bSSatish Balay   ierr = PetscLogFlops(C->n);CHKERRQ(ierr);
541d2276718SHong Zhang   if (sctx.nshift){
542118f5789SBarry Smith     if (info->shiftnz) {
54363ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
544118f5789SBarry Smith     } else if (info->shiftpd) {
54563ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top));CHKERRQ(ierr);
5466cc28720Svictorle     }
5470697b6caSHong Zhang   }
5483a40ed3dSBarry Smith   PetscFunctionReturn(0);
549289bc588SBarry Smith }
5500889b5dcSvictorle 
5510889b5dcSvictorle #undef __FUNCT__
5520889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
553dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5540889b5dcSvictorle {
5550889b5dcSvictorle   PetscFunctionBegin;
5560889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5570889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5580889b5dcSvictorle   PetscFunctionReturn(0);
5590889b5dcSvictorle }
5600889b5dcSvictorle 
5610889b5dcSvictorle 
5620a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5634a2ae208SSatish Balay #undef __FUNCT__
5644a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
565dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
566da3a660dSBarry Smith {
567dfbe8321SBarry Smith   PetscErrorCode ierr;
568416022c9SBarry Smith   Mat            C;
569416022c9SBarry Smith 
5703a40ed3dSBarry Smith   PetscFunctionBegin;
5719e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
572af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
573273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
57452e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5753a40ed3dSBarry Smith   PetscFunctionReturn(0);
576da3a660dSBarry Smith }
5770a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5784a2ae208SSatish Balay #undef __FUNCT__
5794a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
580dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5818c37ef55SBarry Smith {
582416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
583416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
5846849ba73SBarry Smith   PetscErrorCode ierr;
58538baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
58638baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
58787828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
5888c37ef55SBarry Smith 
5893a40ed3dSBarry Smith   PetscFunctionBegin;
5903a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
59191bf9adeSBarry Smith 
5921ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
5931ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
594416022c9SBarry Smith   tmp  = a->solve_work;
5958c37ef55SBarry Smith 
5963b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
5973b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
5988c37ef55SBarry Smith 
5998c37ef55SBarry Smith   /* forward solve the lower triangular */
6008c37ef55SBarry Smith   tmp[0] = b[*r++];
601010a20caSHong Zhang   tmps   = tmp;
6028c37ef55SBarry Smith   for (i=1; i<n; i++) {
603010a20caSHong Zhang     v   = aa + ai[i] ;
604010a20caSHong Zhang     vi  = aj + ai[i] ;
60553e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
60653e7425aSBarry Smith     sum = b[*r++];
607ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6088c37ef55SBarry Smith     tmp[i] = sum;
6098c37ef55SBarry Smith   }
6108c37ef55SBarry Smith 
6118c37ef55SBarry Smith   /* backward solve the upper triangular */
6128c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
613010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
614010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
615416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6168c37ef55SBarry Smith     sum = tmp[i];
617ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
618010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6198c37ef55SBarry Smith   }
6208c37ef55SBarry Smith 
6213b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6223b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6231ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6241ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
625efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
6263a40ed3dSBarry Smith   PetscFunctionReturn(0);
6278c37ef55SBarry Smith }
628026e39d0SSatish Balay 
629930ae53cSBarry Smith /* ----------------------------------------------------------- */
6304a2ae208SSatish Balay #undef __FUNCT__
6314a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
632dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
633930ae53cSBarry Smith {
634930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6356849ba73SBarry Smith   PetscErrorCode ierr;
63638baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
637362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
638aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
63938baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
640362ced78SSatish Balay   PetscScalar    *v,sum;
641d85afc42SSatish Balay #endif
642930ae53cSBarry Smith 
6433a40ed3dSBarry Smith   PetscFunctionBegin;
6443a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
645930ae53cSBarry Smith 
6461ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6471ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
648930ae53cSBarry Smith 
649aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6501c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6511c4feecaSSatish Balay #else
652930ae53cSBarry Smith   /* forward solve the lower triangular */
653930ae53cSBarry Smith   x[0] = b[0];
654930ae53cSBarry Smith   for (i=1; i<n; i++) {
655930ae53cSBarry Smith     ai_i = ai[i];
656930ae53cSBarry Smith     v    = aa + ai_i;
657930ae53cSBarry Smith     vi   = aj + ai_i;
658930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
659930ae53cSBarry Smith     sum  = b[i];
660930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
661930ae53cSBarry Smith     x[i] = sum;
662930ae53cSBarry Smith   }
663930ae53cSBarry Smith 
664930ae53cSBarry Smith   /* backward solve the upper triangular */
665930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
666930ae53cSBarry Smith     adiag_i = adiag[i];
667d708749eSBarry Smith     v       = aa + adiag_i + 1;
668d708749eSBarry Smith     vi      = aj + adiag_i + 1;
669930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
670930ae53cSBarry Smith     sum     = x[i];
671930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
672930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
673930ae53cSBarry Smith   }
6741c4feecaSSatish Balay #endif
675efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
6761ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6771ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6783a40ed3dSBarry Smith   PetscFunctionReturn(0);
679930ae53cSBarry Smith }
680930ae53cSBarry Smith 
6814a2ae208SSatish Balay #undef __FUNCT__
6824a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
683dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
684da3a660dSBarry Smith {
685416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
686416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6876849ba73SBarry Smith   PetscErrorCode ierr;
68838baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
68938baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
69087828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
691da3a660dSBarry Smith 
6923a40ed3dSBarry Smith   PetscFunctionBegin;
69378b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
694da3a660dSBarry Smith 
6951ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6961ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
697416022c9SBarry Smith   tmp  = a->solve_work;
698da3a660dSBarry Smith 
6993b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7003b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
701da3a660dSBarry Smith 
702da3a660dSBarry Smith   /* forward solve the lower triangular */
703da3a660dSBarry Smith   tmp[0] = b[*r++];
704da3a660dSBarry Smith   for (i=1; i<n; i++) {
705010a20caSHong Zhang     v   = aa + ai[i] ;
706010a20caSHong Zhang     vi  = aj + ai[i] ;
707416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
708da3a660dSBarry Smith     sum = b[*r++];
709010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
710da3a660dSBarry Smith     tmp[i] = sum;
711da3a660dSBarry Smith   }
712da3a660dSBarry Smith 
713da3a660dSBarry Smith   /* backward solve the upper triangular */
714da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
715010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
716010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
717416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
718da3a660dSBarry Smith     sum = tmp[i];
719010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
720010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
721da3a660dSBarry Smith     x[*c--] += tmp[i];
722da3a660dSBarry Smith   }
723da3a660dSBarry Smith 
7243b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7253b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7261ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
728efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
729898183ecSLois Curfman McInnes 
7303a40ed3dSBarry Smith   PetscFunctionReturn(0);
731da3a660dSBarry Smith }
732da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7334a2ae208SSatish Balay #undef __FUNCT__
7344a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
735dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
736da3a660dSBarry Smith {
737416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7382235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7396849ba73SBarry Smith   PetscErrorCode ierr;
74038baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
74138baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
74287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
743da3a660dSBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
7451ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7461ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
747416022c9SBarry Smith   tmp  = a->solve_work;
748da3a660dSBarry Smith 
7492235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7502235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
751da3a660dSBarry Smith 
752da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7532235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
754da3a660dSBarry Smith 
755da3a660dSBarry Smith   /* forward solve the U^T */
756da3a660dSBarry Smith   for (i=0; i<n; i++) {
757010a20caSHong Zhang     v   = aa + diag[i] ;
758010a20caSHong Zhang     vi  = aj + diag[i] + 1;
759f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
760c38d4ed2SBarry Smith     s1  = tmp[i];
761ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
762da3a660dSBarry Smith     while (nz--) {
763010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
764da3a660dSBarry Smith     }
765c38d4ed2SBarry Smith     tmp[i] = s1;
766da3a660dSBarry Smith   }
767da3a660dSBarry Smith 
768da3a660dSBarry Smith   /* backward solve the L^T */
769da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
770010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
771010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
772f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
773f1af5d2fSBarry Smith     s1  = tmp[i];
774da3a660dSBarry Smith     while (nz--) {
775010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
776da3a660dSBarry Smith     }
777da3a660dSBarry Smith   }
778da3a660dSBarry Smith 
779da3a660dSBarry Smith   /* copy tmp into x according to permutation */
780da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
781da3a660dSBarry Smith 
7822235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7832235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7841ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7851ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
786da3a660dSBarry Smith 
787efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr);
7883a40ed3dSBarry Smith   PetscFunctionReturn(0);
789da3a660dSBarry Smith }
790da3a660dSBarry Smith 
7914a2ae208SSatish Balay #undef __FUNCT__
7924a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
793dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
794da3a660dSBarry Smith {
795416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7962235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7976849ba73SBarry Smith   PetscErrorCode ierr;
79838baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
79938baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
80087828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8016abc6512SBarry Smith 
8023a40ed3dSBarry Smith   PetscFunctionBegin;
80323bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8046abc6512SBarry Smith 
8051ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8061ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
807416022c9SBarry Smith   tmp = a->solve_work;
8086abc6512SBarry Smith 
8092235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8102235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8116abc6512SBarry Smith 
8126abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8132235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8146abc6512SBarry Smith 
8156abc6512SBarry Smith   /* forward solve the U^T */
8166abc6512SBarry Smith   for (i=0; i<n; i++) {
817010a20caSHong Zhang     v   = aa + diag[i] ;
818010a20caSHong Zhang     vi  = aj + diag[i] + 1;
819f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8206abc6512SBarry Smith     tmp[i] *= *v++;
8216abc6512SBarry Smith     while (nz--) {
822010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8236abc6512SBarry Smith     }
8246abc6512SBarry Smith   }
8256abc6512SBarry Smith 
8266abc6512SBarry Smith   /* backward solve the L^T */
8276abc6512SBarry Smith   for (i=n-1; i>=0; i--){
828010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
829010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
830f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8316abc6512SBarry Smith     while (nz--) {
832010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8336abc6512SBarry Smith     }
8346abc6512SBarry Smith   }
8356abc6512SBarry Smith 
8366abc6512SBarry Smith   /* copy tmp into x according to permutation */
8376abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8386abc6512SBarry Smith 
8392235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8402235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8411ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8421ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8436abc6512SBarry Smith 
844efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
846da3a660dSBarry Smith }
8479e25ed09SBarry Smith /* ----------------------------------------------------------------*/
848dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
84908480c60SBarry Smith 
8504a2ae208SSatish Balay #undef __FUNCT__
8514a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
852dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8539e25ed09SBarry Smith {
854416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8559e25ed09SBarry Smith   IS             isicol;
8566849ba73SBarry Smith   PetscErrorCode ierr;
8575a8e39fbSHong Zhang   PetscInt       *r,*ic,n=A->m,*ai=a->i,*aj=a->j;
8585a8e39fbSHong Zhang   PetscInt       *bi,*cols,nnz,*cols_lvl;
8595a8e39fbSHong Zhang   PetscInt       *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
8605a8e39fbSHong Zhang   PetscInt       i,levels,diagonal_fill;
86177c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
862329f5518SBarry Smith   PetscReal      f;
8635a8e39fbSHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
8645a8e39fbSHong Zhang   PetscBT        lnkbt;
8655a8e39fbSHong Zhang   PetscInt       nzi,*bj,**bj_ptr,**bjlvl_ptr;
8665a8e39fbSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
8675a8e39fbSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
8689e25ed09SBarry Smith 
8693a40ed3dSBarry Smith   PetscFunctionBegin;
8706cf997b0SBarry Smith   f             = info->fill;
87138baddfdSBarry Smith   levels        = (PetscInt)info->levels;
87238baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8734c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
87482bf6240SBarry Smith 
87508480c60SBarry Smith   /* special case that simply copies fill pattern */
876be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
877be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
87886bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8792e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
88008480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
88108480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
88208480c60SBarry Smith     if (!b->diag) {
8837c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
88408480c60SBarry Smith     }
8857c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
88608480c60SBarry Smith     b->row              = isrow;
88708480c60SBarry Smith     b->col              = iscol;
88882bf6240SBarry Smith     b->icol             = isicol;
88987828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
890f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
891c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
892c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
8933a40ed3dSBarry Smith     PetscFunctionReturn(0);
89408480c60SBarry Smith   }
89508480c60SBarry Smith 
896898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
897898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
8989e25ed09SBarry Smith 
8999e25ed09SBarry Smith   /* get new row pointers */
9005a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9015a8e39fbSHong Zhang   bi[0] = 0;
9025a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9035a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9045a8e39fbSHong Zhang   bdiag[0]  = 0;
9056cf997b0SBarry Smith 
9065a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9075a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9085a8e39fbSHong Zhang 
9095a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9105a8e39fbSHong Zhang   nlnk = n + 1;
9115a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9125a8e39fbSHong Zhang 
9135a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
9145a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9155a8e39fbSHong Zhang   current_space = free_space;
9165a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9175a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
9185a8e39fbSHong Zhang 
9195a8e39fbSHong Zhang   for (i=0; i<n; i++) {
9205a8e39fbSHong Zhang     nzi = 0;
9216cf997b0SBarry Smith     /* copy current row into linked list */
9225a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
9235a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
9245a8e39fbSHong Zhang     cols = aj + ai[r[i]];
9255a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
9265a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9275a8e39fbSHong Zhang     nzi += nlnk;
9286cf997b0SBarry Smith 
9296cf997b0SBarry Smith     /* make sure diagonal entry is included */
9305a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
9316cf997b0SBarry Smith       fm = n;
9325a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
9335a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
9345a8e39fbSHong Zhang       lnk[fm]    = i;
9355a8e39fbSHong Zhang       lnk_lvl[i] = 0;
9365a8e39fbSHong Zhang       nzi++; dcount++;
9376cf997b0SBarry Smith     }
9386cf997b0SBarry Smith 
9395a8e39fbSHong Zhang     /* add pivot rows into the active row */
9405a8e39fbSHong Zhang     nzbd = 0;
9415a8e39fbSHong Zhang     prow = lnk[n];
9425a8e39fbSHong Zhang     while (prow < i) {
9435a8e39fbSHong Zhang       nnz      = bdiag[prow];
9445a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
9455a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
9465a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
9475a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
9485a8e39fbSHong Zhang       nzi += nlnk;
9495a8e39fbSHong Zhang       prow = lnk[prow];
9505a8e39fbSHong Zhang       nzbd++;
95156beaf04SBarry Smith     }
9525a8e39fbSHong Zhang     bdiag[i] = nzbd;
9535a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
954ecf371e4SBarry Smith 
9555a8e39fbSHong Zhang     /* if free space is not available, make more free space */
9565a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
9575a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
9585a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
9595a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space_lvl);CHKERRQ(ierr);
9605a8e39fbSHong Zhang       reallocs++;
9615a8e39fbSHong Zhang     }
962ecf371e4SBarry Smith 
9635a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
9645a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
9655a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
9665a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
9675a8e39fbSHong Zhang 
9685a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
9695a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
97077431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
9715a8e39fbSHong Zhang     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i);
9726cf997b0SBarry Smith     }
9735a8e39fbSHong Zhang 
9745a8e39fbSHong Zhang     current_space->array           += nzi;
9755a8e39fbSHong Zhang     current_space->local_used      += nzi;
9765a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
9775a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
9785a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
9795a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
9809e25ed09SBarry Smith   }
9815a8e39fbSHong Zhang 
982898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
983898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
9845a8e39fbSHong Zhang 
9855a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
9865a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
9875a8e39fbSHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
9885a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
9895a8e39fbSHong Zhang   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
9905a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
9919e25ed09SBarry Smith 
99263ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
993f64065bbSBarry Smith   {
9945a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
99563ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr);
99663ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af));CHKERRQ(ierr);
99763ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af));CHKERRQ(ierr);
99863ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr);
999335d9088SBarry Smith     if (diagonal_fill) {
100063ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount));CHKERRQ(ierr);
1001335d9088SBarry Smith     }
1002f64065bbSBarry Smith   }
100363ba0a88SBarry Smith #endif
1004bbb6d6a8SBarry Smith 
10059e25ed09SBarry Smith   /* put together the new matrix */
1006f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1007f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1008ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
100952e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1010416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1011a96a251dSBarry Smith   b->freedata     = PETSC_TRUE;
10127c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10135a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1014b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10155a8e39fbSHong Zhang   b->j          = bj;
10165a8e39fbSHong Zhang   b->i          = bi;
10175a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
10185a8e39fbSHong Zhang   b->diag       = bdiag;
1019416022c9SBarry Smith   b->ilen       = 0;
1020416022c9SBarry Smith   b->imax       = 0;
1021416022c9SBarry Smith   b->row        = isrow;
1022416022c9SBarry Smith   b->col        = iscol;
1023c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1024c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
102582bf6240SBarry Smith   b->icol       = isicol;
102687828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1027416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
10285a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
10295a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
10305a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
10319e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1032418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1033ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
10345a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
103571c2f376SKris Buschelman 
103654e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
103771c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
10384846f1f5SKris Buschelman 
10393a40ed3dSBarry Smith   PetscFunctionReturn(0);
10409e25ed09SBarry Smith }
1041d5d45c9bSBarry Smith 
10423c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1043a6175056SHong Zhang #undef __FUNCT__
1044a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1045af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1046a6175056SHong Zhang {
10472ed133dbSHong Zhang   Mat            C = *B;
10480968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10492ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
10502ed133dbSHong Zhang   IS             ip=b->row;
10512ed133dbSHong Zhang   PetscErrorCode ierr;
10522ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
10532ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
10542ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1055622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1056540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1057540703acSHong Zhang   PetscTruth     shiftpd;
1058540703acSHong Zhang   ChShift_Ctx    sctx;
1059540703acSHong Zhang   PetscInt       newshift;
1060d5d45c9bSBarry Smith 
1061a6175056SHong Zhang   PetscFunctionBegin;
1062540703acSHong Zhang   shiftnz   = info->shiftnz;
1063540703acSHong Zhang   shiftpd   = info->shiftpd;
1064ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1065ee45ca4aSHong Zhang 
10662ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
10672ed133dbSHong Zhang 
10682ed133dbSHong Zhang   /* initialization */
10692ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
10702ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
10712ed133dbSHong Zhang   jl   = il + mbs;
10722ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
10732ed133dbSHong Zhang 
1074540703acSHong Zhang   sctx.shift_amount = 0;
1075540703acSHong Zhang   sctx.nshift       = 0;
10762ed133dbSHong Zhang   do {
1077540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
10782ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
10792ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1080a6175056SHong Zhang     }
10812ed133dbSHong Zhang 
10822ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1083c05c3958SHong Zhang       bval = ba + bi[k];
10842ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
10852ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
10862ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
10872ed133dbSHong Zhang         col = rip[aj[j]];
10882ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
10892ed133dbSHong Zhang           rtmp[col] = aa[j];
1090c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
10912ed133dbSHong Zhang         }
10922ed133dbSHong Zhang       }
109339e3d630SHong Zhang       /* shift the diagonal of the matrix */
1094540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
10952ed133dbSHong Zhang 
10962ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
10972ed133dbSHong Zhang       dk = rtmp[k];
10982ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
10992ed133dbSHong Zhang 
11002ed133dbSHong Zhang       while (i < k){
11012ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11022ed133dbSHong Zhang 
11032ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11042ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11052ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11062ed133dbSHong Zhang         dk += uikdi*ba[ili];
11072ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11082ed133dbSHong Zhang 
11092ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11102ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11112ed133dbSHong Zhang         if (jmin < jmax){
11122ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11132ed133dbSHong Zhang           /* update il and jl for row i */
11142ed133dbSHong Zhang           il[i] = jmin;
11152ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11162ed133dbSHong Zhang         }
11172ed133dbSHong Zhang         i = nexti;
11182ed133dbSHong Zhang       }
11192ed133dbSHong Zhang 
11200697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11210697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11220697b6caSHong Zhang       rs   = 0.0;
11233c9e8598SHong Zhang       jmin = bi[k]+1;
11243c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11253c9e8598SHong Zhang       if (nz){
11263c9e8598SHong Zhang         bcol = bj + jmin;
11273c9e8598SHong Zhang         while (nz--){
112889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
112989f845c8SHong Zhang           bcol++;
11303c9e8598SHong Zhang         }
11313c9e8598SHong Zhang       }
1132540703acSHong Zhang 
1133540703acSHong Zhang       sctx.rs = rs;
1134540703acSHong Zhang       sctx.pv = dk;
11353e8c821fSHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1136540703acSHong Zhang       if (newshift == 1){
1137540703acSHong Zhang         break;    /* sctx.shift_amount is updated */
1138540703acSHong Zhang       } else if (newshift == -1){
11390697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11403c9e8598SHong Zhang       }
11413c9e8598SHong Zhang 
11423c9e8598SHong Zhang       /* copy data into U(k,:) */
114339e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
114439e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
114539e3d630SHong Zhang       if (jmin < jmax) {
114639e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
114739e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
11483c9e8598SHong Zhang         }
114939e3d630SHong Zhang         /* add the k-th row into il and jl */
11503c9e8598SHong Zhang         il[k] = jmin;
11513c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
11523c9e8598SHong Zhang       }
11533c9e8598SHong Zhang     }
1154540703acSHong Zhang   } while (sctx.chshift);
11553c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
11563c9e8598SHong Zhang 
115739e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
11583c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
11593c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
11603c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1161efee365bSSatish Balay   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
1162540703acSHong Zhang   if (sctx.nshift){
1163540703acSHong Zhang     if (shiftnz) {
116463ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1165540703acSHong Zhang     } else if (shiftpd) {
116663ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
11670697b6caSHong Zhang     }
11683c9e8598SHong Zhang   }
11693c9e8598SHong Zhang   PetscFunctionReturn(0);
11703c9e8598SHong Zhang }
1171a6175056SHong Zhang 
1172a6175056SHong Zhang #undef __FUNCT__
1173a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1174dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1175a6175056SHong Zhang {
11760968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1177ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1178ed59401aSHong Zhang   Mat            B;
1179dfbe8321SBarry Smith   PetscErrorCode ierr;
1180653a6975SHong Zhang   PetscTruth     perm_identity;
1181622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1182ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1183391eac55SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
11845a8e39fbSHong Zhang   PetscInt       ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1185ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
11867a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
11877a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
11887a48dd6fSHong Zhang   PetscBT        lnkbt;
1189a6175056SHong Zhang 
1190a6175056SHong Zhang   PetscFunctionBegin;
1191653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
11927a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
11937a48dd6fSHong Zhang 
119439e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
119539e3d630SHong Zhang   ui[0] = 0;
119639e3d630SHong Zhang 
1197622e2028SHong Zhang   /* special case that simply copies fill pattern */
1198622e2028SHong Zhang   if (!levels && perm_identity) {
1199622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1200ed59401aSHong Zhang     for (i=0; i<am; i++) {
120139e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1202ed59401aSHong Zhang     }
120339e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
120439e3d630SHong Zhang     cols = uj;
1205ed59401aSHong Zhang     for (i=0; i<am; i++) {
1206ed59401aSHong Zhang       aj    = a->j + a->diag[i];
120739e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
120839e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1209ed59401aSHong Zhang     }
121039e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12117a48dd6fSHong Zhang     /* initialization */
12125a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12137a48dd6fSHong Zhang 
12147a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12157a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12161393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12177a48dd6fSHong Zhang     il         = jl + am;
12187a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12197a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12207a48dd6fSHong Zhang     for (i=0; i<am; i++){
12217a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12227a48dd6fSHong Zhang     }
12237a48dd6fSHong Zhang 
12247a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12257a48dd6fSHong Zhang     nlnk = am + 1;
12262ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12277a48dd6fSHong Zhang 
12287a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12297a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12307a48dd6fSHong Zhang     current_space = free_space;
12317a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12327a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12337a48dd6fSHong Zhang 
12347a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12357a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12367a48dd6fSHong Zhang       nzk   = 0;
1237622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1238622e2028SHong Zhang       ncols_upper = 0;
1239622e2028SHong Zhang       for (j=0; j<ncols; j++){
12405a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
12415a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
12425a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1243622e2028SHong Zhang           ncols_upper++;
1244622e2028SHong Zhang         }
1245622e2028SHong Zhang       }
12465a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12477a48dd6fSHong Zhang       nzk += nlnk;
12487a48dd6fSHong Zhang 
12497a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
12507a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
12517a48dd6fSHong Zhang 
12527a48dd6fSHong Zhang       while (prow < k){
12537a48dd6fSHong Zhang         nextprow = jl[prow];
12547a48dd6fSHong Zhang 
12557a48dd6fSHong Zhang         /* merge prow into k-th row */
12567a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
12577a48dd6fSHong Zhang         jmax = ui[prow+1];
12587a48dd6fSHong Zhang         ncols = jmax-jmin;
1259ed59401aSHong Zhang         i     = jmin - ui[prow];
1260ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
12615a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
12625a8e39fbSHong Zhang         j     = *(uj - 1);
12635a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
12647a48dd6fSHong Zhang         nzk += nlnk;
12657a48dd6fSHong Zhang 
12667a48dd6fSHong Zhang         /* update il and jl for prow */
12677a48dd6fSHong Zhang         if (jmin < jmax){
12687a48dd6fSHong Zhang           il[prow] = jmin;
12697a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
12707a48dd6fSHong Zhang         }
12717a48dd6fSHong Zhang         prow = nextprow;
12727a48dd6fSHong Zhang       }
12737a48dd6fSHong Zhang 
12747a48dd6fSHong Zhang       /* if free space is not available, make more free space */
12757a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
12767a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
12777a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
12787a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
12797a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
12807a48dd6fSHong Zhang         reallocs++;
12817a48dd6fSHong Zhang       }
12827a48dd6fSHong Zhang 
12837a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
12842ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
12857a48dd6fSHong Zhang 
12867a48dd6fSHong Zhang       /* add the k-th row into il and jl */
128739e3d630SHong Zhang       if (nzk > 1){
12887a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
12897a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
12907a48dd6fSHong Zhang         il[k] = ui[k] + 1;
12917a48dd6fSHong Zhang       }
12927a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
12937a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
12947a48dd6fSHong Zhang 
12957a48dd6fSHong Zhang       current_space->array           += nzk;
12967a48dd6fSHong Zhang       current_space->local_used      += nzk;
12977a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
12987a48dd6fSHong Zhang 
12997a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13007a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13017a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13027a48dd6fSHong Zhang 
13037a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13047a48dd6fSHong Zhang     }
13057a48dd6fSHong Zhang 
130663ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
13077a48dd6fSHong Zhang     if (ai[am] != 0) {
130839e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
130963ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
131063ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
131163ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
13127a48dd6fSHong Zhang     } else {
131363ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr);
13147a48dd6fSHong Zhang     }
131563ba0a88SBarry Smith #endif
13167a48dd6fSHong Zhang 
13177a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13187a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13195a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13207a48dd6fSHong Zhang 
13217a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13227a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13237a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13242ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13257a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13267a48dd6fSHong Zhang 
132739e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
132839e3d630SHong Zhang 
13297a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13307a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1331ed59401aSHong Zhang   B = *fact;
1332ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1333ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13347a48dd6fSHong Zhang 
1335ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13367a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13377a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13387a48dd6fSHong Zhang   b->j    = uj;
13397a48dd6fSHong Zhang   b->i    = ui;
13407a48dd6fSHong Zhang   b->diag = 0;
13417a48dd6fSHong Zhang   b->ilen = 0;
13427a48dd6fSHong Zhang   b->imax = 0;
13437a48dd6fSHong Zhang   b->row  = perm;
13447a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
13457a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13467a48dd6fSHong Zhang   b->icol = perm;
13477a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13487a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
134952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13507a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
13517a48dd6fSHong Zhang 
1352ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1353ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1354ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
13557a48dd6fSHong Zhang   if (ai[am] != 0) {
1356ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
13577a48dd6fSHong Zhang   } else {
1358ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
13597a48dd6fSHong Zhang   }
136039e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1361622e2028SHong Zhang   if (perm_identity){
1362ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1363ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1364622e2028SHong Zhang   }
1365a6175056SHong Zhang   PetscFunctionReturn(0);
1366a6175056SHong Zhang }
1367d5d45c9bSBarry Smith 
1368f76d2b81SHong Zhang #undef __FUNCT__
1369f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1370dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1371f76d2b81SHong Zhang {
1372f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
137310c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
13742ed133dbSHong Zhang   Mat            B;
1375dfbe8321SBarry Smith   PetscErrorCode ierr;
1376f76d2b81SHong Zhang   PetscTruth     perm_identity;
137710c27e3dSHong Zhang   PetscReal      fill = info->fill;
13781393bc97SHong Zhang   PetscInt       *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
137910c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
13802ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
138110c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
138210c27e3dSHong Zhang   PetscBT        lnkbt;
13832ed133dbSHong Zhang   IS             iperm;
1384f76d2b81SHong Zhang 
1385f76d2b81SHong Zhang   PetscFunctionBegin;
13862ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1387f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
13882ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
13892ed133dbSHong Zhang 
1390f76d2b81SHong Zhang   if (!perm_identity){
13912ed133dbSHong Zhang     /* check if perm is symmetric! */
13922ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
13932ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
13941393bc97SHong Zhang     for (i=0; i<am; i++) {
13952ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
13962ed133dbSHong Zhang     }
13972ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
13982ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1399f76d2b81SHong Zhang   }
140010c27e3dSHong Zhang 
140110c27e3dSHong Zhang   /* initialization */
14021393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
140310c27e3dSHong Zhang   ui[0] = 0;
140410c27e3dSHong Zhang 
140510c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14061393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14071393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14081393bc97SHong Zhang   il     = jl + am;
14091393bc97SHong Zhang   cols   = il + am;
14101393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14111393bc97SHong Zhang   for (i=0; i<am; i++){
14121393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1413f76d2b81SHong Zhang   }
141410c27e3dSHong Zhang 
141510c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14161393bc97SHong Zhang   nlnk = am + 1;
14171393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
141810c27e3dSHong Zhang 
14191393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
14201393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
142110c27e3dSHong Zhang   current_space = free_space;
142210c27e3dSHong Zhang 
14231393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
142410c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
142510c27e3dSHong Zhang     nzk   = 0;
14262ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14272ed133dbSHong Zhang     ncols_upper = 0;
14282ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1429622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14302ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14312ed133dbSHong Zhang         cols[ncols_upper] = i;
14322ed133dbSHong Zhang         ncols_upper++;
14332ed133dbSHong Zhang       }
14342ed133dbSHong Zhang     }
14351393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
143610c27e3dSHong Zhang     nzk += nlnk;
143710c27e3dSHong Zhang 
143810c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
143910c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
144010c27e3dSHong Zhang 
144110c27e3dSHong Zhang     while (prow < k){
144210c27e3dSHong Zhang       nextprow = jl[prow];
144310c27e3dSHong Zhang       /* merge prow into k-th row */
14441393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
144510c27e3dSHong Zhang       jmax = ui[prow+1];
144610c27e3dSHong Zhang       ncols = jmax-jmin;
14471393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
14485a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
144910c27e3dSHong Zhang       nzk += nlnk;
145010c27e3dSHong Zhang 
145110c27e3dSHong Zhang       /* update il and jl for prow */
145210c27e3dSHong Zhang       if (jmin < jmax){
145310c27e3dSHong Zhang         il[prow] = jmin;
14542ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
145510c27e3dSHong Zhang       }
145610c27e3dSHong Zhang       prow = nextprow;
145710c27e3dSHong Zhang     }
145810c27e3dSHong Zhang 
145910c27e3dSHong Zhang     /* if free space is not available, make more free space */
146010c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
14611393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
146210c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
146310c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
146410c27e3dSHong Zhang       reallocs++;
146510c27e3dSHong Zhang     }
146610c27e3dSHong Zhang 
146710c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
14681393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
146910c27e3dSHong Zhang 
147010c27e3dSHong Zhang     /* add the k-th row into il and jl */
147110c27e3dSHong Zhang     if (nzk-1 > 0){
14721393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
147310c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
147410c27e3dSHong Zhang       il[k] = ui[k] + 1;
147510c27e3dSHong Zhang     }
14762ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
147710c27e3dSHong Zhang     current_space->array           += nzk;
147810c27e3dSHong Zhang     current_space->local_used      += nzk;
147910c27e3dSHong Zhang     current_space->local_remaining -= nzk;
148010c27e3dSHong Zhang 
148110c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
148210c27e3dSHong Zhang   }
148310c27e3dSHong Zhang 
148463ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
14851393bc97SHong Zhang   if (ai[am] != 0) {
14861393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
148763ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
148863ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
148963ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
149010c27e3dSHong Zhang   } else {
149163ba0a88SBarry Smith      ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr);
149210c27e3dSHong Zhang   }
149363ba0a88SBarry Smith #endif
149410c27e3dSHong Zhang 
149510c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
149610c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
149710c27e3dSHong Zhang 
149810c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
14991393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
150010c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
150110c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
150210c27e3dSHong Zhang 
150310c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15041393bc97SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
15052ed133dbSHong Zhang   B    = *fact;
15062ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1507ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
150810c27e3dSHong Zhang 
15092ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
151010c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
15111393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
151210c27e3dSHong Zhang   b->j    = uj;
151310c27e3dSHong Zhang   b->i    = ui;
151410c27e3dSHong Zhang   b->diag = 0;
151510c27e3dSHong Zhang   b->ilen = 0;
151610c27e3dSHong Zhang   b->imax = 0;
151710c27e3dSHong Zhang   b->row  = perm;
151810c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
151910c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
152010c27e3dSHong Zhang   b->icol = perm;
152110c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15221393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15231393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
15241393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
152510c27e3dSHong Zhang 
15262ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15272ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15282ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15291393bc97SHong Zhang   if (ai[am] != 0) {
15301393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
153110c27e3dSHong Zhang   } else {
15322ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
153310c27e3dSHong Zhang   }
153439e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15352ed133dbSHong Zhang   if (perm_identity){
153610c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
153710c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15382ed133dbSHong Zhang   }
1539f76d2b81SHong Zhang   PetscFunctionReturn(0);
1540f76d2b81SHong Zhang }
1541