xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 356650c2cadef3e29157f31f7bffbb31d8b16dd8)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2289bc588SBarry Smith 
370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
41c4feecaSSatish Balay #include "src/inline/dot.h"
5ed480e8bSBarry Smith #include "src/inline/spops.h"
61393bc97SHong Zhang #include "petscbt.h"
71393bc97SHong Zhang #include "src/mat/utils/freespace.h"
83b2fbd54SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
11dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
123b2fbd54SBarry Smith {
133a40ed3dSBarry Smith   PetscFunctionBegin;
143a40ed3dSBarry Smith 
1529bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
16aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
17d64ed03dSBarry Smith   PetscFunctionReturn(0);
18d64ed03dSBarry Smith #endif
193b2fbd54SBarry Smith }
203b2fbd54SBarry Smith 
2186bacbd2SBarry Smith 
22e34fafa9SBarry Smith #if !defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
23a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
24a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
25a77337e4SBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,MatScalar*,PetscInt*,MatScalar*,PetscInt*,PetscInt*,MatScalar*,PetscInt*);
26e34fafa9SBarry Smith #endif
2786bacbd2SBarry Smith 
284a2ae208SSatish Balay #undef __FUNCT__
294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
3086bacbd2SBarry Smith   /* ------------------------------------------------------------
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith           This interface was contribed by Tony Caola
3386bacbd2SBarry Smith 
3486bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
355bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
365bd2ddc7SBarry Smith      SPARSEKIT2.
375bd2ddc7SBarry Smith 
385bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
395bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
4086bacbd2SBarry Smith 
4186bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4286bacbd2SBarry Smith      help in getting this routine ironed out.
4386bacbd2SBarry Smith 
445bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
455bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
465bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
475bd2ddc7SBarry Smith      Yousef Saad's code.
4886bacbd2SBarry Smith 
4986bacbd2SBarry Smith      ------------------------------------------------------------
5086bacbd2SBarry Smith */
517529efd4SKris Buschelman PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
5286bacbd2SBarry Smith {
5360ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5460ec86bdSBarry Smith   PetscFunctionBegin;
55e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5660ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5760ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5860ec86bdSBarry Smith #else
5986bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
6007d2056aSBarry Smith   IS             iscolf,isicol,isirow;
6186bacbd2SBarry Smith   PetscTruth     reorder;
629cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
63899cda47SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->rmap.n;
6438baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6538baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6638baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6754f21887SBarry Smith   MatScalar      *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6854a8161fSBarry Smith   PetscReal      af;
6986bacbd2SBarry Smith 
7086bacbd2SBarry Smith   PetscFunctionBegin;
7186bacbd2SBarry Smith 
7286bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7338baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7486bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7533259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7638baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7738baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7886bacbd2SBarry Smith 
7986bacbd2SBarry Smith 
8086bacbd2SBarry Smith   /* ------------------------------------------------------------
8186bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8286bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8386bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8486bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8586bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8686bacbd2SBarry Smith      ------------------------------------------------------------
8786bacbd2SBarry Smith   */
8886bacbd2SBarry Smith 
8986bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
9086bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
9186bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9286bacbd2SBarry Smith   reorder = PetscNot(reorder);
9386bacbd2SBarry Smith 
9486bacbd2SBarry Smith 
9586bacbd2SBarry Smith   /* storage for ilu factor */
9638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9738baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
98a77337e4SBarry Smith   ierr = PetscMalloc(jmax*sizeof(MatScalar),&new_a);CHKERRQ(ierr);
9938baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
10086bacbd2SBarry Smith 
10186bacbd2SBarry Smith   /* ------------------------------------------------------------
10286bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10386bacbd2SBarry Smith      ------------------------------------------------------------
10486bacbd2SBarry Smith   */
105b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
106b19713cbSBarry Smith     old_j[i]++;
10786bacbd2SBarry Smith   }
108b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
109b19713cbSBarry Smith     old_i[i]++;
110b19713cbSBarry Smith   };
111010a20caSHong Zhang 
11286bacbd2SBarry Smith 
11386bacbd2SBarry Smith   if (reorder) {
114c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
115c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
116c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
117c0c2c81eSBarry Smith       r[i]  = r[i]+1;
118c0c2c81eSBarry Smith       c[i]  = c[i]+1;
119c0c2c81eSBarry Smith     }
12038baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
12138baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
122a77337e4SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(MatScalar),&old_a2);CHKERRQ(ierr);
1235bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
124c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
125c0c2c81eSBarry Smith       r[i]  = r[i]-1;
126c0c2c81eSBarry Smith       c[i]  = c[i]-1;
127c0c2c81eSBarry Smith     }
128c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
129c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
130b19713cbSBarry Smith     o_a = old_a2;
131b19713cbSBarry Smith     o_j = old_j2;
132b19713cbSBarry Smith     o_i = old_i2;
133b19713cbSBarry Smith   } else {
134b19713cbSBarry Smith     o_a = old_a;
135b19713cbSBarry Smith     o_j = old_j;
136b19713cbSBarry Smith     o_i = old_i;
13786bacbd2SBarry Smith   }
13886bacbd2SBarry Smith 
13986bacbd2SBarry Smith   /* ------------------------------------------------------------
14086bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
14186bacbd2SBarry Smith      ------------------------------------------------------------
14286bacbd2SBarry Smith   */
14386bacbd2SBarry Smith 
14438baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14538baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14687828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14786bacbd2SBarry Smith 
14854a8161fSBarry Smith   SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&info->dtcol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&sierr);
1496849ba73SBarry Smith   if (sierr) {
1506849ba73SBarry Smith     switch (sierr) {
151a83599f4SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
152a83599f4SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %G space allocated %D",info->fill,jmax);
153e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
154e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15577431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15677431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15786bacbd2SBarry Smith     }
15886bacbd2SBarry Smith   }
15986bacbd2SBarry Smith 
16086bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
16186bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16286bacbd2SBarry Smith 
16386bacbd2SBarry Smith   /* ------------------------------------------------------------
16486bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16586bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16686bacbd2SBarry Smith      ------------------------------------------------------------
16786bacbd2SBarry Smith   */
16886bacbd2SBarry Smith 
16987828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
17038baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
17186bacbd2SBarry Smith 
1725bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17386bacbd2SBarry Smith 
17486bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17586bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17686bacbd2SBarry Smith 
17786bacbd2SBarry Smith   if (reorder) {
17886bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17986bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
18086bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
181313ae024SBarry Smith   } else {
182b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
183f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
184b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
185b19713cbSBarry Smith     }
186b19713cbSBarry Smith   }
187b19713cbSBarry Smith 
188b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18986bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
190b19713cbSBarry Smith     old_i[i]--;
191b19713cbSBarry Smith   }
192b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
193b19713cbSBarry Smith     old_j[i]--;
194b19713cbSBarry Smith   }
19586bacbd2SBarry Smith 
196b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
198b19713cbSBarry Smith     new_i[i]--;
199b19713cbSBarry Smith   }
200b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
201b19713cbSBarry Smith     new_j[i]--;
202b19713cbSBarry Smith   }
20386bacbd2SBarry Smith 
20486bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20586bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2064c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2074c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
208c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20986bacbd2SBarry Smith   for(i=0; i<n; i++) {
21086bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
21186bacbd2SBarry Smith   };
21286bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21307d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21486bacbd2SBarry Smith 
215c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
216c0c2c81eSBarry Smith 
217329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21807d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21986bacbd2SBarry Smith 
22086bacbd2SBarry Smith   /*----- put together the new matrix -----*/
22186bacbd2SBarry Smith 
2227adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
223f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
2247adad957SLisandro Dalcin   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
225ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
22686bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22786bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22886bacbd2SBarry Smith 
22986bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
230e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
231e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
23207d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23386bacbd2SBarry Smith   b->a             = new_a;
23486bacbd2SBarry Smith   b->j             = new_j;
23586bacbd2SBarry Smith   b->i             = new_i;
23686bacbd2SBarry Smith   b->ilen          = 0;
23786bacbd2SBarry Smith   b->imax          = 0;
2381f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
239313ae024SBarry Smith   b->row           = isirow;
24086bacbd2SBarry Smith   b->col           = iscolf;
24187828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24286bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24386bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24486bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24586bacbd2SBarry Smith 
246431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
247ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
248ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
250ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
251431e710aSBarry Smith 
2527529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25371c2f376SKris Buschelman 
25486bacbd2SBarry Smith   PetscFunctionReturn(0);
25560ec86bdSBarry Smith #endif
25686bacbd2SBarry Smith }
25786bacbd2SBarry Smith 
2584a2ae208SSatish Balay #undef __FUNCT__
259b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
260dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
261289bc588SBarry Smith {
262416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
263289bc588SBarry Smith   IS                 isicol;
2646849ba73SBarry Smith   PetscErrorCode     ierr;
265899cda47SBarry Smith   PetscInt           *r,*ic,i,n=A->rmap.n,*ai=a->i,*aj=a->j;
2661393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
2671393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2689e046878SBarry Smith   PetscReal          f;
2694875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
270a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2711393bc97SHong Zhang   PetscBT            lnkbt;
272289bc588SBarry Smith 
2733a40ed3dSBarry Smith   PetscFunctionBegin;
274899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2754c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
276ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
277ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
278289bc588SBarry Smith 
279289bc588SBarry Smith   /* get new row pointers */
2801393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2811393bc97SHong Zhang   bi[0] = 0;
2821393bc97SHong Zhang 
2831393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2841393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2851393bc97SHong Zhang   bdiag[0] = 0;
2861393bc97SHong Zhang 
2871393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2881393bc97SHong Zhang   nlnk = n + 1;
2891393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2901393bc97SHong Zhang 
29135baf27eSSatish Balay   ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr);
2921393bc97SHong Zhang 
2931393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
294d3d32019SBarry Smith   f = info->fill;
295a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
2961393bc97SHong Zhang   current_space = free_space;
297289bc588SBarry Smith 
298289bc588SBarry Smith   for (i=0; i<n; i++) {
2991393bc97SHong Zhang     /* copy previous fill into linked list */
300289bc588SBarry Smith     nzi = 0;
3011393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3021393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3031393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
304afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3051393bc97SHong Zhang     nzi += nlnk;
3061393bc97SHong Zhang 
3071393bc97SHong Zhang     /* add pivot rows into linked list */
3081393bc97SHong Zhang     row = lnk[n];
3091393bc97SHong Zhang     while (row < i) {
310d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
311d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
312d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3131393bc97SHong Zhang       nzi += nlnk;
3141393bc97SHong Zhang       row  = lnk[row];
315289bc588SBarry Smith     }
3161393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3171393bc97SHong Zhang     im[i]   = nzi;
3181393bc97SHong Zhang 
3191393bc97SHong Zhang     /* mark bdiag */
3201393bc97SHong Zhang     nzbd = 0;
3211393bc97SHong Zhang     nnz  = nzi;
3221393bc97SHong Zhang     k    = lnk[n];
3231393bc97SHong Zhang     while (nnz-- && k < i){
3241393bc97SHong Zhang       nzbd++;
3251393bc97SHong Zhang       k = lnk[k];
3261393bc97SHong Zhang     }
3271393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3281393bc97SHong Zhang 
3291393bc97SHong Zhang     /* if free space is not available, make more free space */
3301393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3311393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
332a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3331393bc97SHong Zhang       reallocs++;
3341393bc97SHong Zhang     }
3351393bc97SHong Zhang 
3361393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3371393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3381393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3391393bc97SHong Zhang     current_space->array           += nzi;
3401393bc97SHong Zhang     current_space->local_used      += nzi;
3411393bc97SHong Zhang     current_space->local_remaining -= nzi;
342289bc588SBarry Smith   }
3436cf91177SBarry Smith #if defined(PETSC_USE_INFO)
344464e8d44SSatish Balay   if (ai[n] != 0) {
3451393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
346ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
347ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
348ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
349ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
35005bf559cSSatish Balay   } else {
351ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
35205bf559cSSatish Balay   }
35363ba0a88SBarry Smith #endif
3542fb3ed76SBarry Smith 
355898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
356898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3571987afe7SBarry Smith 
3581393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3591393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
360a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3611393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
36235baf27eSSatish Balay   ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr);
363289bc588SBarry Smith 
364289bc588SBarry Smith   /* put together the new matrix */
3657adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
366f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
3677adad957SLisandro Dalcin   ierr = MatSetType(*B,((PetscObject)A)->type_name);CHKERRQ(ierr);
368ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
36952e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
370416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
371e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
372e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
3737c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3741393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3751393bc97SHong Zhang   b->j          = bj;
3761393bc97SHong Zhang   b->i          = bi;
3771393bc97SHong Zhang   b->diag       = bdiag;
378416022c9SBarry Smith   b->ilen       = 0;
379416022c9SBarry Smith   b->imax       = 0;
380416022c9SBarry Smith   b->row        = isrow;
381416022c9SBarry Smith   b->col        = iscol;
382c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
383c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
38482bf6240SBarry Smith   b->icol       = isicol;
38587828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3861393bc97SHong Zhang 
3871393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3881393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3891393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
3907fd98bd6SLois Curfman McInnes 
391329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
392418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
393ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
394703deb49SSatish Balay 
395e93ccc4dSBarry Smith   if (ai[n] != 0) {
3961393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
397eea2913aSSatish Balay   } else {
398eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
399eea2913aSSatish Balay   }
4004846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
40171c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
403289bc588SBarry Smith }
4041393bc97SHong Zhang 
4055b5bc046SBarry Smith /*
4065b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4075b5bc046SBarry Smith */
4085b5bc046SBarry Smith #undef __FUNCT__
4095b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4105b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4115b5bc046SBarry Smith {
4125b5bc046SBarry Smith   PetscErrorCode ierr;
4135b5bc046SBarry Smith   PetscTruth     flg;
4145b5bc046SBarry Smith 
4155b5bc046SBarry Smith   PetscFunctionBegin;
4165b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4175b5bc046SBarry Smith   if (flg) {
4185b5bc046SBarry Smith     PetscViewer viewer;
4195b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4205b5bc046SBarry Smith 
4215b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4227adad957SLisandro Dalcin     ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4235b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4245b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4255b5bc046SBarry Smith   }
4265b5bc046SBarry Smith   PetscFunctionReturn(0);
4275b5bc046SBarry Smith }
4285b5bc046SBarry Smith 
4290a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4304a2ae208SSatish Balay #undef __FUNCT__
4314a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
432af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
433289bc588SBarry Smith {
434416022c9SBarry Smith   Mat            C=*B;
435416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
43682bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4376849ba73SBarry Smith   PetscErrorCode ierr;
438899cda47SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap.n,*bi=b->i,*bj=b->j;
439b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
440d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
44154f21887SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,*rtmps;
44254f21887SBarry Smith   MatScalar      *v,*pv;
4436a7c8fc2SHong Zhang   PetscScalar    d;
4446a7c8fc2SHong Zhang   PetscReal      rs;
445b3bf805bSHong Zhang   LUShift_Ctx    sctx;
44642898933SHong Zhang   PetscInt       newshift,*ddiag;
447289bc588SBarry Smith 
4483a40ed3dSBarry Smith   PetscFunctionBegin;
44978b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
45078b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
45187828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45287828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
453010a20caSHong Zhang   rtmps = rtmp; ics = ic;
454289bc588SBarry Smith 
455ada7143aSSatish Balay   sctx.shift_top  = 0;
456ada7143aSSatish Balay   sctx.nshift_max = 0;
457ada7143aSSatish Balay   sctx.shift_lo   = 0;
458ada7143aSSatish Balay   sctx.shift_hi   = 0;
459ada7143aSSatish Balay 
460118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
461118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4621a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
4639deb78dcSHong Zhang     PetscInt *aai = a->i;
46442898933SHong Zhang     ddiag          = a->diag;
465118f5789SBarry Smith     sctx.shift_top = 0;
4666cc28720Svictorle     for (i=0; i<n; i++) {
4679f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4689f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4691a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
470010a20caSHong Zhang       v  = a->a+aai[i];
471e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
472e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4731a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4741a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4756cc28720Svictorle     }
476c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
477118f5789SBarry Smith     sctx.shift_top    *= 1.1;
478d2276718SHong Zhang     sctx.nshift_max   = 5;
479d2276718SHong Zhang     sctx.shift_lo     = 0.;
480d2276718SHong Zhang     sctx.shift_hi     = 1.;
481a0a356efSHong Zhang   }
482a0a356efSHong Zhang 
483a0a356efSHong Zhang   sctx.shift_amount = 0;
484a0a356efSHong Zhang   sctx.nshift       = 0;
485085256b3SBarry Smith   do {
486d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
487289bc588SBarry Smith     for (i=0; i<n; i++){
488b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
489b3bf805bSHong Zhang       bjtmp = bj + bi[i];
490b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
491289bc588SBarry Smith 
492289bc588SBarry Smith       /* load in initial (unfactored row) */
493416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
494b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
495010a20caSHong Zhang       v     = a->a + a->i[r[i]];
496085256b3SBarry Smith       for (j=0; j<nz; j++) {
497b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
498085256b3SBarry Smith       }
499d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
500289bc588SBarry Smith 
501b3bf805bSHong Zhang       row = *bjtmp++;
502f2582111SSatish Balay       while  (row < i) {
5038c37ef55SBarry Smith         pc = rtmp + row;
5048c37ef55SBarry Smith         if (*pc != 0.0) {
505010a20caSHong Zhang           pv         = b->a + diag_offset[row];
506010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50735aab85fSBarry Smith           multiplier = *pc / *pv++;
5088c37ef55SBarry Smith           *pc        = multiplier;
509b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
510f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
511efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
512289bc588SBarry Smith         }
513b3bf805bSHong Zhang         row = *bjtmp++;
514289bc588SBarry Smith       }
515416022c9SBarry Smith       /* finished row so stick it into b->a */
516b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
517b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
518b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
519b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
520d3d32019SBarry Smith       rs   = 0.0;
521d3d32019SBarry Smith       for (j=0; j<nz; j++) {
522d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
523d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
524d3d32019SBarry Smith       }
525d2276718SHong Zhang 
526b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
527d2276718SHong Zhang       sctx.rs  = rs;
528d2276718SHong Zhang       sctx.pv  = pv[diag];
529c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
5305b5bc046SBarry Smith       if (newshift == 1) break;
53135aab85fSBarry Smith     }
532d2276718SHong Zhang 
533118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5346cc28720Svictorle       /*
5356c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5366cc28720Svictorle        * then try lower shift
5376cc28720Svictorle        */
538d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
539d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
540d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
541d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
542d2276718SHong Zhang       sctx.nshift++;
5436cc28720Svictorle     }
544d2276718SHong Zhang   } while (sctx.lushift);
5450f11f9aeSSatish Balay 
54635aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
54735aab85fSBarry Smith   for (i=0; i<n; i++) {
548010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
54935aab85fSBarry Smith   }
55035aab85fSBarry Smith 
551606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55278b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55378b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
554416022c9SBarry Smith   C->factor = FACTOR_LU;
5557dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
556c456f294SBarry Smith   C->assembled = PETSC_TRUE;
557899cda47SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
558d2276718SHong Zhang   if (sctx.nshift){
559118f5789SBarry Smith     if (info->shiftnz) {
5601e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
561118f5789SBarry Smith     } else if (info->shiftpd) {
5621e2582c4SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
5636cc28720Svictorle     }
5640697b6caSHong Zhang   }
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
566289bc588SBarry Smith }
5670889b5dcSvictorle 
568137fb511SHong Zhang /*
569137fb511SHong Zhang    This routine implements inplace ILU(0) with row or/and column permutations.
570137fb511SHong Zhang    Input:
571137fb511SHong Zhang      A - original matrix
572137fb511SHong Zhang    Output;
573137fb511SHong Zhang      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
574137fb511SHong Zhang          a->j (col index) is permuted by the inverse of colperm, then sorted
575137fb511SHong Zhang          a->a reordered accordingly with a->j
576137fb511SHong Zhang          a->diag (ptr to diagonal elements) is updated.
577137fb511SHong Zhang */
578137fb511SHong Zhang #undef __FUNCT__
579137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm"
580137fb511SHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat A,MatFactorInfo *info,Mat *B)
581137fb511SHong Zhang {
582b051339dSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
583b051339dSHong Zhang   IS             isrow = a->row,isicol = a->icol;
584137fb511SHong Zhang   PetscErrorCode ierr;
585b051339dSHong Zhang   PetscInt       *r,*ic,i,j,n=A->rmap.n,*ai=a->i,*aj=a->j;
586b051339dSHong Zhang   PetscInt       *ajtmp,nz,row,*ics;
587b051339dSHong Zhang   PetscInt       *diag = a->diag,nbdiag,*pj;
588a77337e4SBarry Smith   PetscScalar    *rtmp,*pc,multiplier,d;
589a77337e4SBarry Smith   MatScalar      *v,*pv;
590137fb511SHong Zhang   PetscReal      rs;
591137fb511SHong Zhang   LUShift_Ctx    sctx;
592b051339dSHong Zhang   PetscInt       newshift;
593137fb511SHong Zhang 
594137fb511SHong Zhang   PetscFunctionBegin;
595b051339dSHong Zhang   if (A != *B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");
596137fb511SHong Zhang   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
597137fb511SHong Zhang   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
598137fb511SHong Zhang   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
599137fb511SHong Zhang   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
600b051339dSHong Zhang   ics = ic;
601137fb511SHong Zhang 
602137fb511SHong Zhang   sctx.shift_top  = 0;
603137fb511SHong Zhang   sctx.nshift_max = 0;
604137fb511SHong Zhang   sctx.shift_lo   = 0;
605137fb511SHong Zhang   sctx.shift_hi   = 0;
606137fb511SHong Zhang 
607137fb511SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
608137fb511SHong Zhang   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
609137fb511SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
610137fb511SHong Zhang     sctx.shift_top = 0;
611137fb511SHong Zhang     for (i=0; i<n; i++) {
612137fb511SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
613b051339dSHong Zhang       d  = (a->a)[diag[i]];
614137fb511SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
615b051339dSHong Zhang       v  = a->a+ai[i];
616b051339dSHong Zhang       nz = ai[i+1] - ai[i];
617137fb511SHong Zhang       for (j=0; j<nz; j++)
618137fb511SHong Zhang 	rs += PetscAbsScalar(v[j]);
619137fb511SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
620137fb511SHong Zhang     }
621137fb511SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
622137fb511SHong Zhang     sctx.shift_top    *= 1.1;
623137fb511SHong Zhang     sctx.nshift_max   = 5;
624137fb511SHong Zhang     sctx.shift_lo     = 0.;
625137fb511SHong Zhang     sctx.shift_hi     = 1.;
626137fb511SHong Zhang   }
627137fb511SHong Zhang 
628137fb511SHong Zhang   sctx.shift_amount = 0;
629137fb511SHong Zhang   sctx.nshift       = 0;
630137fb511SHong Zhang   do {
631137fb511SHong Zhang     sctx.lushift = PETSC_FALSE;
632137fb511SHong Zhang     for (i=0; i<n; i++){
633b051339dSHong Zhang       /* load in initial unfactored row */
634b051339dSHong Zhang       nz    = ai[r[i]+1] - ai[r[i]];
635b051339dSHong Zhang       ajtmp = aj + ai[r[i]];
636b051339dSHong Zhang       v     = a->a + ai[r[i]];
637137fb511SHong Zhang       /* sort permuted ajtmp and values v accordingly */
638b051339dSHong Zhang       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
639137fb511SHong Zhang       ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr);
640137fb511SHong Zhang 
641b051339dSHong Zhang       diag[r[i]] = ai[r[i]];
642137fb511SHong Zhang       for (j=0; j<nz; j++) {
643137fb511SHong Zhang         rtmp[ajtmp[j]] = v[j];
644b051339dSHong Zhang         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
645137fb511SHong Zhang       }
646137fb511SHong Zhang       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */
647137fb511SHong Zhang 
648b051339dSHong Zhang       row = *ajtmp++;
649137fb511SHong Zhang       while  (row < i) {
650137fb511SHong Zhang         pc = rtmp + row;
651137fb511SHong Zhang         if (*pc != 0.0) {
652b051339dSHong Zhang           pv         = a->a + diag[r[row]];
653b051339dSHong Zhang           pj         = aj + diag[r[row]] + 1;
654137fb511SHong Zhang 
655137fb511SHong Zhang           multiplier = *pc / *pv++;
656137fb511SHong Zhang           *pc        = multiplier;
657b051339dSHong Zhang           nz         = ai[r[row]+1] - diag[r[row]] - 1;
658b051339dSHong Zhang           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
659137fb511SHong Zhang           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
660137fb511SHong Zhang         }
661b051339dSHong Zhang         row = *ajtmp++;
662137fb511SHong Zhang       }
663b051339dSHong Zhang       /* finished row so overwrite it onto a->a */
664b051339dSHong Zhang       pv   = a->a + ai[r[i]] ;
665b051339dSHong Zhang       pj   = aj + ai[r[i]] ;
666b051339dSHong Zhang       nz   = ai[r[i]+1] - ai[r[i]];
667b051339dSHong Zhang       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */
668137fb511SHong Zhang 
669137fb511SHong Zhang       rs   = 0.0;
670137fb511SHong Zhang       for (j=0; j<nz; j++) {
671b051339dSHong Zhang         pv[j] = rtmp[pj[j]];
672b051339dSHong Zhang         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
673137fb511SHong Zhang       }
674137fb511SHong Zhang 
675137fb511SHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
676137fb511SHong Zhang       sctx.rs  = rs;
677b051339dSHong Zhang       sctx.pv  = pv[nbdiag];
678c6b1f410SHong Zhang       ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr);
679137fb511SHong Zhang       if (newshift == 1) break;
680137fb511SHong Zhang     }
681137fb511SHong Zhang 
682137fb511SHong Zhang     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
683137fb511SHong Zhang       /*
684137fb511SHong Zhang        * if no shift in this attempt & shifting & started shifting & can refine,
685137fb511SHong Zhang        * then try lower shift
686137fb511SHong Zhang        */
687137fb511SHong Zhang       sctx.shift_hi        = info->shift_fraction;
688137fb511SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
689137fb511SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
690137fb511SHong Zhang       sctx.lushift         = PETSC_TRUE;
691137fb511SHong Zhang       sctx.nshift++;
692137fb511SHong Zhang     }
693137fb511SHong Zhang   } while (sctx.lushift);
694137fb511SHong Zhang 
695137fb511SHong Zhang   /* invert diagonal entries for simplier triangular solves */
696137fb511SHong Zhang   for (i=0; i<n; i++) {
697b051339dSHong Zhang     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
698137fb511SHong Zhang   }
699137fb511SHong Zhang 
700137fb511SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
701137fb511SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
702137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
703b051339dSHong Zhang   A->factor     = FACTOR_LU;
704b051339dSHong Zhang   A->ops->solve = MatSolve_SeqAIJ_InplaceWithPerm;
705b051339dSHong Zhang   A->assembled = PETSC_TRUE;
706b051339dSHong Zhang   ierr = PetscLogFlops(A->cmap.n);CHKERRQ(ierr);
707137fb511SHong Zhang   if (sctx.nshift){
708137fb511SHong Zhang     if (info->shiftnz) {
7091e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
710137fb511SHong Zhang     } else if (info->shiftpd) {
7111e2582c4SBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
712137fb511SHong Zhang     }
713137fb511SHong Zhang   }
714137fb511SHong Zhang   PetscFunctionReturn(0);
715137fb511SHong Zhang }
716137fb511SHong Zhang 
7170889b5dcSvictorle #undef __FUNCT__
7180889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
719dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
7200889b5dcSvictorle {
7210889b5dcSvictorle   PetscFunctionBegin;
7220889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
7230889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
7240889b5dcSvictorle   PetscFunctionReturn(0);
7250889b5dcSvictorle }
7260889b5dcSvictorle 
7270889b5dcSvictorle 
7280a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
731dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
732da3a660dSBarry Smith {
733dfbe8321SBarry Smith   PetscErrorCode ierr;
734416022c9SBarry Smith   Mat            C;
735416022c9SBarry Smith 
7363a40ed3dSBarry Smith   PetscFunctionBegin;
7379e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
738af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
739273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
74052e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
7413a40ed3dSBarry Smith   PetscFunctionReturn(0);
742da3a660dSBarry Smith }
7430a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
7444a2ae208SSatish Balay #undef __FUNCT__
7454a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
746dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
7478c37ef55SBarry Smith {
748416022c9SBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
749416022c9SBarry Smith   IS                iscol = a->col,isrow = a->row;
7506849ba73SBarry Smith   PetscErrorCode    ierr;
751899cda47SBarry Smith   PetscInt          *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
75238baddfdSBarry Smith   PetscInt          nz,*rout,*cout;
753dd6ea824SBarry Smith   PetscScalar       *x,*tmp,*tmps,sum;
754d9fead3dSBarry Smith   const PetscScalar *b;
755dd6ea824SBarry Smith   const MatScalar   *aa = a->a,*v;
7568c37ef55SBarry Smith 
7573a40ed3dSBarry Smith   PetscFunctionBegin;
7583a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
75991bf9adeSBarry Smith 
760d9fead3dSBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
7611ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
762416022c9SBarry Smith   tmp  = a->solve_work;
7638c37ef55SBarry Smith 
7643b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7653b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
7668c37ef55SBarry Smith 
7678c37ef55SBarry Smith   /* forward solve the lower triangular */
7688c37ef55SBarry Smith   tmp[0] = b[*r++];
769010a20caSHong Zhang   tmps   = tmp;
7708c37ef55SBarry Smith   for (i=1; i<n; i++) {
771010a20caSHong Zhang     v   = aa + ai[i] ;
772010a20caSHong Zhang     vi  = aj + ai[i] ;
77353e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
77453e7425aSBarry Smith     sum = b[*r++];
775ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
7768c37ef55SBarry Smith     tmp[i] = sum;
7778c37ef55SBarry Smith   }
7788c37ef55SBarry Smith 
7798c37ef55SBarry Smith   /* backward solve the upper triangular */
7808c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
781010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
782010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
783416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
7848c37ef55SBarry Smith     sum = tmp[i];
785ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
786010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
7878c37ef55SBarry Smith   }
7888c37ef55SBarry Smith 
7893b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7903b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
791d9fead3dSBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
7921ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
793899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
7943a40ed3dSBarry Smith   PetscFunctionReturn(0);
7958c37ef55SBarry Smith }
796026e39d0SSatish Balay 
7972bebee5dSHong Zhang #undef __FUNCT__
7982bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
7992bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
8002bebee5dSHong Zhang {
8012bebee5dSHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
8022bebee5dSHong Zhang   IS              iscol = a->col,isrow = a->row;
8032bebee5dSHong Zhang   PetscErrorCode  ierr;
804899cda47SBarry Smith   PetscInt        *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
8052bebee5dSHong Zhang   PetscInt        nz,*rout,*cout,neq;
806dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
807dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
8082bebee5dSHong Zhang 
8092bebee5dSHong Zhang   PetscFunctionBegin;
8102bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
8112bebee5dSHong Zhang 
8122bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
8132bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
8142bebee5dSHong Zhang 
8152bebee5dSHong Zhang   tmp  = a->solve_work;
8162bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8172bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8182bebee5dSHong Zhang 
8192bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
8202bebee5dSHong Zhang     /* forward solve the lower triangular */
8212bebee5dSHong Zhang     tmp[0] = b[r[0]];
8222bebee5dSHong Zhang     tmps   = tmp;
8232bebee5dSHong Zhang     for (i=1; i<n; i++) {
8242bebee5dSHong Zhang       v   = aa + ai[i] ;
8252bebee5dSHong Zhang       vi  = aj + ai[i] ;
8262bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
8272bebee5dSHong Zhang       sum = b[r[i]];
8282bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8292bebee5dSHong Zhang       tmp[i] = sum;
8302bebee5dSHong Zhang     }
8312bebee5dSHong Zhang     /* backward solve the upper triangular */
8322bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
8332bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
8342bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
8352bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
8362bebee5dSHong Zhang       sum = tmp[i];
8372bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
8382bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
8392bebee5dSHong Zhang     }
8402bebee5dSHong Zhang 
8412bebee5dSHong Zhang     b += n;
8422bebee5dSHong Zhang     x += n;
8432bebee5dSHong Zhang   }
8442bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8452bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8462bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
8472bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
8482bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
8492bebee5dSHong Zhang   PetscFunctionReturn(0);
8502bebee5dSHong Zhang }
8512bebee5dSHong Zhang 
852137fb511SHong Zhang #undef __FUNCT__
853137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm"
854137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
855137fb511SHong Zhang {
856137fb511SHong Zhang   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
857137fb511SHong Zhang   IS              iscol = a->col,isrow = a->row;
858137fb511SHong Zhang   PetscErrorCode  ierr;
859137fb511SHong Zhang   PetscInt        *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
860137fb511SHong Zhang   PetscInt        nz,*rout,*cout,row;
861dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,*tmps,sum;
862dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
863137fb511SHong Zhang 
864137fb511SHong Zhang   PetscFunctionBegin;
865137fb511SHong Zhang   if (!n) PetscFunctionReturn(0);
866137fb511SHong Zhang 
867137fb511SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
868137fb511SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
869137fb511SHong Zhang   tmp  = a->solve_work;
870137fb511SHong Zhang 
871137fb511SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
872137fb511SHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
873137fb511SHong Zhang 
874137fb511SHong Zhang   /* forward solve the lower triangular */
875137fb511SHong Zhang   tmp[0] = b[*r++];
876137fb511SHong Zhang   tmps   = tmp;
877137fb511SHong Zhang   for (row=1; row<n; row++) {
878137fb511SHong Zhang     i   = rout[row]; /* permuted row */
879137fb511SHong Zhang     v   = aa + ai[i] ;
880137fb511SHong Zhang     vi  = aj + ai[i] ;
881137fb511SHong Zhang     nz  = a->diag[i] - ai[i];
882137fb511SHong Zhang     sum = b[*r++];
883137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
884137fb511SHong Zhang     tmp[row] = sum;
885137fb511SHong Zhang   }
886137fb511SHong Zhang 
887137fb511SHong Zhang   /* backward solve the upper triangular */
888137fb511SHong Zhang   for (row=n-1; row>=0; row--){
889137fb511SHong Zhang     i   = rout[row]; /* permuted row */
890137fb511SHong Zhang     v   = aa + a->diag[i] + 1;
891137fb511SHong Zhang     vi  = aj + a->diag[i] + 1;
892137fb511SHong Zhang     nz  = ai[i+1] - a->diag[i] - 1;
893137fb511SHong Zhang     sum = tmp[row];
894137fb511SHong Zhang     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
895137fb511SHong Zhang     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
896137fb511SHong Zhang   }
897137fb511SHong Zhang 
898137fb511SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
899137fb511SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
900137fb511SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
901137fb511SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
902137fb511SHong Zhang   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
903137fb511SHong Zhang   PetscFunctionReturn(0);
904137fb511SHong Zhang }
905137fb511SHong Zhang 
906930ae53cSBarry Smith /* ----------------------------------------------------------- */
9074a2ae208SSatish Balay #undef __FUNCT__
9084a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
909dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
910930ae53cSBarry Smith {
911930ae53cSBarry Smith   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
9126849ba73SBarry Smith   PetscErrorCode    ierr;
913*356650c2SBarry Smith   PetscInt          n = A->rmap.n;
914*356650c2SBarry Smith   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
915*356650c2SBarry Smith   PetscScalar       *x;
916*356650c2SBarry Smith   const PetscScalar *b;
917dd6ea824SBarry Smith   const MatScalar   *aa = a->a;
918aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
919*356650c2SBarry Smith   PetscInt          adiag_i,i,nz,ai_i;
920dd6ea824SBarry Smith   const MatScalar   *v;
921dd6ea824SBarry Smith   PetscScalar       sum;
922d85afc42SSatish Balay #endif
923930ae53cSBarry Smith 
9243a40ed3dSBarry Smith   PetscFunctionBegin;
9253a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
926930ae53cSBarry Smith 
927*356650c2SBarry Smith   ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
929930ae53cSBarry Smith 
930aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
9311c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
9321c4feecaSSatish Balay #else
933930ae53cSBarry Smith   /* forward solve the lower triangular */
934930ae53cSBarry Smith   x[0] = b[0];
935930ae53cSBarry Smith   for (i=1; i<n; i++) {
936930ae53cSBarry Smith     ai_i = ai[i];
937930ae53cSBarry Smith     v    = aa + ai_i;
938930ae53cSBarry Smith     vi   = aj + ai_i;
939930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
940930ae53cSBarry Smith     sum  = b[i];
941930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
942930ae53cSBarry Smith     x[i] = sum;
943930ae53cSBarry Smith   }
944930ae53cSBarry Smith 
945930ae53cSBarry Smith   /* backward solve the upper triangular */
946930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
947930ae53cSBarry Smith     adiag_i = adiag[i];
948d708749eSBarry Smith     v       = aa + adiag_i + 1;
949d708749eSBarry Smith     vi      = aj + adiag_i + 1;
950930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
951930ae53cSBarry Smith     sum     = x[i];
952930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
953930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
954930ae53cSBarry Smith   }
9551c4feecaSSatish Balay #endif
956899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
957*356650c2SBarry Smith   ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
9581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9593a40ed3dSBarry Smith   PetscFunctionReturn(0);
960930ae53cSBarry Smith }
961930ae53cSBarry Smith 
9624a2ae208SSatish Balay #undef __FUNCT__
9634a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
964dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
965da3a660dSBarry Smith {
966416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
967416022c9SBarry Smith   IS              iscol = a->col,isrow = a->row;
9686849ba73SBarry Smith   PetscErrorCode  ierr;
969899cda47SBarry Smith   PetscInt        *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
97038baddfdSBarry Smith   PetscInt        nz,*rout,*cout;
971dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,sum;
972dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
973da3a660dSBarry Smith 
9743a40ed3dSBarry Smith   PetscFunctionBegin;
97578b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
976da3a660dSBarry Smith 
9771ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
9781ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
979416022c9SBarry Smith   tmp  = a->solve_work;
980da3a660dSBarry Smith 
9813b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
9823b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
983da3a660dSBarry Smith 
984da3a660dSBarry Smith   /* forward solve the lower triangular */
985da3a660dSBarry Smith   tmp[0] = b[*r++];
986da3a660dSBarry Smith   for (i=1; i<n; i++) {
987010a20caSHong Zhang     v   = aa + ai[i] ;
988010a20caSHong Zhang     vi  = aj + ai[i] ;
989416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
990da3a660dSBarry Smith     sum = b[*r++];
991010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
992da3a660dSBarry Smith     tmp[i] = sum;
993da3a660dSBarry Smith   }
994da3a660dSBarry Smith 
995da3a660dSBarry Smith   /* backward solve the upper triangular */
996da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
997010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
998010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
999416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
1000da3a660dSBarry Smith     sum = tmp[i];
1001010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
1002010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
1003da3a660dSBarry Smith     x[*c--] += tmp[i];
1004da3a660dSBarry Smith   }
1005da3a660dSBarry Smith 
10063b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10073b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10081ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10091ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1010efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
1011898183ecSLois Curfman McInnes 
10123a40ed3dSBarry Smith   PetscFunctionReturn(0);
1013da3a660dSBarry Smith }
1014da3a660dSBarry Smith /* -------------------------------------------------------------------*/
10154a2ae208SSatish Balay #undef __FUNCT__
10164a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
1017dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1018da3a660dSBarry Smith {
1019416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10202235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10216849ba73SBarry Smith   PetscErrorCode  ierr;
1022899cda47SBarry Smith   PetscInt        *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
102338baddfdSBarry Smith   PetscInt        nz,*rout,*cout,*diag = a->diag;
1024dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp,s1;
1025dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
1026da3a660dSBarry Smith 
10273a40ed3dSBarry Smith   PetscFunctionBegin;
10281ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10291ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1030416022c9SBarry Smith   tmp  = a->solve_work;
1031da3a660dSBarry Smith 
10322235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10332235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1034da3a660dSBarry Smith 
1035da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
10362235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
1037da3a660dSBarry Smith 
1038da3a660dSBarry Smith   /* forward solve the U^T */
1039da3a660dSBarry Smith   for (i=0; i<n; i++) {
1040010a20caSHong Zhang     v   = aa + diag[i] ;
1041010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1042f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
1043c38d4ed2SBarry Smith     s1  = tmp[i];
1044ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1045da3a660dSBarry Smith     while (nz--) {
1046010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
1047da3a660dSBarry Smith     }
1048c38d4ed2SBarry Smith     tmp[i] = s1;
1049da3a660dSBarry Smith   }
1050da3a660dSBarry Smith 
1051da3a660dSBarry Smith   /* backward solve the L^T */
1052da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
1053010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1054010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1055f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
1056f1af5d2fSBarry Smith     s1  = tmp[i];
1057da3a660dSBarry Smith     while (nz--) {
1058010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
1059da3a660dSBarry Smith     }
1060da3a660dSBarry Smith   }
1061da3a660dSBarry Smith 
1062da3a660dSBarry Smith   /* copy tmp into x according to permutation */
1063da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
1064da3a660dSBarry Smith 
10652235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
10662235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
10671ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
10681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1069da3a660dSBarry Smith 
1070899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap.n);CHKERRQ(ierr);
10713a40ed3dSBarry Smith   PetscFunctionReturn(0);
1072da3a660dSBarry Smith }
1073da3a660dSBarry Smith 
10744a2ae208SSatish Balay #undef __FUNCT__
10754a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
1076dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1077da3a660dSBarry Smith {
1078416022c9SBarry Smith   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
10792235e667SBarry Smith   IS              iscol = a->col,isrow = a->row;
10806849ba73SBarry Smith   PetscErrorCode  ierr;
1081899cda47SBarry Smith   PetscInt        *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
108238baddfdSBarry Smith   PetscInt        nz,*rout,*cout,*diag = a->diag;
1083dd6ea824SBarry Smith   PetscScalar     *x,*b,*tmp;
1084dd6ea824SBarry Smith   const MatScalar *aa = a->a,*v;
10856abc6512SBarry Smith 
10863a40ed3dSBarry Smith   PetscFunctionBegin;
108723bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
10886abc6512SBarry Smith 
10891ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
10901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1091416022c9SBarry Smith   tmp = a->solve_work;
10926abc6512SBarry Smith 
10932235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
10942235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
10956abc6512SBarry Smith 
10966abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
10972235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
10986abc6512SBarry Smith 
10996abc6512SBarry Smith   /* forward solve the U^T */
11006abc6512SBarry Smith   for (i=0; i<n; i++) {
1101010a20caSHong Zhang     v   = aa + diag[i] ;
1102010a20caSHong Zhang     vi  = aj + diag[i] + 1;
1103f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
11046abc6512SBarry Smith     tmp[i] *= *v++;
11056abc6512SBarry Smith     while (nz--) {
1106010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
11076abc6512SBarry Smith     }
11086abc6512SBarry Smith   }
11096abc6512SBarry Smith 
11106abc6512SBarry Smith   /* backward solve the L^T */
11116abc6512SBarry Smith   for (i=n-1; i>=0; i--){
1112010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
1113010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
1114f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
11156abc6512SBarry Smith     while (nz--) {
1116010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
11176abc6512SBarry Smith     }
11186abc6512SBarry Smith   }
11196abc6512SBarry Smith 
11206abc6512SBarry Smith   /* copy tmp into x according to permutation */
11216abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
11226abc6512SBarry Smith 
11232235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
11242235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
11251ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11261ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11276abc6512SBarry Smith 
1128efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
11293a40ed3dSBarry Smith   PetscFunctionReturn(0);
1130da3a660dSBarry Smith }
11319e25ed09SBarry Smith /* ----------------------------------------------------------------*/
11323c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
113308480c60SBarry Smith 
11344a2ae208SSatish Balay #undef __FUNCT__
11354a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
1136dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
11379e25ed09SBarry Smith {
1138416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
11399e25ed09SBarry Smith   IS                 isicol;
11406849ba73SBarry Smith   PetscErrorCode     ierr;
114109f38230SBarry Smith   PetscInt           *r,*ic,n=A->rmap.n,*ai=a->i,*aj=a->j,d;
11425a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
11435a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
11445a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
114577c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
1146329f5518SBarry Smith   PetscReal          f;
11475a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
11485a8e39fbSHong Zhang   PetscBT            lnkbt;
11495a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1150a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1151a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
115209f38230SBarry Smith   PetscTruth         missing;
11539e25ed09SBarry Smith 
11543a40ed3dSBarry Smith   PetscFunctionBegin;
115509192fe3SBarry Smith   if (A->rmap.n != A->cmap.n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap.n,A->cmap.n);
11566cf997b0SBarry Smith   f             = info->fill;
115738baddfdSBarry Smith   levels        = (PetscInt)info->levels;
115838baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
11594c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
116082bf6240SBarry Smith 
116108480c60SBarry Smith   /* special case that simply copies fill pattern */
1162be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1163be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
116486bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
11652e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
116608480c60SBarry Smith     (*fact)->factor                 = FACTOR_LU;
1167f3a39becSBarry Smith     (*fact)->info.factor_mallocs    = 0;
1168f3a39becSBarry Smith     (*fact)->info.fill_ratio_given  = info->fill;
1169f3a39becSBarry Smith     (*fact)->info.fill_ratio_needed = 1.0;
117008480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
11718ef3462fSBarry Smith     ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
117209f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
117308480c60SBarry Smith     b->row              = isrow;
117408480c60SBarry Smith     b->col              = iscol;
117582bf6240SBarry Smith     b->icol             = isicol;
1176899cda47SBarry Smith     ierr                = PetscMalloc(((*fact)->rmap.n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1177f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
1178c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1179c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
11803c5fc038SBarry Smith     ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
11813a40ed3dSBarry Smith     PetscFunctionReturn(0);
118208480c60SBarry Smith   }
118308480c60SBarry Smith 
1184898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1185898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
11869e25ed09SBarry Smith 
11879e25ed09SBarry Smith   /* get new row pointers */
11885a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
11895a8e39fbSHong Zhang   bi[0] = 0;
11905a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
11915a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
11925a8e39fbSHong Zhang   bdiag[0]  = 0;
11936cf997b0SBarry Smith 
11945a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
11955a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
11965a8e39fbSHong Zhang 
11975a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
11985a8e39fbSHong Zhang   nlnk = n + 1;
11995a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12005a8e39fbSHong Zhang 
12015a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
1202a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
12035a8e39fbSHong Zhang   current_space = free_space;
1204a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
12055a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
12065a8e39fbSHong Zhang 
12075a8e39fbSHong Zhang   for (i=0; i<n; i++) {
12085a8e39fbSHong Zhang     nzi = 0;
12096cf997b0SBarry Smith     /* copy current row into linked list */
12105a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
12115a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
12125a8e39fbSHong Zhang     cols = aj + ai[r[i]];
12135a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
12145a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12155a8e39fbSHong Zhang     nzi += nlnk;
12166cf997b0SBarry Smith 
12176cf997b0SBarry Smith     /* make sure diagonal entry is included */
12185a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
12196cf997b0SBarry Smith       fm = n;
12205a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
12215a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
12225a8e39fbSHong Zhang       lnk[fm]    = i;
12235a8e39fbSHong Zhang       lnk_lvl[i] = 0;
12245a8e39fbSHong Zhang       nzi++; dcount++;
12256cf997b0SBarry Smith     }
12266cf997b0SBarry Smith 
12275a8e39fbSHong Zhang     /* add pivot rows into the active row */
12285a8e39fbSHong Zhang     nzbd = 0;
12295a8e39fbSHong Zhang     prow = lnk[n];
12305a8e39fbSHong Zhang     while (prow < i) {
12315a8e39fbSHong Zhang       nnz      = bdiag[prow];
12325a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
12335a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
12345a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
12355a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
12365a8e39fbSHong Zhang       nzi += nlnk;
12375a8e39fbSHong Zhang       prow = lnk[prow];
12385a8e39fbSHong Zhang       nzbd++;
123956beaf04SBarry Smith     }
12405a8e39fbSHong Zhang     bdiag[i] = nzbd;
12415a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1242ecf371e4SBarry Smith 
12435a8e39fbSHong Zhang     /* if free space is not available, make more free space */
12445a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
12455a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1246a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1247a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
12485a8e39fbSHong Zhang       reallocs++;
12495a8e39fbSHong Zhang     }
1250ecf371e4SBarry Smith 
12515a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
12525a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
12535a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
12545a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
12555a8e39fbSHong Zhang 
12565a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
12575a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
125877431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1259d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
12606cf997b0SBarry Smith     }
12615a8e39fbSHong Zhang 
12625a8e39fbSHong Zhang     current_space->array           += nzi;
12635a8e39fbSHong Zhang     current_space->local_used      += nzi;
12645a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
12655a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
12665a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
12675a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
12689e25ed09SBarry Smith   }
12695a8e39fbSHong Zhang 
1270898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1271898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
12725a8e39fbSHong Zhang 
12735a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
12745a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1275a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
12765a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1277a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
12785a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
12799e25ed09SBarry Smith 
12806cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1281f64065bbSBarry Smith   {
12825a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1283ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1284ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1285ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1286ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1287335d9088SBarry Smith     if (diagonal_fill) {
1288ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1289335d9088SBarry Smith     }
1290f64065bbSBarry Smith   }
129163ba0a88SBarry Smith #endif
1292bbb6d6a8SBarry Smith 
12939e25ed09SBarry Smith   /* put together the new matrix */
12947adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr);
1295f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
12967adad957SLisandro Dalcin   ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr);
1297ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
129852e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1299416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1300e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1301e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
13027c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
13035a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1304b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
13055a8e39fbSHong Zhang   b->j          = bj;
13065a8e39fbSHong Zhang   b->i          = bi;
13075a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
13085a8e39fbSHong Zhang   b->diag       = bdiag;
1309416022c9SBarry Smith   b->ilen       = 0;
1310416022c9SBarry Smith   b->imax       = 0;
1311416022c9SBarry Smith   b->row        = isrow;
1312416022c9SBarry Smith   b->col        = iscol;
1313c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1314c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
131582bf6240SBarry Smith   b->icol       = isicol;
131687828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1317416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
13185a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
13195a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
13205a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
13219e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1322418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1323ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
13245a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
132571c2f376SKris Buschelman 
132654e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
132771c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
13284846f1f5SKris Buschelman 
13293a40ed3dSBarry Smith   PetscFunctionReturn(0);
13309e25ed09SBarry Smith }
1331d5d45c9bSBarry Smith 
13323c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1333a6175056SHong Zhang #undef __FUNCT__
1334a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1335af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1336a6175056SHong Zhang {
13372ed133dbSHong Zhang   Mat            C = *B;
13380968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
13392ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
13409bfd6278SHong Zhang   IS             ip=b->row,iip = b->icol;
13412ed133dbSHong Zhang   PetscErrorCode ierr;
13429bfd6278SHong Zhang   PetscInt       *rip,*riip,i,j,mbs=A->rmap.n,*bi=b->i,*bj=b->j,*bcol;
13432ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
13442ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1345622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1346540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1347fbf22428SSatish Balay   PetscReal      shiftpd;
1348540703acSHong Zhang   ChShift_Ctx    sctx;
1349540703acSHong Zhang   PetscInt       newshift;
1350d5d45c9bSBarry Smith 
1351a6175056SHong Zhang   PetscFunctionBegin;
1352117f1891Stmunson 
1353540703acSHong Zhang   shiftnz   = info->shiftnz;
1354540703acSHong Zhang   shiftpd   = info->shiftpd;
1355ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1356ee45ca4aSHong Zhang 
13572ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
13589bfd6278SHong Zhang   ierr  = ISGetIndices(iip,&riip);CHKERRQ(ierr);
13592ed133dbSHong Zhang 
13602ed133dbSHong Zhang   /* initialization */
13612ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
13622ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
13632ed133dbSHong Zhang   jl   = il + mbs;
13642ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
13652ed133dbSHong Zhang 
1366540703acSHong Zhang   sctx.shift_amount = 0;
1367540703acSHong Zhang   sctx.nshift       = 0;
13682ed133dbSHong Zhang   do {
1369540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
13702ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
13712ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1372a6175056SHong Zhang     }
13732ed133dbSHong Zhang 
13742ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1375c05c3958SHong Zhang       bval = ba + bi[k];
13762ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
13772ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
13782ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
13799bfd6278SHong Zhang         col = riip[aj[j]];
13802ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
13812ed133dbSHong Zhang           rtmp[col] = aa[j];
1382c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
13832ed133dbSHong Zhang         }
13842ed133dbSHong Zhang       }
138539e3d630SHong Zhang       /* shift the diagonal of the matrix */
1386540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
13872ed133dbSHong Zhang 
13882ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
13892ed133dbSHong Zhang       dk = rtmp[k];
13902ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
13912ed133dbSHong Zhang 
13922ed133dbSHong Zhang       while (i < k){
13932ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
13942ed133dbSHong Zhang 
13952ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
13962ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
13972ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
13982ed133dbSHong Zhang         dk += uikdi*ba[ili];
13992ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
14002ed133dbSHong Zhang 
14012ed133dbSHong Zhang         /* add multiple of row i to k-th row */
14022ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
14032ed133dbSHong Zhang         if (jmin < jmax){
14042ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
14052ed133dbSHong Zhang           /* update il and jl for row i */
14062ed133dbSHong Zhang           il[i] = jmin;
14072ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
14082ed133dbSHong Zhang         }
14092ed133dbSHong Zhang         i = nexti;
14102ed133dbSHong Zhang       }
14112ed133dbSHong Zhang 
14120697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
14130697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
14140697b6caSHong Zhang       rs   = 0.0;
14153c9e8598SHong Zhang       jmin = bi[k]+1;
14163c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
14173c9e8598SHong Zhang       bcol = bj + jmin;
14183c9e8598SHong Zhang       while (nz--){
141989f845c8SHong Zhang         rs += PetscAbsScalar(rtmp[*bcol]);
142089f845c8SHong Zhang         bcol++;
14213c9e8598SHong Zhang       }
1422540703acSHong Zhang 
1423540703acSHong Zhang       sctx.rs = rs;
1424540703acSHong Zhang       sctx.pv = dk;
14255b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
1426117f1891Stmunson 
1427117f1891Stmunson       if (newshift == 1) {
1428117f1891Stmunson         if (!sctx.shift_amount) {
1429117f1891Stmunson           sctx.shift_amount = 1e-5;
1430117f1891Stmunson         }
1431117f1891Stmunson         break;
1432117f1891Stmunson       }
14333c9e8598SHong Zhang 
14343c9e8598SHong Zhang       /* copy data into U(k,:) */
143539e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
143639e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
143739e3d630SHong Zhang       if (jmin < jmax) {
143839e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
143939e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
14403c9e8598SHong Zhang         }
144139e3d630SHong Zhang         /* add the k-th row into il and jl */
14423c9e8598SHong Zhang         il[k] = jmin;
14433c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
14443c9e8598SHong Zhang       }
14453c9e8598SHong Zhang     }
1446540703acSHong Zhang   } while (sctx.chshift);
14473c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
14483c9e8598SHong Zhang 
144939e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
14509bfd6278SHong Zhang   ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr);
14513c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
14523c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
14533c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1454899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.n);CHKERRQ(ierr);
1455540703acSHong Zhang   if (sctx.nshift){
1456540703acSHong Zhang     if (shiftnz) {
14571e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1458540703acSHong Zhang     } else if (shiftpd) {
14591e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
14600697b6caSHong Zhang     }
14613c9e8598SHong Zhang   }
14623c9e8598SHong Zhang   PetscFunctionReturn(0);
14633c9e8598SHong Zhang }
1464a6175056SHong Zhang 
1465a6175056SHong Zhang #undef __FUNCT__
1466a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1467dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1468a6175056SHong Zhang {
14690968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1470ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1471ed59401aSHong Zhang   Mat                B;
1472dfbe8321SBarry Smith   PetscErrorCode     ierr;
147358ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
1474b635d36bSHong Zhang   PetscInt           reallocs=0,*rip,*riip,i,*ai=a->i,*aj=a->j,am=A->rmap.n,*ui;
1475ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
147658ebbce7SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,d;
14775a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1478ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1479a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1480a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
14817a48dd6fSHong Zhang   PetscBT            lnkbt;
1482b635d36bSHong Zhang   IS                 iperm;
1483a6175056SHong Zhang 
1484a6175056SHong Zhang   PetscFunctionBegin;
148509192fe3SBarry Smith   if (A->rmap.n != A->cmap.n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap.n,A->cmap.n);
148658ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
148758ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
1488653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1489b635d36bSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14907a48dd6fSHong Zhang 
149139e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
149239e3d630SHong Zhang   ui[0] = 0;
149339e3d630SHong Zhang 
1494b635d36bSHong Zhang   /* ICC(0) without matrix ordering: simply copies fill pattern */
1495622e2028SHong Zhang   if (!levels && perm_identity) {
149658ebbce7SBarry Smith 
1497ed59401aSHong Zhang     for (i=0; i<am; i++) {
149839e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1499ed59401aSHong Zhang     }
150039e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
150139e3d630SHong Zhang     cols = uj;
1502ed59401aSHong Zhang     for (i=0; i<am; i++) {
1503ed59401aSHong Zhang       aj    = a->j + a->diag[i];
150439e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
150539e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1506ed59401aSHong Zhang     }
150739e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
1508b635d36bSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1509b635d36bSHong Zhang     ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1510b635d36bSHong Zhang 
15117a48dd6fSHong Zhang     /* initialization */
15125a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
15137a48dd6fSHong Zhang 
15147a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
15157a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
15161393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
15177a48dd6fSHong Zhang     il         = jl + am;
15187a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
15197a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
15207a48dd6fSHong Zhang     for (i=0; i<am; i++){
15217a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
15227a48dd6fSHong Zhang     }
15237a48dd6fSHong Zhang 
15247a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
15257a48dd6fSHong Zhang     nlnk = am + 1;
15262ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15277a48dd6fSHong Zhang 
15287a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1529a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
15307a48dd6fSHong Zhang     current_space = free_space;
1531a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
15327a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
15337a48dd6fSHong Zhang 
15347a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
15357a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
15367a48dd6fSHong Zhang       nzk   = 0;
1537622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1538b635d36bSHong Zhang       if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
1539622e2028SHong Zhang       ncols_upper = 0;
1540622e2028SHong Zhang       for (j=0; j<ncols; j++){
1541b635d36bSHong Zhang         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
1542b635d36bSHong Zhang         if (riip[i] >= k){ /* only take upper triangular entry */
15435a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1544622e2028SHong Zhang           ncols_upper++;
1545622e2028SHong Zhang         }
1546622e2028SHong Zhang       }
1547b635d36bSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15487a48dd6fSHong Zhang       nzk += nlnk;
15497a48dd6fSHong Zhang 
15507a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
15517a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
15527a48dd6fSHong Zhang 
15537a48dd6fSHong Zhang       while (prow < k){
15547a48dd6fSHong Zhang         nextprow = jl[prow];
15557a48dd6fSHong Zhang 
15567a48dd6fSHong Zhang         /* merge prow into k-th row */
15577a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
15587a48dd6fSHong Zhang         jmax = ui[prow+1];
15597a48dd6fSHong Zhang         ncols = jmax-jmin;
1560ed59401aSHong Zhang         i     = jmin - ui[prow];
1561ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15625a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
15635a8e39fbSHong Zhang         j     = *(uj - 1);
15645a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
15657a48dd6fSHong Zhang         nzk += nlnk;
15667a48dd6fSHong Zhang 
15677a48dd6fSHong Zhang         /* update il and jl for prow */
15687a48dd6fSHong Zhang         if (jmin < jmax){
15697a48dd6fSHong Zhang           il[prow] = jmin;
15707a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
15717a48dd6fSHong Zhang         }
15727a48dd6fSHong Zhang         prow = nextprow;
15737a48dd6fSHong Zhang       }
15747a48dd6fSHong Zhang 
15757a48dd6fSHong Zhang       /* if free space is not available, make more free space */
15767a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
15777a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
15787a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1579a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1580a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
15817a48dd6fSHong Zhang         reallocs++;
15827a48dd6fSHong Zhang       }
15837a48dd6fSHong Zhang 
15847a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
1585b635d36bSHong Zhang       if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
15862ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
15877a48dd6fSHong Zhang 
15887a48dd6fSHong Zhang       /* add the k-th row into il and jl */
158939e3d630SHong Zhang       if (nzk > 1){
15907a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
15917a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
15927a48dd6fSHong Zhang         il[k] = ui[k] + 1;
15937a48dd6fSHong Zhang       }
15947a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
15957a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
15967a48dd6fSHong Zhang 
15977a48dd6fSHong Zhang       current_space->array           += nzk;
15987a48dd6fSHong Zhang       current_space->local_used      += nzk;
15997a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
16007a48dd6fSHong Zhang 
16017a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
16027a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
16037a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
16047a48dd6fSHong Zhang 
16057a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
16067a48dd6fSHong Zhang     }
16077a48dd6fSHong Zhang 
16086cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16097a48dd6fSHong Zhang     if (ai[am] != 0) {
161039e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1611ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1612ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1613ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
16147a48dd6fSHong Zhang     } else {
1615ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16167a48dd6fSHong Zhang     }
161763ba0a88SBarry Smith #endif
16187a48dd6fSHong Zhang 
16197a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1620b635d36bSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16217a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
16225a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
16237a48dd6fSHong Zhang 
16247a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
16257a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1626a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
16272ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1628a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
16297a48dd6fSHong Zhang 
163039e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
163139e3d630SHong Zhang 
16327a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1633f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1634f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
1635ed59401aSHong Zhang   B = *fact;
1636ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1637ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
16387a48dd6fSHong Zhang 
1639ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
16407a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
16417a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
16427a48dd6fSHong Zhang   b->j    = uj;
16437a48dd6fSHong Zhang   b->i    = ui;
16447a48dd6fSHong Zhang   b->diag = 0;
16457a48dd6fSHong Zhang   b->ilen = 0;
16467a48dd6fSHong Zhang   b->imax = 0;
16477a48dd6fSHong Zhang   b->row  = perm;
1648b635d36bSHong Zhang   b->col  = perm;
1649b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1650b635d36bSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1651b635d36bSHong Zhang   b->icol = iperm;
16527a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
16537a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
165452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
16557a48dd6fSHong Zhang   b->maxnz   = b->nz = ui[am];
1656e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1657e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
16587a48dd6fSHong Zhang 
1659ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1660ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1661ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
16627a48dd6fSHong Zhang   if (ai[am] != 0) {
1663ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
16647a48dd6fSHong Zhang   } else {
1665ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
16667a48dd6fSHong Zhang   }
166739e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1668622e2028SHong Zhang   if (perm_identity){
1669ed59401aSHong Zhang     B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1670ed59401aSHong Zhang     B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1671622e2028SHong Zhang   }
1672a6175056SHong Zhang   PetscFunctionReturn(0);
1673a6175056SHong Zhang }
1674d5d45c9bSBarry Smith 
1675f76d2b81SHong Zhang #undef __FUNCT__
1676f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1677dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1678f76d2b81SHong Zhang {
1679f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
168010c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
16812ed133dbSHong Zhang   Mat                B;
1682dfbe8321SBarry Smith   PetscErrorCode     ierr;
1683f76d2b81SHong Zhang   PetscTruth         perm_identity;
168410c27e3dSHong Zhang   PetscReal          fill = info->fill;
1685899cda47SBarry Smith   PetscInt           *rip,*riip,i,am=A->rmap.n,*ai=a->i,*aj=a->j,reallocs=0,prow;
168610c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
16872ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1688a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
168910c27e3dSHong Zhang   PetscBT            lnkbt;
16902ed133dbSHong Zhang   IS                 iperm;
1691f76d2b81SHong Zhang 
1692f76d2b81SHong Zhang   PetscFunctionBegin;
169309192fe3SBarry Smith   if (A->rmap.n != A->cmap.n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap.n,A->cmap.n);
16942ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1695f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
16962ed133dbSHong Zhang   ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
16972ed133dbSHong Zhang   ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
16989bfd6278SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
169910c27e3dSHong Zhang 
170010c27e3dSHong Zhang   /* initialization */
17011393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
170210c27e3dSHong Zhang   ui[0] = 0;
170310c27e3dSHong Zhang 
170410c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
17051393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
17061393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
17071393bc97SHong Zhang   il     = jl + am;
17081393bc97SHong Zhang   cols   = il + am;
17091393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
17101393bc97SHong Zhang   for (i=0; i<am; i++){
17111393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1712f76d2b81SHong Zhang   }
171310c27e3dSHong Zhang 
171410c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
17151393bc97SHong Zhang   nlnk = am + 1;
17161393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
171710c27e3dSHong Zhang 
17181393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1719a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
172010c27e3dSHong Zhang   current_space = free_space;
172110c27e3dSHong Zhang 
17221393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
172310c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
172410c27e3dSHong Zhang     nzk   = 0;
17252ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
17269bfd6278SHong Zhang     if (!ncols) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix");
17272ed133dbSHong Zhang     ncols_upper = 0;
17282ed133dbSHong Zhang     for (j=0; j<ncols; j++){
17299bfd6278SHong Zhang       i = riip[*(aj + ai[rip[k]] + j)];
17302ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
17312ed133dbSHong Zhang         cols[ncols_upper] = i;
17322ed133dbSHong Zhang         ncols_upper++;
17332ed133dbSHong Zhang       }
17342ed133dbSHong Zhang     }
17351393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
173610c27e3dSHong Zhang     nzk += nlnk;
173710c27e3dSHong Zhang 
173810c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
173910c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
174010c27e3dSHong Zhang 
174110c27e3dSHong Zhang     while (prow < k){
174210c27e3dSHong Zhang       nextprow = jl[prow];
174310c27e3dSHong Zhang       /* merge prow into k-th row */
17441393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
174510c27e3dSHong Zhang       jmax = ui[prow+1];
174610c27e3dSHong Zhang       ncols = jmax-jmin;
17471393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
17485a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
174910c27e3dSHong Zhang       nzk += nlnk;
175010c27e3dSHong Zhang 
175110c27e3dSHong Zhang       /* update il and jl for prow */
175210c27e3dSHong Zhang       if (jmin < jmax){
175310c27e3dSHong Zhang         il[prow] = jmin;
17542ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
175510c27e3dSHong Zhang       }
175610c27e3dSHong Zhang       prow = nextprow;
175710c27e3dSHong Zhang     }
175810c27e3dSHong Zhang 
175910c27e3dSHong Zhang     /* if free space is not available, make more free space */
176010c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
17611393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
176210c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1763a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
176410c27e3dSHong Zhang       reallocs++;
176510c27e3dSHong Zhang     }
176610c27e3dSHong Zhang 
176710c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
17681393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
176910c27e3dSHong Zhang 
177010c27e3dSHong Zhang     /* add the k-th row into il and jl */
177110c27e3dSHong Zhang     if (nzk-1 > 0){
17721393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
177310c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
177410c27e3dSHong Zhang       il[k] = ui[k] + 1;
177510c27e3dSHong Zhang     }
17762ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
177710c27e3dSHong Zhang     current_space->array           += nzk;
177810c27e3dSHong Zhang     current_space->local_used      += nzk;
177910c27e3dSHong Zhang     current_space->local_remaining -= nzk;
178010c27e3dSHong Zhang 
178110c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
178210c27e3dSHong Zhang   }
178310c27e3dSHong Zhang 
17846cf91177SBarry Smith #if defined(PETSC_USE_INFO)
17851393bc97SHong Zhang   if (ai[am] != 0) {
17861393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1787ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1788ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1789ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
179010c27e3dSHong Zhang   } else {
1791ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
179210c27e3dSHong Zhang   }
179363ba0a88SBarry Smith #endif
179410c27e3dSHong Zhang 
179510c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
17969bfd6278SHong Zhang   ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
179710c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
179810c27e3dSHong Zhang 
179910c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
18001393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1801a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
180210c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
180310c27e3dSHong Zhang 
180410c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1805f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1806f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
18072ed133dbSHong Zhang   B    = *fact;
18082ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1809ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
181010c27e3dSHong Zhang 
18112ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
181210c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1813e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1814e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
18151393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
181610c27e3dSHong Zhang   b->j    = uj;
181710c27e3dSHong Zhang   b->i    = ui;
181810c27e3dSHong Zhang   b->diag = 0;
181910c27e3dSHong Zhang   b->ilen = 0;
182010c27e3dSHong Zhang   b->imax = 0;
182110c27e3dSHong Zhang   b->row  = perm;
18229bfd6278SHong Zhang   b->col  = perm;
18239bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18249bfd6278SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18259bfd6278SHong Zhang   b->icol = iperm;
182610c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18271393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
18281393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
18291393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
183010c27e3dSHong Zhang 
18312ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
18322ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
18332ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
18341393bc97SHong Zhang   if (ai[am] != 0) {
18351393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
183610c27e3dSHong Zhang   } else {
18372ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
183810c27e3dSHong Zhang   }
183939e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
18402ed133dbSHong Zhang   if (perm_identity){
184110c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
184210c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1843fd531716SHong Zhang     (*fact)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1844fd531716SHong Zhang     (*fact)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1845fd531716SHong Zhang   } else {
1846fd531716SHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1;
1847fd531716SHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
1848fd531716SHong Zhang     (*fact)->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1;
1849fd531716SHong Zhang     (*fact)->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1;
18502ed133dbSHong Zhang   }
1851f76d2b81SHong Zhang   PetscFunctionReturn(0);
1852f76d2b81SHong Zhang }
1853