xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*);
26e34fafa9SBarry Smith #endif
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
63*5d0c19d7SBarry Smith   const PetscInt *c,*r,*ic;
64*5d0c19d7SBarry 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;
6954a8161fSBarry Smith   PetscReal      af;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   PetscFunctionBegin;
7286bacbd2SBarry Smith 
7386bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7438baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7586bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7633259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7738baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7838baddfdSBarry Smith   jmax    = (PetscInt)(info->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);
117*5d0c19d7SBarry Smith     ierr = PetscMalloc2(n,PetscInt,&cc,n,PetscInt,&cr);CHKERRQ(ierr);
118c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
119*5d0c19d7SBarry Smith       cr[i]  = r[i]+1;
120*5d0c19d7SBarry Smith       cc[i]  = c[i]+1;
121c0c2c81eSBarry Smith     }
122c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
123c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
124*5d0c19d7SBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
125*5d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
126*5d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr);
127*5d0c19d7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,cr,cc,&job);
128*5d0c19d7SBarry 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 
14754a8161fSBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1486849ba73SBarry Smith   if (sierr) {
1496849ba73SBarry Smith     switch (sierr) {
150a83599f4SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
151a83599f4SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
152e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
153e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15477431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15577431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15686bacbd2SBarry Smith     }
15786bacbd2SBarry Smith   }
15886bacbd2SBarry Smith 
15986bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16086bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   /* ------------------------------------------------------------
16386bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16486bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16586bacbd2SBarry Smith      ------------------------------------------------------------
16686bacbd2SBarry Smith   */
16786bacbd2SBarry Smith 
16887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16938baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17086bacbd2SBarry Smith 
1715bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17286bacbd2SBarry Smith 
17386bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17486bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   if (reorder) {
17786bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
180313ae024SBarry Smith   } else {
181b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
182f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
183b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
184b19713cbSBarry Smith     }
185b19713cbSBarry Smith   }
186b19713cbSBarry Smith 
187b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18886bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
189b19713cbSBarry Smith     old_i[i]--;
190b19713cbSBarry Smith   }
191b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
192b19713cbSBarry Smith     old_j[i]--;
193b19713cbSBarry Smith   }
19486bacbd2SBarry Smith 
195b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19686bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
197b19713cbSBarry Smith     new_i[i]--;
198b19713cbSBarry Smith   }
199b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
200b19713cbSBarry Smith     new_j[i]--;
201b19713cbSBarry Smith   }
20286bacbd2SBarry Smith 
20386bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20486bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2054c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2064c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
207c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20886bacbd2SBarry Smith   for(i=0; i<n; i++) {
20986bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21086bacbd2SBarry Smith   };
21186bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21207d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21386bacbd2SBarry Smith 
214c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
215c0c2c81eSBarry Smith 
216329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21707d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21886bacbd2SBarry Smith 
21986bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22086bacbd2SBarry Smith 
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;
246ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
247ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
248ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
250431e710aSBarry Smith 
2517529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25271c2f376SKris Buschelman 
25386bacbd2SBarry Smith   PetscFunctionReturn(0);
25460ec86bdSBarry Smith #endif
25586bacbd2SBarry Smith }
25686bacbd2SBarry Smith 
257e631078cSBarry Smith EXTERN_C_BEGIN
2584a2ae208SSatish Balay #undef __FUNCT__
259db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
260db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
261db4efbfdSBarry Smith {
262db4efbfdSBarry Smith   PetscFunctionBegin;
263db4efbfdSBarry Smith   *flg = PETSC_TRUE;
264db4efbfdSBarry Smith   PetscFunctionReturn(0);
265db4efbfdSBarry Smith }
266db4efbfdSBarry Smith EXTERN_C_END
267db4efbfdSBarry Smith 
268db4efbfdSBarry Smith EXTERN_C_BEGIN
269db4efbfdSBarry Smith #undef __FUNCT__
270b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc"
271b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
272b24902e0SBarry Smith {
273d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
274b24902e0SBarry Smith   PetscErrorCode     ierr;
275b24902e0SBarry Smith 
276b24902e0SBarry Smith   PetscFunctionBegin;
277b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
278b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
279b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
280b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
281db4efbfdSBarry Smith     (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
282db4efbfdSBarry Smith     (*B)->ops->ilufactorsymbolic= MatILUFactorSymbolic_SeqAIJ;
283b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
284b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
285b24902e0SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2865c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
2875c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
288b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
289db4efbfdSBarry Smith   (*B)->factor = ftype;
290b24902e0SBarry Smith   PetscFunctionReturn(0);
291b24902e0SBarry Smith }
292e631078cSBarry Smith EXTERN_C_END
293b24902e0SBarry Smith 
294b24902e0SBarry Smith #undef __FUNCT__
295b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
296719d5645SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,MatFactorInfo *info)
297289bc588SBarry Smith {
298416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
299289bc588SBarry Smith   IS                 isicol;
3006849ba73SBarry Smith   PetscErrorCode     ierr;
301*5d0c19d7SBarry Smith   const PetscInt     *r,*ic;
302*5d0c19d7SBarry Smith   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
3031393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
3041393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
3059e046878SBarry Smith   PetscReal          f;
3064875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
307a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
3081393bc97SHong Zhang   PetscBT            lnkbt;
309289bc588SBarry Smith 
3103a40ed3dSBarry Smith   PetscFunctionBegin;
311d0f46423SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
3124c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
313ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
314ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
315289bc588SBarry Smith 
316289bc588SBarry Smith   /* get new row pointers */
3171393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
3181393bc97SHong Zhang   bi[0] = 0;
3191393bc97SHong Zhang 
3201393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
3211393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
3221393bc97SHong Zhang   bdiag[0] = 0;
3231393bc97SHong Zhang 
3241393bc97SHong Zhang   /* linked list for storing column indices of the active row */
3251393bc97SHong Zhang   nlnk = n + 1;
3261393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3271393bc97SHong Zhang 
32835baf27eSSatish Balay   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
3291393bc97SHong Zhang 
3301393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
331d3d32019SBarry Smith   f = info->fill;
332a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3331393bc97SHong Zhang   current_space = free_space;
334289bc588SBarry Smith 
335289bc588SBarry Smith   for (i=0; i<n; i++) {
3361393bc97SHong Zhang     /* copy previous fill into linked list */
337289bc588SBarry Smith     nzi = 0;
3381393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3391393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3401393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
341afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3421393bc97SHong Zhang     nzi += nlnk;
3431393bc97SHong Zhang 
3441393bc97SHong Zhang     /* add pivot rows into linked list */
3451393bc97SHong Zhang     row = lnk[n];
3461393bc97SHong Zhang     while (row < i) {
347d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
348d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
349d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3501393bc97SHong Zhang       nzi += nlnk;
3511393bc97SHong Zhang       row  = lnk[row];
352289bc588SBarry Smith     }
3531393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3541393bc97SHong Zhang     im[i]   = nzi;
3551393bc97SHong Zhang 
3561393bc97SHong Zhang     /* mark bdiag */
3571393bc97SHong Zhang     nzbd = 0;
3581393bc97SHong Zhang     nnz  = nzi;
3591393bc97SHong Zhang     k    = lnk[n];
3601393bc97SHong Zhang     while (nnz-- && k < i){
3611393bc97SHong Zhang       nzbd++;
3621393bc97SHong Zhang       k = lnk[k];
3631393bc97SHong Zhang     }
3641393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3651393bc97SHong Zhang 
3661393bc97SHong Zhang     /* if free space is not available, make more free space */
3671393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3681393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
369a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3701393bc97SHong Zhang       reallocs++;
3711393bc97SHong Zhang     }
3721393bc97SHong Zhang 
3731393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3741393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3751393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3761393bc97SHong Zhang     current_space->array           += nzi;
3771393bc97SHong Zhang     current_space->local_used      += nzi;
3781393bc97SHong Zhang     current_space->local_remaining -= nzi;
379289bc588SBarry Smith   }
3806cf91177SBarry Smith #if defined(PETSC_USE_INFO)
381464e8d44SSatish Balay   if (ai[n] != 0) {
3821393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
383ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
384ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
385ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
386ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
38705bf559cSSatish Balay   } else {
388ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
38905bf559cSSatish Balay   }
39063ba0a88SBarry Smith #endif
3912fb3ed76SBarry Smith 
392898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
393898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3941987afe7SBarry Smith 
3951393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3961393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
397a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3981393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
39935baf27eSSatish Balay   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
400289bc588SBarry Smith 
401289bc588SBarry Smith   /* put together the new matrix */
402719d5645SBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
403719d5645SBarry Smith   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
404719d5645SBarry Smith   b    = (Mat_SeqAIJ*)(B)->data;
405e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
406e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
4077c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
4081393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
4091393bc97SHong Zhang   b->j          = bj;
4101393bc97SHong Zhang   b->i          = bi;
4111393bc97SHong Zhang   b->diag       = bdiag;
412416022c9SBarry Smith   b->ilen       = 0;
413416022c9SBarry Smith   b->imax       = 0;
414416022c9SBarry Smith   b->row        = isrow;
415416022c9SBarry Smith   b->col        = iscol;
416c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
417c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
41882bf6240SBarry Smith   b->icol       = isicol;
41987828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
4201393bc97SHong Zhang 
4211393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
422719d5645SBarry Smith   ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
4231393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4247fd98bd6SLois Curfman McInnes 
425719d5645SBarry Smith   (B)->factor                 =  MAT_FACTOR_LU;
426719d5645SBarry Smith   (B)->info.factor_mallocs    = reallocs;
427719d5645SBarry Smith   (B)->info.fill_ratio_given  = f;
428703deb49SSatish Balay 
429e93ccc4dSBarry Smith   if (ai[n] != 0) {
430719d5645SBarry Smith     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
431eea2913aSSatish Balay   } else {
432719d5645SBarry Smith     (B)->info.fill_ratio_needed = 0.0;
433eea2913aSSatish Balay   }
434719d5645SBarry Smith   (B)->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
435719d5645SBarry Smith   (B)->ops->solve            = MatSolve_SeqAIJ;
436719d5645SBarry Smith   (B)->ops->solvetranspose   = MatSolveTranspose_SeqAIJ;
437db4efbfdSBarry Smith   /* switch to inodes if appropriate */
438719d5645SBarry Smith   ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr);
4393a40ed3dSBarry Smith   PetscFunctionReturn(0);
440289bc588SBarry Smith }
4411393bc97SHong Zhang 
4425b5bc046SBarry Smith /*
4435b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4445b5bc046SBarry Smith */
4455b5bc046SBarry Smith #undef __FUNCT__
4465b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4475b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4485b5bc046SBarry Smith {
4495b5bc046SBarry Smith   PetscErrorCode ierr;
4505b5bc046SBarry Smith   PetscTruth     flg;
4515b5bc046SBarry Smith 
4525b5bc046SBarry Smith   PetscFunctionBegin;
4535b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4545b5bc046SBarry Smith   if (flg) {
4555b5bc046SBarry Smith     PetscViewer viewer;
4565b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4575b5bc046SBarry Smith 
4585b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4597adad957SLisandro Dalcin     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4605b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4615b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4625b5bc046SBarry Smith   }
4635b5bc046SBarry Smith   PetscFunctionReturn(0);
4645b5bc046SBarry Smith }
4655b5bc046SBarry Smith 
466db4efbfdSBarry Smith extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec);
467db4efbfdSBarry Smith 
4680a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4694a2ae208SSatish Balay #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
471719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,MatFactorInfo *info)
472289bc588SBarry Smith {
473719d5645SBarry Smith   Mat            C=B;
474416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
47582bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4766849ba73SBarry Smith   PetscErrorCode ierr;
477*5d0c19d7SBarry Smith   const PetscInt  *r,*ic,*ics;
478*5d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*bi=b->i,*bj=b->j;
479*5d0c19d7SBarry Smith   PetscInt       *ajtmp,*bjtmp,nz,row;
480d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
48154f21887SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,*rtmps;
48254f21887SBarry Smith   MatScalar      *v,*pv;
4836a7c8fc2SHong Zhang   PetscScalar    d;
4846a7c8fc2SHong Zhang   PetscReal      rs;
485b3bf805bSHong Zhang   LUShift_Ctx    sctx;
48642898933SHong Zhang   PetscInt       newshift,*ddiag;
487289bc588SBarry Smith 
4883a40ed3dSBarry Smith   PetscFunctionBegin;
48978b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
49078b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
49187828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
49287828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
493010a20caSHong Zhang   rtmps = rtmp; ics = ic;
494289bc588SBarry Smith 
495ada7143aSSatish Balay   sctx.shift_top  = 0;
496ada7143aSSatish Balay   sctx.nshift_max = 0;
497ada7143aSSatish Balay   sctx.shift_lo   = 0;
498ada7143aSSatish Balay   sctx.shift_hi   = 0;
499ada7143aSSatish Balay 
500118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
501118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
5021a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
5039deb78dcSHong Zhang     PetscInt *aai = a->i;
50442898933SHong Zhang     ddiag          = a->diag;
505118f5789SBarry Smith     sctx.shift_top = 0;
5066cc28720Svictorle     for (i=0; i<n; i++) {
5079f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
5089f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
5091a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
510010a20caSHong Zhang       v  = a->a+aai[i];
511e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
512e24cfe64SHong Zhang       for (j=0; j<nz; j++)
5131a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
5141a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
5156cc28720Svictorle     }
516c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
517118f5789SBarry Smith     sctx.shift_top    *= 1.1;
518d2276718SHong Zhang     sctx.nshift_max   = 5;
519d2276718SHong Zhang     sctx.shift_lo     = 0.;
520d2276718SHong Zhang     sctx.shift_hi     = 1.;
521a0a356efSHong Zhang   }
522a0a356efSHong Zhang 
523a0a356efSHong Zhang   sctx.shift_amount = 0;
524a0a356efSHong Zhang   sctx.nshift       = 0;
525085256b3SBarry Smith   do {
526d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
527289bc588SBarry Smith     for (i=0; i<n; i++){
528b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
529b3bf805bSHong Zhang       bjtmp = bj + bi[i];
530b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
531289bc588SBarry Smith 
532289bc588SBarry Smith       /* load in initial (unfactored row) */
533416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
534b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
535010a20caSHong Zhang       v     = a->a + a->i[r[i]];
536085256b3SBarry Smith       for (j=0; j<nz; j++) {
537b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
538085256b3SBarry Smith       }
539d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
540289bc588SBarry Smith 
541b3bf805bSHong Zhang       row = *bjtmp++;
542f2582111SSatish Balay       while  (row < i) {
5438c37ef55SBarry Smith         pc = rtmp + row;
5448c37ef55SBarry Smith         if (*pc != 0.0) {
545010a20caSHong Zhang           pv         = b->a + diag_offset[row];
546010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
54735aab85fSBarry Smith           multiplier = *pc / *pv++;
5488c37ef55SBarry Smith           *pc        = multiplier;
549b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
550f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
551efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
552289bc588SBarry Smith         }
553b3bf805bSHong Zhang         row = *bjtmp++;
554289bc588SBarry Smith       }
555416022c9SBarry Smith       /* finished row so stick it into b->a */
556b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
557b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
558b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
559b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
560d3d32019SBarry Smith       rs   = 0.0;
561d3d32019SBarry Smith       for (j=0; j<nz; j++) {
562d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
563d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
564d3d32019SBarry Smith       }
565d2276718SHong Zhang 
566b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
567d2276718SHong Zhang       sctx.rs  = rs;
568d2276718SHong Zhang       sctx.pv  = pv[diag];
569c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
5705b5bc046SBarry Smith       if (newshift == 1) break;
57135aab85fSBarry Smith     }
572d2276718SHong Zhang 
573118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5746cc28720Svictorle       /*
5756c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5766cc28720Svictorle        * then try lower shift
5776cc28720Svictorle        */
578d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
579d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
580d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
581d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
582d2276718SHong Zhang       sctx.nshift++;
5836cc28720Svictorle     }
584d2276718SHong Zhang   } while (sctx.lushift);
5850f11f9aeSSatish Balay 
58635aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
58735aab85fSBarry Smith   for (i=0; i<n; i++) {
588010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
58935aab85fSBarry Smith   }
590606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
59178b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
59278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
593db4efbfdSBarry Smith   if (b->inode.use) {
594db4efbfdSBarry Smith     C->ops->solve   = MatSolve_Inode;
595db4efbfdSBarry Smith   } else {
596db4efbfdSBarry Smith     C->ops->solve   = MatSolve_SeqAIJ;
597db4efbfdSBarry Smith   }
598db4efbfdSBarry Smith   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
599db4efbfdSBarry Smith   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
600db4efbfdSBarry Smith   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
601de26497bSBarry Smith   C->ops->matsolve           = MatMatSolve_SeqAIJ;
602c456f294SBarry Smith   C->assembled    = PETSC_TRUE;
6035c9eb25fSBarry Smith   C->preallocated = PETSC_TRUE;
604d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
605d2276718SHong Zhang   if (sctx.nshift){
606118f5789SBarry Smith     if (info->shiftnz) {
6071e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
608118f5789SBarry Smith     } else if (info->shiftpd) {
6091e2582c4SBarry 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,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
6106cc28720Svictorle     }
6110697b6caSHong Zhang   }
6123a40ed3dSBarry Smith   PetscFunctionReturn(0);
613289bc588SBarry Smith }
6140889b5dcSvictorle 
615137fb511SHong Zhang /*
616137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
617137fb511SHong Zhang    Input:
618137fb511SHong Zhang      A - original matrix
619137fb511SHong Zhang    Output;
620137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
621137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
622137fb511SHong Zhang          a->a reordered accordingly with a->j
623137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
624137fb511SHong Zhang */
625137fb511SHong Zhang #undef __FUNCT__
626137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
627719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,MatFactorInfo *info)
628137fb511SHong Zhang {
629b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
630b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
631137fb511SHong Zhang   PetscErrorCode ierr;
632*5d0c19d7SBarry Smith   const PetscInt *r,*ic,*ics;
633*5d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
634*5d0c19d7SBarry Smith   PetscInt       *ajtmp,nz,row;
635b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
636a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
637a77337e4SBarry Smith   MatScalar      *v,*pv;
638137fb511SHong Zhang   PetscReal      rs;
639137fb511SHong Zhang   LUShift_Ctx    sctx;
640b051339dSHong Zhang   PetscInt       newshift;
641137fb511SHong Zhang 
642137fb511SHong Zhang   PetscFunctionBegin;
643719d5645SBarry Smith   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
644137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
645137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
646137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
647137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
648b051339dSHong Zhang   ics = ic;
649137fb511SHong Zhang 
650137fb511SHong Zhang   sctx.shift_top  = 0;
651137fb511SHong Zhang   sctx.nshift_max = 0;
652137fb511SHong Zhang   sctx.shift_lo   = 0;
653137fb511SHong Zhang   sctx.shift_hi   = 0;
654137fb511SHong Zhang 
655137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
656137fb511SHong Zhang   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
657137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
658137fb511SHong Zhang     sctx.shift_top = 0;
659137fb511SHong Zhang     for (i=0; i<n; i++) {
660137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
661b051339dSHong Zhang       d  = (a->a)[diag[i]];
662137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
663b051339dSHong Zhang       v  = a->a+ai[i];
664b051339dSHong Zhang       nz = ai[i+1] - ai[i];
665137fb511SHong Zhang       for (j=0; j<nz; j++)
666137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
667137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
668137fb511SHong Zhang     }
669137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
670137fb511SHong Zhang     sctx.shift_top    *= 1.1;
671137fb511SHong Zhang     sctx.nshift_max   = 5;
672137fb511SHong Zhang     sctx.shift_lo     = 0.;
673137fb511SHong Zhang     sctx.shift_hi     = 1.;
674137fb511SHong Zhang   }
675137fb511SHong Zhang 
676137fb511SHong Zhang   sctx.shift_amount = 0;
677137fb511SHong Zhang   sctx.nshift       = 0;
678137fb511SHong Zhang   do {
679137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
680137fb511SHong Zhang     for (i=0; i<n; i++){
681b051339dSHong Zhang       /* load in initial unfactored row */
682b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
683b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
684b051339dSHong Zhang       v     = a->a + ai[r[i]];
685137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
686b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
687137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
688137fb511SHong Zhang 
689b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
690137fb511SHong Zhang       for (j=0; j<nz; j++) {
691137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
692b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
693137fb511SHong Zhang       }
694137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
695137fb511SHong Zhang 
696b051339dSHong Zhang       row = *ajtmp++;
697137fb511SHong Zhang       while  (row < i) {
698137fb511SHong Zhang         pc = rtmp + row;
699137fb511SHong Zhang         if (*pc != 0.0) {
700b051339dSHong Zhang           pv         = a->a + diag[r[row]];
701b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
702137fb511SHong Zhang 
703137fb511SHong Zhang           multiplier = *pc / *pv++;
704137fb511SHong Zhang           *pc        = multiplier;
705b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
706b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
707137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
708137fb511SHong Zhang         }
709b051339dSHong Zhang         row = *ajtmp++;
710137fb511SHong Zhang       }
711b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
712b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
713b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
714b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
715b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
716137fb511SHong Zhang 
717137fb511SHong Zhang       rs   = 0.0;
718137fb511SHong Zhang       for (j=0; j<nz; j++) {
719b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
720b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
721137fb511SHong Zhang       }
722137fb511SHong Zhang 
723137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
724137fb511SHong Zhang       sctx.rs  = rs;
725b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
726c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
727137fb511SHong Zhang       if (newshift == 1) break;
728137fb511SHong Zhang     }
729137fb511SHong Zhang 
730137fb511SHong Zhang     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
731137fb511SHong Zhang       /*
732137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
733137fb511SHong Zhang        * then try lower shift
734137fb511SHong Zhang        */
735137fb511SHong Zhang       sctx.shift_hi        = info->shift_fraction;
736137fb511SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
737137fb511SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
738137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
739137fb511SHong Zhang       sctx.nshift++;
740137fb511SHong Zhang     }
741137fb511SHong Zhang   } while (sctx.lushift);
742137fb511SHong Zhang 
743137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
744137fb511SHong Zhang   for (i=0; i<n; i++) {
745b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
746137fb511SHong Zhang   }
747137fb511SHong Zhang 
748137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
749137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
750137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
751b051339dSHong Zhang   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
752db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
753db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
754db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
755b051339dSHong Zhang   A->assembled = PETSC_TRUE;
7565c9eb25fSBarry Smith   A->preallocated = PETSC_TRUE;
757d0f46423SBarry Smith   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
758137fb511SHong Zhang   if (sctx.nshift){
759137fb511SHong Zhang     if (info->shiftnz) {
7601e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
761137fb511SHong Zhang     } else if (info->shiftpd) {
7621e2582c4SBarry 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,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
763137fb511SHong Zhang     }
764137fb511SHong Zhang   }
765137fb511SHong Zhang   PetscFunctionReturn(0);
766137fb511SHong Zhang }
767137fb511SHong Zhang 
7680a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7694a2ae208SSatish Balay #undef __FUNCT__
7704a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
771dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
772da3a660dSBarry Smith {
773dfbe8321SBarry Smith   PetscErrorCode ierr;
774416022c9SBarry Smith   Mat            C;
775416022c9SBarry Smith 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
77743244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
778719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
779719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
780db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
781db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
782273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
78352e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
7843a40ed3dSBarry Smith   PetscFunctionReturn(0);
785da3a660dSBarry Smith }
7860a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7874a2ae208SSatish Balay #undef __FUNCT__
7884a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
789dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
7908c37ef55SBarry Smith {
791416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
792416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
7936849ba73SBarry Smith   PetscErrorCode    ierr;
794*5d0c19d7SBarry Smith   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
795*5d0c19d7SBarry Smith   PetscInt          nz;
796*5d0c19d7SBarry Smith   const PetscInt    *rout,*cout,*r,*c;
797dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
798d9fead3dSBarry Smith   const PetscScalar *b;
799dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
8008c37ef55SBarry Smith 
8013a40ed3dSBarry Smith   PetscFunctionBegin;
8023a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
80391bf9adeSBarry Smith 
804d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8051ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
806416022c9SBarry Smith   tmp  = a->solve_work;
8078c37ef55SBarry Smith 
8083b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8093b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
8108c37ef55SBarry Smith 
8118c37ef55SBarry Smith   /* forward solve the lower triangular */
8128c37ef55SBarry Smith   tmp[0] = b[*r++];
813010a20caSHong Zhang   tmps   = tmp;
8148c37ef55SBarry Smith   for (i=1; i<n; i++) {
815010a20caSHong Zhang     v   = aa + ai[i] ;
816010a20caSHong Zhang     vi  = aj + ai[i] ;
81753e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
81853e7425aSBarry Smith     sum = b[*r++];
819ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8208c37ef55SBarry Smith     tmp[i] = sum;
8218c37ef55SBarry Smith   }
8228c37ef55SBarry Smith 
8238c37ef55SBarry Smith   /* backward solve the upper triangular */
8248c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
825010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
826010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
827416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
8288c37ef55SBarry Smith     sum = tmp[i];
829ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
830010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
8318c37ef55SBarry Smith   }
8328c37ef55SBarry Smith 
8333b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8343b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
835d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8361ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
837d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
8383a40ed3dSBarry Smith   PetscFunctionReturn(0);
8398c37ef55SBarry Smith }
840026e39d0SSatish Balay 
8412bebee5dSHong Zhang #undef __FUNCT__
8422bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
8432bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8442bebee5dSHong Zhang {
8452bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8462bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8472bebee5dSHong Zhang   PetscErrorCode  ierr;
848*5d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
849*5d0c19d7SBarry Smith   PetscInt        nz,neq;
850*5d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
851dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
852dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
853db4efbfdSBarry Smith   PetscTruth      bisdense,xisdense;
8542bebee5dSHong Zhang 
8552bebee5dSHong Zhang   PetscFunctionBegin;
8562bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8572bebee5dSHong Zhang 
858db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
859db4efbfdSBarry Smith   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
860db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
861db4efbfdSBarry Smith   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
862db4efbfdSBarry Smith 
8632bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8642bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8652bebee5dSHong Zhang 
8662bebee5dSHong Zhang   tmp  = a->solve_work;
8672bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8682bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8692bebee5dSHong Zhang 
870a36b22ccSSatish Balay   for (neq=0; neq<B->cmap->n; neq++){
8712bebee5dSHong Zhang     /* forward solve the lower triangular */
8722bebee5dSHong Zhang     tmp[0] = b[r[0]];
8732bebee5dSHong Zhang     tmps   = tmp;
8742bebee5dSHong Zhang     for (i=1; i<n; i++) {
8752bebee5dSHong Zhang       v   = aa + ai[i] ;
8762bebee5dSHong Zhang       vi  = aj + ai[i] ;
8772bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8782bebee5dSHong Zhang       sum = b[r[i]];
8792bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8802bebee5dSHong Zhang       tmp[i] = sum;
8812bebee5dSHong Zhang     }
8822bebee5dSHong Zhang     /* backward solve the upper triangular */
8832bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
8842bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
8852bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
8862bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
8872bebee5dSHong Zhang       sum = tmp[i];
8882bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8892bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
8902bebee5dSHong Zhang     }
8912bebee5dSHong Zhang 
8922bebee5dSHong Zhang     b += n;
8932bebee5dSHong Zhang     x += n;
8942bebee5dSHong Zhang   }
8952bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8962bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8972bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
8982bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
8992bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
9002bebee5dSHong Zhang   PetscFunctionReturn(0);
9012bebee5dSHong Zhang }
9022bebee5dSHong Zhang 
903137fb511SHong Zhang #undef __FUNCT__
904137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
905137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
906137fb511SHong Zhang {
907137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
908137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
909137fb511SHong Zhang   PetscErrorCode  ierr;
910*5d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
911*5d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
912*5d0c19d7SBarry Smith   PetscInt        nz,row;
913dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
914dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
915137fb511SHong Zhang 
916137fb511SHong Zhang   PetscFunctionBegin;
917137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
918137fb511SHong Zhang 
919137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
920137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
921137fb511SHong Zhang   tmp  = a->solve_work;
922137fb511SHong Zhang 
923137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
924137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
925137fb511SHong Zhang 
926137fb511SHong Zhang   /* forward solve the lower triangular */
927137fb511SHong Zhang   tmp[0] = b[*r++];
928137fb511SHong Zhang   tmps   = tmp;
929137fb511SHong Zhang   for (row=1; row<n; row++) {
930137fb511SHong Zhang     i   = rout[row]; /* permuted row */
931137fb511SHong Zhang     v   = aa + ai[i] ;
932137fb511SHong Zhang     vi  = aj + ai[i] ;
933137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
934137fb511SHong Zhang     sum = b[*r++];
935137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
936137fb511SHong Zhang     tmp[row] = sum;
937137fb511SHong Zhang   }
938137fb511SHong Zhang 
939137fb511SHong Zhang   /* backward solve the upper triangular */
940137fb511SHong Zhang   for (row=n-1; row>=0; row--){
941137fb511SHong Zhang     i   = rout[row]; /* permuted row */
942137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
943137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
944137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
945137fb511SHong Zhang     sum = tmp[row];
946137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
947137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
948137fb511SHong Zhang   }
949137fb511SHong Zhang 
950137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
951137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
952137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
953137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
954d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
955137fb511SHong Zhang   PetscFunctionReturn(0);
956137fb511SHong Zhang }
957137fb511SHong Zhang 
958930ae53cSBarry Smith /* ----------------------------------------------------------- */
9594a2ae208SSatish Balay #undef __FUNCT__
9604a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
961dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
962930ae53cSBarry Smith {
963930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9646849ba73SBarry Smith   PetscErrorCode    ierr;
965d0f46423SBarry Smith   PetscInt          n = A->rmap->n;
966356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
967356650c2SBarry Smith   PetscScalar       *x;
968356650c2SBarry Smith   const PetscScalar *b;
969dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
970aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
971356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
972dd6ea824SBarry Smith   const MatScalar   *v;
973dd6ea824SBarry Smith   PetscScalar       sum;
974d85afc42SSatish Balay #endif
975930ae53cSBarry Smith 
9763a40ed3dSBarry Smith   PetscFunctionBegin;
9773a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
978930ae53cSBarry Smith 
979356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9801ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
981930ae53cSBarry Smith 
982aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
9831c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
9841c4feecaSSatish Balay #else
985930ae53cSBarry Smith   /* forward solve the lower triangular */
986930ae53cSBarry Smith   x[0] = b[0];
987930ae53cSBarry Smith   for (i=1; i<n; i++) {
988930ae53cSBarry Smith     ai_i = ai[i];
989930ae53cSBarry Smith     v    = aa + ai_i;
990930ae53cSBarry Smith     vi   = aj + ai_i;
991930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
992930ae53cSBarry Smith     sum  = b[i];
993930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
994930ae53cSBarry Smith     x[i] = sum;
995930ae53cSBarry Smith   }
996930ae53cSBarry Smith 
997930ae53cSBarry Smith   /* backward solve the upper triangular */
998930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
999930ae53cSBarry Smith     adiag_i = adiag[i];
1000d708749eSBarry Smith     v       = aa + adiag_i + 1;
1001d708749eSBarry Smith     vi      = aj + adiag_i + 1;
1002930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
1003930ae53cSBarry Smith     sum     = x[i];
1004930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1005930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
1006930ae53cSBarry Smith   }
10071c4feecaSSatish Balay #endif
1008d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
1009356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
10101ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10113a40ed3dSBarry Smith   PetscFunctionReturn(0);
1012930ae53cSBarry Smith }
1013930ae53cSBarry Smith 
10144a2ae208SSatish Balay #undef __FUNCT__
10154a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
1016dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1017da3a660dSBarry Smith {
1018416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1019416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
10206849ba73SBarry Smith   PetscErrorCode  ierr;
1021*5d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1022*5d0c19d7SBarry Smith   PetscInt        nz;
1023*5d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
1024dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
1025dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1026da3a660dSBarry Smith 
10273a40ed3dSBarry Smith   PetscFunctionBegin;
102878b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
1029da3a660dSBarry Smith 
10301ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10311ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1032416022c9SBarry Smith   tmp  = a->solve_work;
1033da3a660dSBarry Smith 
10343b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10353b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1036da3a660dSBarry Smith 
1037da3a660dSBarry Smith   /* forward solve the lower triangular */
1038da3a660dSBarry Smith   tmp[0] = b[*r++];
1039da3a660dSBarry Smith   for (i=1; i<n; i++) {
1040010a20caSHong Zhang     v   = aa + ai[i] ;
1041010a20caSHong Zhang     vi  = aj + ai[i] ;
1042416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
1043da3a660dSBarry Smith     sum = b[*r++];
1044010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1045da3a660dSBarry Smith     tmp[i] = sum;
1046da3a660dSBarry Smith   }
1047da3a660dSBarry Smith 
1048da3a660dSBarry Smith   /* backward solve the upper triangular */
1049da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1050010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
1051010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
1052416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1053da3a660dSBarry Smith     sum = tmp[i];
1054010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1055010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1056da3a660dSBarry Smith     x[*c--] += tmp[i];
1057da3a660dSBarry Smith   }
1058da3a660dSBarry Smith 
10593b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10603b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10611ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1063efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1064898183ecSLois Curfman McInnes 
10653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1066da3a660dSBarry Smith }
1067da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10684a2ae208SSatish Balay #undef __FUNCT__
10694a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1070dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1071da3a660dSBarry Smith {
1072416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10732235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10746849ba73SBarry Smith   PetscErrorCode  ierr;
1075*5d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
1076*5d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1077*5d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1078dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1079dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1080da3a660dSBarry Smith 
10813a40ed3dSBarry Smith   PetscFunctionBegin;
10821ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10831ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1084416022c9SBarry Smith   tmp  = a->solve_work;
1085da3a660dSBarry Smith 
10862235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10872235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1088da3a660dSBarry Smith 
1089da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
10902235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1091da3a660dSBarry Smith 
1092da3a660dSBarry Smith   /* forward solve the U^T */
1093da3a660dSBarry Smith   for (i=0; i<n; i++) {
1094010a20caSHong Zhang     v   = aa + diag[i] ;
1095010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1096f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1097c38d4ed2SBarry Smith     s1  = tmp[i];
1098ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1099da3a660dSBarry Smith     while (nz--) {
1100010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1101da3a660dSBarry Smith     }
1102c38d4ed2SBarry Smith     tmp[i] = s1;
1103da3a660dSBarry Smith   }
1104da3a660dSBarry Smith 
1105da3a660dSBarry Smith   /* backward solve the L^T */
1106da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1107010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1108010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1109f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1110f1af5d2fSBarry Smith     s1  = tmp[i];
1111da3a660dSBarry Smith     while (nz--) {
1112010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1113da3a660dSBarry Smith     }
1114da3a660dSBarry Smith   }
1115da3a660dSBarry Smith 
1116da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1117da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1118da3a660dSBarry Smith 
11192235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11202235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11211ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11221ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1123da3a660dSBarry Smith 
1124d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap->n);CHKERRQ(ierr);
11253a40ed3dSBarry Smith   PetscFunctionReturn(0);
1126da3a660dSBarry Smith }
1127da3a660dSBarry Smith 
11284a2ae208SSatish Balay #undef __FUNCT__
11294a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1130dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1131da3a660dSBarry Smith {
1132416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
11332235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
11346849ba73SBarry Smith   PetscErrorCode  ierr;
1135*5d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
1136*5d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1137*5d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1138dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1139dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
11406abc6512SBarry Smith 
11413a40ed3dSBarry Smith   PetscFunctionBegin;
114223bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
11436abc6512SBarry Smith 
11441ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1146416022c9SBarry Smith   tmp = a->solve_work;
11476abc6512SBarry Smith 
11482235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11492235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11506abc6512SBarry Smith 
11516abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
11522235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
11536abc6512SBarry Smith 
11546abc6512SBarry Smith   /* forward solve the U^T */
11556abc6512SBarry Smith   for (i=0; i<n; i++) {
1156010a20caSHong Zhang     v   = aa + diag[i] ;
1157010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1158f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11596abc6512SBarry Smith     tmp[i] *= *v++;
11606abc6512SBarry Smith     while (nz--) {
1161010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11626abc6512SBarry Smith     }
11636abc6512SBarry Smith   }
11646abc6512SBarry Smith 
11656abc6512SBarry Smith   /* backward solve the L^T */
11666abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1167010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1168010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1169f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11706abc6512SBarry Smith     while (nz--) {
1171010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11726abc6512SBarry Smith     }
11736abc6512SBarry Smith   }
11746abc6512SBarry Smith 
11756abc6512SBarry Smith   /* copy tmp into x according to permutation */
11766abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11776abc6512SBarry Smith 
11782235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11792235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11801ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11811ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11826abc6512SBarry Smith 
1183efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
11843a40ed3dSBarry Smith   PetscFunctionReturn(0);
1185da3a660dSBarry Smith }
11869e25ed09SBarry Smith /* ----------------------------------------------------------------*/
11873c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1188719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption);
118908480c60SBarry Smith 
11904a2ae208SSatish Balay #undef __FUNCT__
11914a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
1192719d5645SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,MatFactorInfo *info)
11939e25ed09SBarry Smith {
1194416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
11959e25ed09SBarry Smith   IS                 isicol;
11966849ba73SBarry Smith   PetscErrorCode     ierr;
1197*5d0c19d7SBarry Smith   const PetscInt     *r,*ic;
1198*5d0c19d7SBarry Smith   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
11995a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
12005a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
12015a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
120277c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1203329f5518SBarry Smith   PetscReal          f;
12045a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12055a8e39fbSHong Zhang   PetscBT            lnkbt;
12065a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1207a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1208a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
120909f38230SBarry Smith   PetscTruth         missing;
12109e25ed09SBarry Smith 
12113a40ed3dSBarry Smith   PetscFunctionBegin;
1212d0f46423SBarry 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);
12136cf997b0SBarry Smith   f             = info->fill;
121438baddfdSBarry Smith   levels        = (PetscInt)info->levels;
121538baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
12164c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
121782bf6240SBarry Smith 
121808480c60SBarry Smith   /* special case that simply copies fill pattern */
1219be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1220be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
122186bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
1222719d5645SBarry Smith     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr);
1223719d5645SBarry Smith     fact->factor = MAT_FACTOR_ILU;
1224719d5645SBarry Smith     (fact)->info.factor_mallocs    = 0;
1225719d5645SBarry Smith     (fact)->info.fill_ratio_given  = info->fill;
1226719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 1.0;
1227719d5645SBarry Smith     b               = (Mat_SeqAIJ*)(fact)->data;
12288ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
122909f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
123008480c60SBarry Smith     b->row              = isrow;
123108480c60SBarry Smith     b->col              = iscol;
123282bf6240SBarry Smith     b->icol             = isicol;
1233719d5645SBarry Smith     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1234c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1235c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1236719d5645SBarry Smith     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1237719d5645SBarry Smith     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
12383a40ed3dSBarry Smith     PetscFunctionReturn(0);
123908480c60SBarry Smith   }
124008480c60SBarry Smith 
1241898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1242898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
12439e25ed09SBarry Smith 
12449e25ed09SBarry Smith   /* get new row pointers */
12455a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12465a8e39fbSHong Zhang   bi[0] = 0;
12475a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
12485a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12495a8e39fbSHong Zhang   bdiag[0]  = 0;
12506cf997b0SBarry Smith 
12515a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
12525a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
12535a8e39fbSHong Zhang 
12545a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
12555a8e39fbSHong Zhang   nlnk = n + 1;
12565a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12575a8e39fbSHong Zhang 
12585a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1259a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12605a8e39fbSHong Zhang   current_space = free_space;
1261a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12625a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12635a8e39fbSHong Zhang 
12645a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12655a8e39fbSHong Zhang     nzi = 0;
12666cf997b0SBarry Smith     /* copy current row into linked list */
12675a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
12685a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
12695a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12705a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12715a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12725a8e39fbSHong Zhang     nzi += nlnk;
12736cf997b0SBarry Smith 
12746cf997b0SBarry Smith     /* make sure diagonal entry is included */
12755a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12766cf997b0SBarry Smith       fm = n;
12775a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12785a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12795a8e39fbSHong Zhang       lnk[fm]    = i;
12805a8e39fbSHong Zhang       lnk_lvl[i] = 0;
12815a8e39fbSHong Zhang       nzi++; dcount++;
12826cf997b0SBarry Smith     }
12836cf997b0SBarry Smith 
12845a8e39fbSHong Zhang     /* add pivot rows into the active row */
12855a8e39fbSHong Zhang     nzbd = 0;
12865a8e39fbSHong Zhang     prow = lnk[n];
12875a8e39fbSHong Zhang     while (prow < i) {
12885a8e39fbSHong Zhang       nnz      = bdiag[prow];
12895a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
12905a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
12915a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
12925a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
12935a8e39fbSHong Zhang       nzi += nlnk;
12945a8e39fbSHong Zhang       prow = lnk[prow];
12955a8e39fbSHong Zhang       nzbd++;
129656beaf04SBarry Smith     }
12975a8e39fbSHong Zhang     bdiag[i] = nzbd;
12985a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1299ecf371e4SBarry Smith 
13005a8e39fbSHong Zhang     /* if free space is not available, make more free space */
13015a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
13025a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1303a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1304a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
13055a8e39fbSHong Zhang       reallocs++;
13065a8e39fbSHong Zhang     }
1307ecf371e4SBarry Smith 
13085a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
13095a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13105a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
13115a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
13125a8e39fbSHong Zhang 
13135a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
13145a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
131577431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1316d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
13176cf997b0SBarry Smith     }
13185a8e39fbSHong Zhang 
13195a8e39fbSHong Zhang     current_space->array           += nzi;
13205a8e39fbSHong Zhang     current_space->local_used      += nzi;
13215a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
13225a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
13235a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
13245a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
13259e25ed09SBarry Smith   }
13265a8e39fbSHong Zhang 
1327898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1328898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
13295a8e39fbSHong Zhang 
13305a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
13315a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1332a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13335a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1334a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13355a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
13369e25ed09SBarry Smith 
13376cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1338f64065bbSBarry Smith   {
13395a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1340ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1341ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1342ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1343ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1344335d9088SBarry Smith     if (diagonal_fill) {
1345ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1346335d9088SBarry Smith     }
1347f64065bbSBarry Smith   }
134863ba0a88SBarry Smith #endif
1349bbb6d6a8SBarry Smith 
13509e25ed09SBarry Smith   /* put together the new matrix */
1351719d5645SBarry Smith   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1352719d5645SBarry Smith   b = (Mat_SeqAIJ*)(fact)->data;
1353e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1354e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13557c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13565a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1357b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13585a8e39fbSHong Zhang   b->j          = bj;
13595a8e39fbSHong Zhang   b->i          = bi;
13605a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13615a8e39fbSHong Zhang   b->diag       = bdiag;
1362416022c9SBarry Smith   b->ilen       = 0;
1363416022c9SBarry Smith   b->imax       = 0;
1364416022c9SBarry Smith   b->row        = isrow;
1365416022c9SBarry Smith   b->col        = iscol;
1366c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1367c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
136882bf6240SBarry Smith   b->icol       = isicol;
136987828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1370416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13715a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
1372719d5645SBarry Smith   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13735a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
1374719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1375719d5645SBarry Smith   (fact)->info.fill_ratio_given  = f;
1376719d5645SBarry Smith   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1377719d5645SBarry Smith   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1378719d5645SBarry Smith   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
13793a40ed3dSBarry Smith   PetscFunctionReturn(0);
13809e25ed09SBarry Smith }
1381d5d45c9bSBarry Smith 
13823c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1383a6175056SHong Zhang #undef __FUNCT__
1384a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1385719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,MatFactorInfo *info)
1386a6175056SHong Zhang {
1387719d5645SBarry Smith   Mat            C = B;
13880968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
13892ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
13909bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
13912ed133dbSHong Zhang   PetscErrorCode ierr;
1392*5d0c19d7SBarry Smith   const PetscInt *rip,*riip;
1393*5d0c19d7SBarry Smith   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
13942ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
13952ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1396622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1397540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1398fbf22428SSatish Balay   PetscReal      shiftpd;
1399540703acSHong Zhang   ChShift_Ctx    sctx;
1400540703acSHong Zhang   PetscInt       newshift;
1401db4efbfdSBarry Smith   PetscTruth     perm_identity;
1402d5d45c9bSBarry Smith 
1403a6175056SHong Zhang   PetscFunctionBegin;
1404117f1891Stmunson 
1405540703acSHong Zhang   shiftnz   = info->shiftnz;
1406540703acSHong Zhang   shiftpd   = info->shiftpd;
1407ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1408ee45ca4aSHong Zhang 
14092ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
14109bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
14112ed133dbSHong Zhang 
14122ed133dbSHong Zhang   /* initialization */
14132ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
14142ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
14152ed133dbSHong Zhang   jl   = il + mbs;
14162ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
14172ed133dbSHong Zhang 
1418540703acSHong Zhang   sctx.shift_amount = 0;
1419540703acSHong Zhang   sctx.nshift       = 0;
14202ed133dbSHong Zhang   do {
1421540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
14222ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14232ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1424a6175056SHong Zhang     }
14252ed133dbSHong Zhang 
14262ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1427c05c3958SHong Zhang       bval = ba + bi[k];
14282ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
14292ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
14302ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
14319bfd6278SHong Zhang         col = riip[aj[j]];
14322ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
14332ed133dbSHong Zhang           rtmp[col] = aa[j];
1434c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
14352ed133dbSHong Zhang         }
14362ed133dbSHong Zhang       }
143739e3d630SHong Zhang       /* shift the diagonal of the matrix */
1438540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
14392ed133dbSHong Zhang 
14402ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
14412ed133dbSHong Zhang       dk = rtmp[k];
14422ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
14432ed133dbSHong Zhang 
14442ed133dbSHong Zhang       while (i < k){
14452ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
14462ed133dbSHong Zhang 
14472ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
14482ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
14492ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
14502ed133dbSHong Zhang         dk += uikdi*ba[ili];
14512ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14522ed133dbSHong Zhang 
14532ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14542ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14552ed133dbSHong Zhang         if (jmin < jmax){
14562ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14572ed133dbSHong Zhang           /* update il and jl for row i */
14582ed133dbSHong Zhang           il[i] = jmin;
14592ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14602ed133dbSHong Zhang         }
14612ed133dbSHong Zhang         i = nexti;
14622ed133dbSHong Zhang       }
14632ed133dbSHong Zhang 
14640697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14650697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14660697b6caSHong Zhang       rs   = 0.0;
14673c9e8598SHong Zhang       jmin = bi[k]+1;
14683c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14693c9e8598SHong Zhang       bcol = bj + jmin;
14703c9e8598SHong Zhang       while (nz--){
147189f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
147289f845c8SHong Zhang         bcol++;
14733c9e8598SHong Zhang       }
1474540703acSHong Zhang 
1475540703acSHong Zhang       sctx.rs = rs;
1476540703acSHong Zhang       sctx.pv = dk;
14775b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1478117f1891Stmunson 
1479117f1891Stmunson       if (newshift == 1) {
1480117f1891Stmunson         if (!sctx.shift_amount) {
1481117f1891Stmunson           sctx.shift_amount = 1e-5;
1482117f1891Stmunson         }
1483117f1891Stmunson         break;
1484117f1891Stmunson       }
14853c9e8598SHong Zhang 
14863c9e8598SHong Zhang       /* copy data into U(k,:) */
148739e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
148839e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
148939e3d630SHong Zhang       if (jmin < jmax) {
149039e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
149139e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
14923c9e8598SHong Zhang         }
149339e3d630SHong Zhang         /* add the k-th row into il and jl */
14943c9e8598SHong Zhang         il[k] = jmin;
14953c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
14963c9e8598SHong Zhang       }
14973c9e8598SHong Zhang     }
1498540703acSHong Zhang   } while (sctx.chshift);
14993c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
15003c9e8598SHong Zhang 
150139e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
15029bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1503db4efbfdSBarry Smith 
1504db4efbfdSBarry Smith   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1505db4efbfdSBarry Smith   if (perm_identity){
1506719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1507719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1508719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1509719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1510db4efbfdSBarry Smith   } else {
1511719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1512719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1513719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1514719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1515db4efbfdSBarry Smith   }
1516db4efbfdSBarry Smith 
15173c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
15183c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1519d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1520540703acSHong Zhang   if (sctx.nshift){
1521540703acSHong Zhang     if (shiftnz) {
15221e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1523540703acSHong Zhang     } else if (shiftpd) {
15241e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
15250697b6caSHong Zhang     }
15263c9e8598SHong Zhang   }
15273c9e8598SHong Zhang   PetscFunctionReturn(0);
15283c9e8598SHong Zhang }
1529a6175056SHong Zhang 
1530a6175056SHong Zhang #undef __FUNCT__
1531a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1532719d5645SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info)
1533a6175056SHong Zhang {
15340968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1535ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1536dfbe8321SBarry Smith   PetscErrorCode     ierr;
153758ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
1538*5d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
1539*5d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
1540ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
154158ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
15425a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1543ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1544a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1545a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
15467a48dd6fSHong Zhang   PetscBT            lnkbt;
1547b635d36bSHong Zhang   IS                 iperm;
1548a6175056SHong Zhang 
1549a6175056SHong Zhang   PetscFunctionBegin;
1550d0f46423SBarry 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);
155158ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
155258ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1553653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1554b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
15557a48dd6fSHong Zhang 
155639e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
155739e3d630SHong Zhang   ui[0] = 0;
155839e3d630SHong Zhang 
1559b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1560622e2028SHong Zhang   if (!levels && perm_identity) {
156158ebbce7SBarry Smith 
1562ed59401aSHong Zhang     for (i=0; i<am; i++) {
156339e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1564ed59401aSHong Zhang     }
156539e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
156639e3d630SHong Zhang     cols = uj;
1567ed59401aSHong Zhang     for (i=0; i<am; i++) {
1568ed59401aSHong Zhang       aj    = a->j + a->diag[i];
156939e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
157039e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1571ed59401aSHong Zhang     }
157239e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1573b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1574b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1575b635d36bSHong Zhang 
15767a48dd6fSHong Zhang     /* initialization */
15775a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15787a48dd6fSHong Zhang 
15797a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
15807a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
15811393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
15827a48dd6fSHong Zhang     il         = jl + am;
15837a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
15847a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
15857a48dd6fSHong Zhang     for (i=0; i<am; i++){
15867a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
15877a48dd6fSHong Zhang     }
15887a48dd6fSHong Zhang 
15897a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
15907a48dd6fSHong Zhang     nlnk = am + 1;
15912ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15927a48dd6fSHong Zhang 
15937a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1594a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
15957a48dd6fSHong Zhang     current_space = free_space;
1596a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
15977a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
15987a48dd6fSHong Zhang 
15997a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
16007a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
16017a48dd6fSHong Zhang       nzk   = 0;
1602622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1603b635d36bSHong Zhang       if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
1604622e2028SHong Zhang       ncols_upper = 0;
1605622e2028SHong Zhang       for (j=0; j<ncols; j++){
1606b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1607b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
16085a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1609622e2028SHong Zhang           ncols_upper++;
1610622e2028SHong Zhang         }
1611622e2028SHong Zhang       }
1612b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16137a48dd6fSHong Zhang       nzk += nlnk;
16147a48dd6fSHong Zhang 
16157a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
16167a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
16177a48dd6fSHong Zhang 
16187a48dd6fSHong Zhang       while (prow < k){
16197a48dd6fSHong Zhang         nextprow = jl[prow];
16207a48dd6fSHong Zhang 
16217a48dd6fSHong Zhang         /* merge prow into k-th row */
16227a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
16237a48dd6fSHong Zhang         jmax = ui[prow+1];
16247a48dd6fSHong Zhang         ncols = jmax-jmin;
1625ed59401aSHong Zhang         i     = jmin - ui[prow];
1626ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
16275a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
16285a8e39fbSHong Zhang         j     = *(uj - 1);
16295a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
16307a48dd6fSHong Zhang         nzk += nlnk;
16317a48dd6fSHong Zhang 
16327a48dd6fSHong Zhang         /* update il and jl for prow */
16337a48dd6fSHong Zhang         if (jmin < jmax){
16347a48dd6fSHong Zhang           il[prow] = jmin;
16357a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
16367a48dd6fSHong Zhang         }
16377a48dd6fSHong Zhang         prow = nextprow;
16387a48dd6fSHong Zhang       }
16397a48dd6fSHong Zhang 
16407a48dd6fSHong Zhang       /* if free space is not available, make more free space */
16417a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
16427a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
16437a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1644a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1645a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
16467a48dd6fSHong Zhang         reallocs++;
16477a48dd6fSHong Zhang       }
16487a48dd6fSHong Zhang 
16497a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1650b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
16512ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
16527a48dd6fSHong Zhang 
16537a48dd6fSHong Zhang       /* add the k-th row into il and jl */
165439e3d630SHong Zhang       if (nzk > 1){
16557a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
16567a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
16577a48dd6fSHong Zhang         il[k] = ui[k] + 1;
16587a48dd6fSHong Zhang       }
16597a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
16607a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
16617a48dd6fSHong Zhang 
16627a48dd6fSHong Zhang       current_space->array           += nzk;
16637a48dd6fSHong Zhang       current_space->local_used      += nzk;
16647a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16657a48dd6fSHong Zhang 
16667a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16677a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16687a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16697a48dd6fSHong Zhang 
16707a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16717a48dd6fSHong Zhang     }
16727a48dd6fSHong Zhang 
16736cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16747a48dd6fSHong Zhang     if (ai[am] != 0) {
167539e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1676ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1677ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1678ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16797a48dd6fSHong Zhang     } else {
1680ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16817a48dd6fSHong Zhang     }
168263ba0a88SBarry Smith #endif
16837a48dd6fSHong Zhang 
16847a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1685b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16867a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
16875a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
16887a48dd6fSHong Zhang 
16897a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
16907a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1691a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
16922ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1693a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
16947a48dd6fSHong Zhang 
169539e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
169639e3d630SHong Zhang 
16977a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
16987a48dd6fSHong Zhang 
1699719d5645SBarry Smith   b    = (Mat_SeqSBAIJ*)(fact)->data;
17007a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
17017a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
17027a48dd6fSHong Zhang   b->j    = uj;
17037a48dd6fSHong Zhang   b->i    = ui;
17047a48dd6fSHong Zhang   b->diag = 0;
17057a48dd6fSHong Zhang   b->ilen = 0;
17067a48dd6fSHong Zhang   b->imax = 0;
17077a48dd6fSHong Zhang   b->row  = perm;
1708b635d36bSHong Zhang   b->col  = perm;
1709b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1710b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1711b635d36bSHong Zhang   b->icol = iperm;
17127a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
17137a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1714719d5645SBarry Smith   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
17157a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1716e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1717e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
17187a48dd6fSHong Zhang 
1719719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1720719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
17217a48dd6fSHong Zhang   if (ai[am] != 0) {
1722719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
17237a48dd6fSHong Zhang   } else {
1724719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
17257a48dd6fSHong Zhang   }
1726719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1727a6175056SHong Zhang   PetscFunctionReturn(0);
1728a6175056SHong Zhang }
1729d5d45c9bSBarry Smith 
1730f76d2b81SHong Zhang #undef __FUNCT__
1731f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1732719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info)
1733f76d2b81SHong Zhang {
1734f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
173510c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
1736dfbe8321SBarry Smith   PetscErrorCode     ierr;
1737f76d2b81SHong Zhang   PetscTruth         perm_identity;
173810c27e3dSHong Zhang   PetscReal          fill = info->fill;
1739*5d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
1740*5d0c19d7SBarry Smith   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
174110c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
17422ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1743a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
174410c27e3dSHong Zhang   PetscBT            lnkbt;
17452ed133dbSHong Zhang   IS                 iperm;
1746f76d2b81SHong Zhang 
1747f76d2b81SHong Zhang   PetscFunctionBegin;
1748d0f46423SBarry 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);
17492ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1750f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
17512ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
17522ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
17539bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
175410c27e3dSHong Zhang 
175510c27e3dSHong Zhang   /* initialization */
17561393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
175710c27e3dSHong Zhang   ui[0] = 0;
175810c27e3dSHong Zhang 
175910c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17601393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17611393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17621393bc97SHong Zhang   il     = jl + am;
17631393bc97SHong Zhang   cols   = il + am;
17641393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17651393bc97SHong Zhang   for (i=0; i<am; i++){
17661393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1767f76d2b81SHong Zhang   }
176810c27e3dSHong Zhang 
176910c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17701393bc97SHong Zhang   nlnk = am + 1;
17711393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
177210c27e3dSHong Zhang 
17731393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1774a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
177510c27e3dSHong Zhang   current_space = free_space;
177610c27e3dSHong Zhang 
17771393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
177810c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
177910c27e3dSHong Zhang     nzk   = 0;
17802ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
17819bfd6278SHong Zhang     if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
17822ed133dbSHong Zhang     ncols_upper = 0;
17832ed133dbSHong Zhang     for (j=0; j<ncols; j++){
17849bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
17852ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
17862ed133dbSHong Zhang         cols[ncols_upper] = i;
17872ed133dbSHong Zhang         ncols_upper++;
17882ed133dbSHong Zhang       }
17892ed133dbSHong Zhang     }
17901393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
179110c27e3dSHong Zhang     nzk += nlnk;
179210c27e3dSHong Zhang 
179310c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
179410c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
179510c27e3dSHong Zhang 
179610c27e3dSHong Zhang     while (prow < k){
179710c27e3dSHong Zhang       nextprow = jl[prow];
179810c27e3dSHong Zhang       /* merge prow into k-th row */
17991393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
180010c27e3dSHong Zhang       jmax = ui[prow+1];
180110c27e3dSHong Zhang       ncols = jmax-jmin;
18021393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
18035a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
180410c27e3dSHong Zhang       nzk += nlnk;
180510c27e3dSHong Zhang 
180610c27e3dSHong Zhang       /* update il and jl for prow */
180710c27e3dSHong Zhang       if (jmin < jmax){
180810c27e3dSHong Zhang         il[prow] = jmin;
18092ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
181010c27e3dSHong Zhang       }
181110c27e3dSHong Zhang       prow = nextprow;
181210c27e3dSHong Zhang     }
181310c27e3dSHong Zhang 
181410c27e3dSHong Zhang     /* if free space is not available, make more free space */
181510c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
18161393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
181710c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1818a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
181910c27e3dSHong Zhang       reallocs++;
182010c27e3dSHong Zhang     }
182110c27e3dSHong Zhang 
182210c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
18231393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
182410c27e3dSHong Zhang 
182510c27e3dSHong Zhang     /* add the k-th row into il and jl */
182610c27e3dSHong Zhang     if (nzk-1 > 0){
18271393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
182810c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
182910c27e3dSHong Zhang       il[k] = ui[k] + 1;
183010c27e3dSHong Zhang     }
18312ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
183210c27e3dSHong Zhang     current_space->array           += nzk;
183310c27e3dSHong Zhang     current_space->local_used      += nzk;
183410c27e3dSHong Zhang     current_space->local_remaining -= nzk;
183510c27e3dSHong Zhang 
183610c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
183710c27e3dSHong Zhang   }
183810c27e3dSHong Zhang 
18396cf91177SBarry Smith #if defined(PETSC_USE_INFO)
18401393bc97SHong Zhang   if (ai[am] != 0) {
18411393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1842ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1843ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1844ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
184510c27e3dSHong Zhang   } else {
1846ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
184710c27e3dSHong Zhang   }
184863ba0a88SBarry Smith #endif
184910c27e3dSHong Zhang 
185010c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
18519bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
185210c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
185310c27e3dSHong Zhang 
185410c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18551393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1856a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
185710c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
185810c27e3dSHong Zhang 
185910c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
186010c27e3dSHong Zhang 
1861719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
186210c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1863e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1864e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18651393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
186610c27e3dSHong Zhang   b->j    = uj;
186710c27e3dSHong Zhang   b->i    = ui;
186810c27e3dSHong Zhang   b->diag = 0;
186910c27e3dSHong Zhang   b->ilen = 0;
187010c27e3dSHong Zhang   b->imax = 0;
187110c27e3dSHong Zhang   b->row  = perm;
18729bfd6278SHong Zhang   b->col  = perm;
18739bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18749bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18759bfd6278SHong Zhang   b->icol = iperm;
187610c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18771393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1878719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18791393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
188010c27e3dSHong Zhang 
1881719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1882719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
18831393bc97SHong Zhang   if (ai[am] != 0) {
1884719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
188510c27e3dSHong Zhang   } else {
1886719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
188710c27e3dSHong Zhang   }
1888719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1889f76d2b81SHong Zhang   PetscFunctionReturn(0);
1890f76d2b81SHong Zhang }
1891