xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 35bd34fa3344fb53ae05f23bf99eddd2f9bb15e6)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*);
26e34fafa9SBarry Smith #endif
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 */
510481f469SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
635d0c19d7SBarry Smith   const PetscInt *c,*r,*ic;
645d0c19d7SBarry Smith   PetscInt       i,n = A->rmap->n,*cc,*cr;
6538baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6638baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6738baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6854f21887SBarry Smith   MatScalar      *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
690481f469SBarry Smith   PetscReal      af,dt,dtcount,dtcol,fill;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   PetscFunctionBegin;
7286bacbd2SBarry Smith 
730481f469SBarry Smith   if (info->dt == PETSC_DEFAULT)      dt      = .005; else dt = info->dt;
740481f469SBarry Smith   if (info->dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax);  else dtcount = info->dtcount;
750481f469SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   dtcol   = .01; else dtcol = info->dtcol;
760481f469SBarry Smith   if (info->fill == PETSC_DEFAULT)    fill    = ((double)(n*(info->dtcount+1)))/a->nz; else fill = info->fill;
770481f469SBarry Smith   lfill   = (PetscInt)(dtcount/2.0);
780481f469SBarry Smith   jmax    = (PetscInt)(fill*a->nz);
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith 
8186bacbd2SBarry Smith   /* ------------------------------------------------------------
8286bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8386bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8486bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8586bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8686bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8786bacbd2SBarry Smith      ------------------------------------------------------------
8886bacbd2SBarry Smith   */
8986bacbd2SBarry Smith 
9086bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9186bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9286bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9386bacbd2SBarry Smith   reorder = PetscNot(reorder);
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith 
9686bacbd2SBarry Smith   /* storage for ilu factor */
9738baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9838baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
99a77337e4SBarry Smith   ierr = PetscMalloc(jmax*sizeof(MatScalar),&new_a);CHKERRQ(ierr);
10038baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10186bacbd2SBarry Smith 
10286bacbd2SBarry Smith   /* ------------------------------------------------------------
10386bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10486bacbd2SBarry Smith      ------------------------------------------------------------
10586bacbd2SBarry Smith   */
106b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
107b19713cbSBarry Smith     old_j[i]++;
10886bacbd2SBarry Smith   }
109b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
110b19713cbSBarry Smith     old_i[i]++;
111b19713cbSBarry Smith   };
112010a20caSHong Zhang 
11386bacbd2SBarry Smith 
11486bacbd2SBarry Smith   if (reorder) {
115c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
116c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1175d0c19d7SBarry Smith     ierr = PetscMalloc2(n,PetscInt,&cc,n,PetscInt,&cr);CHKERRQ(ierr);
118c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
1195d0c19d7SBarry Smith       cr[i]  = r[i]+1;
1205d0c19d7SBarry Smith       cc[i]  = c[i]+1;
121c0c2c81eSBarry Smith     }
122c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
123c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1245d0c19d7SBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
1255d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
1265d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr);
1275d0c19d7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,cr,cc,&job);
1285d0c19d7SBarry Smith     ierr = PetscFree2(cc,cr);CHKERRQ(ierr);
129b19713cbSBarry Smith     o_a = old_a2;
130b19713cbSBarry Smith     o_j = old_j2;
131b19713cbSBarry Smith     o_i = old_i2;
132b19713cbSBarry Smith   } else {
133b19713cbSBarry Smith     o_a = old_a;
134b19713cbSBarry Smith     o_j = old_j;
135b19713cbSBarry Smith     o_i = old_i;
13686bacbd2SBarry Smith   }
13786bacbd2SBarry Smith 
13886bacbd2SBarry Smith   /* ------------------------------------------------------------
13986bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14086bacbd2SBarry Smith      ------------------------------------------------------------
14186bacbd2SBarry Smith   */
14286bacbd2SBarry Smith 
14338baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14587828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14686bacbd2SBarry Smith 
1470481f469SBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)dt,&dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1486849ba73SBarry Smith   if (sierr) {
1496849ba73SBarry Smith     switch (sierr) {
1500481f469SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger fill current fill %G space allocated %D",fill,jmax);
1510481f469SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger fill current fill %G space allocated %D",fill,jmax);
152e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
153e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
1540481f469SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal fill value %D",jmax);
15577431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15686bacbd2SBarry Smith     }
15786bacbd2SBarry Smith   }
15886bacbd2SBarry Smith 
15986bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16086bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   /* ------------------------------------------------------------
16386bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16486bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16586bacbd2SBarry Smith      ------------------------------------------------------------
16686bacbd2SBarry Smith   */
16786bacbd2SBarry Smith 
16887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16938baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17086bacbd2SBarry Smith 
1715bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17286bacbd2SBarry Smith 
17386bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17486bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   if (reorder) {
17786bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
180313ae024SBarry Smith   } else {
181b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
182f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
183b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
184b19713cbSBarry Smith     }
185b19713cbSBarry Smith   }
186b19713cbSBarry Smith 
187b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18886bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
189b19713cbSBarry Smith     old_i[i]--;
190b19713cbSBarry Smith   }
191b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
192b19713cbSBarry Smith     old_j[i]--;
193b19713cbSBarry Smith   }
19486bacbd2SBarry Smith 
195b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19686bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
197b19713cbSBarry Smith     new_i[i]--;
198b19713cbSBarry Smith   }
199b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
200b19713cbSBarry Smith     new_j[i]--;
201b19713cbSBarry Smith   }
20286bacbd2SBarry Smith 
20386bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20486bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2054c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2064c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
207c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20886bacbd2SBarry Smith   for(i=0; i<n; i++) {
20986bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21086bacbd2SBarry Smith   };
21186bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21207d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21386bacbd2SBarry Smith 
214c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
215c0c2c81eSBarry Smith 
216329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21707d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21886bacbd2SBarry Smith 
21986bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22086bacbd2SBarry Smith 
2217adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
222f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
2237adad957SLisandro Dalcin   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
224ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2255c9eb25fSBarry Smith   (*fact)->factor    = MAT_FACTOR_LU;
22686bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22786bacbd2SBarry Smith 
22886bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
229e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
230e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
23107d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23286bacbd2SBarry Smith   b->a             = new_a;
23386bacbd2SBarry Smith   b->j             = new_j;
23486bacbd2SBarry Smith   b->i             = new_i;
23586bacbd2SBarry Smith   b->ilen          = 0;
23686bacbd2SBarry Smith   b->imax          = 0;
2371f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
238313ae024SBarry Smith   b->row           = isirow;
23986bacbd2SBarry Smith   b->col           = iscolf;
24087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24186bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24286bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24386bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24486bacbd2SBarry Smith 
245431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
2460481f469SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",fill,af);CHKERRQ(ierr);
247ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
248ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
250431e710aSBarry Smith 
2517529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25271c2f376SKris Buschelman 
25386bacbd2SBarry Smith   PetscFunctionReturn(0);
25460ec86bdSBarry Smith #endif
25586bacbd2SBarry Smith }
25686bacbd2SBarry Smith 
257e631078cSBarry Smith EXTERN_C_BEGIN
2584a2ae208SSatish Balay #undef __FUNCT__
259db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
260db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
261db4efbfdSBarry Smith {
262db4efbfdSBarry Smith   PetscFunctionBegin;
263db4efbfdSBarry Smith   *flg = PETSC_TRUE;
264db4efbfdSBarry Smith   PetscFunctionReturn(0);
265db4efbfdSBarry Smith }
266db4efbfdSBarry Smith EXTERN_C_END
267db4efbfdSBarry Smith 
268db4efbfdSBarry Smith EXTERN_C_BEGIN
269db4efbfdSBarry Smith #undef __FUNCT__
270b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc"
271b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
272b24902e0SBarry Smith {
273d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
274b24902e0SBarry Smith   PetscErrorCode     ierr;
275b24902e0SBarry Smith 
276b24902e0SBarry Smith   PetscFunctionBegin;
277b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
278b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
279b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
280b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
281*35bd34faSBarry Smith     if (ftype == MAT_FACTOR_ILU) {
282db4efbfdSBarry Smith       (*B)->ops->ilufactorsymbolic= MatILUFactorSymbolic_SeqAIJ;
283*35bd34faSBarry Smith     } else {
284*35bd34faSBarry Smith       (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
285*35bd34faSBarry Smith     }
286b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
287b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
288b24902e0SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
289*35bd34faSBarry Smith     if (ftype == MAT_FACTOR_ICC) {
2905c9eb25fSBarry Smith       (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
291*35bd34faSBarry Smith     } else {
2925c9eb25fSBarry Smith       (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
293*35bd34faSBarry Smith     }
294b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
295db4efbfdSBarry Smith   (*B)->factor = ftype;
296b24902e0SBarry Smith   PetscFunctionReturn(0);
297b24902e0SBarry Smith }
298e631078cSBarry Smith EXTERN_C_END
299b24902e0SBarry Smith 
300b24902e0SBarry Smith #undef __FUNCT__
301b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
3020481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
303289bc588SBarry Smith {
304416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
305289bc588SBarry Smith   IS                 isicol;
3066849ba73SBarry Smith   PetscErrorCode     ierr;
3075d0c19d7SBarry Smith   const PetscInt     *r,*ic;
3085d0c19d7SBarry Smith   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
3091393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
3101393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
3119e046878SBarry Smith   PetscReal          f;
3124875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
313a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3141393bc97SHong Zhang   PetscBT            lnkbt;
315289bc588SBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
317d0f46423SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
3184c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
319ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
320ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
321289bc588SBarry Smith 
322289bc588SBarry Smith   /* get new row pointers */
3231393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
3241393bc97SHong Zhang   bi[0] = 0;
3251393bc97SHong Zhang 
3261393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
3271393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
3281393bc97SHong Zhang   bdiag[0] = 0;
3291393bc97SHong Zhang 
3301393bc97SHong Zhang   /* linked list for storing column indices of the active row */
3311393bc97SHong Zhang   nlnk = n + 1;
3321393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3331393bc97SHong Zhang 
33435baf27eSSatish Balay   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
3351393bc97SHong Zhang 
3361393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
337d3d32019SBarry Smith   f = info->fill;
338a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3391393bc97SHong Zhang   current_space = free_space;
340289bc588SBarry Smith 
341289bc588SBarry Smith   for (i=0; i<n; i++) {
3421393bc97SHong Zhang     /* copy previous fill into linked list */
343289bc588SBarry Smith     nzi = 0;
3441393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3451393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3461393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
347afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3481393bc97SHong Zhang     nzi += nlnk;
3491393bc97SHong Zhang 
3501393bc97SHong Zhang     /* add pivot rows into linked list */
3511393bc97SHong Zhang     row = lnk[n];
3521393bc97SHong Zhang     while (row < i) {
353d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
354d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
355d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3561393bc97SHong Zhang       nzi += nlnk;
3571393bc97SHong Zhang       row  = lnk[row];
358289bc588SBarry Smith     }
3591393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3601393bc97SHong Zhang     im[i]   = nzi;
3611393bc97SHong Zhang 
3621393bc97SHong Zhang     /* mark bdiag */
3631393bc97SHong Zhang     nzbd = 0;
3641393bc97SHong Zhang     nnz  = nzi;
3651393bc97SHong Zhang     k    = lnk[n];
3661393bc97SHong Zhang     while (nnz-- && k < i){
3671393bc97SHong Zhang       nzbd++;
3681393bc97SHong Zhang       k = lnk[k];
3691393bc97SHong Zhang     }
3701393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3711393bc97SHong Zhang 
3721393bc97SHong Zhang     /* if free space is not available, make more free space */
3731393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3741393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
375a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3761393bc97SHong Zhang       reallocs++;
3771393bc97SHong Zhang     }
3781393bc97SHong Zhang 
3791393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3801393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3811393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3821393bc97SHong Zhang     current_space->array           += nzi;
3831393bc97SHong Zhang     current_space->local_used      += nzi;
3841393bc97SHong Zhang     current_space->local_remaining -= nzi;
385289bc588SBarry Smith   }
3866cf91177SBarry Smith #if defined(PETSC_USE_INFO)
387464e8d44SSatish Balay   if (ai[n] != 0) {
3881393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
389ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
390ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
391ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
392ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
39305bf559cSSatish Balay   } else {
394ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
39505bf559cSSatish Balay   }
39663ba0a88SBarry Smith #endif
3972fb3ed76SBarry Smith 
398898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
399898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
4001987afe7SBarry Smith 
4011393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
4021393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
403a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
4041393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
40535baf27eSSatish Balay   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
406289bc588SBarry Smith 
407289bc588SBarry Smith   /* put together the new matrix */
408719d5645SBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
409719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
410719d5645SBarry Smith   b    = (Mat_SeqAIJ*)(B)->data;
411e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
412e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
4137c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
4141393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
4151393bc97SHong Zhang   b->j          = bj;
4161393bc97SHong Zhang   b->i          = bi;
4171393bc97SHong Zhang   b->diag       = bdiag;
418416022c9SBarry Smith   b->ilen       = 0;
419416022c9SBarry Smith   b->imax       = 0;
420416022c9SBarry Smith   b->row        = isrow;
421416022c9SBarry Smith   b->col        = iscol;
422c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
423c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
42482bf6240SBarry Smith   b->icol       = isicol;
42587828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
4261393bc97SHong Zhang 
4271393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
428719d5645SBarry Smith   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
4291393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4307fd98bd6SLois Curfman McInnes 
431719d5645SBarry Smith   (B)->factor                 =  MAT_FACTOR_LU;
432719d5645SBarry Smith   (B)->info.factor_mallocs    = reallocs;
433719d5645SBarry Smith   (B)->info.fill_ratio_given  = f;
434703deb49SSatish Balay 
435e93ccc4dSBarry Smith   if (ai[n] != 0) {
436719d5645SBarry Smith     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
437eea2913aSSatish Balay   } else {
438719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
439eea2913aSSatish Balay   }
440719d5645SBarry Smith   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
441719d5645SBarry Smith   (B)->ops->solve            = MatSolve_SeqAIJ;
442719d5645SBarry Smith   (B)->ops->solvetranspose   = MatSolveTranspose_SeqAIJ;
443db4efbfdSBarry Smith   /* switch to inodes if appropriate */
444719d5645SBarry Smith   ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr);
4453a40ed3dSBarry Smith   PetscFunctionReturn(0);
446289bc588SBarry Smith }
4471393bc97SHong Zhang 
4485b5bc046SBarry Smith /*
4495b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4505b5bc046SBarry Smith */
4515b5bc046SBarry Smith #undef __FUNCT__
4525b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4535b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4545b5bc046SBarry Smith {
4555b5bc046SBarry Smith   PetscErrorCode ierr;
4565b5bc046SBarry Smith   PetscTruth     flg;
4575b5bc046SBarry Smith 
4585b5bc046SBarry Smith   PetscFunctionBegin;
4595b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4605b5bc046SBarry Smith   if (flg) {
4615b5bc046SBarry Smith     PetscViewer viewer;
4625b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4635b5bc046SBarry Smith 
4645b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4657adad957SLisandro Dalcin     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4665b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4675b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4685b5bc046SBarry Smith   }
4695b5bc046SBarry Smith   PetscFunctionReturn(0);
4705b5bc046SBarry Smith }
4715b5bc046SBarry Smith 
472db4efbfdSBarry Smith extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec);
473db4efbfdSBarry Smith 
4740a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4754a2ae208SSatish Balay #undef __FUNCT__
4764a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
4770481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
478289bc588SBarry Smith {
479719d5645SBarry Smith   Mat            C=B;
480416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
48182bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4826849ba73SBarry Smith   PetscErrorCode ierr;
4835d0c19d7SBarry Smith   const PetscInt  *r,*ic,*ics;
4845d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*bi=b->i,*bj=b->j;
4855d0c19d7SBarry Smith   PetscInt       *ajtmp,*bjtmp,nz,row;
486d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
48754f21887SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,*rtmps;
48854f21887SBarry Smith   MatScalar      *v,*pv;
4896a7c8fc2SHong Zhang   PetscScalar    d;
4906a7c8fc2SHong Zhang   PetscReal      rs;
491b3bf805bSHong Zhang   LUShift_Ctx    sctx;
49242898933SHong Zhang   PetscInt       newshift,*ddiag;
493289bc588SBarry Smith 
4943a40ed3dSBarry Smith   PetscFunctionBegin;
49578b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
49678b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
49787828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
49887828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
499010a20caSHong Zhang   rtmps = rtmp; ics = ic;
500289bc588SBarry Smith 
501ada7143aSSatish Balay   sctx.shift_top      = 0;
502ada7143aSSatish Balay   sctx.nshift_max     = 0;
503ada7143aSSatish Balay   sctx.shift_lo       = 0;
504ada7143aSSatish Balay   sctx.shift_hi       = 0;
5050481f469SBarry Smith   sctx.shift_fraction = 0;
506ada7143aSSatish Balay 
507118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
5081a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
5099deb78dcSHong Zhang     PetscInt *aai = a->i;
51042898933SHong Zhang     ddiag          = a->diag;
511118f5789SBarry Smith     sctx.shift_top = 0;
5126cc28720Svictorle     for (i=0; i<n; i++) {
5139f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
5149f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
5151a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
516010a20caSHong Zhang       v  = a->a+aai[i];
517e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
518e24cfe64SHong Zhang       for (j=0; j<nz; j++)
5191a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
5201a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
5216cc28720Svictorle     }
522c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
523118f5789SBarry Smith     sctx.shift_top    *= 1.1;
524d2276718SHong Zhang     sctx.nshift_max   = 5;
525d2276718SHong Zhang     sctx.shift_lo     = 0.;
526d2276718SHong Zhang     sctx.shift_hi     = 1.;
527a0a356efSHong Zhang   }
528a0a356efSHong Zhang 
529a0a356efSHong Zhang   sctx.shift_amount = 0;
530a0a356efSHong Zhang   sctx.nshift       = 0;
531085256b3SBarry Smith   do {
532d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
533289bc588SBarry Smith     for (i=0; i<n; i++){
534b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
535b3bf805bSHong Zhang       bjtmp = bj + bi[i];
536b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
537289bc588SBarry Smith 
538289bc588SBarry Smith       /* load in initial (unfactored row) */
539416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
540b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
541010a20caSHong Zhang       v     = a->a + a->i[r[i]];
542085256b3SBarry Smith       for (j=0; j<nz; j++) {
543b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
544085256b3SBarry Smith       }
545d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
546289bc588SBarry Smith 
547b3bf805bSHong Zhang       row = *bjtmp++;
548f2582111SSatish Balay       while  (row < i) {
5498c37ef55SBarry Smith         pc = rtmp + row;
5508c37ef55SBarry Smith         if (*pc != 0.0) {
551010a20caSHong Zhang           pv         = b->a + diag_offset[row];
552010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
55335aab85fSBarry Smith           multiplier = *pc / *pv++;
5548c37ef55SBarry Smith           *pc        = multiplier;
555b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
556f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
557efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
558289bc588SBarry Smith         }
559b3bf805bSHong Zhang         row = *bjtmp++;
560289bc588SBarry Smith       }
561416022c9SBarry Smith       /* finished row so stick it into b->a */
562b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
563b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
564b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
565b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
566d3d32019SBarry Smith       rs   = 0.0;
567d3d32019SBarry Smith       for (j=0; j<nz; j++) {
568d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
569d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
570d3d32019SBarry Smith       }
571d2276718SHong Zhang 
572b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
573d2276718SHong Zhang       sctx.rs  = rs;
574d2276718SHong Zhang       sctx.pv  = pv[diag];
575c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
5765b5bc046SBarry Smith       if (newshift == 1) break;
57735aab85fSBarry Smith     }
578d2276718SHong Zhang 
5790481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5806cc28720Svictorle       /*
5816c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5826cc28720Svictorle        * then try lower shift
5836cc28720Svictorle        */
5840481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
5850481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
5860481f469SBarry Smith       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
587d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
588d2276718SHong Zhang       sctx.nshift++;
5896cc28720Svictorle     }
590d2276718SHong Zhang   } while (sctx.lushift);
5910f11f9aeSSatish Balay 
59235aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
59335aab85fSBarry Smith   for (i=0; i<n; i++) {
594010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
59535aab85fSBarry Smith   }
596606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
59778b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
59878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
599db4efbfdSBarry Smith   if (b->inode.use) {
600db4efbfdSBarry Smith     C->ops->solve   = MatSolve_Inode;
601db4efbfdSBarry Smith   } else {
6028d8c7c43SBarry Smith     PetscTruth row_identity, col_identity;
6038d8c7c43SBarry Smith     ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
6048d8c7c43SBarry Smith     ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
6058d8c7c43SBarry Smith     if (row_identity && col_identity) {
6068d8c7c43SBarry Smith       C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering;
6078d8c7c43SBarry Smith     } else {
608db4efbfdSBarry Smith       C->ops->solve   = MatSolve_SeqAIJ;
609db4efbfdSBarry Smith     }
6108d8c7c43SBarry Smith   }
611db4efbfdSBarry Smith   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
612db4efbfdSBarry Smith   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
613db4efbfdSBarry Smith   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
614de26497bSBarry Smith   C->ops->matsolve           = MatMatSolve_SeqAIJ;
615c456f294SBarry Smith   C->assembled    = PETSC_TRUE;
6165c9eb25fSBarry Smith   C->preallocated = PETSC_TRUE;
617d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
618d2276718SHong Zhang   if (sctx.nshift){
619e738e422SBarry Smith      if (info->shiftpd) {
6200481f469SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
621e738e422SBarry Smith     } else if (info->shiftnz) {
622e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
6236cc28720Svictorle     }
6240697b6caSHong Zhang   }
6253a40ed3dSBarry Smith   PetscFunctionReturn(0);
626289bc588SBarry Smith }
6270889b5dcSvictorle 
628137fb511SHong Zhang /*
629137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
630137fb511SHong Zhang    Input:
631137fb511SHong Zhang      A - original matrix
632137fb511SHong Zhang    Output;
633137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
634137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
635137fb511SHong Zhang          a->a reordered accordingly with a->j
636137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
637137fb511SHong Zhang */
638137fb511SHong Zhang #undef __FUNCT__
639137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
6400481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
641137fb511SHong Zhang {
642b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
643b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
644137fb511SHong Zhang   PetscErrorCode ierr;
6455d0c19d7SBarry Smith   const PetscInt *r,*ic,*ics;
6465d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
6475d0c19d7SBarry Smith   PetscInt       *ajtmp,nz,row;
648b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
649a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
650a77337e4SBarry Smith   MatScalar      *v,*pv;
651137fb511SHong Zhang   PetscReal      rs;
652137fb511SHong Zhang   LUShift_Ctx    sctx;
653b051339dSHong Zhang   PetscInt       newshift;
654137fb511SHong Zhang 
655137fb511SHong Zhang   PetscFunctionBegin;
656719d5645SBarry Smith   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
657137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
658137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
659137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
660137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
661b051339dSHong Zhang   ics = ic;
662137fb511SHong Zhang 
663137fb511SHong Zhang   sctx.shift_top      = 0;
664137fb511SHong Zhang   sctx.nshift_max     = 0;
665137fb511SHong Zhang   sctx.shift_lo       = 0;
666137fb511SHong Zhang   sctx.shift_hi       = 0;
6670481f469SBarry Smith   sctx.shift_fraction = 0;
668137fb511SHong Zhang 
669137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
670137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
671137fb511SHong Zhang     sctx.shift_top = 0;
672137fb511SHong Zhang     for (i=0; i<n; i++) {
673137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
674b051339dSHong Zhang       d  = (a->a)[diag[i]];
675137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
676b051339dSHong Zhang       v  = a->a+ai[i];
677b051339dSHong Zhang       nz = ai[i+1] - ai[i];
678137fb511SHong Zhang       for (j=0; j<nz; j++)
679137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
680137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
681137fb511SHong Zhang     }
682137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
683137fb511SHong Zhang     sctx.shift_top    *= 1.1;
684137fb511SHong Zhang     sctx.nshift_max   = 5;
685137fb511SHong Zhang     sctx.shift_lo     = 0.;
686137fb511SHong Zhang     sctx.shift_hi     = 1.;
687137fb511SHong Zhang   }
688137fb511SHong Zhang 
689137fb511SHong Zhang   sctx.shift_amount = 0;
690137fb511SHong Zhang   sctx.nshift       = 0;
691137fb511SHong Zhang   do {
692137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
693137fb511SHong Zhang     for (i=0; i<n; i++){
694b051339dSHong Zhang       /* load in initial unfactored row */
695b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
696b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
697b051339dSHong Zhang       v     = a->a + ai[r[i]];
698137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
699b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
700137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
701137fb511SHong Zhang 
702b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
703137fb511SHong Zhang       for (j=0; j<nz; j++) {
704137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
705b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
706137fb511SHong Zhang       }
707137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
708137fb511SHong Zhang 
709b051339dSHong Zhang       row = *ajtmp++;
710137fb511SHong Zhang       while  (row < i) {
711137fb511SHong Zhang         pc = rtmp + row;
712137fb511SHong Zhang         if (*pc != 0.0) {
713b051339dSHong Zhang           pv         = a->a + diag[r[row]];
714b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
715137fb511SHong Zhang 
716137fb511SHong Zhang           multiplier = *pc / *pv++;
717137fb511SHong Zhang           *pc        = multiplier;
718b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
719b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
720137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
721137fb511SHong Zhang         }
722b051339dSHong Zhang         row = *ajtmp++;
723137fb511SHong Zhang       }
724b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
725b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
726b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
727b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
728b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
729137fb511SHong Zhang 
730137fb511SHong Zhang       rs   = 0.0;
731137fb511SHong Zhang       for (j=0; j<nz; j++) {
732b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
733b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
734137fb511SHong Zhang       }
735137fb511SHong Zhang 
736137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
737137fb511SHong Zhang       sctx.rs  = rs;
738b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
739c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
740137fb511SHong Zhang       if (newshift == 1) break;
741137fb511SHong Zhang     }
742137fb511SHong Zhang 
7430481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
744137fb511SHong Zhang       /*
745137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
746137fb511SHong Zhang        * then try lower shift
747137fb511SHong Zhang        */
7480481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
7490481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
7500481f469SBarry Smith       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
751137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
752137fb511SHong Zhang       sctx.nshift++;
753137fb511SHong Zhang     }
754137fb511SHong Zhang   } while (sctx.lushift);
755137fb511SHong Zhang 
756137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
757137fb511SHong Zhang   for (i=0; i<n; i++) {
758b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
759137fb511SHong Zhang   }
760137fb511SHong Zhang 
761137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
762137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
763137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
764b051339dSHong Zhang   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
765db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
766db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
767db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
768b051339dSHong Zhang   A->assembled = PETSC_TRUE;
7695c9eb25fSBarry Smith   A->preallocated = PETSC_TRUE;
770d0f46423SBarry Smith   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
771137fb511SHong Zhang   if (sctx.nshift){
772e738e422SBarry Smith     if (info->shiftpd) {
7730481f469SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
774e738e422SBarry Smith     } else if (info->shiftnz) {
775e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
776137fb511SHong Zhang     }
777137fb511SHong Zhang   }
778137fb511SHong Zhang   PetscFunctionReturn(0);
779137fb511SHong Zhang }
780137fb511SHong Zhang 
7810a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
7840481f469SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
785da3a660dSBarry Smith {
786dfbe8321SBarry Smith   PetscErrorCode ierr;
787416022c9SBarry Smith   Mat            C;
788416022c9SBarry Smith 
7893a40ed3dSBarry Smith   PetscFunctionBegin;
79043244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
791719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
792719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
793db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
794db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
795273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
79652e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
7973a40ed3dSBarry Smith   PetscFunctionReturn(0);
798da3a660dSBarry Smith }
7990a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
8004a2ae208SSatish Balay #undef __FUNCT__
8014a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
802dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
8038c37ef55SBarry Smith {
804416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
805416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
8066849ba73SBarry Smith   PetscErrorCode    ierr;
8075d0c19d7SBarry Smith   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8085d0c19d7SBarry Smith   PetscInt          nz;
8095d0c19d7SBarry Smith   const PetscInt    *rout,*cout,*r,*c;
810dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
811d9fead3dSBarry Smith   const PetscScalar *b;
812dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
8138c37ef55SBarry Smith 
8143a40ed3dSBarry Smith   PetscFunctionBegin;
8153a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
81691bf9adeSBarry Smith 
817d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
819416022c9SBarry Smith   tmp  = a->solve_work;
8208c37ef55SBarry Smith 
8213b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8223b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
8238c37ef55SBarry Smith 
8248c37ef55SBarry Smith   /* forward solve the lower triangular */
8258c37ef55SBarry Smith   tmp[0] = b[*r++];
826010a20caSHong Zhang   tmps   = tmp;
8278c37ef55SBarry Smith   for (i=1; i<n; i++) {
828010a20caSHong Zhang     v   = aa + ai[i] ;
829010a20caSHong Zhang     vi  = aj + ai[i] ;
83053e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
83153e7425aSBarry Smith     sum = b[*r++];
832ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8338c37ef55SBarry Smith     tmp[i] = sum;
8348c37ef55SBarry Smith   }
8358c37ef55SBarry Smith 
8368c37ef55SBarry Smith   /* backward solve the upper triangular */
8378c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
838010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
839010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
840416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
8418c37ef55SBarry Smith     sum = tmp[i];
842ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
843010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
8448c37ef55SBarry Smith   }
8458c37ef55SBarry Smith 
8463b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8473b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
848d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8491ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
850d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
8513a40ed3dSBarry Smith   PetscFunctionReturn(0);
8528c37ef55SBarry Smith }
853026e39d0SSatish Balay 
8542bebee5dSHong Zhang #undef __FUNCT__
8552bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
8562bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8572bebee5dSHong Zhang {
8582bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8592bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8602bebee5dSHong Zhang   PetscErrorCode  ierr;
8615d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8625d0c19d7SBarry Smith   PetscInt        nz,neq;
8635d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
864dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
865dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
866db4efbfdSBarry Smith   PetscTruth      bisdense,xisdense;
8672bebee5dSHong Zhang 
8682bebee5dSHong Zhang   PetscFunctionBegin;
8692bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8702bebee5dSHong Zhang 
871db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
872db4efbfdSBarry Smith   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
873db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
874db4efbfdSBarry Smith   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
875db4efbfdSBarry Smith 
8762bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8772bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8782bebee5dSHong Zhang 
8792bebee5dSHong Zhang   tmp  = a->solve_work;
8802bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8812bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8822bebee5dSHong Zhang 
883a36b22ccSSatish Balay   for (neq=0; neq<B->cmap->n; neq++){
8842bebee5dSHong Zhang     /* forward solve the lower triangular */
8852bebee5dSHong Zhang     tmp[0] = b[r[0]];
8862bebee5dSHong Zhang     tmps   = tmp;
8872bebee5dSHong Zhang     for (i=1; i<n; i++) {
8882bebee5dSHong Zhang       v   = aa + ai[i] ;
8892bebee5dSHong Zhang       vi  = aj + ai[i] ;
8902bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8912bebee5dSHong Zhang       sum = b[r[i]];
8922bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8932bebee5dSHong Zhang       tmp[i] = sum;
8942bebee5dSHong Zhang     }
8952bebee5dSHong Zhang     /* backward solve the upper triangular */
8962bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
8972bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
8982bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
8992bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
9002bebee5dSHong Zhang       sum = tmp[i];
9012bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
9022bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
9032bebee5dSHong Zhang     }
9042bebee5dSHong Zhang 
9052bebee5dSHong Zhang     b += n;
9062bebee5dSHong Zhang     x += n;
9072bebee5dSHong Zhang   }
9082bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
9092bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
9102bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
9112bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
9122bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
9132bebee5dSHong Zhang   PetscFunctionReturn(0);
9142bebee5dSHong Zhang }
9152bebee5dSHong Zhang 
916137fb511SHong Zhang #undef __FUNCT__
917137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
918137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
919137fb511SHong Zhang {
920137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
921137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
922137fb511SHong Zhang   PetscErrorCode  ierr;
9235d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
9245d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
9255d0c19d7SBarry Smith   PetscInt        nz,row;
926dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
927dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
928137fb511SHong Zhang 
929137fb511SHong Zhang   PetscFunctionBegin;
930137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
931137fb511SHong Zhang 
932137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
933137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
934137fb511SHong Zhang   tmp  = a->solve_work;
935137fb511SHong Zhang 
936137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
937137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
938137fb511SHong Zhang 
939137fb511SHong Zhang   /* forward solve the lower triangular */
940137fb511SHong Zhang   tmp[0] = b[*r++];
941137fb511SHong Zhang   tmps   = tmp;
942137fb511SHong Zhang   for (row=1; row<n; row++) {
943137fb511SHong Zhang     i   = rout[row]; /* permuted row */
944137fb511SHong Zhang     v   = aa + ai[i] ;
945137fb511SHong Zhang     vi  = aj + ai[i] ;
946137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
947137fb511SHong Zhang     sum = b[*r++];
948137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
949137fb511SHong Zhang     tmp[row] = sum;
950137fb511SHong Zhang   }
951137fb511SHong Zhang 
952137fb511SHong Zhang   /* backward solve the upper triangular */
953137fb511SHong Zhang   for (row=n-1; row>=0; row--){
954137fb511SHong Zhang     i   = rout[row]; /* permuted row */
955137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
956137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
957137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
958137fb511SHong Zhang     sum = tmp[row];
959137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
960137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
961137fb511SHong Zhang   }
962137fb511SHong Zhang 
963137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
964137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
965137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
966137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
967d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
968137fb511SHong Zhang   PetscFunctionReturn(0);
969137fb511SHong Zhang }
970137fb511SHong Zhang 
971930ae53cSBarry Smith /* ----------------------------------------------------------- */
9724a2ae208SSatish Balay #undef __FUNCT__
9734a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
974dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
975930ae53cSBarry Smith {
976930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9776849ba73SBarry Smith   PetscErrorCode    ierr;
978d0f46423SBarry Smith   PetscInt          n = A->rmap->n;
979356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
980356650c2SBarry Smith   PetscScalar       *x;
981356650c2SBarry Smith   const PetscScalar *b;
982dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
983aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
984356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
985dd6ea824SBarry Smith   const MatScalar   *v;
986dd6ea824SBarry Smith   PetscScalar       sum;
987d85afc42SSatish Balay #endif
988930ae53cSBarry Smith 
9893a40ed3dSBarry Smith   PetscFunctionBegin;
9903a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
991930ae53cSBarry Smith 
992356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9931ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
994930ae53cSBarry Smith 
995aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
9961c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
9971c4feecaSSatish Balay #else
998930ae53cSBarry Smith   /* forward solve the lower triangular */
999930ae53cSBarry Smith   x[0] = b[0];
1000930ae53cSBarry Smith   for (i=1; i<n; i++) {
1001930ae53cSBarry Smith     ai_i = ai[i];
1002930ae53cSBarry Smith     v    = aa + ai_i;
1003930ae53cSBarry Smith     vi   = aj + ai_i;
1004930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
1005930ae53cSBarry Smith     sum  = b[i];
1006930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1007930ae53cSBarry Smith     x[i] = sum;
1008930ae53cSBarry Smith   }
1009930ae53cSBarry Smith 
1010930ae53cSBarry Smith   /* backward solve the upper triangular */
1011930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
1012930ae53cSBarry Smith     adiag_i = adiag[i];
1013d708749eSBarry Smith     v       = aa + adiag_i + 1;
1014d708749eSBarry Smith     vi      = aj + adiag_i + 1;
1015930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
1016930ae53cSBarry Smith     sum     = x[i];
1017930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1018930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
1019930ae53cSBarry Smith   }
10201c4feecaSSatish Balay #endif
1021d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
1022356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
10231ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1025930ae53cSBarry Smith }
1026930ae53cSBarry Smith 
10274a2ae208SSatish Balay #undef __FUNCT__
10284a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
1029dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1030da3a660dSBarry Smith {
1031416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1032416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
10336849ba73SBarry Smith   PetscErrorCode  ierr;
10345d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10355d0c19d7SBarry Smith   PetscInt        nz;
10365d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
1037dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
1038dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1039da3a660dSBarry Smith 
10403a40ed3dSBarry Smith   PetscFunctionBegin;
104178b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
1042da3a660dSBarry Smith 
10431ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10441ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1045416022c9SBarry Smith   tmp  = a->solve_work;
1046da3a660dSBarry Smith 
10473b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10483b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1049da3a660dSBarry Smith 
1050da3a660dSBarry Smith   /* forward solve the lower triangular */
1051da3a660dSBarry Smith   tmp[0] = b[*r++];
1052da3a660dSBarry Smith   for (i=1; i<n; i++) {
1053010a20caSHong Zhang     v   = aa + ai[i] ;
1054010a20caSHong Zhang     vi  = aj + ai[i] ;
1055416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
1056da3a660dSBarry Smith     sum = b[*r++];
1057010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1058da3a660dSBarry Smith     tmp[i] = sum;
1059da3a660dSBarry Smith   }
1060da3a660dSBarry Smith 
1061da3a660dSBarry Smith   /* backward solve the upper triangular */
1062da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1063010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
1064010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
1065416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1066da3a660dSBarry Smith     sum = tmp[i];
1067010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1068010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1069da3a660dSBarry Smith     x[*c--] += tmp[i];
1070da3a660dSBarry Smith   }
1071da3a660dSBarry Smith 
10723b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10733b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10741ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10751ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1076efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1077898183ecSLois Curfman McInnes 
10783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1079da3a660dSBarry Smith }
1080da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10814a2ae208SSatish Balay #undef __FUNCT__
10824a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1083dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1084da3a660dSBarry Smith {
1085416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10862235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10876849ba73SBarry Smith   PetscErrorCode  ierr;
10885d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
10895d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10905d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1091dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1092dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1093da3a660dSBarry Smith 
10943a40ed3dSBarry Smith   PetscFunctionBegin;
10951ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10961ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1097416022c9SBarry Smith   tmp  = a->solve_work;
1098da3a660dSBarry Smith 
10992235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11002235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1101da3a660dSBarry Smith 
1102da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
11032235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1104da3a660dSBarry Smith 
1105da3a660dSBarry Smith   /* forward solve the U^T */
1106da3a660dSBarry Smith   for (i=0; i<n; i++) {
1107010a20caSHong Zhang     v   = aa + diag[i] ;
1108010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1109f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1110c38d4ed2SBarry Smith     s1  = tmp[i];
1111ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1112da3a660dSBarry Smith     while (nz--) {
1113010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1114da3a660dSBarry Smith     }
1115c38d4ed2SBarry Smith     tmp[i] = s1;
1116da3a660dSBarry Smith   }
1117da3a660dSBarry Smith 
1118da3a660dSBarry Smith   /* backward solve the L^T */
1119da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1120010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1121010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1122f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1123f1af5d2fSBarry Smith     s1  = tmp[i];
1124da3a660dSBarry Smith     while (nz--) {
1125010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1126da3a660dSBarry Smith     }
1127da3a660dSBarry Smith   }
1128da3a660dSBarry Smith 
1129da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1130da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1131da3a660dSBarry Smith 
11322235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11332235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11341ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1136da3a660dSBarry Smith 
1137d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap->n);CHKERRQ(ierr);
11383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1139da3a660dSBarry Smith }
1140da3a660dSBarry Smith 
11414a2ae208SSatish Balay #undef __FUNCT__
11424a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1143dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1144da3a660dSBarry Smith {
1145416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
11462235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
11476849ba73SBarry Smith   PetscErrorCode  ierr;
11485d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
11495d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
11505d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1151dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1152dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
11536abc6512SBarry Smith 
11543a40ed3dSBarry Smith   PetscFunctionBegin;
115523bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
11566abc6512SBarry Smith 
11571ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1159416022c9SBarry Smith   tmp = a->solve_work;
11606abc6512SBarry Smith 
11612235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11622235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11636abc6512SBarry Smith 
11646abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
11652235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
11666abc6512SBarry Smith 
11676abc6512SBarry Smith   /* forward solve the U^T */
11686abc6512SBarry Smith   for (i=0; i<n; i++) {
1169010a20caSHong Zhang     v   = aa + diag[i] ;
1170010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1171f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11726abc6512SBarry Smith     tmp[i] *= *v++;
11736abc6512SBarry Smith     while (nz--) {
1174010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11756abc6512SBarry Smith     }
11766abc6512SBarry Smith   }
11776abc6512SBarry Smith 
11786abc6512SBarry Smith   /* backward solve the L^T */
11796abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1180010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1181010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1182f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11836abc6512SBarry Smith     while (nz--) {
1184010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11856abc6512SBarry Smith     }
11866abc6512SBarry Smith   }
11876abc6512SBarry Smith 
11886abc6512SBarry Smith   /* copy tmp into x according to permutation */
11896abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11906abc6512SBarry Smith 
11912235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11922235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11931ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11956abc6512SBarry Smith 
1196efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
11973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1198da3a660dSBarry Smith }
11999e25ed09SBarry Smith /* ----------------------------------------------------------------*/
12003c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1201719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption);
120208480c60SBarry Smith 
12034a2ae208SSatish Balay #undef __FUNCT__
12044a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
12050481f469SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
12069e25ed09SBarry Smith {
1207416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
12089e25ed09SBarry Smith   IS                 isicol;
12096849ba73SBarry Smith   PetscErrorCode     ierr;
12105d0c19d7SBarry Smith   const PetscInt     *r,*ic;
12115d0c19d7SBarry Smith   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
12125a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
12135a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
12145a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
121577c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1216329f5518SBarry Smith   PetscReal          f;
12175a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12185a8e39fbSHong Zhang   PetscBT            lnkbt;
12195a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1220a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1221a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
122209f38230SBarry Smith   PetscTruth         missing;
12239e25ed09SBarry Smith 
12243a40ed3dSBarry Smith   PetscFunctionBegin;
1225d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
12266cf997b0SBarry Smith   f             = info->fill;
122738baddfdSBarry Smith   levels        = (PetscInt)info->levels;
122838baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
12294c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
123082bf6240SBarry Smith 
123108480c60SBarry Smith   /* special case that simply copies fill pattern */
1232be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1233be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
123486bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
1235719d5645SBarry Smith     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr);
1236719d5645SBarry Smith     fact->factor = MAT_FACTOR_ILU;
1237719d5645SBarry Smith     (fact)->info.factor_mallocs    = 0;
1238719d5645SBarry Smith     (fact)->info.fill_ratio_given  = info->fill;
1239719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 1.0;
1240719d5645SBarry Smith     b               = (Mat_SeqAIJ*)(fact)->data;
12418ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
124209f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
124308480c60SBarry Smith     b->row              = isrow;
124408480c60SBarry Smith     b->col              = iscol;
124582bf6240SBarry Smith     b->icol             = isicol;
1246719d5645SBarry Smith     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1247c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1248c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1249719d5645SBarry Smith     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1250719d5645SBarry Smith     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
12513a40ed3dSBarry Smith     PetscFunctionReturn(0);
125208480c60SBarry Smith   }
125308480c60SBarry Smith 
1254898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1255898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
12569e25ed09SBarry Smith 
12579e25ed09SBarry Smith   /* get new row pointers */
12585a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12595a8e39fbSHong Zhang   bi[0] = 0;
12605a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
12615a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12625a8e39fbSHong Zhang   bdiag[0]  = 0;
12636cf997b0SBarry Smith 
12645a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
12655a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
12665a8e39fbSHong Zhang 
12675a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
12685a8e39fbSHong Zhang   nlnk = n + 1;
12695a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12705a8e39fbSHong Zhang 
12715a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1272a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12735a8e39fbSHong Zhang   current_space = free_space;
1274a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12755a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12765a8e39fbSHong Zhang 
12775a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12785a8e39fbSHong Zhang     nzi = 0;
12796cf997b0SBarry Smith     /* copy current row into linked list */
12805a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
12815a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
12825a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12835a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12845a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12855a8e39fbSHong Zhang     nzi += nlnk;
12866cf997b0SBarry Smith 
12876cf997b0SBarry Smith     /* make sure diagonal entry is included */
12885a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12896cf997b0SBarry Smith       fm = n;
12905a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12915a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12925a8e39fbSHong Zhang       lnk[fm]    = i;
12935a8e39fbSHong Zhang       lnk_lvl[i] = 0;
12945a8e39fbSHong Zhang       nzi++; dcount++;
12956cf997b0SBarry Smith     }
12966cf997b0SBarry Smith 
12975a8e39fbSHong Zhang     /* add pivot rows into the active row */
12985a8e39fbSHong Zhang     nzbd = 0;
12995a8e39fbSHong Zhang     prow = lnk[n];
13005a8e39fbSHong Zhang     while (prow < i) {
13015a8e39fbSHong Zhang       nnz      = bdiag[prow];
13025a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
13035a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
13045a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
13055a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
13065a8e39fbSHong Zhang       nzi += nlnk;
13075a8e39fbSHong Zhang       prow = lnk[prow];
13085a8e39fbSHong Zhang       nzbd++;
130956beaf04SBarry Smith     }
13105a8e39fbSHong Zhang     bdiag[i] = nzbd;
13115a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1312ecf371e4SBarry Smith 
13135a8e39fbSHong Zhang     /* if free space is not available, make more free space */
13145a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
13155a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1316a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1317a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
13185a8e39fbSHong Zhang       reallocs++;
13195a8e39fbSHong Zhang     }
1320ecf371e4SBarry Smith 
13215a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
13225a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13235a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
13245a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
13255a8e39fbSHong Zhang 
13265a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
13275a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
132877431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1329d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
13306cf997b0SBarry Smith     }
13315a8e39fbSHong Zhang 
13325a8e39fbSHong Zhang     current_space->array           += nzi;
13335a8e39fbSHong Zhang     current_space->local_used      += nzi;
13345a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
13355a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
13365a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
13375a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
13389e25ed09SBarry Smith   }
13395a8e39fbSHong Zhang 
1340898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1341898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
13425a8e39fbSHong Zhang 
13435a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
13445a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1345a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13465a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1347a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13485a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
13499e25ed09SBarry Smith 
13506cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1351f64065bbSBarry Smith   {
13525a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1353ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1354ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1355ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1356ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1357335d9088SBarry Smith     if (diagonal_fill) {
1358ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1359335d9088SBarry Smith     }
1360f64065bbSBarry Smith   }
136163ba0a88SBarry Smith #endif
1362bbb6d6a8SBarry Smith 
13639e25ed09SBarry Smith   /* put together the new matrix */
1364719d5645SBarry Smith   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1365719d5645SBarry Smith   b = (Mat_SeqAIJ*)(fact)->data;
1366e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1367e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13687c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13695a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1370b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13715a8e39fbSHong Zhang   b->j          = bj;
13725a8e39fbSHong Zhang   b->i          = bi;
13735a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13745a8e39fbSHong Zhang   b->diag       = bdiag;
1375416022c9SBarry Smith   b->ilen       = 0;
1376416022c9SBarry Smith   b->imax       = 0;
1377416022c9SBarry Smith   b->row        = isrow;
1378416022c9SBarry Smith   b->col        = iscol;
1379c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1380c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
138182bf6240SBarry Smith   b->icol       = isicol;
138287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1383416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13845a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
1385719d5645SBarry Smith   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13865a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
1387719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1388719d5645SBarry Smith   (fact)->info.fill_ratio_given  = f;
1389719d5645SBarry Smith   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1390719d5645SBarry Smith   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1391719d5645SBarry Smith   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
13923a40ed3dSBarry Smith   PetscFunctionReturn(0);
13939e25ed09SBarry Smith }
1394d5d45c9bSBarry Smith 
13953c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1396a6175056SHong Zhang #undef __FUNCT__
1397a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
13980481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1399a6175056SHong Zhang {
1400719d5645SBarry Smith   Mat            C = B;
14010968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
14022ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
14039bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
14042ed133dbSHong Zhang   PetscErrorCode ierr;
14055d0c19d7SBarry Smith   const PetscInt *rip,*riip;
14065d0c19d7SBarry Smith   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
14072ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
14082ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1409622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1410540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1411fbf22428SSatish Balay   PetscReal      shiftpd;
1412540703acSHong Zhang   ChShift_Ctx    sctx;
1413540703acSHong Zhang   PetscInt       newshift;
1414db4efbfdSBarry Smith   PetscTruth     perm_identity;
1415d5d45c9bSBarry Smith 
1416a6175056SHong Zhang   PetscFunctionBegin;
1417117f1891Stmunson 
1418540703acSHong Zhang   shiftnz   = info->shiftnz;
1419540703acSHong Zhang   shiftpd   = info->shiftpd;
1420ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1421ee45ca4aSHong Zhang 
14222ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
14239bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
14242ed133dbSHong Zhang 
14252ed133dbSHong Zhang   /* initialization */
14262ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
14272ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
14282ed133dbSHong Zhang   jl   = il + mbs;
14292ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
14302ed133dbSHong Zhang 
1431540703acSHong Zhang   sctx.shift_amount = 0;
1432540703acSHong Zhang   sctx.nshift       = 0;
14332ed133dbSHong Zhang   do {
1434540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
14352ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14362ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1437a6175056SHong Zhang     }
14382ed133dbSHong Zhang 
14392ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1440c05c3958SHong Zhang       bval = ba + bi[k];
14412ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
14422ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
14432ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
14449bfd6278SHong Zhang         col = riip[aj[j]];
14452ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
14462ed133dbSHong Zhang           rtmp[col] = aa[j];
1447c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
14482ed133dbSHong Zhang         }
14492ed133dbSHong Zhang       }
145039e3d630SHong Zhang       /* shift the diagonal of the matrix */
1451540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
14522ed133dbSHong Zhang 
14532ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
14542ed133dbSHong Zhang       dk = rtmp[k];
14552ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
14562ed133dbSHong Zhang 
14572ed133dbSHong Zhang       while (i < k){
14582ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
14592ed133dbSHong Zhang 
14602ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
14612ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
14622ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
14632ed133dbSHong Zhang         dk += uikdi*ba[ili];
14642ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14652ed133dbSHong Zhang 
14662ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14672ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14682ed133dbSHong Zhang         if (jmin < jmax){
14692ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14702ed133dbSHong Zhang           /* update il and jl for row i */
14712ed133dbSHong Zhang           il[i] = jmin;
14722ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14732ed133dbSHong Zhang         }
14742ed133dbSHong Zhang         i = nexti;
14752ed133dbSHong Zhang       }
14762ed133dbSHong Zhang 
14770697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14780697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14790697b6caSHong Zhang       rs   = 0.0;
14803c9e8598SHong Zhang       jmin = bi[k]+1;
14813c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14823c9e8598SHong Zhang       bcol = bj + jmin;
14833c9e8598SHong Zhang       while (nz--){
148489f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
148589f845c8SHong Zhang         bcol++;
14863c9e8598SHong Zhang       }
1487540703acSHong Zhang 
1488540703acSHong Zhang       sctx.rs = rs;
1489540703acSHong Zhang       sctx.pv = dk;
14905b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1491117f1891Stmunson 
1492117f1891Stmunson       if (newshift == 1) {
1493117f1891Stmunson         if (!sctx.shift_amount) {
1494117f1891Stmunson           sctx.shift_amount = 1e-5;
1495117f1891Stmunson         }
1496117f1891Stmunson         break;
1497117f1891Stmunson       }
14983c9e8598SHong Zhang 
14993c9e8598SHong Zhang       /* copy data into U(k,:) */
150039e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
150139e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
150239e3d630SHong Zhang       if (jmin < jmax) {
150339e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
150439e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
15053c9e8598SHong Zhang         }
150639e3d630SHong Zhang         /* add the k-th row into il and jl */
15073c9e8598SHong Zhang         il[k] = jmin;
15083c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
15093c9e8598SHong Zhang       }
15103c9e8598SHong Zhang     }
1511540703acSHong Zhang   } while (sctx.chshift);
15123c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
15133c9e8598SHong Zhang 
151439e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
15159bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1516db4efbfdSBarry Smith 
1517db4efbfdSBarry Smith   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1518db4efbfdSBarry Smith   if (perm_identity){
1519719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1520719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1521719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1522719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1523db4efbfdSBarry Smith   } else {
1524719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1525719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1526719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1527719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1528db4efbfdSBarry Smith   }
1529db4efbfdSBarry Smith 
15303c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
15313c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1532d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1533540703acSHong Zhang   if (sctx.nshift){
1534540703acSHong Zhang     if (shiftnz) {
15351e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1536540703acSHong Zhang     } else if (shiftpd) {
15371e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
15380697b6caSHong Zhang     }
15393c9e8598SHong Zhang   }
15403c9e8598SHong Zhang   PetscFunctionReturn(0);
15413c9e8598SHong Zhang }
1542a6175056SHong Zhang 
1543a6175056SHong Zhang #undef __FUNCT__
1544a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
15450481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1546a6175056SHong Zhang {
15470968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1548ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1549dfbe8321SBarry Smith   PetscErrorCode     ierr;
155058ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
15515d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
15525d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
1553ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
155458ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
15555a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1556ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1557a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1558a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
15597a48dd6fSHong Zhang   PetscBT            lnkbt;
1560b635d36bSHong Zhang   IS                 iperm;
1561a6175056SHong Zhang 
1562a6175056SHong Zhang   PetscFunctionBegin;
1563d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
156458ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
156558ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1566653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1567b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
15687a48dd6fSHong Zhang 
156939e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
157039e3d630SHong Zhang   ui[0] = 0;
157139e3d630SHong Zhang 
1572b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1573622e2028SHong Zhang   if (!levels && perm_identity) {
157458ebbce7SBarry Smith 
1575ed59401aSHong Zhang     for (i=0; i<am; i++) {
157639e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1577ed59401aSHong Zhang     }
157839e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
157939e3d630SHong Zhang     cols = uj;
1580ed59401aSHong Zhang     for (i=0; i<am; i++) {
1581ed59401aSHong Zhang       aj    = a->j + a->diag[i];
158239e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
158339e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1584ed59401aSHong Zhang     }
158539e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1586b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1587b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1588b635d36bSHong Zhang 
15897a48dd6fSHong Zhang     /* initialization */
15905a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15917a48dd6fSHong Zhang 
15927a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
15937a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
15941393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
15957a48dd6fSHong Zhang     il         = jl + am;
15967a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
15977a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
15987a48dd6fSHong Zhang     for (i=0; i<am; i++){
15997a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
16007a48dd6fSHong Zhang     }
16017a48dd6fSHong Zhang 
16027a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
16037a48dd6fSHong Zhang     nlnk = am + 1;
16042ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16057a48dd6fSHong Zhang 
16067a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1607a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
16087a48dd6fSHong Zhang     current_space = free_space;
1609a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
16107a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
16117a48dd6fSHong Zhang 
16127a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
16137a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
16147a48dd6fSHong Zhang       nzk   = 0;
1615622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1616b635d36bSHong Zhang       if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
1617622e2028SHong Zhang       ncols_upper = 0;
1618622e2028SHong Zhang       for (j=0; j<ncols; j++){
1619b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1620b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
16215a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1622622e2028SHong Zhang           ncols_upper++;
1623622e2028SHong Zhang         }
1624622e2028SHong Zhang       }
1625b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16267a48dd6fSHong Zhang       nzk += nlnk;
16277a48dd6fSHong Zhang 
16287a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
16297a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
16307a48dd6fSHong Zhang 
16317a48dd6fSHong Zhang       while (prow < k){
16327a48dd6fSHong Zhang         nextprow = jl[prow];
16337a48dd6fSHong Zhang 
16347a48dd6fSHong Zhang         /* merge prow into k-th row */
16357a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
16367a48dd6fSHong Zhang         jmax = ui[prow+1];
16377a48dd6fSHong Zhang         ncols = jmax-jmin;
1638ed59401aSHong Zhang         i     = jmin - ui[prow];
1639ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
16405a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
16415a8e39fbSHong Zhang         j     = *(uj - 1);
16425a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
16437a48dd6fSHong Zhang         nzk += nlnk;
16447a48dd6fSHong Zhang 
16457a48dd6fSHong Zhang         /* update il and jl for prow */
16467a48dd6fSHong Zhang         if (jmin < jmax){
16477a48dd6fSHong Zhang           il[prow] = jmin;
16487a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
16497a48dd6fSHong Zhang         }
16507a48dd6fSHong Zhang         prow = nextprow;
16517a48dd6fSHong Zhang       }
16527a48dd6fSHong Zhang 
16537a48dd6fSHong Zhang       /* if free space is not available, make more free space */
16547a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
16557a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
16567a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1657a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1658a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
16597a48dd6fSHong Zhang         reallocs++;
16607a48dd6fSHong Zhang       }
16617a48dd6fSHong Zhang 
16627a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1663b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
16642ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
16657a48dd6fSHong Zhang 
16667a48dd6fSHong Zhang       /* add the k-th row into il and jl */
166739e3d630SHong Zhang       if (nzk > 1){
16687a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
16697a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
16707a48dd6fSHong Zhang         il[k] = ui[k] + 1;
16717a48dd6fSHong Zhang       }
16727a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
16737a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
16747a48dd6fSHong Zhang 
16757a48dd6fSHong Zhang       current_space->array           += nzk;
16767a48dd6fSHong Zhang       current_space->local_used      += nzk;
16777a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16787a48dd6fSHong Zhang 
16797a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16807a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16817a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16827a48dd6fSHong Zhang 
16837a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16847a48dd6fSHong Zhang     }
16857a48dd6fSHong Zhang 
16866cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16877a48dd6fSHong Zhang     if (ai[am] != 0) {
168839e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1689ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1690ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1691ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16927a48dd6fSHong Zhang     } else {
1693ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16947a48dd6fSHong Zhang     }
169563ba0a88SBarry Smith #endif
16967a48dd6fSHong Zhang 
16977a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1698b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16997a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
17005a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
17017a48dd6fSHong Zhang 
17027a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
17037a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1704a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
17052ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1706a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
17077a48dd6fSHong Zhang 
170839e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
170939e3d630SHong Zhang 
17107a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
17117a48dd6fSHong Zhang 
1712719d5645SBarry Smith   b    = (Mat_SeqSBAIJ*)(fact)->data;
17137a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
17147a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
17157a48dd6fSHong Zhang   b->j    = uj;
17167a48dd6fSHong Zhang   b->i    = ui;
17177a48dd6fSHong Zhang   b->diag = 0;
17187a48dd6fSHong Zhang   b->ilen = 0;
17197a48dd6fSHong Zhang   b->imax = 0;
17207a48dd6fSHong Zhang   b->row  = perm;
1721b635d36bSHong Zhang   b->col  = perm;
1722b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1723b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1724b635d36bSHong Zhang   b->icol = iperm;
17257a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
17267a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1727719d5645SBarry Smith   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
17287a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1729e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1730e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
17317a48dd6fSHong Zhang 
1732719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1733719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
17347a48dd6fSHong Zhang   if (ai[am] != 0) {
1735719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
17367a48dd6fSHong Zhang   } else {
1737719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
17387a48dd6fSHong Zhang   }
1739719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1740a6175056SHong Zhang   PetscFunctionReturn(0);
1741a6175056SHong Zhang }
1742d5d45c9bSBarry Smith 
1743f76d2b81SHong Zhang #undef __FUNCT__
1744f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
17450481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1746f76d2b81SHong Zhang {
1747f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
174810c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
1749dfbe8321SBarry Smith   PetscErrorCode     ierr;
1750f76d2b81SHong Zhang   PetscTruth         perm_identity;
175110c27e3dSHong Zhang   PetscReal          fill = info->fill;
17525d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
17535d0c19d7SBarry Smith   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
175410c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
17552ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1756a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
175710c27e3dSHong Zhang   PetscBT            lnkbt;
17582ed133dbSHong Zhang   IS                 iperm;
1759f76d2b81SHong Zhang 
1760f76d2b81SHong Zhang   PetscFunctionBegin;
1761d0f46423SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
17622ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1763f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
17642ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
17652ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
17669bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
176710c27e3dSHong Zhang 
176810c27e3dSHong Zhang   /* initialization */
17691393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
177010c27e3dSHong Zhang   ui[0] = 0;
177110c27e3dSHong Zhang 
177210c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17731393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17741393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17751393bc97SHong Zhang   il     = jl + am;
17761393bc97SHong Zhang   cols   = il + am;
17771393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17781393bc97SHong Zhang   for (i=0; i<am; i++){
17791393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1780f76d2b81SHong Zhang   }
178110c27e3dSHong Zhang 
178210c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17831393bc97SHong Zhang   nlnk = am + 1;
17841393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
178510c27e3dSHong Zhang 
17861393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1787a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
178810c27e3dSHong Zhang   current_space = free_space;
178910c27e3dSHong Zhang 
17901393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
179110c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
179210c27e3dSHong Zhang     nzk   = 0;
17932ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
17949bfd6278SHong Zhang     if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
17952ed133dbSHong Zhang     ncols_upper = 0;
17962ed133dbSHong Zhang     for (j=0; j<ncols; j++){
17979bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
17982ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
17992ed133dbSHong Zhang         cols[ncols_upper] = i;
18002ed133dbSHong Zhang         ncols_upper++;
18012ed133dbSHong Zhang       }
18022ed133dbSHong Zhang     }
18031393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
180410c27e3dSHong Zhang     nzk += nlnk;
180510c27e3dSHong Zhang 
180610c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
180710c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
180810c27e3dSHong Zhang 
180910c27e3dSHong Zhang     while (prow < k){
181010c27e3dSHong Zhang       nextprow = jl[prow];
181110c27e3dSHong Zhang       /* merge prow into k-th row */
18121393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
181310c27e3dSHong Zhang       jmax = ui[prow+1];
181410c27e3dSHong Zhang       ncols = jmax-jmin;
18151393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
18165a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
181710c27e3dSHong Zhang       nzk += nlnk;
181810c27e3dSHong Zhang 
181910c27e3dSHong Zhang       /* update il and jl for prow */
182010c27e3dSHong Zhang       if (jmin < jmax){
182110c27e3dSHong Zhang         il[prow] = jmin;
18222ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
182310c27e3dSHong Zhang       }
182410c27e3dSHong Zhang       prow = nextprow;
182510c27e3dSHong Zhang     }
182610c27e3dSHong Zhang 
182710c27e3dSHong Zhang     /* if free space is not available, make more free space */
182810c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
18291393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
183010c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1831a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
183210c27e3dSHong Zhang       reallocs++;
183310c27e3dSHong Zhang     }
183410c27e3dSHong Zhang 
183510c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
18361393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
183710c27e3dSHong Zhang 
183810c27e3dSHong Zhang     /* add the k-th row into il and jl */
183910c27e3dSHong Zhang     if (nzk-1 > 0){
18401393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
184110c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
184210c27e3dSHong Zhang       il[k] = ui[k] + 1;
184310c27e3dSHong Zhang     }
18442ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
184510c27e3dSHong Zhang     current_space->array           += nzk;
184610c27e3dSHong Zhang     current_space->local_used      += nzk;
184710c27e3dSHong Zhang     current_space->local_remaining -= nzk;
184810c27e3dSHong Zhang 
184910c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
185010c27e3dSHong Zhang   }
185110c27e3dSHong Zhang 
18526cf91177SBarry Smith #if defined(PETSC_USE_INFO)
18531393bc97SHong Zhang   if (ai[am] != 0) {
18541393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1855ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1856ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1857ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
185810c27e3dSHong Zhang   } else {
1859ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
186010c27e3dSHong Zhang   }
186163ba0a88SBarry Smith #endif
186210c27e3dSHong Zhang 
186310c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
18649bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
186510c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
186610c27e3dSHong Zhang 
186710c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18681393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1869a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
187010c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
187110c27e3dSHong Zhang 
187210c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
187310c27e3dSHong Zhang 
1874719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
187510c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1876e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1877e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18781393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
187910c27e3dSHong Zhang   b->j    = uj;
188010c27e3dSHong Zhang   b->i    = ui;
188110c27e3dSHong Zhang   b->diag = 0;
188210c27e3dSHong Zhang   b->ilen = 0;
188310c27e3dSHong Zhang   b->imax = 0;
188410c27e3dSHong Zhang   b->row  = perm;
18859bfd6278SHong Zhang   b->col  = perm;
18869bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18879bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18889bfd6278SHong Zhang   b->icol = iperm;
188910c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18901393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1891719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18921393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
189310c27e3dSHong Zhang 
1894719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1895719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
18961393bc97SHong Zhang   if (ai[am] != 0) {
1897719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
189810c27e3dSHong Zhang   } else {
1899719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
190010c27e3dSHong Zhang   }
1901719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1902f76d2b81SHong Zhang   PetscFunctionReturn(0);
1903f76d2b81SHong Zhang }
1904