xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 3b4a8b6d5dd2aa52bbea25c73643b7fba0754cce)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
37c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
47c4f633dSBarry Smith #include "../src/inline/dot.h"
57c4f633dSBarry Smith #include "../src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
77c4f633dSBarry Smith #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;
7623395705SBarry Smith   if (info->fill == PETSC_DEFAULT)    fill    = ((double)(n*(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 
25123395705SBarry Smith   if (reorder) (*fact)->ops->solve = MatSolve_SeqAIJ;
25223395705SBarry Smith   else         (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
25323395705SBarry Smith   (*fact)->ops->solveadd           = MatSolveAdd_SeqAIJ;
25423395705SBarry Smith   (*fact)->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
25523395705SBarry Smith   (*fact)->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
25623395705SBarry Smith 
2577529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25871c2f376SKris Buschelman 
25986bacbd2SBarry Smith   PetscFunctionReturn(0);
26060ec86bdSBarry Smith #endif
26186bacbd2SBarry Smith }
26286bacbd2SBarry Smith 
263e631078cSBarry Smith EXTERN_C_BEGIN
2644a2ae208SSatish Balay #undef __FUNCT__
265db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
266db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
267db4efbfdSBarry Smith {
268db4efbfdSBarry Smith   PetscFunctionBegin;
269db4efbfdSBarry Smith   *flg = PETSC_TRUE;
270db4efbfdSBarry Smith   PetscFunctionReturn(0);
271db4efbfdSBarry Smith }
272db4efbfdSBarry Smith EXTERN_C_END
273db4efbfdSBarry Smith 
274db4efbfdSBarry Smith EXTERN_C_BEGIN
275db4efbfdSBarry Smith #undef __FUNCT__
276b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc"
277b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
278b24902e0SBarry Smith {
279d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
280b24902e0SBarry Smith   PetscErrorCode     ierr;
281b24902e0SBarry Smith 
282b24902e0SBarry Smith   PetscFunctionBegin;
283b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
284b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
285b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
286b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
28735bd34faSBarry Smith     if (ftype == MAT_FACTOR_ILU) {
288db4efbfdSBarry Smith       (*B)->ops->ilufactorsymbolic= MatILUFactorSymbolic_SeqAIJ;
28935bd34faSBarry Smith     } else {
29035bd34faSBarry Smith       (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
29135bd34faSBarry Smith     }
292b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
293b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
294b24902e0SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
29535bd34faSBarry Smith     if (ftype == MAT_FACTOR_ICC) {
2965c9eb25fSBarry Smith       (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
29735bd34faSBarry Smith     } else {
2985c9eb25fSBarry Smith       (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
29935bd34faSBarry Smith     }
300b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
301db4efbfdSBarry Smith   (*B)->factor = ftype;
302b24902e0SBarry Smith   PetscFunctionReturn(0);
303b24902e0SBarry Smith }
304e631078cSBarry Smith EXTERN_C_END
305b24902e0SBarry Smith 
306b24902e0SBarry Smith #undef __FUNCT__
307b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
3080481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
309289bc588SBarry Smith {
310416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
311289bc588SBarry Smith   IS                 isicol;
3126849ba73SBarry Smith   PetscErrorCode     ierr;
3135d0c19d7SBarry Smith   const PetscInt     *r,*ic;
3145d0c19d7SBarry Smith   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
3151393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
3161393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
3179e046878SBarry Smith   PetscReal          f;
3184875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
319a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3201393bc97SHong Zhang   PetscBT            lnkbt;
321289bc588SBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
323d0f46423SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
3244c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
325ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
326ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
327289bc588SBarry Smith 
328289bc588SBarry Smith   /* get new row pointers */
3291393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
3301393bc97SHong Zhang   bi[0] = 0;
3311393bc97SHong Zhang 
3321393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
3331393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
3341393bc97SHong Zhang   bdiag[0] = 0;
3351393bc97SHong Zhang 
3361393bc97SHong Zhang   /* linked list for storing column indices of the active row */
3371393bc97SHong Zhang   nlnk = n + 1;
3381393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3391393bc97SHong Zhang 
34035baf27eSSatish Balay   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
3411393bc97SHong Zhang 
3421393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
343d3d32019SBarry Smith   f = info->fill;
344a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3451393bc97SHong Zhang   current_space = free_space;
346289bc588SBarry Smith 
347289bc588SBarry Smith   for (i=0; i<n; i++) {
3481393bc97SHong Zhang     /* copy previous fill into linked list */
349289bc588SBarry Smith     nzi = 0;
3501393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
351*3b4a8b6dSBarry Smith     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
3521393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
353afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3541393bc97SHong Zhang     nzi += nlnk;
3551393bc97SHong Zhang 
3561393bc97SHong Zhang     /* add pivot rows into linked list */
3571393bc97SHong Zhang     row = lnk[n];
3581393bc97SHong Zhang     while (row < i) {
359d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
360d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
361d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3621393bc97SHong Zhang       nzi += nlnk;
3631393bc97SHong Zhang       row  = lnk[row];
364289bc588SBarry Smith     }
3651393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3661393bc97SHong Zhang     im[i]   = nzi;
3671393bc97SHong Zhang 
3681393bc97SHong Zhang     /* mark bdiag */
3691393bc97SHong Zhang     nzbd = 0;
3701393bc97SHong Zhang     nnz  = nzi;
3711393bc97SHong Zhang     k    = lnk[n];
3721393bc97SHong Zhang     while (nnz-- && k < i){
3731393bc97SHong Zhang       nzbd++;
3741393bc97SHong Zhang       k = lnk[k];
3751393bc97SHong Zhang     }
3761393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3771393bc97SHong Zhang 
3781393bc97SHong Zhang     /* if free space is not available, make more free space */
3791393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3801393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
381a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3821393bc97SHong Zhang       reallocs++;
3831393bc97SHong Zhang     }
3841393bc97SHong Zhang 
3851393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3861393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3871393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3881393bc97SHong Zhang     current_space->array           += nzi;
3891393bc97SHong Zhang     current_space->local_used      += nzi;
3901393bc97SHong Zhang     current_space->local_remaining -= nzi;
391289bc588SBarry Smith   }
3926cf91177SBarry Smith #if defined(PETSC_USE_INFO)
393464e8d44SSatish Balay   if (ai[n] != 0) {
3941393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
395ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
396ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
397ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
398ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
39905bf559cSSatish Balay   } else {
400ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
40105bf559cSSatish Balay   }
40263ba0a88SBarry Smith #endif
4032fb3ed76SBarry Smith 
404898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
405898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
4061987afe7SBarry Smith 
4071393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
4081393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
409a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
4101393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
41135baf27eSSatish Balay   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
412289bc588SBarry Smith 
413289bc588SBarry Smith   /* put together the new matrix */
414719d5645SBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
415719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
416719d5645SBarry Smith   b    = (Mat_SeqAIJ*)(B)->data;
417e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
418e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
4197c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
4201393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
4211393bc97SHong Zhang   b->j          = bj;
4221393bc97SHong Zhang   b->i          = bi;
4231393bc97SHong Zhang   b->diag       = bdiag;
424416022c9SBarry Smith   b->ilen       = 0;
425416022c9SBarry Smith   b->imax       = 0;
426416022c9SBarry Smith   b->row        = isrow;
427416022c9SBarry Smith   b->col        = iscol;
428c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
429c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
43082bf6240SBarry Smith   b->icol       = isicol;
43187828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
4321393bc97SHong Zhang 
4331393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
434719d5645SBarry Smith   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
4351393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4367fd98bd6SLois Curfman McInnes 
437719d5645SBarry Smith   (B)->factor                 =  MAT_FACTOR_LU;
438719d5645SBarry Smith   (B)->info.factor_mallocs    = reallocs;
439719d5645SBarry Smith   (B)->info.fill_ratio_given  = f;
440703deb49SSatish Balay 
441e93ccc4dSBarry Smith   if (ai[n] != 0) {
442719d5645SBarry Smith     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
443eea2913aSSatish Balay   } else {
444719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
445eea2913aSSatish Balay   }
446719d5645SBarry Smith   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
447719d5645SBarry Smith   (B)->ops->solve            = MatSolve_SeqAIJ;
448719d5645SBarry Smith   (B)->ops->solvetranspose   = MatSolveTranspose_SeqAIJ;
449db4efbfdSBarry Smith   /* switch to inodes if appropriate */
450719d5645SBarry Smith   ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr);
4513a40ed3dSBarry Smith   PetscFunctionReturn(0);
452289bc588SBarry Smith }
4531393bc97SHong Zhang 
4545b5bc046SBarry Smith /*
4555b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4565b5bc046SBarry Smith */
4575b5bc046SBarry Smith #undef __FUNCT__
4585b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4595b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4605b5bc046SBarry Smith {
4615b5bc046SBarry Smith   PetscErrorCode ierr;
4625b5bc046SBarry Smith   PetscTruth     flg;
4635b5bc046SBarry Smith 
4645b5bc046SBarry Smith   PetscFunctionBegin;
4655b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4665b5bc046SBarry Smith   if (flg) {
4675b5bc046SBarry Smith     PetscViewer viewer;
4685b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4695b5bc046SBarry Smith 
4705b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4717adad957SLisandro Dalcin     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4725b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4735b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4745b5bc046SBarry Smith   }
4755b5bc046SBarry Smith   PetscFunctionReturn(0);
4765b5bc046SBarry Smith }
4775b5bc046SBarry Smith 
478db4efbfdSBarry Smith extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec);
479db4efbfdSBarry Smith 
4800a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4814a2ae208SSatish Balay #undef __FUNCT__
4824a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
4830481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
484289bc588SBarry Smith {
485719d5645SBarry Smith   Mat            C=B;
486416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
48782bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4886849ba73SBarry Smith   PetscErrorCode ierr;
4895d0c19d7SBarry Smith   const PetscInt  *r,*ic,*ics;
4905d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*bi=b->i,*bj=b->j;
4915d0c19d7SBarry Smith   PetscInt       *ajtmp,*bjtmp,nz,row;
492d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
49354f21887SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,*rtmps;
49454f21887SBarry Smith   MatScalar      *v,*pv;
4956a7c8fc2SHong Zhang   PetscScalar    d;
4966a7c8fc2SHong Zhang   PetscReal      rs;
497b3bf805bSHong Zhang   LUShift_Ctx    sctx;
49842898933SHong Zhang   PetscInt       newshift,*ddiag;
499289bc588SBarry Smith 
5003a40ed3dSBarry Smith   PetscFunctionBegin;
50178b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
50278b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
50387828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
50487828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
505010a20caSHong Zhang   rtmps = rtmp; ics = ic;
506289bc588SBarry Smith 
507ada7143aSSatish Balay   sctx.shift_top      = 0;
508ada7143aSSatish Balay   sctx.nshift_max     = 0;
509ada7143aSSatish Balay   sctx.shift_lo       = 0;
510ada7143aSSatish Balay   sctx.shift_hi       = 0;
5110481f469SBarry Smith   sctx.shift_fraction = 0;
512ada7143aSSatish Balay 
513118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
5141a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
5159deb78dcSHong Zhang     PetscInt *aai = a->i;
51642898933SHong Zhang     ddiag          = a->diag;
517118f5789SBarry Smith     sctx.shift_top = 0;
5186cc28720Svictorle     for (i=0; i<n; i++) {
5199f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
5209f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
5211a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
522010a20caSHong Zhang       v  = a->a+aai[i];
523e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
524e24cfe64SHong Zhang       for (j=0; j<nz; j++)
5251a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
5261a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
5276cc28720Svictorle     }
528c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
529118f5789SBarry Smith     sctx.shift_top    *= 1.1;
530d2276718SHong Zhang     sctx.nshift_max   = 5;
531d2276718SHong Zhang     sctx.shift_lo     = 0.;
532d2276718SHong Zhang     sctx.shift_hi     = 1.;
533a0a356efSHong Zhang   }
534a0a356efSHong Zhang 
535a0a356efSHong Zhang   sctx.shift_amount = 0;
536a0a356efSHong Zhang   sctx.nshift       = 0;
537085256b3SBarry Smith   do {
538d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
539289bc588SBarry Smith     for (i=0; i<n; i++){
540b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
541b3bf805bSHong Zhang       bjtmp = bj + bi[i];
542b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
543289bc588SBarry Smith 
544289bc588SBarry Smith       /* load in initial (unfactored row) */
545416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
546b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
547010a20caSHong Zhang       v     = a->a + a->i[r[i]];
548085256b3SBarry Smith       for (j=0; j<nz; j++) {
549b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
550085256b3SBarry Smith       }
551d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
552289bc588SBarry Smith 
553b3bf805bSHong Zhang       row = *bjtmp++;
554f2582111SSatish Balay       while  (row < i) {
5558c37ef55SBarry Smith         pc = rtmp + row;
5568c37ef55SBarry Smith         if (*pc != 0.0) {
557010a20caSHong Zhang           pv         = b->a + diag_offset[row];
558010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
55935aab85fSBarry Smith           multiplier = *pc / *pv++;
5608c37ef55SBarry Smith           *pc        = multiplier;
561b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
562f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
563efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
564289bc588SBarry Smith         }
565b3bf805bSHong Zhang         row = *bjtmp++;
566289bc588SBarry Smith       }
567416022c9SBarry Smith       /* finished row so stick it into b->a */
568b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
569b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
570b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
571b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
572d3d32019SBarry Smith       rs   = 0.0;
573d3d32019SBarry Smith       for (j=0; j<nz; j++) {
574d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
575d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
576d3d32019SBarry Smith       }
577d2276718SHong Zhang 
578b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
579d2276718SHong Zhang       sctx.rs  = rs;
580d2276718SHong Zhang       sctx.pv  = pv[diag];
581c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
5825b5bc046SBarry Smith       if (newshift == 1) break;
58335aab85fSBarry Smith     }
584d2276718SHong Zhang 
5850481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5866cc28720Svictorle       /*
5876c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5886cc28720Svictorle        * then try lower shift
5896cc28720Svictorle        */
5900481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
5910481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
5920481f469SBarry Smith       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
593d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
594d2276718SHong Zhang       sctx.nshift++;
5956cc28720Svictorle     }
596d2276718SHong Zhang   } while (sctx.lushift);
5970f11f9aeSSatish Balay 
59835aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
59935aab85fSBarry Smith   for (i=0; i<n; i++) {
600010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
60135aab85fSBarry Smith   }
602606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
60378b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
60478b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
605db4efbfdSBarry Smith   if (b->inode.use) {
606db4efbfdSBarry Smith     C->ops->solve   = MatSolve_Inode;
607db4efbfdSBarry Smith   } else {
6088d8c7c43SBarry Smith     PetscTruth row_identity, col_identity;
6098d8c7c43SBarry Smith     ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
6108d8c7c43SBarry Smith     ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
6118d8c7c43SBarry Smith     if (row_identity && col_identity) {
6128d8c7c43SBarry Smith       C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering;
6138d8c7c43SBarry Smith     } else {
614db4efbfdSBarry Smith       C->ops->solve   = MatSolve_SeqAIJ;
615db4efbfdSBarry Smith     }
6168d8c7c43SBarry Smith   }
617db4efbfdSBarry Smith   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
618db4efbfdSBarry Smith   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
619db4efbfdSBarry Smith   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
620de26497bSBarry Smith   C->ops->matsolve           = MatMatSolve_SeqAIJ;
621c456f294SBarry Smith   C->assembled    = PETSC_TRUE;
6225c9eb25fSBarry Smith   C->preallocated = PETSC_TRUE;
623d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
624d2276718SHong Zhang   if (sctx.nshift){
625e738e422SBarry Smith      if (info->shiftpd) {
6260481f469SBarry 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);
627e738e422SBarry Smith     } else if (info->shiftnz) {
628e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
6296cc28720Svictorle     }
6300697b6caSHong Zhang   }
6313a40ed3dSBarry Smith   PetscFunctionReturn(0);
632289bc588SBarry Smith }
6330889b5dcSvictorle 
634137fb511SHong Zhang /*
635137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
636137fb511SHong Zhang    Input:
637137fb511SHong Zhang      A - original matrix
638137fb511SHong Zhang    Output;
639137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
640137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
641137fb511SHong Zhang          a->a reordered accordingly with a->j
642137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
643137fb511SHong Zhang */
644137fb511SHong Zhang #undef __FUNCT__
645137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
6460481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
647137fb511SHong Zhang {
648b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
649b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
650137fb511SHong Zhang   PetscErrorCode ierr;
6515d0c19d7SBarry Smith   const PetscInt *r,*ic,*ics;
6525d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
6535d0c19d7SBarry Smith   PetscInt       *ajtmp,nz,row;
654b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
655a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
656a77337e4SBarry Smith   MatScalar      *v,*pv;
657137fb511SHong Zhang   PetscReal      rs;
658137fb511SHong Zhang   LUShift_Ctx    sctx;
659b051339dSHong Zhang   PetscInt       newshift;
660137fb511SHong Zhang 
661137fb511SHong Zhang   PetscFunctionBegin;
662719d5645SBarry Smith   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
663137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
664137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
665137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
666137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
667b051339dSHong Zhang   ics = ic;
668137fb511SHong Zhang 
669137fb511SHong Zhang   sctx.shift_top      = 0;
670137fb511SHong Zhang   sctx.nshift_max     = 0;
671137fb511SHong Zhang   sctx.shift_lo       = 0;
672137fb511SHong Zhang   sctx.shift_hi       = 0;
6730481f469SBarry Smith   sctx.shift_fraction = 0;
674137fb511SHong Zhang 
675137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
676137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
677137fb511SHong Zhang     sctx.shift_top = 0;
678137fb511SHong Zhang     for (i=0; i<n; i++) {
679137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
680b051339dSHong Zhang       d  = (a->a)[diag[i]];
681137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
682b051339dSHong Zhang       v  = a->a+ai[i];
683b051339dSHong Zhang       nz = ai[i+1] - ai[i];
684137fb511SHong Zhang       for (j=0; j<nz; j++)
685137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
686137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
687137fb511SHong Zhang     }
688137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
689137fb511SHong Zhang     sctx.shift_top    *= 1.1;
690137fb511SHong Zhang     sctx.nshift_max   = 5;
691137fb511SHong Zhang     sctx.shift_lo     = 0.;
692137fb511SHong Zhang     sctx.shift_hi     = 1.;
693137fb511SHong Zhang   }
694137fb511SHong Zhang 
695137fb511SHong Zhang   sctx.shift_amount = 0;
696137fb511SHong Zhang   sctx.nshift       = 0;
697137fb511SHong Zhang   do {
698137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
699137fb511SHong Zhang     for (i=0; i<n; i++){
700b051339dSHong Zhang       /* load in initial unfactored row */
701b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
702b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
703b051339dSHong Zhang       v     = a->a + ai[r[i]];
704137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
705b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
706137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
707137fb511SHong Zhang 
708b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
709137fb511SHong Zhang       for (j=0; j<nz; j++) {
710137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
711b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
712137fb511SHong Zhang       }
713137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
714137fb511SHong Zhang 
715b051339dSHong Zhang       row = *ajtmp++;
716137fb511SHong Zhang       while  (row < i) {
717137fb511SHong Zhang         pc = rtmp + row;
718137fb511SHong Zhang         if (*pc != 0.0) {
719b051339dSHong Zhang           pv         = a->a + diag[r[row]];
720b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
721137fb511SHong Zhang 
722137fb511SHong Zhang           multiplier = *pc / *pv++;
723137fb511SHong Zhang           *pc        = multiplier;
724b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
725b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
726137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
727137fb511SHong Zhang         }
728b051339dSHong Zhang         row = *ajtmp++;
729137fb511SHong Zhang       }
730b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
731b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
732b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
733b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
734b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
735137fb511SHong Zhang 
736137fb511SHong Zhang       rs   = 0.0;
737137fb511SHong Zhang       for (j=0; j<nz; j++) {
738b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
739b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
740137fb511SHong Zhang       }
741137fb511SHong Zhang 
742137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
743137fb511SHong Zhang       sctx.rs  = rs;
744b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
745c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
746137fb511SHong Zhang       if (newshift == 1) break;
747137fb511SHong Zhang     }
748137fb511SHong Zhang 
7490481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
750137fb511SHong Zhang       /*
751137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
752137fb511SHong Zhang        * then try lower shift
753137fb511SHong Zhang        */
7540481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
7550481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
7560481f469SBarry Smith       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
757137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
758137fb511SHong Zhang       sctx.nshift++;
759137fb511SHong Zhang     }
760137fb511SHong Zhang   } while (sctx.lushift);
761137fb511SHong Zhang 
762137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
763137fb511SHong Zhang   for (i=0; i<n; i++) {
764b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
765137fb511SHong Zhang   }
766137fb511SHong Zhang 
767137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
768137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
769137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
770b051339dSHong Zhang   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
771db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
772db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
773db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
774b051339dSHong Zhang   A->assembled = PETSC_TRUE;
7755c9eb25fSBarry Smith   A->preallocated = PETSC_TRUE;
776d0f46423SBarry Smith   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
777137fb511SHong Zhang   if (sctx.nshift){
778e738e422SBarry Smith     if (info->shiftpd) {
7790481f469SBarry 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);
780e738e422SBarry Smith     } else if (info->shiftnz) {
781e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
782137fb511SHong Zhang     }
783137fb511SHong Zhang   }
784137fb511SHong Zhang   PetscFunctionReturn(0);
785137fb511SHong Zhang }
786137fb511SHong Zhang 
7870a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7884a2ae208SSatish Balay #undef __FUNCT__
7894a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
7900481f469SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
791da3a660dSBarry Smith {
792dfbe8321SBarry Smith   PetscErrorCode ierr;
793416022c9SBarry Smith   Mat            C;
794416022c9SBarry Smith 
7953a40ed3dSBarry Smith   PetscFunctionBegin;
79643244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
797719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
798719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
799db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
800db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
801273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
80252e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
8033a40ed3dSBarry Smith   PetscFunctionReturn(0);
804da3a660dSBarry Smith }
8050a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
8064a2ae208SSatish Balay #undef __FUNCT__
8074a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
808dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
8098c37ef55SBarry Smith {
810416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
811416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
8126849ba73SBarry Smith   PetscErrorCode    ierr;
8135d0c19d7SBarry Smith   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8145d0c19d7SBarry Smith   PetscInt          nz;
8155d0c19d7SBarry Smith   const PetscInt    *rout,*cout,*r,*c;
816dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
817d9fead3dSBarry Smith   const PetscScalar *b;
818dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
8198c37ef55SBarry Smith 
8203a40ed3dSBarry Smith   PetscFunctionBegin;
8213a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
82291bf9adeSBarry Smith 
823d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8241ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
825416022c9SBarry Smith   tmp  = a->solve_work;
8268c37ef55SBarry Smith 
8273b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8283b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
8298c37ef55SBarry Smith 
8308c37ef55SBarry Smith   /* forward solve the lower triangular */
8318c37ef55SBarry Smith   tmp[0] = b[*r++];
832010a20caSHong Zhang   tmps   = tmp;
8338c37ef55SBarry Smith   for (i=1; i<n; i++) {
834010a20caSHong Zhang     v   = aa + ai[i] ;
835010a20caSHong Zhang     vi  = aj + ai[i] ;
83653e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
83753e7425aSBarry Smith     sum = b[*r++];
838ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8398c37ef55SBarry Smith     tmp[i] = sum;
8408c37ef55SBarry Smith   }
8418c37ef55SBarry Smith 
8428c37ef55SBarry Smith   /* backward solve the upper triangular */
8438c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
844010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
845010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
846416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
8478c37ef55SBarry Smith     sum = tmp[i];
848ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
849010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
8508c37ef55SBarry Smith   }
8518c37ef55SBarry Smith 
8523b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8533b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
854d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8551ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
856d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
8588c37ef55SBarry Smith }
859026e39d0SSatish Balay 
8602bebee5dSHong Zhang #undef __FUNCT__
8612bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
8622bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8632bebee5dSHong Zhang {
8642bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8652bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8662bebee5dSHong Zhang   PetscErrorCode  ierr;
8675d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8685d0c19d7SBarry Smith   PetscInt        nz,neq;
8695d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
870dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
871dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
872db4efbfdSBarry Smith   PetscTruth      bisdense,xisdense;
8732bebee5dSHong Zhang 
8742bebee5dSHong Zhang   PetscFunctionBegin;
8752bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8762bebee5dSHong Zhang 
877db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
878db4efbfdSBarry Smith   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
879db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
880db4efbfdSBarry Smith   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
881db4efbfdSBarry Smith 
8822bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8832bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8842bebee5dSHong Zhang 
8852bebee5dSHong Zhang   tmp  = a->solve_work;
8862bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8872bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8882bebee5dSHong Zhang 
889a36b22ccSSatish Balay   for (neq=0; neq<B->cmap->n; neq++){
8902bebee5dSHong Zhang     /* forward solve the lower triangular */
8912bebee5dSHong Zhang     tmp[0] = b[r[0]];
8922bebee5dSHong Zhang     tmps   = tmp;
8932bebee5dSHong Zhang     for (i=1; i<n; i++) {
8942bebee5dSHong Zhang       v   = aa + ai[i] ;
8952bebee5dSHong Zhang       vi  = aj + ai[i] ;
8962bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8972bebee5dSHong Zhang       sum = b[r[i]];
8982bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8992bebee5dSHong Zhang       tmp[i] = sum;
9002bebee5dSHong Zhang     }
9012bebee5dSHong Zhang     /* backward solve the upper triangular */
9022bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
9032bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
9042bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
9052bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
9062bebee5dSHong Zhang       sum = tmp[i];
9072bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
9082bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
9092bebee5dSHong Zhang     }
9102bebee5dSHong Zhang 
9112bebee5dSHong Zhang     b += n;
9122bebee5dSHong Zhang     x += n;
9132bebee5dSHong Zhang   }
9142bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
9152bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
9162bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
9172bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
9183639cd1fSBarry Smith   ierr = PetscLogFlops(B->cmap->n*(2*a->nz - n));CHKERRQ(ierr);
9192bebee5dSHong Zhang   PetscFunctionReturn(0);
9202bebee5dSHong Zhang }
9212bebee5dSHong Zhang 
922137fb511SHong Zhang #undef __FUNCT__
923137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
924137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
925137fb511SHong Zhang {
926137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
927137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
928137fb511SHong Zhang   PetscErrorCode  ierr;
9295d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
9305d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
9315d0c19d7SBarry Smith   PetscInt        nz,row;
932dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
933dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
934137fb511SHong Zhang 
935137fb511SHong Zhang   PetscFunctionBegin;
936137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
937137fb511SHong Zhang 
938137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
939137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
940137fb511SHong Zhang   tmp  = a->solve_work;
941137fb511SHong Zhang 
942137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
943137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
944137fb511SHong Zhang 
945137fb511SHong Zhang   /* forward solve the lower triangular */
946137fb511SHong Zhang   tmp[0] = b[*r++];
947137fb511SHong Zhang   tmps   = tmp;
948137fb511SHong Zhang   for (row=1; row<n; row++) {
949137fb511SHong Zhang     i   = rout[row]; /* permuted row */
950137fb511SHong Zhang     v   = aa + ai[i] ;
951137fb511SHong Zhang     vi  = aj + ai[i] ;
952137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
953137fb511SHong Zhang     sum = b[*r++];
954137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
955137fb511SHong Zhang     tmp[row] = sum;
956137fb511SHong Zhang   }
957137fb511SHong Zhang 
958137fb511SHong Zhang   /* backward solve the upper triangular */
959137fb511SHong Zhang   for (row=n-1; row>=0; row--){
960137fb511SHong Zhang     i   = rout[row]; /* permuted row */
961137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
962137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
963137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
964137fb511SHong Zhang     sum = tmp[row];
965137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
966137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
967137fb511SHong Zhang   }
968137fb511SHong Zhang 
969137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
970137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
971137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
972137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
973d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
974137fb511SHong Zhang   PetscFunctionReturn(0);
975137fb511SHong Zhang }
976137fb511SHong Zhang 
977930ae53cSBarry Smith /* ----------------------------------------------------------- */
9784a2ae208SSatish Balay #undef __FUNCT__
9794a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
980dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
981930ae53cSBarry Smith {
982930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9836849ba73SBarry Smith   PetscErrorCode    ierr;
984d0f46423SBarry Smith   PetscInt          n = A->rmap->n;
985356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
986356650c2SBarry Smith   PetscScalar       *x;
987356650c2SBarry Smith   const PetscScalar *b;
988dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
989aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
990356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
991dd6ea824SBarry Smith   const MatScalar   *v;
992dd6ea824SBarry Smith   PetscScalar       sum;
993d85afc42SSatish Balay #endif
994930ae53cSBarry Smith 
9953a40ed3dSBarry Smith   PetscFunctionBegin;
9963a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
997930ae53cSBarry Smith 
998356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9991ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1000930ae53cSBarry Smith 
1001aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
10021c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
10031c4feecaSSatish Balay #else
1004930ae53cSBarry Smith   /* forward solve the lower triangular */
1005930ae53cSBarry Smith   x[0] = b[0];
1006930ae53cSBarry Smith   for (i=1; i<n; i++) {
1007930ae53cSBarry Smith     ai_i = ai[i];
1008930ae53cSBarry Smith     v    = aa + ai_i;
1009930ae53cSBarry Smith     vi   = aj + ai_i;
1010930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
1011930ae53cSBarry Smith     sum  = b[i];
1012930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1013930ae53cSBarry Smith     x[i] = sum;
1014930ae53cSBarry Smith   }
1015930ae53cSBarry Smith 
1016930ae53cSBarry Smith   /* backward solve the upper triangular */
1017930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
1018930ae53cSBarry Smith     adiag_i = adiag[i];
1019d708749eSBarry Smith     v       = aa + adiag_i + 1;
1020d708749eSBarry Smith     vi      = aj + adiag_i + 1;
1021930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
1022930ae53cSBarry Smith     sum     = x[i];
1023930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1024930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
1025930ae53cSBarry Smith   }
10261c4feecaSSatish Balay #endif
1027d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
1028356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
10291ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1031930ae53cSBarry Smith }
1032930ae53cSBarry Smith 
10334a2ae208SSatish Balay #undef __FUNCT__
10344a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
1035dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1036da3a660dSBarry Smith {
1037416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1038416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
10396849ba73SBarry Smith   PetscErrorCode  ierr;
10405d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10415d0c19d7SBarry Smith   PetscInt        nz;
10425d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
1043dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
1044dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1045da3a660dSBarry Smith 
10463a40ed3dSBarry Smith   PetscFunctionBegin;
104778b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
1048da3a660dSBarry Smith 
10491ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10501ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1051416022c9SBarry Smith   tmp  = a->solve_work;
1052da3a660dSBarry Smith 
10533b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10543b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1055da3a660dSBarry Smith 
1056da3a660dSBarry Smith   /* forward solve the lower triangular */
1057da3a660dSBarry Smith   tmp[0] = b[*r++];
1058da3a660dSBarry Smith   for (i=1; i<n; i++) {
1059010a20caSHong Zhang     v   = aa + ai[i] ;
1060010a20caSHong Zhang     vi  = aj + ai[i] ;
1061416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
1062da3a660dSBarry Smith     sum = b[*r++];
1063010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1064da3a660dSBarry Smith     tmp[i] = sum;
1065da3a660dSBarry Smith   }
1066da3a660dSBarry Smith 
1067da3a660dSBarry Smith   /* backward solve the upper triangular */
1068da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1069010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
1070010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
1071416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1072da3a660dSBarry Smith     sum = tmp[i];
1073010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1074010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1075da3a660dSBarry Smith     x[*c--] += tmp[i];
1076da3a660dSBarry Smith   }
1077da3a660dSBarry Smith 
10783b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10793b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10801ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10811ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1082efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1083898183ecSLois Curfman McInnes 
10843a40ed3dSBarry Smith   PetscFunctionReturn(0);
1085da3a660dSBarry Smith }
1086da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10874a2ae208SSatish Balay #undef __FUNCT__
10884a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1089dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1090da3a660dSBarry Smith {
1091416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10922235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10936849ba73SBarry Smith   PetscErrorCode  ierr;
10945d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
10955d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10965d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1097dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1098dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1099da3a660dSBarry Smith 
11003a40ed3dSBarry Smith   PetscFunctionBegin;
11011ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11021ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1103416022c9SBarry Smith   tmp  = a->solve_work;
1104da3a660dSBarry Smith 
11052235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11062235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1107da3a660dSBarry Smith 
1108da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
11092235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1110da3a660dSBarry Smith 
1111da3a660dSBarry Smith   /* forward solve the U^T */
1112da3a660dSBarry Smith   for (i=0; i<n; i++) {
1113010a20caSHong Zhang     v   = aa + diag[i] ;
1114010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1115f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1116c38d4ed2SBarry Smith     s1  = tmp[i];
1117ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1118da3a660dSBarry Smith     while (nz--) {
1119010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1120da3a660dSBarry Smith     }
1121c38d4ed2SBarry Smith     tmp[i] = s1;
1122da3a660dSBarry Smith   }
1123da3a660dSBarry Smith 
1124da3a660dSBarry Smith   /* backward solve the L^T */
1125da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1126010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1127010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1128f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1129f1af5d2fSBarry Smith     s1  = tmp[i];
1130da3a660dSBarry Smith     while (nz--) {
1131010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1132da3a660dSBarry Smith     }
1133da3a660dSBarry Smith   }
1134da3a660dSBarry Smith 
1135da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1136da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1137da3a660dSBarry Smith 
11382235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11392235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11401ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11411ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1142da3a660dSBarry Smith 
1143d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap->n);CHKERRQ(ierr);
11443a40ed3dSBarry Smith   PetscFunctionReturn(0);
1145da3a660dSBarry Smith }
1146da3a660dSBarry Smith 
11474a2ae208SSatish Balay #undef __FUNCT__
11484a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1149dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1150da3a660dSBarry Smith {
1151416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
11522235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
11536849ba73SBarry Smith   PetscErrorCode  ierr;
11545d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
11555d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
11565d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1157dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1158dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
11596abc6512SBarry Smith 
11603a40ed3dSBarry Smith   PetscFunctionBegin;
116123bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
11626abc6512SBarry Smith 
11631ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11641ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1165416022c9SBarry Smith   tmp = a->solve_work;
11666abc6512SBarry Smith 
11672235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11682235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11696abc6512SBarry Smith 
11706abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
11712235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
11726abc6512SBarry Smith 
11736abc6512SBarry Smith   /* forward solve the U^T */
11746abc6512SBarry Smith   for (i=0; i<n; i++) {
1175010a20caSHong Zhang     v   = aa + diag[i] ;
1176010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1177f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11786abc6512SBarry Smith     tmp[i] *= *v++;
11796abc6512SBarry Smith     while (nz--) {
1180010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11816abc6512SBarry Smith     }
11826abc6512SBarry Smith   }
11836abc6512SBarry Smith 
11846abc6512SBarry Smith   /* backward solve the L^T */
11856abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1186010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1187010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1188f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11896abc6512SBarry Smith     while (nz--) {
1190010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11916abc6512SBarry Smith     }
11926abc6512SBarry Smith   }
11936abc6512SBarry Smith 
11946abc6512SBarry Smith   /* copy tmp into x according to permutation */
11956abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11966abc6512SBarry Smith 
11972235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11982235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11991ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
12001ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
12016abc6512SBarry Smith 
1202efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
12033a40ed3dSBarry Smith   PetscFunctionReturn(0);
1204da3a660dSBarry Smith }
12059e25ed09SBarry Smith /* ----------------------------------------------------------------*/
12063c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1207719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption);
120808480c60SBarry Smith 
12094a2ae208SSatish Balay #undef __FUNCT__
12104a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
12110481f469SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
12129e25ed09SBarry Smith {
1213416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
12149e25ed09SBarry Smith   IS                 isicol;
12156849ba73SBarry Smith   PetscErrorCode     ierr;
12165d0c19d7SBarry Smith   const PetscInt     *r,*ic;
12175d0c19d7SBarry Smith   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
12185a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
12195a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
12205a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
122177c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1222329f5518SBarry Smith   PetscReal          f;
12235a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12245a8e39fbSHong Zhang   PetscBT            lnkbt;
12255a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1226a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1227a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
122809f38230SBarry Smith   PetscTruth         missing;
12299e25ed09SBarry Smith 
12303a40ed3dSBarry Smith   PetscFunctionBegin;
1231d0f46423SBarry 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);
12326cf997b0SBarry Smith   f             = info->fill;
123338baddfdSBarry Smith   levels        = (PetscInt)info->levels;
123438baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
12354c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
123682bf6240SBarry Smith 
123708480c60SBarry Smith   /* special case that simply copies fill pattern */
1238be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1239be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
124086bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
1241719d5645SBarry Smith     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr);
1242719d5645SBarry Smith     fact->factor = MAT_FACTOR_ILU;
1243719d5645SBarry Smith     (fact)->info.factor_mallocs    = 0;
1244719d5645SBarry Smith     (fact)->info.fill_ratio_given  = info->fill;
1245719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 1.0;
1246719d5645SBarry Smith     b               = (Mat_SeqAIJ*)(fact)->data;
12478ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
124809f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
124908480c60SBarry Smith     b->row              = isrow;
125008480c60SBarry Smith     b->col              = iscol;
125182bf6240SBarry Smith     b->icol             = isicol;
1252719d5645SBarry Smith     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1253c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1254c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1255719d5645SBarry Smith     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1256719d5645SBarry Smith     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
12573a40ed3dSBarry Smith     PetscFunctionReturn(0);
125808480c60SBarry Smith   }
125908480c60SBarry Smith 
1260898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1261898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
12629e25ed09SBarry Smith 
12639e25ed09SBarry Smith   /* get new row pointers */
12645a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12655a8e39fbSHong Zhang   bi[0] = 0;
12665a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
12675a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12685a8e39fbSHong Zhang   bdiag[0]  = 0;
12696cf997b0SBarry Smith 
12705a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
12715a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
12725a8e39fbSHong Zhang 
12735a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
12745a8e39fbSHong Zhang   nlnk = n + 1;
12755a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12765a8e39fbSHong Zhang 
12775a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1278a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12795a8e39fbSHong Zhang   current_space = free_space;
1280a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12815a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12825a8e39fbSHong Zhang 
12835a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12845a8e39fbSHong Zhang     nzi = 0;
12856cf997b0SBarry Smith     /* copy current row into linked list */
12865a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
1287*3b4a8b6dSBarry Smith     if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
12885a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12895a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12905a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12915a8e39fbSHong Zhang     nzi += nlnk;
12926cf997b0SBarry Smith 
12936cf997b0SBarry Smith     /* make sure diagonal entry is included */
12945a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12956cf997b0SBarry Smith       fm = n;
12965a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12975a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12985a8e39fbSHong Zhang       lnk[fm]    = i;
12995a8e39fbSHong Zhang       lnk_lvl[i] = 0;
13005a8e39fbSHong Zhang       nzi++; dcount++;
13016cf997b0SBarry Smith     }
13026cf997b0SBarry Smith 
13035a8e39fbSHong Zhang     /* add pivot rows into the active row */
13045a8e39fbSHong Zhang     nzbd = 0;
13055a8e39fbSHong Zhang     prow = lnk[n];
13065a8e39fbSHong Zhang     while (prow < i) {
13075a8e39fbSHong Zhang       nnz      = bdiag[prow];
13085a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
13095a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
13105a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
13115a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
13125a8e39fbSHong Zhang       nzi += nlnk;
13135a8e39fbSHong Zhang       prow = lnk[prow];
13145a8e39fbSHong Zhang       nzbd++;
131556beaf04SBarry Smith     }
13165a8e39fbSHong Zhang     bdiag[i] = nzbd;
13175a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1318ecf371e4SBarry Smith 
13195a8e39fbSHong Zhang     /* if free space is not available, make more free space */
13205a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
13215a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1322a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1323a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
13245a8e39fbSHong Zhang       reallocs++;
13255a8e39fbSHong Zhang     }
1326ecf371e4SBarry Smith 
13275a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
13285a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13295a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
13305a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
13315a8e39fbSHong Zhang 
13325a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
13335a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
133477431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1335d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
13366cf997b0SBarry Smith     }
13375a8e39fbSHong Zhang 
13385a8e39fbSHong Zhang     current_space->array           += nzi;
13395a8e39fbSHong Zhang     current_space->local_used      += nzi;
13405a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
13415a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
13425a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
13435a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
13449e25ed09SBarry Smith   }
13455a8e39fbSHong Zhang 
1346898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1347898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
13485a8e39fbSHong Zhang 
13495a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
13505a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1351a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13525a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1353a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13545a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
13559e25ed09SBarry Smith 
13566cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1357f64065bbSBarry Smith   {
13585a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1359ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1360ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1361ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1362ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1363335d9088SBarry Smith     if (diagonal_fill) {
1364ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1365335d9088SBarry Smith     }
1366f64065bbSBarry Smith   }
136763ba0a88SBarry Smith #endif
1368bbb6d6a8SBarry Smith 
13699e25ed09SBarry Smith   /* put together the new matrix */
137007dbf95fSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1371719d5645SBarry Smith   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1372719d5645SBarry Smith   b = (Mat_SeqAIJ*)(fact)->data;
1373e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1374e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13757c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13765a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1377b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13785a8e39fbSHong Zhang   b->j          = bj;
13795a8e39fbSHong Zhang   b->i          = bi;
13805a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13815a8e39fbSHong Zhang   b->diag       = bdiag;
1382416022c9SBarry Smith   b->ilen       = 0;
1383416022c9SBarry Smith   b->imax       = 0;
1384416022c9SBarry Smith   b->row        = isrow;
1385416022c9SBarry Smith   b->col        = iscol;
1386c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1387c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
138882bf6240SBarry Smith   b->icol       = isicol;
138987828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1390416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13915a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
1392719d5645SBarry Smith   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13935a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
1394719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1395719d5645SBarry Smith   (fact)->info.fill_ratio_given  = f;
1396719d5645SBarry Smith   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1397719d5645SBarry Smith   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1398719d5645SBarry Smith   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
13993a40ed3dSBarry Smith   PetscFunctionReturn(0);
14009e25ed09SBarry Smith }
1401d5d45c9bSBarry Smith 
14027c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
1403a6175056SHong Zhang #undef __FUNCT__
1404a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
14050481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1406a6175056SHong Zhang {
1407719d5645SBarry Smith   Mat            C = B;
14080968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
14092ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
14109bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
14112ed133dbSHong Zhang   PetscErrorCode ierr;
14125d0c19d7SBarry Smith   const PetscInt *rip,*riip;
14135d0c19d7SBarry Smith   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
14142ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
14152ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1416622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1417540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1418fbf22428SSatish Balay   PetscReal      shiftpd;
1419540703acSHong Zhang   ChShift_Ctx    sctx;
1420540703acSHong Zhang   PetscInt       newshift;
1421db4efbfdSBarry Smith   PetscTruth     perm_identity;
1422d5d45c9bSBarry Smith 
1423a6175056SHong Zhang   PetscFunctionBegin;
1424117f1891Stmunson 
1425540703acSHong Zhang   shiftnz   = info->shiftnz;
1426540703acSHong Zhang   shiftpd   = info->shiftpd;
1427ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1428ee45ca4aSHong Zhang 
14292ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
14309bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
14312ed133dbSHong Zhang 
14322ed133dbSHong Zhang   /* initialization */
14332ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
14342ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
14352ed133dbSHong Zhang   jl   = il + mbs;
14362ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
14372ed133dbSHong Zhang 
1438540703acSHong Zhang   sctx.shift_amount = 0;
1439540703acSHong Zhang   sctx.nshift       = 0;
14402ed133dbSHong Zhang   do {
1441540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
14422ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14432ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1444a6175056SHong Zhang     }
14452ed133dbSHong Zhang 
14462ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1447c05c3958SHong Zhang       bval = ba + bi[k];
14482ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
14492ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
14502ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
14519bfd6278SHong Zhang         col = riip[aj[j]];
14522ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
14532ed133dbSHong Zhang           rtmp[col] = aa[j];
1454c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
14552ed133dbSHong Zhang         }
14562ed133dbSHong Zhang       }
145739e3d630SHong Zhang       /* shift the diagonal of the matrix */
1458540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
14592ed133dbSHong Zhang 
14602ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
14612ed133dbSHong Zhang       dk = rtmp[k];
14622ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
14632ed133dbSHong Zhang 
14642ed133dbSHong Zhang       while (i < k){
14652ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
14662ed133dbSHong Zhang 
14672ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
14682ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
14692ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
14702ed133dbSHong Zhang         dk += uikdi*ba[ili];
14712ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14722ed133dbSHong Zhang 
14732ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14742ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14752ed133dbSHong Zhang         if (jmin < jmax){
14762ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14772ed133dbSHong Zhang           /* update il and jl for row i */
14782ed133dbSHong Zhang           il[i] = jmin;
14792ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14802ed133dbSHong Zhang         }
14812ed133dbSHong Zhang         i = nexti;
14822ed133dbSHong Zhang       }
14832ed133dbSHong Zhang 
14840697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14850697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14860697b6caSHong Zhang       rs   = 0.0;
14873c9e8598SHong Zhang       jmin = bi[k]+1;
14883c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14893c9e8598SHong Zhang       bcol = bj + jmin;
14903c9e8598SHong Zhang       while (nz--){
149189f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
149289f845c8SHong Zhang         bcol++;
14933c9e8598SHong Zhang       }
1494540703acSHong Zhang 
1495540703acSHong Zhang       sctx.rs = rs;
1496540703acSHong Zhang       sctx.pv = dk;
14975b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1498117f1891Stmunson 
1499117f1891Stmunson       if (newshift == 1) {
1500117f1891Stmunson         if (!sctx.shift_amount) {
1501117f1891Stmunson           sctx.shift_amount = 1e-5;
1502117f1891Stmunson         }
1503117f1891Stmunson         break;
1504117f1891Stmunson       }
15053c9e8598SHong Zhang 
15063c9e8598SHong Zhang       /* copy data into U(k,:) */
150739e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
150839e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
150939e3d630SHong Zhang       if (jmin < jmax) {
151039e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
151139e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
15123c9e8598SHong Zhang         }
151339e3d630SHong Zhang         /* add the k-th row into il and jl */
15143c9e8598SHong Zhang         il[k] = jmin;
15153c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
15163c9e8598SHong Zhang       }
15173c9e8598SHong Zhang     }
1518540703acSHong Zhang   } while (sctx.chshift);
15193c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
15203c9e8598SHong Zhang 
152139e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
15229bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1523db4efbfdSBarry Smith 
1524db4efbfdSBarry Smith   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1525db4efbfdSBarry Smith   if (perm_identity){
1526719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1527719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1528719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1529719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1530db4efbfdSBarry Smith   } else {
1531719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1532719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1533719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1534719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1535db4efbfdSBarry Smith   }
1536db4efbfdSBarry Smith 
15373c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
15383c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1539d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1540540703acSHong Zhang   if (sctx.nshift){
1541540703acSHong Zhang     if (shiftnz) {
15421e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1543540703acSHong Zhang     } else if (shiftpd) {
15441e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
15450697b6caSHong Zhang     }
15463c9e8598SHong Zhang   }
15473c9e8598SHong Zhang   PetscFunctionReturn(0);
15483c9e8598SHong Zhang }
1549a6175056SHong Zhang 
1550a6175056SHong Zhang #undef __FUNCT__
1551a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
15520481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1553a6175056SHong Zhang {
15540968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1555ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1556dfbe8321SBarry Smith   PetscErrorCode     ierr;
155758ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
15585d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
15595d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
1560ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
156158ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
15625a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1563ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1564a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1565a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
15667a48dd6fSHong Zhang   PetscBT            lnkbt;
1567b635d36bSHong Zhang   IS                 iperm;
1568a6175056SHong Zhang 
1569a6175056SHong Zhang   PetscFunctionBegin;
1570d0f46423SBarry 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);
157158ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
157258ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1573653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1574b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
15757a48dd6fSHong Zhang 
157639e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
157739e3d630SHong Zhang   ui[0] = 0;
157839e3d630SHong Zhang 
1579b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1580622e2028SHong Zhang   if (!levels && perm_identity) {
158158ebbce7SBarry Smith 
1582ed59401aSHong Zhang     for (i=0; i<am; i++) {
158339e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1584ed59401aSHong Zhang     }
158539e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
158639e3d630SHong Zhang     cols = uj;
1587ed59401aSHong Zhang     for (i=0; i<am; i++) {
1588ed59401aSHong Zhang       aj    = a->j + a->diag[i];
158939e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
159039e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1591ed59401aSHong Zhang     }
159239e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1593b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1594b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1595b635d36bSHong Zhang 
15967a48dd6fSHong Zhang     /* initialization */
15975a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15987a48dd6fSHong Zhang 
15997a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
16007a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
16011393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
16027a48dd6fSHong Zhang     il         = jl + am;
16037a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
16047a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
16057a48dd6fSHong Zhang     for (i=0; i<am; i++){
16067a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
16077a48dd6fSHong Zhang     }
16087a48dd6fSHong Zhang 
16097a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
16107a48dd6fSHong Zhang     nlnk = am + 1;
16112ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16127a48dd6fSHong Zhang 
16137a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1614a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
16157a48dd6fSHong Zhang     current_space = free_space;
1616a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
16177a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
16187a48dd6fSHong Zhang 
16197a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
16207a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
16217a48dd6fSHong Zhang       nzk   = 0;
1622622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1623*3b4a8b6dSBarry Smith       if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
1624622e2028SHong Zhang       ncols_upper = 0;
1625622e2028SHong Zhang       for (j=0; j<ncols; j++){
1626b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1627b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
16285a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1629622e2028SHong Zhang           ncols_upper++;
1630622e2028SHong Zhang         }
1631622e2028SHong Zhang       }
1632b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16337a48dd6fSHong Zhang       nzk += nlnk;
16347a48dd6fSHong Zhang 
16357a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
16367a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
16377a48dd6fSHong Zhang 
16387a48dd6fSHong Zhang       while (prow < k){
16397a48dd6fSHong Zhang         nextprow = jl[prow];
16407a48dd6fSHong Zhang 
16417a48dd6fSHong Zhang         /* merge prow into k-th row */
16427a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
16437a48dd6fSHong Zhang         jmax = ui[prow+1];
16447a48dd6fSHong Zhang         ncols = jmax-jmin;
1645ed59401aSHong Zhang         i     = jmin - ui[prow];
1646ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
16475a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
16485a8e39fbSHong Zhang         j     = *(uj - 1);
16495a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
16507a48dd6fSHong Zhang         nzk += nlnk;
16517a48dd6fSHong Zhang 
16527a48dd6fSHong Zhang         /* update il and jl for prow */
16537a48dd6fSHong Zhang         if (jmin < jmax){
16547a48dd6fSHong Zhang           il[prow] = jmin;
16557a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
16567a48dd6fSHong Zhang         }
16577a48dd6fSHong Zhang         prow = nextprow;
16587a48dd6fSHong Zhang       }
16597a48dd6fSHong Zhang 
16607a48dd6fSHong Zhang       /* if free space is not available, make more free space */
16617a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
16627a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
16637a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1664a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1665a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
16667a48dd6fSHong Zhang         reallocs++;
16677a48dd6fSHong Zhang       }
16687a48dd6fSHong Zhang 
16697a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1670b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
16712ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
16727a48dd6fSHong Zhang 
16737a48dd6fSHong Zhang       /* add the k-th row into il and jl */
167439e3d630SHong Zhang       if (nzk > 1){
16757a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
16767a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
16777a48dd6fSHong Zhang         il[k] = ui[k] + 1;
16787a48dd6fSHong Zhang       }
16797a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
16807a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
16817a48dd6fSHong Zhang 
16827a48dd6fSHong Zhang       current_space->array           += nzk;
16837a48dd6fSHong Zhang       current_space->local_used      += nzk;
16847a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16857a48dd6fSHong Zhang 
16867a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16877a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16887a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16897a48dd6fSHong Zhang 
16907a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16917a48dd6fSHong Zhang     }
16927a48dd6fSHong Zhang 
16936cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16947a48dd6fSHong Zhang     if (ai[am] != 0) {
169539e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1696ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1697ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1698ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16997a48dd6fSHong Zhang     } else {
1700ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
17017a48dd6fSHong Zhang     }
170263ba0a88SBarry Smith #endif
17037a48dd6fSHong Zhang 
17047a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1705b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
17067a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
17075a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
17087a48dd6fSHong Zhang 
17097a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
17107a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1711a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
17122ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1713a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
17147a48dd6fSHong Zhang 
171539e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
171639e3d630SHong Zhang 
17177a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
17187a48dd6fSHong Zhang 
1719719d5645SBarry Smith   b    = (Mat_SeqSBAIJ*)(fact)->data;
17207a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
17217a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
17227a48dd6fSHong Zhang   b->j    = uj;
17237a48dd6fSHong Zhang   b->i    = ui;
17247a48dd6fSHong Zhang   b->diag = 0;
17257a48dd6fSHong Zhang   b->ilen = 0;
17267a48dd6fSHong Zhang   b->imax = 0;
17277a48dd6fSHong Zhang   b->row  = perm;
1728b635d36bSHong Zhang   b->col  = perm;
1729b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1730b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1731b635d36bSHong Zhang   b->icol = iperm;
17327a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
17337a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1734719d5645SBarry Smith   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
17357a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1736e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1737e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
17387a48dd6fSHong Zhang 
1739719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1740719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
17417a48dd6fSHong Zhang   if (ai[am] != 0) {
1742719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
17437a48dd6fSHong Zhang   } else {
1744719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
17457a48dd6fSHong Zhang   }
1746719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1747a6175056SHong Zhang   PetscFunctionReturn(0);
1748a6175056SHong Zhang }
1749d5d45c9bSBarry Smith 
1750f76d2b81SHong Zhang #undef __FUNCT__
1751f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
17520481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1753f76d2b81SHong Zhang {
1754f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
175510c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
1756dfbe8321SBarry Smith   PetscErrorCode     ierr;
1757f76d2b81SHong Zhang   PetscTruth         perm_identity;
175810c27e3dSHong Zhang   PetscReal          fill = info->fill;
17595d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
17605d0c19d7SBarry Smith   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
176110c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
17622ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1763a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
176410c27e3dSHong Zhang   PetscBT            lnkbt;
17652ed133dbSHong Zhang   IS                 iperm;
1766f76d2b81SHong Zhang 
1767f76d2b81SHong Zhang   PetscFunctionBegin;
1768d0f46423SBarry 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);
17692ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1770f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
17712ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
17722ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
17739bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
177410c27e3dSHong Zhang 
177510c27e3dSHong Zhang   /* initialization */
17761393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
177710c27e3dSHong Zhang   ui[0] = 0;
177810c27e3dSHong Zhang 
177910c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17801393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17811393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17821393bc97SHong Zhang   il     = jl + am;
17831393bc97SHong Zhang   cols   = il + am;
17841393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17851393bc97SHong Zhang   for (i=0; i<am; i++){
17861393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1787f76d2b81SHong Zhang   }
178810c27e3dSHong Zhang 
178910c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17901393bc97SHong Zhang   nlnk = am + 1;
17911393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
179210c27e3dSHong Zhang 
17931393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1794a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
179510c27e3dSHong Zhang   current_space = free_space;
179610c27e3dSHong Zhang 
17971393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
179810c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
179910c27e3dSHong Zhang     nzk   = 0;
18002ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
1801*3b4a8b6dSBarry Smith     if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
18022ed133dbSHong Zhang     ncols_upper = 0;
18032ed133dbSHong Zhang     for (j=0; j<ncols; j++){
18049bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
18052ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
18062ed133dbSHong Zhang         cols[ncols_upper] = i;
18072ed133dbSHong Zhang         ncols_upper++;
18082ed133dbSHong Zhang       }
18092ed133dbSHong Zhang     }
18101393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
181110c27e3dSHong Zhang     nzk += nlnk;
181210c27e3dSHong Zhang 
181310c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
181410c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
181510c27e3dSHong Zhang 
181610c27e3dSHong Zhang     while (prow < k){
181710c27e3dSHong Zhang       nextprow = jl[prow];
181810c27e3dSHong Zhang       /* merge prow into k-th row */
18191393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
182010c27e3dSHong Zhang       jmax = ui[prow+1];
182110c27e3dSHong Zhang       ncols = jmax-jmin;
18221393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
18235a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
182410c27e3dSHong Zhang       nzk += nlnk;
182510c27e3dSHong Zhang 
182610c27e3dSHong Zhang       /* update il and jl for prow */
182710c27e3dSHong Zhang       if (jmin < jmax){
182810c27e3dSHong Zhang         il[prow] = jmin;
18292ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
183010c27e3dSHong Zhang       }
183110c27e3dSHong Zhang       prow = nextprow;
183210c27e3dSHong Zhang     }
183310c27e3dSHong Zhang 
183410c27e3dSHong Zhang     /* if free space is not available, make more free space */
183510c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
18361393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
183710c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1838a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
183910c27e3dSHong Zhang       reallocs++;
184010c27e3dSHong Zhang     }
184110c27e3dSHong Zhang 
184210c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
18431393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
184410c27e3dSHong Zhang 
184510c27e3dSHong Zhang     /* add the k-th row into il and jl */
184610c27e3dSHong Zhang     if (nzk-1 > 0){
18471393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
184810c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
184910c27e3dSHong Zhang       il[k] = ui[k] + 1;
185010c27e3dSHong Zhang     }
18512ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
185210c27e3dSHong Zhang     current_space->array           += nzk;
185310c27e3dSHong Zhang     current_space->local_used      += nzk;
185410c27e3dSHong Zhang     current_space->local_remaining -= nzk;
185510c27e3dSHong Zhang 
185610c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
185710c27e3dSHong Zhang   }
185810c27e3dSHong Zhang 
18596cf91177SBarry Smith #if defined(PETSC_USE_INFO)
18601393bc97SHong Zhang   if (ai[am] != 0) {
18611393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1862ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1863ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1864ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
186510c27e3dSHong Zhang   } else {
1866ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
186710c27e3dSHong Zhang   }
186863ba0a88SBarry Smith #endif
186910c27e3dSHong Zhang 
187010c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
18719bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
187210c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
187310c27e3dSHong Zhang 
187410c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18751393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1876a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
187710c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
187810c27e3dSHong Zhang 
187910c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
188010c27e3dSHong Zhang 
1881719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
188210c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1883e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1884e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18851393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
188610c27e3dSHong Zhang   b->j    = uj;
188710c27e3dSHong Zhang   b->i    = ui;
188810c27e3dSHong Zhang   b->diag = 0;
188910c27e3dSHong Zhang   b->ilen = 0;
189010c27e3dSHong Zhang   b->imax = 0;
189110c27e3dSHong Zhang   b->row  = perm;
18929bfd6278SHong Zhang   b->col  = perm;
18939bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18949bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18959bfd6278SHong Zhang   b->icol = iperm;
189610c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18971393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1898719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18991393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
190010c27e3dSHong Zhang 
1901719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1902719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
19031393bc97SHong Zhang   if (ai[am] != 0) {
1904719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
190510c27e3dSHong Zhang   } else {
1906719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
190710c27e3dSHong Zhang   }
1908719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1909f76d2b81SHong Zhang   PetscFunctionReturn(0);
1910f76d2b81SHong Zhang }
1911