xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision d0f46423f772d871f92d805d83c5bf0c050e33b4)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*);
26e34fafa9SBarry Smith #endif
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
63*d0f46423SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->rmap->n;
6438baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6538baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6638baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6754f21887SBarry Smith   MatScalar      *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6854a8161fSBarry Smith   PetscReal      af;
6986bacbd2SBarry Smith 
7086bacbd2SBarry Smith   PetscFunctionBegin;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7338baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7486bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7533259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7638baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7738baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7886bacbd2SBarry Smith 
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith   /* ------------------------------------------------------------
8186bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8286bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8386bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8486bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8586bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8686bacbd2SBarry Smith      ------------------------------------------------------------
8786bacbd2SBarry Smith   */
8886bacbd2SBarry Smith 
8986bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9086bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9186bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9286bacbd2SBarry Smith   reorder = PetscNot(reorder);
9386bacbd2SBarry Smith 
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith   /* storage for ilu factor */
9638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9738baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
98a77337e4SBarry Smith   ierr = PetscMalloc(jmax*sizeof(MatScalar),&new_a);CHKERRQ(ierr);
9938baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10086bacbd2SBarry Smith 
10186bacbd2SBarry Smith   /* ------------------------------------------------------------
10286bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10386bacbd2SBarry Smith      ------------------------------------------------------------
10486bacbd2SBarry Smith   */
105b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
106b19713cbSBarry Smith     old_j[i]++;
10786bacbd2SBarry Smith   }
108b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
109b19713cbSBarry Smith     old_i[i]++;
110b19713cbSBarry Smith   };
111010a20caSHong Zhang 
11286bacbd2SBarry Smith 
11386bacbd2SBarry Smith   if (reorder) {
114c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
115c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
116c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
117c0c2c81eSBarry Smith       r[i]  = r[i]+1;
118c0c2c81eSBarry Smith       c[i]  = c[i]+1;
119c0c2c81eSBarry Smith     }
12038baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12138baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
122a77337e4SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr);
1235bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
124c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
125c0c2c81eSBarry Smith       r[i]  = r[i]-1;
126c0c2c81eSBarry Smith       c[i]  = c[i]-1;
127c0c2c81eSBarry Smith     }
128c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
129c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
130b19713cbSBarry Smith     o_a = old_a2;
131b19713cbSBarry Smith     o_j = old_j2;
132b19713cbSBarry Smith     o_i = old_i2;
133b19713cbSBarry Smith   } else {
134b19713cbSBarry Smith     o_a = old_a;
135b19713cbSBarry Smith     o_j = old_j;
136b19713cbSBarry Smith     o_i = old_i;
13786bacbd2SBarry Smith   }
13886bacbd2SBarry Smith 
13986bacbd2SBarry Smith   /* ------------------------------------------------------------
14086bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14186bacbd2SBarry Smith      ------------------------------------------------------------
14286bacbd2SBarry Smith   */
14386bacbd2SBarry Smith 
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14538baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14687828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14786bacbd2SBarry Smith 
14854a8161fSBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1496849ba73SBarry Smith   if (sierr) {
1506849ba73SBarry Smith     switch (sierr) {
151a83599f4SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
152a83599f4SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
153e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
154e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15577431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15677431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15786bacbd2SBarry Smith     }
15886bacbd2SBarry Smith   }
15986bacbd2SBarry Smith 
16086bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16186bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16286bacbd2SBarry Smith 
16386bacbd2SBarry Smith   /* ------------------------------------------------------------
16486bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16586bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16686bacbd2SBarry Smith      ------------------------------------------------------------
16786bacbd2SBarry Smith   */
16886bacbd2SBarry Smith 
16987828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
17038baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17186bacbd2SBarry Smith 
1725bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17386bacbd2SBarry Smith 
17486bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17586bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17686bacbd2SBarry Smith 
17786bacbd2SBarry Smith   if (reorder) {
17886bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18086bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
181313ae024SBarry Smith   } else {
182b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
183f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
184b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
185b19713cbSBarry Smith     }
186b19713cbSBarry Smith   }
187b19713cbSBarry Smith 
188b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
190b19713cbSBarry Smith     old_i[i]--;
191b19713cbSBarry Smith   }
192b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
193b19713cbSBarry Smith     old_j[i]--;
194b19713cbSBarry Smith   }
19586bacbd2SBarry Smith 
196b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
198b19713cbSBarry Smith     new_i[i]--;
199b19713cbSBarry Smith   }
200b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
201b19713cbSBarry Smith     new_j[i]--;
202b19713cbSBarry Smith   }
20386bacbd2SBarry Smith 
20486bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20586bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2064c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2074c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
208c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20986bacbd2SBarry Smith   for(i=0; i<n; i++) {
21086bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21186bacbd2SBarry Smith   };
21286bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21307d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21486bacbd2SBarry Smith 
215c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
216c0c2c81eSBarry Smith 
217329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21807d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21986bacbd2SBarry Smith 
22086bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22186bacbd2SBarry Smith 
2227adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
223f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
2247adad957SLisandro Dalcin   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
225ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2265c9eb25fSBarry Smith   (*fact)->factor    = MAT_FACTOR_LU;
22786bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22886bacbd2SBarry Smith 
22986bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
230e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
231e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
23207d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23386bacbd2SBarry Smith   b->a             = new_a;
23486bacbd2SBarry Smith   b->j             = new_j;
23586bacbd2SBarry Smith   b->i             = new_i;
23686bacbd2SBarry Smith   b->ilen          = 0;
23786bacbd2SBarry Smith   b->imax          = 0;
2381f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
239313ae024SBarry Smith   b->row           = isirow;
24086bacbd2SBarry Smith   b->col           = iscolf;
24187828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24286bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24386bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24486bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24586bacbd2SBarry Smith 
246431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
247ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
248ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
250ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
251431e710aSBarry Smith 
2527529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25371c2f376SKris Buschelman 
25486bacbd2SBarry Smith   PetscFunctionReturn(0);
25560ec86bdSBarry Smith #endif
25686bacbd2SBarry Smith }
25786bacbd2SBarry Smith 
258e631078cSBarry Smith EXTERN_C_BEGIN
2594a2ae208SSatish Balay #undef __FUNCT__
260b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc"
261b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
262b24902e0SBarry Smith {
263*d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
264b24902e0SBarry Smith   PetscErrorCode     ierr;
265b24902e0SBarry Smith 
266b24902e0SBarry Smith   PetscFunctionBegin;
267b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
268b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
269b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
270b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr);
271b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
272b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
273b24902e0SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2745c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
2755c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
2765c9eb25fSBarry Smith     (*B)->ops->choleskyfactornumeric  = MatCholeskyFactorNumeric_SeqAIJ;
277b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
278b24902e0SBarry Smith   PetscFunctionReturn(0);
279b24902e0SBarry Smith }
280e631078cSBarry Smith EXTERN_C_END
281b24902e0SBarry Smith 
282b24902e0SBarry Smith #undef __FUNCT__
283b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
284dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
285289bc588SBarry Smith {
286416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
287289bc588SBarry Smith   IS                 isicol;
2886849ba73SBarry Smith   PetscErrorCode     ierr;
289*d0f46423SBarry Smith   PetscInt           *r,*ic,i,n=A->rmap->n,*ai=a->i,*aj=a->j;
2901393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
2911393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2929e046878SBarry Smith   PetscReal          f;
2934875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
294a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2951393bc97SHong Zhang   PetscBT            lnkbt;
296289bc588SBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
298*d0f46423SBarry Smith   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2994c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
300ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
301ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
302289bc588SBarry Smith 
303289bc588SBarry Smith   /* get new row pointers */
3041393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
3051393bc97SHong Zhang   bi[0] = 0;
3061393bc97SHong Zhang 
3071393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
3081393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
3091393bc97SHong Zhang   bdiag[0] = 0;
3101393bc97SHong Zhang 
3111393bc97SHong Zhang   /* linked list for storing column indices of the active row */
3121393bc97SHong Zhang   nlnk = n + 1;
3131393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3141393bc97SHong Zhang 
31535baf27eSSatish Balay   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
3161393bc97SHong Zhang 
3171393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
318d3d32019SBarry Smith   f = info->fill;
319a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3201393bc97SHong Zhang   current_space = free_space;
321289bc588SBarry Smith 
322289bc588SBarry Smith   for (i=0; i<n; i++) {
3231393bc97SHong Zhang     /* copy previous fill into linked list */
324289bc588SBarry Smith     nzi = 0;
3251393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3261393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3271393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
328afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3291393bc97SHong Zhang     nzi += nlnk;
3301393bc97SHong Zhang 
3311393bc97SHong Zhang     /* add pivot rows into linked list */
3321393bc97SHong Zhang     row = lnk[n];
3331393bc97SHong Zhang     while (row < i) {
334d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
335d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
336d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3371393bc97SHong Zhang       nzi += nlnk;
3381393bc97SHong Zhang       row  = lnk[row];
339289bc588SBarry Smith     }
3401393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3411393bc97SHong Zhang     im[i]   = nzi;
3421393bc97SHong Zhang 
3431393bc97SHong Zhang     /* mark bdiag */
3441393bc97SHong Zhang     nzbd = 0;
3451393bc97SHong Zhang     nnz  = nzi;
3461393bc97SHong Zhang     k    = lnk[n];
3471393bc97SHong Zhang     while (nnz-- && k < i){
3481393bc97SHong Zhang       nzbd++;
3491393bc97SHong Zhang       k = lnk[k];
3501393bc97SHong Zhang     }
3511393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3521393bc97SHong Zhang 
3531393bc97SHong Zhang     /* if free space is not available, make more free space */
3541393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3551393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
356a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3571393bc97SHong Zhang       reallocs++;
3581393bc97SHong Zhang     }
3591393bc97SHong Zhang 
3601393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3611393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3621393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3631393bc97SHong Zhang     current_space->array           += nzi;
3641393bc97SHong Zhang     current_space->local_used      += nzi;
3651393bc97SHong Zhang     current_space->local_remaining -= nzi;
366289bc588SBarry Smith   }
3676cf91177SBarry Smith #if defined(PETSC_USE_INFO)
368464e8d44SSatish Balay   if (ai[n] != 0) {
3691393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
370ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
371ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
372ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
373ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
37405bf559cSSatish Balay   } else {
375ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
37605bf559cSSatish Balay   }
37763ba0a88SBarry Smith #endif
3782fb3ed76SBarry Smith 
379898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
380898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3811987afe7SBarry Smith 
3821393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3831393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
384a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3851393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
38635baf27eSSatish Balay   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
387289bc588SBarry Smith 
388289bc588SBarry Smith   /* put together the new matrix */
3895c9eb25fSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
39052e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
391416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
392e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
393e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
3947c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3951393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3961393bc97SHong Zhang   b->j          = bj;
3971393bc97SHong Zhang   b->i          = bi;
3981393bc97SHong Zhang   b->diag       = bdiag;
399416022c9SBarry Smith   b->ilen       = 0;
400416022c9SBarry Smith   b->imax       = 0;
401416022c9SBarry Smith   b->row        = isrow;
402416022c9SBarry Smith   b->col        = iscol;
403c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
404c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40582bf6240SBarry Smith   b->icol       = isicol;
40687828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
4071393bc97SHong Zhang 
4081393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
4091393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
4101393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
4117fd98bd6SLois Curfman McInnes 
4125c9eb25fSBarry Smith   (*B)->factor                 =  MAT_FACTOR_LU;
413418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
414ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
415703deb49SSatish Balay 
416e93ccc4dSBarry Smith   if (ai[n] != 0) {
4171393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
418eea2913aSSatish Balay   } else {
419eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
420eea2913aSSatish Balay   }
4214846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
42271c2f376SKris Buschelman   (*B)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
424289bc588SBarry Smith }
4251393bc97SHong Zhang 
4265b5bc046SBarry Smith /*
4275b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4285b5bc046SBarry Smith */
4295b5bc046SBarry Smith #undef __FUNCT__
4305b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4315b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4325b5bc046SBarry Smith {
4335b5bc046SBarry Smith   PetscErrorCode ierr;
4345b5bc046SBarry Smith   PetscTruth     flg;
4355b5bc046SBarry Smith 
4365b5bc046SBarry Smith   PetscFunctionBegin;
4375b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4385b5bc046SBarry Smith   if (flg) {
4395b5bc046SBarry Smith     PetscViewer viewer;
4405b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4415b5bc046SBarry Smith 
4425b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4437adad957SLisandro Dalcin     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4445b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4455b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4465b5bc046SBarry Smith   }
4475b5bc046SBarry Smith   PetscFunctionReturn(0);
4485b5bc046SBarry Smith }
4495b5bc046SBarry Smith 
4500a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4514a2ae208SSatish Balay #undef __FUNCT__
4524a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
453af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
454289bc588SBarry Smith {
455416022c9SBarry Smith   Mat            C=*B;
456416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
45782bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4586849ba73SBarry Smith   PetscErrorCode ierr;
459*d0f46423SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap->n,*bi=b->i,*bj=b->j;
460b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
461d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
46254f21887SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,*rtmps;
46354f21887SBarry Smith   MatScalar      *v,*pv;
4646a7c8fc2SHong Zhang   PetscScalar    d;
4656a7c8fc2SHong Zhang   PetscReal      rs;
466b3bf805bSHong Zhang   LUShift_Ctx    sctx;
46742898933SHong Zhang   PetscInt       newshift,*ddiag;
468289bc588SBarry Smith 
4693a40ed3dSBarry Smith   PetscFunctionBegin;
47078b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
47178b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
47287828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
47387828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
474010a20caSHong Zhang   rtmps = rtmp; ics = ic;
475289bc588SBarry Smith 
476ada7143aSSatish Balay   sctx.shift_top  = 0;
477ada7143aSSatish Balay   sctx.nshift_max = 0;
478ada7143aSSatish Balay   sctx.shift_lo   = 0;
479ada7143aSSatish Balay   sctx.shift_hi   = 0;
480ada7143aSSatish Balay 
481118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
482118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4831a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
4849deb78dcSHong Zhang     PetscInt *aai = a->i;
48542898933SHong Zhang     ddiag          = a->diag;
486118f5789SBarry Smith     sctx.shift_top = 0;
4876cc28720Svictorle     for (i=0; i<n; i++) {
4889f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4899f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4901a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
491010a20caSHong Zhang       v  = a->a+aai[i];
492e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
493e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4941a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4951a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4966cc28720Svictorle     }
497c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
498118f5789SBarry Smith     sctx.shift_top    *= 1.1;
499d2276718SHong Zhang     sctx.nshift_max   = 5;
500d2276718SHong Zhang     sctx.shift_lo     = 0.;
501d2276718SHong Zhang     sctx.shift_hi     = 1.;
502a0a356efSHong Zhang   }
503a0a356efSHong Zhang 
504a0a356efSHong Zhang   sctx.shift_amount = 0;
505a0a356efSHong Zhang   sctx.nshift       = 0;
506085256b3SBarry Smith   do {
507d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
508289bc588SBarry Smith     for (i=0; i<n; i++){
509b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
510b3bf805bSHong Zhang       bjtmp = bj + bi[i];
511b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
512289bc588SBarry Smith 
513289bc588SBarry Smith       /* load in initial (unfactored row) */
514416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
515b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
516010a20caSHong Zhang       v     = a->a + a->i[r[i]];
517085256b3SBarry Smith       for (j=0; j<nz; j++) {
518b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
519085256b3SBarry Smith       }
520d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
521289bc588SBarry Smith 
522b3bf805bSHong Zhang       row = *bjtmp++;
523f2582111SSatish Balay       while  (row < i) {
5248c37ef55SBarry Smith         pc = rtmp + row;
5258c37ef55SBarry Smith         if (*pc != 0.0) {
526010a20caSHong Zhang           pv         = b->a + diag_offset[row];
527010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
52835aab85fSBarry Smith           multiplier = *pc / *pv++;
5298c37ef55SBarry Smith           *pc        = multiplier;
530b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
531f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
532efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
533289bc588SBarry Smith         }
534b3bf805bSHong Zhang         row = *bjtmp++;
535289bc588SBarry Smith       }
536416022c9SBarry Smith       /* finished row so stick it into b->a */
537b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
538b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
539b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
540b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
541d3d32019SBarry Smith       rs   = 0.0;
542d3d32019SBarry Smith       for (j=0; j<nz; j++) {
543d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
544d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
545d3d32019SBarry Smith       }
546d2276718SHong Zhang 
547b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
548d2276718SHong Zhang       sctx.rs  = rs;
549d2276718SHong Zhang       sctx.pv  = pv[diag];
550c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
5515b5bc046SBarry Smith       if (newshift == 1) break;
55235aab85fSBarry Smith     }
553d2276718SHong Zhang 
554118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5556cc28720Svictorle       /*
5566c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5576cc28720Svictorle        * then try lower shift
5586cc28720Svictorle        */
559d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
560d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
561d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
562d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
563d2276718SHong Zhang       sctx.nshift++;
5646cc28720Svictorle     }
565d2276718SHong Zhang   } while (sctx.lushift);
5660f11f9aeSSatish Balay 
56735aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
56835aab85fSBarry Smith   for (i=0; i<n; i++) {
569010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
57035aab85fSBarry Smith   }
571606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
57278b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
57378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
5745c9eb25fSBarry Smith   C->factor       = MAT_FACTOR_LU;
575c456f294SBarry Smith   C->assembled    = PETSC_TRUE;
5765c9eb25fSBarry Smith   C->preallocated = PETSC_TRUE;
577*d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
578d2276718SHong Zhang   if (sctx.nshift){
579118f5789SBarry Smith     if (info->shiftnz) {
5801e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
581118f5789SBarry Smith     } else if (info->shiftpd) {
5821e2582c4SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
5836cc28720Svictorle     }
5840697b6caSHong Zhang   }
5853a40ed3dSBarry Smith   PetscFunctionReturn(0);
586289bc588SBarry Smith }
5870889b5dcSvictorle 
588137fb511SHong Zhang /*
589137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
590137fb511SHong Zhang    Input:
591137fb511SHong Zhang      A - original matrix
592137fb511SHong Zhang    Output;
593137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
594137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
595137fb511SHong Zhang          a->a reordered accordingly with a->j
596137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
597137fb511SHong Zhang */
598137fb511SHong Zhang #undef __FUNCT__
599137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
600137fb511SHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat A,MatFactorInfo *info,Mat *B)
601137fb511SHong Zhang {
602b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
603b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
604137fb511SHong Zhang   PetscErrorCode ierr;
605*d0f46423SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
606b051339dSHong Zhang   PetscInt       *ajtmp,nz,row,*ics;
607b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
608a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
609a77337e4SBarry Smith   MatScalar      *v,*pv;
610137fb511SHong Zhang   PetscReal      rs;
611137fb511SHong Zhang   LUShift_Ctx    sctx;
612b051339dSHong Zhang   PetscInt       newshift;
613137fb511SHong Zhang 
614137fb511SHong Zhang   PetscFunctionBegin;
615b051339dSHong Zhang   if (A != *B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
616137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
617137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
618137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
619137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
620b051339dSHong Zhang   ics = ic;
621137fb511SHong Zhang 
622137fb511SHong Zhang   sctx.shift_top  = 0;
623137fb511SHong Zhang   sctx.nshift_max = 0;
624137fb511SHong Zhang   sctx.shift_lo   = 0;
625137fb511SHong Zhang   sctx.shift_hi   = 0;
626137fb511SHong Zhang 
627137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
628137fb511SHong Zhang   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
629137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
630137fb511SHong Zhang     sctx.shift_top = 0;
631137fb511SHong Zhang     for (i=0; i<n; i++) {
632137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
633b051339dSHong Zhang       d  = (a->a)[diag[i]];
634137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
635b051339dSHong Zhang       v  = a->a+ai[i];
636b051339dSHong Zhang       nz = ai[i+1] - ai[i];
637137fb511SHong Zhang       for (j=0; j<nz; j++)
638137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
639137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
640137fb511SHong Zhang     }
641137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
642137fb511SHong Zhang     sctx.shift_top    *= 1.1;
643137fb511SHong Zhang     sctx.nshift_max   = 5;
644137fb511SHong Zhang     sctx.shift_lo     = 0.;
645137fb511SHong Zhang     sctx.shift_hi     = 1.;
646137fb511SHong Zhang   }
647137fb511SHong Zhang 
648137fb511SHong Zhang   sctx.shift_amount = 0;
649137fb511SHong Zhang   sctx.nshift       = 0;
650137fb511SHong Zhang   do {
651137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
652137fb511SHong Zhang     for (i=0; i<n; i++){
653b051339dSHong Zhang       /* load in initial unfactored row */
654b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
655b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
656b051339dSHong Zhang       v     = a->a + ai[r[i]];
657137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
658b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
659137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
660137fb511SHong Zhang 
661b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
662137fb511SHong Zhang       for (j=0; j<nz; j++) {
663137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
664b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
665137fb511SHong Zhang       }
666137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
667137fb511SHong Zhang 
668b051339dSHong Zhang       row = *ajtmp++;
669137fb511SHong Zhang       while  (row < i) {
670137fb511SHong Zhang         pc = rtmp + row;
671137fb511SHong Zhang         if (*pc != 0.0) {
672b051339dSHong Zhang           pv         = a->a + diag[r[row]];
673b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
674137fb511SHong Zhang 
675137fb511SHong Zhang           multiplier = *pc / *pv++;
676137fb511SHong Zhang           *pc        = multiplier;
677b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
678b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
679137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
680137fb511SHong Zhang         }
681b051339dSHong Zhang         row = *ajtmp++;
682137fb511SHong Zhang       }
683b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
684b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
685b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
686b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
687b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
688137fb511SHong Zhang 
689137fb511SHong Zhang       rs   = 0.0;
690137fb511SHong Zhang       for (j=0; j<nz; j++) {
691b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
692b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
693137fb511SHong Zhang       }
694137fb511SHong Zhang 
695137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
696137fb511SHong Zhang       sctx.rs  = rs;
697b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
698c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
699137fb511SHong Zhang       if (newshift == 1) break;
700137fb511SHong Zhang     }
701137fb511SHong Zhang 
702137fb511SHong Zhang     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
703137fb511SHong Zhang       /*
704137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
705137fb511SHong Zhang        * then try lower shift
706137fb511SHong Zhang        */
707137fb511SHong Zhang       sctx.shift_hi        = info->shift_fraction;
708137fb511SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
709137fb511SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
710137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
711137fb511SHong Zhang       sctx.nshift++;
712137fb511SHong Zhang     }
713137fb511SHong Zhang   } while (sctx.lushift);
714137fb511SHong Zhang 
715137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
716137fb511SHong Zhang   for (i=0; i<n; i++) {
717b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
718137fb511SHong Zhang   }
719137fb511SHong Zhang 
720137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
721137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
722137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
7235c9eb25fSBarry Smith   A->factor     = MAT_FACTOR_LU;
724b051339dSHong Zhang   A->ops->solve = MatSolve_SeqAIJ_InplaceWithPerm;
725b051339dSHong Zhang   A->assembled = PETSC_TRUE;
7265c9eb25fSBarry Smith   A->preallocated = PETSC_TRUE;
727*d0f46423SBarry Smith   ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr);
728137fb511SHong Zhang   if (sctx.nshift){
729137fb511SHong Zhang     if (info->shiftnz) {
7301e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
731137fb511SHong Zhang     } else if (info->shiftpd) {
7321e2582c4SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
733137fb511SHong Zhang     }
734137fb511SHong Zhang   }
735137fb511SHong Zhang   PetscFunctionReturn(0);
736137fb511SHong Zhang }
737137fb511SHong Zhang 
7380a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7394a2ae208SSatish Balay #undef __FUNCT__
7404a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
741dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
742da3a660dSBarry Smith {
743dfbe8321SBarry Smith   PetscErrorCode ierr;
744416022c9SBarry Smith   Mat            C;
745416022c9SBarry Smith 
7463a40ed3dSBarry Smith   PetscFunctionBegin;
74743244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
7489e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
749af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
750273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
75152e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
7523a40ed3dSBarry Smith   PetscFunctionReturn(0);
753da3a660dSBarry Smith }
7540a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7554a2ae208SSatish Balay #undef __FUNCT__
7564a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
757dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
7588c37ef55SBarry Smith {
759416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
760416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
7616849ba73SBarry Smith   PetscErrorCode    ierr;
762*d0f46423SBarry Smith   PetscInt          *r,*c,i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
76338baddfdSBarry Smith   PetscInt          nz,*rout,*cout;
764dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
765d9fead3dSBarry Smith   const PetscScalar *b;
766dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
7678c37ef55SBarry Smith 
7683a40ed3dSBarry Smith   PetscFunctionBegin;
7693a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
77091bf9adeSBarry Smith 
771d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
7721ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
773416022c9SBarry Smith   tmp  = a->solve_work;
7748c37ef55SBarry Smith 
7753b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7763b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
7778c37ef55SBarry Smith 
7788c37ef55SBarry Smith   /* forward solve the lower triangular */
7798c37ef55SBarry Smith   tmp[0] = b[*r++];
780010a20caSHong Zhang   tmps   = tmp;
7818c37ef55SBarry Smith   for (i=1; i<n; i++) {
782010a20caSHong Zhang     v   = aa + ai[i] ;
783010a20caSHong Zhang     vi  = aj + ai[i] ;
78453e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
78553e7425aSBarry Smith     sum = b[*r++];
786ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
7878c37ef55SBarry Smith     tmp[i] = sum;
7888c37ef55SBarry Smith   }
7898c37ef55SBarry Smith 
7908c37ef55SBarry Smith   /* backward solve the upper triangular */
7918c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
792010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
793010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
794416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
7958c37ef55SBarry Smith     sum = tmp[i];
796ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
797010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
7988c37ef55SBarry Smith   }
7998c37ef55SBarry Smith 
8003b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8013b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
802d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
8031ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
804*d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
8053a40ed3dSBarry Smith   PetscFunctionReturn(0);
8068c37ef55SBarry Smith }
807026e39d0SSatish Balay 
8082bebee5dSHong Zhang #undef __FUNCT__
8092bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
8102bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8112bebee5dSHong Zhang {
8122bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8132bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8142bebee5dSHong Zhang   PetscErrorCode  ierr;
815*d0f46423SBarry Smith   PetscInt        *r,*c,i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
8162bebee5dSHong Zhang   PetscInt        nz,*rout,*cout,neq;
817dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
818dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
8192bebee5dSHong Zhang 
8202bebee5dSHong Zhang   PetscFunctionBegin;
8212bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8222bebee5dSHong Zhang 
8232bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8242bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8252bebee5dSHong Zhang 
8262bebee5dSHong Zhang   tmp  = a->solve_work;
8272bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8282bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8292bebee5dSHong Zhang 
8302bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
8312bebee5dSHong Zhang     /* forward solve the lower triangular */
8322bebee5dSHong Zhang     tmp[0] = b[r[0]];
8332bebee5dSHong Zhang     tmps   = tmp;
8342bebee5dSHong Zhang     for (i=1; i<n; i++) {
8352bebee5dSHong Zhang       v   = aa + ai[i] ;
8362bebee5dSHong Zhang       vi  = aj + ai[i] ;
8372bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8382bebee5dSHong Zhang       sum = b[r[i]];
8392bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8402bebee5dSHong Zhang       tmp[i] = sum;
8412bebee5dSHong Zhang     }
8422bebee5dSHong Zhang     /* backward solve the upper triangular */
8432bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
8442bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
8452bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
8462bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
8472bebee5dSHong Zhang       sum = tmp[i];
8482bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8492bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
8502bebee5dSHong Zhang     }
8512bebee5dSHong Zhang 
8522bebee5dSHong Zhang     b += n;
8532bebee5dSHong Zhang     x += n;
8542bebee5dSHong Zhang   }
8552bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8562bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8572bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
8582bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
8592bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
8602bebee5dSHong Zhang   PetscFunctionReturn(0);
8612bebee5dSHong Zhang }
8622bebee5dSHong Zhang 
863137fb511SHong Zhang #undef __FUNCT__
864137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
865137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
866137fb511SHong Zhang {
867137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
868137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
869137fb511SHong Zhang   PetscErrorCode  ierr;
870*d0f46423SBarry Smith   PetscInt        *r,*c,i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
871137fb511SHong Zhang   PetscInt        nz,*rout,*cout,row;
872dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
873dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
874137fb511SHong Zhang 
875137fb511SHong Zhang   PetscFunctionBegin;
876137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
877137fb511SHong Zhang 
878137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
879137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
880137fb511SHong Zhang   tmp  = a->solve_work;
881137fb511SHong Zhang 
882137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
883137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
884137fb511SHong Zhang 
885137fb511SHong Zhang   /* forward solve the lower triangular */
886137fb511SHong Zhang   tmp[0] = b[*r++];
887137fb511SHong Zhang   tmps   = tmp;
888137fb511SHong Zhang   for (row=1; row<n; row++) {
889137fb511SHong Zhang     i   = rout[row]; /* permuted row */
890137fb511SHong Zhang     v   = aa + ai[i] ;
891137fb511SHong Zhang     vi  = aj + ai[i] ;
892137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
893137fb511SHong Zhang     sum = b[*r++];
894137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
895137fb511SHong Zhang     tmp[row] = sum;
896137fb511SHong Zhang   }
897137fb511SHong Zhang 
898137fb511SHong Zhang   /* backward solve the upper triangular */
899137fb511SHong Zhang   for (row=n-1; row>=0; row--){
900137fb511SHong Zhang     i   = rout[row]; /* permuted row */
901137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
902137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
903137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
904137fb511SHong Zhang     sum = tmp[row];
905137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
906137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
907137fb511SHong Zhang   }
908137fb511SHong Zhang 
909137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
910137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
911137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
912137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
913*d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
914137fb511SHong Zhang   PetscFunctionReturn(0);
915137fb511SHong Zhang }
916137fb511SHong Zhang 
917930ae53cSBarry Smith /* ----------------------------------------------------------- */
9184a2ae208SSatish Balay #undef __FUNCT__
9194a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
920dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
921930ae53cSBarry Smith {
922930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9236849ba73SBarry Smith   PetscErrorCode    ierr;
924*d0f46423SBarry Smith   PetscInt          n = A->rmap->n;
925356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
926356650c2SBarry Smith   PetscScalar       *x;
927356650c2SBarry Smith   const PetscScalar *b;
928dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
929aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
930356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
931dd6ea824SBarry Smith   const MatScalar   *v;
932dd6ea824SBarry Smith   PetscScalar       sum;
933d85afc42SSatish Balay #endif
934930ae53cSBarry Smith 
9353a40ed3dSBarry Smith   PetscFunctionBegin;
9363a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
937930ae53cSBarry Smith 
938356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
940930ae53cSBarry Smith 
941aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
9421c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
9431c4feecaSSatish Balay #else
944930ae53cSBarry Smith   /* forward solve the lower triangular */
945930ae53cSBarry Smith   x[0] = b[0];
946930ae53cSBarry Smith   for (i=1; i<n; i++) {
947930ae53cSBarry Smith     ai_i = ai[i];
948930ae53cSBarry Smith     v    = aa + ai_i;
949930ae53cSBarry Smith     vi   = aj + ai_i;
950930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
951930ae53cSBarry Smith     sum  = b[i];
952930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
953930ae53cSBarry Smith     x[i] = sum;
954930ae53cSBarry Smith   }
955930ae53cSBarry Smith 
956930ae53cSBarry Smith   /* backward solve the upper triangular */
957930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
958930ae53cSBarry Smith     adiag_i = adiag[i];
959d708749eSBarry Smith     v       = aa + adiag_i + 1;
960d708749eSBarry Smith     vi      = aj + adiag_i + 1;
961930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
962930ae53cSBarry Smith     sum     = x[i];
963930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
964930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
965930ae53cSBarry Smith   }
9661c4feecaSSatish Balay #endif
967*d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr);
968356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9691ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9703a40ed3dSBarry Smith   PetscFunctionReturn(0);
971930ae53cSBarry Smith }
972930ae53cSBarry Smith 
9734a2ae208SSatish Balay #undef __FUNCT__
9744a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
975dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
976da3a660dSBarry Smith {
977416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
978416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
9796849ba73SBarry Smith   PetscErrorCode  ierr;
980*d0f46423SBarry Smith   PetscInt        *r,*c,i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
98138baddfdSBarry Smith   PetscInt        nz,*rout,*cout;
982dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
983dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
984da3a660dSBarry Smith 
9853a40ed3dSBarry Smith   PetscFunctionBegin;
98678b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
987da3a660dSBarry Smith 
9881ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
9891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
990416022c9SBarry Smith   tmp  = a->solve_work;
991da3a660dSBarry Smith 
9923b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
9933b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
994da3a660dSBarry Smith 
995da3a660dSBarry Smith   /* forward solve the lower triangular */
996da3a660dSBarry Smith   tmp[0] = b[*r++];
997da3a660dSBarry Smith   for (i=1; i<n; i++) {
998010a20caSHong Zhang     v   = aa + ai[i] ;
999010a20caSHong Zhang     vi  = aj + ai[i] ;
1000416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
1001da3a660dSBarry Smith     sum = b[*r++];
1002010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1003da3a660dSBarry Smith     tmp[i] = sum;
1004da3a660dSBarry Smith   }
1005da3a660dSBarry Smith 
1006da3a660dSBarry Smith   /* backward solve the upper triangular */
1007da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1008010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
1009010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
1010416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1011da3a660dSBarry Smith     sum = tmp[i];
1012010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1013010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1014da3a660dSBarry Smith     x[*c--] += tmp[i];
1015da3a660dSBarry Smith   }
1016da3a660dSBarry Smith 
10173b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10183b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10191ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10201ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1021efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1022898183ecSLois Curfman McInnes 
10233a40ed3dSBarry Smith   PetscFunctionReturn(0);
1024da3a660dSBarry Smith }
1025da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10264a2ae208SSatish Balay #undef __FUNCT__
10274a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1028dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1029da3a660dSBarry Smith {
1030416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10312235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10326849ba73SBarry Smith   PetscErrorCode  ierr;
1033*d0f46423SBarry Smith   PetscInt        *r,*c,i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
103438baddfdSBarry Smith   PetscInt        nz,*rout,*cout,*diag = a->diag;
1035dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1036dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1037da3a660dSBarry Smith 
10383a40ed3dSBarry Smith   PetscFunctionBegin;
10391ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1041416022c9SBarry Smith   tmp  = a->solve_work;
1042da3a660dSBarry Smith 
10432235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10442235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1045da3a660dSBarry Smith 
1046da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
10472235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1048da3a660dSBarry Smith 
1049da3a660dSBarry Smith   /* forward solve the U^T */
1050da3a660dSBarry Smith   for (i=0; i<n; i++) {
1051010a20caSHong Zhang     v   = aa + diag[i] ;
1052010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1053f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1054c38d4ed2SBarry Smith     s1  = tmp[i];
1055ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1056da3a660dSBarry Smith     while (nz--) {
1057010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1058da3a660dSBarry Smith     }
1059c38d4ed2SBarry Smith     tmp[i] = s1;
1060da3a660dSBarry Smith   }
1061da3a660dSBarry Smith 
1062da3a660dSBarry Smith   /* backward solve the L^T */
1063da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1064010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1065010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1066f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1067f1af5d2fSBarry Smith     s1  = tmp[i];
1068da3a660dSBarry Smith     while (nz--) {
1069010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1070da3a660dSBarry Smith     }
1071da3a660dSBarry Smith   }
1072da3a660dSBarry Smith 
1073da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1074da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1075da3a660dSBarry Smith 
10762235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10772235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10781ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10791ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1080da3a660dSBarry Smith 
1081*d0f46423SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap->n);CHKERRQ(ierr);
10823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1083da3a660dSBarry Smith }
1084da3a660dSBarry Smith 
10854a2ae208SSatish Balay #undef __FUNCT__
10864a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1087dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1088da3a660dSBarry Smith {
1089416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10902235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10916849ba73SBarry Smith   PetscErrorCode  ierr;
1092*d0f46423SBarry Smith   PetscInt        *r,*c,i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
109338baddfdSBarry Smith   PetscInt        nz,*rout,*cout,*diag = a->diag;
1094dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1095dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
10966abc6512SBarry Smith 
10973a40ed3dSBarry Smith   PetscFunctionBegin;
109823bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
10996abc6512SBarry Smith 
11001ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1102416022c9SBarry Smith   tmp = a->solve_work;
11036abc6512SBarry Smith 
11042235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
11052235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11066abc6512SBarry Smith 
11076abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
11082235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
11096abc6512SBarry Smith 
11106abc6512SBarry Smith   /* forward solve the U^T */
11116abc6512SBarry Smith   for (i=0; i<n; i++) {
1112010a20caSHong Zhang     v   = aa + diag[i] ;
1113010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1114f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11156abc6512SBarry Smith     tmp[i] *= *v++;
11166abc6512SBarry Smith     while (nz--) {
1117010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11186abc6512SBarry Smith     }
11196abc6512SBarry Smith   }
11206abc6512SBarry Smith 
11216abc6512SBarry Smith   /* backward solve the L^T */
11226abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1123010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1124010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1125f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11266abc6512SBarry Smith     while (nz--) {
1127010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11286abc6512SBarry Smith     }
11296abc6512SBarry Smith   }
11306abc6512SBarry Smith 
11316abc6512SBarry Smith   /* copy tmp into x according to permutation */
11326abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11336abc6512SBarry Smith 
11342235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11352235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11361ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11386abc6512SBarry Smith 
1139efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
11403a40ed3dSBarry Smith   PetscFunctionReturn(0);
1141da3a660dSBarry Smith }
11429e25ed09SBarry Smith /* ----------------------------------------------------------------*/
11433c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
1144b24902e0SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
114508480c60SBarry Smith 
11464a2ae208SSatish Balay #undef __FUNCT__
11474a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
1148dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
11499e25ed09SBarry Smith {
1150416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
11519e25ed09SBarry Smith   IS                 isicol;
11526849ba73SBarry Smith   PetscErrorCode     ierr;
1153*d0f46423SBarry Smith   PetscInt           *r,*ic,n=A->rmap->n,*ai=a->i,*aj=a->j,d;
11545a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
11555a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
11565a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
115777c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1158329f5518SBarry Smith   PetscReal          f;
11595a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
11605a8e39fbSHong Zhang   PetscBT            lnkbt;
11615a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1162a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1163a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
116409f38230SBarry Smith   PetscTruth         missing;
11659e25ed09SBarry Smith 
11663a40ed3dSBarry Smith   PetscFunctionBegin;
1167*d0f46423SBarry 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);
11686cf997b0SBarry Smith   f             = info->fill;
116938baddfdSBarry Smith   levels        = (PetscInt)info->levels;
117038baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
11714c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
117282bf6240SBarry Smith 
117308480c60SBarry Smith   /* special case that simply copies fill pattern */
1174be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1175be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
117686bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
1177b24902e0SBarry Smith     ierr = MatDuplicateNoCreate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
11785c9eb25fSBarry Smith     (*fact)->factor                 = MAT_FACTOR_LU;
1179f3a39becSBarry Smith     (*fact)->info.factor_mallocs    = 0;
1180f3a39becSBarry Smith     (*fact)->info.fill_ratio_given  = info->fill;
1181f3a39becSBarry Smith     (*fact)->info.fill_ratio_needed = 1.0;
118208480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
11838ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
118409f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
118508480c60SBarry Smith     b->row              = isrow;
118608480c60SBarry Smith     b->col              = iscol;
118782bf6240SBarry Smith     b->icol             = isicol;
1188*d0f46423SBarry Smith     ierr                = PetscMalloc(((*fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1189f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
1190c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1191c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
11923c5fc038SBarry Smith     ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
11933a40ed3dSBarry Smith     PetscFunctionReturn(0);
119408480c60SBarry Smith   }
119508480c60SBarry Smith 
1196898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1197898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
11989e25ed09SBarry Smith 
11999e25ed09SBarry Smith   /* get new row pointers */
12005a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12015a8e39fbSHong Zhang   bi[0] = 0;
12025a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
12035a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12045a8e39fbSHong Zhang   bdiag[0]  = 0;
12056cf997b0SBarry Smith 
12065a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
12075a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
12085a8e39fbSHong Zhang 
12095a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
12105a8e39fbSHong Zhang   nlnk = n + 1;
12115a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12125a8e39fbSHong Zhang 
12135a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1214a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12155a8e39fbSHong Zhang   current_space = free_space;
1216a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12175a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12185a8e39fbSHong Zhang 
12195a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12205a8e39fbSHong Zhang     nzi = 0;
12216cf997b0SBarry Smith     /* copy current row into linked list */
12225a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
12235a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
12245a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12255a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12265a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12275a8e39fbSHong Zhang     nzi += nlnk;
12286cf997b0SBarry Smith 
12296cf997b0SBarry Smith     /* make sure diagonal entry is included */
12305a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12316cf997b0SBarry Smith       fm = n;
12325a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12335a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12345a8e39fbSHong Zhang       lnk[fm]    = i;
12355a8e39fbSHong Zhang       lnk_lvl[i] = 0;
12365a8e39fbSHong Zhang       nzi++; dcount++;
12376cf997b0SBarry Smith     }
12386cf997b0SBarry Smith 
12395a8e39fbSHong Zhang     /* add pivot rows into the active row */
12405a8e39fbSHong Zhang     nzbd = 0;
12415a8e39fbSHong Zhang     prow = lnk[n];
12425a8e39fbSHong Zhang     while (prow < i) {
12435a8e39fbSHong Zhang       nnz      = bdiag[prow];
12445a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
12455a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
12465a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
12475a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
12485a8e39fbSHong Zhang       nzi += nlnk;
12495a8e39fbSHong Zhang       prow = lnk[prow];
12505a8e39fbSHong Zhang       nzbd++;
125156beaf04SBarry Smith     }
12525a8e39fbSHong Zhang     bdiag[i] = nzbd;
12535a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1254ecf371e4SBarry Smith 
12555a8e39fbSHong Zhang     /* if free space is not available, make more free space */
12565a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
12575a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1258a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1259a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
12605a8e39fbSHong Zhang       reallocs++;
12615a8e39fbSHong Zhang     }
1262ecf371e4SBarry Smith 
12635a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
12645a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
12655a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
12665a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
12675a8e39fbSHong Zhang 
12685a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
12695a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
127077431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1271d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
12726cf997b0SBarry Smith     }
12735a8e39fbSHong Zhang 
12745a8e39fbSHong Zhang     current_space->array           += nzi;
12755a8e39fbSHong Zhang     current_space->local_used      += nzi;
12765a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
12775a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
12785a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
12795a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
12809e25ed09SBarry Smith   }
12815a8e39fbSHong Zhang 
1282898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1283898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
12845a8e39fbSHong Zhang 
12855a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
12865a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1287a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
12885a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1289a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
12905a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
12919e25ed09SBarry Smith 
12926cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1293f64065bbSBarry Smith   {
12945a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1295ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1296ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1297ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1298ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1299335d9088SBarry Smith     if (diagonal_fill) {
1300ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1301335d9088SBarry Smith     }
1302f64065bbSBarry Smith   }
130363ba0a88SBarry Smith #endif
1304bbb6d6a8SBarry Smith 
13059e25ed09SBarry Smith   /* put together the new matrix */
130652e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1307416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1308e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1309e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13107c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13115a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1312b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13135a8e39fbSHong Zhang   b->j          = bj;
13145a8e39fbSHong Zhang   b->i          = bi;
13155a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13165a8e39fbSHong Zhang   b->diag       = bdiag;
1317416022c9SBarry Smith   b->ilen       = 0;
1318416022c9SBarry Smith   b->imax       = 0;
1319416022c9SBarry Smith   b->row        = isrow;
1320416022c9SBarry Smith   b->col        = iscol;
1321c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1322c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
132382bf6240SBarry Smith   b->icol       = isicol;
132487828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1325416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13265a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
13275a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13285a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
13295c9eb25fSBarry Smith   (*fact)->factor = MAT_FACTOR_LU;
1330418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1331ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
13325a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
133371c2f376SKris Buschelman 
133454e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
133571c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
13364846f1f5SKris Buschelman 
13373a40ed3dSBarry Smith   PetscFunctionReturn(0);
13389e25ed09SBarry Smith }
1339d5d45c9bSBarry Smith 
13403c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1341a6175056SHong Zhang #undef __FUNCT__
1342a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1343af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1344a6175056SHong Zhang {
13452ed133dbSHong Zhang   Mat            C = *B;
13460968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
13472ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
13489bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
13492ed133dbSHong Zhang   PetscErrorCode ierr;
1350*d0f46423SBarry Smith   PetscInt       *rip,*riip,i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol;
13512ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
13522ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1353622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1354540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1355fbf22428SSatish Balay   PetscReal      shiftpd;
1356540703acSHong Zhang   ChShift_Ctx    sctx;
1357540703acSHong Zhang   PetscInt       newshift;
1358d5d45c9bSBarry Smith 
1359a6175056SHong Zhang   PetscFunctionBegin;
1360117f1891Stmunson 
1361540703acSHong Zhang   shiftnz   = info->shiftnz;
1362540703acSHong Zhang   shiftpd   = info->shiftpd;
1363ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1364ee45ca4aSHong Zhang 
13652ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
13669bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
13672ed133dbSHong Zhang 
13682ed133dbSHong Zhang   /* initialization */
13692ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
13702ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
13712ed133dbSHong Zhang   jl   = il + mbs;
13722ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
13732ed133dbSHong Zhang 
1374540703acSHong Zhang   sctx.shift_amount = 0;
1375540703acSHong Zhang   sctx.nshift       = 0;
13762ed133dbSHong Zhang   do {
1377540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
13782ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
13792ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1380a6175056SHong Zhang     }
13812ed133dbSHong Zhang 
13822ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1383c05c3958SHong Zhang       bval = ba + bi[k];
13842ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
13852ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
13862ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
13879bfd6278SHong Zhang         col = riip[aj[j]];
13882ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
13892ed133dbSHong Zhang           rtmp[col] = aa[j];
1390c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
13912ed133dbSHong Zhang         }
13922ed133dbSHong Zhang       }
139339e3d630SHong Zhang       /* shift the diagonal of the matrix */
1394540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
13952ed133dbSHong Zhang 
13962ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
13972ed133dbSHong Zhang       dk = rtmp[k];
13982ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
13992ed133dbSHong Zhang 
14002ed133dbSHong Zhang       while (i < k){
14012ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
14022ed133dbSHong Zhang 
14032ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
14042ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
14052ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
14062ed133dbSHong Zhang         dk += uikdi*ba[ili];
14072ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14082ed133dbSHong Zhang 
14092ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14102ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14112ed133dbSHong Zhang         if (jmin < jmax){
14122ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14132ed133dbSHong Zhang           /* update il and jl for row i */
14142ed133dbSHong Zhang           il[i] = jmin;
14152ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14162ed133dbSHong Zhang         }
14172ed133dbSHong Zhang         i = nexti;
14182ed133dbSHong Zhang       }
14192ed133dbSHong Zhang 
14200697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14210697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14220697b6caSHong Zhang       rs   = 0.0;
14233c9e8598SHong Zhang       jmin = bi[k]+1;
14243c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14253c9e8598SHong Zhang       bcol = bj + jmin;
14263c9e8598SHong Zhang       while (nz--){
142789f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
142889f845c8SHong Zhang         bcol++;
14293c9e8598SHong Zhang       }
1430540703acSHong Zhang 
1431540703acSHong Zhang       sctx.rs = rs;
1432540703acSHong Zhang       sctx.pv = dk;
14335b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1434117f1891Stmunson 
1435117f1891Stmunson       if (newshift == 1) {
1436117f1891Stmunson         if (!sctx.shift_amount) {
1437117f1891Stmunson           sctx.shift_amount = 1e-5;
1438117f1891Stmunson         }
1439117f1891Stmunson         break;
1440117f1891Stmunson       }
14413c9e8598SHong Zhang 
14423c9e8598SHong Zhang       /* copy data into U(k,:) */
144339e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
144439e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
144539e3d630SHong Zhang       if (jmin < jmax) {
144639e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
144739e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
14483c9e8598SHong Zhang         }
144939e3d630SHong Zhang         /* add the k-th row into il and jl */
14503c9e8598SHong Zhang         il[k] = jmin;
14513c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
14523c9e8598SHong Zhang       }
14533c9e8598SHong Zhang     }
1454540703acSHong Zhang   } while (sctx.chshift);
14553c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
14563c9e8598SHong Zhang 
145739e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
14589bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
14595c9eb25fSBarry Smith   C->factor       = MAT_FACTOR_CHOLESKY;
14603c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
14613c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1462*d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr);
1463540703acSHong Zhang   if (sctx.nshift){
1464540703acSHong Zhang     if (shiftnz) {
14651e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1466540703acSHong Zhang     } else if (shiftpd) {
14671e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
14680697b6caSHong Zhang     }
14693c9e8598SHong Zhang   }
14703c9e8598SHong Zhang   PetscFunctionReturn(0);
14713c9e8598SHong Zhang }
1472a6175056SHong Zhang 
1473a6175056SHong Zhang #undef __FUNCT__
1474a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1475dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1476a6175056SHong Zhang {
14770968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1478ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1479dfbe8321SBarry Smith   PetscErrorCode     ierr;
148058ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
1481*d0f46423SBarry Smith   PetscInt           reallocs=0,*rip,*riip,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui;
1482ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
148358ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
14845a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1485ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1486a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1487a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
14887a48dd6fSHong Zhang   PetscBT            lnkbt;
1489b635d36bSHong Zhang   IS                 iperm;
1490a6175056SHong Zhang 
1491a6175056SHong Zhang   PetscFunctionBegin;
1492*d0f46423SBarry 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);
149358ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
149458ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1495653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1496b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14977a48dd6fSHong Zhang 
149839e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
149939e3d630SHong Zhang   ui[0] = 0;
150039e3d630SHong Zhang 
1501b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1502622e2028SHong Zhang   if (!levels && perm_identity) {
150358ebbce7SBarry Smith 
1504ed59401aSHong Zhang     for (i=0; i<am; i++) {
150539e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1506ed59401aSHong Zhang     }
150739e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
150839e3d630SHong Zhang     cols = uj;
1509ed59401aSHong Zhang     for (i=0; i<am; i++) {
1510ed59401aSHong Zhang       aj    = a->j + a->diag[i];
151139e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
151239e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1513ed59401aSHong Zhang     }
151439e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1515b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1516b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1517b635d36bSHong Zhang 
15187a48dd6fSHong Zhang     /* initialization */
15195a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15207a48dd6fSHong Zhang 
15217a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
15227a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
15231393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
15247a48dd6fSHong Zhang     il         = jl + am;
15257a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
15267a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
15277a48dd6fSHong Zhang     for (i=0; i<am; i++){
15287a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
15297a48dd6fSHong Zhang     }
15307a48dd6fSHong Zhang 
15317a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
15327a48dd6fSHong Zhang     nlnk = am + 1;
15332ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15347a48dd6fSHong Zhang 
15357a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1536a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
15377a48dd6fSHong Zhang     current_space = free_space;
1538a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
15397a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
15407a48dd6fSHong Zhang 
15417a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
15427a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
15437a48dd6fSHong Zhang       nzk   = 0;
1544622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1545b635d36bSHong Zhang       if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
1546622e2028SHong Zhang       ncols_upper = 0;
1547622e2028SHong Zhang       for (j=0; j<ncols; j++){
1548b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1549b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
15505a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1551622e2028SHong Zhang           ncols_upper++;
1552622e2028SHong Zhang         }
1553622e2028SHong Zhang       }
1554b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15557a48dd6fSHong Zhang       nzk += nlnk;
15567a48dd6fSHong Zhang 
15577a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
15587a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
15597a48dd6fSHong Zhang 
15607a48dd6fSHong Zhang       while (prow < k){
15617a48dd6fSHong Zhang         nextprow = jl[prow];
15627a48dd6fSHong Zhang 
15637a48dd6fSHong Zhang         /* merge prow into k-th row */
15647a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
15657a48dd6fSHong Zhang         jmax = ui[prow+1];
15667a48dd6fSHong Zhang         ncols = jmax-jmin;
1567ed59401aSHong Zhang         i     = jmin - ui[prow];
1568ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15695a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
15705a8e39fbSHong Zhang         j     = *(uj - 1);
15715a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
15727a48dd6fSHong Zhang         nzk += nlnk;
15737a48dd6fSHong Zhang 
15747a48dd6fSHong Zhang         /* update il and jl for prow */
15757a48dd6fSHong Zhang         if (jmin < jmax){
15767a48dd6fSHong Zhang           il[prow] = jmin;
15777a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
15787a48dd6fSHong Zhang         }
15797a48dd6fSHong Zhang         prow = nextprow;
15807a48dd6fSHong Zhang       }
15817a48dd6fSHong Zhang 
15827a48dd6fSHong Zhang       /* if free space is not available, make more free space */
15837a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
15847a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
15857a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1586a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1587a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
15887a48dd6fSHong Zhang         reallocs++;
15897a48dd6fSHong Zhang       }
15907a48dd6fSHong Zhang 
15917a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1592b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
15932ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
15947a48dd6fSHong Zhang 
15957a48dd6fSHong Zhang       /* add the k-th row into il and jl */
159639e3d630SHong Zhang       if (nzk > 1){
15977a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
15987a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
15997a48dd6fSHong Zhang         il[k] = ui[k] + 1;
16007a48dd6fSHong Zhang       }
16017a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
16027a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
16037a48dd6fSHong Zhang 
16047a48dd6fSHong Zhang       current_space->array           += nzk;
16057a48dd6fSHong Zhang       current_space->local_used      += nzk;
16067a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16077a48dd6fSHong Zhang 
16087a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16097a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16107a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16117a48dd6fSHong Zhang 
16127a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16137a48dd6fSHong Zhang     }
16147a48dd6fSHong Zhang 
16156cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16167a48dd6fSHong Zhang     if (ai[am] != 0) {
161739e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1618ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1619ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1620ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16217a48dd6fSHong Zhang     } else {
1622ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16237a48dd6fSHong Zhang     }
162463ba0a88SBarry Smith #endif
16257a48dd6fSHong Zhang 
16267a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1627b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16287a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
16295a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
16307a48dd6fSHong Zhang 
16317a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
16327a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1633a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
16342ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1635a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
16367a48dd6fSHong Zhang 
163739e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
163839e3d630SHong Zhang 
16397a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
16407a48dd6fSHong Zhang 
1641b24902e0SBarry Smith   b    = (Mat_SeqSBAIJ*)(*fact)->data;
16427a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
16437a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
16447a48dd6fSHong Zhang   b->j    = uj;
16457a48dd6fSHong Zhang   b->i    = ui;
16467a48dd6fSHong Zhang   b->diag = 0;
16477a48dd6fSHong Zhang   b->ilen = 0;
16487a48dd6fSHong Zhang   b->imax = 0;
16497a48dd6fSHong Zhang   b->row  = perm;
1650b635d36bSHong Zhang   b->col  = perm;
1651b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1652b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1653b635d36bSHong Zhang   b->icol = iperm;
16547a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
16557a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1656b24902e0SBarry Smith   ierr = PetscLogObjectMemory((*fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
16577a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1658e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1659e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
16607a48dd6fSHong Zhang 
16615c9eb25fSBarry Smith   (*fact)->factor                 = MAT_FACTOR_CHOLESKY;
1662b24902e0SBarry Smith   (*fact)->info.factor_mallocs    = reallocs;
1663b24902e0SBarry Smith   (*fact)->info.fill_ratio_given  = fill;
16647a48dd6fSHong Zhang   if (ai[am] != 0) {
1665b24902e0SBarry Smith     (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
16667a48dd6fSHong Zhang   } else {
1667b24902e0SBarry Smith     (*fact)->info.fill_ratio_needed = 0.0;
16687a48dd6fSHong Zhang   }
166939e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1670622e2028SHong Zhang   if (perm_identity){
1671b24902e0SBarry Smith     (*fact)->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1672b24902e0SBarry Smith     (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1673622e2028SHong Zhang   }
1674a6175056SHong Zhang   PetscFunctionReturn(0);
1675a6175056SHong Zhang }
1676d5d45c9bSBarry Smith 
1677f76d2b81SHong Zhang #undef __FUNCT__
1678f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1679dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1680f76d2b81SHong Zhang {
1681f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
168210c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
1683dfbe8321SBarry Smith   PetscErrorCode     ierr;
1684f76d2b81SHong Zhang   PetscTruth         perm_identity;
168510c27e3dSHong Zhang   PetscReal          fill = info->fill;
1686*d0f46423SBarry Smith   PetscInt           *rip,*riip,i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
168710c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
16882ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1689a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
169010c27e3dSHong Zhang   PetscBT            lnkbt;
16912ed133dbSHong Zhang   IS                 iperm;
1692f76d2b81SHong Zhang 
1693f76d2b81SHong Zhang   PetscFunctionBegin;
1694*d0f46423SBarry 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);
16952ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1696f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
16972ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
16982ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
16999bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
170010c27e3dSHong Zhang 
170110c27e3dSHong Zhang   /* initialization */
17021393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
170310c27e3dSHong Zhang   ui[0] = 0;
170410c27e3dSHong Zhang 
170510c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17061393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17071393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17081393bc97SHong Zhang   il     = jl + am;
17091393bc97SHong Zhang   cols   = il + am;
17101393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17111393bc97SHong Zhang   for (i=0; i<am; i++){
17121393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1713f76d2b81SHong Zhang   }
171410c27e3dSHong Zhang 
171510c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17161393bc97SHong Zhang   nlnk = am + 1;
17171393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
171810c27e3dSHong Zhang 
17191393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1720a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
172110c27e3dSHong Zhang   current_space = free_space;
172210c27e3dSHong Zhang 
17231393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
172410c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
172510c27e3dSHong Zhang     nzk   = 0;
17262ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
17279bfd6278SHong Zhang     if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
17282ed133dbSHong Zhang     ncols_upper = 0;
17292ed133dbSHong Zhang     for (j=0; j<ncols; j++){
17309bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
17312ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
17322ed133dbSHong Zhang         cols[ncols_upper] = i;
17332ed133dbSHong Zhang         ncols_upper++;
17342ed133dbSHong Zhang       }
17352ed133dbSHong Zhang     }
17361393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
173710c27e3dSHong Zhang     nzk += nlnk;
173810c27e3dSHong Zhang 
173910c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
174010c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
174110c27e3dSHong Zhang 
174210c27e3dSHong Zhang     while (prow < k){
174310c27e3dSHong Zhang       nextprow = jl[prow];
174410c27e3dSHong Zhang       /* merge prow into k-th row */
17451393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
174610c27e3dSHong Zhang       jmax = ui[prow+1];
174710c27e3dSHong Zhang       ncols = jmax-jmin;
17481393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
17495a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
175010c27e3dSHong Zhang       nzk += nlnk;
175110c27e3dSHong Zhang 
175210c27e3dSHong Zhang       /* update il and jl for prow */
175310c27e3dSHong Zhang       if (jmin < jmax){
175410c27e3dSHong Zhang         il[prow] = jmin;
17552ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
175610c27e3dSHong Zhang       }
175710c27e3dSHong Zhang       prow = nextprow;
175810c27e3dSHong Zhang     }
175910c27e3dSHong Zhang 
176010c27e3dSHong Zhang     /* if free space is not available, make more free space */
176110c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
17621393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
176310c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1764a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
176510c27e3dSHong Zhang       reallocs++;
176610c27e3dSHong Zhang     }
176710c27e3dSHong Zhang 
176810c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
17691393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
177010c27e3dSHong Zhang 
177110c27e3dSHong Zhang     /* add the k-th row into il and jl */
177210c27e3dSHong Zhang     if (nzk-1 > 0){
17731393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
177410c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
177510c27e3dSHong Zhang       il[k] = ui[k] + 1;
177610c27e3dSHong Zhang     }
17772ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
177810c27e3dSHong Zhang     current_space->array           += nzk;
177910c27e3dSHong Zhang     current_space->local_used      += nzk;
178010c27e3dSHong Zhang     current_space->local_remaining -= nzk;
178110c27e3dSHong Zhang 
178210c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
178310c27e3dSHong Zhang   }
178410c27e3dSHong Zhang 
17856cf91177SBarry Smith #if defined(PETSC_USE_INFO)
17861393bc97SHong Zhang   if (ai[am] != 0) {
17871393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1788ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1789ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1790ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
179110c27e3dSHong Zhang   } else {
1792ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
179310c27e3dSHong Zhang   }
179463ba0a88SBarry Smith #endif
179510c27e3dSHong Zhang 
179610c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
17979bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
179810c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
179910c27e3dSHong Zhang 
180010c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18011393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1802a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
180310c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
180410c27e3dSHong Zhang 
180510c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
180610c27e3dSHong Zhang 
1807b24902e0SBarry Smith   b = (Mat_SeqSBAIJ*)(*fact)->data;
180810c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1809e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1810e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18111393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
181210c27e3dSHong Zhang   b->j    = uj;
181310c27e3dSHong Zhang   b->i    = ui;
181410c27e3dSHong Zhang   b->diag = 0;
181510c27e3dSHong Zhang   b->ilen = 0;
181610c27e3dSHong Zhang   b->imax = 0;
181710c27e3dSHong Zhang   b->row  = perm;
18189bfd6278SHong Zhang   b->col  = perm;
18199bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18209bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18219bfd6278SHong Zhang   b->icol = iperm;
182210c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18231393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1824b24902e0SBarry Smith   ierr    = PetscLogObjectMemory(*fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18251393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
182610c27e3dSHong Zhang 
18275c9eb25fSBarry Smith   (*fact)->factor                 = MAT_FACTOR_CHOLESKY;
1828b24902e0SBarry Smith   (*fact)->info.factor_mallocs    = reallocs;
1829b24902e0SBarry Smith   (*fact)->info.fill_ratio_given  = fill;
18301393bc97SHong Zhang   if (ai[am] != 0) {
1831b24902e0SBarry Smith     (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
183210c27e3dSHong Zhang   } else {
1833b24902e0SBarry Smith     (*fact)->info.fill_ratio_needed = 0.0;
183410c27e3dSHong Zhang   }
183539e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
18362ed133dbSHong Zhang   if (perm_identity){
183710c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
183810c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1839fd531716SHong Zhang     (*fact)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1840fd531716SHong Zhang     (*fact)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1841fd531716SHong Zhang   } else {
1842fd531716SHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1;
1843fd531716SHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1844fd531716SHong Zhang     (*fact)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1845fd531716SHong Zhang     (*fact)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
18462ed133dbSHong Zhang   }
1847f76d2b81SHong Zhang   PetscFunctionReturn(0);
1848f76d2b81SHong Zhang }
1849