xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 09f382308e660e6aeafeafa6047bab56e09352c9)
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)
2338baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
249cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2538baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,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;
6787828ca2SBarry Smith   PetscScalar    *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);
9887828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&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);
12287828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&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 
222f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
223f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
224f204ca49SKris Buschelman   ierr = MatSetType(*fact,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;
23286bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
23307d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
23486bacbd2SBarry Smith   b->a             = new_a;
23586bacbd2SBarry Smith   b->j             = new_j;
23686bacbd2SBarry Smith   b->i             = new_i;
23786bacbd2SBarry Smith   b->ilen          = 0;
23886bacbd2SBarry Smith   b->imax          = 0;
2391f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
240313ae024SBarry Smith   b->row           = isirow;
24186bacbd2SBarry Smith   b->col           = iscolf;
24287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24386bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24486bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24586bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24686bacbd2SBarry Smith 
247431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
248ae15b995SBarry Smith   ierr = PetscInfo2(A,"Fill ratio:given %G needed %G\n",info->fill,af);CHKERRQ(ierr);
249ae15b995SBarry Smith   ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
250ae15b995SBarry Smith   ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
251ae15b995SBarry Smith   ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
252431e710aSBarry Smith 
2537529efd4SKris Buschelman   ierr = MatILUDTFactor_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
25471c2f376SKris Buschelman 
25586bacbd2SBarry Smith   PetscFunctionReturn(0);
25660ec86bdSBarry Smith #endif
25786bacbd2SBarry Smith }
25886bacbd2SBarry Smith 
2594a2ae208SSatish Balay #undef __FUNCT__
260b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
261dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
262289bc588SBarry Smith {
263416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
264289bc588SBarry Smith   IS                 isicol;
2656849ba73SBarry Smith   PetscErrorCode     ierr;
266899cda47SBarry Smith   PetscInt           *r,*ic,i,n=A->rmap.n,*ai=a->i,*aj=a->j;
2671393bc97SHong Zhang   PetscInt           *bi,*bj,*ajtmp;
2681393bc97SHong Zhang   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
2699e046878SBarry Smith   PetscReal          f;
2704875ba38SHong Zhang   PetscInt           nlnk,*lnk,k,**bi_ptr;
271a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
2721393bc97SHong Zhang   PetscBT            lnkbt;
273289bc588SBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
275899cda47SBarry Smith   if (A->rmap.N != A->cmap.N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2764c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
277ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
278ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
279289bc588SBarry Smith 
280289bc588SBarry Smith   /* get new row pointers */
2811393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
2821393bc97SHong Zhang   bi[0] = 0;
2831393bc97SHong Zhang 
2841393bc97SHong Zhang   /* bdiag is location of diagonal in factor */
2851393bc97SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
2861393bc97SHong Zhang   bdiag[0] = 0;
2871393bc97SHong Zhang 
2881393bc97SHong Zhang   /* linked list for storing column indices of the active row */
2891393bc97SHong Zhang   nlnk = n + 1;
2901393bc97SHong Zhang   ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
2911393bc97SHong Zhang 
2924875ba38SHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt)+n*sizeof(PetscInt**),&im);CHKERRQ(ierr);
2931393bc97SHong Zhang   bi_ptr = (PetscInt**)(im + n);
2941393bc97SHong Zhang 
2951393bc97SHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
296d3d32019SBarry Smith   f = info->fill;
297a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
2981393bc97SHong Zhang   current_space = free_space;
299289bc588SBarry Smith 
300289bc588SBarry Smith   for (i=0; i<n; i++) {
3011393bc97SHong Zhang     /* copy previous fill into linked list */
302289bc588SBarry Smith     nzi = 0;
3031393bc97SHong Zhang     nnz = ai[r[i]+1] - ai[r[i]];
3041393bc97SHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
3051393bc97SHong Zhang     ajtmp = aj + ai[r[i]];
306afcc9904SHong Zhang     ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr);
3071393bc97SHong Zhang     nzi += nlnk;
3081393bc97SHong Zhang 
3091393bc97SHong Zhang     /* add pivot rows into linked list */
3101393bc97SHong Zhang     row = lnk[n];
3111393bc97SHong Zhang     while (row < i) {
312d1170560SHong Zhang       nzbd    = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */
313d1170560SHong Zhang       ajtmp   = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
314d1170560SHong Zhang       ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr);
3151393bc97SHong Zhang       nzi += nlnk;
3161393bc97SHong Zhang       row  = lnk[row];
317289bc588SBarry Smith     }
3181393bc97SHong Zhang     bi[i+1] = bi[i] + nzi;
3191393bc97SHong Zhang     im[i]   = nzi;
3201393bc97SHong Zhang 
3211393bc97SHong Zhang     /* mark bdiag */
3221393bc97SHong Zhang     nzbd = 0;
3231393bc97SHong Zhang     nnz  = nzi;
3241393bc97SHong Zhang     k    = lnk[n];
3251393bc97SHong Zhang     while (nnz-- && k < i){
3261393bc97SHong Zhang       nzbd++;
3271393bc97SHong Zhang       k = lnk[k];
3281393bc97SHong Zhang     }
3291393bc97SHong Zhang     bdiag[i] = bi[i] + nzbd;
3301393bc97SHong Zhang 
3311393bc97SHong Zhang     /* if free space is not available, make more free space */
3321393bc97SHong Zhang     if (current_space->local_remaining<nzi) {
3331393bc97SHong Zhang       nnz = (n - i)*nzi; /* estimated and max additional space needed */
334a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
3351393bc97SHong Zhang       reallocs++;
3361393bc97SHong Zhang     }
3371393bc97SHong Zhang 
3381393bc97SHong Zhang     /* copy data into free space, then initialize lnk */
3391393bc97SHong Zhang     ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
3401393bc97SHong Zhang     bi_ptr[i] = current_space->array;
3411393bc97SHong Zhang     current_space->array           += nzi;
3421393bc97SHong Zhang     current_space->local_used      += nzi;
3431393bc97SHong Zhang     current_space->local_remaining -= nzi;
344289bc588SBarry Smith   }
3456cf91177SBarry Smith #if defined(PETSC_USE_INFO)
346464e8d44SSatish Balay   if (ai[n] != 0) {
3471393bc97SHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
348ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
349ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
350ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
351ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
35205bf559cSSatish Balay   } else {
353ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
35405bf559cSSatish Balay   }
35563ba0a88SBarry Smith #endif
3562fb3ed76SBarry Smith 
357898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
358898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3591987afe7SBarry Smith 
3601393bc97SHong Zhang   /* destroy list of free space and other temporary array(s) */
3611393bc97SHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
362a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
3631393bc97SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
3644875ba38SHong Zhang   ierr = PetscFree(im);CHKERRQ(ierr);
365289bc588SBarry Smith 
366289bc588SBarry Smith   /* put together the new matrix */
367f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,B);CHKERRQ(ierr);
368f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
369f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
370ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
37152e6d16bSBarry Smith   ierr = PetscLogObjectParent(*B,isicol);CHKERRQ(ierr);
372416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
373e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
374e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
3757c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
3761393bc97SHong Zhang   ierr          = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
3771393bc97SHong Zhang   b->j          = bj;
3781393bc97SHong Zhang   b->i          = bi;
3791393bc97SHong Zhang   b->diag       = bdiag;
380416022c9SBarry Smith   b->ilen       = 0;
381416022c9SBarry Smith   b->imax       = 0;
382416022c9SBarry Smith   b->row        = isrow;
383416022c9SBarry Smith   b->col        = iscol;
384c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
385c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
38682bf6240SBarry Smith   b->icol       = isicol;
38787828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3881393bc97SHong Zhang 
3891393bc97SHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
3901393bc97SHong Zhang   ierr = PetscLogObjectMemory(*B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
3911393bc97SHong Zhang   b->maxnz = b->nz = bi[n] ;
3927fd98bd6SLois Curfman McInnes 
393329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
394418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
395ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
396703deb49SSatish Balay 
397e93ccc4dSBarry Smith   if (ai[n] != 0) {
3981393bc97SHong Zhang     (*B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
399eea2913aSSatish Balay   } else {
400eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
401eea2913aSSatish Balay   }
4024846f1f5SKris Buschelman   ierr = MatLUFactorSymbolic_Inode(A,isrow,iscol,info,B);CHKERRQ(ierr);
40371c2f376SKris Buschelman   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
405289bc588SBarry Smith }
4061393bc97SHong Zhang 
4075b5bc046SBarry Smith /*
4085b5bc046SBarry Smith     Trouble in factorization, should we dump the original matrix?
4095b5bc046SBarry Smith */
4105b5bc046SBarry Smith #undef __FUNCT__
4115b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix"
4125b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A)
4135b5bc046SBarry Smith {
4145b5bc046SBarry Smith   PetscErrorCode ierr;
4155b5bc046SBarry Smith   PetscTruth     flg;
4165b5bc046SBarry Smith 
4175b5bc046SBarry Smith   PetscFunctionBegin;
4185b5bc046SBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_factor_dump_on_error",&flg);CHKERRQ(ierr);
4195b5bc046SBarry Smith   if (flg) {
4205b5bc046SBarry Smith     PetscViewer viewer;
4215b5bc046SBarry Smith     char        filename[PETSC_MAX_PATH_LEN];
4225b5bc046SBarry Smith 
4235b5bc046SBarry Smith     ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr);
4245b5bc046SBarry Smith     ierr = PetscViewerBinaryOpen(A->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4255b5bc046SBarry Smith     ierr = MatView(A,viewer);CHKERRQ(ierr);
4265b5bc046SBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
4275b5bc046SBarry Smith   }
4285b5bc046SBarry Smith   PetscFunctionReturn(0);
4295b5bc046SBarry Smith }
4305b5bc046SBarry Smith 
4310a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
4324a2ae208SSatish Balay #undef __FUNCT__
4334a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
434af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
435289bc588SBarry Smith {
436416022c9SBarry Smith   Mat            C=*B;
437416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
43882bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4396849ba73SBarry Smith   PetscErrorCode ierr;
440899cda47SBarry Smith   PetscInt       *r,*ic,i,j,n=A->rmap.n,*bi=b->i,*bj=b->j;
441b3bf805bSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,row,*ics;
442d2276718SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj;
4435b5bc046SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4446a7c8fc2SHong Zhang   PetscScalar    d;
4456a7c8fc2SHong Zhang   PetscReal      rs;
446b3bf805bSHong Zhang   LUShift_Ctx    sctx;
44742898933SHong Zhang   PetscInt       newshift,*ddiag;
448289bc588SBarry Smith 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
45078b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
45178b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
45287828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
45387828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
454010a20caSHong Zhang   rtmps = rtmp; ics = ic;
455289bc588SBarry Smith 
456ada7143aSSatish Balay   sctx.shift_top  = 0;
457ada7143aSSatish Balay   sctx.nshift_max = 0;
458ada7143aSSatish Balay   sctx.shift_lo   = 0;
459ada7143aSSatish Balay   sctx.shift_hi   = 0;
460ada7143aSSatish Balay 
461118f5789SBarry Smith   /* if both shift schemes are chosen by user, only use info->shiftpd */
462118f5789SBarry Smith   if (info->shiftpd && info->shiftnz) info->shiftnz = 0.0;
4631a890ab1SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
4649deb78dcSHong Zhang     PetscInt *aai = a->i;
46542898933SHong Zhang     ddiag          = a->diag;
466118f5789SBarry Smith     sctx.shift_top = 0;
4676cc28720Svictorle     for (i=0; i<n; i++) {
4689f95998fSHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
4699f95998fSHong Zhang       d  = (a->a)[ddiag[i]];
4701a890ab1SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
471010a20caSHong Zhang       v  = a->a+aai[i];
472e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
473e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4741a890ab1SHong Zhang 	rs += PetscAbsScalar(v[j]);
4751a890ab1SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
4766cc28720Svictorle     }
477c3af1dc1SHong Zhang     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
478118f5789SBarry Smith     sctx.shift_top    *= 1.1;
479d2276718SHong Zhang     sctx.nshift_max   = 5;
480d2276718SHong Zhang     sctx.shift_lo     = 0.;
481d2276718SHong Zhang     sctx.shift_hi     = 1.;
482a0a356efSHong Zhang   }
483a0a356efSHong Zhang 
484a0a356efSHong Zhang   sctx.shift_amount = 0;
485a0a356efSHong Zhang   sctx.nshift       = 0;
486085256b3SBarry Smith   do {
487d2276718SHong Zhang     sctx.lushift = PETSC_FALSE;
488289bc588SBarry Smith     for (i=0; i<n; i++){
489b3bf805bSHong Zhang       nz    = bi[i+1] - bi[i];
490b3bf805bSHong Zhang       bjtmp = bj + bi[i];
491b3bf805bSHong Zhang       for  (j=0; j<nz; j++) rtmps[bjtmp[j]] = 0.0;
492289bc588SBarry Smith 
493289bc588SBarry Smith       /* load in initial (unfactored row) */
494416022c9SBarry Smith       nz    = a->i[r[i]+1] - a->i[r[i]];
495b3bf805bSHong Zhang       ajtmp = a->j + a->i[r[i]];
496010a20caSHong Zhang       v     = a->a + a->i[r[i]];
497085256b3SBarry Smith       for (j=0; j<nz; j++) {
498b3bf805bSHong Zhang         rtmp[ics[ajtmp[j]]] = v[j];
499085256b3SBarry Smith       }
500d2276718SHong Zhang       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */
501289bc588SBarry Smith 
502b3bf805bSHong Zhang       row = *bjtmp++;
503f2582111SSatish Balay       while  (row < i) {
5048c37ef55SBarry Smith         pc = rtmp + row;
5058c37ef55SBarry Smith         if (*pc != 0.0) {
506010a20caSHong Zhang           pv         = b->a + diag_offset[row];
507010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50835aab85fSBarry Smith           multiplier = *pc / *pv++;
5098c37ef55SBarry Smith           *pc        = multiplier;
510b3bf805bSHong Zhang           nz         = bi[row+1] - diag_offset[row] - 1;
511f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
512efee365bSSatish Balay           ierr = PetscLogFlops(2*nz);CHKERRQ(ierr);
513289bc588SBarry Smith         }
514b3bf805bSHong Zhang         row = *bjtmp++;
515289bc588SBarry Smith       }
516416022c9SBarry Smith       /* finished row so stick it into b->a */
517b3bf805bSHong Zhang       pv   = b->a + bi[i] ;
518b3bf805bSHong Zhang       pj   = b->j + bi[i] ;
519b3bf805bSHong Zhang       nz   = bi[i+1] - bi[i];
520b3bf805bSHong Zhang       diag = diag_offset[i] - bi[i];
521d3d32019SBarry Smith       rs   = 0.0;
522d3d32019SBarry Smith       for (j=0; j<nz; j++) {
523d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
524d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
525d3d32019SBarry Smith       }
526d2276718SHong Zhang 
527b3bf805bSHong Zhang       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
528d2276718SHong Zhang       sctx.rs  = rs;
529d2276718SHong Zhang       sctx.pv  = pv[diag];
5305b5bc046SBarry Smith       ierr = MatLUCheckShift_inline(info,sctx,i,a->diag,newshift);CHKERRQ(ierr);
5315b5bc046SBarry Smith       if (newshift == 1) break;
53235aab85fSBarry Smith     }
533d2276718SHong Zhang 
534118f5789SBarry Smith     if (info->shiftpd && !sctx.lushift && info->shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
5356cc28720Svictorle       /*
5366c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5376cc28720Svictorle        * then try lower shift
5386cc28720Svictorle        */
539d2276718SHong Zhang       sctx.shift_hi        = info->shift_fraction;
540d2276718SHong Zhang       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
541d2276718SHong Zhang       sctx.shift_amount    = info->shift_fraction * sctx.shift_top;
542d2276718SHong Zhang       sctx.lushift         = PETSC_TRUE;
543d2276718SHong Zhang       sctx.nshift++;
5446cc28720Svictorle     }
545d2276718SHong Zhang   } while (sctx.lushift);
5460f11f9aeSSatish Balay 
54735aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
54835aab85fSBarry Smith   for (i=0; i<n; i++) {
549010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55035aab85fSBarry Smith   }
55135aab85fSBarry Smith 
552606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
55378b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
55478b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
555416022c9SBarry Smith   C->factor = FACTOR_LU;
5567dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
557c456f294SBarry Smith   C->assembled = PETSC_TRUE;
558899cda47SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
559d2276718SHong Zhang   if (sctx.nshift){
560118f5789SBarry Smith     if (info->shiftnz) {
561ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
562118f5789SBarry Smith     } else if (info->shiftpd) {
563ae15b995SBarry Smith       ierr = PetscInfo4(0,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,info->shift_fraction,sctx.shift_top);CHKERRQ(ierr);
5646cc28720Svictorle     }
5650697b6caSHong Zhang   }
5663a40ed3dSBarry Smith   PetscFunctionReturn(0);
567289bc588SBarry Smith }
5680889b5dcSvictorle 
5690889b5dcSvictorle #undef __FUNCT__
5700889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
571dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5720889b5dcSvictorle {
5730889b5dcSvictorle   PetscFunctionBegin;
5740889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5750889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5760889b5dcSvictorle   PetscFunctionReturn(0);
5770889b5dcSvictorle }
5780889b5dcSvictorle 
5790889b5dcSvictorle 
5800a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5814a2ae208SSatish Balay #undef __FUNCT__
5824a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
583dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
584da3a660dSBarry Smith {
585dfbe8321SBarry Smith   PetscErrorCode ierr;
586416022c9SBarry Smith   Mat            C;
587416022c9SBarry Smith 
5883a40ed3dSBarry Smith   PetscFunctionBegin;
5899e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
590af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
591273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
59252e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr);
5933a40ed3dSBarry Smith   PetscFunctionReturn(0);
594da3a660dSBarry Smith }
5950a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5964a2ae208SSatish Balay #undef __FUNCT__
5974a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
598dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
5998c37ef55SBarry Smith {
600416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
601416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6026849ba73SBarry Smith   PetscErrorCode ierr;
603899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
60438baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
60587828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6068c37ef55SBarry Smith 
6073a40ed3dSBarry Smith   PetscFunctionBegin;
6083a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
60991bf9adeSBarry Smith 
6101ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
612416022c9SBarry Smith   tmp  = a->solve_work;
6138c37ef55SBarry Smith 
6143b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6153b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6168c37ef55SBarry Smith 
6178c37ef55SBarry Smith   /* forward solve the lower triangular */
6188c37ef55SBarry Smith   tmp[0] = b[*r++];
619010a20caSHong Zhang   tmps   = tmp;
6208c37ef55SBarry Smith   for (i=1; i<n; i++) {
621010a20caSHong Zhang     v   = aa + ai[i] ;
622010a20caSHong Zhang     vi  = aj + ai[i] ;
62353e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
62453e7425aSBarry Smith     sum = b[*r++];
625ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6268c37ef55SBarry Smith     tmp[i] = sum;
6278c37ef55SBarry Smith   }
6288c37ef55SBarry Smith 
6298c37ef55SBarry Smith   /* backward solve the upper triangular */
6308c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
631010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
632010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
633416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6348c37ef55SBarry Smith     sum = tmp[i];
635ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
636010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6378c37ef55SBarry Smith   }
6388c37ef55SBarry Smith 
6393b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6403b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6411ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6421ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
643899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
6443a40ed3dSBarry Smith   PetscFunctionReturn(0);
6458c37ef55SBarry Smith }
646026e39d0SSatish Balay 
6472bebee5dSHong Zhang #undef __FUNCT__
6482bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ"
6492bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
6502bebee5dSHong Zhang {
6512bebee5dSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6522bebee5dSHong Zhang   IS             iscol = a->col,isrow = a->row;
6532bebee5dSHong Zhang   PetscErrorCode ierr;
654899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
6552bebee5dSHong Zhang   PetscInt       nz,*rout,*cout,neq;
6562bebee5dSHong Zhang   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6572bebee5dSHong Zhang 
6582bebee5dSHong Zhang   PetscFunctionBegin;
6592bebee5dSHong Zhang   if (!n) PetscFunctionReturn(0);
6602bebee5dSHong Zhang 
6612bebee5dSHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
6622bebee5dSHong Zhang   ierr = MatGetArray(X,&x);CHKERRQ(ierr);
6632bebee5dSHong Zhang 
6642bebee5dSHong Zhang   tmp  = a->solve_work;
6652bebee5dSHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6662bebee5dSHong Zhang   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
6672bebee5dSHong Zhang 
6682bebee5dSHong Zhang   for (neq=0; neq<n; neq++){
6692bebee5dSHong Zhang     /* forward solve the lower triangular */
6702bebee5dSHong Zhang     tmp[0] = b[r[0]];
6712bebee5dSHong Zhang     tmps   = tmp;
6722bebee5dSHong Zhang     for (i=1; i<n; i++) {
6732bebee5dSHong Zhang       v   = aa + ai[i] ;
6742bebee5dSHong Zhang       vi  = aj + ai[i] ;
6752bebee5dSHong Zhang       nz  = a->diag[i] - ai[i];
6762bebee5dSHong Zhang       sum = b[r[i]];
6772bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6782bebee5dSHong Zhang       tmp[i] = sum;
6792bebee5dSHong Zhang     }
6802bebee5dSHong Zhang     /* backward solve the upper triangular */
6812bebee5dSHong Zhang     for (i=n-1; i>=0; i--){
6822bebee5dSHong Zhang       v   = aa + a->diag[i] + 1;
6832bebee5dSHong Zhang       vi  = aj + a->diag[i] + 1;
6842bebee5dSHong Zhang       nz  = ai[i+1] - a->diag[i] - 1;
6852bebee5dSHong Zhang       sum = tmp[i];
6862bebee5dSHong Zhang       SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6872bebee5dSHong Zhang       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
6882bebee5dSHong Zhang     }
6892bebee5dSHong Zhang 
6902bebee5dSHong Zhang     b += n;
6912bebee5dSHong Zhang     x += n;
6922bebee5dSHong Zhang   }
6932bebee5dSHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6942bebee5dSHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6952bebee5dSHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
6962bebee5dSHong Zhang   ierr = MatRestoreArray(X,&x);CHKERRQ(ierr);
6972bebee5dSHong Zhang   ierr = PetscLogFlops(n*(2*a->nz - n));CHKERRQ(ierr);
6982bebee5dSHong Zhang   PetscFunctionReturn(0);
6992bebee5dSHong Zhang }
7002bebee5dSHong Zhang 
701930ae53cSBarry Smith /* ----------------------------------------------------------- */
7024a2ae208SSatish Balay #undef __FUNCT__
7034a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
704dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
705930ae53cSBarry Smith {
706930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7076849ba73SBarry Smith   PetscErrorCode ierr;
708899cda47SBarry Smith   PetscInt       n = A->rmap.n,*ai = a->i,*aj = a->j,*adiag = a->diag;
709362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
710aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
71138baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
712362ced78SSatish Balay   PetscScalar    *v,sum;
713d85afc42SSatish Balay #endif
714930ae53cSBarry Smith 
7153a40ed3dSBarry Smith   PetscFunctionBegin;
7163a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
717930ae53cSBarry Smith 
7181ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7191ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
720930ae53cSBarry Smith 
721aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
7221c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
7231c4feecaSSatish Balay #else
724930ae53cSBarry Smith   /* forward solve the lower triangular */
725930ae53cSBarry Smith   x[0] = b[0];
726930ae53cSBarry Smith   for (i=1; i<n; i++) {
727930ae53cSBarry Smith     ai_i = ai[i];
728930ae53cSBarry Smith     v    = aa + ai_i;
729930ae53cSBarry Smith     vi   = aj + ai_i;
730930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
731930ae53cSBarry Smith     sum  = b[i];
732930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
733930ae53cSBarry Smith     x[i] = sum;
734930ae53cSBarry Smith   }
735930ae53cSBarry Smith 
736930ae53cSBarry Smith   /* backward solve the upper triangular */
737930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
738930ae53cSBarry Smith     adiag_i = adiag[i];
739d708749eSBarry Smith     v       = aa + adiag_i + 1;
740d708749eSBarry Smith     vi      = aj + adiag_i + 1;
741930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
742930ae53cSBarry Smith     sum     = x[i];
743930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
744930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
745930ae53cSBarry Smith   }
7461c4feecaSSatish Balay #endif
747899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz - A->cmap.n);CHKERRQ(ierr);
7481ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7491ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7503a40ed3dSBarry Smith   PetscFunctionReturn(0);
751930ae53cSBarry Smith }
752930ae53cSBarry Smith 
7534a2ae208SSatish Balay #undef __FUNCT__
7544a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
755dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
756da3a660dSBarry Smith {
757416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
758416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7596849ba73SBarry Smith   PetscErrorCode ierr;
760899cda47SBarry Smith   PetscInt       *r,*c,i, n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
76138baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
76287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
763da3a660dSBarry Smith 
7643a40ed3dSBarry Smith   PetscFunctionBegin;
76578b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
766da3a660dSBarry Smith 
7671ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7681ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
769416022c9SBarry Smith   tmp  = a->solve_work;
770da3a660dSBarry Smith 
7713b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7723b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
773da3a660dSBarry Smith 
774da3a660dSBarry Smith   /* forward solve the lower triangular */
775da3a660dSBarry Smith   tmp[0] = b[*r++];
776da3a660dSBarry Smith   for (i=1; i<n; i++) {
777010a20caSHong Zhang     v   = aa + ai[i] ;
778010a20caSHong Zhang     vi  = aj + ai[i] ;
779416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
780da3a660dSBarry Smith     sum = b[*r++];
781010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
782da3a660dSBarry Smith     tmp[i] = sum;
783da3a660dSBarry Smith   }
784da3a660dSBarry Smith 
785da3a660dSBarry Smith   /* backward solve the upper triangular */
786da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
787010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
788010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
789416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
790da3a660dSBarry Smith     sum = tmp[i];
791010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
792010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
793da3a660dSBarry Smith     x[*c--] += tmp[i];
794da3a660dSBarry Smith   }
795da3a660dSBarry Smith 
7963b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7973b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7981ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7991ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
800efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
801898183ecSLois Curfman McInnes 
8023a40ed3dSBarry Smith   PetscFunctionReturn(0);
803da3a660dSBarry Smith }
804da3a660dSBarry Smith /* -------------------------------------------------------------------*/
8054a2ae208SSatish Balay #undef __FUNCT__
8064a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
807dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
808da3a660dSBarry Smith {
809416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8102235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8116849ba73SBarry Smith   PetscErrorCode ierr;
812899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
81338baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
81487828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
815da3a660dSBarry Smith 
8163a40ed3dSBarry Smith   PetscFunctionBegin;
8171ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
819416022c9SBarry Smith   tmp  = a->solve_work;
820da3a660dSBarry Smith 
8212235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8222235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
823da3a660dSBarry Smith 
824da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
8252235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
826da3a660dSBarry Smith 
827da3a660dSBarry Smith   /* forward solve the U^T */
828da3a660dSBarry Smith   for (i=0; i<n; i++) {
829010a20caSHong Zhang     v   = aa + diag[i] ;
830010a20caSHong Zhang     vi  = aj + diag[i] + 1;
831f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
832c38d4ed2SBarry Smith     s1  = tmp[i];
833ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
834da3a660dSBarry Smith     while (nz--) {
835010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
836da3a660dSBarry Smith     }
837c38d4ed2SBarry Smith     tmp[i] = s1;
838da3a660dSBarry Smith   }
839da3a660dSBarry Smith 
840da3a660dSBarry Smith   /* backward solve the L^T */
841da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
842010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
843010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
844f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
845f1af5d2fSBarry Smith     s1  = tmp[i];
846da3a660dSBarry Smith     while (nz--) {
847010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
848da3a660dSBarry Smith     }
849da3a660dSBarry Smith   }
850da3a660dSBarry Smith 
851da3a660dSBarry Smith   /* copy tmp into x according to permutation */
852da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
853da3a660dSBarry Smith 
8542235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8552235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8561ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8571ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
858da3a660dSBarry Smith 
859899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz-A->cmap.n);CHKERRQ(ierr);
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861da3a660dSBarry Smith }
862da3a660dSBarry Smith 
8634a2ae208SSatish Balay #undef __FUNCT__
8644a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
865dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
866da3a660dSBarry Smith {
867416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8682235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8696849ba73SBarry Smith   PetscErrorCode ierr;
870899cda47SBarry Smith   PetscInt       *r,*c,i,n = A->rmap.n,*vi,*ai = a->i,*aj = a->j;
87138baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
87287828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8736abc6512SBarry Smith 
8743a40ed3dSBarry Smith   PetscFunctionBegin;
87523bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8766abc6512SBarry Smith 
8771ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8781ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
879416022c9SBarry Smith   tmp = a->solve_work;
8806abc6512SBarry Smith 
8812235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8822235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8836abc6512SBarry Smith 
8846abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8852235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8866abc6512SBarry Smith 
8876abc6512SBarry Smith   /* forward solve the U^T */
8886abc6512SBarry Smith   for (i=0; i<n; i++) {
889010a20caSHong Zhang     v   = aa + diag[i] ;
890010a20caSHong Zhang     vi  = aj + diag[i] + 1;
891f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8926abc6512SBarry Smith     tmp[i] *= *v++;
8936abc6512SBarry Smith     while (nz--) {
894010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8956abc6512SBarry Smith     }
8966abc6512SBarry Smith   }
8976abc6512SBarry Smith 
8986abc6512SBarry Smith   /* backward solve the L^T */
8996abc6512SBarry Smith   for (i=n-1; i>=0; i--){
900010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
901010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
902f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
9036abc6512SBarry Smith     while (nz--) {
904010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
9056abc6512SBarry Smith     }
9066abc6512SBarry Smith   }
9076abc6512SBarry Smith 
9086abc6512SBarry Smith   /* copy tmp into x according to permutation */
9096abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
9106abc6512SBarry Smith 
9112235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
9122235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
9131ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
9141ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9156abc6512SBarry Smith 
916efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
9173a40ed3dSBarry Smith   PetscFunctionReturn(0);
918da3a660dSBarry Smith }
9199e25ed09SBarry Smith /* ----------------------------------------------------------------*/
920*09f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
92108480c60SBarry Smith 
9224a2ae208SSatish Balay #undef __FUNCT__
9234a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
924dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
9259e25ed09SBarry Smith {
926416022c9SBarry Smith   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
9279e25ed09SBarry Smith   IS                 isicol;
9286849ba73SBarry Smith   PetscErrorCode     ierr;
929*09f38230SBarry Smith   PetscInt           *r,*ic,n=A->rmap.n,*ai=a->i,*aj=a->j,d;
9305a8e39fbSHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
9315a8e39fbSHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,len, reallocs=0,dcount=0;
9325a8e39fbSHong Zhang   PetscInt           i,levels,diagonal_fill;
93377c4ece6SBarry Smith   PetscTruth         col_identity,row_identity;
934329f5518SBarry Smith   PetscReal          f;
9355a8e39fbSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
9365a8e39fbSHong Zhang   PetscBT            lnkbt;
9375a8e39fbSHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
938a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
939a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
940*09f38230SBarry Smith   PetscTruth         missing;
9419e25ed09SBarry Smith 
9423a40ed3dSBarry Smith   PetscFunctionBegin;
9436cf997b0SBarry Smith   f             = info->fill;
94438baddfdSBarry Smith   levels        = (PetscInt)info->levels;
94538baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9464c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
94782bf6240SBarry Smith 
94808480c60SBarry Smith   /* special case that simply copies fill pattern */
949be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
950be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
95186bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9522e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
95308480c60SBarry Smith     (*fact)->factor                 = FACTOR_LU;
954f3a39becSBarry Smith     (*fact)->info.factor_mallocs    = 0;
955f3a39becSBarry Smith     (*fact)->info.fill_ratio_given  = info->fill;
956f3a39becSBarry Smith     (*fact)->info.fill_ratio_needed = 1.0;
95708480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
958*09f38230SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact,&missing,&d);CHKERRQ(ierr);
959*09f38230SBarry Smith     if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
96008480c60SBarry Smith     b->row              = isrow;
96108480c60SBarry Smith     b->col              = iscol;
96282bf6240SBarry Smith     b->icol             = isicol;
963899cda47SBarry Smith     ierr                = PetscMalloc(((*fact)->rmap.n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
964f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
965c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
966c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9673a40ed3dSBarry Smith     PetscFunctionReturn(0);
96808480c60SBarry Smith   }
96908480c60SBarry Smith 
970898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
971898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9729e25ed09SBarry Smith 
9739e25ed09SBarry Smith   /* get new row pointers */
9745a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
9755a8e39fbSHong Zhang   bi[0] = 0;
9765a8e39fbSHong Zhang   /* bdiag is location of diagonal in factor */
9775a8e39fbSHong Zhang   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
9785a8e39fbSHong Zhang   bdiag[0]  = 0;
9796cf997b0SBarry Smith 
9805a8e39fbSHong Zhang   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr);
9815a8e39fbSHong Zhang   bjlvl_ptr = (PetscInt**)(bj_ptr + n);
9825a8e39fbSHong Zhang 
9835a8e39fbSHong Zhang   /* create a linked list for storing column indices of the active row */
9845a8e39fbSHong Zhang   nlnk = n + 1;
9855a8e39fbSHong Zhang   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
9865a8e39fbSHong Zhang 
9875a8e39fbSHong Zhang   /* initial FreeSpace size is f*(ai[n]+1) */
988a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
9895a8e39fbSHong Zhang   current_space = free_space;
990a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
9915a8e39fbSHong Zhang   current_space_lvl = free_space_lvl;
9925a8e39fbSHong Zhang 
9935a8e39fbSHong Zhang   for (i=0; i<n; i++) {
9945a8e39fbSHong Zhang     nzi = 0;
9956cf997b0SBarry Smith     /* copy current row into linked list */
9965a8e39fbSHong Zhang     nnz  = ai[r[i]+1] - ai[r[i]];
9975a8e39fbSHong Zhang     if (!nnz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
9985a8e39fbSHong Zhang     cols = aj + ai[r[i]];
9995a8e39fbSHong Zhang     lnk[i] = -1; /* marker to indicate if diagonal exists */
10005a8e39fbSHong Zhang     ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
10015a8e39fbSHong Zhang     nzi += nlnk;
10026cf997b0SBarry Smith 
10036cf997b0SBarry Smith     /* make sure diagonal entry is included */
10045a8e39fbSHong Zhang     if (diagonal_fill && lnk[i] == -1) {
10056cf997b0SBarry Smith       fm = n;
10065a8e39fbSHong Zhang       while (lnk[fm] < i) fm = lnk[fm];
10075a8e39fbSHong Zhang       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
10085a8e39fbSHong Zhang       lnk[fm]    = i;
10095a8e39fbSHong Zhang       lnk_lvl[i] = 0;
10105a8e39fbSHong Zhang       nzi++; dcount++;
10116cf997b0SBarry Smith     }
10126cf997b0SBarry Smith 
10135a8e39fbSHong Zhang     /* add pivot rows into the active row */
10145a8e39fbSHong Zhang     nzbd = 0;
10155a8e39fbSHong Zhang     prow = lnk[n];
10165a8e39fbSHong Zhang     while (prow < i) {
10175a8e39fbSHong Zhang       nnz      = bdiag[prow];
10185a8e39fbSHong Zhang       cols     = bj_ptr[prow] + nnz + 1;
10195a8e39fbSHong Zhang       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
10205a8e39fbSHong Zhang       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
10215a8e39fbSHong Zhang       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
10225a8e39fbSHong Zhang       nzi += nlnk;
10235a8e39fbSHong Zhang       prow = lnk[prow];
10245a8e39fbSHong Zhang       nzbd++;
102556beaf04SBarry Smith     }
10265a8e39fbSHong Zhang     bdiag[i] = nzbd;
10275a8e39fbSHong Zhang     bi[i+1]  = bi[i] + nzi;
1028ecf371e4SBarry Smith 
10295a8e39fbSHong Zhang     /* if free space is not available, make more free space */
10305a8e39fbSHong Zhang     if (current_space->local_remaining<nzi) {
10315a8e39fbSHong Zhang       nnz = nzi*(n - i); /* estimated and max additional space needed */
1032a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
1033a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
10345a8e39fbSHong Zhang       reallocs++;
10355a8e39fbSHong Zhang     }
1036ecf371e4SBarry Smith 
10375a8e39fbSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
10385a8e39fbSHong Zhang     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
10395a8e39fbSHong Zhang     bj_ptr[i]    = current_space->array;
10405a8e39fbSHong Zhang     bjlvl_ptr[i] = current_space_lvl->array;
10415a8e39fbSHong Zhang 
10425a8e39fbSHong Zhang     /* make sure the active row i has diagonal entry */
10435a8e39fbSHong Zhang     if (*(bj_ptr[i]+bdiag[i]) != i) {
104477431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
1045d7ee0231SBarry Smith     try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
10466cf997b0SBarry Smith     }
10475a8e39fbSHong Zhang 
10485a8e39fbSHong Zhang     current_space->array           += nzi;
10495a8e39fbSHong Zhang     current_space->local_used      += nzi;
10505a8e39fbSHong Zhang     current_space->local_remaining -= nzi;
10515a8e39fbSHong Zhang     current_space_lvl->array           += nzi;
10525a8e39fbSHong Zhang     current_space_lvl->local_used      += nzi;
10535a8e39fbSHong Zhang     current_space_lvl->local_remaining -= nzi;
10549e25ed09SBarry Smith   }
10555a8e39fbSHong Zhang 
1056898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1057898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
10585a8e39fbSHong Zhang 
10595a8e39fbSHong Zhang   /* destroy list of free space and other temporary arrays */
10605a8e39fbSHong Zhang   ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1061a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
10625a8e39fbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1063a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
10645a8e39fbSHong Zhang   ierr = PetscFree(bj_ptr);CHKERRQ(ierr);
10659e25ed09SBarry Smith 
10666cf91177SBarry Smith #if defined(PETSC_USE_INFO)
1067f64065bbSBarry Smith   {
10685a8e39fbSHong Zhang     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1069ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
1070ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1071ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr);
1072ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
1073335d9088SBarry Smith     if (diagonal_fill) {
1074ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr);
1075335d9088SBarry Smith     }
1076f64065bbSBarry Smith   }
107763ba0a88SBarry Smith #endif
1078bbb6d6a8SBarry Smith 
10799e25ed09SBarry Smith   /* put together the new matrix */
1080f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,fact);CHKERRQ(ierr);
1081f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,n,n,n,n);CHKERRQ(ierr);
1082f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1083ab93d7beSBarry Smith   ierr = MatSeqAIJSetPreallocation_SeqAIJ(*fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
108452e6d16bSBarry Smith   ierr = PetscLogObjectParent(*fact,isicol);CHKERRQ(ierr);
1085416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1086e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1087e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
10887c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
10895a8e39fbSHong Zhang   len = (bi[n] )*sizeof(PetscScalar);
1090b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
10915a8e39fbSHong Zhang   b->j          = bj;
10925a8e39fbSHong Zhang   b->i          = bi;
10935a8e39fbSHong Zhang   for (i=0; i<n; i++) bdiag[i] += bi[i];
10945a8e39fbSHong Zhang   b->diag       = bdiag;
1095416022c9SBarry Smith   b->ilen       = 0;
1096416022c9SBarry Smith   b->imax       = 0;
1097416022c9SBarry Smith   b->row        = isrow;
1098416022c9SBarry Smith   b->col        = iscol;
1099c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1100c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
110182bf6240SBarry Smith   b->icol       = isicol;
110287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1103416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
11045a8e39fbSHong Zhang      Allocate bdiag, solve_work, new a, new j */
11055a8e39fbSHong Zhang   ierr = PetscLogObjectMemory(*fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr);
11065a8e39fbSHong Zhang   b->maxnz             = b->nz = bi[n] ;
11079e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
1108418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1109ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
11105a8e39fbSHong Zhang   (*fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
111171c2f376SKris Buschelman 
111254e71613SHong Zhang   ierr = MatILUFactorSymbolic_Inode(A,isrow,iscol,info,fact);CHKERRQ(ierr);
111371c2f376SKris Buschelman   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
11144846f1f5SKris Buschelman 
11153a40ed3dSBarry Smith   PetscFunctionReturn(0);
11169e25ed09SBarry Smith }
1117d5d45c9bSBarry Smith 
11183c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1119a6175056SHong Zhang #undef __FUNCT__
1120a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1121af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1122a6175056SHong Zhang {
11232ed133dbSHong Zhang   Mat            C = *B;
11240968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11252ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11262ed133dbSHong Zhang   IS             ip=b->row;
11272ed133dbSHong Zhang   PetscErrorCode ierr;
1128899cda47SBarry Smith   PetscInt       *rip,i,j,mbs=A->rmap.n,*bi=b->i,*bj=b->j,*bcol;
11292ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11302ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1131622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1132540703acSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1133fbf22428SSatish Balay   PetscReal      shiftpd;
1134540703acSHong Zhang   ChShift_Ctx    sctx;
1135540703acSHong Zhang   PetscInt       newshift;
1136d5d45c9bSBarry Smith 
1137a6175056SHong Zhang   PetscFunctionBegin;
1138540703acSHong Zhang   shiftnz   = info->shiftnz;
1139540703acSHong Zhang   shiftpd   = info->shiftpd;
1140ee45ca4aSHong Zhang   zeropivot = info->zeropivot;
1141ee45ca4aSHong Zhang 
11422ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11432ed133dbSHong Zhang 
11442ed133dbSHong Zhang   /* initialization */
11452ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11462ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11472ed133dbSHong Zhang   jl   = il + mbs;
11482ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11492ed133dbSHong Zhang 
1150540703acSHong Zhang   sctx.shift_amount = 0;
1151540703acSHong Zhang   sctx.nshift       = 0;
11522ed133dbSHong Zhang   do {
1153540703acSHong Zhang     sctx.chshift = PETSC_FALSE;
11542ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11552ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1156a6175056SHong Zhang     }
11572ed133dbSHong Zhang 
11582ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1159c05c3958SHong Zhang       bval = ba + bi[k];
11602ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11612ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11622ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11632ed133dbSHong Zhang         col = rip[aj[j]];
11642ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11652ed133dbSHong Zhang           rtmp[col] = aa[j];
1166c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11672ed133dbSHong Zhang         }
11682ed133dbSHong Zhang       }
116939e3d630SHong Zhang       /* shift the diagonal of the matrix */
1170540703acSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
11712ed133dbSHong Zhang 
11722ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11732ed133dbSHong Zhang       dk = rtmp[k];
11742ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11752ed133dbSHong Zhang 
11762ed133dbSHong Zhang       while (i < k){
11772ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11782ed133dbSHong Zhang 
11792ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11802ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11812ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11822ed133dbSHong Zhang         dk += uikdi*ba[ili];
11832ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11842ed133dbSHong Zhang 
11852ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11862ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11872ed133dbSHong Zhang         if (jmin < jmax){
11882ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11892ed133dbSHong Zhang           /* update il and jl for row i */
11902ed133dbSHong Zhang           il[i] = jmin;
11912ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11922ed133dbSHong Zhang         }
11932ed133dbSHong Zhang         i = nexti;
11942ed133dbSHong Zhang       }
11952ed133dbSHong Zhang 
11960697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11970697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11980697b6caSHong Zhang       rs   = 0.0;
11993c9e8598SHong Zhang       jmin = bi[k]+1;
12003c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
12013c9e8598SHong Zhang       if (nz){
12023c9e8598SHong Zhang         bcol = bj + jmin;
12033c9e8598SHong Zhang         while (nz--){
120489f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
120589f845c8SHong Zhang           bcol++;
12063c9e8598SHong Zhang         }
12073c9e8598SHong Zhang       }
1208540703acSHong Zhang 
1209540703acSHong Zhang       sctx.rs = rs;
1210540703acSHong Zhang       sctx.pv = dk;
12115b5bc046SBarry Smith       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
12125b5bc046SBarry Smith       if (newshift == 1) break;
12133c9e8598SHong Zhang 
12143c9e8598SHong Zhang       /* copy data into U(k,:) */
121539e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
121639e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
121739e3d630SHong Zhang       if (jmin < jmax) {
121839e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
121939e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12203c9e8598SHong Zhang         }
122139e3d630SHong Zhang         /* add the k-th row into il and jl */
12223c9e8598SHong Zhang         il[k] = jmin;
12233c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12243c9e8598SHong Zhang       }
12253c9e8598SHong Zhang     }
1226540703acSHong Zhang   } while (sctx.chshift);
12273c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12283c9e8598SHong Zhang 
122939e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12303c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12313c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12323c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
1233899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.n);CHKERRQ(ierr);
1234540703acSHong Zhang   if (sctx.nshift){
1235540703acSHong Zhang     if (shiftnz) {
1236ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1237540703acSHong Zhang     } else if (shiftpd) {
1238ae15b995SBarry Smith       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
12390697b6caSHong Zhang     }
12403c9e8598SHong Zhang   }
12413c9e8598SHong Zhang   PetscFunctionReturn(0);
12423c9e8598SHong Zhang }
1243a6175056SHong Zhang 
1244a6175056SHong Zhang #undef __FUNCT__
1245a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1246dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1247a6175056SHong Zhang {
12480968510aSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1249ed59401aSHong Zhang   Mat_SeqSBAIJ       *b;
1250ed59401aSHong Zhang   Mat                B;
1251dfbe8321SBarry Smith   PetscErrorCode     ierr;
1252653a6975SHong Zhang   PetscTruth         perm_identity;
1253899cda47SBarry Smith   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->rmap.n,*ui;
1254ed59401aSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1255391eac55SHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL;
12565a8e39fbSHong Zhang   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
1257ed59401aSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
1258a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1259a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12607a48dd6fSHong Zhang   PetscBT            lnkbt;
1261a6175056SHong Zhang 
1262a6175056SHong Zhang   PetscFunctionBegin;
1263653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12647a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12657a48dd6fSHong Zhang 
126639e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
126739e3d630SHong Zhang   ui[0] = 0;
126839e3d630SHong Zhang 
1269622e2028SHong Zhang   /* special case that simply copies fill pattern */
1270622e2028SHong Zhang   if (!levels && perm_identity) {
1271ed59401aSHong Zhang     for (i=0; i<am; i++) {
127239e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1273ed59401aSHong Zhang     }
127439e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
127539e3d630SHong Zhang     cols = uj;
1276ed59401aSHong Zhang     for (i=0; i<am; i++) {
1277ed59401aSHong Zhang       aj    = a->j + a->diag[i];
127839e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
127939e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1280ed59401aSHong Zhang     }
128139e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12827a48dd6fSHong Zhang     /* initialization */
12835a8e39fbSHong Zhang     ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
12847a48dd6fSHong Zhang 
12857a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12867a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12871393bc97SHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
12887a48dd6fSHong Zhang     il         = jl + am;
12897a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12907a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12917a48dd6fSHong Zhang     for (i=0; i<am; i++){
12927a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12937a48dd6fSHong Zhang     }
12947a48dd6fSHong Zhang 
12957a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12967a48dd6fSHong Zhang     nlnk = am + 1;
12972ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12987a48dd6fSHong Zhang 
12997a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
1300a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
13017a48dd6fSHong Zhang     current_space = free_space;
1302a1a86e44SBarry Smith     ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
13037a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
13047a48dd6fSHong Zhang 
13057a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
13067a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
13077a48dd6fSHong Zhang       nzk   = 0;
1308622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1309622e2028SHong Zhang       ncols_upper = 0;
1310622e2028SHong Zhang       for (j=0; j<ncols; j++){
13115a8e39fbSHong Zhang         i = *(aj + ai[rip[k]] + j);
13125a8e39fbSHong Zhang         if (rip[i] >= k){ /* only take upper triangular entry */
13135a8e39fbSHong Zhang           ajtmp[ncols_upper] = i;
1314622e2028SHong Zhang           ncols_upper++;
1315622e2028SHong Zhang         }
1316622e2028SHong Zhang       }
13175a8e39fbSHong Zhang       ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13187a48dd6fSHong Zhang       nzk += nlnk;
13197a48dd6fSHong Zhang 
13207a48dd6fSHong Zhang       /* update lnk by computing fill-in for each pivot row to be merged in */
13217a48dd6fSHong Zhang       prow = jl[k]; /* 1st pivot row */
13227a48dd6fSHong Zhang 
13237a48dd6fSHong Zhang       while (prow < k){
13247a48dd6fSHong Zhang         nextprow = jl[prow];
13257a48dd6fSHong Zhang 
13267a48dd6fSHong Zhang         /* merge prow into k-th row */
13277a48dd6fSHong Zhang         jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
13287a48dd6fSHong Zhang         jmax = ui[prow+1];
13297a48dd6fSHong Zhang         ncols = jmax-jmin;
1330ed59401aSHong Zhang         i     = jmin - ui[prow];
1331ed59401aSHong Zhang         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
13325a8e39fbSHong Zhang         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
13335a8e39fbSHong Zhang         j     = *(uj - 1);
13345a8e39fbSHong Zhang         ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr);
13357a48dd6fSHong Zhang         nzk += nlnk;
13367a48dd6fSHong Zhang 
13377a48dd6fSHong Zhang         /* update il and jl for prow */
13387a48dd6fSHong Zhang         if (jmin < jmax){
13397a48dd6fSHong Zhang           il[prow] = jmin;
13407a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13417a48dd6fSHong Zhang         }
13427a48dd6fSHong Zhang         prow = nextprow;
13437a48dd6fSHong Zhang       }
13447a48dd6fSHong Zhang 
13457a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13467a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13477a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13487a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1349a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1350a1a86e44SBarry Smith         ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
13517a48dd6fSHong Zhang         reallocs++;
13527a48dd6fSHong Zhang       }
13537a48dd6fSHong Zhang 
13547a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13552ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13567a48dd6fSHong Zhang 
13577a48dd6fSHong Zhang       /* add the k-th row into il and jl */
135839e3d630SHong Zhang       if (nzk > 1){
13597a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13607a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13617a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13627a48dd6fSHong Zhang       }
13637a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13647a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13657a48dd6fSHong Zhang 
13667a48dd6fSHong Zhang       current_space->array           += nzk;
13677a48dd6fSHong Zhang       current_space->local_used      += nzk;
13687a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13697a48dd6fSHong Zhang 
13707a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13717a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13727a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13737a48dd6fSHong Zhang 
13747a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13757a48dd6fSHong Zhang     }
13767a48dd6fSHong Zhang 
13776cf91177SBarry Smith #if defined(PETSC_USE_INFO)
13787a48dd6fSHong Zhang     if (ai[am] != 0) {
137939e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
1380ae15b995SBarry Smith       ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1381ae15b995SBarry Smith       ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1382ae15b995SBarry Smith       ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
13837a48dd6fSHong Zhang     } else {
1384ae15b995SBarry Smith       ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
13857a48dd6fSHong Zhang     }
138663ba0a88SBarry Smith #endif
13877a48dd6fSHong Zhang 
13887a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13897a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13905a8e39fbSHong Zhang     ierr = PetscFree(ajtmp);CHKERRQ(ierr);
13917a48dd6fSHong Zhang 
13927a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13937a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1394a1a86e44SBarry Smith     ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13952ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1396a1a86e44SBarry Smith     ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
13977a48dd6fSHong Zhang 
139839e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
139939e3d630SHong Zhang 
14007a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1401f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1402f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
1403ed59401aSHong Zhang   B = *fact;
1404ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1405ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
14067a48dd6fSHong Zhang 
1407ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
14087a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
14097a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
14107a48dd6fSHong Zhang   b->j    = uj;
14117a48dd6fSHong Zhang   b->i    = ui;
14127a48dd6fSHong Zhang   b->diag = 0;
14137a48dd6fSHong Zhang   b->ilen = 0;
14147a48dd6fSHong Zhang   b->imax = 0;
14157a48dd6fSHong Zhang   b->row  = perm;
14167a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
14177a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14187a48dd6fSHong Zhang   b->icol = perm;
14197a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
14207a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
142152e6d16bSBarry Smith   ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
14227a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
1423e6b907acSBarry Smith   b->free_a  = PETSC_TRUE;
1424e6b907acSBarry Smith   b->free_ij = PETSC_TRUE;
14257a48dd6fSHong Zhang 
1426ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1427ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1428ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14297a48dd6fSHong Zhang   if (ai[am] != 0) {
1430ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14317a48dd6fSHong Zhang   } else {
1432ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14337a48dd6fSHong Zhang   }
143439e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1435622e2028SHong Zhang   if (perm_identity){
1436ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1437ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1438622e2028SHong Zhang   }
1439a6175056SHong Zhang   PetscFunctionReturn(0);
1440a6175056SHong Zhang }
1441d5d45c9bSBarry Smith 
1442f76d2b81SHong Zhang #undef __FUNCT__
1443f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1444dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1445f76d2b81SHong Zhang {
1446f76d2b81SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
144710c27e3dSHong Zhang   Mat_SeqSBAIJ       *b;
14482ed133dbSHong Zhang   Mat                B;
1449dfbe8321SBarry Smith   PetscErrorCode     ierr;
1450f76d2b81SHong Zhang   PetscTruth         perm_identity;
145110c27e3dSHong Zhang   PetscReal          fill = info->fill;
1452899cda47SBarry Smith   PetscInt           *rip,*riip,i,am=A->rmap.n,*ai=a->i,*aj=a->j,reallocs=0,prow;
145310c27e3dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14542ed133dbSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1455a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
145610c27e3dSHong Zhang   PetscBT            lnkbt;
14572ed133dbSHong Zhang   IS                 iperm;
1458f76d2b81SHong Zhang 
1459f76d2b81SHong Zhang   PetscFunctionBegin;
14602ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1461f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14622ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14632ed133dbSHong Zhang 
1464f76d2b81SHong Zhang   if (!perm_identity){
14652ed133dbSHong Zhang     /* check if perm is symmetric! */
14662ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14672ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14681393bc97SHong Zhang     for (i=0; i<am; i++) {
14692ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14702ed133dbSHong Zhang     }
14712ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14722ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1473f76d2b81SHong Zhang   }
147410c27e3dSHong Zhang 
147510c27e3dSHong Zhang   /* initialization */
14761393bc97SHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
147710c27e3dSHong Zhang   ui[0] = 0;
147810c27e3dSHong Zhang 
147910c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14801393bc97SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14811393bc97SHong Zhang   ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr);
14821393bc97SHong Zhang   il     = jl + am;
14831393bc97SHong Zhang   cols   = il + am;
14841393bc97SHong Zhang   ui_ptr = (PetscInt**)(cols + am);
14851393bc97SHong Zhang   for (i=0; i<am; i++){
14861393bc97SHong Zhang     jl[i] = am; il[i] = 0;
1487f76d2b81SHong Zhang   }
148810c27e3dSHong Zhang 
148910c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14901393bc97SHong Zhang   nlnk = am + 1;
14911393bc97SHong Zhang   ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
149210c27e3dSHong Zhang 
14931393bc97SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
1494a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
149510c27e3dSHong Zhang   current_space = free_space;
149610c27e3dSHong Zhang 
14971393bc97SHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
149810c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
149910c27e3dSHong Zhang     nzk   = 0;
15002ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
15012ed133dbSHong Zhang     ncols_upper = 0;
15022ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1503622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
15042ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
15052ed133dbSHong Zhang         cols[ncols_upper] = i;
15062ed133dbSHong Zhang         ncols_upper++;
15072ed133dbSHong Zhang       }
15082ed133dbSHong Zhang     }
15091393bc97SHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
151010c27e3dSHong Zhang     nzk += nlnk;
151110c27e3dSHong Zhang 
151210c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
151310c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
151410c27e3dSHong Zhang 
151510c27e3dSHong Zhang     while (prow < k){
151610c27e3dSHong Zhang       nextprow = jl[prow];
151710c27e3dSHong Zhang       /* merge prow into k-th row */
15181393bc97SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
151910c27e3dSHong Zhang       jmax = ui[prow+1];
152010c27e3dSHong Zhang       ncols = jmax-jmin;
15211393bc97SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
15225a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr);
152310c27e3dSHong Zhang       nzk += nlnk;
152410c27e3dSHong Zhang 
152510c27e3dSHong Zhang       /* update il and jl for prow */
152610c27e3dSHong Zhang       if (jmin < jmax){
152710c27e3dSHong Zhang         il[prow] = jmin;
15282ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
152910c27e3dSHong Zhang       }
153010c27e3dSHong Zhang       prow = nextprow;
153110c27e3dSHong Zhang     }
153210c27e3dSHong Zhang 
153310c27e3dSHong Zhang     /* if free space is not available, make more free space */
153410c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15351393bc97SHong Zhang       i = am - k + 1; /* num of unfactored rows */
153610c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1537a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
153810c27e3dSHong Zhang       reallocs++;
153910c27e3dSHong Zhang     }
154010c27e3dSHong Zhang 
154110c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15421393bc97SHong Zhang     ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
154310c27e3dSHong Zhang 
154410c27e3dSHong Zhang     /* add the k-th row into il and jl */
154510c27e3dSHong Zhang     if (nzk-1 > 0){
15461393bc97SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
154710c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
154810c27e3dSHong Zhang       il[k] = ui[k] + 1;
154910c27e3dSHong Zhang     }
15502ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
155110c27e3dSHong Zhang     current_space->array           += nzk;
155210c27e3dSHong Zhang     current_space->local_used      += nzk;
155310c27e3dSHong Zhang     current_space->local_remaining -= nzk;
155410c27e3dSHong Zhang 
155510c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
155610c27e3dSHong Zhang   }
155710c27e3dSHong Zhang 
15586cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15591393bc97SHong Zhang   if (ai[am] != 0) {
15601393bc97SHong Zhang     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
1561ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1562ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1563ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
156410c27e3dSHong Zhang   } else {
1565ae15b995SBarry Smith      ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
156610c27e3dSHong Zhang   }
156763ba0a88SBarry Smith #endif
156810c27e3dSHong Zhang 
156910c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
157010c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
157110c27e3dSHong Zhang 
157210c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15731393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1574a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
157510c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
157610c27e3dSHong Zhang 
157710c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1578f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
1579f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
15802ed133dbSHong Zhang   B    = *fact;
15812ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1582ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
158310c27e3dSHong Zhang 
15842ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
158510c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
1586e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1587e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
15881393bc97SHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
158910c27e3dSHong Zhang   b->j    = uj;
159010c27e3dSHong Zhang   b->i    = ui;
159110c27e3dSHong Zhang   b->diag = 0;
159210c27e3dSHong Zhang   b->ilen = 0;
159310c27e3dSHong Zhang   b->imax = 0;
159410c27e3dSHong Zhang   b->row  = perm;
159510c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
159610c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
159710c27e3dSHong Zhang   b->icol = perm;
159810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15991393bc97SHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
16001393bc97SHong Zhang   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
16011393bc97SHong Zhang   b->maxnz = b->nz = ui[am];
160210c27e3dSHong Zhang 
16032ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
16042ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
16052ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
16061393bc97SHong Zhang   if (ai[am] != 0) {
16071393bc97SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
160810c27e3dSHong Zhang   } else {
16092ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
161010c27e3dSHong Zhang   }
161139e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
16122ed133dbSHong Zhang   if (perm_identity){
161310c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
161410c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
16152ed133dbSHong Zhang   }
1616f76d2b81SHong Zhang   PetscFunctionReturn(0);
1617f76d2b81SHong Zhang }
1618