xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 2bebee5db09b2907f7718f484d1164e7f8ef8d84)
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;
659cc1f4e3SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->m;
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) {
15377431f27SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
15477431f27SBarry 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;
25163ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af));CHKERRQ(ierr);
25263ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af));CHKERRQ(ierr);
25363ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af));CHKERRQ(ierr);
25463ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatILUDTFactor_SeqAIJ: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;
26938baddfdSBarry Smith   PetscInt           *r,*ic,i,n=A->m,*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;
27829bbc08cSBarry Smith   if (A->M != A->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   }
34863ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
349464e8d44SSatish Balay   if (ai[n] != 0) {
3501393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
35163ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr);
3520c451bc4SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %G or use \n",af));CHKERRQ(ierr);
35363ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af));CHKERRQ(ierr);
35463ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr);
35505bf559cSSatish Balay   } else {
35663ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatLUFactorSymbolic_SeqAIJ: 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;
418b3bf805bSHong Zhang   PetscInt       *r,*ic,i,j,n=A->m,*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 
4346cc28720Svictorle   if (!a->diag) {
4350cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4360cf777b8SBarry Smith   }
437118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
438118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4391a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
44038baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
441118f5789SBarry Smith     sctx.shift_top = 0;
4426cc28720Svictorle     for (i=0; i<n; i++) {
4439f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4449f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4451a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
446010a20caSHong Zhang       v  = a->a+aai[i];
447e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
448e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4491a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4501a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4516cc28720Svictorle     }
452118f5789SBarry Smith     if (sctx.shift_top == 0.0) sctx.shift_top += 1.e-12;
453118f5789SBarry Smith     sctx.shift_top    *= 1.1;
454d2276718SHong Zhang     sctx.nshift_max   = 5;
455d2276718SHong Zhang     sctx.shift_lo     = 0.;
456d2276718SHong Zhang     sctx.shift_hi     = 1.;
457a0a356efSHong Zhang   }
458a0a356efSHong Zhang 
459a0a356efSHong Zhang   sctx.shift_amount = 0;
460a0a356efSHong Zhang   sctx.nshift       = 0;
461085256b3SBarry Smith   do {
462d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
463289bc588SBarry Smith     for (i=0; i<n; i++){
464b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
465b3bf805bSHong Zhang       bjtmp = bj + bi[i];
466b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
467289bc588SBarry Smith 
468289bc588SBarry Smith       /* load in initial (unfactored row) */
469416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
470b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
471010a20caSHong Zhang       v     = a->a + a->i[r[i]];
472085256b3SBarry Smith       for (j=0; j<nz; j++) {
473b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
474085256b3SBarry Smith       }
475d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
476289bc588SBarry Smith 
477b3bf805bSHong Zhang       row = *bjtmp++;
478f2582111SSatish Balay       while  (row < i) {
4798c37ef55SBarry Smith         pc = rtmp + row;
4808c37ef55SBarry Smith         if (*pc != 0.0) {
481010a20caSHong Zhang           pv         = b->a + diag_offset[row];
482010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
48335aab85fSBarry Smith           multiplier = *pc / *pv++;
4848c37ef55SBarry Smith           *pc        = multiplier;
485b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
486f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
487efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
488289bc588SBarry Smith         }
489b3bf805bSHong Zhang         row = *bjtmp++;
490289bc588SBarry Smith       }
491416022c9SBarry Smith       /* finished row so stick it into b->a */
492b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
493b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
494b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
495b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
496d3d32019SBarry Smith       rs   = 0.0;
497d3d32019SBarry Smith       for (j=0; j<nz; j++) {
498d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
499d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
500d3d32019SBarry Smith       }
501d2276718SHong Zhang 
502b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
503d2276718SHong Zhang       sctx.rs  = rs;
504d2276718SHong Zhang       sctx.pv  = pv[diag];
5053e8c821fSHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
506d2276718SHong Zhang       if (newshift == 1){
507b3bf805bSHong Zhang         break;    /* sctx.shift_amount is updated */
508d2276718SHong Zhang       } else if (newshift == -1){
509d2276718SHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(sctx.pv),info->zeropivot,rs);
510d708749eSBarry Smith       }
51135aab85fSBarry Smith     }
512d2276718SHong Zhang 
513118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5146cc28720Svictorle       /*
5156c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5166cc28720Svictorle        * then try lower shift
5176cc28720Svictorle        */
518d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
519d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
520d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
521d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
522d2276718SHong Zhang       sctx.nshift++;
5236cc28720Svictorle     }
524d2276718SHong Zhang   } while (sctx.lushift);
5250f11f9aeSSatish Balay 
52635aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
52735aab85fSBarry Smith   for (i=0; i<n; i++) {
528010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
52935aab85fSBarry Smith   }
53035aab85fSBarry Smith 
531606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
53278b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
53378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
534416022c9SBarry Smith   C->factor = FACTOR_LU;
5357dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
536c456f294SBarry Smith   C->assembled = PETSC_TRUE;
537efee365bSSatish Balay   ierr = PetscLogFlops(C->n);CHKERRQ(ierr);
538d2276718SHong Zhang   if (sctx.nshift){
539118f5789SBarry Smith     if (info->shiftnz) {
54063ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
541118f5789SBarry Smith     } else if (info->shiftpd) {
54263ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatLUFactorNumeric_SeqAIJ: 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);
5436cc28720Svictorle     }
5440697b6caSHong Zhang   }
5453a40ed3dSBarry Smith   PetscFunctionReturn(0);
546289bc588SBarry Smith }
5470889b5dcSvictorle 
5480889b5dcSvictorle #undef __FUNCT__
5490889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
550dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5510889b5dcSvictorle {
5520889b5dcSvictorle   PetscFunctionBegin;
5530889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5540889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5550889b5dcSvictorle   PetscFunctionReturn(0);
5560889b5dcSvictorle }
5570889b5dcSvictorle 
5580889b5dcSvictorle 
5590a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5604a2ae208SSatish Balay #undef __FUNCT__
5614a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
562dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
563da3a660dSBarry Smith {
564dfbe8321SBarry Smith   PetscErrorCode ierr;
565416022c9SBarry Smith   Mat            C;
566416022c9SBarry Smith 
5673a40ed3dSBarry Smith   PetscFunctionBegin;
5689e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
569af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
570273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
57152e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5723a40ed3dSBarry Smith   PetscFunctionReturn(0);
573da3a660dSBarry Smith }
5740a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5754a2ae208SSatish Balay #undef __FUNCT__
5764a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
577dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5788c37ef55SBarry Smith {
579416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
580416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
5816849ba73SBarry Smith   PetscErrorCode ierr;
58238baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
58338baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
58487828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
5858c37ef55SBarry Smith 
5863a40ed3dSBarry Smith   PetscFunctionBegin;
5873a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
58891bf9adeSBarry Smith 
5891ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
5901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
591416022c9SBarry Smith   tmp  = a->solve_work;
5928c37ef55SBarry Smith 
5933b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
5943b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
5958c37ef55SBarry Smith 
5968c37ef55SBarry Smith   /* forward solve the lower triangular */
5978c37ef55SBarry Smith   tmp[0] = b[*r++];
598010a20caSHong Zhang   tmps   = tmp;
5998c37ef55SBarry Smith   for (i=1; i<n; i++) {
600010a20caSHong Zhang     v   = aa + ai[i] ;
601010a20caSHong Zhang     vi  = aj + ai[i] ;
60253e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
60353e7425aSBarry Smith     sum = b[*r++];
604ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6058c37ef55SBarry Smith     tmp[i] = sum;
6068c37ef55SBarry Smith   }
6078c37ef55SBarry Smith 
6088c37ef55SBarry Smith   /* backward solve the upper triangular */
6098c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
610010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
611010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
612416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6138c37ef55SBarry Smith     sum = tmp[i];
614ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
615010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6168c37ef55SBarry Smith   }
6178c37ef55SBarry Smith 
6183b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6193b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6201ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6211ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
622efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
6233a40ed3dSBarry Smith   PetscFunctionReturn(0);
6248c37ef55SBarry Smith }
625026e39d0SSatish Balay 
626*2bebee5dSHong Zhang #undef __FUNCT__
627*2bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
628*2bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
629*2bebee5dSHong Zhang {
630*2bebee5dSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
631*2bebee5dSHong Zhang   IS             iscol = a->col,isrow = a->row;
632*2bebee5dSHong Zhang   PetscErrorCode ierr;
633*2bebee5dSHong Zhang   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
634*2bebee5dSHong Zhang   PetscInt       nz,*rout,*cout,neq;
635*2bebee5dSHong Zhang   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
636*2bebee5dSHong Zhang 
637*2bebee5dSHong Zhang   PetscFunctionBegin;
638*2bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
639*2bebee5dSHong Zhang 
640*2bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
641*2bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
642*2bebee5dSHong Zhang 
643*2bebee5dSHong Zhang   tmp  = a->solve_work;
644*2bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
645*2bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
646*2bebee5dSHong Zhang 
647*2bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
648*2bebee5dSHong Zhang     /* forward solve the lower triangular */
649*2bebee5dSHong Zhang     tmp[0] = b[r[0]];
650*2bebee5dSHong Zhang     tmps   = tmp;
651*2bebee5dSHong Zhang     for (i=1; i<n; i++) {
652*2bebee5dSHong Zhang       v   = aa + ai[i] ;
653*2bebee5dSHong Zhang       vi  = aj + ai[i] ;
654*2bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
655*2bebee5dSHong Zhang       sum = b[r[i]];
656*2bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
657*2bebee5dSHong Zhang       tmp[i] = sum;
658*2bebee5dSHong Zhang     }
659*2bebee5dSHong Zhang     /* backward solve the upper triangular */
660*2bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
661*2bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
662*2bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
663*2bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
664*2bebee5dSHong Zhang       sum = tmp[i];
665*2bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
666*2bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
667*2bebee5dSHong Zhang     }
668*2bebee5dSHong Zhang 
669*2bebee5dSHong Zhang     b += n;
670*2bebee5dSHong Zhang     x += n;
671*2bebee5dSHong Zhang   }
672*2bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
673*2bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
674*2bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
675*2bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
676*2bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
677*2bebee5dSHong Zhang   PetscFunctionReturn(0);
678*2bebee5dSHong Zhang }
679*2bebee5dSHong Zhang 
680930ae53cSBarry Smith /* ----------------------------------------------------------- */
6814a2ae208SSatish Balay #undef __FUNCT__
6824a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
683dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
684930ae53cSBarry Smith {
685930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6866849ba73SBarry Smith   PetscErrorCode ierr;
68738baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
688362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
689aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
69038baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
691362ced78SSatish Balay   PetscScalar    *v,sum;
692d85afc42SSatish Balay #endif
693930ae53cSBarry Smith 
6943a40ed3dSBarry Smith   PetscFunctionBegin;
6953a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
696930ae53cSBarry Smith 
6971ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
699930ae53cSBarry Smith 
700aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
7011c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
7021c4feecaSSatish Balay #else
703930ae53cSBarry Smith   /* forward solve the lower triangular */
704930ae53cSBarry Smith   x[0] = b[0];
705930ae53cSBarry Smith   for (i=1; i<n; i++) {
706930ae53cSBarry Smith     ai_i = ai[i];
707930ae53cSBarry Smith     v    = aa + ai_i;
708930ae53cSBarry Smith     vi   = aj + ai_i;
709930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
710930ae53cSBarry Smith     sum  = b[i];
711930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
712930ae53cSBarry Smith     x[i] = sum;
713930ae53cSBarry Smith   }
714930ae53cSBarry Smith 
715930ae53cSBarry Smith   /* backward solve the upper triangular */
716930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
717930ae53cSBarry Smith     adiag_i = adiag[i];
718d708749eSBarry Smith     v       = aa + adiag_i + 1;
719d708749eSBarry Smith     vi      = aj + adiag_i + 1;
720930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
721930ae53cSBarry Smith     sum     = x[i];
722930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
723930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
724930ae53cSBarry Smith   }
7251c4feecaSSatish Balay #endif
726efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - A->n);CHKERRQ(ierr);
7271ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7281ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
730930ae53cSBarry Smith }
731930ae53cSBarry Smith 
7324a2ae208SSatish Balay #undef __FUNCT__
7334a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
734dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
735da3a660dSBarry Smith {
736416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
737416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7386849ba73SBarry Smith   PetscErrorCode ierr;
73938baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
74038baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
74187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
742da3a660dSBarry Smith 
7433a40ed3dSBarry Smith   PetscFunctionBegin;
74478b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
745da3a660dSBarry Smith 
7461ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7471ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
748416022c9SBarry Smith   tmp  = a->solve_work;
749da3a660dSBarry Smith 
7503b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7513b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
752da3a660dSBarry Smith 
753da3a660dSBarry Smith   /* forward solve the lower triangular */
754da3a660dSBarry Smith   tmp[0] = b[*r++];
755da3a660dSBarry Smith   for (i=1; i<n; i++) {
756010a20caSHong Zhang     v   = aa + ai[i] ;
757010a20caSHong Zhang     vi  = aj + ai[i] ;
758416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
759da3a660dSBarry Smith     sum = b[*r++];
760010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
761da3a660dSBarry Smith     tmp[i] = sum;
762da3a660dSBarry Smith   }
763da3a660dSBarry Smith 
764da3a660dSBarry Smith   /* backward solve the upper triangular */
765da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
766010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
767010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
768416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
769da3a660dSBarry Smith     sum = tmp[i];
770010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
771010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
772da3a660dSBarry Smith     x[*c--] += tmp[i];
773da3a660dSBarry Smith   }
774da3a660dSBarry Smith 
7753b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7763b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7771ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7781ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
779efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
780898183ecSLois Curfman McInnes 
7813a40ed3dSBarry Smith   PetscFunctionReturn(0);
782da3a660dSBarry Smith }
783da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7844a2ae208SSatish Balay #undef __FUNCT__
7854a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
786dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
787da3a660dSBarry Smith {
788416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7892235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7906849ba73SBarry Smith   PetscErrorCode ierr;
79138baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
79238baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
79387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
794da3a660dSBarry Smith 
7953a40ed3dSBarry Smith   PetscFunctionBegin;
7961ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7971ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
798416022c9SBarry Smith   tmp  = a->solve_work;
799da3a660dSBarry Smith 
8002235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8012235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
802da3a660dSBarry Smith 
803da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
8042235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
805da3a660dSBarry Smith 
806da3a660dSBarry Smith   /* forward solve the U^T */
807da3a660dSBarry Smith   for (i=0; i<n; i++) {
808010a20caSHong Zhang     v   = aa + diag[i] ;
809010a20caSHong Zhang     vi  = aj + diag[i] + 1;
810f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
811c38d4ed2SBarry Smith     s1  = tmp[i];
812ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
813da3a660dSBarry Smith     while (nz--) {
814010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
815da3a660dSBarry Smith     }
816c38d4ed2SBarry Smith     tmp[i] = s1;
817da3a660dSBarry Smith   }
818da3a660dSBarry Smith 
819da3a660dSBarry Smith   /* backward solve the L^T */
820da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
821010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
822010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
823f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
824f1af5d2fSBarry Smith     s1  = tmp[i];
825da3a660dSBarry Smith     while (nz--) {
826010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
827da3a660dSBarry Smith     }
828da3a660dSBarry Smith   }
829da3a660dSBarry Smith 
830da3a660dSBarry Smith   /* copy tmp into x according to permutation */
831da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
832da3a660dSBarry Smith 
8332235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8342235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8351ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8361ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
837da3a660dSBarry Smith 
838efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz-A->n);CHKERRQ(ierr);
8393a40ed3dSBarry Smith   PetscFunctionReturn(0);
840da3a660dSBarry Smith }
841da3a660dSBarry Smith 
8424a2ae208SSatish Balay #undef __FUNCT__
8434a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
844dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
845da3a660dSBarry Smith {
846416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8472235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8486849ba73SBarry Smith   PetscErrorCode ierr;
84938baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
85038baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
85187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8526abc6512SBarry Smith 
8533a40ed3dSBarry Smith   PetscFunctionBegin;
85423bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8556abc6512SBarry Smith 
8561ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8571ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
858416022c9SBarry Smith   tmp = a->solve_work;
8596abc6512SBarry Smith 
8602235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8612235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8626abc6512SBarry Smith 
8636abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8642235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8656abc6512SBarry Smith 
8666abc6512SBarry Smith   /* forward solve the U^T */
8676abc6512SBarry Smith   for (i=0; i<n; i++) {
868010a20caSHong Zhang     v   = aa + diag[i] ;
869010a20caSHong Zhang     vi  = aj + diag[i] + 1;
870f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8716abc6512SBarry Smith     tmp[i] *= *v++;
8726abc6512SBarry Smith     while (nz--) {
873010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8746abc6512SBarry Smith     }
8756abc6512SBarry Smith   }
8766abc6512SBarry Smith 
8776abc6512SBarry Smith   /* backward solve the L^T */
8786abc6512SBarry Smith   for (i=n-1; i>=0; i--){
879010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
880010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
881f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8826abc6512SBarry Smith     while (nz--) {
883010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8846abc6512SBarry Smith     }
8856abc6512SBarry Smith   }
8866abc6512SBarry Smith 
8876abc6512SBarry Smith   /* copy tmp into x according to permutation */
8886abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8896abc6512SBarry Smith 
8902235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8912235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8921ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8931ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8946abc6512SBarry Smith 
895efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
8963a40ed3dSBarry Smith   PetscFunctionReturn(0);
897da3a660dSBarry Smith }
8989e25ed09SBarry Smith /* ----------------------------------------------------------------*/
899dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
90008480c60SBarry Smith 
9014a2ae208SSatish Balay #undef __FUNCT__
9024a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
903dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
9049e25ed09SBarry Smith {
905416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
9069e25ed09SBarry Smith   IS                 isicol;
9076849ba73SBarry Smith   PetscErrorCode     ierr;
9085a8e39fbSHong Zhang   PetscInt           *r,*ic,n=A->m,*ai=a->i,*aj=a->j;
9095a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
9105a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
9115a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
91277c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
913329f5518SBarry Smith   PetscReal          f;
9145a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
9155a8e39fbSHong Zhang   PetscBT            lnkbt;
9165a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
917a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
918a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
9199e25ed09SBarry Smith 
9203a40ed3dSBarry Smith   PetscFunctionBegin;
9216cf997b0SBarry Smith   f             = info->fill;
92238baddfdSBarry Smith   levels        = (PetscInt)info->levels;
92338baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9244c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
92582bf6240SBarry Smith 
92608480c60SBarry Smith   /* special case that simply copies fill pattern */
927be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
928be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
92986bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9302e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
93108480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
93208480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
93308480c60SBarry Smith     if (!b->diag) {
9347c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
93508480c60SBarry Smith     }
9367c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
93708480c60SBarry Smith     b->row              = isrow;
93808480c60SBarry Smith     b->col              = iscol;
93982bf6240SBarry Smith     b->icol             = isicol;
94087828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
941f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
942c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
943c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9443a40ed3dSBarry Smith     PetscFunctionReturn(0);
94508480c60SBarry Smith   }
94608480c60SBarry Smith 
947898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
948898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9499e25ed09SBarry Smith 
9509e25ed09SBarry Smith   /* get new row pointers */
9515a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9525a8e39fbSHong Zhang   bi[0] = 0;
9535a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9545a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9555a8e39fbSHong Zhang   bdiag[0]  = 0;
9566cf997b0SBarry Smith 
9575a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9585a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9595a8e39fbSHong Zhang 
9605a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9615a8e39fbSHong Zhang   nlnk = n + 1;
9625a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9635a8e39fbSHong Zhang 
9645a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
965a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9665a8e39fbSHong Zhang   current_space = free_space;
967a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9685a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
9695a8e39fbSHong Zhang 
9705a8e39fbSHong Zhang   for (i=0; i<n; i++) {
9715a8e39fbSHong Zhang     nzi = 0;
9726cf997b0SBarry Smith     /* copy current row into linked list */
9735a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
9745a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
9755a8e39fbSHong Zhang     cols = aj + ai[r[i]];
9765a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
9775a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9785a8e39fbSHong Zhang     nzi += nlnk;
9796cf997b0SBarry Smith 
9806cf997b0SBarry Smith     /* make sure diagonal entry is included */
9815a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
9826cf997b0SBarry Smith       fm = n;
9835a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
9845a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
9855a8e39fbSHong Zhang       lnk[fm]    = i;
9865a8e39fbSHong Zhang       lnk_lvl[i] = 0;
9875a8e39fbSHong Zhang       nzi++; dcount++;
9886cf997b0SBarry Smith     }
9896cf997b0SBarry Smith 
9905a8e39fbSHong Zhang     /* add pivot rows into the active row */
9915a8e39fbSHong Zhang     nzbd = 0;
9925a8e39fbSHong Zhang     prow = lnk[n];
9935a8e39fbSHong Zhang     while (prow < i) {
9945a8e39fbSHong Zhang       nnz      = bdiag[prow];
9955a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
9965a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
9975a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
9985a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
9995a8e39fbSHong Zhang       nzi += nlnk;
10005a8e39fbSHong Zhang       prow = lnk[prow];
10015a8e39fbSHong Zhang       nzbd++;
100256beaf04SBarry Smith     }
10035a8e39fbSHong Zhang     bdiag[i] = nzbd;
10045a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1005ecf371e4SBarry Smith 
10065a8e39fbSHong Zhang     /* if free space is not available, make more free space */
10075a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
10085a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1009a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1010a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
10115a8e39fbSHong Zhang       reallocs++;
10125a8e39fbSHong Zhang     }
1013ecf371e4SBarry Smith 
10145a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
10155a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
10165a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
10175a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
10185a8e39fbSHong Zhang 
10195a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
10205a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
102177431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10225a8e39fbSHong Zhang     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",i);
10236cf997b0SBarry Smith     }
10245a8e39fbSHong Zhang 
10255a8e39fbSHong Zhang     current_space->array           += nzi;
10265a8e39fbSHong Zhang     current_space->local_used      += nzi;
10275a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
10285a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
10295a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
10305a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
10319e25ed09SBarry Smith   }
10325a8e39fbSHong Zhang 
1033898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1034898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
10355a8e39fbSHong Zhang 
10365a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
10375a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1038a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
10395a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1040a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
10415a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10429e25ed09SBarry Smith 
104363ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
1044f64065bbSBarry Smith   {
10455a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
104663ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr);
104763ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af));CHKERRQ(ierr);
104863ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af));CHKERRQ(ierr);
104963ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"));CHKERRQ(ierr);
1050335d9088SBarry Smith     if (diagonal_fill) {
105163ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount));CHKERRQ(ierr);
1052335d9088SBarry Smith     }
1053f64065bbSBarry Smith   }
105463ba0a88SBarry Smith #endif
1055bbb6d6a8SBarry Smith 
10569e25ed09SBarry Smith   /* put together the new matrix */
1057f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
1058f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
1059f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1060ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
106152e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1062416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1063a96a251dSBarry Smith   b->freedata     = PETSC_TRUE;
10647c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10655a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1066b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10675a8e39fbSHong Zhang   b->j          = bj;
10685a8e39fbSHong Zhang   b->i          = bi;
10695a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
10705a8e39fbSHong Zhang   b->diag       = bdiag;
1071416022c9SBarry Smith   b->ilen       = 0;
1072416022c9SBarry Smith   b->imax       = 0;
1073416022c9SBarry Smith   b->row        = isrow;
1074416022c9SBarry Smith   b->col        = iscol;
1075c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1076c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
107782bf6240SBarry Smith   b->icol       = isicol;
107887828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1079416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
10805a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
10815a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
10825a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
10839e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1084418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1085ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
10865a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
108771c2f376SKris Buschelman 
108854e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
108971c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
10904846f1f5SKris Buschelman 
10913a40ed3dSBarry Smith   PetscFunctionReturn(0);
10929e25ed09SBarry Smith }
1093d5d45c9bSBarry Smith 
10943c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1095a6175056SHong Zhang #undef __FUNCT__
1096a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1097af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1098a6175056SHong Zhang {
10992ed133dbSHong Zhang   Mat            C = *B;
11000968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11012ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11022ed133dbSHong Zhang   IS             ip=b->row;
11032ed133dbSHong Zhang   PetscErrorCode ierr;
11042ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
11052ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11062ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1107622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1108540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1109fbf22428SSatish Balay   PetscReal      shiftpd;
1110540703acSHong Zhang   ChShift_Ctx    sctx;
1111540703acSHong Zhang   PetscInt       newshift;
1112d5d45c9bSBarry Smith 
1113a6175056SHong Zhang   PetscFunctionBegin;
1114540703acSHong Zhang   shiftnz   = info->shiftnz;
1115540703acSHong Zhang   shiftpd   = info->shiftpd;
1116ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1117ee45ca4aSHong Zhang 
11182ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11192ed133dbSHong Zhang 
11202ed133dbSHong Zhang   /* initialization */
11212ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11222ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11232ed133dbSHong Zhang   jl   = il + mbs;
11242ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11252ed133dbSHong Zhang 
1126540703acSHong Zhang   sctx.shift_amount = 0;
1127540703acSHong Zhang   sctx.nshift       = 0;
11282ed133dbSHong Zhang   do {
1129540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
11302ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11312ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1132a6175056SHong Zhang     }
11332ed133dbSHong Zhang 
11342ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1135c05c3958SHong Zhang       bval = ba + bi[k];
11362ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11372ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11382ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11392ed133dbSHong Zhang         col = rip[aj[j]];
11402ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11412ed133dbSHong Zhang           rtmp[col] = aa[j];
1142c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11432ed133dbSHong Zhang         }
11442ed133dbSHong Zhang       }
114539e3d630SHong Zhang       /* shift the diagonal of the matrix */
1146540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11472ed133dbSHong Zhang 
11482ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11492ed133dbSHong Zhang       dk = rtmp[k];
11502ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11512ed133dbSHong Zhang 
11522ed133dbSHong Zhang       while (i < k){
11532ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11542ed133dbSHong Zhang 
11552ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11562ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11572ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11582ed133dbSHong Zhang         dk += uikdi*ba[ili];
11592ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11602ed133dbSHong Zhang 
11612ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11622ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11632ed133dbSHong Zhang         if (jmin < jmax){
11642ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11652ed133dbSHong Zhang           /* update il and jl for row i */
11662ed133dbSHong Zhang           il[i] = jmin;
11672ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11682ed133dbSHong Zhang         }
11692ed133dbSHong Zhang         i = nexti;
11702ed133dbSHong Zhang       }
11712ed133dbSHong Zhang 
11720697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11730697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11740697b6caSHong Zhang       rs   = 0.0;
11753c9e8598SHong Zhang       jmin = bi[k]+1;
11763c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11773c9e8598SHong Zhang       if (nz){
11783c9e8598SHong Zhang         bcol = bj + jmin;
11793c9e8598SHong Zhang         while (nz--){
118089f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
118189f845c8SHong Zhang           bcol++;
11823c9e8598SHong Zhang         }
11833c9e8598SHong Zhang       }
1184540703acSHong Zhang 
1185540703acSHong Zhang       sctx.rs = rs;
1186540703acSHong Zhang       sctx.pv = dk;
11873e8c821fSHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1188540703acSHong Zhang       if (newshift == 1){
1189540703acSHong Zhang         break;    /* sctx.shift_amount is updated */
1190540703acSHong Zhang       } else if (newshift == -1){
11910697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
11923c9e8598SHong Zhang       }
11933c9e8598SHong Zhang 
11943c9e8598SHong Zhang       /* copy data into U(k,:) */
119539e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
119639e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
119739e3d630SHong Zhang       if (jmin < jmax) {
119839e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
119939e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12003c9e8598SHong Zhang         }
120139e3d630SHong Zhang         /* add the k-th row into il and jl */
12023c9e8598SHong Zhang         il[k] = jmin;
12033c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12043c9e8598SHong Zhang       }
12053c9e8598SHong Zhang     }
1206540703acSHong Zhang   } while (sctx.chshift);
12073c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12083c9e8598SHong Zhang 
120939e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12103c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12113c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12123c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1213efee365bSSatish Balay   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
1214540703acSHong Zhang   if (sctx.nshift){
1215540703acSHong Zhang     if (shiftnz) {
121663ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1217540703acSHong Zhang     } else if (shiftpd) {
121863ba0a88SBarry Smith       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqAIJ: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
12190697b6caSHong Zhang     }
12203c9e8598SHong Zhang   }
12213c9e8598SHong Zhang   PetscFunctionReturn(0);
12223c9e8598SHong Zhang }
1223a6175056SHong Zhang 
1224a6175056SHong Zhang #undef __FUNCT__
1225a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1226dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1227a6175056SHong Zhang {
12280968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1229ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1230ed59401aSHong Zhang   Mat                B;
1231dfbe8321SBarry Smith   PetscErrorCode     ierr;
1232653a6975SHong Zhang   PetscTruth         perm_identity;
1233622e2028SHong Zhang   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1234ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1235391eac55SHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12365a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1237ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1238a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1239a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12407a48dd6fSHong Zhang   PetscBT            lnkbt;
1241a6175056SHong Zhang 
1242a6175056SHong Zhang   PetscFunctionBegin;
1243653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12447a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12457a48dd6fSHong Zhang 
124639e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
124739e3d630SHong Zhang   ui[0] = 0;
124839e3d630SHong Zhang 
1249622e2028SHong Zhang   /* special case that simply copies fill pattern */
1250622e2028SHong Zhang   if (!levels && perm_identity) {
1251622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1252ed59401aSHong Zhang     for (i=0; i<am; i++) {
125339e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1254ed59401aSHong Zhang     }
125539e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
125639e3d630SHong Zhang     cols = uj;
1257ed59401aSHong Zhang     for (i=0; i<am; i++) {
1258ed59401aSHong Zhang       aj    = a->j + a->diag[i];
125939e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
126039e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1261ed59401aSHong Zhang     }
126239e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12637a48dd6fSHong Zhang     /* initialization */
12645a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12657a48dd6fSHong Zhang 
12667a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12677a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12681393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12697a48dd6fSHong Zhang     il         = jl + am;
12707a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12717a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12727a48dd6fSHong Zhang     for (i=0; i<am; i++){
12737a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12747a48dd6fSHong Zhang     }
12757a48dd6fSHong Zhang 
12767a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12777a48dd6fSHong Zhang     nlnk = am + 1;
12782ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12797a48dd6fSHong Zhang 
12807a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1281a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12827a48dd6fSHong Zhang     current_space = free_space;
1283a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
12847a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
12857a48dd6fSHong Zhang 
12867a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
12877a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
12887a48dd6fSHong Zhang       nzk   = 0;
1289622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1290622e2028SHong Zhang       ncols_upper = 0;
1291622e2028SHong Zhang       for (j=0; j<ncols; j++){
12925a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
12935a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
12945a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1295622e2028SHong Zhang           ncols_upper++;
1296622e2028SHong Zhang         }
1297622e2028SHong Zhang       }
12985a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12997a48dd6fSHong Zhang       nzk += nlnk;
13007a48dd6fSHong Zhang 
13017a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13027a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13037a48dd6fSHong Zhang 
13047a48dd6fSHong Zhang       while (prow < k){
13057a48dd6fSHong Zhang         nextprow = jl[prow];
13067a48dd6fSHong Zhang 
13077a48dd6fSHong Zhang         /* merge prow into k-th row */
13087a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13097a48dd6fSHong Zhang         jmax = ui[prow+1];
13107a48dd6fSHong Zhang         ncols = jmax-jmin;
1311ed59401aSHong Zhang         i     = jmin - ui[prow];
1312ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
13135a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
13145a8e39fbSHong Zhang         j     = *(uj - 1);
13155a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
13167a48dd6fSHong Zhang         nzk += nlnk;
13177a48dd6fSHong Zhang 
13187a48dd6fSHong Zhang         /* update il and jl for prow */
13197a48dd6fSHong Zhang         if (jmin < jmax){
13207a48dd6fSHong Zhang           il[prow] = jmin;
13217a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13227a48dd6fSHong Zhang         }
13237a48dd6fSHong Zhang         prow = nextprow;
13247a48dd6fSHong Zhang       }
13257a48dd6fSHong Zhang 
13267a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13277a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13287a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13297a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1330a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1331a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
13327a48dd6fSHong Zhang         reallocs++;
13337a48dd6fSHong Zhang       }
13347a48dd6fSHong Zhang 
13357a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13362ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13377a48dd6fSHong Zhang 
13387a48dd6fSHong Zhang       /* add the k-th row into il and jl */
133939e3d630SHong Zhang       if (nzk > 1){
13407a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13417a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13427a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13437a48dd6fSHong Zhang       }
13447a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13457a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13467a48dd6fSHong Zhang 
13477a48dd6fSHong Zhang       current_space->array           += nzk;
13487a48dd6fSHong Zhang       current_space->local_used      += nzk;
13497a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13507a48dd6fSHong Zhang 
13517a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13527a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13537a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13547a48dd6fSHong Zhang 
13557a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13567a48dd6fSHong Zhang     }
13577a48dd6fSHong Zhang 
135863ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
13597a48dd6fSHong Zhang     if (ai[am] != 0) {
136039e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
136163ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
136263ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
136363ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
13647a48dd6fSHong Zhang     } else {
136563ba0a88SBarry Smith       ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr);
13667a48dd6fSHong Zhang     }
136763ba0a88SBarry Smith #endif
13687a48dd6fSHong Zhang 
13697a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13707a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13715a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13727a48dd6fSHong Zhang 
13737a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13747a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1375a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13762ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1377a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13787a48dd6fSHong Zhang 
137939e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
138039e3d630SHong Zhang 
13817a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1382f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1383f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
1384ed59401aSHong Zhang   B = *fact;
1385ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1386ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13877a48dd6fSHong Zhang 
1388ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
13897a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
13907a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
13917a48dd6fSHong Zhang   b->j    = uj;
13927a48dd6fSHong Zhang   b->i    = ui;
13937a48dd6fSHong Zhang   b->diag = 0;
13947a48dd6fSHong Zhang   b->ilen = 0;
13957a48dd6fSHong Zhang   b->imax = 0;
13967a48dd6fSHong Zhang   b->row  = perm;
13977a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
13987a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
13997a48dd6fSHong Zhang   b->icol = perm;
14007a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14017a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
140252e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
14037a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
14047a48dd6fSHong Zhang 
1405ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1406ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1407ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14087a48dd6fSHong Zhang   if (ai[am] != 0) {
1409ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14107a48dd6fSHong Zhang   } else {
1411ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14127a48dd6fSHong Zhang   }
141339e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1414622e2028SHong Zhang   if (perm_identity){
1415ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1416ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1417622e2028SHong Zhang   }
1418a6175056SHong Zhang   PetscFunctionReturn(0);
1419a6175056SHong Zhang }
1420d5d45c9bSBarry Smith 
1421f76d2b81SHong Zhang #undef __FUNCT__
1422f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1423dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1424f76d2b81SHong Zhang {
1425f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
142610c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
14272ed133dbSHong Zhang   Mat                B;
1428dfbe8321SBarry Smith   PetscErrorCode     ierr;
1429f76d2b81SHong Zhang   PetscTruth         perm_identity;
143010c27e3dSHong Zhang   PetscReal          fill = info->fill;
14311393bc97SHong Zhang   PetscInt           *rip,*riip,i,am=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
143210c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14332ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1434a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
143510c27e3dSHong Zhang   PetscBT            lnkbt;
14362ed133dbSHong Zhang   IS                 iperm;
1437f76d2b81SHong Zhang 
1438f76d2b81SHong Zhang   PetscFunctionBegin;
14392ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1440f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14412ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14422ed133dbSHong Zhang 
1443f76d2b81SHong Zhang   if (!perm_identity){
14442ed133dbSHong Zhang     /* check if perm is symmetric! */
14452ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14462ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14471393bc97SHong Zhang     for (i=0; i<am; i++) {
14482ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14492ed133dbSHong Zhang     }
14502ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14512ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1452f76d2b81SHong Zhang   }
145310c27e3dSHong Zhang 
145410c27e3dSHong Zhang   /* initialization */
14551393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
145610c27e3dSHong Zhang   ui[0] = 0;
145710c27e3dSHong Zhang 
145810c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14591393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14601393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14611393bc97SHong Zhang   il     = jl + am;
14621393bc97SHong Zhang   cols   = il + am;
14631393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14641393bc97SHong Zhang   for (i=0; i<am; i++){
14651393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1466f76d2b81SHong Zhang   }
146710c27e3dSHong Zhang 
146810c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14691393bc97SHong Zhang   nlnk = am + 1;
14701393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
147110c27e3dSHong Zhang 
14721393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1473a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
147410c27e3dSHong Zhang   current_space = free_space;
147510c27e3dSHong Zhang 
14761393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
147710c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
147810c27e3dSHong Zhang     nzk   = 0;
14792ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14802ed133dbSHong Zhang     ncols_upper = 0;
14812ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1482622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
14832ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
14842ed133dbSHong Zhang         cols[ncols_upper] = i;
14852ed133dbSHong Zhang         ncols_upper++;
14862ed133dbSHong Zhang       }
14872ed133dbSHong Zhang     }
14881393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
148910c27e3dSHong Zhang     nzk += nlnk;
149010c27e3dSHong Zhang 
149110c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
149210c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
149310c27e3dSHong Zhang 
149410c27e3dSHong Zhang     while (prow < k){
149510c27e3dSHong Zhang       nextprow = jl[prow];
149610c27e3dSHong Zhang       /* merge prow into k-th row */
14971393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
149810c27e3dSHong Zhang       jmax = ui[prow+1];
149910c27e3dSHong Zhang       ncols = jmax-jmin;
15001393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15015a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150210c27e3dSHong Zhang       nzk += nlnk;
150310c27e3dSHong Zhang 
150410c27e3dSHong Zhang       /* update il and jl for prow */
150510c27e3dSHong Zhang       if (jmin < jmax){
150610c27e3dSHong Zhang         il[prow] = jmin;
15072ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
150810c27e3dSHong Zhang       }
150910c27e3dSHong Zhang       prow = nextprow;
151010c27e3dSHong Zhang     }
151110c27e3dSHong Zhang 
151210c27e3dSHong Zhang     /* if free space is not available, make more free space */
151310c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15141393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
151510c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1516a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
151710c27e3dSHong Zhang       reallocs++;
151810c27e3dSHong Zhang     }
151910c27e3dSHong Zhang 
152010c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15211393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
152210c27e3dSHong Zhang 
152310c27e3dSHong Zhang     /* add the k-th row into il and jl */
152410c27e3dSHong Zhang     if (nzk-1 > 0){
15251393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
152610c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
152710c27e3dSHong Zhang       il[k] = ui[k] + 1;
152810c27e3dSHong Zhang     }
15292ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
153010c27e3dSHong Zhang     current_space->array           += nzk;
153110c27e3dSHong Zhang     current_space->local_used      += nzk;
153210c27e3dSHong Zhang     current_space->local_remaining -= nzk;
153310c27e3dSHong Zhang 
153410c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
153510c27e3dSHong Zhang   }
153610c27e3dSHong Zhang 
153763ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG)
15381393bc97SHong Zhang   if (ai[am] != 0) {
15391393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
154063ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
154163ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
154263ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
154310c27e3dSHong Zhang   } else {
154463ba0a88SBarry Smith      ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"));CHKERRQ(ierr);
154510c27e3dSHong Zhang   }
154663ba0a88SBarry Smith #endif
154710c27e3dSHong Zhang 
154810c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
154910c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
155010c27e3dSHong Zhang 
155110c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15521393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1553a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
155410c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
155510c27e3dSHong Zhang 
155610c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1557f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1558f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
15592ed133dbSHong Zhang   B    = *fact;
15602ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1561ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
156210c27e3dSHong Zhang 
15632ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
156410c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
15651393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
156610c27e3dSHong Zhang   b->j    = uj;
156710c27e3dSHong Zhang   b->i    = ui;
156810c27e3dSHong Zhang   b->diag = 0;
156910c27e3dSHong Zhang   b->ilen = 0;
157010c27e3dSHong Zhang   b->imax = 0;
157110c27e3dSHong Zhang   b->row  = perm;
157210c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
157310c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
157410c27e3dSHong Zhang   b->icol = perm;
157510c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15761393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15771393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
15781393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
157910c27e3dSHong Zhang 
15802ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
15812ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
15822ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
15831393bc97SHong Zhang   if (ai[am] != 0) {
15841393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
158510c27e3dSHong Zhang   } else {
15862ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
158710c27e3dSHong Zhang   }
158839e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
15892ed133dbSHong Zhang   if (perm_identity){
159010c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
159110c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15922ed133dbSHong Zhang   }
1593f76d2b81SHong Zhang   PetscFunctionReturn(0);
1594f76d2b81SHong Zhang }
1595