xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision e6b907ac9e74227fac7d27fb51c2636d4aaee38b)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
2386bacbd2SBarry Smith 
24e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
269cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2738baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
28e34fafa9SBarry Smith #endif
2986bacbd2SBarry Smith 
304a2ae208SSatish Balay #undef __FUNCT__
314a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3286bacbd2SBarry Smith   /* ------------------------------------------------------------
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith           This interface was contribed by Tony Caola
3586bacbd2SBarry Smith 
3686bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
375bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
385bd2ddc7SBarry Smith      SPARSEKIT2.
395bd2ddc7SBarry Smith 
405bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
415bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4286bacbd2SBarry Smith 
4386bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4486bacbd2SBarry Smith      help in getting this routine ironed out.
4586bacbd2SBarry Smith 
465bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
475bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
485bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
495bd2ddc7SBarry Smith      Yousef Saad's code.
5086bacbd2SBarry Smith 
5186bacbd2SBarry Smith      ------------------------------------------------------------
5286bacbd2SBarry Smith */
537529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5486bacbd2SBarry Smith {
5560ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5660ec86bdSBarry Smith   PetscFunctionBegin;
57e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5860ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5960ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
6060ec86bdSBarry Smith #else
6186bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6207d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6386bacbd2SBarry Smith   PetscTruth     reorder;
649cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
65899cda47SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->rmap.n;
6638baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6738baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6838baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6987828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
7054a8161fSBarry Smith   PetscReal      af;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   PetscFunctionBegin;
7386bacbd2SBarry Smith 
7486bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7538baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7686bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7733259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7838baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7938baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
8086bacbd2SBarry Smith 
8186bacbd2SBarry Smith 
8286bacbd2SBarry Smith   /* ------------------------------------------------------------
8386bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8486bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8586bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8686bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8786bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8886bacbd2SBarry Smith      ------------------------------------------------------------
8986bacbd2SBarry Smith   */
9086bacbd2SBarry Smith 
9186bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9286bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9386bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9486bacbd2SBarry Smith   reorder = PetscNot(reorder);
9586bacbd2SBarry Smith 
9686bacbd2SBarry Smith 
9786bacbd2SBarry Smith   /* storage for ilu factor */
9838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9938baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
10087828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
10138baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10286bacbd2SBarry Smith 
10386bacbd2SBarry Smith   /* ------------------------------------------------------------
10486bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10586bacbd2SBarry Smith      ------------------------------------------------------------
10686bacbd2SBarry Smith   */
107b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
108b19713cbSBarry Smith     old_j[i]++;
10986bacbd2SBarry Smith   }
110b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
111b19713cbSBarry Smith     old_i[i]++;
112b19713cbSBarry Smith   };
113010a20caSHong Zhang 
11486bacbd2SBarry Smith 
11586bacbd2SBarry Smith   if (reorder) {
116c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
117c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
118c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
119c0c2c81eSBarry Smith       r[i]  = r[i]+1;
120c0c2c81eSBarry Smith       c[i]  = c[i]+1;
121c0c2c81eSBarry Smith     }
12238baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12338baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12487828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1255bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
126c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
127c0c2c81eSBarry Smith       r[i]  = r[i]-1;
128c0c2c81eSBarry Smith       c[i]  = c[i]-1;
129c0c2c81eSBarry Smith     }
130c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
131c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
132b19713cbSBarry Smith     o_a = old_a2;
133b19713cbSBarry Smith     o_j = old_j2;
134b19713cbSBarry Smith     o_i = old_i2;
135b19713cbSBarry Smith   } else {
136b19713cbSBarry Smith     o_a = old_a;
137b19713cbSBarry Smith     o_j = old_j;
138b19713cbSBarry Smith     o_i = old_i;
13986bacbd2SBarry Smith   }
14086bacbd2SBarry Smith 
14186bacbd2SBarry Smith   /* ------------------------------------------------------------
14286bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14386bacbd2SBarry Smith      ------------------------------------------------------------
14486bacbd2SBarry Smith   */
14586bacbd2SBarry Smith 
14638baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14738baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14887828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14986bacbd2SBarry Smith 
15054a8161fSBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1516849ba73SBarry Smith   if (sierr) {
1526849ba73SBarry Smith     switch (sierr) {
153a83599f4SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
154a83599f4SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
155e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
156e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15777431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15877431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15986bacbd2SBarry Smith     }
16086bacbd2SBarry Smith   }
16186bacbd2SBarry Smith 
16286bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16386bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16486bacbd2SBarry Smith 
16586bacbd2SBarry Smith   /* ------------------------------------------------------------
16686bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16786bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16886bacbd2SBarry Smith      ------------------------------------------------------------
16986bacbd2SBarry Smith   */
17086bacbd2SBarry Smith 
17187828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
17238baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17386bacbd2SBarry Smith 
1745bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17586bacbd2SBarry Smith 
17686bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17786bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17886bacbd2SBarry Smith 
17986bacbd2SBarry Smith   if (reorder) {
18086bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
18186bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18286bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
183313ae024SBarry Smith   } else {
184b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
185f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
186b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
187b19713cbSBarry Smith     }
188b19713cbSBarry Smith   }
189b19713cbSBarry Smith 
190b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
19186bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
192b19713cbSBarry Smith     old_i[i]--;
193b19713cbSBarry Smith   }
194b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
195b19713cbSBarry Smith     old_j[i]--;
196b19713cbSBarry Smith   }
19786bacbd2SBarry Smith 
198b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
200b19713cbSBarry Smith     new_i[i]--;
201b19713cbSBarry Smith   }
202b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
203b19713cbSBarry Smith     new_j[i]--;
204b19713cbSBarry Smith   }
20586bacbd2SBarry Smith 
20686bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20786bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2084c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2094c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
210c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
21186bacbd2SBarry Smith   for(i=0; i<n; i++) {
21286bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21386bacbd2SBarry Smith   };
21486bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21507d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21686bacbd2SBarry Smith 
217c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
218c0c2c81eSBarry Smith 
219329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
22007d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
22186bacbd2SBarry Smith 
22286bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22386bacbd2SBarry Smith 
224f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
225f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
226f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
227ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
22886bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22986bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
23086bacbd2SBarry Smith 
23186bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
232*e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
233*e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
23486bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23507d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23686bacbd2SBarry Smith   b->a             = new_a;
23786bacbd2SBarry Smith   b->j             = new_j;
23886bacbd2SBarry Smith   b->i             = new_i;
23986bacbd2SBarry Smith   b->ilen          = 0;
24086bacbd2SBarry Smith   b->imax          = 0;
2411f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
242313ae024SBarry Smith   b->row           = isirow;
24386bacbd2SBarry Smith   b->col           = iscolf;
24487828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24586bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24686bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24786bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24886bacbd2SBarry Smith 
24986bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
25086bacbd2SBarry Smith 
251431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
252ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
253ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
254ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
255ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
256431e710aSBarry Smith 
2577529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25871c2f376SKris Buschelman 
25986bacbd2SBarry Smith   PetscFunctionReturn(0);
26060ec86bdSBarry Smith #endif
26186bacbd2SBarry Smith }
26286bacbd2SBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
264b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
265dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
266289bc588SBarry Smith {
267416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
268289bc588SBarry Smith   IS                 isicol;
2696849ba73SBarry Smith   PetscErrorCode     ierr;
270899cda47SBarry Smith   PetscInt           *r,*ic,i,n=A->rmap.n,*ai=a->i,*aj=a->j;
2711393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
2721393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2739e046878SBarry Smith   PetscReal          f;
2744875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
275a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2761393bc97SHong Zhang   PetscBT            lnkbt;
277289bc588SBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionBegin;
279899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2804c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
281ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
282ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
283289bc588SBarry Smith 
284289bc588SBarry Smith   /* get new row pointers */
2851393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2861393bc97SHong Zhang   bi[0] = 0;
2871393bc97SHong Zhang 
2881393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2891393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2901393bc97SHong Zhang   bdiag[0] = 0;
2911393bc97SHong Zhang 
2921393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2931393bc97SHong Zhang   nlnk = n + 1;
2941393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2951393bc97SHong Zhang 
2964875ba38SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&im);CHKERRQ(ierr);
2971393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2981393bc97SHong Zhang 
2991393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
300d3d32019SBarry Smith   f = info->fill;
301a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
3021393bc97SHong Zhang   current_space = free_space;
303289bc588SBarry Smith 
304289bc588SBarry Smith   for (i=0; i<n; i++) {
3051393bc97SHong Zhang     /* copy previous fill into linked list */
306289bc588SBarry Smith     nzi = 0;
3071393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3081393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3091393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
310afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3111393bc97SHong Zhang     nzi += nlnk;
3121393bc97SHong Zhang 
3131393bc97SHong Zhang     /* add pivot rows into linked list */
3141393bc97SHong Zhang     row = lnk[n];
3151393bc97SHong Zhang     while (row < i) {
316d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
317d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
318d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3191393bc97SHong Zhang       nzi += nlnk;
3201393bc97SHong Zhang       row  = lnk[row];
321289bc588SBarry Smith     }
3221393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3231393bc97SHong Zhang     im[i]   = nzi;
3241393bc97SHong Zhang 
3251393bc97SHong Zhang     /* mark bdiag */
3261393bc97SHong Zhang     nzbd = 0;
3271393bc97SHong Zhang     nnz  = nzi;
3281393bc97SHong Zhang     k    = lnk[n];
3291393bc97SHong Zhang     while (nnz-- && k < i){
3301393bc97SHong Zhang       nzbd++;
3311393bc97SHong Zhang       k = lnk[k];
3321393bc97SHong Zhang     }
3331393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3341393bc97SHong Zhang 
3351393bc97SHong Zhang     /* if free space is not available, make more free space */
3361393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3371393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
338a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3391393bc97SHong Zhang       reallocs++;
3401393bc97SHong Zhang     }
3411393bc97SHong Zhang 
3421393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3431393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3441393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3451393bc97SHong Zhang     current_space->array           += nzi;
3461393bc97SHong Zhang     current_space->local_used      += nzi;
3471393bc97SHong Zhang     current_space->local_remaining -= nzi;
348289bc588SBarry Smith   }
3496cf91177SBarry Smith #if defined(PETSC_USE_INFO)
350464e8d44SSatish Balay   if (ai[n] != 0) {
3511393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
352ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
353ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
354ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
355ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
35605bf559cSSatish Balay   } else {
357ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
35805bf559cSSatish Balay   }
35963ba0a88SBarry Smith #endif
3602fb3ed76SBarry Smith 
361898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
362898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3631987afe7SBarry Smith 
3641393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3651393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
366a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3671393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3684875ba38SHong Zhang   ierr = PetscFree(im);CHKERRQ(ierr);
369289bc588SBarry Smith 
370289bc588SBarry Smith   /* put together the new matrix */
371f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,B);CHKERRQ(ierr);
372f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
373f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
374ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
37552e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
376416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
377*e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
378*e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
3797c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3801393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3811393bc97SHong Zhang   b->j          = bj;
3821393bc97SHong Zhang   b->i          = bi;
3831393bc97SHong Zhang   b->diag       = bdiag;
384416022c9SBarry Smith   b->ilen       = 0;
385416022c9SBarry Smith   b->imax       = 0;
386416022c9SBarry Smith   b->row        = isrow;
387416022c9SBarry Smith   b->col        = iscol;
388c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
389c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
39082bf6240SBarry Smith   b->icol       = isicol;
39187828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3921393bc97SHong Zhang 
3931393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3941393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3951393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
3967fd98bd6SLois Curfman McInnes 
397329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
398418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
399ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
400703deb49SSatish Balay 
401e93ccc4dSBarry Smith   if (ai[n] != 0) {
4021393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
403eea2913aSSatish Balay   } else {
404eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
405eea2913aSSatish Balay   }
4064846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
40771c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4083a40ed3dSBarry Smith   PetscFunctionReturn(0);
409289bc588SBarry Smith }
4101393bc97SHong Zhang 
4115b5bc046SBarry Smith /*
4125b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4135b5bc046SBarry Smith */
4145b5bc046SBarry Smith #undef __FUNCT__
4155b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4165b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4175b5bc046SBarry Smith {
4185b5bc046SBarry Smith   PetscErrorCode ierr;
4195b5bc046SBarry Smith   PetscTruth     flg;
4205b5bc046SBarry Smith 
4215b5bc046SBarry Smith   PetscFunctionBegin;
4225b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4235b5bc046SBarry Smith   if (flg) {
4245b5bc046SBarry Smith     PetscViewer viewer;
4255b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4265b5bc046SBarry Smith 
4275b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4285b5bc046SBarry Smith     ierr = PetscViewerBinaryOpen(A->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4295b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4305b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4315b5bc046SBarry Smith   }
4325b5bc046SBarry Smith   PetscFunctionReturn(0);
4335b5bc046SBarry Smith }
4345b5bc046SBarry Smith 
4350a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4364a2ae208SSatish Balay #undef __FUNCT__
4374a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
438af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
439289bc588SBarry Smith {
440416022c9SBarry Smith   Mat            C=*B;
441416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
44282bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4436849ba73SBarry Smith   PetscErrorCode ierr;
444899cda47SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap.n,*bi=b->i,*bj=b->j;
445b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
446d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
4475b5bc046SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4486a7c8fc2SHong Zhang   PetscScalar    d;
4496a7c8fc2SHong Zhang   PetscReal      rs;
450b3bf805bSHong Zhang   LUShift_Ctx    sctx;
45142898933SHong Zhang   PetscInt       newshift,*ddiag;
452289bc588SBarry Smith 
4533a40ed3dSBarry Smith   PetscFunctionBegin;
45478b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
45578b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
45687828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45787828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
458010a20caSHong Zhang   rtmps = rtmp; ics = ic;
459289bc588SBarry Smith 
460ada7143aSSatish Balay   sctx.shift_top  = 0;
461ada7143aSSatish Balay   sctx.nshift_max = 0;
462ada7143aSSatish Balay   sctx.shift_lo   = 0;
463ada7143aSSatish Balay   sctx.shift_hi   = 0;
464ada7143aSSatish Balay 
46542898933SHong Zhang   if (!a->diag) {
4660cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4670cf777b8SBarry Smith   }
468118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
469118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4701a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
4719deb78dcSHong Zhang     PetscInt *aai = a->i;
47242898933SHong Zhang     ddiag          = a->diag;
473118f5789SBarry Smith     sctx.shift_top = 0;
4746cc28720Svictorle     for (i=0; i<n; i++) {
4759f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4769f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4771a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
478010a20caSHong Zhang       v  = a->a+aai[i];
479e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
480e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4811a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4821a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4836cc28720Svictorle     }
484c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
485118f5789SBarry Smith     sctx.shift_top    *= 1.1;
486d2276718SHong Zhang     sctx.nshift_max   = 5;
487d2276718SHong Zhang     sctx.shift_lo     = 0.;
488d2276718SHong Zhang     sctx.shift_hi     = 1.;
489a0a356efSHong Zhang   }
490a0a356efSHong Zhang 
491a0a356efSHong Zhang   sctx.shift_amount = 0;
492a0a356efSHong Zhang   sctx.nshift       = 0;
493085256b3SBarry Smith   do {
494d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
495289bc588SBarry Smith     for (i=0; i<n; i++){
496b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
497b3bf805bSHong Zhang       bjtmp = bj + bi[i];
498b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
499289bc588SBarry Smith 
500289bc588SBarry Smith       /* load in initial (unfactored row) */
501416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
502b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
503010a20caSHong Zhang       v     = a->a + a->i[r[i]];
504085256b3SBarry Smith       for (j=0; j<nz; j++) {
505b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
506085256b3SBarry Smith       }
507d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
508289bc588SBarry Smith 
509b3bf805bSHong Zhang       row = *bjtmp++;
510f2582111SSatish Balay       while  (row < i) {
5118c37ef55SBarry Smith         pc = rtmp + row;
5128c37ef55SBarry Smith         if (*pc != 0.0) {
513010a20caSHong Zhang           pv         = b->a + diag_offset[row];
514010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
51535aab85fSBarry Smith           multiplier = *pc / *pv++;
5168c37ef55SBarry Smith           *pc        = multiplier;
517b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
518f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
519efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
520289bc588SBarry Smith         }
521b3bf805bSHong Zhang         row = *bjtmp++;
522289bc588SBarry Smith       }
523416022c9SBarry Smith       /* finished row so stick it into b->a */
524b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
525b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
526b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
527b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
528d3d32019SBarry Smith       rs   = 0.0;
529d3d32019SBarry Smith       for (j=0; j<nz; j++) {
530d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
531d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
532d3d32019SBarry Smith       }
533d2276718SHong Zhang 
534b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
535d2276718SHong Zhang       sctx.rs  = rs;
536d2276718SHong Zhang       sctx.pv  = pv[diag];
5375b5bc046SBarry Smith       ierr = MatLUCheckShift_inline(info,sctx,i,a->diag,newshift);CHKERRQ(ierr);
5385b5bc046SBarry Smith       if (newshift == 1) break;
53935aab85fSBarry Smith     }
540d2276718SHong Zhang 
541118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5426cc28720Svictorle       /*
5436c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5446cc28720Svictorle        * then try lower shift
5456cc28720Svictorle        */
546d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
547d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
548d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
549d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
550d2276718SHong Zhang       sctx.nshift++;
5516cc28720Svictorle     }
552d2276718SHong Zhang   } while (sctx.lushift);
5530f11f9aeSSatish Balay 
55435aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55535aab85fSBarry Smith   for (i=0; i<n; i++) {
556010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55735aab85fSBarry Smith   }
55835aab85fSBarry Smith 
559606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
56078b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
56178b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
562416022c9SBarry Smith   C->factor = FACTOR_LU;
5637dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
564c456f294SBarry Smith   C->assembled = PETSC_TRUE;
565899cda47SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
566d2276718SHong Zhang   if (sctx.nshift){
567118f5789SBarry Smith     if (info->shiftnz) {
568ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
569118f5789SBarry Smith     } else if (info->shiftpd) {
570ae15b995SBarry Smith       ierr = PetscInfo4(0,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
5716cc28720Svictorle     }
5720697b6caSHong Zhang   }
5733a40ed3dSBarry Smith   PetscFunctionReturn(0);
574289bc588SBarry Smith }
5750889b5dcSvictorle 
5760889b5dcSvictorle #undef __FUNCT__
5770889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
578dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5790889b5dcSvictorle {
5800889b5dcSvictorle   PetscFunctionBegin;
5810889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5820889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5830889b5dcSvictorle   PetscFunctionReturn(0);
5840889b5dcSvictorle }
5850889b5dcSvictorle 
5860889b5dcSvictorle 
5870a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5884a2ae208SSatish Balay #undef __FUNCT__
5894a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
590dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
591da3a660dSBarry Smith {
592dfbe8321SBarry Smith   PetscErrorCode ierr;
593416022c9SBarry Smith   Mat            C;
594416022c9SBarry Smith 
5953a40ed3dSBarry Smith   PetscFunctionBegin;
5969e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
597af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
598273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
59952e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
6003a40ed3dSBarry Smith   PetscFunctionReturn(0);
601da3a660dSBarry Smith }
6020a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
6034a2ae208SSatish Balay #undef __FUNCT__
6044a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
605dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6068c37ef55SBarry Smith {
607416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
608416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6096849ba73SBarry Smith   PetscErrorCode ierr;
610899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
61138baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
61287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6138c37ef55SBarry Smith 
6143a40ed3dSBarry Smith   PetscFunctionBegin;
6153a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
61691bf9adeSBarry Smith 
6171ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
619416022c9SBarry Smith   tmp  = a->solve_work;
6208c37ef55SBarry Smith 
6213b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6223b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6238c37ef55SBarry Smith 
6248c37ef55SBarry Smith   /* forward solve the lower triangular */
6258c37ef55SBarry Smith   tmp[0] = b[*r++];
626010a20caSHong Zhang   tmps   = tmp;
6278c37ef55SBarry Smith   for (i=1; i<n; i++) {
628010a20caSHong Zhang     v   = aa + ai[i] ;
629010a20caSHong Zhang     vi  = aj + ai[i] ;
63053e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
63153e7425aSBarry Smith     sum = b[*r++];
632ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6338c37ef55SBarry Smith     tmp[i] = sum;
6348c37ef55SBarry Smith   }
6358c37ef55SBarry Smith 
6368c37ef55SBarry Smith   /* backward solve the upper triangular */
6378c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
638010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
639010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
640416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6418c37ef55SBarry Smith     sum = tmp[i];
642ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
643010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6448c37ef55SBarry Smith   }
6458c37ef55SBarry Smith 
6463b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6473b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6481ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6491ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
650899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
6528c37ef55SBarry Smith }
653026e39d0SSatish Balay 
6542bebee5dSHong Zhang #undef __FUNCT__
6552bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
6562bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
6572bebee5dSHong Zhang {
6582bebee5dSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6592bebee5dSHong Zhang   IS             iscol = a->col,isrow = a->row;
6602bebee5dSHong Zhang   PetscErrorCode ierr;
661899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
6622bebee5dSHong Zhang   PetscInt       nz,*rout,*cout,neq;
6632bebee5dSHong Zhang   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6642bebee5dSHong Zhang 
6652bebee5dSHong Zhang   PetscFunctionBegin;
6662bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
6672bebee5dSHong Zhang 
6682bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
6692bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
6702bebee5dSHong Zhang 
6712bebee5dSHong Zhang   tmp  = a->solve_work;
6722bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6732bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
6742bebee5dSHong Zhang 
6752bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
6762bebee5dSHong Zhang     /* forward solve the lower triangular */
6772bebee5dSHong Zhang     tmp[0] = b[r[0]];
6782bebee5dSHong Zhang     tmps   = tmp;
6792bebee5dSHong Zhang     for (i=1; i<n; i++) {
6802bebee5dSHong Zhang       v   = aa + ai[i] ;
6812bebee5dSHong Zhang       vi  = aj + ai[i] ;
6822bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
6832bebee5dSHong Zhang       sum = b[r[i]];
6842bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6852bebee5dSHong Zhang       tmp[i] = sum;
6862bebee5dSHong Zhang     }
6872bebee5dSHong Zhang     /* backward solve the upper triangular */
6882bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
6892bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
6902bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
6912bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
6922bebee5dSHong Zhang       sum = tmp[i];
6932bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6942bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
6952bebee5dSHong Zhang     }
6962bebee5dSHong Zhang 
6972bebee5dSHong Zhang     b += n;
6982bebee5dSHong Zhang     x += n;
6992bebee5dSHong Zhang   }
7002bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7012bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7022bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
7032bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
7042bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
7052bebee5dSHong Zhang   PetscFunctionReturn(0);
7062bebee5dSHong Zhang }
7072bebee5dSHong Zhang 
708930ae53cSBarry Smith /* ----------------------------------------------------------- */
7094a2ae208SSatish Balay #undef __FUNCT__
7104a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
711dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
712930ae53cSBarry Smith {
713930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7146849ba73SBarry Smith   PetscErrorCode ierr;
715899cda47SBarry Smith   PetscInt       n = A->rmap.n,*ai = a->i,*aj = a->j,*adiag = a->diag;
716362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
717aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
71838baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
719362ced78SSatish Balay   PetscScalar    *v,sum;
720d85afc42SSatish Balay #endif
721930ae53cSBarry Smith 
7223a40ed3dSBarry Smith   PetscFunctionBegin;
7233a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
724930ae53cSBarry Smith 
7251ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
727930ae53cSBarry Smith 
728aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
7291c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
7301c4feecaSSatish Balay #else
731930ae53cSBarry Smith   /* forward solve the lower triangular */
732930ae53cSBarry Smith   x[0] = b[0];
733930ae53cSBarry Smith   for (i=1; i<n; i++) {
734930ae53cSBarry Smith     ai_i = ai[i];
735930ae53cSBarry Smith     v    = aa + ai_i;
736930ae53cSBarry Smith     vi   = aj + ai_i;
737930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
738930ae53cSBarry Smith     sum  = b[i];
739930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
740930ae53cSBarry Smith     x[i] = sum;
741930ae53cSBarry Smith   }
742930ae53cSBarry Smith 
743930ae53cSBarry Smith   /* backward solve the upper triangular */
744930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
745930ae53cSBarry Smith     adiag_i = adiag[i];
746d708749eSBarry Smith     v       = aa + adiag_i + 1;
747d708749eSBarry Smith     vi      = aj + adiag_i + 1;
748930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
749930ae53cSBarry Smith     sum     = x[i];
750930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
751930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
752930ae53cSBarry Smith   }
7531c4feecaSSatish Balay #endif
754899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
7551ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7561ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7573a40ed3dSBarry Smith   PetscFunctionReturn(0);
758930ae53cSBarry Smith }
759930ae53cSBarry Smith 
7604a2ae208SSatish Balay #undef __FUNCT__
7614a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
762dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
763da3a660dSBarry Smith {
764416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
765416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7666849ba73SBarry Smith   PetscErrorCode ierr;
767899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
76838baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
76987828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
770da3a660dSBarry Smith 
7713a40ed3dSBarry Smith   PetscFunctionBegin;
77278b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
773da3a660dSBarry Smith 
7741ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7751ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
776416022c9SBarry Smith   tmp  = a->solve_work;
777da3a660dSBarry Smith 
7783b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7793b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
780da3a660dSBarry Smith 
781da3a660dSBarry Smith   /* forward solve the lower triangular */
782da3a660dSBarry Smith   tmp[0] = b[*r++];
783da3a660dSBarry Smith   for (i=1; i<n; i++) {
784010a20caSHong Zhang     v   = aa + ai[i] ;
785010a20caSHong Zhang     vi  = aj + ai[i] ;
786416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
787da3a660dSBarry Smith     sum = b[*r++];
788010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
789da3a660dSBarry Smith     tmp[i] = sum;
790da3a660dSBarry Smith   }
791da3a660dSBarry Smith 
792da3a660dSBarry Smith   /* backward solve the upper triangular */
793da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
794010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
795010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
796416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
797da3a660dSBarry Smith     sum = tmp[i];
798010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
799010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
800da3a660dSBarry Smith     x[*c--] += tmp[i];
801da3a660dSBarry Smith   }
802da3a660dSBarry Smith 
8033b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8043b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8051ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
807efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
808898183ecSLois Curfman McInnes 
8093a40ed3dSBarry Smith   PetscFunctionReturn(0);
810da3a660dSBarry Smith }
811da3a660dSBarry Smith /* -------------------------------------------------------------------*/
8124a2ae208SSatish Balay #undef __FUNCT__
8134a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
814dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
815da3a660dSBarry Smith {
816416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8172235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8186849ba73SBarry Smith   PetscErrorCode ierr;
819899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
82038baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
82187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
822da3a660dSBarry Smith 
8233a40ed3dSBarry Smith   PetscFunctionBegin;
8241ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8251ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
826416022c9SBarry Smith   tmp  = a->solve_work;
827da3a660dSBarry Smith 
8282235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8292235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
830da3a660dSBarry Smith 
831da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
8322235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
833da3a660dSBarry Smith 
834da3a660dSBarry Smith   /* forward solve the U^T */
835da3a660dSBarry Smith   for (i=0; i<n; i++) {
836010a20caSHong Zhang     v   = aa + diag[i] ;
837010a20caSHong Zhang     vi  = aj + diag[i] + 1;
838f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
839c38d4ed2SBarry Smith     s1  = tmp[i];
840ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
841da3a660dSBarry Smith     while (nz--) {
842010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
843da3a660dSBarry Smith     }
844c38d4ed2SBarry Smith     tmp[i] = s1;
845da3a660dSBarry Smith   }
846da3a660dSBarry Smith 
847da3a660dSBarry Smith   /* backward solve the L^T */
848da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
849010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
850010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
851f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
852f1af5d2fSBarry Smith     s1  = tmp[i];
853da3a660dSBarry Smith     while (nz--) {
854010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
855da3a660dSBarry Smith     }
856da3a660dSBarry Smith   }
857da3a660dSBarry Smith 
858da3a660dSBarry Smith   /* copy tmp into x according to permutation */
859da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
860da3a660dSBarry Smith 
8612235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8622235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8631ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
865da3a660dSBarry Smith 
866899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap.n);CHKERRQ(ierr);
8673a40ed3dSBarry Smith   PetscFunctionReturn(0);
868da3a660dSBarry Smith }
869da3a660dSBarry Smith 
8704a2ae208SSatish Balay #undef __FUNCT__
8714a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
872dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
873da3a660dSBarry Smith {
874416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8752235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8766849ba73SBarry Smith   PetscErrorCode ierr;
877899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
87838baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
87987828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8806abc6512SBarry Smith 
8813a40ed3dSBarry Smith   PetscFunctionBegin;
88223bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8836abc6512SBarry Smith 
8841ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8851ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
886416022c9SBarry Smith   tmp = a->solve_work;
8876abc6512SBarry Smith 
8882235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8892235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8906abc6512SBarry Smith 
8916abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8922235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8936abc6512SBarry Smith 
8946abc6512SBarry Smith   /* forward solve the U^T */
8956abc6512SBarry Smith   for (i=0; i<n; i++) {
896010a20caSHong Zhang     v   = aa + diag[i] ;
897010a20caSHong Zhang     vi  = aj + diag[i] + 1;
898f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8996abc6512SBarry Smith     tmp[i] *= *v++;
9006abc6512SBarry Smith     while (nz--) {
901010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
9026abc6512SBarry Smith     }
9036abc6512SBarry Smith   }
9046abc6512SBarry Smith 
9056abc6512SBarry Smith   /* backward solve the L^T */
9066abc6512SBarry Smith   for (i=n-1; i>=0; i--){
907010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
908010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
909f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
9106abc6512SBarry Smith     while (nz--) {
911010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
9126abc6512SBarry Smith     }
9136abc6512SBarry Smith   }
9146abc6512SBarry Smith 
9156abc6512SBarry Smith   /* copy tmp into x according to permutation */
9166abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
9176abc6512SBarry Smith 
9182235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
9192235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
9201ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
9211ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9226abc6512SBarry Smith 
923efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
9243a40ed3dSBarry Smith   PetscFunctionReturn(0);
925da3a660dSBarry Smith }
9269e25ed09SBarry Smith /* ----------------------------------------------------------------*/
927dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
92808480c60SBarry Smith 
9294a2ae208SSatish Balay #undef __FUNCT__
9304a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
931dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
9329e25ed09SBarry Smith {
933416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
9349e25ed09SBarry Smith   IS                 isicol;
9356849ba73SBarry Smith   PetscErrorCode     ierr;
936899cda47SBarry Smith   PetscInt           *r,*ic,n=A->rmap.n,*ai=a->i,*aj=a->j;
9375a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
9385a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
9395a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
94077c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
941329f5518SBarry Smith   PetscReal          f;
9425a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
9435a8e39fbSHong Zhang   PetscBT            lnkbt;
9445a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
945a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
946a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
9479e25ed09SBarry Smith 
9483a40ed3dSBarry Smith   PetscFunctionBegin;
9496cf997b0SBarry Smith   f             = info->fill;
95038baddfdSBarry Smith   levels        = (PetscInt)info->levels;
95138baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9524c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
95382bf6240SBarry Smith 
95408480c60SBarry Smith   /* special case that simply copies fill pattern */
955be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
956be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
95786bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9582e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
95908480c60SBarry Smith     (*fact)->factor                 = FACTOR_LU;
960f3a39becSBarry Smith     (*fact)->info.factor_mallocs    = 0;
961f3a39becSBarry Smith     (*fact)->info.fill_ratio_given  = info->fill;
962f3a39becSBarry Smith     (*fact)->info.fill_ratio_needed = 1.0;
96308480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
96408480c60SBarry Smith     if (!b->diag) {
9657c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
96608480c60SBarry Smith     }
9677c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
96808480c60SBarry Smith     b->row              = isrow;
96908480c60SBarry Smith     b->col              = iscol;
97082bf6240SBarry Smith     b->icol             = isicol;
971899cda47SBarry Smith     ierr                = PetscMalloc(((*fact)->rmap.n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
972f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
973c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
974c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9753a40ed3dSBarry Smith     PetscFunctionReturn(0);
97608480c60SBarry Smith   }
97708480c60SBarry Smith 
978898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
979898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9809e25ed09SBarry Smith 
9819e25ed09SBarry Smith   /* get new row pointers */
9825a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9835a8e39fbSHong Zhang   bi[0] = 0;
9845a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9855a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9865a8e39fbSHong Zhang   bdiag[0]  = 0;
9876cf997b0SBarry Smith 
9885a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9895a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9905a8e39fbSHong Zhang 
9915a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9925a8e39fbSHong Zhang   nlnk = n + 1;
9935a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9945a8e39fbSHong Zhang 
9955a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
996a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9975a8e39fbSHong Zhang   current_space = free_space;
998a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9995a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
10005a8e39fbSHong Zhang 
10015a8e39fbSHong Zhang   for (i=0; i<n; i++) {
10025a8e39fbSHong Zhang     nzi = 0;
10036cf997b0SBarry Smith     /* copy current row into linked list */
10045a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
10055a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
10065a8e39fbSHong Zhang     cols = aj + ai[r[i]];
10075a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
10085a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
10095a8e39fbSHong Zhang     nzi += nlnk;
10106cf997b0SBarry Smith 
10116cf997b0SBarry Smith     /* make sure diagonal entry is included */
10125a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
10136cf997b0SBarry Smith       fm = n;
10145a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
10155a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
10165a8e39fbSHong Zhang       lnk[fm]    = i;
10175a8e39fbSHong Zhang       lnk_lvl[i] = 0;
10185a8e39fbSHong Zhang       nzi++; dcount++;
10196cf997b0SBarry Smith     }
10206cf997b0SBarry Smith 
10215a8e39fbSHong Zhang     /* add pivot rows into the active row */
10225a8e39fbSHong Zhang     nzbd = 0;
10235a8e39fbSHong Zhang     prow = lnk[n];
10245a8e39fbSHong Zhang     while (prow < i) {
10255a8e39fbSHong Zhang       nnz      = bdiag[prow];
10265a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
10275a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
10285a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
10295a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
10305a8e39fbSHong Zhang       nzi += nlnk;
10315a8e39fbSHong Zhang       prow = lnk[prow];
10325a8e39fbSHong Zhang       nzbd++;
103356beaf04SBarry Smith     }
10345a8e39fbSHong Zhang     bdiag[i] = nzbd;
10355a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1036ecf371e4SBarry Smith 
10375a8e39fbSHong Zhang     /* if free space is not available, make more free space */
10385a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
10395a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1040a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1041a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
10425a8e39fbSHong Zhang       reallocs++;
10435a8e39fbSHong Zhang     }
1044ecf371e4SBarry Smith 
10455a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
10465a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
10475a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
10485a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
10495a8e39fbSHong Zhang 
10505a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
10515a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
105277431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1053d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
10546cf997b0SBarry Smith     }
10555a8e39fbSHong Zhang 
10565a8e39fbSHong Zhang     current_space->array           += nzi;
10575a8e39fbSHong Zhang     current_space->local_used      += nzi;
10585a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
10595a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
10605a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
10615a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
10629e25ed09SBarry Smith   }
10635a8e39fbSHong Zhang 
1064898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1065898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
10665a8e39fbSHong Zhang 
10675a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
10685a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1069a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
10705a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1071a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
10725a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10739e25ed09SBarry Smith 
10746cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1075f64065bbSBarry Smith   {
10765a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1077ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1078ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1079ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1080ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1081335d9088SBarry Smith     if (diagonal_fill) {
1082ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1083335d9088SBarry Smith     }
1084f64065bbSBarry Smith   }
108563ba0a88SBarry Smith #endif
1086bbb6d6a8SBarry Smith 
10879e25ed09SBarry Smith   /* put together the new matrix */
1088f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
1089f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
1090f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1091ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
109252e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1093416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1094*e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1095*e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
10967c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10975a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1098b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10995a8e39fbSHong Zhang   b->j          = bj;
11005a8e39fbSHong Zhang   b->i          = bi;
11015a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
11025a8e39fbSHong Zhang   b->diag       = bdiag;
1103416022c9SBarry Smith   b->ilen       = 0;
1104416022c9SBarry Smith   b->imax       = 0;
1105416022c9SBarry Smith   b->row        = isrow;
1106416022c9SBarry Smith   b->col        = iscol;
1107c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1108c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
110982bf6240SBarry Smith   b->icol       = isicol;
111087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1111416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
11125a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
11135a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
11145a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
11159e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1116418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1117ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
11185a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
111971c2f376SKris Buschelman 
112054e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
112171c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
11224846f1f5SKris Buschelman 
11233a40ed3dSBarry Smith   PetscFunctionReturn(0);
11249e25ed09SBarry Smith }
1125d5d45c9bSBarry Smith 
11263c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1127a6175056SHong Zhang #undef __FUNCT__
1128a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1129af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1130a6175056SHong Zhang {
11312ed133dbSHong Zhang   Mat            C = *B;
11320968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11332ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11342ed133dbSHong Zhang   IS             ip=b->row;
11352ed133dbSHong Zhang   PetscErrorCode ierr;
1136899cda47SBarry Smith   PetscInt       *rip,i,j,mbs=A->rmap.n,*bi=b->i,*bj=b->j,*bcol;
11372ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11382ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1139622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1140540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1141fbf22428SSatish Balay   PetscReal      shiftpd;
1142540703acSHong Zhang   ChShift_Ctx    sctx;
1143540703acSHong Zhang   PetscInt       newshift;
1144d5d45c9bSBarry Smith 
1145a6175056SHong Zhang   PetscFunctionBegin;
1146540703acSHong Zhang   shiftnz   = info->shiftnz;
1147540703acSHong Zhang   shiftpd   = info->shiftpd;
1148ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1149ee45ca4aSHong Zhang 
11502ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11512ed133dbSHong Zhang 
11522ed133dbSHong Zhang   /* initialization */
11532ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11542ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11552ed133dbSHong Zhang   jl   = il + mbs;
11562ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11572ed133dbSHong Zhang 
1158540703acSHong Zhang   sctx.shift_amount = 0;
1159540703acSHong Zhang   sctx.nshift       = 0;
11602ed133dbSHong Zhang   do {
1161540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
11622ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11632ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1164a6175056SHong Zhang     }
11652ed133dbSHong Zhang 
11662ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1167c05c3958SHong Zhang       bval = ba + bi[k];
11682ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11692ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11702ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11712ed133dbSHong Zhang         col = rip[aj[j]];
11722ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11732ed133dbSHong Zhang           rtmp[col] = aa[j];
1174c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11752ed133dbSHong Zhang         }
11762ed133dbSHong Zhang       }
117739e3d630SHong Zhang       /* shift the diagonal of the matrix */
1178540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11792ed133dbSHong Zhang 
11802ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11812ed133dbSHong Zhang       dk = rtmp[k];
11822ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11832ed133dbSHong Zhang 
11842ed133dbSHong Zhang       while (i < k){
11852ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11862ed133dbSHong Zhang 
11872ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11882ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11892ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11902ed133dbSHong Zhang         dk += uikdi*ba[ili];
11912ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11922ed133dbSHong Zhang 
11932ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11942ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11952ed133dbSHong Zhang         if (jmin < jmax){
11962ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11972ed133dbSHong Zhang           /* update il and jl for row i */
11982ed133dbSHong Zhang           il[i] = jmin;
11992ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
12002ed133dbSHong Zhang         }
12012ed133dbSHong Zhang         i = nexti;
12022ed133dbSHong Zhang       }
12032ed133dbSHong Zhang 
12040697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
12050697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12060697b6caSHong Zhang       rs   = 0.0;
12073c9e8598SHong Zhang       jmin = bi[k]+1;
12083c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
12093c9e8598SHong Zhang       if (nz){
12103c9e8598SHong Zhang         bcol = bj + jmin;
12113c9e8598SHong Zhang         while (nz--){
121289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
121389f845c8SHong Zhang           bcol++;
12143c9e8598SHong Zhang         }
12153c9e8598SHong Zhang       }
1216540703acSHong Zhang 
1217540703acSHong Zhang       sctx.rs = rs;
1218540703acSHong Zhang       sctx.pv = dk;
12195b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
12205b5bc046SBarry Smith       if (newshift == 1) break;
12213c9e8598SHong Zhang 
12223c9e8598SHong Zhang       /* copy data into U(k,:) */
122339e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
122439e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
122539e3d630SHong Zhang       if (jmin < jmax) {
122639e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
122739e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12283c9e8598SHong Zhang         }
122939e3d630SHong Zhang         /* add the k-th row into il and jl */
12303c9e8598SHong Zhang         il[k] = jmin;
12313c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12323c9e8598SHong Zhang       }
12333c9e8598SHong Zhang     }
1234540703acSHong Zhang   } while (sctx.chshift);
12353c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12363c9e8598SHong Zhang 
123739e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12383c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12393c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12403c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1241899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.n);CHKERRQ(ierr);
1242540703acSHong Zhang   if (sctx.nshift){
1243540703acSHong Zhang     if (shiftnz) {
1244ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1245540703acSHong Zhang     } else if (shiftpd) {
1246ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
12470697b6caSHong Zhang     }
12483c9e8598SHong Zhang   }
12493c9e8598SHong Zhang   PetscFunctionReturn(0);
12503c9e8598SHong Zhang }
1251a6175056SHong Zhang 
1252a6175056SHong Zhang #undef __FUNCT__
1253a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1254dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1255a6175056SHong Zhang {
12560968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1257ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1258ed59401aSHong Zhang   Mat                B;
1259dfbe8321SBarry Smith   PetscErrorCode     ierr;
1260653a6975SHong Zhang   PetscTruth         perm_identity;
1261899cda47SBarry Smith   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->rmap.n,*ui;
1262ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1263391eac55SHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12645a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1265ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1266a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1267a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12687a48dd6fSHong Zhang   PetscBT            lnkbt;
1269a6175056SHong Zhang 
1270a6175056SHong Zhang   PetscFunctionBegin;
1271653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12727a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12737a48dd6fSHong Zhang 
127439e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
127539e3d630SHong Zhang   ui[0] = 0;
127639e3d630SHong Zhang 
1277622e2028SHong Zhang   /* special case that simply copies fill pattern */
1278622e2028SHong Zhang   if (!levels && perm_identity) {
1279622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1280ed59401aSHong Zhang     for (i=0; i<am; i++) {
128139e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1282ed59401aSHong Zhang     }
128339e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
128439e3d630SHong Zhang     cols = uj;
1285ed59401aSHong Zhang     for (i=0; i<am; i++) {
1286ed59401aSHong Zhang       aj    = a->j + a->diag[i];
128739e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
128839e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1289ed59401aSHong Zhang     }
129039e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12917a48dd6fSHong Zhang     /* initialization */
12925a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12937a48dd6fSHong Zhang 
12947a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12957a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12961393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12977a48dd6fSHong Zhang     il         = jl + am;
12987a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12997a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
13007a48dd6fSHong Zhang     for (i=0; i<am; i++){
13017a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
13027a48dd6fSHong Zhang     }
13037a48dd6fSHong Zhang 
13047a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
13057a48dd6fSHong Zhang     nlnk = am + 1;
13062ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13077a48dd6fSHong Zhang 
13087a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1309a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
13107a48dd6fSHong Zhang     current_space = free_space;
1311a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
13127a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
13137a48dd6fSHong Zhang 
13147a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
13157a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
13167a48dd6fSHong Zhang       nzk   = 0;
1317622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1318622e2028SHong Zhang       ncols_upper = 0;
1319622e2028SHong Zhang       for (j=0; j<ncols; j++){
13205a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
13215a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
13225a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1323622e2028SHong Zhang           ncols_upper++;
1324622e2028SHong Zhang         }
1325622e2028SHong Zhang       }
13265a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13277a48dd6fSHong Zhang       nzk += nlnk;
13287a48dd6fSHong Zhang 
13297a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13307a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13317a48dd6fSHong Zhang 
13327a48dd6fSHong Zhang       while (prow < k){
13337a48dd6fSHong Zhang         nextprow = jl[prow];
13347a48dd6fSHong Zhang 
13357a48dd6fSHong Zhang         /* merge prow into k-th row */
13367a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13377a48dd6fSHong Zhang         jmax = ui[prow+1];
13387a48dd6fSHong Zhang         ncols = jmax-jmin;
1339ed59401aSHong Zhang         i     = jmin - ui[prow];
1340ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
13415a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
13425a8e39fbSHong Zhang         j     = *(uj - 1);
13435a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
13447a48dd6fSHong Zhang         nzk += nlnk;
13457a48dd6fSHong Zhang 
13467a48dd6fSHong Zhang         /* update il and jl for prow */
13477a48dd6fSHong Zhang         if (jmin < jmax){
13487a48dd6fSHong Zhang           il[prow] = jmin;
13497a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13507a48dd6fSHong Zhang         }
13517a48dd6fSHong Zhang         prow = nextprow;
13527a48dd6fSHong Zhang       }
13537a48dd6fSHong Zhang 
13547a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13557a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13567a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13577a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1358a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1359a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
13607a48dd6fSHong Zhang         reallocs++;
13617a48dd6fSHong Zhang       }
13627a48dd6fSHong Zhang 
13637a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13642ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13657a48dd6fSHong Zhang 
13667a48dd6fSHong Zhang       /* add the k-th row into il and jl */
136739e3d630SHong Zhang       if (nzk > 1){
13687a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13697a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13707a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13717a48dd6fSHong Zhang       }
13727a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13737a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13747a48dd6fSHong Zhang 
13757a48dd6fSHong Zhang       current_space->array           += nzk;
13767a48dd6fSHong Zhang       current_space->local_used      += nzk;
13777a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13787a48dd6fSHong Zhang 
13797a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13807a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13817a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13827a48dd6fSHong Zhang 
13837a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13847a48dd6fSHong Zhang     }
13857a48dd6fSHong Zhang 
13866cf91177SBarry Smith #if defined(PETSC_USE_INFO)
13877a48dd6fSHong Zhang     if (ai[am] != 0) {
138839e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1389ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1390ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1391ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
13927a48dd6fSHong Zhang     } else {
1393ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
13947a48dd6fSHong Zhang     }
139563ba0a88SBarry Smith #endif
13967a48dd6fSHong Zhang 
13977a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13987a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13995a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
14007a48dd6fSHong Zhang 
14017a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
14027a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1403a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
14042ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1405a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
14067a48dd6fSHong Zhang 
140739e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
140839e3d630SHong Zhang 
14097a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1410f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1411f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
1412ed59401aSHong Zhang   B = *fact;
1413ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1414ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
14157a48dd6fSHong Zhang 
1416ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
14177a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
14187a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
14197a48dd6fSHong Zhang   b->j    = uj;
14207a48dd6fSHong Zhang   b->i    = ui;
14217a48dd6fSHong Zhang   b->diag = 0;
14227a48dd6fSHong Zhang   b->ilen = 0;
14237a48dd6fSHong Zhang   b->imax = 0;
14247a48dd6fSHong Zhang   b->row  = perm;
14257a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
14267a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14277a48dd6fSHong Zhang   b->icol = perm;
14287a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14297a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
143052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
14317a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
1432*e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1433*e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
14347a48dd6fSHong Zhang 
1435ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1436ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1437ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14387a48dd6fSHong Zhang   if (ai[am] != 0) {
1439ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14407a48dd6fSHong Zhang   } else {
1441ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14427a48dd6fSHong Zhang   }
144339e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1444622e2028SHong Zhang   if (perm_identity){
1445ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1446ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1447622e2028SHong Zhang   }
1448a6175056SHong Zhang   PetscFunctionReturn(0);
1449a6175056SHong Zhang }
1450d5d45c9bSBarry Smith 
1451f76d2b81SHong Zhang #undef __FUNCT__
1452f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1453dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1454f76d2b81SHong Zhang {
1455f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
145610c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
14572ed133dbSHong Zhang   Mat                B;
1458dfbe8321SBarry Smith   PetscErrorCode     ierr;
1459f76d2b81SHong Zhang   PetscTruth         perm_identity;
146010c27e3dSHong Zhang   PetscReal          fill = info->fill;
1461899cda47SBarry Smith   PetscInt           *rip,*riip,i,am=A->rmap.n,*ai=a->i,*aj=a->j,reallocs=0,prow;
146210c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14632ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1464a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
146510c27e3dSHong Zhang   PetscBT            lnkbt;
14662ed133dbSHong Zhang   IS                 iperm;
1467f76d2b81SHong Zhang 
1468f76d2b81SHong Zhang   PetscFunctionBegin;
14692ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1470f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14712ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14722ed133dbSHong Zhang 
1473f76d2b81SHong Zhang   if (!perm_identity){
14742ed133dbSHong Zhang     /* check if perm is symmetric! */
14752ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14762ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14771393bc97SHong Zhang     for (i=0; i<am; i++) {
14782ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14792ed133dbSHong Zhang     }
14802ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14812ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1482f76d2b81SHong Zhang   }
148310c27e3dSHong Zhang 
148410c27e3dSHong Zhang   /* initialization */
14851393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
148610c27e3dSHong Zhang   ui[0] = 0;
148710c27e3dSHong Zhang 
148810c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14891393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14901393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14911393bc97SHong Zhang   il     = jl + am;
14921393bc97SHong Zhang   cols   = il + am;
14931393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14941393bc97SHong Zhang   for (i=0; i<am; i++){
14951393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1496f76d2b81SHong Zhang   }
149710c27e3dSHong Zhang 
149810c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14991393bc97SHong Zhang   nlnk = am + 1;
15001393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150110c27e3dSHong Zhang 
15021393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1503a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
150410c27e3dSHong Zhang   current_space = free_space;
150510c27e3dSHong Zhang 
15061393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
150710c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
150810c27e3dSHong Zhang     nzk   = 0;
15092ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
15102ed133dbSHong Zhang     ncols_upper = 0;
15112ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1512622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
15132ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
15142ed133dbSHong Zhang         cols[ncols_upper] = i;
15152ed133dbSHong Zhang         ncols_upper++;
15162ed133dbSHong Zhang       }
15172ed133dbSHong Zhang     }
15181393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
151910c27e3dSHong Zhang     nzk += nlnk;
152010c27e3dSHong Zhang 
152110c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
152210c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
152310c27e3dSHong Zhang 
152410c27e3dSHong Zhang     while (prow < k){
152510c27e3dSHong Zhang       nextprow = jl[prow];
152610c27e3dSHong Zhang       /* merge prow into k-th row */
15271393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
152810c27e3dSHong Zhang       jmax = ui[prow+1];
152910c27e3dSHong Zhang       ncols = jmax-jmin;
15301393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15315a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
153210c27e3dSHong Zhang       nzk += nlnk;
153310c27e3dSHong Zhang 
153410c27e3dSHong Zhang       /* update il and jl for prow */
153510c27e3dSHong Zhang       if (jmin < jmax){
153610c27e3dSHong Zhang         il[prow] = jmin;
15372ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
153810c27e3dSHong Zhang       }
153910c27e3dSHong Zhang       prow = nextprow;
154010c27e3dSHong Zhang     }
154110c27e3dSHong Zhang 
154210c27e3dSHong Zhang     /* if free space is not available, make more free space */
154310c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15441393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
154510c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1546a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
154710c27e3dSHong Zhang       reallocs++;
154810c27e3dSHong Zhang     }
154910c27e3dSHong Zhang 
155010c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15511393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
155210c27e3dSHong Zhang 
155310c27e3dSHong Zhang     /* add the k-th row into il and jl */
155410c27e3dSHong Zhang     if (nzk-1 > 0){
15551393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
155610c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
155710c27e3dSHong Zhang       il[k] = ui[k] + 1;
155810c27e3dSHong Zhang     }
15592ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
156010c27e3dSHong Zhang     current_space->array           += nzk;
156110c27e3dSHong Zhang     current_space->local_used      += nzk;
156210c27e3dSHong Zhang     current_space->local_remaining -= nzk;
156310c27e3dSHong Zhang 
156410c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
156510c27e3dSHong Zhang   }
156610c27e3dSHong Zhang 
15676cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15681393bc97SHong Zhang   if (ai[am] != 0) {
15691393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1570ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1571ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1572ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
157310c27e3dSHong Zhang   } else {
1574ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
157510c27e3dSHong Zhang   }
157663ba0a88SBarry Smith #endif
157710c27e3dSHong Zhang 
157810c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
157910c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
158010c27e3dSHong Zhang 
158110c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15821393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1583a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
158410c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
158510c27e3dSHong Zhang 
158610c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1587f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1588f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
15892ed133dbSHong Zhang   B    = *fact;
15902ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1591ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
159210c27e3dSHong Zhang 
15932ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
159410c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1595*e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1596*e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
15971393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
159810c27e3dSHong Zhang   b->j    = uj;
159910c27e3dSHong Zhang   b->i    = ui;
160010c27e3dSHong Zhang   b->diag = 0;
160110c27e3dSHong Zhang   b->ilen = 0;
160210c27e3dSHong Zhang   b->imax = 0;
160310c27e3dSHong Zhang   b->row  = perm;
160410c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
160510c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
160610c27e3dSHong Zhang   b->icol = perm;
160710c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
16081393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
16091393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
16101393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
161110c27e3dSHong Zhang 
16122ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
16132ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
16142ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
16151393bc97SHong Zhang   if (ai[am] != 0) {
16161393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
161710c27e3dSHong Zhang   } else {
16182ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
161910c27e3dSHong Zhang   }
162039e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
16212ed133dbSHong Zhang   if (perm_identity){
162210c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
162310c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
16242ed133dbSHong Zhang   }
1625f76d2b81SHong Zhang   PetscFunctionReturn(0);
1626f76d2b81SHong Zhang }
1627