xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision c3af1dc1987bf6e39eadad1cce0fe90c4e2d5aeb)
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 
22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2386bacbd2SBarry Smith 
24e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
269cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2738baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
28e34fafa9SBarry Smith #endif
2986bacbd2SBarry Smith 
304a2ae208SSatish Balay #undef __FUNCT__
314a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3286bacbd2SBarry Smith   /* ------------------------------------------------------------
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith           This interface was contribed by Tony Caola
3586bacbd2SBarry Smith 
3686bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
375bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
385bd2ddc7SBarry Smith      SPARSEKIT2.
395bd2ddc7SBarry Smith 
405bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
415bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4286bacbd2SBarry Smith 
4386bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4486bacbd2SBarry Smith      help in getting this routine ironed out.
4586bacbd2SBarry Smith 
465bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
475bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
485bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
495bd2ddc7SBarry Smith      Yousef Saad's code.
5086bacbd2SBarry Smith 
5186bacbd2SBarry Smith      ------------------------------------------------------------
5286bacbd2SBarry Smith */
537529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5486bacbd2SBarry Smith {
5560ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5660ec86bdSBarry Smith   PetscFunctionBegin;
57e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5860ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5960ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
6060ec86bdSBarry Smith #else
6186bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6207d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6386bacbd2SBarry Smith   PetscTruth     reorder;
649cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
65899cda47SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->rmap.n;
6638baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6738baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6838baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6987828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
7054a8161fSBarry Smith   PetscReal      af;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   PetscFunctionBegin;
7386bacbd2SBarry Smith 
7486bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7538baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7686bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7733259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7838baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7938baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
8086bacbd2SBarry Smith 
8186bacbd2SBarry Smith 
8286bacbd2SBarry Smith   /* ------------------------------------------------------------
8386bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8486bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8586bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8686bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8786bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8886bacbd2SBarry Smith      ------------------------------------------------------------
8986bacbd2SBarry Smith   */
9086bacbd2SBarry Smith 
9186bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9286bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9386bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9486bacbd2SBarry Smith   reorder = PetscNot(reorder);
9586bacbd2SBarry Smith 
9686bacbd2SBarry Smith 
9786bacbd2SBarry Smith   /* storage for ilu factor */
9838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9938baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
10087828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
10138baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10286bacbd2SBarry Smith 
10386bacbd2SBarry Smith   /* ------------------------------------------------------------
10486bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10586bacbd2SBarry Smith      ------------------------------------------------------------
10686bacbd2SBarry Smith   */
107b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
108b19713cbSBarry Smith     old_j[i]++;
10986bacbd2SBarry Smith   }
110b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
111b19713cbSBarry Smith     old_i[i]++;
112b19713cbSBarry Smith   };
113010a20caSHong Zhang 
11486bacbd2SBarry Smith 
11586bacbd2SBarry Smith   if (reorder) {
116c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
117c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
118c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
119c0c2c81eSBarry Smith       r[i]  = r[i]+1;
120c0c2c81eSBarry Smith       c[i]  = c[i]+1;
121c0c2c81eSBarry Smith     }
12238baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12338baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12487828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1255bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
126c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
127c0c2c81eSBarry Smith       r[i]  = r[i]-1;
128c0c2c81eSBarry Smith       c[i]  = c[i]-1;
129c0c2c81eSBarry Smith     }
130c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
131c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
132b19713cbSBarry Smith     o_a = old_a2;
133b19713cbSBarry Smith     o_j = old_j2;
134b19713cbSBarry Smith     o_i = old_i2;
135b19713cbSBarry Smith   } else {
136b19713cbSBarry Smith     o_a = old_a;
137b19713cbSBarry Smith     o_j = old_j;
138b19713cbSBarry Smith     o_i = old_i;
13986bacbd2SBarry Smith   }
14086bacbd2SBarry Smith 
14186bacbd2SBarry Smith   /* ------------------------------------------------------------
14286bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14386bacbd2SBarry Smith      ------------------------------------------------------------
14486bacbd2SBarry Smith   */
14586bacbd2SBarry Smith 
14638baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14738baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14986bacbd2SBarry Smith 
15054a8161fSBarry 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);
1516849ba73SBarry Smith   if (sierr) {
1526849ba73SBarry Smith     switch (sierr) {
153a83599f4SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
154a83599f4SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
155e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
156e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15777431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15877431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15986bacbd2SBarry Smith     }
16086bacbd2SBarry Smith   }
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16386bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16486bacbd2SBarry Smith 
16586bacbd2SBarry Smith   /* ------------------------------------------------------------
16686bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16786bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16886bacbd2SBarry Smith      ------------------------------------------------------------
16986bacbd2SBarry Smith   */
17086bacbd2SBarry Smith 
17187828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
17238baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17386bacbd2SBarry Smith 
1745bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17786bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17886bacbd2SBarry Smith 
17986bacbd2SBarry Smith   if (reorder) {
18086bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
18186bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18286bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
183313ae024SBarry Smith   } else {
184b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
185f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
186b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
187b19713cbSBarry Smith     }
188b19713cbSBarry Smith   }
189b19713cbSBarry Smith 
190b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
19186bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
192b19713cbSBarry Smith     old_i[i]--;
193b19713cbSBarry Smith   }
194b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
195b19713cbSBarry Smith     old_j[i]--;
196b19713cbSBarry Smith   }
19786bacbd2SBarry Smith 
198b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
200b19713cbSBarry Smith     new_i[i]--;
201b19713cbSBarry Smith   }
202b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
203b19713cbSBarry Smith     new_j[i]--;
204b19713cbSBarry Smith   }
20586bacbd2SBarry Smith 
20686bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20786bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2084c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2094c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
210c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21186bacbd2SBarry Smith   for(i=0; i<n; i++) {
21286bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21386bacbd2SBarry Smith   };
21486bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21507d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21686bacbd2SBarry Smith 
217c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
218c0c2c81eSBarry Smith 
219329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
22007d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
22186bacbd2SBarry Smith 
22286bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22386bacbd2SBarry Smith 
224f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
225f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
226f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
227ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
22886bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22986bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
23086bacbd2SBarry Smith 
23186bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
232a96a251dSBarry Smith   b->freedata      = PETSC_TRUE;
23386bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23407d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23586bacbd2SBarry Smith   b->a             = new_a;
23686bacbd2SBarry Smith   b->j             = new_j;
23786bacbd2SBarry Smith   b->i             = new_i;
23886bacbd2SBarry Smith   b->ilen          = 0;
23986bacbd2SBarry Smith   b->imax          = 0;
2401f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
241313ae024SBarry Smith   b->row           = isirow;
24286bacbd2SBarry Smith   b->col           = iscolf;
24387828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24486bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24586bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24686bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24786bacbd2SBarry Smith 
24886bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24986bacbd2SBarry Smith 
250431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
251ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
252ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
253ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
254ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
255431e710aSBarry Smith 
2567529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25771c2f376SKris Buschelman 
25886bacbd2SBarry Smith   PetscFunctionReturn(0);
25960ec86bdSBarry Smith #endif
26086bacbd2SBarry Smith }
26186bacbd2SBarry Smith 
2624a2ae208SSatish Balay #undef __FUNCT__
263b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
264dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
265289bc588SBarry Smith {
266416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
267289bc588SBarry Smith   IS                 isicol;
2686849ba73SBarry Smith   PetscErrorCode     ierr;
269899cda47SBarry Smith   PetscInt           *r,*ic,i,n=A->rmap.n,*ai=a->i,*aj=a->j;
2701393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
2711393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2729e046878SBarry Smith   PetscReal          f;
2734875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
274a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2751393bc97SHong Zhang   PetscBT            lnkbt;
276289bc588SBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
278899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2794c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
280ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
281ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
282289bc588SBarry Smith 
283289bc588SBarry Smith   /* get new row pointers */
2841393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2851393bc97SHong Zhang   bi[0] = 0;
2861393bc97SHong Zhang 
2871393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2881393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2891393bc97SHong Zhang   bdiag[0] = 0;
2901393bc97SHong Zhang 
2911393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2921393bc97SHong Zhang   nlnk = n + 1;
2931393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2941393bc97SHong Zhang 
2954875ba38SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&im);CHKERRQ(ierr);
2961393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2971393bc97SHong Zhang 
2981393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
299d3d32019SBarry Smith   f = info->fill;
300a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3011393bc97SHong Zhang   current_space = free_space;
302289bc588SBarry Smith 
303289bc588SBarry Smith   for (i=0; i<n; i++) {
3041393bc97SHong Zhang     /* copy previous fill into linked list */
305289bc588SBarry Smith     nzi = 0;
3061393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3071393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3081393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
309afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3101393bc97SHong Zhang     nzi += nlnk;
3111393bc97SHong Zhang 
3121393bc97SHong Zhang     /* add pivot rows into linked list */
3131393bc97SHong Zhang     row = lnk[n];
3141393bc97SHong Zhang     while (row < i) {
315d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
316d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
317d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3181393bc97SHong Zhang       nzi += nlnk;
3191393bc97SHong Zhang       row  = lnk[row];
320289bc588SBarry Smith     }
3211393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3221393bc97SHong Zhang     im[i]   = nzi;
3231393bc97SHong Zhang 
3241393bc97SHong Zhang     /* mark bdiag */
3251393bc97SHong Zhang     nzbd = 0;
3261393bc97SHong Zhang     nnz  = nzi;
3271393bc97SHong Zhang     k    = lnk[n];
3281393bc97SHong Zhang     while (nnz-- && k < i){
3291393bc97SHong Zhang       nzbd++;
3301393bc97SHong Zhang       k = lnk[k];
3311393bc97SHong Zhang     }
3321393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3331393bc97SHong Zhang 
3341393bc97SHong Zhang     /* if free space is not available, make more free space */
3351393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3361393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
337a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3381393bc97SHong Zhang       reallocs++;
3391393bc97SHong Zhang     }
3401393bc97SHong Zhang 
3411393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3421393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3431393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3441393bc97SHong Zhang     current_space->array           += nzi;
3451393bc97SHong Zhang     current_space->local_used      += nzi;
3461393bc97SHong Zhang     current_space->local_remaining -= nzi;
347289bc588SBarry Smith   }
3486cf91177SBarry Smith #if defined(PETSC_USE_INFO)
349464e8d44SSatish Balay   if (ai[n] != 0) {
3501393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
351ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
352ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
353ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
354ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
35505bf559cSSatish Balay   } else {
356ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
35705bf559cSSatish Balay   }
35863ba0a88SBarry Smith #endif
3592fb3ed76SBarry Smith 
360898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
361898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3621987afe7SBarry Smith 
3631393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3641393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
365a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3661393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3674875ba38SHong Zhang   ierr = PetscFree(im);CHKERRQ(ierr);
368289bc588SBarry Smith 
369289bc588SBarry Smith   /* put together the new matrix */
370f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,B);CHKERRQ(ierr);
371f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
372f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
373ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
37452e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
375416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
376a96a251dSBarry Smith   b->freedata     = PETSC_TRUE;
3777c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3781393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3791393bc97SHong Zhang   b->j          = bj;
3801393bc97SHong Zhang   b->i          = bi;
3811393bc97SHong Zhang   b->diag       = bdiag;
382416022c9SBarry Smith   b->ilen       = 0;
383416022c9SBarry Smith   b->imax       = 0;
384416022c9SBarry Smith   b->row        = isrow;
385416022c9SBarry Smith   b->col        = iscol;
386c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
387c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
38882bf6240SBarry Smith   b->icol       = isicol;
38987828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3901393bc97SHong Zhang 
3911393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3921393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3931393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
3947fd98bd6SLois Curfman McInnes 
395329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
396418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
397ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
398703deb49SSatish Balay 
399e93ccc4dSBarry Smith   if (ai[n] != 0) {
4001393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
401eea2913aSSatish Balay   } else {
402eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
403eea2913aSSatish Balay   }
4044846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
40571c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4063a40ed3dSBarry Smith   PetscFunctionReturn(0);
407289bc588SBarry Smith }
4081393bc97SHong Zhang 
4090a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4104a2ae208SSatish Balay #undef __FUNCT__
4114a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
412af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
413289bc588SBarry Smith {
414416022c9SBarry Smith   Mat            C=*B;
415416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
41682bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4176849ba73SBarry Smith   PetscErrorCode ierr;
418899cda47SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap.n,*bi=b->i,*bj=b->j;
419b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
420d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
42187828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4226a7c8fc2SHong Zhang   PetscScalar    d;
4236a7c8fc2SHong Zhang   PetscReal      rs;
424b3bf805bSHong Zhang   LUShift_Ctx    sctx;
425d2276718SHong Zhang   PetscInt       newshift;
426289bc588SBarry Smith 
4273a40ed3dSBarry Smith   PetscFunctionBegin;
42878b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
42978b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
43087828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
43187828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
432010a20caSHong Zhang   rtmps = rtmp; ics = ic;
433289bc588SBarry Smith 
434ada7143aSSatish Balay   sctx.shift_top  = 0;
435ada7143aSSatish Balay   sctx.nshift_max = 0;
436ada7143aSSatish Balay   sctx.shift_lo   = 0;
437ada7143aSSatish Balay   sctx.shift_hi   = 0;
438ada7143aSSatish Balay 
4396cc28720Svictorle   if (!a->diag) {
4400cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4410cf777b8SBarry Smith   }
442118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
443118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4441a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
44538baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
446118f5789SBarry Smith     sctx.shift_top = 0;
4476cc28720Svictorle     for (i=0; i<n; i++) {
4489f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4499f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4501a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
451010a20caSHong Zhang       v  = a->a+aai[i];
452e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
453e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4541a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4551a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4566cc28720Svictorle     }
457*c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
458118f5789SBarry Smith     sctx.shift_top    *= 1.1;
459d2276718SHong Zhang     sctx.nshift_max   = 5;
460d2276718SHong Zhang     sctx.shift_lo     = 0.;
461d2276718SHong Zhang     sctx.shift_hi     = 1.;
462a0a356efSHong Zhang   }
463a0a356efSHong Zhang 
464a0a356efSHong Zhang   sctx.shift_amount = 0;
465a0a356efSHong Zhang   sctx.nshift       = 0;
466085256b3SBarry Smith   do {
467d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
468289bc588SBarry Smith     for (i=0; i<n; i++){
469b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
470b3bf805bSHong Zhang       bjtmp = bj + bi[i];
471b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
472289bc588SBarry Smith 
473289bc588SBarry Smith       /* load in initial (unfactored row) */
474416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
475b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
476010a20caSHong Zhang       v     = a->a + a->i[r[i]];
477085256b3SBarry Smith       for (j=0; j<nz; j++) {
478b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
479085256b3SBarry Smith       }
480d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
481289bc588SBarry Smith 
482b3bf805bSHong Zhang       row = *bjtmp++;
483f2582111SSatish Balay       while  (row < i) {
4848c37ef55SBarry Smith         pc = rtmp + row;
4858c37ef55SBarry Smith         if (*pc != 0.0) {
486010a20caSHong Zhang           pv         = b->a + diag_offset[row];
487010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
48835aab85fSBarry Smith           multiplier = *pc / *pv++;
4898c37ef55SBarry Smith           *pc        = multiplier;
490b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
491f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
492efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
493289bc588SBarry Smith         }
494b3bf805bSHong Zhang         row = *bjtmp++;
495289bc588SBarry Smith       }
496416022c9SBarry Smith       /* finished row so stick it into b->a */
497b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
498b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
499b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
500b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
501d3d32019SBarry Smith       rs   = 0.0;
502d3d32019SBarry Smith       for (j=0; j<nz; j++) {
503d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
504d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
505d3d32019SBarry Smith       }
506d2276718SHong Zhang 
507b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
508d2276718SHong Zhang       sctx.rs  = rs;
509d2276718SHong Zhang       sctx.pv  = pv[diag];
5103e8c821fSHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
511d2276718SHong Zhang       if (newshift == 1){
512b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
513d2276718SHong Zhang       } else if (newshift == -1){
514a83599f4SBarry Smith         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %G tolerance %G * rs %G",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
515d708749eSBarry Smith       }
51635aab85fSBarry Smith     }
517d2276718SHong Zhang 
518118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5196cc28720Svictorle       /*
5206c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5216cc28720Svictorle        * then try lower shift
5226cc28720Svictorle        */
523d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
524d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
525d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
526d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
527d2276718SHong Zhang       sctx.nshift++;
5286cc28720Svictorle     }
529d2276718SHong Zhang   } while (sctx.lushift);
5300f11f9aeSSatish Balay 
53135aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
53235aab85fSBarry Smith   for (i=0; i<n; i++) {
533010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
53435aab85fSBarry Smith   }
53535aab85fSBarry Smith 
536606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
53778b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
53878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
539416022c9SBarry Smith   C->factor = FACTOR_LU;
5407dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
541c456f294SBarry Smith   C->assembled = PETSC_TRUE;
542899cda47SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
543d2276718SHong Zhang   if (sctx.nshift){
544118f5789SBarry Smith     if (info->shiftnz) {
545ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
546118f5789SBarry Smith     } else if (info->shiftpd) {
547ae15b995SBarry Smith       ierr = PetscInfo4(0,"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);
5486cc28720Svictorle     }
5490697b6caSHong Zhang   }
5503a40ed3dSBarry Smith   PetscFunctionReturn(0);
551289bc588SBarry Smith }
5520889b5dcSvictorle 
5530889b5dcSvictorle #undef __FUNCT__
5540889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
555dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5560889b5dcSvictorle {
5570889b5dcSvictorle   PetscFunctionBegin;
5580889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5590889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5600889b5dcSvictorle   PetscFunctionReturn(0);
5610889b5dcSvictorle }
5620889b5dcSvictorle 
5630889b5dcSvictorle 
5640a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5654a2ae208SSatish Balay #undef __FUNCT__
5664a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
567dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
568da3a660dSBarry Smith {
569dfbe8321SBarry Smith   PetscErrorCode ierr;
570416022c9SBarry Smith   Mat            C;
571416022c9SBarry Smith 
5723a40ed3dSBarry Smith   PetscFunctionBegin;
5739e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
574af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
575273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
57652e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5773a40ed3dSBarry Smith   PetscFunctionReturn(0);
578da3a660dSBarry Smith }
5790a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
582dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5838c37ef55SBarry Smith {
584416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
585416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
5866849ba73SBarry Smith   PetscErrorCode ierr;
587899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
58838baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
58987828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
5908c37ef55SBarry Smith 
5913a40ed3dSBarry Smith   PetscFunctionBegin;
5923a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
59391bf9adeSBarry Smith 
5941ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
5951ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
596416022c9SBarry Smith   tmp  = a->solve_work;
5978c37ef55SBarry Smith 
5983b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
5993b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6008c37ef55SBarry Smith 
6018c37ef55SBarry Smith   /* forward solve the lower triangular */
6028c37ef55SBarry Smith   tmp[0] = b[*r++];
603010a20caSHong Zhang   tmps   = tmp;
6048c37ef55SBarry Smith   for (i=1; i<n; i++) {
605010a20caSHong Zhang     v   = aa + ai[i] ;
606010a20caSHong Zhang     vi  = aj + ai[i] ;
60753e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
60853e7425aSBarry Smith     sum = b[*r++];
609ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6108c37ef55SBarry Smith     tmp[i] = sum;
6118c37ef55SBarry Smith   }
6128c37ef55SBarry Smith 
6138c37ef55SBarry Smith   /* backward solve the upper triangular */
6148c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
615010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
616010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
617416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6188c37ef55SBarry Smith     sum = tmp[i];
619ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
620010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6218c37ef55SBarry Smith   }
6228c37ef55SBarry Smith 
6233b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6243b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6251ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6261ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
627899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
6283a40ed3dSBarry Smith   PetscFunctionReturn(0);
6298c37ef55SBarry Smith }
630026e39d0SSatish Balay 
6312bebee5dSHong Zhang #undef __FUNCT__
6322bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
6332bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
6342bebee5dSHong Zhang {
6352bebee5dSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6362bebee5dSHong Zhang   IS             iscol = a->col,isrow = a->row;
6372bebee5dSHong Zhang   PetscErrorCode ierr;
638899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
6392bebee5dSHong Zhang   PetscInt       nz,*rout,*cout,neq;
6402bebee5dSHong Zhang   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6412bebee5dSHong Zhang 
6422bebee5dSHong Zhang   PetscFunctionBegin;
6432bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
6442bebee5dSHong Zhang 
6452bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
6462bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
6472bebee5dSHong Zhang 
6482bebee5dSHong Zhang   tmp  = a->solve_work;
6492bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6502bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
6512bebee5dSHong Zhang 
6522bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
6532bebee5dSHong Zhang     /* forward solve the lower triangular */
6542bebee5dSHong Zhang     tmp[0] = b[r[0]];
6552bebee5dSHong Zhang     tmps   = tmp;
6562bebee5dSHong Zhang     for (i=1; i<n; i++) {
6572bebee5dSHong Zhang       v   = aa + ai[i] ;
6582bebee5dSHong Zhang       vi  = aj + ai[i] ;
6592bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
6602bebee5dSHong Zhang       sum = b[r[i]];
6612bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6622bebee5dSHong Zhang       tmp[i] = sum;
6632bebee5dSHong Zhang     }
6642bebee5dSHong Zhang     /* backward solve the upper triangular */
6652bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
6662bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
6672bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
6682bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
6692bebee5dSHong Zhang       sum = tmp[i];
6702bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6712bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
6722bebee5dSHong Zhang     }
6732bebee5dSHong Zhang 
6742bebee5dSHong Zhang     b += n;
6752bebee5dSHong Zhang     x += n;
6762bebee5dSHong Zhang   }
6772bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6782bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6792bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
6802bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
6812bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
6822bebee5dSHong Zhang   PetscFunctionReturn(0);
6832bebee5dSHong Zhang }
6842bebee5dSHong Zhang 
685930ae53cSBarry Smith /* ----------------------------------------------------------- */
6864a2ae208SSatish Balay #undef __FUNCT__
6874a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
688dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
689930ae53cSBarry Smith {
690930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6916849ba73SBarry Smith   PetscErrorCode ierr;
692899cda47SBarry Smith   PetscInt       n = A->rmap.n,*ai = a->i,*aj = a->j,*adiag = a->diag;
693362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
694aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
69538baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
696362ced78SSatish Balay   PetscScalar    *v,sum;
697d85afc42SSatish Balay #endif
698930ae53cSBarry Smith 
6993a40ed3dSBarry Smith   PetscFunctionBegin;
7003a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
701930ae53cSBarry Smith 
7021ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7031ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
704930ae53cSBarry Smith 
705aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
7061c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
7071c4feecaSSatish Balay #else
708930ae53cSBarry Smith   /* forward solve the lower triangular */
709930ae53cSBarry Smith   x[0] = b[0];
710930ae53cSBarry Smith   for (i=1; i<n; i++) {
711930ae53cSBarry Smith     ai_i = ai[i];
712930ae53cSBarry Smith     v    = aa + ai_i;
713930ae53cSBarry Smith     vi   = aj + ai_i;
714930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
715930ae53cSBarry Smith     sum  = b[i];
716930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
717930ae53cSBarry Smith     x[i] = sum;
718930ae53cSBarry Smith   }
719930ae53cSBarry Smith 
720930ae53cSBarry Smith   /* backward solve the upper triangular */
721930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
722930ae53cSBarry Smith     adiag_i = adiag[i];
723d708749eSBarry Smith     v       = aa + adiag_i + 1;
724d708749eSBarry Smith     vi      = aj + adiag_i + 1;
725930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
726930ae53cSBarry Smith     sum     = x[i];
727930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
728930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
729930ae53cSBarry Smith   }
7301c4feecaSSatish Balay #endif
731899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
7321ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7331ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
735930ae53cSBarry Smith }
736930ae53cSBarry Smith 
7374a2ae208SSatish Balay #undef __FUNCT__
7384a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
739dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
740da3a660dSBarry Smith {
741416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
742416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7436849ba73SBarry Smith   PetscErrorCode ierr;
744899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
74538baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
74687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
747da3a660dSBarry Smith 
7483a40ed3dSBarry Smith   PetscFunctionBegin;
74978b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
750da3a660dSBarry Smith 
7511ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
753416022c9SBarry Smith   tmp  = a->solve_work;
754da3a660dSBarry Smith 
7553b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7563b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
757da3a660dSBarry Smith 
758da3a660dSBarry Smith   /* forward solve the lower triangular */
759da3a660dSBarry Smith   tmp[0] = b[*r++];
760da3a660dSBarry Smith   for (i=1; i<n; i++) {
761010a20caSHong Zhang     v   = aa + ai[i] ;
762010a20caSHong Zhang     vi  = aj + ai[i] ;
763416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
764da3a660dSBarry Smith     sum = b[*r++];
765010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
766da3a660dSBarry Smith     tmp[i] = sum;
767da3a660dSBarry Smith   }
768da3a660dSBarry Smith 
769da3a660dSBarry Smith   /* backward solve the upper triangular */
770da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
771010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
772010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
773416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
774da3a660dSBarry Smith     sum = tmp[i];
775010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
776010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
777da3a660dSBarry Smith     x[*c--] += tmp[i];
778da3a660dSBarry Smith   }
779da3a660dSBarry Smith 
7803b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7813b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7821ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7831ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
784efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
785898183ecSLois Curfman McInnes 
7863a40ed3dSBarry Smith   PetscFunctionReturn(0);
787da3a660dSBarry Smith }
788da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7894a2ae208SSatish Balay #undef __FUNCT__
7904a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
791dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
792da3a660dSBarry Smith {
793416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7942235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7956849ba73SBarry Smith   PetscErrorCode ierr;
796899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
79738baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
79887828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
799da3a660dSBarry Smith 
8003a40ed3dSBarry Smith   PetscFunctionBegin;
8011ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8021ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
803416022c9SBarry Smith   tmp  = a->solve_work;
804da3a660dSBarry Smith 
8052235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8062235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
807da3a660dSBarry Smith 
808da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
8092235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
810da3a660dSBarry Smith 
811da3a660dSBarry Smith   /* forward solve the U^T */
812da3a660dSBarry Smith   for (i=0; i<n; i++) {
813010a20caSHong Zhang     v   = aa + diag[i] ;
814010a20caSHong Zhang     vi  = aj + diag[i] + 1;
815f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
816c38d4ed2SBarry Smith     s1  = tmp[i];
817ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
818da3a660dSBarry Smith     while (nz--) {
819010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
820da3a660dSBarry Smith     }
821c38d4ed2SBarry Smith     tmp[i] = s1;
822da3a660dSBarry Smith   }
823da3a660dSBarry Smith 
824da3a660dSBarry Smith   /* backward solve the L^T */
825da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
826010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
827010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
828f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
829f1af5d2fSBarry Smith     s1  = tmp[i];
830da3a660dSBarry Smith     while (nz--) {
831010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
832da3a660dSBarry Smith     }
833da3a660dSBarry Smith   }
834da3a660dSBarry Smith 
835da3a660dSBarry Smith   /* copy tmp into x according to permutation */
836da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
837da3a660dSBarry Smith 
8382235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8392235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8401ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8411ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
842da3a660dSBarry Smith 
843899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap.n);CHKERRQ(ierr);
8443a40ed3dSBarry Smith   PetscFunctionReturn(0);
845da3a660dSBarry Smith }
846da3a660dSBarry Smith 
8474a2ae208SSatish Balay #undef __FUNCT__
8484a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
849dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
850da3a660dSBarry Smith {
851416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8522235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8536849ba73SBarry Smith   PetscErrorCode ierr;
854899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
85538baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
85687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8576abc6512SBarry Smith 
8583a40ed3dSBarry Smith   PetscFunctionBegin;
85923bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8606abc6512SBarry Smith 
8611ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8621ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
863416022c9SBarry Smith   tmp = a->solve_work;
8646abc6512SBarry Smith 
8652235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8662235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8676abc6512SBarry Smith 
8686abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8692235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8706abc6512SBarry Smith 
8716abc6512SBarry Smith   /* forward solve the U^T */
8726abc6512SBarry Smith   for (i=0; i<n; i++) {
873010a20caSHong Zhang     v   = aa + diag[i] ;
874010a20caSHong Zhang     vi  = aj + diag[i] + 1;
875f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8766abc6512SBarry Smith     tmp[i] *= *v++;
8776abc6512SBarry Smith     while (nz--) {
878010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8796abc6512SBarry Smith     }
8806abc6512SBarry Smith   }
8816abc6512SBarry Smith 
8826abc6512SBarry Smith   /* backward solve the L^T */
8836abc6512SBarry Smith   for (i=n-1; i>=0; i--){
884010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
885010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
886f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8876abc6512SBarry Smith     while (nz--) {
888010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8896abc6512SBarry Smith     }
8906abc6512SBarry Smith   }
8916abc6512SBarry Smith 
8926abc6512SBarry Smith   /* copy tmp into x according to permutation */
8936abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8946abc6512SBarry Smith 
8952235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8962235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8971ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8981ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8996abc6512SBarry Smith 
900efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
9013a40ed3dSBarry Smith   PetscFunctionReturn(0);
902da3a660dSBarry Smith }
9039e25ed09SBarry Smith /* ----------------------------------------------------------------*/
904dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
90508480c60SBarry Smith 
9064a2ae208SSatish Balay #undef __FUNCT__
9074a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
908dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
9099e25ed09SBarry Smith {
910416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
9119e25ed09SBarry Smith   IS                 isicol;
9126849ba73SBarry Smith   PetscErrorCode     ierr;
913899cda47SBarry Smith   PetscInt           *r,*ic,n=A->rmap.n,*ai=a->i,*aj=a->j;
9145a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
9155a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
9165a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
91777c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
918329f5518SBarry Smith   PetscReal          f;
9195a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
9205a8e39fbSHong Zhang   PetscBT            lnkbt;
9215a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
922a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
923a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
9249e25ed09SBarry Smith 
9253a40ed3dSBarry Smith   PetscFunctionBegin;
9266cf997b0SBarry Smith   f             = info->fill;
92738baddfdSBarry Smith   levels        = (PetscInt)info->levels;
92838baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9294c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
93082bf6240SBarry Smith 
93108480c60SBarry Smith   /* special case that simply copies fill pattern */
932be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
933be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
93486bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9352e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
93608480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
93708480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
93808480c60SBarry Smith     if (!b->diag) {
9397c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
94008480c60SBarry Smith     }
9417c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
94208480c60SBarry Smith     b->row              = isrow;
94308480c60SBarry Smith     b->col              = iscol;
94482bf6240SBarry Smith     b->icol             = isicol;
945899cda47SBarry Smith     ierr                = PetscMalloc(((*fact)->rmap.n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
946f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
947c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
948c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9493a40ed3dSBarry Smith     PetscFunctionReturn(0);
95008480c60SBarry Smith   }
95108480c60SBarry Smith 
952898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
953898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9549e25ed09SBarry Smith 
9559e25ed09SBarry Smith   /* get new row pointers */
9565a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9575a8e39fbSHong Zhang   bi[0] = 0;
9585a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9595a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9605a8e39fbSHong Zhang   bdiag[0]  = 0;
9616cf997b0SBarry Smith 
9625a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9635a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9645a8e39fbSHong Zhang 
9655a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9665a8e39fbSHong Zhang   nlnk = n + 1;
9675a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9685a8e39fbSHong Zhang 
9695a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
970a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9715a8e39fbSHong Zhang   current_space = free_space;
972a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9735a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
9745a8e39fbSHong Zhang 
9755a8e39fbSHong Zhang   for (i=0; i<n; i++) {
9765a8e39fbSHong Zhang     nzi = 0;
9776cf997b0SBarry Smith     /* copy current row into linked list */
9785a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
9795a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
9805a8e39fbSHong Zhang     cols = aj + ai[r[i]];
9815a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
9825a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9835a8e39fbSHong Zhang     nzi += nlnk;
9846cf997b0SBarry Smith 
9856cf997b0SBarry Smith     /* make sure diagonal entry is included */
9865a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
9876cf997b0SBarry Smith       fm = n;
9885a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
9895a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
9905a8e39fbSHong Zhang       lnk[fm]    = i;
9915a8e39fbSHong Zhang       lnk_lvl[i] = 0;
9925a8e39fbSHong Zhang       nzi++; dcount++;
9936cf997b0SBarry Smith     }
9946cf997b0SBarry Smith 
9955a8e39fbSHong Zhang     /* add pivot rows into the active row */
9965a8e39fbSHong Zhang     nzbd = 0;
9975a8e39fbSHong Zhang     prow = lnk[n];
9985a8e39fbSHong Zhang     while (prow < i) {
9995a8e39fbSHong Zhang       nnz      = bdiag[prow];
10005a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
10015a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
10025a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
10035a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
10045a8e39fbSHong Zhang       nzi += nlnk;
10055a8e39fbSHong Zhang       prow = lnk[prow];
10065a8e39fbSHong Zhang       nzbd++;
100756beaf04SBarry Smith     }
10085a8e39fbSHong Zhang     bdiag[i] = nzbd;
10095a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1010ecf371e4SBarry Smith 
10115a8e39fbSHong Zhang     /* if free space is not available, make more free space */
10125a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
10135a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1014a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1015a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
10165a8e39fbSHong Zhang       reallocs++;
10175a8e39fbSHong Zhang     }
1018ecf371e4SBarry Smith 
10195a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
10205a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
10215a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
10225a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
10235a8e39fbSHong Zhang 
10245a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
10255a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
102677431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10275a8e39fbSHong Zhang     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i);
10286cf997b0SBarry Smith     }
10295a8e39fbSHong Zhang 
10305a8e39fbSHong Zhang     current_space->array           += nzi;
10315a8e39fbSHong Zhang     current_space->local_used      += nzi;
10325a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
10335a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
10345a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
10355a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
10369e25ed09SBarry Smith   }
10375a8e39fbSHong Zhang 
1038898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1039898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
10405a8e39fbSHong Zhang 
10415a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
10425a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1043a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
10445a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1045a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
10465a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10479e25ed09SBarry Smith 
10486cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1049f64065bbSBarry Smith   {
10505a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1051ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1052ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1053ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1054ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1055335d9088SBarry Smith     if (diagonal_fill) {
1056ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1057335d9088SBarry Smith     }
1058f64065bbSBarry Smith   }
105963ba0a88SBarry Smith #endif
1060bbb6d6a8SBarry Smith 
10619e25ed09SBarry Smith   /* put together the new matrix */
1062f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
1063f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
1064f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1065ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
106652e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1067416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1068a96a251dSBarry Smith   b->freedata     = PETSC_TRUE;
10697c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10705a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1071b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10725a8e39fbSHong Zhang   b->j          = bj;
10735a8e39fbSHong Zhang   b->i          = bi;
10745a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
10755a8e39fbSHong Zhang   b->diag       = bdiag;
1076416022c9SBarry Smith   b->ilen       = 0;
1077416022c9SBarry Smith   b->imax       = 0;
1078416022c9SBarry Smith   b->row        = isrow;
1079416022c9SBarry Smith   b->col        = iscol;
1080c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1081c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
108282bf6240SBarry Smith   b->icol       = isicol;
108387828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1084416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
10855a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
10865a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
10875a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
10889e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1089418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1090ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
10915a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
109271c2f376SKris Buschelman 
109354e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
109471c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
10954846f1f5SKris Buschelman 
10963a40ed3dSBarry Smith   PetscFunctionReturn(0);
10979e25ed09SBarry Smith }
1098d5d45c9bSBarry Smith 
10993c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1100a6175056SHong Zhang #undef __FUNCT__
1101a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1102af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1103a6175056SHong Zhang {
11042ed133dbSHong Zhang   Mat            C = *B;
11050968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11062ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11072ed133dbSHong Zhang   IS             ip=b->row;
11082ed133dbSHong Zhang   PetscErrorCode ierr;
1109899cda47SBarry Smith   PetscInt       *rip,i,j,mbs=A->rmap.n,*bi=b->i,*bj=b->j,*bcol;
11102ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11112ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1112622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1113540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1114fbf22428SSatish Balay   PetscReal      shiftpd;
1115540703acSHong Zhang   ChShift_Ctx    sctx;
1116540703acSHong Zhang   PetscInt       newshift;
1117d5d45c9bSBarry Smith 
1118a6175056SHong Zhang   PetscFunctionBegin;
1119540703acSHong Zhang   shiftnz   = info->shiftnz;
1120540703acSHong Zhang   shiftpd   = info->shiftpd;
1121ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1122ee45ca4aSHong Zhang 
11232ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11242ed133dbSHong Zhang 
11252ed133dbSHong Zhang   /* initialization */
11262ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11272ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11282ed133dbSHong Zhang   jl   = il + mbs;
11292ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11302ed133dbSHong Zhang 
1131540703acSHong Zhang   sctx.shift_amount = 0;
1132540703acSHong Zhang   sctx.nshift       = 0;
11332ed133dbSHong Zhang   do {
1134540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
11352ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11362ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1137a6175056SHong Zhang     }
11382ed133dbSHong Zhang 
11392ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1140c05c3958SHong Zhang       bval = ba + bi[k];
11412ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11422ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11432ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11442ed133dbSHong Zhang         col = rip[aj[j]];
11452ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11462ed133dbSHong Zhang           rtmp[col] = aa[j];
1147c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11482ed133dbSHong Zhang         }
11492ed133dbSHong Zhang       }
115039e3d630SHong Zhang       /* shift the diagonal of the matrix */
1151540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11522ed133dbSHong Zhang 
11532ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11542ed133dbSHong Zhang       dk = rtmp[k];
11552ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11562ed133dbSHong Zhang 
11572ed133dbSHong Zhang       while (i < k){
11582ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11592ed133dbSHong Zhang 
11602ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11612ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11622ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11632ed133dbSHong Zhang         dk += uikdi*ba[ili];
11642ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11652ed133dbSHong Zhang 
11662ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11672ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11682ed133dbSHong Zhang         if (jmin < jmax){
11692ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11702ed133dbSHong Zhang           /* update il and jl for row i */
11712ed133dbSHong Zhang           il[i] = jmin;
11722ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11732ed133dbSHong Zhang         }
11742ed133dbSHong Zhang         i = nexti;
11752ed133dbSHong Zhang       }
11762ed133dbSHong Zhang 
11770697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11780697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11790697b6caSHong Zhang       rs   = 0.0;
11803c9e8598SHong Zhang       jmin = bi[k]+1;
11813c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11823c9e8598SHong Zhang       if (nz){
11833c9e8598SHong Zhang         bcol = bj + jmin;
11843c9e8598SHong Zhang         while (nz--){
118589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
118689f845c8SHong Zhang           bcol++;
11873c9e8598SHong Zhang         }
11883c9e8598SHong Zhang       }
1189540703acSHong Zhang 
1190540703acSHong Zhang       sctx.rs = rs;
1191540703acSHong Zhang       sctx.pv = dk;
11923e8c821fSHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1193540703acSHong Zhang       if (newshift == 1){
1194540703acSHong Zhang         break;    /* sctx.shift_amount is updated */
1195540703acSHong Zhang       } else if (newshift == -1){
1196a83599f4SBarry Smith         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %G tolerance %G * rs %G",k,PetscAbsScalar(dk),zeropivot,rs);
11973c9e8598SHong Zhang       }
11983c9e8598SHong Zhang 
11993c9e8598SHong Zhang       /* copy data into U(k,:) */
120039e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
120139e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
120239e3d630SHong Zhang       if (jmin < jmax) {
120339e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
120439e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12053c9e8598SHong Zhang         }
120639e3d630SHong Zhang         /* add the k-th row into il and jl */
12073c9e8598SHong Zhang         il[k] = jmin;
12083c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12093c9e8598SHong Zhang       }
12103c9e8598SHong Zhang     }
1211540703acSHong Zhang   } while (sctx.chshift);
12123c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12133c9e8598SHong Zhang 
121439e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12153c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12163c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12173c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1218899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.n);CHKERRQ(ierr);
1219540703acSHong Zhang   if (sctx.nshift){
1220540703acSHong Zhang     if (shiftnz) {
1221ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1222540703acSHong Zhang     } else if (shiftpd) {
1223ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
12240697b6caSHong Zhang     }
12253c9e8598SHong Zhang   }
12263c9e8598SHong Zhang   PetscFunctionReturn(0);
12273c9e8598SHong Zhang }
1228a6175056SHong Zhang 
1229a6175056SHong Zhang #undef __FUNCT__
1230a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1231dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1232a6175056SHong Zhang {
12330968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1234ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1235ed59401aSHong Zhang   Mat                B;
1236dfbe8321SBarry Smith   PetscErrorCode     ierr;
1237653a6975SHong Zhang   PetscTruth         perm_identity;
1238899cda47SBarry Smith   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->rmap.n,*ui;
1239ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1240391eac55SHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12415a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1242ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1243a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1244a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12457a48dd6fSHong Zhang   PetscBT            lnkbt;
1246a6175056SHong Zhang 
1247a6175056SHong Zhang   PetscFunctionBegin;
1248653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12497a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12507a48dd6fSHong Zhang 
125139e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
125239e3d630SHong Zhang   ui[0] = 0;
125339e3d630SHong Zhang 
1254622e2028SHong Zhang   /* special case that simply copies fill pattern */
1255622e2028SHong Zhang   if (!levels && perm_identity) {
1256622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1257ed59401aSHong Zhang     for (i=0; i<am; i++) {
125839e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1259ed59401aSHong Zhang     }
126039e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
126139e3d630SHong Zhang     cols = uj;
1262ed59401aSHong Zhang     for (i=0; i<am; i++) {
1263ed59401aSHong Zhang       aj    = a->j + a->diag[i];
126439e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
126539e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1266ed59401aSHong Zhang     }
126739e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12687a48dd6fSHong Zhang     /* initialization */
12695a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12707a48dd6fSHong Zhang 
12717a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12727a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12731393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12747a48dd6fSHong Zhang     il         = jl + am;
12757a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12767a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12777a48dd6fSHong Zhang     for (i=0; i<am; i++){
12787a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12797a48dd6fSHong Zhang     }
12807a48dd6fSHong Zhang 
12817a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12827a48dd6fSHong Zhang     nlnk = am + 1;
12832ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12847a48dd6fSHong Zhang 
12857a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1286a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12877a48dd6fSHong Zhang     current_space = free_space;
1288a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12897a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12907a48dd6fSHong Zhang 
12917a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12927a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12937a48dd6fSHong Zhang       nzk   = 0;
1294622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1295622e2028SHong Zhang       ncols_upper = 0;
1296622e2028SHong Zhang       for (j=0; j<ncols; j++){
12975a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
12985a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
12995a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1300622e2028SHong Zhang           ncols_upper++;
1301622e2028SHong Zhang         }
1302622e2028SHong Zhang       }
13035a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13047a48dd6fSHong Zhang       nzk += nlnk;
13057a48dd6fSHong Zhang 
13067a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13077a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13087a48dd6fSHong Zhang 
13097a48dd6fSHong Zhang       while (prow < k){
13107a48dd6fSHong Zhang         nextprow = jl[prow];
13117a48dd6fSHong Zhang 
13127a48dd6fSHong Zhang         /* merge prow into k-th row */
13137a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13147a48dd6fSHong Zhang         jmax = ui[prow+1];
13157a48dd6fSHong Zhang         ncols = jmax-jmin;
1316ed59401aSHong Zhang         i     = jmin - ui[prow];
1317ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
13185a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
13195a8e39fbSHong Zhang         j     = *(uj - 1);
13205a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
13217a48dd6fSHong Zhang         nzk += nlnk;
13227a48dd6fSHong Zhang 
13237a48dd6fSHong Zhang         /* update il and jl for prow */
13247a48dd6fSHong Zhang         if (jmin < jmax){
13257a48dd6fSHong Zhang           il[prow] = jmin;
13267a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13277a48dd6fSHong Zhang         }
13287a48dd6fSHong Zhang         prow = nextprow;
13297a48dd6fSHong Zhang       }
13307a48dd6fSHong Zhang 
13317a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13327a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13337a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13347a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1335a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1336a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
13377a48dd6fSHong Zhang         reallocs++;
13387a48dd6fSHong Zhang       }
13397a48dd6fSHong Zhang 
13407a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13412ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13427a48dd6fSHong Zhang 
13437a48dd6fSHong Zhang       /* add the k-th row into il and jl */
134439e3d630SHong Zhang       if (nzk > 1){
13457a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13467a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13477a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13487a48dd6fSHong Zhang       }
13497a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13507a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13517a48dd6fSHong Zhang 
13527a48dd6fSHong Zhang       current_space->array           += nzk;
13537a48dd6fSHong Zhang       current_space->local_used      += nzk;
13547a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13557a48dd6fSHong Zhang 
13567a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13577a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13587a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13597a48dd6fSHong Zhang 
13607a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13617a48dd6fSHong Zhang     }
13627a48dd6fSHong Zhang 
13636cf91177SBarry Smith #if defined(PETSC_USE_INFO)
13647a48dd6fSHong Zhang     if (ai[am] != 0) {
136539e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1366ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1367ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1368ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
13697a48dd6fSHong Zhang     } else {
1370ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
13717a48dd6fSHong Zhang     }
137263ba0a88SBarry Smith #endif
13737a48dd6fSHong Zhang 
13747a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13757a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13765a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13777a48dd6fSHong Zhang 
13787a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13797a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1380a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13812ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1382a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13837a48dd6fSHong Zhang 
138439e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
138539e3d630SHong Zhang 
13867a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1387f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1388f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
1389ed59401aSHong Zhang   B = *fact;
1390ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1391ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13927a48dd6fSHong Zhang 
1393ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13947a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13957a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13967a48dd6fSHong Zhang   b->j    = uj;
13977a48dd6fSHong Zhang   b->i    = ui;
13987a48dd6fSHong Zhang   b->diag = 0;
13997a48dd6fSHong Zhang   b->ilen = 0;
14007a48dd6fSHong Zhang   b->imax = 0;
14017a48dd6fSHong Zhang   b->row  = perm;
14027a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
14037a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14047a48dd6fSHong Zhang   b->icol = perm;
14057a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14067a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
140752e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
14087a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
14097a48dd6fSHong Zhang 
1410ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1411ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1412ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14137a48dd6fSHong Zhang   if (ai[am] != 0) {
1414ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14157a48dd6fSHong Zhang   } else {
1416ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14177a48dd6fSHong Zhang   }
141839e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1419622e2028SHong Zhang   if (perm_identity){
1420ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1421ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1422622e2028SHong Zhang   }
1423a6175056SHong Zhang   PetscFunctionReturn(0);
1424a6175056SHong Zhang }
1425d5d45c9bSBarry Smith 
1426f76d2b81SHong Zhang #undef __FUNCT__
1427f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1428dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1429f76d2b81SHong Zhang {
1430f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
143110c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
14322ed133dbSHong Zhang   Mat                B;
1433dfbe8321SBarry Smith   PetscErrorCode     ierr;
1434f76d2b81SHong Zhang   PetscTruth         perm_identity;
143510c27e3dSHong Zhang   PetscReal          fill = info->fill;
1436899cda47SBarry Smith   PetscInt           *rip,*riip,i,am=A->rmap.n,*ai=a->i,*aj=a->j,reallocs=0,prow;
143710c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14382ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1439a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
144010c27e3dSHong Zhang   PetscBT            lnkbt;
14412ed133dbSHong Zhang   IS                 iperm;
1442f76d2b81SHong Zhang 
1443f76d2b81SHong Zhang   PetscFunctionBegin;
14442ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1445f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14462ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14472ed133dbSHong Zhang 
1448f76d2b81SHong Zhang   if (!perm_identity){
14492ed133dbSHong Zhang     /* check if perm is symmetric! */
14502ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14512ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14521393bc97SHong Zhang     for (i=0; i<am; i++) {
14532ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14542ed133dbSHong Zhang     }
14552ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14562ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1457f76d2b81SHong Zhang   }
145810c27e3dSHong Zhang 
145910c27e3dSHong Zhang   /* initialization */
14601393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
146110c27e3dSHong Zhang   ui[0] = 0;
146210c27e3dSHong Zhang 
146310c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14641393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14651393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14661393bc97SHong Zhang   il     = jl + am;
14671393bc97SHong Zhang   cols   = il + am;
14681393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14691393bc97SHong Zhang   for (i=0; i<am; i++){
14701393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1471f76d2b81SHong Zhang   }
147210c27e3dSHong Zhang 
147310c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14741393bc97SHong Zhang   nlnk = am + 1;
14751393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
147610c27e3dSHong Zhang 
14771393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1478a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
147910c27e3dSHong Zhang   current_space = free_space;
148010c27e3dSHong Zhang 
14811393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
148210c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
148310c27e3dSHong Zhang     nzk   = 0;
14842ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14852ed133dbSHong Zhang     ncols_upper = 0;
14862ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1487622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14882ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14892ed133dbSHong Zhang         cols[ncols_upper] = i;
14902ed133dbSHong Zhang         ncols_upper++;
14912ed133dbSHong Zhang       }
14922ed133dbSHong Zhang     }
14931393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
149410c27e3dSHong Zhang     nzk += nlnk;
149510c27e3dSHong Zhang 
149610c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
149710c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
149810c27e3dSHong Zhang 
149910c27e3dSHong Zhang     while (prow < k){
150010c27e3dSHong Zhang       nextprow = jl[prow];
150110c27e3dSHong Zhang       /* merge prow into k-th row */
15021393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
150310c27e3dSHong Zhang       jmax = ui[prow+1];
150410c27e3dSHong Zhang       ncols = jmax-jmin;
15051393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15065a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150710c27e3dSHong Zhang       nzk += nlnk;
150810c27e3dSHong Zhang 
150910c27e3dSHong Zhang       /* update il and jl for prow */
151010c27e3dSHong Zhang       if (jmin < jmax){
151110c27e3dSHong Zhang         il[prow] = jmin;
15122ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
151310c27e3dSHong Zhang       }
151410c27e3dSHong Zhang       prow = nextprow;
151510c27e3dSHong Zhang     }
151610c27e3dSHong Zhang 
151710c27e3dSHong Zhang     /* if free space is not available, make more free space */
151810c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15191393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
152010c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1521a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
152210c27e3dSHong Zhang       reallocs++;
152310c27e3dSHong Zhang     }
152410c27e3dSHong Zhang 
152510c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15261393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
152710c27e3dSHong Zhang 
152810c27e3dSHong Zhang     /* add the k-th row into il and jl */
152910c27e3dSHong Zhang     if (nzk-1 > 0){
15301393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
153110c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
153210c27e3dSHong Zhang       il[k] = ui[k] + 1;
153310c27e3dSHong Zhang     }
15342ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
153510c27e3dSHong Zhang     current_space->array           += nzk;
153610c27e3dSHong Zhang     current_space->local_used      += nzk;
153710c27e3dSHong Zhang     current_space->local_remaining -= nzk;
153810c27e3dSHong Zhang 
153910c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
154010c27e3dSHong Zhang   }
154110c27e3dSHong Zhang 
15426cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15431393bc97SHong Zhang   if (ai[am] != 0) {
15441393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1545ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1546ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1547ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
154810c27e3dSHong Zhang   } else {
1549ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
155010c27e3dSHong Zhang   }
155163ba0a88SBarry Smith #endif
155210c27e3dSHong Zhang 
155310c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
155410c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
155510c27e3dSHong Zhang 
155610c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15571393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1558a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
155910c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
156010c27e3dSHong Zhang 
156110c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1562f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1563f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
15642ed133dbSHong Zhang   B    = *fact;
15652ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1566ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
156710c27e3dSHong Zhang 
15682ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
156910c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
15701393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
157110c27e3dSHong Zhang   b->j    = uj;
157210c27e3dSHong Zhang   b->i    = ui;
157310c27e3dSHong Zhang   b->diag = 0;
157410c27e3dSHong Zhang   b->ilen = 0;
157510c27e3dSHong Zhang   b->imax = 0;
157610c27e3dSHong Zhang   b->row  = perm;
157710c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
157810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
157910c27e3dSHong Zhang   b->icol = perm;
158010c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15811393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15821393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
15831393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
158410c27e3dSHong Zhang 
15852ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15862ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15872ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15881393bc97SHong Zhang   if (ai[am] != 0) {
15891393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
159010c27e3dSHong Zhang   } else {
15912ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
159210c27e3dSHong Zhang   }
159339e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15942ed133dbSHong Zhang   if (perm_identity){
159510c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
159610c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15972ed133dbSHong Zhang   }
1598f76d2b81SHong Zhang   PetscFunctionReturn(0);
1599f76d2b81SHong Zhang }
1600