xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 5a8e39fb7ca75afe734b11fd5263a37b1f94dc5e)
1289bc588SBarry Smith 
270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
31c4feecaSSatish Balay #include "src/inline/dot.h"
4ed480e8bSBarry Smith #include "src/inline/spops.h"
51393bc97SHong Zhang #include "petscbt.h"
61393bc97SHong Zhang #include "src/mat/utils/freespace.h"
73b2fbd54SBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
10dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
113b2fbd54SBarry Smith {
123a40ed3dSBarry Smith   PetscFunctionBegin;
133a40ed3dSBarry Smith 
1429bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
15aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
16d64ed03dSBarry Smith   PetscFunctionReturn(0);
17d64ed03dSBarry Smith #endif
183b2fbd54SBarry Smith }
193b2fbd54SBarry Smith 
2086bacbd2SBarry Smith 
21dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
22dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
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 */
51dfbe8321SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
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);
224f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,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;
22986bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
23086bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23107d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
232b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23386bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
23486bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
23586bacbd2SBarry Smith   b->a             = new_a;
23686bacbd2SBarry Smith   b->j             = new_j;
23786bacbd2SBarry Smith   b->i             = new_i;
23886bacbd2SBarry Smith   b->ilen          = 0;
23986bacbd2SBarry Smith   b->imax          = 0;
2401f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
241313ae024SBarry Smith   b->row           = isirow;
24286bacbd2SBarry Smith   b->col           = iscolf;
24387828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24486bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24586bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24686bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24786bacbd2SBarry Smith 
24886bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24986bacbd2SBarry Smith 
25086bacbd2SBarry Smith   /* check out for identical nodes. If found, use inode functions */
2513a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
25286bacbd2SBarry Smith 
253431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
254b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
255b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
256b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
257b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
258431e710aSBarry Smith 
25986bacbd2SBarry Smith   PetscFunctionReturn(0);
26060ec86bdSBarry Smith #endif
26186bacbd2SBarry Smith }
26286bacbd2SBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
264b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
265dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
266289bc588SBarry Smith {
267416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
268289bc588SBarry Smith   IS             isicol;
2696849ba73SBarry Smith   PetscErrorCode ierr;
27038baddfdSBarry Smith   PetscInt       *r,*ic,i,n=A->m,*ai=a->i,*aj=a->j;
2711393bc97SHong Zhang   PetscInt       *bi,*bj,*ajtmp;
2721393bc97SHong Zhang   PetscInt       *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2739e046878SBarry Smith   PetscReal      f;
2741393bc97SHong Zhang   PetscInt       nlnk,*lnk,k,*cols,**bi_ptr;
2751393bc97SHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
2761393bc97SHong Zhang   PetscBT        lnkbt;
277289bc588SBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionBegin;
27929bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2804c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
281ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
282ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
283289bc588SBarry Smith 
284289bc588SBarry Smith   /* get new row pointers */
2851393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2861393bc97SHong Zhang   bi[0] = 0;
2871393bc97SHong Zhang 
2881393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2891393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2901393bc97SHong Zhang   bdiag[0] = 0;
2911393bc97SHong Zhang 
2921393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2931393bc97SHong Zhang   nlnk = n + 1;
2941393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2951393bc97SHong Zhang 
2961393bc97SHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&cols);CHKERRQ(ierr);
2971393bc97SHong Zhang   im     = cols + n;
2981393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2991393bc97SHong Zhang 
3001393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
301d3d32019SBarry Smith   f = info->fill;
3021393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3031393bc97SHong Zhang   current_space = free_space;
304289bc588SBarry Smith 
305289bc588SBarry Smith   for (i=0; i<n; i++) {
3061393bc97SHong Zhang     /* copy previous fill into linked list */
307289bc588SBarry Smith     nzi = 0;
3081393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3091393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3101393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
3111393bc97SHong Zhang     for (k=0; k<nnz; k++) cols[k] = ic[*(ajtmp+k)]; /* note: cols is not sorted when iscol!=indentity */
3121393bc97SHong Zhang     ierr = PetscLLAdd(nnz,cols,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3131393bc97SHong Zhang     nzi += nlnk;
3141393bc97SHong Zhang 
3151393bc97SHong Zhang     /* add pivot rows into linked list */
3161393bc97SHong Zhang     row = lnk[n];
3171393bc97SHong Zhang     while (row < i) {
3181393bc97SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1;
3191393bc97SHong Zhang       ajtmp   = bi_ptr[row] + nzbd;
3201393bc97SHong Zhang       nnz     = im[row] - nzbd; /* num of columns with row<indices<=i */
3211393bc97SHong Zhang       im[row] = nzbd;
3221393bc97SHong Zhang       ierr = PetscLLAddSortedLU(nnz,ajtmp,row,nlnk,lnk,lnkbt,i,nzbd);CHKERRQ(ierr);
3231393bc97SHong Zhang       nzi     += nlnk;
3241393bc97SHong Zhang       im[row] += nzbd;  /* update im[row]: num of cols with index<=i */
3251393bc97SHong Zhang 
3261393bc97SHong Zhang       row = lnk[row];
327289bc588SBarry Smith     }
3281393bc97SHong Zhang 
3291393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3301393bc97SHong Zhang     im[i]   = nzi;
3311393bc97SHong Zhang 
3321393bc97SHong Zhang     /* mark bdiag */
3331393bc97SHong Zhang     nzbd = 0;
3341393bc97SHong Zhang     nnz  = nzi;
3351393bc97SHong Zhang     k    = lnk[n];
3361393bc97SHong Zhang     while (nnz-- && k < i){
3371393bc97SHong Zhang       nzbd++;
3381393bc97SHong Zhang       k = lnk[k];
3391393bc97SHong Zhang     }
3401393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3411393bc97SHong Zhang 
3421393bc97SHong Zhang     /* if free space is not available, make more free space */
3431393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3441393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
3451393bc97SHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
3461393bc97SHong Zhang       reallocs++;
3471393bc97SHong Zhang     }
3481393bc97SHong Zhang 
3491393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3501393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3511393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3521393bc97SHong Zhang     current_space->array           += nzi;
3531393bc97SHong Zhang     current_space->local_used      += nzi;
3541393bc97SHong Zhang     current_space->local_remaining -= nzi;
355289bc588SBarry Smith   }
356464e8d44SSatish Balay   if (ai[n] != 0) {
3571393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
358418422e8SSatish Balay     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
359b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
360b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
361b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
36205bf559cSSatish Balay   } else {
363b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
36405bf559cSSatish Balay   }
3652fb3ed76SBarry Smith 
366898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
367898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3681987afe7SBarry Smith 
3691393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3701393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
3711393bc97SHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3721393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3731393bc97SHong Zhang   ierr = PetscFree(cols);CHKERRQ(ierr);
374289bc588SBarry Smith 
375289bc588SBarry Smith   /* put together the new matrix */
376f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr);
377f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
378f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr);
37952e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
380416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
381606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3827c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
383e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
384606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
385606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
3861393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3871393bc97SHong Zhang   b->j          = bj;
3881393bc97SHong Zhang   b->i          = bi;
3891393bc97SHong Zhang   b->diag       = bdiag;
390416022c9SBarry Smith   b->ilen       = 0;
391416022c9SBarry Smith   b->imax       = 0;
392416022c9SBarry Smith   b->row        = isrow;
393416022c9SBarry Smith   b->col        = iscol;
394c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
395c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
39682bf6240SBarry Smith   b->icol       = isicol;
39787828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3981393bc97SHong Zhang 
3991393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
4001393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
4011393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4027fd98bd6SLois Curfman McInnes 
403329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
404418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
405ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4063a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4077dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
408703deb49SSatish Balay 
409e93ccc4dSBarry Smith   if (ai[n] != 0) {
4101393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
411eea2913aSSatish Balay   } else {
412eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
413eea2913aSSatish Balay   }
4143a40ed3dSBarry Smith   PetscFunctionReturn(0);
415289bc588SBarry Smith }
4161393bc97SHong Zhang 
4170a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
418dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
41941c01911SSatish Balay 
4204a2ae208SSatish Balay #undef __FUNCT__
4214a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
422af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
423289bc588SBarry Smith {
424416022c9SBarry Smith   Mat            C=*B;
425416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
42682bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4276849ba73SBarry Smith   PetscErrorCode ierr;
428b3bf805bSHong Zhang   PetscInt       *r,*ic,i,j,n=A->m,*bi=b->i,*bj=b->j;
429b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
430d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
43187828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
432540703acSHong Zhang   PetscReal      zeropivot,rs,d,shift_nz;
433d2276718SHong Zhang   PetscReal      row_shift,shift_top=0.;
434d2276718SHong Zhang   PetscTruth     shift_pd;
435b3bf805bSHong Zhang   LUShift_Ctx    sctx;
436d2276718SHong Zhang   PetscInt       newshift;
437289bc588SBarry Smith 
4383a40ed3dSBarry Smith   PetscFunctionBegin;
439540703acSHong Zhang   shift_nz  = info->shiftnz;
4400a29876aSHong Zhang   shift_pd  = info->shiftpd;
4410a29876aSHong Zhang   zeropivot = info->zeropivot;
4420a29876aSHong Zhang 
44378b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44478b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44587828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44687828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
447010a20caSHong Zhang   rtmps = rtmp; ics = ic;
448289bc588SBarry Smith 
4496cc28720Svictorle   if (!a->diag) {
4500cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4510cf777b8SBarry Smith   }
452e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
453540703acSHong Zhang   if (shift_pd && shift_nz) shift_nz = 0.0;
454e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45538baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4560cf777b8SBarry Smith     shift_top = 0;
4576cc28720Svictorle     for (i=0; i<n; i++) {
458010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4596cc28720Svictorle       /* calculate amt of shift needed for this row */
4606c037d1bSvictorle       if (d<=0) {
4613aeef889SHong Zhang         row_shift = 0;
4623aeef889SHong Zhang       } else {
4633aeef889SHong Zhang         row_shift = -2*d;
4643aeef889SHong Zhang       }
465010a20caSHong Zhang       v  = a->a+aai[i];
466e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
467e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4686cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4696cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4706cc28720Svictorle     }
4716cc28720Svictorle   }
4726cc28720Svictorle 
473d2276718SHong Zhang   sctx.shift_amount = 0;
474d2276718SHong Zhang   sctx.shift_top    = shift_top;
475d2276718SHong Zhang   sctx.nshift       = 0;
476d2276718SHong Zhang   sctx.nshift_max   = 5;
477d2276718SHong Zhang   sctx.shift_lo     = 0.;
478d2276718SHong Zhang   sctx.shift_hi     = 1.;
479085256b3SBarry Smith   do {
480d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
481289bc588SBarry Smith     for (i=0; i<n; i++){
482b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
483b3bf805bSHong Zhang       bjtmp = bj + bi[i];
484b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
485289bc588SBarry Smith 
486289bc588SBarry Smith       /* load in initial (unfactored row) */
487416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
488b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
489010a20caSHong Zhang       v     = a->a + a->i[r[i]];
490085256b3SBarry Smith       for (j=0; j<nz; j++) {
491b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
492085256b3SBarry Smith       }
493d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
494289bc588SBarry Smith 
495b3bf805bSHong Zhang       row = *bjtmp++;
496f2582111SSatish Balay       while  (row < i) {
4978c37ef55SBarry Smith         pc = rtmp + row;
4988c37ef55SBarry Smith         if (*pc != 0.0) {
499010a20caSHong Zhang           pv         = b->a + diag_offset[row];
500010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50135aab85fSBarry Smith           multiplier = *pc / *pv++;
5028c37ef55SBarry Smith           *pc        = multiplier;
503b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
504f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
505b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
506289bc588SBarry Smith         }
507b3bf805bSHong Zhang         row = *bjtmp++;
508289bc588SBarry Smith       }
509416022c9SBarry Smith       /* finished row so stick it into b->a */
510b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
511b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
512b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
513b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
514d3d32019SBarry Smith       rs   = 0.0;
515d3d32019SBarry Smith       for (j=0; j<nz; j++) {
516d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
517d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
518d3d32019SBarry Smith       }
519d2276718SHong Zhang 
520b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
521d2276718SHong Zhang       sctx.rs  = rs;
522d2276718SHong Zhang       sctx.pv  = pv[diag];
523540703acSHong Zhang       ierr = Mat_LUCheckShift(info,&sctx,&newshift);CHKERRQ(ierr);
524d2276718SHong Zhang       if (newshift == 1){
525b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
526d2276718SHong Zhang       } else if (newshift == -1){
527d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
528d708749eSBarry Smith       }
52935aab85fSBarry Smith     }
530d2276718SHong Zhang 
531d2276718SHong Zhang     if (shift_pd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5326cc28720Svictorle       /*
5336c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5346cc28720Svictorle        * then try lower shift
5356cc28720Svictorle        */
536d2276718SHong Zhang       sctx.shift_hi       = info->shift_fraction;
537d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
538d2276718SHong Zhang       sctx.shift_amount   = info->shift_fraction * sctx.shift_top;
539d2276718SHong Zhang       sctx.lushift        = PETSC_TRUE;
540d2276718SHong Zhang       sctx.nshift++;
5416cc28720Svictorle     }
542d2276718SHong Zhang   } while (sctx.lushift);
5430f11f9aeSSatish Balay 
54435aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
54535aab85fSBarry Smith   for (i=0; i<n; i++) {
546010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
54735aab85fSBarry Smith   }
54835aab85fSBarry Smith 
549606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55078b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
552416022c9SBarry Smith   C->factor = FACTOR_LU;
5537dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
554c456f294SBarry Smith   C->assembled = PETSC_TRUE;
555b0a32e0cSBarry Smith   PetscLogFlops(C->n);
556d2276718SHong Zhang   if (sctx.nshift){
557540703acSHong Zhang     if (shift_nz) {
558540703acSHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
5590697b6caSHong Zhang     } else if (shift_pd) {
560d2276718SHong Zhang       b->lu_shift_fraction = info->shift_fraction;
561d2276718SHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %D\n",info->shift_fraction,shift_top,sctx.nshift);
5626cc28720Svictorle     }
5630697b6caSHong Zhang   }
5643a40ed3dSBarry Smith   PetscFunctionReturn(0);
565289bc588SBarry Smith }
5660889b5dcSvictorle 
5670889b5dcSvictorle #undef __FUNCT__
5680889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
569dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5700889b5dcSvictorle {
5710889b5dcSvictorle   PetscFunctionBegin;
5720889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5730889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5740889b5dcSvictorle   PetscFunctionReturn(0);
5750889b5dcSvictorle }
5760889b5dcSvictorle 
5770889b5dcSvictorle 
5780a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5794a2ae208SSatish Balay #undef __FUNCT__
5804a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
581dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
582da3a660dSBarry Smith {
583dfbe8321SBarry Smith   PetscErrorCode ierr;
584416022c9SBarry Smith   Mat            C;
585416022c9SBarry Smith 
5863a40ed3dSBarry Smith   PetscFunctionBegin;
5879e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
588af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
589273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
59052e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5913a40ed3dSBarry Smith   PetscFunctionReturn(0);
592da3a660dSBarry Smith }
5930a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5944a2ae208SSatish Balay #undef __FUNCT__
5954a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
596dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5978c37ef55SBarry Smith {
598416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
599416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6006849ba73SBarry Smith   PetscErrorCode ierr;
60138baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
60238baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
60387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6048c37ef55SBarry Smith 
6053a40ed3dSBarry Smith   PetscFunctionBegin;
6063a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
60791bf9adeSBarry Smith 
6081ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6091ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
610416022c9SBarry Smith   tmp  = a->solve_work;
6118c37ef55SBarry Smith 
6123b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6133b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6148c37ef55SBarry Smith 
6158c37ef55SBarry Smith   /* forward solve the lower triangular */
6168c37ef55SBarry Smith   tmp[0] = b[*r++];
617010a20caSHong Zhang   tmps   = tmp;
6188c37ef55SBarry Smith   for (i=1; i<n; i++) {
619010a20caSHong Zhang     v   = aa + ai[i] ;
620010a20caSHong Zhang     vi  = aj + ai[i] ;
62153e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
62253e7425aSBarry Smith     sum = b[*r++];
623ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6248c37ef55SBarry Smith     tmp[i] = sum;
6258c37ef55SBarry Smith   }
6268c37ef55SBarry Smith 
6278c37ef55SBarry Smith   /* backward solve the upper triangular */
6288c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
629010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
630010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
631416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6328c37ef55SBarry Smith     sum = tmp[i];
633ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
634010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6358c37ef55SBarry Smith   }
6368c37ef55SBarry Smith 
6373b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6383b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6391ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
641b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6423a40ed3dSBarry Smith   PetscFunctionReturn(0);
6438c37ef55SBarry Smith }
644026e39d0SSatish Balay 
645930ae53cSBarry Smith /* ----------------------------------------------------------- */
6464a2ae208SSatish Balay #undef __FUNCT__
6474a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
648dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
649930ae53cSBarry Smith {
650930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6516849ba73SBarry Smith   PetscErrorCode ierr;
65238baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
653362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
654aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
65538baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
656362ced78SSatish Balay   PetscScalar    *v,sum;
657d85afc42SSatish Balay #endif
658930ae53cSBarry Smith 
6593a40ed3dSBarry Smith   PetscFunctionBegin;
6603a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
661930ae53cSBarry Smith 
6621ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6631ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
664930ae53cSBarry Smith 
665aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6661c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6671c4feecaSSatish Balay #else
668930ae53cSBarry Smith   /* forward solve the lower triangular */
669930ae53cSBarry Smith   x[0] = b[0];
670930ae53cSBarry Smith   for (i=1; i<n; i++) {
671930ae53cSBarry Smith     ai_i = ai[i];
672930ae53cSBarry Smith     v    = aa + ai_i;
673930ae53cSBarry Smith     vi   = aj + ai_i;
674930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
675930ae53cSBarry Smith     sum  = b[i];
676930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
677930ae53cSBarry Smith     x[i] = sum;
678930ae53cSBarry Smith   }
679930ae53cSBarry Smith 
680930ae53cSBarry Smith   /* backward solve the upper triangular */
681930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
682930ae53cSBarry Smith     adiag_i = adiag[i];
683d708749eSBarry Smith     v       = aa + adiag_i + 1;
684d708749eSBarry Smith     vi      = aj + adiag_i + 1;
685930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
686930ae53cSBarry Smith     sum     = x[i];
687930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
688930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
689930ae53cSBarry Smith   }
6901c4feecaSSatish Balay #endif
691b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6921ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6931ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
695930ae53cSBarry Smith }
696930ae53cSBarry Smith 
6974a2ae208SSatish Balay #undef __FUNCT__
6984a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
699dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
700da3a660dSBarry Smith {
701416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
702416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7036849ba73SBarry Smith   PetscErrorCode ierr;
70438baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
70538baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
70687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
707da3a660dSBarry Smith 
7083a40ed3dSBarry Smith   PetscFunctionBegin;
70978b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
710da3a660dSBarry Smith 
7111ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7121ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
713416022c9SBarry Smith   tmp  = a->solve_work;
714da3a660dSBarry Smith 
7153b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7163b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
717da3a660dSBarry Smith 
718da3a660dSBarry Smith   /* forward solve the lower triangular */
719da3a660dSBarry Smith   tmp[0] = b[*r++];
720da3a660dSBarry Smith   for (i=1; i<n; i++) {
721010a20caSHong Zhang     v   = aa + ai[i] ;
722010a20caSHong Zhang     vi  = aj + ai[i] ;
723416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
724da3a660dSBarry Smith     sum = b[*r++];
725010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
726da3a660dSBarry Smith     tmp[i] = sum;
727da3a660dSBarry Smith   }
728da3a660dSBarry Smith 
729da3a660dSBarry Smith   /* backward solve the upper triangular */
730da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
731010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
732010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
733416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
734da3a660dSBarry Smith     sum = tmp[i];
735010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
736010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
737da3a660dSBarry Smith     x[*c--] += tmp[i];
738da3a660dSBarry Smith   }
739da3a660dSBarry Smith 
7403b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7413b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7421ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7431ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
744b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
745898183ecSLois Curfman McInnes 
7463a40ed3dSBarry Smith   PetscFunctionReturn(0);
747da3a660dSBarry Smith }
748da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7494a2ae208SSatish Balay #undef __FUNCT__
7504a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
751dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
752da3a660dSBarry Smith {
753416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7542235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7556849ba73SBarry Smith   PetscErrorCode ierr;
75638baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
75738baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
75887828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
759da3a660dSBarry Smith 
7603a40ed3dSBarry Smith   PetscFunctionBegin;
7611ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7621ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
763416022c9SBarry Smith   tmp  = a->solve_work;
764da3a660dSBarry Smith 
7652235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7662235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
767da3a660dSBarry Smith 
768da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7692235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
770da3a660dSBarry Smith 
771da3a660dSBarry Smith   /* forward solve the U^T */
772da3a660dSBarry Smith   for (i=0; i<n; i++) {
773010a20caSHong Zhang     v   = aa + diag[i] ;
774010a20caSHong Zhang     vi  = aj + diag[i] + 1;
775f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
776c38d4ed2SBarry Smith     s1  = tmp[i];
777ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
778da3a660dSBarry Smith     while (nz--) {
779010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
780da3a660dSBarry Smith     }
781c38d4ed2SBarry Smith     tmp[i] = s1;
782da3a660dSBarry Smith   }
783da3a660dSBarry Smith 
784da3a660dSBarry Smith   /* backward solve the L^T */
785da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
786010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
787010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
788f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
789f1af5d2fSBarry Smith     s1  = tmp[i];
790da3a660dSBarry Smith     while (nz--) {
791010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
792da3a660dSBarry Smith     }
793da3a660dSBarry Smith   }
794da3a660dSBarry Smith 
795da3a660dSBarry Smith   /* copy tmp into x according to permutation */
796da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
797da3a660dSBarry Smith 
7982235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7992235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8001ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8011ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
802da3a660dSBarry Smith 
803b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8043a40ed3dSBarry Smith   PetscFunctionReturn(0);
805da3a660dSBarry Smith }
806da3a660dSBarry Smith 
8074a2ae208SSatish Balay #undef __FUNCT__
8084a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
809dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
810da3a660dSBarry Smith {
811416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8122235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8136849ba73SBarry Smith   PetscErrorCode ierr;
81438baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
81538baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
81687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8176abc6512SBarry Smith 
8183a40ed3dSBarry Smith   PetscFunctionBegin;
81923bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8206abc6512SBarry Smith 
8211ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
823416022c9SBarry Smith   tmp = a->solve_work;
8246abc6512SBarry Smith 
8252235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8262235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8276abc6512SBarry Smith 
8286abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8292235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8306abc6512SBarry Smith 
8316abc6512SBarry Smith   /* forward solve the U^T */
8326abc6512SBarry Smith   for (i=0; i<n; i++) {
833010a20caSHong Zhang     v   = aa + diag[i] ;
834010a20caSHong Zhang     vi  = aj + diag[i] + 1;
835f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8366abc6512SBarry Smith     tmp[i] *= *v++;
8376abc6512SBarry Smith     while (nz--) {
838010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8396abc6512SBarry Smith     }
8406abc6512SBarry Smith   }
8416abc6512SBarry Smith 
8426abc6512SBarry Smith   /* backward solve the L^T */
8436abc6512SBarry Smith   for (i=n-1; i>=0; i--){
844010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
845010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
846f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8476abc6512SBarry Smith     while (nz--) {
848010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8496abc6512SBarry Smith     }
8506abc6512SBarry Smith   }
8516abc6512SBarry Smith 
8526abc6512SBarry Smith   /* copy tmp into x according to permutation */
8536abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8546abc6512SBarry Smith 
8552235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8562235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8571ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8596abc6512SBarry Smith 
860b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
862da3a660dSBarry Smith }
8639e25ed09SBarry Smith /* ----------------------------------------------------------------*/
864dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
86508480c60SBarry Smith 
8664a2ae208SSatish Balay #undef __FUNCT__
8674a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
868dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8699e25ed09SBarry Smith {
870416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8719e25ed09SBarry Smith   IS             isicol;
8726849ba73SBarry Smith   PetscErrorCode ierr;
873*5a8e39fbSHong Zhang   PetscInt       *r,*ic,n=A->m,*ai=a->i,*aj=a->j;
874*5a8e39fbSHong Zhang   PetscInt       *bi,*cols,nnz,*cols_lvl;
875*5a8e39fbSHong Zhang   PetscInt       *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
876*5a8e39fbSHong Zhang   PetscInt       i,levels,diagonal_fill;
87777c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
878329f5518SBarry Smith   PetscReal      f;
879*5a8e39fbSHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
880*5a8e39fbSHong Zhang   PetscBT        lnkbt;
881*5a8e39fbSHong Zhang   PetscInt       nzi,*bj,**bj_ptr,**bjlvl_ptr;
882*5a8e39fbSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
883*5a8e39fbSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
8849e25ed09SBarry Smith 
8853a40ed3dSBarry Smith   PetscFunctionBegin;
8866cf997b0SBarry Smith   f             = info->fill;
88738baddfdSBarry Smith   levels        = (PetscInt)info->levels;
88838baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8894c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
89082bf6240SBarry Smith 
89108480c60SBarry Smith   /* special case that simply copies fill pattern */
892be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
893be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
89486bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
8952e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
89608480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
89708480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
89808480c60SBarry Smith     if (!b->diag) {
8997c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90008480c60SBarry Smith     }
9017c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90208480c60SBarry Smith     b->row              = isrow;
90308480c60SBarry Smith     b->col              = iscol;
90482bf6240SBarry Smith     b->icol             = isicol;
90587828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
906f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
907c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
908c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9093a40ed3dSBarry Smith     PetscFunctionReturn(0);
91008480c60SBarry Smith   }
91108480c60SBarry Smith 
912898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
913898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9149e25ed09SBarry Smith 
9159e25ed09SBarry Smith   /* get new row pointers */
916*5a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
917*5a8e39fbSHong Zhang   bi[0] = 0;
918*5a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
919*5a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
920*5a8e39fbSHong Zhang   bdiag[0]  = 0;
9216cf997b0SBarry Smith 
922*5a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
923*5a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
924*5a8e39fbSHong Zhang 
925*5a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
926*5a8e39fbSHong Zhang   nlnk = n + 1;
927*5a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
928*5a8e39fbSHong Zhang 
929*5a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
930*5a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
931*5a8e39fbSHong Zhang   current_space = free_space;
932*5a8e39fbSHong Zhang   ierr = GetMoreSpace((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
933*5a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
934*5a8e39fbSHong Zhang 
935*5a8e39fbSHong Zhang   for (i=0; i<n; i++) {
936*5a8e39fbSHong Zhang     nzi = 0;
9376cf997b0SBarry Smith     /* copy current row into linked list */
938*5a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
939*5a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
940*5a8e39fbSHong Zhang     cols = aj + ai[r[i]];
941*5a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
942*5a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
943*5a8e39fbSHong Zhang     nzi += nlnk;
9446cf997b0SBarry Smith 
9456cf997b0SBarry Smith     /* make sure diagonal entry is included */
946*5a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
9476cf997b0SBarry Smith       fm = n;
948*5a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
949*5a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
950*5a8e39fbSHong Zhang       lnk[fm]    = i;
951*5a8e39fbSHong Zhang       lnk_lvl[i] = 0;
952*5a8e39fbSHong Zhang       nzi++; dcount++;
9536cf997b0SBarry Smith     }
9546cf997b0SBarry Smith 
955*5a8e39fbSHong Zhang     /* add pivot rows into the active row */
956*5a8e39fbSHong Zhang     nzbd = 0;
957*5a8e39fbSHong Zhang     prow = lnk[n];
958*5a8e39fbSHong Zhang     while (prow < i) {
959*5a8e39fbSHong Zhang       nnz      = bdiag[prow];
960*5a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
961*5a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
962*5a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
963*5a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
964*5a8e39fbSHong Zhang       nzi += nlnk;
965*5a8e39fbSHong Zhang       prow = lnk[prow];
966*5a8e39fbSHong Zhang       nzbd++;
96756beaf04SBarry Smith     }
968*5a8e39fbSHong Zhang     bdiag[i] = nzbd;
969*5a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
970ecf371e4SBarry Smith 
971*5a8e39fbSHong Zhang     /* if free space is not available, make more free space */
972*5a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
973*5a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
974*5a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space);CHKERRQ(ierr);
975*5a8e39fbSHong Zhang       ierr = GetMoreSpace(nnz,&current_space_lvl);CHKERRQ(ierr);
976*5a8e39fbSHong Zhang       reallocs++;
977*5a8e39fbSHong Zhang     }
978ecf371e4SBarry Smith 
979*5a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
980*5a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
981*5a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
982*5a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
983*5a8e39fbSHong Zhang 
984*5a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
985*5a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
98677431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
987*5a8e39fbSHong Zhang     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i);
9886cf997b0SBarry Smith     }
989*5a8e39fbSHong Zhang 
990*5a8e39fbSHong Zhang     current_space->array           += nzi;
991*5a8e39fbSHong Zhang     current_space->local_used      += nzi;
992*5a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
993*5a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
994*5a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
995*5a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
9969e25ed09SBarry Smith   }
997*5a8e39fbSHong Zhang 
998898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
999898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1000*5a8e39fbSHong Zhang 
1001*5a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
1002*5a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1003*5a8e39fbSHong Zhang   ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
1004*5a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1005*5a8e39fbSHong Zhang   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
1006*5a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10079e25ed09SBarry Smith 
1008f64065bbSBarry Smith   {
1009*5a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1010418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1011c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1012c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1013b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1014335d9088SBarry Smith     if (diagonal_fill) {
101577431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1016335d9088SBarry Smith     }
1017f64065bbSBarry Smith   }
1018bbb6d6a8SBarry Smith 
10199e25ed09SBarry Smith   /* put together the new matrix */
1020f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1021f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1022f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
102352e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1024416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1025606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10267c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1027*5a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
10289e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1029606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1030606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1031b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1032*5a8e39fbSHong Zhang   b->j          = bj;
1033*5a8e39fbSHong Zhang   b->i          = bi;
1034*5a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
1035*5a8e39fbSHong Zhang   b->diag       = bdiag;
1036416022c9SBarry Smith   b->ilen       = 0;
1037416022c9SBarry Smith   b->imax       = 0;
1038416022c9SBarry Smith   b->row        = isrow;
1039416022c9SBarry Smith   b->col        = iscol;
1040c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1041c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
104282bf6240SBarry Smith   b->icol       = isicol;
104387828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1044416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
1045*5a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
1046*5a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
1047*5a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
10489e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10493a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10503a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1051418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1052ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1053*5a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
10543a40ed3dSBarry Smith   PetscFunctionReturn(0);
10559e25ed09SBarry Smith }
1056d5d45c9bSBarry Smith 
10573c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1058a6175056SHong Zhang #undef __FUNCT__
1059a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1060af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1061a6175056SHong Zhang {
10622ed133dbSHong Zhang   Mat            C = *B;
10630968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
10642ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
10652ed133dbSHong Zhang   IS             ip=b->row;
10662ed133dbSHong Zhang   PetscErrorCode ierr;
10672ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
10682ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
10692ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1070622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1071540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1072540703acSHong Zhang   PetscTruth     shiftpd;
1073540703acSHong Zhang   ChShift_Ctx    sctx;
1074540703acSHong Zhang   PetscInt       newshift;
1075d5d45c9bSBarry Smith 
1076a6175056SHong Zhang   PetscFunctionBegin;
1077540703acSHong Zhang   shiftnz   = info->shiftnz;
1078540703acSHong Zhang   shiftpd   = info->shiftpd;
1079ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1080ee45ca4aSHong Zhang 
10812ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
10822ed133dbSHong Zhang 
10832ed133dbSHong Zhang   /* initialization */
10842ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
10852ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
10862ed133dbSHong Zhang   jl   = il + mbs;
10872ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
10882ed133dbSHong Zhang 
1089540703acSHong Zhang   sctx.shift_amount = 0;
1090540703acSHong Zhang   sctx.nshift       = 0;
10912ed133dbSHong Zhang   do {
1092540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
10932ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
10942ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1095a6175056SHong Zhang     }
10962ed133dbSHong Zhang 
10972ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1098c05c3958SHong Zhang       bval = ba + bi[k];
10992ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11002ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11012ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11022ed133dbSHong Zhang         col = rip[aj[j]];
11032ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11042ed133dbSHong Zhang           rtmp[col] = aa[j];
1105c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11062ed133dbSHong Zhang         }
11072ed133dbSHong Zhang       }
110839e3d630SHong Zhang       /* shift the diagonal of the matrix */
1109540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11102ed133dbSHong Zhang 
11112ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11122ed133dbSHong Zhang       dk = rtmp[k];
11132ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11142ed133dbSHong Zhang 
11152ed133dbSHong Zhang       while (i < k){
11162ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11172ed133dbSHong Zhang 
11182ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11192ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11202ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11212ed133dbSHong Zhang         dk += uikdi*ba[ili];
11222ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11232ed133dbSHong Zhang 
11242ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11252ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11262ed133dbSHong Zhang         if (jmin < jmax){
11272ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11282ed133dbSHong Zhang           /* update il and jl for row i */
11292ed133dbSHong Zhang           il[i] = jmin;
11302ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11312ed133dbSHong Zhang         }
11322ed133dbSHong Zhang         i = nexti;
11332ed133dbSHong Zhang       }
11342ed133dbSHong Zhang 
11350697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11360697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11370697b6caSHong Zhang       rs   = 0.0;
11383c9e8598SHong Zhang       jmin = bi[k]+1;
11393c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11403c9e8598SHong Zhang       if (nz){
11413c9e8598SHong Zhang         bcol = bj + jmin;
11423c9e8598SHong Zhang         while (nz--){
114389f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
114489f845c8SHong Zhang           bcol++;
11453c9e8598SHong Zhang         }
11463c9e8598SHong Zhang       }
1147540703acSHong Zhang 
1148540703acSHong Zhang       sctx.rs = rs;
1149540703acSHong Zhang       sctx.pv = dk;
1150540703acSHong Zhang       ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr);
1151540703acSHong Zhang       if (newshift == 1){
1152540703acSHong Zhang         break;    /* sctx.shift_amount is updated */
1153540703acSHong Zhang       } else if (newshift == -1){
11540697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11553c9e8598SHong Zhang       }
11563c9e8598SHong Zhang 
11573c9e8598SHong Zhang       /* copy data into U(k,:) */
115839e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
115939e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
116039e3d630SHong Zhang       if (jmin < jmax) {
116139e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
116239e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
11633c9e8598SHong Zhang         }
116439e3d630SHong Zhang         /* add the k-th row into il and jl */
11653c9e8598SHong Zhang         il[k] = jmin;
11663c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
11673c9e8598SHong Zhang       }
11683c9e8598SHong Zhang     }
1169540703acSHong Zhang   } while (sctx.chshift);
11703c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
11713c9e8598SHong Zhang 
117239e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
11733c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
11743c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
11753c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
11763c9e8598SHong Zhang   PetscLogFlops(C->m);
1177540703acSHong Zhang   if (sctx.nshift){
1178540703acSHong Zhang     if (shiftnz) {
1179540703acSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
1180540703acSHong Zhang     } else if (shiftpd) {
1181540703acSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount);
11820697b6caSHong Zhang     }
11833c9e8598SHong Zhang   }
11843c9e8598SHong Zhang   PetscFunctionReturn(0);
11853c9e8598SHong Zhang }
1186a6175056SHong Zhang 
1187a6175056SHong Zhang #undef __FUNCT__
1188a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1189dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1190a6175056SHong Zhang {
11910968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1192ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1193ed59401aSHong Zhang   Mat            B;
1194dfbe8321SBarry Smith   PetscErrorCode ierr;
1195653a6975SHong Zhang   PetscTruth     perm_identity;
1196622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1197ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1198391eac55SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL;
1199*5a8e39fbSHong Zhang   PetscInt       ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1200ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
12017a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
12027a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12037a48dd6fSHong Zhang   PetscBT        lnkbt;
1204a6175056SHong Zhang 
1205a6175056SHong Zhang   PetscFunctionBegin;
1206653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12077a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12087a48dd6fSHong Zhang 
120939e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
121039e3d630SHong Zhang   ui[0] = 0;
121139e3d630SHong Zhang 
1212622e2028SHong Zhang   /* special case that simply copies fill pattern */
1213622e2028SHong Zhang   if (!levels && perm_identity) {
1214622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1215ed59401aSHong Zhang     for (i=0; i<am; i++) {
121639e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1217ed59401aSHong Zhang     }
121839e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
121939e3d630SHong Zhang     cols = uj;
1220ed59401aSHong Zhang     for (i=0; i<am; i++) {
1221ed59401aSHong Zhang       aj    = a->j + a->diag[i];
122239e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
122339e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1224ed59401aSHong Zhang     }
122539e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12267a48dd6fSHong Zhang     /* initialization */
1227*5a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12287a48dd6fSHong Zhang 
12297a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12307a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12311393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12327a48dd6fSHong Zhang     il         = jl + am;
12337a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12347a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12357a48dd6fSHong Zhang     for (i=0; i<am; i++){
12367a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12377a48dd6fSHong Zhang     }
12387a48dd6fSHong Zhang 
12397a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12407a48dd6fSHong Zhang     nlnk = am + 1;
12412ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12427a48dd6fSHong Zhang 
12437a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12447a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12457a48dd6fSHong Zhang     current_space = free_space;
12467a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12477a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12487a48dd6fSHong Zhang 
12497a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12507a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12517a48dd6fSHong Zhang       nzk   = 0;
1252622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1253622e2028SHong Zhang       ncols_upper = 0;
1254622e2028SHong Zhang       for (j=0; j<ncols; j++){
1255*5a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
1256*5a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
1257*5a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1258622e2028SHong Zhang           ncols_upper++;
1259622e2028SHong Zhang         }
1260622e2028SHong Zhang       }
1261*5a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12627a48dd6fSHong Zhang       nzk += nlnk;
12637a48dd6fSHong Zhang 
12647a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
12657a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
12667a48dd6fSHong Zhang 
12677a48dd6fSHong Zhang       while (prow < k){
12687a48dd6fSHong Zhang         nextprow = jl[prow];
12697a48dd6fSHong Zhang 
12707a48dd6fSHong Zhang         /* merge prow into k-th row */
12717a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
12727a48dd6fSHong Zhang         jmax = ui[prow+1];
12737a48dd6fSHong Zhang         ncols = jmax-jmin;
1274ed59401aSHong Zhang         i     = jmin - ui[prow];
1275ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1276*5a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
1277*5a8e39fbSHong Zhang         j     = *(uj - 1);
1278*5a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
12797a48dd6fSHong Zhang         nzk += nlnk;
12807a48dd6fSHong Zhang 
12817a48dd6fSHong Zhang         /* update il and jl for prow */
12827a48dd6fSHong Zhang         if (jmin < jmax){
12837a48dd6fSHong Zhang           il[prow] = jmin;
12847a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
12857a48dd6fSHong Zhang         }
12867a48dd6fSHong Zhang         prow = nextprow;
12877a48dd6fSHong Zhang       }
12887a48dd6fSHong Zhang 
12897a48dd6fSHong Zhang       /* if free space is not available, make more free space */
12907a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
12917a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
12927a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
12937a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
12947a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
12957a48dd6fSHong Zhang         reallocs++;
12967a48dd6fSHong Zhang       }
12977a48dd6fSHong Zhang 
12987a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
12992ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13007a48dd6fSHong Zhang 
13017a48dd6fSHong Zhang       /* add the k-th row into il and jl */
130239e3d630SHong Zhang       if (nzk > 1){
13037a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13047a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13057a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13067a48dd6fSHong Zhang       }
13077a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13087a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13097a48dd6fSHong Zhang 
13107a48dd6fSHong Zhang       current_space->array           += nzk;
13117a48dd6fSHong Zhang       current_space->local_used      += nzk;
13127a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13137a48dd6fSHong Zhang 
13147a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13157a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13167a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13177a48dd6fSHong Zhang 
13187a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13197a48dd6fSHong Zhang     }
13207a48dd6fSHong Zhang 
13217a48dd6fSHong Zhang     if (ai[am] != 0) {
132239e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
13237a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
13247a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
13257a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
13267a48dd6fSHong Zhang     } else {
1327ed59401aSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
13287a48dd6fSHong Zhang     }
13297a48dd6fSHong Zhang 
13307a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13317a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
1332*5a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13337a48dd6fSHong Zhang 
13347a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13357a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13367a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13372ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13387a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13397a48dd6fSHong Zhang 
134039e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
134139e3d630SHong Zhang 
13427a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13437a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1344ed59401aSHong Zhang   B = *fact;
1345ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1346ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
13477a48dd6fSHong Zhang 
1348ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13497a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
13507a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13517a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
13527a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
13537a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
13547a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13557a48dd6fSHong Zhang   b->j    = uj;
13567a48dd6fSHong Zhang   b->i    = ui;
13577a48dd6fSHong Zhang   b->diag = 0;
13587a48dd6fSHong Zhang   b->ilen = 0;
13597a48dd6fSHong Zhang   b->imax = 0;
13607a48dd6fSHong Zhang   b->row  = perm;
13617a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
13627a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13637a48dd6fSHong Zhang   b->icol = perm;
13647a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13657a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
136652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13677a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
13687a48dd6fSHong Zhang 
1369ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1370ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1371ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
13727a48dd6fSHong Zhang   if (ai[am] != 0) {
1373ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
13747a48dd6fSHong Zhang   } else {
1375ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
13767a48dd6fSHong Zhang   }
137739e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1378622e2028SHong Zhang   if (perm_identity){
1379ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1380ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1381622e2028SHong Zhang   }
1382a6175056SHong Zhang   PetscFunctionReturn(0);
1383a6175056SHong Zhang }
1384d5d45c9bSBarry Smith 
1385f76d2b81SHong Zhang #undef __FUNCT__
1386f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1387dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1388f76d2b81SHong Zhang {
1389f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
139010c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
13912ed133dbSHong Zhang   Mat            B;
1392dfbe8321SBarry Smith   PetscErrorCode ierr;
1393f76d2b81SHong Zhang   PetscTruth     perm_identity;
139410c27e3dSHong Zhang   PetscReal      fill = info->fill;
13951393bc97SHong Zhang   PetscInt       *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
139610c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
13972ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
139810c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
139910c27e3dSHong Zhang   PetscBT        lnkbt;
14002ed133dbSHong Zhang   IS             iperm;
1401f76d2b81SHong Zhang 
1402f76d2b81SHong Zhang   PetscFunctionBegin;
14032ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1404f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14052ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14062ed133dbSHong Zhang 
1407f76d2b81SHong Zhang   if (!perm_identity){
14082ed133dbSHong Zhang     /* check if perm is symmetric! */
14092ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14102ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14111393bc97SHong Zhang     for (i=0; i<am; i++) {
14122ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14132ed133dbSHong Zhang     }
14142ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14152ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1416f76d2b81SHong Zhang   }
141710c27e3dSHong Zhang 
141810c27e3dSHong Zhang   /* initialization */
14191393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
142010c27e3dSHong Zhang   ui[0] = 0;
142110c27e3dSHong Zhang 
142210c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14231393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14241393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14251393bc97SHong Zhang   il     = jl + am;
14261393bc97SHong Zhang   cols   = il + am;
14271393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14281393bc97SHong Zhang   for (i=0; i<am; i++){
14291393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1430f76d2b81SHong Zhang   }
143110c27e3dSHong Zhang 
143210c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14331393bc97SHong Zhang   nlnk = am + 1;
14341393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
143510c27e3dSHong Zhang 
14361393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
14371393bc97SHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
143810c27e3dSHong Zhang   current_space = free_space;
143910c27e3dSHong Zhang 
14401393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
144110c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
144210c27e3dSHong Zhang     nzk   = 0;
14432ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14442ed133dbSHong Zhang     ncols_upper = 0;
14452ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1446622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14472ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14482ed133dbSHong Zhang         cols[ncols_upper] = i;
14492ed133dbSHong Zhang         ncols_upper++;
14502ed133dbSHong Zhang       }
14512ed133dbSHong Zhang     }
14521393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
145310c27e3dSHong Zhang     nzk += nlnk;
145410c27e3dSHong Zhang 
145510c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
145610c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
145710c27e3dSHong Zhang 
145810c27e3dSHong Zhang     while (prow < k){
145910c27e3dSHong Zhang       nextprow = jl[prow];
146010c27e3dSHong Zhang       /* merge prow into k-th row */
14611393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
146210c27e3dSHong Zhang       jmax = ui[prow+1];
146310c27e3dSHong Zhang       ncols = jmax-jmin;
14641393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1465*5a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
146610c27e3dSHong Zhang       nzk += nlnk;
146710c27e3dSHong Zhang 
146810c27e3dSHong Zhang       /* update il and jl for prow */
146910c27e3dSHong Zhang       if (jmin < jmax){
147010c27e3dSHong Zhang         il[prow] = jmin;
14712ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
147210c27e3dSHong Zhang       }
147310c27e3dSHong Zhang       prow = nextprow;
147410c27e3dSHong Zhang     }
147510c27e3dSHong Zhang 
147610c27e3dSHong Zhang     /* if free space is not available, make more free space */
147710c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
14781393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
147910c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
148010c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
148110c27e3dSHong Zhang       reallocs++;
148210c27e3dSHong Zhang     }
148310c27e3dSHong Zhang 
148410c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
14851393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
148610c27e3dSHong Zhang 
148710c27e3dSHong Zhang     /* add the k-th row into il and jl */
148810c27e3dSHong Zhang     if (nzk-1 > 0){
14891393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
149010c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
149110c27e3dSHong Zhang       il[k] = ui[k] + 1;
149210c27e3dSHong Zhang     }
14932ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
149410c27e3dSHong Zhang     current_space->array           += nzk;
149510c27e3dSHong Zhang     current_space->local_used      += nzk;
149610c27e3dSHong Zhang     current_space->local_remaining -= nzk;
149710c27e3dSHong Zhang 
149810c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
149910c27e3dSHong Zhang   }
150010c27e3dSHong Zhang 
15011393bc97SHong Zhang   if (ai[am] != 0) {
15021393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
150378910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
150478910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
150578910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
150610c27e3dSHong Zhang   } else {
150778910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
150810c27e3dSHong Zhang   }
150910c27e3dSHong Zhang 
151010c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
151110c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
151210c27e3dSHong Zhang 
151310c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15141393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
151510c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
151610c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
151710c27e3dSHong Zhang 
151810c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15191393bc97SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
15202ed133dbSHong Zhang   B    = *fact;
15212ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
15222ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
152310c27e3dSHong Zhang 
15242ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
152510c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
152610c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
152710c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
152810c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
152910c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15301393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
153110c27e3dSHong Zhang   b->j    = uj;
153210c27e3dSHong Zhang   b->i    = ui;
153310c27e3dSHong Zhang   b->diag = 0;
153410c27e3dSHong Zhang   b->ilen = 0;
153510c27e3dSHong Zhang   b->imax = 0;
153610c27e3dSHong Zhang   b->row  = perm;
153710c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
153810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
153910c27e3dSHong Zhang   b->icol = perm;
154010c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15411393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15421393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
15431393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
154410c27e3dSHong Zhang 
15452ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15462ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15472ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15481393bc97SHong Zhang   if (ai[am] != 0) {
15491393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
155010c27e3dSHong Zhang   } else {
15512ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
155210c27e3dSHong Zhang   }
155339e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15542ed133dbSHong Zhang   if (perm_identity){
155510c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
155610c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15572ed133dbSHong Zhang   }
1558f76d2b81SHong Zhang   PetscFunctionReturn(0);
1559f76d2b81SHong Zhang }
1560