xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 8d8c7c43fb9794c62ea67ae59b1fd03e296ce825)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*);
26e34fafa9SBarry Smith #endif
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
510481f469SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
635d0c19d7SBarry Smith   const PetscInt *c,*r,*ic;
645d0c19d7SBarry Smith   PetscInt       i,n = A->rmap->n,*cc,*cr;
6538baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6638baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6738baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6854f21887SBarry Smith   MatScalar      *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
690481f469SBarry Smith   PetscReal      af,dt,dtcount,dtcol,fill;
7086bacbd2SBarry Smith 
7186bacbd2SBarry Smith   PetscFunctionBegin;
7286bacbd2SBarry Smith 
730481f469SBarry Smith   if (info->dt == PETSC_DEFAULT)      dt      = .005; else dt = info->dt;
740481f469SBarry Smith   if (info->dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax);  else dtcount = info->dtcount;
750481f469SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   dtcol   = .01; else dtcol = info->dtcol;
760481f469SBarry Smith   if (info->fill == PETSC_DEFAULT)    fill    = ((double)(n*(info->dtcount+1)))/a->nz; else fill = info->fill;
770481f469SBarry Smith   lfill   = (PetscInt)(dtcount/2.0);
780481f469SBarry Smith   jmax    = (PetscInt)(fill*a->nz);
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith 
8186bacbd2SBarry Smith   /* ------------------------------------------------------------
8286bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8386bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8486bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8586bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8686bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8786bacbd2SBarry Smith      ------------------------------------------------------------
8886bacbd2SBarry Smith   */
8986bacbd2SBarry Smith 
9086bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9186bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9286bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9386bacbd2SBarry Smith   reorder = PetscNot(reorder);
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith 
9686bacbd2SBarry Smith   /* storage for ilu factor */
9738baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9838baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
99a77337e4SBarry Smith   ierr = PetscMalloc(jmax*sizeof(MatScalar),&new_a);CHKERRQ(ierr);
10038baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10186bacbd2SBarry Smith 
10286bacbd2SBarry Smith   /* ------------------------------------------------------------
10386bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10486bacbd2SBarry Smith      ------------------------------------------------------------
10586bacbd2SBarry Smith   */
106b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
107b19713cbSBarry Smith     old_j[i]++;
10886bacbd2SBarry Smith   }
109b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
110b19713cbSBarry Smith     old_i[i]++;
111b19713cbSBarry Smith   };
112010a20caSHong Zhang 
11386bacbd2SBarry Smith 
11486bacbd2SBarry Smith   if (reorder) {
115c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
116c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1175d0c19d7SBarry Smith     ierr = PetscMalloc2(n,PetscInt,&cc,n,PetscInt,&cr);CHKERRQ(ierr);
118c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
1195d0c19d7SBarry Smith       cr[i]  = r[i]+1;
1205d0c19d7SBarry Smith       cc[i]  = c[i]+1;
121c0c2c81eSBarry Smith     }
122c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
123c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1245d0c19d7SBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
1255d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
1265d0c19d7SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr);
1275d0c19d7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,cr,cc,&job);
1285d0c19d7SBarry Smith     ierr = PetscFree2(cc,cr);CHKERRQ(ierr);
129b19713cbSBarry Smith     o_a = old_a2;
130b19713cbSBarry Smith     o_j = old_j2;
131b19713cbSBarry Smith     o_i = old_i2;
132b19713cbSBarry Smith   } else {
133b19713cbSBarry Smith     o_a = old_a;
134b19713cbSBarry Smith     o_j = old_j;
135b19713cbSBarry Smith     o_i = old_i;
13686bacbd2SBarry Smith   }
13786bacbd2SBarry Smith 
13886bacbd2SBarry Smith   /* ------------------------------------------------------------
13986bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14086bacbd2SBarry Smith      ------------------------------------------------------------
14186bacbd2SBarry Smith   */
14286bacbd2SBarry Smith 
14338baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14587828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14686bacbd2SBarry Smith 
1470481f469SBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)dt,&dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1486849ba73SBarry Smith   if (sierr) {
1496849ba73SBarry Smith     switch (sierr) {
1500481f469SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger fill current fill %G space allocated %D",fill,jmax);
1510481f469SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger fill current fill %G space allocated %D",fill,jmax);
152e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
153e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
1540481f469SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal fill value %D",jmax);
15577431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15686bacbd2SBarry Smith     }
15786bacbd2SBarry Smith   }
15886bacbd2SBarry Smith 
15986bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16086bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   /* ------------------------------------------------------------
16386bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16486bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16586bacbd2SBarry Smith      ------------------------------------------------------------
16686bacbd2SBarry Smith   */
16786bacbd2SBarry Smith 
16887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16938baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17086bacbd2SBarry Smith 
1715bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17286bacbd2SBarry Smith 
17386bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17486bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   if (reorder) {
17786bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
180313ae024SBarry Smith   } else {
181b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
182f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
183b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
184b19713cbSBarry Smith     }
185b19713cbSBarry Smith   }
186b19713cbSBarry Smith 
187b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18886bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
189b19713cbSBarry Smith     old_i[i]--;
190b19713cbSBarry Smith   }
191b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
192b19713cbSBarry Smith     old_j[i]--;
193b19713cbSBarry Smith   }
19486bacbd2SBarry Smith 
195b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19686bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
197b19713cbSBarry Smith     new_i[i]--;
198b19713cbSBarry Smith   }
199b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
200b19713cbSBarry Smith     new_j[i]--;
201b19713cbSBarry Smith   }
20286bacbd2SBarry Smith 
20386bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20486bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2054c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2064c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
207c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20886bacbd2SBarry Smith   for(i=0; i<n; i++) {
20986bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21086bacbd2SBarry Smith   };
21186bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21207d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21386bacbd2SBarry Smith 
214c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
215c0c2c81eSBarry Smith 
216329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21707d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21886bacbd2SBarry Smith 
21986bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22086bacbd2SBarry Smith 
2217adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
222f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
2237adad957SLisandro Dalcin   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
224ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2255c9eb25fSBarry Smith   (*fact)->factor    = MAT_FACTOR_LU;
22686bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22786bacbd2SBarry Smith 
22886bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
229e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
230e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
23107d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23286bacbd2SBarry Smith   b->a             = new_a;
23386bacbd2SBarry Smith   b->j             = new_j;
23486bacbd2SBarry Smith   b->i             = new_i;
23586bacbd2SBarry Smith   b->ilen          = 0;
23686bacbd2SBarry Smith   b->imax          = 0;
2371f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
238313ae024SBarry Smith   b->row           = isirow;
23986bacbd2SBarry Smith   b->col           = iscolf;
24087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24186bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24286bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24386bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24486bacbd2SBarry Smith 
245431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
2460481f469SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",fill,af);CHKERRQ(ierr);
247ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
248ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
250431e710aSBarry Smith 
2517529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25271c2f376SKris Buschelman 
25386bacbd2SBarry Smith   PetscFunctionReturn(0);
25460ec86bdSBarry Smith #endif
25586bacbd2SBarry Smith }
25686bacbd2SBarry Smith 
257e631078cSBarry Smith EXTERN_C_BEGIN
2584a2ae208SSatish Balay #undef __FUNCT__
259db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
260db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
261db4efbfdSBarry Smith {
262db4efbfdSBarry Smith   PetscFunctionBegin;
263db4efbfdSBarry Smith   *flg = PETSC_TRUE;
264db4efbfdSBarry Smith   PetscFunctionReturn(0);
265db4efbfdSBarry Smith }
266db4efbfdSBarry Smith EXTERN_C_END
267db4efbfdSBarry Smith 
268db4efbfdSBarry Smith EXTERN_C_BEGIN
269db4efbfdSBarry Smith #undef __FUNCT__
270b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc"
271b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
272b24902e0SBarry Smith {
273d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
274b24902e0SBarry Smith   PetscErrorCode     ierr;
275b24902e0SBarry Smith 
276b24902e0SBarry Smith   PetscFunctionBegin;
277b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
278b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
279b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
280b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
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"
2960481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
297289bc588SBarry Smith {
298416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
299289bc588SBarry Smith   IS                 isicol;
3006849ba73SBarry Smith   PetscErrorCode     ierr;
3015d0c19d7SBarry Smith   const PetscInt     *r,*ic;
3025d0c19d7SBarry 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"
4710481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const 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;
4775d0c19d7SBarry Smith   const PetscInt  *r,*ic,*ics;
4785d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*bi=b->i,*bj=b->j;
4795d0c19d7SBarry 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;
4990481f469SBarry Smith   sctx.shift_fraction = 0;
500ada7143aSSatish Balay 
501118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
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 
5730481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.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        */
5780481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
5790481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
5800481f469SBarry Smith       sctx.shift_amount    = sctx.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 {
596*8d8c7c43SBarry Smith     PetscTruth row_identity, col_identity;
597*8d8c7c43SBarry Smith     ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
598*8d8c7c43SBarry Smith     ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
599*8d8c7c43SBarry Smith     if (row_identity && col_identity) {
600*8d8c7c43SBarry Smith       C->ops->solve   = MatSolve_SeqAIJ_NaturalOrdering;
601*8d8c7c43SBarry Smith     } else {
602db4efbfdSBarry Smith       C->ops->solve   = MatSolve_SeqAIJ;
603db4efbfdSBarry Smith     }
604*8d8c7c43SBarry Smith   }
605db4efbfdSBarry Smith   C->ops->solveadd           = MatSolveAdd_SeqAIJ;
606db4efbfdSBarry Smith   C->ops->solvetranspose     = MatSolveTranspose_SeqAIJ;
607db4efbfdSBarry Smith   C->ops->solvetransposeadd  = MatSolveTransposeAdd_SeqAIJ;
608de26497bSBarry Smith   C->ops->matsolve           = MatMatSolve_SeqAIJ;
609c456f294SBarry Smith   C->assembled    = PETSC_TRUE;
6105c9eb25fSBarry Smith   C->preallocated = PETSC_TRUE;
611d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
612d2276718SHong Zhang   if (sctx.nshift){
613e738e422SBarry Smith      if (info->shiftpd) {
6140481f469SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
615e738e422SBarry Smith     } else if (info->shiftnz) {
616e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
6176cc28720Svictorle     }
6180697b6caSHong Zhang   }
6193a40ed3dSBarry Smith   PetscFunctionReturn(0);
620289bc588SBarry Smith }
6210889b5dcSvictorle 
622137fb511SHong Zhang /*
623137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
624137fb511SHong Zhang    Input:
625137fb511SHong Zhang      A - original matrix
626137fb511SHong Zhang    Output;
627137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
628137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
629137fb511SHong Zhang          a->a reordered accordingly with a->j
630137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
631137fb511SHong Zhang */
632137fb511SHong Zhang #undef __FUNCT__
633137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
6340481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
635137fb511SHong Zhang {
636b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
637b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
638137fb511SHong Zhang   PetscErrorCode ierr;
6395d0c19d7SBarry Smith   const PetscInt *r,*ic,*ics;
6405d0c19d7SBarry Smith   PetscInt       i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
6415d0c19d7SBarry Smith   PetscInt       *ajtmp,nz,row;
642b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
643a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
644a77337e4SBarry Smith   MatScalar      *v,*pv;
645137fb511SHong Zhang   PetscReal      rs;
646137fb511SHong Zhang   LUShift_Ctx    sctx;
647b051339dSHong Zhang   PetscInt       newshift;
648137fb511SHong Zhang 
649137fb511SHong Zhang   PetscFunctionBegin;
650719d5645SBarry Smith   if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
651137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
652137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
653137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
654137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
655b051339dSHong Zhang   ics = ic;
656137fb511SHong Zhang 
657137fb511SHong Zhang   sctx.shift_top      = 0;
658137fb511SHong Zhang   sctx.nshift_max     = 0;
659137fb511SHong Zhang   sctx.shift_lo       = 0;
660137fb511SHong Zhang   sctx.shift_hi       = 0;
6610481f469SBarry Smith   sctx.shift_fraction = 0;
662137fb511SHong Zhang 
663137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
664137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
665137fb511SHong Zhang     sctx.shift_top = 0;
666137fb511SHong Zhang     for (i=0; i<n; i++) {
667137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
668b051339dSHong Zhang       d  = (a->a)[diag[i]];
669137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
670b051339dSHong Zhang       v  = a->a+ai[i];
671b051339dSHong Zhang       nz = ai[i+1] - ai[i];
672137fb511SHong Zhang       for (j=0; j<nz; j++)
673137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
674137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
675137fb511SHong Zhang     }
676137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
677137fb511SHong Zhang     sctx.shift_top    *= 1.1;
678137fb511SHong Zhang     sctx.nshift_max   = 5;
679137fb511SHong Zhang     sctx.shift_lo     = 0.;
680137fb511SHong Zhang     sctx.shift_hi     = 1.;
681137fb511SHong Zhang   }
682137fb511SHong Zhang 
683137fb511SHong Zhang   sctx.shift_amount = 0;
684137fb511SHong Zhang   sctx.nshift       = 0;
685137fb511SHong Zhang   do {
686137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
687137fb511SHong Zhang     for (i=0; i<n; i++){
688b051339dSHong Zhang       /* load in initial unfactored row */
689b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
690b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
691b051339dSHong Zhang       v     = a->a + ai[r[i]];
692137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
693b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
694137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
695137fb511SHong Zhang 
696b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
697137fb511SHong Zhang       for (j=0; j<nz; j++) {
698137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
699b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
700137fb511SHong Zhang       }
701137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
702137fb511SHong Zhang 
703b051339dSHong Zhang       row = *ajtmp++;
704137fb511SHong Zhang       while  (row < i) {
705137fb511SHong Zhang         pc = rtmp + row;
706137fb511SHong Zhang         if (*pc != 0.0) {
707b051339dSHong Zhang           pv         = a->a + diag[r[row]];
708b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
709137fb511SHong Zhang 
710137fb511SHong Zhang           multiplier = *pc / *pv++;
711137fb511SHong Zhang           *pc        = multiplier;
712b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
713b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
714137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
715137fb511SHong Zhang         }
716b051339dSHong Zhang         row = *ajtmp++;
717137fb511SHong Zhang       }
718b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
719b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
720b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
721b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
722b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
723137fb511SHong Zhang 
724137fb511SHong Zhang       rs   = 0.0;
725137fb511SHong Zhang       for (j=0; j<nz; j++) {
726b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
727b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
728137fb511SHong Zhang       }
729137fb511SHong Zhang 
730137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
731137fb511SHong Zhang       sctx.rs  = rs;
732b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
733c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
734137fb511SHong Zhang       if (newshift == 1) break;
735137fb511SHong Zhang     }
736137fb511SHong Zhang 
7370481f469SBarry Smith     if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
738137fb511SHong Zhang       /*
739137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
740137fb511SHong Zhang        * then try lower shift
741137fb511SHong Zhang        */
7420481f469SBarry Smith       sctx.shift_hi        = sctx.shift_fraction;
7430481f469SBarry Smith       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
7440481f469SBarry Smith       sctx.shift_amount    = sctx.shift_fraction * sctx.shift_top;
745137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
746137fb511SHong Zhang       sctx.nshift++;
747137fb511SHong Zhang     }
748137fb511SHong Zhang   } while (sctx.lushift);
749137fb511SHong Zhang 
750137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
751137fb511SHong Zhang   for (i=0; i<n; i++) {
752b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
753137fb511SHong Zhang   }
754137fb511SHong Zhang 
755137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
756137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
757137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
758b051339dSHong Zhang   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
759db4efbfdSBarry Smith   A->ops->solveadd          = MatSolveAdd_SeqAIJ;
760db4efbfdSBarry Smith   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
761db4efbfdSBarry Smith   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
762b051339dSHong Zhang   A->assembled = PETSC_TRUE;
7635c9eb25fSBarry Smith   A->preallocated = PETSC_TRUE;
764d0f46423SBarry Smith   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
765137fb511SHong Zhang   if (sctx.nshift){
766e738e422SBarry Smith     if (info->shiftpd) {
7670481f469SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
768e738e422SBarry Smith     } else if (info->shiftnz) {
769e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
770137fb511SHong Zhang     }
771137fb511SHong Zhang   }
772137fb511SHong Zhang   PetscFunctionReturn(0);
773137fb511SHong Zhang }
774137fb511SHong Zhang 
7750a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7764a2ae208SSatish Balay #undef __FUNCT__
7774a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
7780481f469SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
779da3a660dSBarry Smith {
780dfbe8321SBarry Smith   PetscErrorCode ierr;
781416022c9SBarry Smith   Mat            C;
782416022c9SBarry Smith 
7833a40ed3dSBarry Smith   PetscFunctionBegin;
78443244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
785719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
786719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
787db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
788db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
789273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
79052e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
7913a40ed3dSBarry Smith   PetscFunctionReturn(0);
792da3a660dSBarry Smith }
7930a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7944a2ae208SSatish Balay #undef __FUNCT__
7954a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
796dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
7978c37ef55SBarry Smith {
798416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
799416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
8006849ba73SBarry Smith   PetscErrorCode    ierr;
8015d0c19d7SBarry Smith   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8025d0c19d7SBarry Smith   PetscInt          nz;
8035d0c19d7SBarry Smith   const PetscInt    *rout,*cout,*r,*c;
804dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
805d9fead3dSBarry Smith   const PetscScalar *b;
806dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
8078c37ef55SBarry Smith 
8083a40ed3dSBarry Smith   PetscFunctionBegin;
8093a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
81091bf9adeSBarry Smith 
811d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8121ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
813416022c9SBarry Smith   tmp  = a->solve_work;
8148c37ef55SBarry Smith 
8153b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8163b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
8178c37ef55SBarry Smith 
8188c37ef55SBarry Smith   /* forward solve the lower triangular */
8198c37ef55SBarry Smith   tmp[0] = b[*r++];
820010a20caSHong Zhang   tmps   = tmp;
8218c37ef55SBarry Smith   for (i=1; i<n; i++) {
822010a20caSHong Zhang     v   = aa + ai[i] ;
823010a20caSHong Zhang     vi  = aj + ai[i] ;
82453e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
82553e7425aSBarry Smith     sum = b[*r++];
826ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8278c37ef55SBarry Smith     tmp[i] = sum;
8288c37ef55SBarry Smith   }
8298c37ef55SBarry Smith 
8308c37ef55SBarry Smith   /* backward solve the upper triangular */
8318c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
832010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
833010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
834416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
8358c37ef55SBarry Smith     sum = tmp[i];
836ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
837010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
8388c37ef55SBarry Smith   }
8398c37ef55SBarry Smith 
8403b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8413b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
842d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8431ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
844d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
8468c37ef55SBarry Smith }
847026e39d0SSatish Balay 
8482bebee5dSHong Zhang #undef __FUNCT__
8492bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
8502bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8512bebee5dSHong Zhang {
8522bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8532bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8542bebee5dSHong Zhang   PetscErrorCode  ierr;
8555d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8565d0c19d7SBarry Smith   PetscInt        nz,neq;
8575d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
858dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
859dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
860db4efbfdSBarry Smith   PetscTruth      bisdense,xisdense;
8612bebee5dSHong Zhang 
8622bebee5dSHong Zhang   PetscFunctionBegin;
8632bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8642bebee5dSHong Zhang 
865db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr);
866db4efbfdSBarry Smith   if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
867db4efbfdSBarry Smith   ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr);
868db4efbfdSBarry Smith   if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");
869db4efbfdSBarry Smith 
8702bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8712bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8722bebee5dSHong Zhang 
8732bebee5dSHong Zhang   tmp  = a->solve_work;
8742bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8752bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8762bebee5dSHong Zhang 
877a36b22ccSSatish Balay   for (neq=0; neq<B->cmap->n; neq++){
8782bebee5dSHong Zhang     /* forward solve the lower triangular */
8792bebee5dSHong Zhang     tmp[0] = b[r[0]];
8802bebee5dSHong Zhang     tmps   = tmp;
8812bebee5dSHong Zhang     for (i=1; i<n; i++) {
8822bebee5dSHong Zhang       v   = aa + ai[i] ;
8832bebee5dSHong Zhang       vi  = aj + ai[i] ;
8842bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8852bebee5dSHong Zhang       sum = b[r[i]];
8862bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8872bebee5dSHong Zhang       tmp[i] = sum;
8882bebee5dSHong Zhang     }
8892bebee5dSHong Zhang     /* backward solve the upper triangular */
8902bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
8912bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
8922bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
8932bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
8942bebee5dSHong Zhang       sum = tmp[i];
8952bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8962bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
8972bebee5dSHong Zhang     }
8982bebee5dSHong Zhang 
8992bebee5dSHong Zhang     b += n;
9002bebee5dSHong Zhang     x += n;
9012bebee5dSHong Zhang   }
9022bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
9032bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
9042bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
9052bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
9062bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
9072bebee5dSHong Zhang   PetscFunctionReturn(0);
9082bebee5dSHong Zhang }
9092bebee5dSHong Zhang 
910137fb511SHong Zhang #undef __FUNCT__
911137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
912137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
913137fb511SHong Zhang {
914137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
915137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
916137fb511SHong Zhang   PetscErrorCode  ierr;
9175d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
9185d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
9195d0c19d7SBarry Smith   PetscInt        nz,row;
920dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
921dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
922137fb511SHong Zhang 
923137fb511SHong Zhang   PetscFunctionBegin;
924137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
925137fb511SHong Zhang 
926137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
927137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
928137fb511SHong Zhang   tmp  = a->solve_work;
929137fb511SHong Zhang 
930137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
931137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
932137fb511SHong Zhang 
933137fb511SHong Zhang   /* forward solve the lower triangular */
934137fb511SHong Zhang   tmp[0] = b[*r++];
935137fb511SHong Zhang   tmps   = tmp;
936137fb511SHong Zhang   for (row=1; row<n; row++) {
937137fb511SHong Zhang     i   = rout[row]; /* permuted row */
938137fb511SHong Zhang     v   = aa + ai[i] ;
939137fb511SHong Zhang     vi  = aj + ai[i] ;
940137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
941137fb511SHong Zhang     sum = b[*r++];
942137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
943137fb511SHong Zhang     tmp[row] = sum;
944137fb511SHong Zhang   }
945137fb511SHong Zhang 
946137fb511SHong Zhang   /* backward solve the upper triangular */
947137fb511SHong Zhang   for (row=n-1; row>=0; row--){
948137fb511SHong Zhang     i   = rout[row]; /* permuted row */
949137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
950137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
951137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
952137fb511SHong Zhang     sum = tmp[row];
953137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
954137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
955137fb511SHong Zhang   }
956137fb511SHong Zhang 
957137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
958137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
959137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
960137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
961d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
962137fb511SHong Zhang   PetscFunctionReturn(0);
963137fb511SHong Zhang }
964137fb511SHong Zhang 
965930ae53cSBarry Smith /* ----------------------------------------------------------- */
9664a2ae208SSatish Balay #undef __FUNCT__
9674a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
968dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
969930ae53cSBarry Smith {
970930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9716849ba73SBarry Smith   PetscErrorCode    ierr;
972d0f46423SBarry Smith   PetscInt          n = A->rmap->n;
973356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
974356650c2SBarry Smith   PetscScalar       *x;
975356650c2SBarry Smith   const PetscScalar *b;
976dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
977aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
978356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
979dd6ea824SBarry Smith   const MatScalar   *v;
980dd6ea824SBarry Smith   PetscScalar       sum;
981d85afc42SSatish Balay #endif
982930ae53cSBarry Smith 
9833a40ed3dSBarry Smith   PetscFunctionBegin;
9843a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
985930ae53cSBarry Smith 
986356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
988930ae53cSBarry Smith 
989aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
9901c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
9911c4feecaSSatish Balay #else
992930ae53cSBarry Smith   /* forward solve the lower triangular */
993930ae53cSBarry Smith   x[0] = b[0];
994930ae53cSBarry Smith   for (i=1; i<n; i++) {
995930ae53cSBarry Smith     ai_i = ai[i];
996930ae53cSBarry Smith     v    = aa + ai_i;
997930ae53cSBarry Smith     vi   = aj + ai_i;
998930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
999930ae53cSBarry Smith     sum  = b[i];
1000930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1001930ae53cSBarry Smith     x[i] = sum;
1002930ae53cSBarry Smith   }
1003930ae53cSBarry Smith 
1004930ae53cSBarry Smith   /* backward solve the upper triangular */
1005930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
1006930ae53cSBarry Smith     adiag_i = adiag[i];
1007d708749eSBarry Smith     v       = aa + adiag_i + 1;
1008d708749eSBarry Smith     vi      = aj + adiag_i + 1;
1009930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
1010930ae53cSBarry Smith     sum     = x[i];
1011930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
1012930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
1013930ae53cSBarry Smith   }
10141c4feecaSSatish Balay #endif
1015d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
1016356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
10171ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10183a40ed3dSBarry Smith   PetscFunctionReturn(0);
1019930ae53cSBarry Smith }
1020930ae53cSBarry Smith 
10214a2ae208SSatish Balay #undef __FUNCT__
10224a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
1023dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1024da3a660dSBarry Smith {
1025416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
1026416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
10276849ba73SBarry Smith   PetscErrorCode  ierr;
10285d0c19d7SBarry Smith   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10295d0c19d7SBarry Smith   PetscInt        nz;
10305d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
1031dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
1032dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1033da3a660dSBarry Smith 
10343a40ed3dSBarry Smith   PetscFunctionBegin;
103578b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
1036da3a660dSBarry Smith 
10371ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1039416022c9SBarry Smith   tmp  = a->solve_work;
1040da3a660dSBarry Smith 
10413b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10423b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1043da3a660dSBarry Smith 
1044da3a660dSBarry Smith   /* forward solve the lower triangular */
1045da3a660dSBarry Smith   tmp[0] = b[*r++];
1046da3a660dSBarry Smith   for (i=1; i<n; i++) {
1047010a20caSHong Zhang     v   = aa + ai[i] ;
1048010a20caSHong Zhang     vi  = aj + ai[i] ;
1049416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
1050da3a660dSBarry Smith     sum = b[*r++];
1051010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1052da3a660dSBarry Smith     tmp[i] = sum;
1053da3a660dSBarry Smith   }
1054da3a660dSBarry Smith 
1055da3a660dSBarry Smith   /* backward solve the upper triangular */
1056da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1057010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
1058010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
1059416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1060da3a660dSBarry Smith     sum = tmp[i];
1061010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1062010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1063da3a660dSBarry Smith     x[*c--] += tmp[i];
1064da3a660dSBarry Smith   }
1065da3a660dSBarry Smith 
10663b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10673b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10681ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10691ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1070efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1071898183ecSLois Curfman McInnes 
10723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1073da3a660dSBarry Smith }
1074da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10754a2ae208SSatish Balay #undef __FUNCT__
10764a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1077dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1078da3a660dSBarry Smith {
1079416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10802235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10816849ba73SBarry Smith   PetscErrorCode  ierr;
10825d0c19d7SBarry Smith   const PetscInt  *rout,*cout,*r,*c;
10835d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
10845d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1085dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1086dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1087da3a660dSBarry Smith 
10883a40ed3dSBarry Smith   PetscFunctionBegin;
10891ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1091416022c9SBarry Smith   tmp  = a->solve_work;
1092da3a660dSBarry Smith 
10932235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10942235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1095da3a660dSBarry Smith 
1096da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
10972235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1098da3a660dSBarry Smith 
1099da3a660dSBarry Smith   /* forward solve the U^T */
1100da3a660dSBarry Smith   for (i=0; i<n; i++) {
1101010a20caSHong Zhang     v   = aa + diag[i] ;
1102010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1103f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1104c38d4ed2SBarry Smith     s1  = tmp[i];
1105ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1106da3a660dSBarry Smith     while (nz--) {
1107010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1108da3a660dSBarry Smith     }
1109c38d4ed2SBarry Smith     tmp[i] = s1;
1110da3a660dSBarry Smith   }
1111da3a660dSBarry Smith 
1112da3a660dSBarry Smith   /* backward solve the L^T */
1113da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1114010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1115010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1116f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1117f1af5d2fSBarry Smith     s1  = tmp[i];
1118da3a660dSBarry Smith     while (nz--) {
1119010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1120da3a660dSBarry Smith     }
1121da3a660dSBarry Smith   }
1122da3a660dSBarry Smith 
1123da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1124da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1125da3a660dSBarry Smith 
11262235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11272235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11281ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11291ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1130da3a660dSBarry Smith 
1131d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap->n);CHKERRQ(ierr);
11323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1133da3a660dSBarry Smith }
1134da3a660dSBarry Smith 
11354a2ae208SSatish Balay #undef __FUNCT__
11364a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1137dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1138da3a660dSBarry Smith {
1139416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
11402235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
11416849ba73SBarry Smith   PetscErrorCode  ierr;
11425d0c19d7SBarry Smith   const PetscInt  *r,*c,*rout,*cout;
11435d0c19d7SBarry Smith   PetscInt        i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
11445d0c19d7SBarry Smith   PetscInt        nz,*diag = a->diag;
1145dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1146dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
11476abc6512SBarry Smith 
11483a40ed3dSBarry Smith   PetscFunctionBegin;
114923bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
11506abc6512SBarry Smith 
11511ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1153416022c9SBarry Smith   tmp = a->solve_work;
11546abc6512SBarry Smith 
11552235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11562235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11576abc6512SBarry Smith 
11586abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
11592235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
11606abc6512SBarry Smith 
11616abc6512SBarry Smith   /* forward solve the U^T */
11626abc6512SBarry Smith   for (i=0; i<n; i++) {
1163010a20caSHong Zhang     v   = aa + diag[i] ;
1164010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1165f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11666abc6512SBarry Smith     tmp[i] *= *v++;
11676abc6512SBarry Smith     while (nz--) {
1168010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11696abc6512SBarry Smith     }
11706abc6512SBarry Smith   }
11716abc6512SBarry Smith 
11726abc6512SBarry Smith   /* backward solve the L^T */
11736abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1174010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1175010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1176f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11776abc6512SBarry Smith     while (nz--) {
1178010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11796abc6512SBarry Smith     }
11806abc6512SBarry Smith   }
11816abc6512SBarry Smith 
11826abc6512SBarry Smith   /* copy tmp into x according to permutation */
11836abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11846abc6512SBarry Smith 
11852235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11862235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11871ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11881ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11896abc6512SBarry Smith 
1190efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
11913a40ed3dSBarry Smith   PetscFunctionReturn(0);
1192da3a660dSBarry Smith }
11939e25ed09SBarry Smith /* ----------------------------------------------------------------*/
11943c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1195719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption);
119608480c60SBarry Smith 
11974a2ae208SSatish Balay #undef __FUNCT__
11984a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
11990481f469SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
12009e25ed09SBarry Smith {
1201416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
12029e25ed09SBarry Smith   IS                 isicol;
12036849ba73SBarry Smith   PetscErrorCode     ierr;
12045d0c19d7SBarry Smith   const PetscInt     *r,*ic;
12055d0c19d7SBarry Smith   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j,d;
12065a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
12075a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
12085a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
120977c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1210329f5518SBarry Smith   PetscReal          f;
12115a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12125a8e39fbSHong Zhang   PetscBT            lnkbt;
12135a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1214a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1215a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
121609f38230SBarry Smith   PetscTruth         missing;
12179e25ed09SBarry Smith 
12183a40ed3dSBarry Smith   PetscFunctionBegin;
1219d0f46423SBarry 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);
12206cf997b0SBarry Smith   f             = info->fill;
122138baddfdSBarry Smith   levels        = (PetscInt)info->levels;
122238baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
12234c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
122482bf6240SBarry Smith 
122508480c60SBarry Smith   /* special case that simply copies fill pattern */
1226be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1227be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
122886bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
1229719d5645SBarry Smith     ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr);
1230719d5645SBarry Smith     fact->factor = MAT_FACTOR_ILU;
1231719d5645SBarry Smith     (fact)->info.factor_mallocs    = 0;
1232719d5645SBarry Smith     (fact)->info.fill_ratio_given  = info->fill;
1233719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 1.0;
1234719d5645SBarry Smith     b               = (Mat_SeqAIJ*)(fact)->data;
12358ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
123609f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
123708480c60SBarry Smith     b->row              = isrow;
123808480c60SBarry Smith     b->col              = iscol;
123982bf6240SBarry Smith     b->icol             = isicol;
1240719d5645SBarry Smith     ierr                = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1241c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1242c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1243719d5645SBarry Smith     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1244719d5645SBarry Smith     ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
12453a40ed3dSBarry Smith     PetscFunctionReturn(0);
124608480c60SBarry Smith   }
124708480c60SBarry Smith 
1248898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1249898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
12509e25ed09SBarry Smith 
12519e25ed09SBarry Smith   /* get new row pointers */
12525a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12535a8e39fbSHong Zhang   bi[0] = 0;
12545a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
12555a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12565a8e39fbSHong Zhang   bdiag[0]  = 0;
12576cf997b0SBarry Smith 
12585a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
12595a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
12605a8e39fbSHong Zhang 
12615a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
12625a8e39fbSHong Zhang   nlnk = n + 1;
12635a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12645a8e39fbSHong Zhang 
12655a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1266a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12675a8e39fbSHong Zhang   current_space = free_space;
1268a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12695a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12705a8e39fbSHong Zhang 
12715a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12725a8e39fbSHong Zhang     nzi = 0;
12736cf997b0SBarry Smith     /* copy current row into linked list */
12745a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
12755a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
12765a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12775a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12785a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12795a8e39fbSHong Zhang     nzi += nlnk;
12806cf997b0SBarry Smith 
12816cf997b0SBarry Smith     /* make sure diagonal entry is included */
12825a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12836cf997b0SBarry Smith       fm = n;
12845a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12855a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12865a8e39fbSHong Zhang       lnk[fm]    = i;
12875a8e39fbSHong Zhang       lnk_lvl[i] = 0;
12885a8e39fbSHong Zhang       nzi++; dcount++;
12896cf997b0SBarry Smith     }
12906cf997b0SBarry Smith 
12915a8e39fbSHong Zhang     /* add pivot rows into the active row */
12925a8e39fbSHong Zhang     nzbd = 0;
12935a8e39fbSHong Zhang     prow = lnk[n];
12945a8e39fbSHong Zhang     while (prow < i) {
12955a8e39fbSHong Zhang       nnz      = bdiag[prow];
12965a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
12975a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
12985a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
12995a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
13005a8e39fbSHong Zhang       nzi += nlnk;
13015a8e39fbSHong Zhang       prow = lnk[prow];
13025a8e39fbSHong Zhang       nzbd++;
130356beaf04SBarry Smith     }
13045a8e39fbSHong Zhang     bdiag[i] = nzbd;
13055a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1306ecf371e4SBarry Smith 
13075a8e39fbSHong Zhang     /* if free space is not available, make more free space */
13085a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
13095a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1310a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1311a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
13125a8e39fbSHong Zhang       reallocs++;
13135a8e39fbSHong Zhang     }
1314ecf371e4SBarry Smith 
13155a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
13165a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13175a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
13185a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
13195a8e39fbSHong Zhang 
13205a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
13215a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
132277431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1323d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
13246cf997b0SBarry Smith     }
13255a8e39fbSHong Zhang 
13265a8e39fbSHong Zhang     current_space->array           += nzi;
13275a8e39fbSHong Zhang     current_space->local_used      += nzi;
13285a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
13295a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
13305a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
13315a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
13329e25ed09SBarry Smith   }
13335a8e39fbSHong Zhang 
1334898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1335898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
13365a8e39fbSHong Zhang 
13375a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
13385a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1339a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
13405a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1341a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13425a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
13439e25ed09SBarry Smith 
13446cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1345f64065bbSBarry Smith   {
13465a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1347ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1348ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1349ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1350ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1351335d9088SBarry Smith     if (diagonal_fill) {
1352ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1353335d9088SBarry Smith     }
1354f64065bbSBarry Smith   }
135563ba0a88SBarry Smith #endif
1356bbb6d6a8SBarry Smith 
13579e25ed09SBarry Smith   /* put together the new matrix */
1358719d5645SBarry Smith   ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr);
1359719d5645SBarry Smith   b = (Mat_SeqAIJ*)(fact)->data;
1360e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1361e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13627c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13635a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1364b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13655a8e39fbSHong Zhang   b->j          = bj;
13665a8e39fbSHong Zhang   b->i          = bi;
13675a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13685a8e39fbSHong Zhang   b->diag       = bdiag;
1369416022c9SBarry Smith   b->ilen       = 0;
1370416022c9SBarry Smith   b->imax       = 0;
1371416022c9SBarry Smith   b->row        = isrow;
1372416022c9SBarry Smith   b->col        = iscol;
1373c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1374c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
137582bf6240SBarry Smith   b->icol       = isicol;
137687828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1377416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13785a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
1379719d5645SBarry Smith   ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13805a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
1381719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1382719d5645SBarry Smith   (fact)->info.fill_ratio_given  = f;
1383719d5645SBarry Smith   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1384719d5645SBarry Smith   (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ;
1385719d5645SBarry Smith   ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr);
13863a40ed3dSBarry Smith   PetscFunctionReturn(0);
13879e25ed09SBarry Smith }
1388d5d45c9bSBarry Smith 
13893c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1390a6175056SHong Zhang #undef __FUNCT__
1391a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
13920481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
1393a6175056SHong Zhang {
1394719d5645SBarry Smith   Mat            C = B;
13950968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
13962ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
13979bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
13982ed133dbSHong Zhang   PetscErrorCode ierr;
13995d0c19d7SBarry Smith   const PetscInt *rip,*riip;
14005d0c19d7SBarry Smith   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
14012ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
14022ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1403622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1404540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1405fbf22428SSatish Balay   PetscReal      shiftpd;
1406540703acSHong Zhang   ChShift_Ctx    sctx;
1407540703acSHong Zhang   PetscInt       newshift;
1408db4efbfdSBarry Smith   PetscTruth     perm_identity;
1409d5d45c9bSBarry Smith 
1410a6175056SHong Zhang   PetscFunctionBegin;
1411117f1891Stmunson 
1412540703acSHong Zhang   shiftnz   = info->shiftnz;
1413540703acSHong Zhang   shiftpd   = info->shiftpd;
1414ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1415ee45ca4aSHong Zhang 
14162ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
14179bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
14182ed133dbSHong Zhang 
14192ed133dbSHong Zhang   /* initialization */
14202ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
14212ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
14222ed133dbSHong Zhang   jl   = il + mbs;
14232ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
14242ed133dbSHong Zhang 
1425540703acSHong Zhang   sctx.shift_amount = 0;
1426540703acSHong Zhang   sctx.nshift       = 0;
14272ed133dbSHong Zhang   do {
1428540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
14292ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14302ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1431a6175056SHong Zhang     }
14322ed133dbSHong Zhang 
14332ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1434c05c3958SHong Zhang       bval = ba + bi[k];
14352ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
14362ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
14372ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
14389bfd6278SHong Zhang         col = riip[aj[j]];
14392ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
14402ed133dbSHong Zhang           rtmp[col] = aa[j];
1441c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
14422ed133dbSHong Zhang         }
14432ed133dbSHong Zhang       }
144439e3d630SHong Zhang       /* shift the diagonal of the matrix */
1445540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
14462ed133dbSHong Zhang 
14472ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
14482ed133dbSHong Zhang       dk = rtmp[k];
14492ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
14502ed133dbSHong Zhang 
14512ed133dbSHong Zhang       while (i < k){
14522ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
14532ed133dbSHong Zhang 
14542ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
14552ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
14562ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
14572ed133dbSHong Zhang         dk += uikdi*ba[ili];
14582ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14592ed133dbSHong Zhang 
14602ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14612ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14622ed133dbSHong Zhang         if (jmin < jmax){
14632ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14642ed133dbSHong Zhang           /* update il and jl for row i */
14652ed133dbSHong Zhang           il[i] = jmin;
14662ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14672ed133dbSHong Zhang         }
14682ed133dbSHong Zhang         i = nexti;
14692ed133dbSHong Zhang       }
14702ed133dbSHong Zhang 
14710697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14720697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14730697b6caSHong Zhang       rs   = 0.0;
14743c9e8598SHong Zhang       jmin = bi[k]+1;
14753c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14763c9e8598SHong Zhang       bcol = bj + jmin;
14773c9e8598SHong Zhang       while (nz--){
147889f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
147989f845c8SHong Zhang         bcol++;
14803c9e8598SHong Zhang       }
1481540703acSHong Zhang 
1482540703acSHong Zhang       sctx.rs = rs;
1483540703acSHong Zhang       sctx.pv = dk;
14845b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1485117f1891Stmunson 
1486117f1891Stmunson       if (newshift == 1) {
1487117f1891Stmunson         if (!sctx.shift_amount) {
1488117f1891Stmunson           sctx.shift_amount = 1e-5;
1489117f1891Stmunson         }
1490117f1891Stmunson         break;
1491117f1891Stmunson       }
14923c9e8598SHong Zhang 
14933c9e8598SHong Zhang       /* copy data into U(k,:) */
149439e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
149539e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
149639e3d630SHong Zhang       if (jmin < jmax) {
149739e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
149839e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
14993c9e8598SHong Zhang         }
150039e3d630SHong Zhang         /* add the k-th row into il and jl */
15013c9e8598SHong Zhang         il[k] = jmin;
15023c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
15033c9e8598SHong Zhang       }
15043c9e8598SHong Zhang     }
1505540703acSHong Zhang   } while (sctx.chshift);
15063c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
15073c9e8598SHong Zhang 
150839e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
15099bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
1510db4efbfdSBarry Smith 
1511db4efbfdSBarry Smith   ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr);
1512db4efbfdSBarry Smith   if (perm_identity){
1513719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1514719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1515719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1516719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1517db4efbfdSBarry Smith   } else {
1518719d5645SBarry Smith     (B)->ops->solve           = MatSolve_SeqSBAIJ_1;
1519719d5645SBarry Smith     (B)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1520719d5645SBarry Smith     (B)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1521719d5645SBarry Smith     (B)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
1522db4efbfdSBarry Smith   }
1523db4efbfdSBarry Smith 
15243c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
15253c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1526d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1527540703acSHong Zhang   if (sctx.nshift){
1528540703acSHong Zhang     if (shiftnz) {
15291e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1530540703acSHong Zhang     } else if (shiftpd) {
15311e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
15320697b6caSHong Zhang     }
15333c9e8598SHong Zhang   }
15343c9e8598SHong Zhang   PetscFunctionReturn(0);
15353c9e8598SHong Zhang }
1536a6175056SHong Zhang 
1537a6175056SHong Zhang #undef __FUNCT__
1538a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
15390481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1540a6175056SHong Zhang {
15410968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1542ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1543dfbe8321SBarry Smith   PetscErrorCode     ierr;
154458ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
15455d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
15465d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
1547ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
154858ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
15495a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1550ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1551a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1552a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
15537a48dd6fSHong Zhang   PetscBT            lnkbt;
1554b635d36bSHong Zhang   IS                 iperm;
1555a6175056SHong Zhang 
1556a6175056SHong Zhang   PetscFunctionBegin;
1557d0f46423SBarry 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);
155858ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
155958ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1560653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1561b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
15627a48dd6fSHong Zhang 
156339e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
156439e3d630SHong Zhang   ui[0] = 0;
156539e3d630SHong Zhang 
1566b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1567622e2028SHong Zhang   if (!levels && perm_identity) {
156858ebbce7SBarry Smith 
1569ed59401aSHong Zhang     for (i=0; i<am; i++) {
157039e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1571ed59401aSHong Zhang     }
157239e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
157339e3d630SHong Zhang     cols = uj;
1574ed59401aSHong Zhang     for (i=0; i<am; i++) {
1575ed59401aSHong Zhang       aj    = a->j + a->diag[i];
157639e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
157739e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1578ed59401aSHong Zhang     }
157939e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1580b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1581b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1582b635d36bSHong Zhang 
15837a48dd6fSHong Zhang     /* initialization */
15845a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15857a48dd6fSHong Zhang 
15867a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
15877a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
15881393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
15897a48dd6fSHong Zhang     il         = jl + am;
15907a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
15917a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
15927a48dd6fSHong Zhang     for (i=0; i<am; i++){
15937a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
15947a48dd6fSHong Zhang     }
15957a48dd6fSHong Zhang 
15967a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
15977a48dd6fSHong Zhang     nlnk = am + 1;
15982ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15997a48dd6fSHong Zhang 
16007a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1601a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
16027a48dd6fSHong Zhang     current_space = free_space;
1603a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
16047a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
16057a48dd6fSHong Zhang 
16067a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
16077a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
16087a48dd6fSHong Zhang       nzk   = 0;
1609622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1610b635d36bSHong Zhang       if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
1611622e2028SHong Zhang       ncols_upper = 0;
1612622e2028SHong Zhang       for (j=0; j<ncols; j++){
1613b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1614b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
16155a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1616622e2028SHong Zhang           ncols_upper++;
1617622e2028SHong Zhang         }
1618622e2028SHong Zhang       }
1619b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
16207a48dd6fSHong Zhang       nzk += nlnk;
16217a48dd6fSHong Zhang 
16227a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
16237a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
16247a48dd6fSHong Zhang 
16257a48dd6fSHong Zhang       while (prow < k){
16267a48dd6fSHong Zhang         nextprow = jl[prow];
16277a48dd6fSHong Zhang 
16287a48dd6fSHong Zhang         /* merge prow into k-th row */
16297a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
16307a48dd6fSHong Zhang         jmax = ui[prow+1];
16317a48dd6fSHong Zhang         ncols = jmax-jmin;
1632ed59401aSHong Zhang         i     = jmin - ui[prow];
1633ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
16345a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
16355a8e39fbSHong Zhang         j     = *(uj - 1);
16365a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
16377a48dd6fSHong Zhang         nzk += nlnk;
16387a48dd6fSHong Zhang 
16397a48dd6fSHong Zhang         /* update il and jl for prow */
16407a48dd6fSHong Zhang         if (jmin < jmax){
16417a48dd6fSHong Zhang           il[prow] = jmin;
16427a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
16437a48dd6fSHong Zhang         }
16447a48dd6fSHong Zhang         prow = nextprow;
16457a48dd6fSHong Zhang       }
16467a48dd6fSHong Zhang 
16477a48dd6fSHong Zhang       /* if free space is not available, make more free space */
16487a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
16497a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
16507a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1651a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1652a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
16537a48dd6fSHong Zhang         reallocs++;
16547a48dd6fSHong Zhang       }
16557a48dd6fSHong Zhang 
16567a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1657b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
16582ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
16597a48dd6fSHong Zhang 
16607a48dd6fSHong Zhang       /* add the k-th row into il and jl */
166139e3d630SHong Zhang       if (nzk > 1){
16627a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
16637a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
16647a48dd6fSHong Zhang         il[k] = ui[k] + 1;
16657a48dd6fSHong Zhang       }
16667a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
16677a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
16687a48dd6fSHong Zhang 
16697a48dd6fSHong Zhang       current_space->array           += nzk;
16707a48dd6fSHong Zhang       current_space->local_used      += nzk;
16717a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16727a48dd6fSHong Zhang 
16737a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16747a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16757a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16767a48dd6fSHong Zhang 
16777a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16787a48dd6fSHong Zhang     }
16797a48dd6fSHong Zhang 
16806cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16817a48dd6fSHong Zhang     if (ai[am] != 0) {
168239e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1683ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1684ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1685ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16867a48dd6fSHong Zhang     } else {
1687ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16887a48dd6fSHong Zhang     }
168963ba0a88SBarry Smith #endif
16907a48dd6fSHong Zhang 
16917a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1692b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16937a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
16945a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
16957a48dd6fSHong Zhang 
16967a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
16977a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1698a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
16992ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1700a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
17017a48dd6fSHong Zhang 
170239e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
170339e3d630SHong Zhang 
17047a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
17057a48dd6fSHong Zhang 
1706719d5645SBarry Smith   b    = (Mat_SeqSBAIJ*)(fact)->data;
17077a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
17087a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
17097a48dd6fSHong Zhang   b->j    = uj;
17107a48dd6fSHong Zhang   b->i    = ui;
17117a48dd6fSHong Zhang   b->diag = 0;
17127a48dd6fSHong Zhang   b->ilen = 0;
17137a48dd6fSHong Zhang   b->imax = 0;
17147a48dd6fSHong Zhang   b->row  = perm;
1715b635d36bSHong Zhang   b->col  = perm;
1716b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1717b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1718b635d36bSHong Zhang   b->icol = iperm;
17197a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
17207a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1721719d5645SBarry Smith   ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
17227a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1723e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1724e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
17257a48dd6fSHong Zhang 
1726719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1727719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
17287a48dd6fSHong Zhang   if (ai[am] != 0) {
1729719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
17307a48dd6fSHong Zhang   } else {
1731719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
17327a48dd6fSHong Zhang   }
1733719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1734a6175056SHong Zhang   PetscFunctionReturn(0);
1735a6175056SHong Zhang }
1736d5d45c9bSBarry Smith 
1737f76d2b81SHong Zhang #undef __FUNCT__
1738f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
17390481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1740f76d2b81SHong Zhang {
1741f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
174210c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
1743dfbe8321SBarry Smith   PetscErrorCode     ierr;
1744f76d2b81SHong Zhang   PetscTruth         perm_identity;
174510c27e3dSHong Zhang   PetscReal          fill = info->fill;
17465d0c19d7SBarry Smith   const PetscInt     *rip,*riip;
17475d0c19d7SBarry Smith   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
174810c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
17492ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1750a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
175110c27e3dSHong Zhang   PetscBT            lnkbt;
17522ed133dbSHong Zhang   IS                 iperm;
1753f76d2b81SHong Zhang 
1754f76d2b81SHong Zhang   PetscFunctionBegin;
1755d0f46423SBarry 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);
17562ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1757f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
17582ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
17592ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
17609bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
176110c27e3dSHong Zhang 
176210c27e3dSHong Zhang   /* initialization */
17631393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
176410c27e3dSHong Zhang   ui[0] = 0;
176510c27e3dSHong Zhang 
176610c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17671393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17681393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17691393bc97SHong Zhang   il     = jl + am;
17701393bc97SHong Zhang   cols   = il + am;
17711393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17721393bc97SHong Zhang   for (i=0; i<am; i++){
17731393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1774f76d2b81SHong Zhang   }
177510c27e3dSHong Zhang 
177610c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17771393bc97SHong Zhang   nlnk = am + 1;
17781393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
177910c27e3dSHong Zhang 
17801393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1781a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
178210c27e3dSHong Zhang   current_space = free_space;
178310c27e3dSHong Zhang 
17841393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
178510c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
178610c27e3dSHong Zhang     nzk   = 0;
17872ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
17889bfd6278SHong Zhang     if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
17892ed133dbSHong Zhang     ncols_upper = 0;
17902ed133dbSHong Zhang     for (j=0; j<ncols; j++){
17919bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
17922ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
17932ed133dbSHong Zhang         cols[ncols_upper] = i;
17942ed133dbSHong Zhang         ncols_upper++;
17952ed133dbSHong Zhang       }
17962ed133dbSHong Zhang     }
17971393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
179810c27e3dSHong Zhang     nzk += nlnk;
179910c27e3dSHong Zhang 
180010c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
180110c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
180210c27e3dSHong Zhang 
180310c27e3dSHong Zhang     while (prow < k){
180410c27e3dSHong Zhang       nextprow = jl[prow];
180510c27e3dSHong Zhang       /* merge prow into k-th row */
18061393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
180710c27e3dSHong Zhang       jmax = ui[prow+1];
180810c27e3dSHong Zhang       ncols = jmax-jmin;
18091393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
18105a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
181110c27e3dSHong Zhang       nzk += nlnk;
181210c27e3dSHong Zhang 
181310c27e3dSHong Zhang       /* update il and jl for prow */
181410c27e3dSHong Zhang       if (jmin < jmax){
181510c27e3dSHong Zhang         il[prow] = jmin;
18162ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
181710c27e3dSHong Zhang       }
181810c27e3dSHong Zhang       prow = nextprow;
181910c27e3dSHong Zhang     }
182010c27e3dSHong Zhang 
182110c27e3dSHong Zhang     /* if free space is not available, make more free space */
182210c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
18231393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
182410c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1825a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
182610c27e3dSHong Zhang       reallocs++;
182710c27e3dSHong Zhang     }
182810c27e3dSHong Zhang 
182910c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
18301393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
183110c27e3dSHong Zhang 
183210c27e3dSHong Zhang     /* add the k-th row into il and jl */
183310c27e3dSHong Zhang     if (nzk-1 > 0){
18341393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
183510c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
183610c27e3dSHong Zhang       il[k] = ui[k] + 1;
183710c27e3dSHong Zhang     }
18382ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
183910c27e3dSHong Zhang     current_space->array           += nzk;
184010c27e3dSHong Zhang     current_space->local_used      += nzk;
184110c27e3dSHong Zhang     current_space->local_remaining -= nzk;
184210c27e3dSHong Zhang 
184310c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
184410c27e3dSHong Zhang   }
184510c27e3dSHong Zhang 
18466cf91177SBarry Smith #if defined(PETSC_USE_INFO)
18471393bc97SHong Zhang   if (ai[am] != 0) {
18481393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1849ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1850ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1851ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
185210c27e3dSHong Zhang   } else {
1853ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
185410c27e3dSHong Zhang   }
185563ba0a88SBarry Smith #endif
185610c27e3dSHong Zhang 
185710c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
18589bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
185910c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
186010c27e3dSHong Zhang 
186110c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18621393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1863a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
186410c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
186510c27e3dSHong Zhang 
186610c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
186710c27e3dSHong Zhang 
1868719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
186910c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1870e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1871e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18721393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
187310c27e3dSHong Zhang   b->j    = uj;
187410c27e3dSHong Zhang   b->i    = ui;
187510c27e3dSHong Zhang   b->diag = 0;
187610c27e3dSHong Zhang   b->ilen = 0;
187710c27e3dSHong Zhang   b->imax = 0;
187810c27e3dSHong Zhang   b->row  = perm;
18799bfd6278SHong Zhang   b->col  = perm;
18809bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18819bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18829bfd6278SHong Zhang   b->icol = iperm;
188310c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18841393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1885719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18861393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
188710c27e3dSHong Zhang 
1888719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
1889719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
18901393bc97SHong Zhang   if (ai[am] != 0) {
1891719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
189210c27e3dSHong Zhang   } else {
1893719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
189410c27e3dSHong Zhang   }
1895719d5645SBarry Smith   (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1896f76d2b81SHong Zhang   PetscFunctionReturn(0);
1897f76d2b81SHong Zhang }
1898