xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision ee45ca4afdde1af4d1deda7cd4dc1a4a63a3ea97)
1289bc588SBarry Smith 
270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h"
31c4feecaSSatish Balay #include "src/inline/dot.h"
4ed480e8bSBarry Smith #include "src/inline/spops.h"
53b2fbd54SBarry Smith 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ"
8dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
93b2fbd54SBarry Smith {
103a40ed3dSBarry Smith   PetscFunctionBegin;
113a40ed3dSBarry Smith 
1229bbc08cSBarry Smith   SETERRQ(PETSC_ERR_SUP,"Code not written");
13aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG)
14d64ed03dSBarry Smith   PetscFunctionReturn(0);
15d64ed03dSBarry Smith #endif
163b2fbd54SBarry Smith }
173b2fbd54SBarry Smith 
1886bacbd2SBarry Smith 
19dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
20dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
2186bacbd2SBarry Smith 
2238baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2dperm(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
239cc1f4e3SBarry Smith EXTERN PetscErrorCode SPARSEKIT2ilutp(PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscReal,PetscReal*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscErrorCode*);
2438baddfdSBarry Smith EXTERN PetscErrorCode SPARSEKIT2msrcsr(PetscInt*,PetscScalar*,PetscInt*,PetscScalar*,PetscInt*,PetscInt*,PetscScalar*,PetscInt*);
2586bacbd2SBarry Smith 
264a2ae208SSatish Balay #undef __FUNCT__
274a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ"
2886bacbd2SBarry Smith   /* ------------------------------------------------------------
2986bacbd2SBarry Smith 
3086bacbd2SBarry Smith           This interface was contribed by Tony Caola
3186bacbd2SBarry Smith 
3286bacbd2SBarry Smith      This routine is an interface to the pivoting drop-tolerance
335bd2ddc7SBarry Smith      ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of
345bd2ddc7SBarry Smith      SPARSEKIT2.
355bd2ddc7SBarry Smith 
365bd2ddc7SBarry Smith      The SPARSEKIT2 routines used here are covered by the GNU
375bd2ddc7SBarry Smith      copyright; see the file gnu in this directory.
3886bacbd2SBarry Smith 
3986bacbd2SBarry Smith      Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their
4086bacbd2SBarry Smith      help in getting this routine ironed out.
4186bacbd2SBarry Smith 
425bd2ddc7SBarry Smith      The major drawback to this routine is that if info->fill is
435bd2ddc7SBarry Smith      not large enough it fails rather than allocating more space;
445bd2ddc7SBarry Smith      this can be fixed by hacking/improving the f2c version of
455bd2ddc7SBarry Smith      Yousef Saad's code.
4686bacbd2SBarry Smith 
4786bacbd2SBarry Smith      ------------------------------------------------------------
4886bacbd2SBarry Smith */
49dfbe8321SBarry Smith PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact)
5086bacbd2SBarry Smith {
5160ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE)
5260ec86bdSBarry Smith   PetscFunctionBegin;
53e005ede5SBarry Smith   SETERRQ(PETSC_ERR_SUP_SYS,"This distribution does not include GNU Copyright code\n\
5460ec86bdSBarry Smith   You can obtain the drop tolerance routines by installing PETSc from\n\
5560ec86bdSBarry Smith   www.mcs.anl.gov/petsc\n");
5660ec86bdSBarry Smith #else
5786bacbd2SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
5807d2056aSBarry Smith   IS             iscolf,isicol,isirow;
5986bacbd2SBarry Smith   PetscTruth     reorder;
609cc1f4e3SBarry Smith   PetscErrorCode ierr,sierr;
619cc1f4e3SBarry Smith   PetscInt       *c,*r,*ic,i,n = A->m;
6238baddfdSBarry Smith   PetscInt       *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j;
6338baddfdSBarry Smith   PetscInt       *ordcol,*iwk,*iperm,*jw;
6438baddfdSBarry Smith   PetscInt       jmax,lfill,job,*o_i,*o_j;
6587828ca2SBarry Smith   PetscScalar    *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a;
6654a8161fSBarry Smith   PetscReal      af;
6786bacbd2SBarry Smith 
6886bacbd2SBarry Smith   PetscFunctionBegin;
6986bacbd2SBarry Smith 
7086bacbd2SBarry Smith   if (info->dt == PETSC_DEFAULT)      info->dt      = .005;
7138baddfdSBarry Smith   if (info->dtcount == PETSC_DEFAULT) info->dtcount = (PetscInt)(1.5*a->rmax);
7286bacbd2SBarry Smith   if (info->dtcol == PETSC_DEFAULT)   info->dtcol   = .01;
7333259db3SBarry Smith   if (info->fill == PETSC_DEFAULT)    info->fill    = ((double)(n*(info->dtcount+1)))/a->nz;
7438baddfdSBarry Smith   lfill   = (PetscInt)(info->dtcount/2.0);
7538baddfdSBarry Smith   jmax    = (PetscInt)(info->fill*a->nz);
7686bacbd2SBarry Smith 
7786bacbd2SBarry Smith 
7886bacbd2SBarry Smith   /* ------------------------------------------------------------
7986bacbd2SBarry Smith      If reorder=.TRUE., then the original matrix has to be
8086bacbd2SBarry Smith      reordered to reflect the user selected ordering scheme, and
8186bacbd2SBarry Smith      then de-reordered so it is in it's original format.
8286bacbd2SBarry Smith      Because Saad's dperm() is NOT in place, we have to copy
8386bacbd2SBarry Smith      the original matrix and allocate more storage. . .
8486bacbd2SBarry Smith      ------------------------------------------------------------
8586bacbd2SBarry Smith   */
8686bacbd2SBarry Smith 
8786bacbd2SBarry Smith   /* set reorder to true if either isrow or iscol is not identity */
8886bacbd2SBarry Smith   ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr);
8986bacbd2SBarry Smith   if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);}
9086bacbd2SBarry Smith   reorder = PetscNot(reorder);
9186bacbd2SBarry Smith 
9286bacbd2SBarry Smith 
9386bacbd2SBarry Smith   /* storage for ilu factor */
9438baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&new_i);CHKERRQ(ierr);
9538baddfdSBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscInt),&new_j);CHKERRQ(ierr);
9687828ca2SBarry Smith   ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr);
9738baddfdSBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ordcol);CHKERRQ(ierr);
9886bacbd2SBarry Smith 
9986bacbd2SBarry Smith   /* ------------------------------------------------------------
10086bacbd2SBarry Smith      Make sure that everything is Fortran formatted (1-Based)
10186bacbd2SBarry Smith      ------------------------------------------------------------
10286bacbd2SBarry Smith   */
103b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
104b19713cbSBarry Smith     old_j[i]++;
10586bacbd2SBarry Smith   }
106b19713cbSBarry Smith   for(i=0;i<n+1;i++) {
107b19713cbSBarry Smith     old_i[i]++;
108b19713cbSBarry Smith   };
109010a20caSHong Zhang 
11086bacbd2SBarry Smith 
11186bacbd2SBarry Smith   if (reorder) {
112c0c2c81eSBarry Smith     ierr = ISGetIndices(iscol,&c);CHKERRQ(ierr);
113c0c2c81eSBarry Smith     ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
114c0c2c81eSBarry Smith     for(i=0;i<n;i++) {
115c0c2c81eSBarry Smith       r[i]  = r[i]+1;
116c0c2c81eSBarry Smith       c[i]  = c[i]+1;
117c0c2c81eSBarry Smith     }
11838baddfdSBarry Smith     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&old_i2);CHKERRQ(ierr);
11938baddfdSBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscInt),&old_j2);CHKERRQ(ierr);
12087828ca2SBarry Smith     ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr);
1215bd2ddc7SBarry Smith     job  = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job);
122c0c2c81eSBarry Smith     for (i=0;i<n;i++) {
123c0c2c81eSBarry Smith       r[i]  = r[i]-1;
124c0c2c81eSBarry Smith       c[i]  = c[i]-1;
125c0c2c81eSBarry Smith     }
126c0c2c81eSBarry Smith     ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr);
127c0c2c81eSBarry Smith     ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
128b19713cbSBarry Smith     o_a = old_a2;
129b19713cbSBarry Smith     o_j = old_j2;
130b19713cbSBarry Smith     o_i = old_i2;
131b19713cbSBarry Smith   } else {
132b19713cbSBarry Smith     o_a = old_a;
133b19713cbSBarry Smith     o_j = old_j;
134b19713cbSBarry Smith     o_i = old_i;
13586bacbd2SBarry Smith   }
13686bacbd2SBarry Smith 
13786bacbd2SBarry Smith   /* ------------------------------------------------------------
13886bacbd2SBarry Smith      Call Saad's ilutp() routine to generate the factorization
13986bacbd2SBarry Smith      ------------------------------------------------------------
14086bacbd2SBarry Smith   */
14186bacbd2SBarry Smith 
14238baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&iperm);CHKERRQ(ierr);
14338baddfdSBarry Smith   ierr = PetscMalloc(2*n*sizeof(PetscInt),&jw);CHKERRQ(ierr);
14487828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr);
14586bacbd2SBarry Smith 
14654a8161fSBarry 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);
1476849ba73SBarry Smith   if (sierr) {
1486849ba73SBarry Smith     switch (sierr) {
14977431f27SBarry Smith       case -3: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
15077431f27SBarry Smith       case -2: SETERRQ2(PETSC_ERR_LIB,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %D",info->fill,jmax);
151e005ede5SBarry Smith       case -5: SETERRQ(PETSC_ERR_LIB,"ilutp(), zero row encountered");
152e005ede5SBarry Smith       case -1: SETERRQ(PETSC_ERR_LIB,"ilutp(), input matrix may be wrong");
15377431f27SBarry Smith       case -4: SETERRQ1(PETSC_ERR_LIB,"ilutp(), illegal info->fill value %D",jmax);
15477431f27SBarry Smith       default: SETERRQ1(PETSC_ERR_LIB,"ilutp(), zero pivot detected on row %D",sierr);
15586bacbd2SBarry Smith     }
15686bacbd2SBarry Smith   }
15786bacbd2SBarry Smith 
15886bacbd2SBarry Smith   ierr = PetscFree(w);CHKERRQ(ierr);
15986bacbd2SBarry Smith   ierr = PetscFree(jw);CHKERRQ(ierr);
16086bacbd2SBarry Smith 
16186bacbd2SBarry Smith   /* ------------------------------------------------------------
16286bacbd2SBarry Smith      Saad's routine gives the result in Modified Sparse Row (msr)
16386bacbd2SBarry Smith      Convert to Compressed Sparse Row format (csr)
16486bacbd2SBarry Smith      ------------------------------------------------------------
16586bacbd2SBarry Smith   */
16686bacbd2SBarry Smith 
16787828ca2SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr);
16838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&iwk);CHKERRQ(ierr);
16986bacbd2SBarry Smith 
1705bd2ddc7SBarry Smith   SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk);
17186bacbd2SBarry Smith 
17286bacbd2SBarry Smith   ierr = PetscFree(iwk);CHKERRQ(ierr);
17386bacbd2SBarry Smith   ierr = PetscFree(wk);CHKERRQ(ierr);
17486bacbd2SBarry Smith 
17586bacbd2SBarry Smith   if (reorder) {
17686bacbd2SBarry Smith     ierr = PetscFree(old_a2);CHKERRQ(ierr);
17786bacbd2SBarry Smith     ierr = PetscFree(old_j2);CHKERRQ(ierr);
17886bacbd2SBarry Smith     ierr = PetscFree(old_i2);CHKERRQ(ierr);
179313ae024SBarry Smith   } else {
180b19713cbSBarry Smith     /* fix permutation of old_j that the factorization introduced */
181f170005cSBarry Smith     for (i=old_i[0]; i<old_i[n]; i++) {
182b19713cbSBarry Smith       old_j[i-1] = iperm[old_j[i-1]-1];
183b19713cbSBarry Smith     }
184b19713cbSBarry Smith   }
185b19713cbSBarry Smith 
186b801d0f9SBarry Smith   /* get rid of the shift to indices starting at 1 */
18786bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
188b19713cbSBarry Smith     old_i[i]--;
189b19713cbSBarry Smith   }
190b19713cbSBarry Smith   for (i=old_i[0];i<old_i[n];i++) {
191b19713cbSBarry Smith     old_j[i]--;
192b19713cbSBarry Smith   }
19386bacbd2SBarry Smith 
194b19713cbSBarry Smith   /* Make the factored matrix 0-based */
19586bacbd2SBarry Smith   for (i=0; i<n+1; i++) {
196b19713cbSBarry Smith     new_i[i]--;
197b19713cbSBarry Smith   }
198b19713cbSBarry Smith   for (i=new_i[0];i<new_i[n];i++) {
199b19713cbSBarry Smith     new_j[i]--;
200b19713cbSBarry Smith   }
20186bacbd2SBarry Smith 
20286bacbd2SBarry Smith   /*-- due to the pivoting, we need to reorder iscol to correctly --*/
20386bacbd2SBarry Smith   /*-- permute the right-hand-side and solution vectors           --*/
2044c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
2054c49b128SBarry Smith   ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr);
206c0c2c81eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
20786bacbd2SBarry Smith   for(i=0; i<n; i++) {
20886bacbd2SBarry Smith     ordcol[i] = ic[iperm[i]-1];
20986bacbd2SBarry Smith   };
21086bacbd2SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
21107d2056aSBarry Smith   ierr = ISDestroy(isicol);CHKERRQ(ierr);
21286bacbd2SBarry Smith 
213c0c2c81eSBarry Smith   ierr = PetscFree(iperm);CHKERRQ(ierr);
214c0c2c81eSBarry Smith 
215329f5518SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr);
21607d2056aSBarry Smith   ierr = PetscFree(ordcol);CHKERRQ(ierr);
21786bacbd2SBarry Smith 
21886bacbd2SBarry Smith   /*----- put together the new matrix -----*/
21986bacbd2SBarry Smith 
220f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
221f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
222f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
22386bacbd2SBarry Smith   (*fact)->factor    = FACTOR_LU;
22486bacbd2SBarry Smith   (*fact)->assembled = PETSC_TRUE;
22586bacbd2SBarry Smith 
22686bacbd2SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
22786bacbd2SBarry Smith   ierr = PetscFree(b->imax);CHKERRQ(ierr);
22886bacbd2SBarry Smith   b->sorted        = PETSC_FALSE;
22907d2056aSBarry Smith   b->singlemalloc  = PETSC_FALSE;
230b19713cbSBarry Smith   /* the next line frees the default space generated by the MatCreate() */
23186bacbd2SBarry Smith   ierr             = PetscFree(b->a);CHKERRQ(ierr);
23286bacbd2SBarry Smith   ierr             = PetscFree(b->ilen);CHKERRQ(ierr);
23386bacbd2SBarry Smith   b->a             = new_a;
23486bacbd2SBarry Smith   b->j             = new_j;
23586bacbd2SBarry Smith   b->i             = new_i;
23686bacbd2SBarry Smith   b->ilen          = 0;
23786bacbd2SBarry Smith   b->imax          = 0;
2381f9e874aSBarry Smith   /*  I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */
239313ae024SBarry Smith   b->row           = isirow;
24086bacbd2SBarry Smith   b->col           = iscolf;
24187828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
24286bacbd2SBarry Smith   b->maxnz = b->nz = new_i[n];
24386bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
24486bacbd2SBarry Smith   (*fact)->info.factor_mallocs = 0;
24586bacbd2SBarry Smith 
24686bacbd2SBarry Smith   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
24786bacbd2SBarry Smith 
24886bacbd2SBarry Smith   /* check out for identical nodes. If found, use inode functions */
2493a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
25086bacbd2SBarry Smith 
251431e710aSBarry Smith   af = ((double)b->nz)/((double)a->nz) + .001;
252b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af);
253b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af);
254b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af);
255b0a32e0cSBarry Smith   PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n");
256431e710aSBarry Smith 
25786bacbd2SBarry Smith   PetscFunctionReturn(0);
25860ec86bdSBarry Smith #endif
25986bacbd2SBarry Smith }
26086bacbd2SBarry Smith 
261289bc588SBarry Smith /*
262289bc588SBarry Smith     Factorization code for AIJ format.
263289bc588SBarry Smith */
2644a2ae208SSatish Balay #undef __FUNCT__
265b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ"
266dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B)
267289bc588SBarry Smith {
268416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
269289bc588SBarry Smith   IS             isicol;
2706849ba73SBarry Smith   PetscErrorCode ierr;
27138baddfdSBarry Smith   PetscInt       *r,*ic,i,n = A->m,*ai = a->i,*aj = a->j;
27238baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*ajtmp,nz;
273418422e8SSatish Balay   PetscInt       *idnew,idx,row,m,fm,nnz,nzi,reallocs = 0,nzbd,*im;
2749e046878SBarry Smith   PetscReal      f;
275289bc588SBarry Smith 
2763a40ed3dSBarry Smith   PetscFunctionBegin;
27729bbc08cSBarry Smith   if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
2784c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
279ac121b3dSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
280ac121b3dSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
281289bc588SBarry Smith 
282289bc588SBarry Smith   /* get new row pointers */
28338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
284010a20caSHong Zhang   ainew[0] = 0;
285289bc588SBarry Smith   /* don't know how many column pointers are needed so estimate */
286d3d32019SBarry Smith   f = info->fill;
28738baddfdSBarry Smith   jmax  = (PetscInt)(f*ai[n]+1);
28838baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
289289bc588SBarry Smith   /* fill is a linked list of nonzeros in active row */
29038baddfdSBarry Smith   ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
2912fb3ed76SBarry Smith   im   = fill + n + 1;
292289bc588SBarry Smith   /* idnew is location of diagonal in factor */
29338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&idnew);CHKERRQ(ierr);
294010a20caSHong Zhang   idnew[0] = 0;
295289bc588SBarry Smith 
296289bc588SBarry Smith   for (i=0; i<n; i++) {
297289bc588SBarry Smith     /* first copy previous fill into linked list */
298289bc588SBarry Smith     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
29929bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
300010a20caSHong Zhang     ajtmp   = aj + ai[r[i]];
301289bc588SBarry Smith     fill[n] = n;
302289bc588SBarry Smith     while (nz--) {
303289bc588SBarry Smith       fm  = n;
304010a20caSHong Zhang       idx = ic[*ajtmp++];
305289bc588SBarry Smith       do {
306289bc588SBarry Smith         m  = fm;
307289bc588SBarry Smith         fm = fill[m];
308289bc588SBarry Smith       } while (fm < idx);
309289bc588SBarry Smith       fill[m]   = idx;
310289bc588SBarry Smith       fill[idx] = fm;
311289bc588SBarry Smith     }
312289bc588SBarry Smith     row = fill[n];
313289bc588SBarry Smith     while (row < i) {
314010a20caSHong Zhang       ajtmp = ajnew + idnew[row] + 1;
3152fb3ed76SBarry Smith       nzbd  = 1 + idnew[row] - ainew[row];
3162fb3ed76SBarry Smith       nz    = im[row] - nzbd;
317289bc588SBarry Smith       fm    = row;
3182fb3ed76SBarry Smith       while (nz-- > 0) {
319010a20caSHong Zhang         idx = *ajtmp++ ;
3202fb3ed76SBarry Smith         nzbd++;
3212fb3ed76SBarry Smith         if (idx == i) im[row] = nzbd;
322289bc588SBarry Smith         do {
323289bc588SBarry Smith           m  = fm;
324289bc588SBarry Smith           fm = fill[m];
325289bc588SBarry Smith         } while (fm < idx);
326289bc588SBarry Smith         if (fm != idx) {
327289bc588SBarry Smith           fill[m]   = idx;
328289bc588SBarry Smith           fill[idx] = fm;
329289bc588SBarry Smith           fm        = idx;
330289bc588SBarry Smith           nnz++;
331289bc588SBarry Smith         }
332289bc588SBarry Smith       }
333289bc588SBarry Smith       row = fill[row];
334289bc588SBarry Smith     }
335289bc588SBarry Smith     /* copy new filled row into permanent storage */
336289bc588SBarry Smith     ainew[i+1] = ainew[i] + nnz;
337496e697eSBarry Smith     if (ainew[i+1] > jmax) {
338ecf371e4SBarry Smith       /* estimate how much additional space we will need */
339ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
340ecf371e4SBarry Smith       /* just double the memory each time */
34138baddfdSBarry Smith       PetscInt maxadd = jmax;
342ecf371e4SBarry Smith       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
343bbb6d6a8SBarry Smith       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
3442fb3ed76SBarry Smith       jmax += maxadd;
345ecf371e4SBarry Smith 
346ecf371e4SBarry Smith       /* allocate a longer ajnew */
34738baddfdSBarry Smith       ierr  = PetscMalloc(jmax*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
34838baddfdSBarry Smith       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i])*sizeof(PetscInt));CHKERRQ(ierr);
349606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
350289bc588SBarry Smith       ajnew = ajtmp;
351418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
352289bc588SBarry Smith     }
353010a20caSHong Zhang     ajtmp = ajnew + ainew[i];
354289bc588SBarry Smith     fm    = fill[n];
355289bc588SBarry Smith     nzi   = 0;
3562fb3ed76SBarry Smith     im[i] = nnz;
357289bc588SBarry Smith     while (nnz--) {
358289bc588SBarry Smith       if (fm < i) nzi++;
359010a20caSHong Zhang       *ajtmp++ = fm ;
360289bc588SBarry Smith       fm       = fill[fm];
361289bc588SBarry Smith     }
362289bc588SBarry Smith     idnew[i] = ainew[i] + nzi;
363289bc588SBarry Smith   }
364464e8d44SSatish Balay   if (ai[n] != 0) {
365329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
366418422e8SSatish Balay     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
367b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
368b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
369b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
37005bf559cSSatish Balay   } else {
371b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
37205bf559cSSatish Balay   }
3732fb3ed76SBarry Smith 
374898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
375898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3761987afe7SBarry Smith 
377606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
378289bc588SBarry Smith 
379289bc588SBarry Smith   /* put together the new matrix */
380f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr);
381f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
382f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr);
383b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
384416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
385606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3867c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
387e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
388606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
389606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
390010a20caSHong Zhang   ierr          = PetscMalloc((ainew[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
391416022c9SBarry Smith   b->j          = ajnew;
392416022c9SBarry Smith   b->i          = ainew;
393416022c9SBarry Smith   b->diag       = idnew;
394416022c9SBarry Smith   b->ilen       = 0;
395416022c9SBarry Smith   b->imax       = 0;
396416022c9SBarry Smith   b->row        = isrow;
397416022c9SBarry Smith   b->col        = iscol;
398c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
399c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40082bf6240SBarry Smith   b->icol       = isicol;
40187828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
402416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4037fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
40438baddfdSBarry Smith   PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));
405010a20caSHong Zhang   b->maxnz = b->nz = ainew[n] ;
4067fd98bd6SLois Curfman McInnes 
407329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
408418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
409ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4103a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4117dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
412703deb49SSatish Balay 
413e93ccc4dSBarry Smith   if (ai[n] != 0) {
414329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
415eea2913aSSatish Balay   } else {
416eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
417eea2913aSSatish Balay   }
4183a40ed3dSBarry Smith   PetscFunctionReturn(0);
419289bc588SBarry Smith }
4200a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
421dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
42241c01911SSatish Balay 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
425af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
426289bc588SBarry Smith {
427416022c9SBarry Smith   Mat            C=*B;
428416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
42982bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4306849ba73SBarry Smith   PetscErrorCode ierr;
43138baddfdSBarry Smith   PetscInt       *r,*ic,i,j,n=A->m,*ai=b->i,*aj=b->j;
43238baddfdSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ics;
433e24cfe64SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj,nshift=0;
43487828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
435*ee45ca4aSHong Zhang   PetscReal      zeropivot,rs,d,shift_nonzero;
4360cf777b8SBarry Smith   PetscReal      row_shift,shift_fraction,shift_amount,shift_lo=0.,shift_hi=1.,shift_top=0.;
437*ee45ca4aSHong Zhang   PetscTruth     lushift,shift_pd;
438289bc588SBarry Smith 
4393a40ed3dSBarry Smith   PetscFunctionBegin;
4400a29876aSHong Zhang   shift_nonzero  = info->shiftnz;
4410a29876aSHong Zhang   shift_pd       = info->shiftpd;
4420a29876aSHong Zhang   zeropivot      = info->zeropivot;
4430a29876aSHong Zhang   shift_fraction = info->shift_fraction;
4440a29876aSHong Zhang 
44578b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44678b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44787828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44887828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
449010a20caSHong Zhang   rtmps = rtmp; ics = ic;
450289bc588SBarry Smith 
4516cc28720Svictorle   if (!a->diag) {
4520cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4530cf777b8SBarry Smith   }
454e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
455e24cfe64SHong Zhang   if (shift_pd && shift_nonzero) shift_nonzero = 0.0;
456e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45738baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4580cf777b8SBarry Smith     shift_top = 0;
4596cc28720Svictorle     for (i=0; i<n; i++) {
460010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4616cc28720Svictorle       /* calculate amt of shift needed for this row */
4626c037d1bSvictorle       if (d<=0) {
4633aeef889SHong Zhang         row_shift = 0;
4643aeef889SHong Zhang       } else {
4653aeef889SHong Zhang         row_shift = -2*d;
4663aeef889SHong Zhang       }
467010a20caSHong Zhang       v  = a->a+aai[i];
468e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
469e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4706cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4716cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4726cc28720Svictorle     }
4736cc28720Svictorle   }
4746cc28720Svictorle 
4750a29876aSHong Zhang   /* shift_fraction = b->lu_shift_fraction; */
4760697b6caSHong Zhang   shift_amount   = 0;
477085256b3SBarry Smith   do {
4786cc28720Svictorle     lushift = PETSC_FALSE;
479289bc588SBarry Smith     for (i=0; i<n; i++){
480289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
481010a20caSHong Zhang       ajtmp = aj + ai[i];
48248e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
483289bc588SBarry Smith 
484289bc588SBarry Smith       /* load in initial (unfactored row) */
485416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
486010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
487010a20caSHong Zhang       v        = a->a + a->i[r[i]];
488085256b3SBarry Smith       for (j=0; j<nz; j++) {
489085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
490085256b3SBarry Smith       }
491e24cfe64SHong Zhang       rtmp[ics[r[i]]] += shift_amount; /* shift the diagonal of the matrix */
492289bc588SBarry Smith 
493010a20caSHong Zhang       row = *ajtmp++;
494f2582111SSatish Balay       while  (row < i) {
4958c37ef55SBarry Smith         pc = rtmp + row;
4968c37ef55SBarry Smith         if (*pc != 0.0) {
497010a20caSHong Zhang           pv         = b->a + diag_offset[row];
498010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
49935aab85fSBarry Smith           multiplier = *pc / *pv++;
5008c37ef55SBarry Smith           *pc        = multiplier;
501f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
502f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
503b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
504289bc588SBarry Smith         }
505010a20caSHong Zhang         row = *ajtmp++;
506289bc588SBarry Smith       }
507416022c9SBarry Smith       /* finished row so stick it into b->a */
508010a20caSHong Zhang       pv   = b->a + ai[i] ;
509010a20caSHong Zhang       pj   = b->j + ai[i] ;
5108c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
51135aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
512d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
513d3d32019SBarry Smith       rs   = 0.0;
514d3d32019SBarry Smith       for (j=0; j<nz; j++) {
515d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
516d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
517d3d32019SBarry Smith       }
518e24cfe64SHong Zhang       /* shift the diagonals when zero pivot is detected */
5196cc28720Svictorle #define MAX_NSHIFT 5
520e24cfe64SHong Zhang       if (PetscAbsScalar(pv[diag]) <= zeropivot*rs && shift_nonzero){
521e24cfe64SHong Zhang         /* force |diag(*B)| > zeropivot*rs */
522e24cfe64SHong Zhang         if (!nshift){
523e24cfe64SHong Zhang           shift_amount = shift_nonzero;
524e24cfe64SHong Zhang         } else {
525e24cfe64SHong Zhang           shift_amount *= 2.0;
526e24cfe64SHong Zhang         }
527e24cfe64SHong Zhang         lushift = PETSC_TRUE;
528e24cfe64SHong Zhang         nshift++;
529e24cfe64SHong Zhang         break;
530e24cfe64SHong Zhang       } else if (PetscRealPart(pv[diag]) <= zeropivot*rs && shift_pd){
531e24cfe64SHong Zhang         /* force *B to be diagonally dominant */
5326cc28720Svictorle 	if (nshift>MAX_NSHIFT) {
533e005ede5SBarry Smith 	  SETERRQ(PETSC_ERR_CONV_FAILED,"Unable to determine shift to enforce positive definite preconditioner");
5346cc28720Svictorle 	} else if (nshift==MAX_NSHIFT) {
5356cc28720Svictorle 	  shift_fraction = shift_hi;
5366c037d1bSvictorle 	  lushift        = PETSC_FALSE;
5376cc28720Svictorle 	} else {
5386cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5396c037d1bSvictorle 	  lushift  = PETSC_TRUE;
5406cc28720Svictorle 	}
5416cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5420cf777b8SBarry Smith         nshift++;
5430cf777b8SBarry Smith         break;
544e24cfe64SHong Zhang       } else if (PetscAbsScalar(pv[diag]) <= zeropivot*rs){
54577431f27SBarry Smith         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
546d708749eSBarry Smith       }
54735aab85fSBarry Smith     }
548e24cfe64SHong Zhang     if (shift_pd && !lushift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5496cc28720Svictorle       /*
5506c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5516cc28720Svictorle        * then try lower shift
5526cc28720Svictorle        */
5530cf777b8SBarry Smith       shift_hi       = shift_fraction;
5540cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5556cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5560cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5570cf777b8SBarry Smith       nshift++;
5586cc28720Svictorle     }
559e24cfe64SHong Zhang   } while (lushift);
5600f11f9aeSSatish Balay 
56135aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
56235aab85fSBarry Smith   for (i=0; i<n; i++) {
563010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
56435aab85fSBarry Smith   }
56535aab85fSBarry Smith 
566606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
56778b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
56878b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
569416022c9SBarry Smith   C->factor = FACTOR_LU;
5707dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
571c456f294SBarry Smith   C->assembled = PETSC_TRUE;
572b0a32e0cSBarry Smith   PetscLogFlops(C->n);
5730697b6caSHong Zhang   if (nshift){
5740697b6caSHong Zhang     if (shift_nonzero) {
5750697b6caSHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",nshift,shift_amount);
5760697b6caSHong Zhang     } else if (shift_pd) {
5776cc28720Svictorle       b->lu_shift_fraction = shift_fraction;
57877431f27SBarry Smith       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %D\n",shift_fraction,shift_top,nshift);
5796cc28720Svictorle     }
5800697b6caSHong Zhang   }
5813a40ed3dSBarry Smith   PetscFunctionReturn(0);
582289bc588SBarry Smith }
5830889b5dcSvictorle 
5840889b5dcSvictorle #undef __FUNCT__
5850889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
586dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5870889b5dcSvictorle {
5880889b5dcSvictorle   PetscFunctionBegin;
5890889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5900889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5910889b5dcSvictorle   PetscFunctionReturn(0);
5920889b5dcSvictorle }
5930889b5dcSvictorle 
5940889b5dcSvictorle 
5950a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5964a2ae208SSatish Balay #undef __FUNCT__
5974a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
598dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
599da3a660dSBarry Smith {
600dfbe8321SBarry Smith   PetscErrorCode ierr;
601416022c9SBarry Smith   Mat            C;
602416022c9SBarry Smith 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
6049e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
605af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
606273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
607440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
6083a40ed3dSBarry Smith   PetscFunctionReturn(0);
609da3a660dSBarry Smith }
6100a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
6114a2ae208SSatish Balay #undef __FUNCT__
6124a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
613dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6148c37ef55SBarry Smith {
615416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
616416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6176849ba73SBarry Smith   PetscErrorCode ierr;
61838baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
61938baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
62087828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6218c37ef55SBarry Smith 
6223a40ed3dSBarry Smith   PetscFunctionBegin;
6233a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
62491bf9adeSBarry Smith 
6251ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
627416022c9SBarry Smith   tmp  = a->solve_work;
6288c37ef55SBarry Smith 
6293b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6303b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6318c37ef55SBarry Smith 
6328c37ef55SBarry Smith   /* forward solve the lower triangular */
6338c37ef55SBarry Smith   tmp[0] = b[*r++];
634010a20caSHong Zhang   tmps   = tmp;
6358c37ef55SBarry Smith   for (i=1; i<n; i++) {
636010a20caSHong Zhang     v   = aa + ai[i] ;
637010a20caSHong Zhang     vi  = aj + ai[i] ;
63853e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
63953e7425aSBarry Smith     sum = b[*r++];
640ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6418c37ef55SBarry Smith     tmp[i] = sum;
6428c37ef55SBarry Smith   }
6438c37ef55SBarry Smith 
6448c37ef55SBarry Smith   /* backward solve the upper triangular */
6458c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
646010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
647010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
648416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6498c37ef55SBarry Smith     sum = tmp[i];
650ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
651010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6528c37ef55SBarry Smith   }
6538c37ef55SBarry Smith 
6543b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6553b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6561ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6571ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
658b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6593a40ed3dSBarry Smith   PetscFunctionReturn(0);
6608c37ef55SBarry Smith }
661026e39d0SSatish Balay 
662930ae53cSBarry Smith /* ----------------------------------------------------------- */
6634a2ae208SSatish Balay #undef __FUNCT__
6644a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
665dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
666930ae53cSBarry Smith {
667930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6686849ba73SBarry Smith   PetscErrorCode ierr;
66938baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
670362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
671aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
67238baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
673362ced78SSatish Balay   PetscScalar    *v,sum;
674d85afc42SSatish Balay #endif
675930ae53cSBarry Smith 
6763a40ed3dSBarry Smith   PetscFunctionBegin;
6773a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
678930ae53cSBarry Smith 
6791ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6801ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
681930ae53cSBarry Smith 
682aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6831c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6841c4feecaSSatish Balay #else
685930ae53cSBarry Smith   /* forward solve the lower triangular */
686930ae53cSBarry Smith   x[0] = b[0];
687930ae53cSBarry Smith   for (i=1; i<n; i++) {
688930ae53cSBarry Smith     ai_i = ai[i];
689930ae53cSBarry Smith     v    = aa + ai_i;
690930ae53cSBarry Smith     vi   = aj + ai_i;
691930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
692930ae53cSBarry Smith     sum  = b[i];
693930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
694930ae53cSBarry Smith     x[i] = sum;
695930ae53cSBarry Smith   }
696930ae53cSBarry Smith 
697930ae53cSBarry Smith   /* backward solve the upper triangular */
698930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
699930ae53cSBarry Smith     adiag_i = adiag[i];
700d708749eSBarry Smith     v       = aa + adiag_i + 1;
701d708749eSBarry Smith     vi      = aj + adiag_i + 1;
702930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
703930ae53cSBarry Smith     sum     = x[i];
704930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
705930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
706930ae53cSBarry Smith   }
7071c4feecaSSatish Balay #endif
708b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
7091ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7101ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7113a40ed3dSBarry Smith   PetscFunctionReturn(0);
712930ae53cSBarry Smith }
713930ae53cSBarry Smith 
7144a2ae208SSatish Balay #undef __FUNCT__
7154a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
716dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
717da3a660dSBarry Smith {
718416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
719416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7206849ba73SBarry Smith   PetscErrorCode ierr;
72138baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
72238baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
72387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
724da3a660dSBarry Smith 
7253a40ed3dSBarry Smith   PetscFunctionBegin;
72678b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
727da3a660dSBarry Smith 
7281ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7291ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
730416022c9SBarry Smith   tmp  = a->solve_work;
731da3a660dSBarry Smith 
7323b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7333b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
734da3a660dSBarry Smith 
735da3a660dSBarry Smith   /* forward solve the lower triangular */
736da3a660dSBarry Smith   tmp[0] = b[*r++];
737da3a660dSBarry Smith   for (i=1; i<n; i++) {
738010a20caSHong Zhang     v   = aa + ai[i] ;
739010a20caSHong Zhang     vi  = aj + ai[i] ;
740416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
741da3a660dSBarry Smith     sum = b[*r++];
742010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
743da3a660dSBarry Smith     tmp[i] = sum;
744da3a660dSBarry Smith   }
745da3a660dSBarry Smith 
746da3a660dSBarry Smith   /* backward solve the upper triangular */
747da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
748010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
749010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
750416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
751da3a660dSBarry Smith     sum = tmp[i];
752010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
753010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
754da3a660dSBarry Smith     x[*c--] += tmp[i];
755da3a660dSBarry Smith   }
756da3a660dSBarry Smith 
7573b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7583b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7591ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7601ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
761b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
762898183ecSLois Curfman McInnes 
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
764da3a660dSBarry Smith }
765da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7664a2ae208SSatish Balay #undef __FUNCT__
7674a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
768dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
769da3a660dSBarry Smith {
770416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7712235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7726849ba73SBarry Smith   PetscErrorCode ierr;
77338baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
77438baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
77587828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
776da3a660dSBarry Smith 
7773a40ed3dSBarry Smith   PetscFunctionBegin;
7781ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7791ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
780416022c9SBarry Smith   tmp  = a->solve_work;
781da3a660dSBarry Smith 
7822235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7832235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
784da3a660dSBarry Smith 
785da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7862235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
787da3a660dSBarry Smith 
788da3a660dSBarry Smith   /* forward solve the U^T */
789da3a660dSBarry Smith   for (i=0; i<n; i++) {
790010a20caSHong Zhang     v   = aa + diag[i] ;
791010a20caSHong Zhang     vi  = aj + diag[i] + 1;
792f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
793c38d4ed2SBarry Smith     s1  = tmp[i];
794ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
795da3a660dSBarry Smith     while (nz--) {
796010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
797da3a660dSBarry Smith     }
798c38d4ed2SBarry Smith     tmp[i] = s1;
799da3a660dSBarry Smith   }
800da3a660dSBarry Smith 
801da3a660dSBarry Smith   /* backward solve the L^T */
802da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
803010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
804010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
805f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
806f1af5d2fSBarry Smith     s1  = tmp[i];
807da3a660dSBarry Smith     while (nz--) {
808010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
809da3a660dSBarry Smith     }
810da3a660dSBarry Smith   }
811da3a660dSBarry Smith 
812da3a660dSBarry Smith   /* copy tmp into x according to permutation */
813da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
814da3a660dSBarry Smith 
8152235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8162235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8171ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8181ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
819da3a660dSBarry Smith 
820b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8213a40ed3dSBarry Smith   PetscFunctionReturn(0);
822da3a660dSBarry Smith }
823da3a660dSBarry Smith 
8244a2ae208SSatish Balay #undef __FUNCT__
8254a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
826dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
827da3a660dSBarry Smith {
828416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8292235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8306849ba73SBarry Smith   PetscErrorCode ierr;
83138baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
83238baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
83387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8346abc6512SBarry Smith 
8353a40ed3dSBarry Smith   PetscFunctionBegin;
83623bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8376abc6512SBarry Smith 
8381ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
840416022c9SBarry Smith   tmp = a->solve_work;
8416abc6512SBarry Smith 
8422235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8432235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8446abc6512SBarry Smith 
8456abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8462235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8476abc6512SBarry Smith 
8486abc6512SBarry Smith   /* forward solve the U^T */
8496abc6512SBarry Smith   for (i=0; i<n; i++) {
850010a20caSHong Zhang     v   = aa + diag[i] ;
851010a20caSHong Zhang     vi  = aj + diag[i] + 1;
852f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8536abc6512SBarry Smith     tmp[i] *= *v++;
8546abc6512SBarry Smith     while (nz--) {
855010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8566abc6512SBarry Smith     }
8576abc6512SBarry Smith   }
8586abc6512SBarry Smith 
8596abc6512SBarry Smith   /* backward solve the L^T */
8606abc6512SBarry Smith   for (i=n-1; i>=0; i--){
861010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
862010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
863f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8646abc6512SBarry Smith     while (nz--) {
865010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8666abc6512SBarry Smith     }
8676abc6512SBarry Smith   }
8686abc6512SBarry Smith 
8696abc6512SBarry Smith   /* copy tmp into x according to permutation */
8706abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8716abc6512SBarry Smith 
8722235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8732235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8741ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8751ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8766abc6512SBarry Smith 
877b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8783a40ed3dSBarry Smith   PetscFunctionReturn(0);
879da3a660dSBarry Smith }
8809e25ed09SBarry Smith /* ----------------------------------------------------------------*/
881dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
88208480c60SBarry Smith 
8834a2ae208SSatish Balay #undef __FUNCT__
8844a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
885dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8869e25ed09SBarry Smith {
887416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8889e25ed09SBarry Smith   IS             isicol;
8896849ba73SBarry Smith   PetscErrorCode ierr;
89038baddfdSBarry Smith   PetscInt       *r,*ic,prow,n = A->m,*ai = a->i,*aj = a->j;
89138baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
892418422e8SSatish Balay   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,len, reallocs = 0,dcount = 0;
89338baddfdSBarry Smith   PetscInt       incrlev,nnz,i,levels,diagonal_fill;
89477c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
895329f5518SBarry Smith   PetscReal      f;
8969e25ed09SBarry Smith 
8973a40ed3dSBarry Smith   PetscFunctionBegin;
8986cf997b0SBarry Smith   f             = info->fill;
89938baddfdSBarry Smith   levels        = (PetscInt)info->levels;
90038baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9014c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
90282bf6240SBarry Smith 
90308480c60SBarry Smith   /* special case that simply copies fill pattern */
904be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
905be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
90686bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9072e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
90808480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
90908480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
91008480c60SBarry Smith     if (!b->diag) {
9117c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
91208480c60SBarry Smith     }
9137c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
91408480c60SBarry Smith     b->row              = isrow;
91508480c60SBarry Smith     b->col              = iscol;
91682bf6240SBarry Smith     b->icol             = isicol;
91787828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
918f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
919c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
920c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9213a40ed3dSBarry Smith     PetscFunctionReturn(0);
92208480c60SBarry Smith   }
92308480c60SBarry Smith 
924898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
925898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9269e25ed09SBarry Smith 
9279e25ed09SBarry Smith   /* get new row pointers */
92838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
929010a20caSHong Zhang   ainew[0] = 0;
9309e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
93138baddfdSBarry Smith   jmax = (PetscInt)(f*(ai[n]+1));
93238baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
9339e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
93438baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajfill);CHKERRQ(ierr);
9359e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
93638baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
93756beaf04SBarry Smith   /* im is level for each filled value */
93838baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
93956beaf04SBarry Smith   /* dloc is location of diagonal in factor */
94038baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&dloc);CHKERRQ(ierr);
94156beaf04SBarry Smith   dloc[0]  = 0;
94256beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9436cf997b0SBarry Smith 
9446cf997b0SBarry Smith     /* copy current row into linked list */
94556beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
94629bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
947010a20caSHong Zhang     xi      = aj + ai[r[prow]];
9489e25ed09SBarry Smith     fill[n] = n;
949435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9509e25ed09SBarry Smith     while (nz--) {
9519e25ed09SBarry Smith       fm  = n;
952010a20caSHong Zhang       idx = ic[*xi++ ];
9539e25ed09SBarry Smith       do {
9549e25ed09SBarry Smith         m  = fm;
9559e25ed09SBarry Smith         fm = fill[m];
9569e25ed09SBarry Smith       } while (fm < idx);
9579e25ed09SBarry Smith       fill[m]   = idx;
9589e25ed09SBarry Smith       fill[idx] = fm;
95956beaf04SBarry Smith       im[idx]   = 0;
9609e25ed09SBarry Smith     }
9616cf997b0SBarry Smith 
9626cf997b0SBarry Smith     /* make sure diagonal entry is included */
963435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9646cf997b0SBarry Smith       fm = n;
965435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
966435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
967435faa5fSBarry Smith       fill[fm]   = prow;
9686cf997b0SBarry Smith       im[prow]   = 0;
9696cf997b0SBarry Smith       nzf++;
970335d9088SBarry Smith       dcount++;
9716cf997b0SBarry Smith     }
9726cf997b0SBarry Smith 
97356beaf04SBarry Smith     nzi = 0;
9749e25ed09SBarry Smith     row = fill[n];
97556beaf04SBarry Smith     while (row < prow) {
97656beaf04SBarry Smith       incrlev = im[row] + 1;
97756beaf04SBarry Smith       nz      = dloc[row];
978010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
979010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
980416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
98156beaf04SBarry Smith       fm      = row;
98256beaf04SBarry Smith       while (nnz-- > 0) {
983010a20caSHong Zhang         idx = *xi++ ;
98456beaf04SBarry Smith         if (*flev + incrlev > levels) {
98556beaf04SBarry Smith           flev++;
98656beaf04SBarry Smith           continue;
98756beaf04SBarry Smith         }
98856beaf04SBarry Smith         do {
9899e25ed09SBarry Smith           m  = fm;
9909e25ed09SBarry Smith           fm = fill[m];
99156beaf04SBarry Smith         } while (fm < idx);
9929e25ed09SBarry Smith         if (fm != idx) {
99356beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9949e25ed09SBarry Smith           fill[m]   = idx;
9959e25ed09SBarry Smith           fill[idx] = fm;
9969e25ed09SBarry Smith           fm        = idx;
99756beaf04SBarry Smith           nzf++;
998ecf371e4SBarry Smith         } else {
99956beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
10009e25ed09SBarry Smith         }
100156beaf04SBarry Smith         flev++;
10029e25ed09SBarry Smith       }
10039e25ed09SBarry Smith       row = fill[row];
100456beaf04SBarry Smith       nzi++;
10059e25ed09SBarry Smith     }
10069e25ed09SBarry Smith     /* copy new filled row into permanent storage */
100756beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
1008010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
1009ecf371e4SBarry Smith 
1010ecf371e4SBarry Smith       /* estimate how much additional space we will need */
1011ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
1012ecf371e4SBarry Smith       /* just double the memory each time */
101338baddfdSBarry Smith       /*  maxadd = (PetscInt)((f*(ai[n]+!shift)*(n-prow+5))/n); */
101438baddfdSBarry Smith       PetscInt maxadd = jmax;
101556beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1016bbb6d6a8SBarry Smith       jmax += maxadd;
1017ecf371e4SBarry Smith 
1018ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
101938baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
102038baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1021606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
102256beaf04SBarry Smith       ajnew  = xi;
102338baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
102438baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1025606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
102656beaf04SBarry Smith       ajfill = xi;
1027418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
10289e25ed09SBarry Smith     }
1029010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1030010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
103156beaf04SBarry Smith     dloc[prow]  = nzi;
10329e25ed09SBarry Smith     fm          = fill[n];
103356beaf04SBarry Smith     while (nzf--) {
1034010a20caSHong Zhang       *xi++   = fm ;
103556beaf04SBarry Smith       *flev++ = im[fm];
10369e25ed09SBarry Smith       fm      = fill[fm];
10379e25ed09SBarry Smith     }
10386cf997b0SBarry Smith     /* make sure row has diagonal entry */
1039010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
104077431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10416cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10426cf997b0SBarry Smith     }
10439e25ed09SBarry Smith   }
1044606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1045898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1046898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1047606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1048606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10499e25ed09SBarry Smith 
1050f64065bbSBarry Smith   {
1051329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1052418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1053c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1054c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1055b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1056335d9088SBarry Smith     if (diagonal_fill) {
105777431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1058335d9088SBarry Smith     }
1059f64065bbSBarry Smith   }
1060bbb6d6a8SBarry Smith 
10619e25ed09SBarry Smith   /* put together the new matrix */
1062f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1063f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1064f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
1065b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1066416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1067606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10687c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1069010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10709e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1071606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1072606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1073b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1074416022c9SBarry Smith   b->j          = ajnew;
1075416022c9SBarry Smith   b->i          = ainew;
107656beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1077416022c9SBarry Smith   b->diag       = dloc;
1078416022c9SBarry Smith   b->ilen       = 0;
1079416022c9SBarry Smith   b->imax       = 0;
1080416022c9SBarry Smith   b->row        = isrow;
1081416022c9SBarry Smith   b->col        = iscol;
1082c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1083c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
108482bf6240SBarry Smith   b->icol       = isicol;
108587828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1086416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
108756beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
108838baddfdSBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
1089010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
10909e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10913a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10923a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1093418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1094ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1095329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10963a40ed3dSBarry Smith   PetscFunctionReturn(0);
10979e25ed09SBarry Smith }
1098d5d45c9bSBarry Smith 
10993c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1100a6175056SHong Zhang #undef __FUNCT__
1101a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1102af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,MatFactorInfo *info,Mat *B)
1103a6175056SHong Zhang {
11042ed133dbSHong Zhang   Mat            C = *B;
11050968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11062ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11072ed133dbSHong Zhang   IS             ip=b->row;
11082ed133dbSHong Zhang   PetscErrorCode ierr;
11092ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
11102ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11112ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1112622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1113*ee45ca4aSHong Zhang   PetscReal      zeropivot,shift_amount,rs,shift_nonzero;
1114*ee45ca4aSHong Zhang   PetscTruth     chshift,shift_pd;
111539e3d630SHong Zhang   PetscInt       nshift=0;
1116d5d45c9bSBarry Smith 
1117a6175056SHong Zhang   PetscFunctionBegin;
1118*ee45ca4aSHong Zhang   shift_nonzero  = info->shiftnz;
1119*ee45ca4aSHong Zhang   shift_pd       = info->shiftpd;
1120*ee45ca4aSHong Zhang   zeropivot      = info->zeropivot;
1121*ee45ca4aSHong Zhang 
11222ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11232ed133dbSHong Zhang 
11242ed133dbSHong Zhang   /* initialization */
11252ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11262ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11272ed133dbSHong Zhang   jl   = il + mbs;
11282ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11292ed133dbSHong Zhang 
11302ed133dbSHong Zhang   shift_amount = 0;
11312ed133dbSHong Zhang   do {
11322ed133dbSHong Zhang     chshift = PETSC_FALSE;
11332ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11342ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1135a6175056SHong Zhang     }
11362ed133dbSHong Zhang 
11372ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1138c05c3958SHong Zhang       bval = ba + bi[k];
11392ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11402ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11412ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11422ed133dbSHong Zhang         col = rip[aj[j]];
11432ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11442ed133dbSHong Zhang           rtmp[col] = aa[j];
1145c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11462ed133dbSHong Zhang         }
11472ed133dbSHong Zhang       }
114839e3d630SHong Zhang       /* shift the diagonal of the matrix */
114939e3d630SHong Zhang       if (nshift) rtmp[k] += shift_amount;
11502ed133dbSHong Zhang 
11512ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11522ed133dbSHong Zhang       dk = rtmp[k];
11532ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11542ed133dbSHong Zhang 
11552ed133dbSHong Zhang       while (i < k){
11562ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11572ed133dbSHong Zhang 
11582ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11592ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11602ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11612ed133dbSHong Zhang         dk += uikdi*ba[ili];
11622ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11632ed133dbSHong Zhang 
11642ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11652ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11662ed133dbSHong Zhang         if (jmin < jmax){
11672ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11682ed133dbSHong Zhang           /* update il and jl for row i */
11692ed133dbSHong Zhang           il[i] = jmin;
11702ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11712ed133dbSHong Zhang         }
11722ed133dbSHong Zhang         i = nexti;
11732ed133dbSHong Zhang       }
11742ed133dbSHong Zhang 
11750697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
11760697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
11770697b6caSHong Zhang       rs   = 0.0;
11783c9e8598SHong Zhang       jmin = bi[k]+1;
11793c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
11803c9e8598SHong Zhang       if (nz){
11813c9e8598SHong Zhang         bcol = bj + jmin;
11823c9e8598SHong Zhang         while (nz--){
11833c9e8598SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol++]);
11843c9e8598SHong Zhang         }
11853c9e8598SHong Zhang       }
11860697b6caSHong Zhang       if (PetscAbsScalar(dk) <= zeropivot*rs && shift_nonzero){
11870697b6caSHong Zhang         /* force |diag(*B)| > zeropivot*rs */
11880697b6caSHong Zhang         if (!nshift){
11890697b6caSHong Zhang           shift_amount = shift_nonzero;
11900697b6caSHong Zhang         } else {
11910697b6caSHong Zhang           shift_amount *= 2.0;
11920697b6caSHong Zhang         }
11930697b6caSHong Zhang         chshift = PETSC_TRUE;
11940697b6caSHong Zhang         nshift++;
11950697b6caSHong Zhang         break;
11960697b6caSHong Zhang       } else if (PetscRealPart(dk) <= zeropivot*rs && shift_pd){
11970697b6caSHong Zhang         /* calculate a shift that would make this row diagonally dominant */
11980697b6caSHong Zhang 	shift_amount = PetscMax(rs+PetscAbs(PetscRealPart(dk)),1.1*shift_amount);
11993c9e8598SHong Zhang 	chshift      = PETSC_TRUE;
12003c9e8598SHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
12013c9e8598SHong Zhang 	   we increase the shift until it converges. There is no guarantee that
12023c9e8598SHong Zhang 	   this algorithm converges faster or slower, or is better or worse
12033c9e8598SHong Zhang 	   than the ILU algorithm. */
12043c9e8598SHong Zhang 	nshift++;
12053c9e8598SHong Zhang 	break;
12060697b6caSHong Zhang       } else if (PetscAbsScalar(dk) <= zeropivot*rs){
12070697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
12083c9e8598SHong Zhang       }
12093c9e8598SHong Zhang 
12103c9e8598SHong Zhang       /* copy data into U(k,:) */
121139e3d630SHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
121239e3d630SHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
121339e3d630SHong Zhang       if (jmin < jmax) {
121439e3d630SHong Zhang         for (j=jmin; j<jmax; j++){
121539e3d630SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12163c9e8598SHong Zhang         }
121739e3d630SHong Zhang         /* add the k-th row into il and jl */
12183c9e8598SHong Zhang         il[k] = jmin;
12193c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12203c9e8598SHong Zhang       }
12213c9e8598SHong Zhang     }
12220697b6caSHong Zhang   } while (chshift);
12233c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12243c9e8598SHong Zhang 
122539e3d630SHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12263c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
12273c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
12283c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
12293c9e8598SHong Zhang   PetscLogFlops(C->m);
12303c9e8598SHong Zhang   if (nshift){
12310697b6caSHong Zhang     if (shift_nonzero) {
123239e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_nonzero tries %D, shift_amount %g\n",nshift,shift_amount);
12330697b6caSHong Zhang     } else if (shift_pd) {
123439e3d630SHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ: number of shift_pd tries %D, shift_amount %g\n",nshift,shift_amount);
12350697b6caSHong Zhang     }
12363c9e8598SHong Zhang   }
12373c9e8598SHong Zhang   PetscFunctionReturn(0);
12383c9e8598SHong Zhang }
1239a6175056SHong Zhang 
12407a48dd6fSHong Zhang #include "petscbt.h"
12417a48dd6fSHong Zhang #include "src/mat/utils/freespace.h"
1242a6175056SHong Zhang #undef __FUNCT__
1243a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1244dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1245a6175056SHong Zhang {
12460968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1247ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1248ed59401aSHong Zhang   Mat            B;
1249dfbe8321SBarry Smith   PetscErrorCode ierr;
1250653a6975SHong Zhang   PetscTruth     perm_identity;
1251622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1252ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1253622e2028SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
1254ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
12557a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
12567a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
12577a48dd6fSHong Zhang   PetscBT        lnkbt;
1258a6175056SHong Zhang 
1259a6175056SHong Zhang   PetscFunctionBegin;
1260653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
12617a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
12627a48dd6fSHong Zhang 
126339e3d630SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
126439e3d630SHong Zhang   ui[0] = 0;
126539e3d630SHong Zhang 
1266622e2028SHong Zhang   /* special case that simply copies fill pattern */
1267622e2028SHong Zhang   if (!levels && perm_identity) {
1268622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1269ed59401aSHong Zhang     for (i=0; i<am; i++) {
127039e3d630SHong Zhang       ui[i+1] = ui[i] + ai[i+1] - a->diag[i];
1271ed59401aSHong Zhang     }
127239e3d630SHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
127339e3d630SHong Zhang     cols = uj;
1274ed59401aSHong Zhang     for (i=0; i<am; i++) {
1275ed59401aSHong Zhang       aj    = a->j + a->diag[i];
127639e3d630SHong Zhang       ncols = ui[i+1] - ui[i];
127739e3d630SHong Zhang       for (j=0; j<ncols; j++) *cols++ = *aj++;
1278ed59401aSHong Zhang     }
127939e3d630SHong Zhang   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
12807a48dd6fSHong Zhang     /* initialization */
1281622e2028SHong Zhang     ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
12827a48dd6fSHong Zhang 
12837a48dd6fSHong Zhang     /* jl: linked list for storing indices of the pivot rows
12847a48dd6fSHong Zhang        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
12857a48dd6fSHong Zhang     ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
12867a48dd6fSHong Zhang     il         = jl + am;
12877a48dd6fSHong Zhang     uj_ptr     = (PetscInt**)(il + am);
12887a48dd6fSHong Zhang     uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
12897a48dd6fSHong Zhang     for (i=0; i<am; i++){
12907a48dd6fSHong Zhang       jl[i] = am; il[i] = 0;
12917a48dd6fSHong Zhang     }
12927a48dd6fSHong Zhang 
12937a48dd6fSHong Zhang     /* create and initialize a linked list for storing column indices of the active row k */
12947a48dd6fSHong Zhang     nlnk = am + 1;
12952ed133dbSHong Zhang     ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
12967a48dd6fSHong Zhang 
12977a48dd6fSHong Zhang     /* initial FreeSpace size is fill*(ai[am]+1) */
12987a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
12997a48dd6fSHong Zhang     current_space = free_space;
13007a48dd6fSHong Zhang     ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
13017a48dd6fSHong Zhang     current_space_lvl = free_space_lvl;
13027a48dd6fSHong Zhang 
13037a48dd6fSHong Zhang     for (k=0; k<am; k++){  /* for each active row k */
13047a48dd6fSHong Zhang       /* initialize lnk by the column indices of row rip[k] of A */
13057a48dd6fSHong Zhang       nzk   = 0;
1306622e2028SHong Zhang       ncols = ai[rip[k]+1] - ai[rip[k]];
1307622e2028SHong Zhang       ncols_upper = 0;
1308622e2028SHong Zhang       cols        = cols_lvl + am;
1309622e2028SHong Zhang       for (j=0; j<ncols; j++){
1310622e2028SHong Zhang         i = rip[*(aj + ai[rip[k]] + j)];
1311622e2028SHong Zhang         if (i >= k){ /* only take upper triangular entry */
1312622e2028SHong Zhang           cols[ncols_upper] = i;
1313622e2028SHong Zhang           cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
1314622e2028SHong Zhang           ncols_upper++;
1315622e2028SHong Zhang         }
1316622e2028SHong Zhang       }
1317622e2028SHong Zhang       ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,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) */
1332ed59401aSHong Zhang         for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
13332ed133dbSHong Zhang         ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
13347a48dd6fSHong Zhang         nzk += nlnk;
13357a48dd6fSHong Zhang 
13367a48dd6fSHong Zhang         /* update il and jl for prow */
13377a48dd6fSHong Zhang         if (jmin < jmax){
13387a48dd6fSHong Zhang           il[prow] = jmin;
13397a48dd6fSHong Zhang           j = *cols; jl[prow] = jl[j]; jl[j] = prow;
13407a48dd6fSHong Zhang         }
13417a48dd6fSHong Zhang         prow = nextprow;
13427a48dd6fSHong Zhang       }
13437a48dd6fSHong Zhang 
13447a48dd6fSHong Zhang       /* if free space is not available, make more free space */
13457a48dd6fSHong Zhang       if (current_space->local_remaining<nzk) {
13467a48dd6fSHong Zhang         i = am - k + 1; /* num of unfactored rows */
13477a48dd6fSHong Zhang         i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
13487a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
13497a48dd6fSHong Zhang         ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
13507a48dd6fSHong Zhang         reallocs++;
13517a48dd6fSHong Zhang       }
13527a48dd6fSHong Zhang 
13537a48dd6fSHong Zhang       /* copy data into free_space and free_space_lvl, then initialize lnk */
13542ed133dbSHong Zhang       ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
13557a48dd6fSHong Zhang 
13567a48dd6fSHong Zhang       /* add the k-th row into il and jl */
135739e3d630SHong Zhang       if (nzk > 1){
13587a48dd6fSHong Zhang         i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
13597a48dd6fSHong Zhang         jl[k] = jl[i]; jl[i] = k;
13607a48dd6fSHong Zhang         il[k] = ui[k] + 1;
13617a48dd6fSHong Zhang       }
13627a48dd6fSHong Zhang       uj_ptr[k]     = current_space->array;
13637a48dd6fSHong Zhang       uj_lvl_ptr[k] = current_space_lvl->array;
13647a48dd6fSHong Zhang 
13657a48dd6fSHong Zhang       current_space->array           += nzk;
13667a48dd6fSHong Zhang       current_space->local_used      += nzk;
13677a48dd6fSHong Zhang       current_space->local_remaining -= nzk;
13687a48dd6fSHong Zhang 
13697a48dd6fSHong Zhang       current_space_lvl->array           += nzk;
13707a48dd6fSHong Zhang       current_space_lvl->local_used      += nzk;
13717a48dd6fSHong Zhang       current_space_lvl->local_remaining -= nzk;
13727a48dd6fSHong Zhang 
13737a48dd6fSHong Zhang       ui[k+1] = ui[k] + nzk;
13747a48dd6fSHong Zhang     }
13757a48dd6fSHong Zhang 
13767a48dd6fSHong Zhang     if (ai[am] != 0) {
137739e3d630SHong Zhang       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
13787a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
13797a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
13807a48dd6fSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
13817a48dd6fSHong Zhang     } else {
1382ed59401aSHong Zhang       PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
13837a48dd6fSHong Zhang     }
13847a48dd6fSHong Zhang 
13857a48dd6fSHong Zhang     ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
13867a48dd6fSHong Zhang     ierr = PetscFree(jl);CHKERRQ(ierr);
13877a48dd6fSHong Zhang     ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
13887a48dd6fSHong Zhang 
13897a48dd6fSHong Zhang     /* destroy list of free space and other temporary array(s) */
13907a48dd6fSHong Zhang     ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
13917a48dd6fSHong Zhang     ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
13922ed133dbSHong Zhang     ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
13937a48dd6fSHong Zhang     ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
13947a48dd6fSHong Zhang 
139539e3d630SHong Zhang   } /* end of case: levels>0 || (levels=0 && !perm_identity) */
139639e3d630SHong Zhang 
13977a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
13987a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1399ed59401aSHong Zhang   B = *fact;
1400ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1401ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
14027a48dd6fSHong Zhang 
1403ed59401aSHong Zhang   b    = (Mat_SeqSBAIJ*)B->data;
14047a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
14057a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
14067a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
14077a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
14087a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
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);
1421ed59401aSHong Zhang   PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
14227a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
14237a48dd6fSHong Zhang 
1424ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1425ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1426ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
14277a48dd6fSHong Zhang   if (ai[am] != 0) {
1428ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
14297a48dd6fSHong Zhang   } else {
1430ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
14317a48dd6fSHong Zhang   }
143239e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1433622e2028SHong Zhang   if (perm_identity){
1434ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1435ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1436622e2028SHong Zhang   }
1437a6175056SHong Zhang   PetscFunctionReturn(0);
1438a6175056SHong Zhang }
1439d5d45c9bSBarry Smith 
1440f76d2b81SHong Zhang #undef __FUNCT__
1441f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1442dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1443f76d2b81SHong Zhang {
1444f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
144510c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
14462ed133dbSHong Zhang   Mat            B;
1447dfbe8321SBarry Smith   PetscErrorCode ierr;
1448f76d2b81SHong Zhang   PetscTruth     perm_identity;
144910c27e3dSHong Zhang   PetscReal      fill = info->fill;
14502ed133dbSHong Zhang   PetscInt       *rip,*riip,i,mbs=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
145110c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
14522ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
145310c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
145410c27e3dSHong Zhang   PetscBT        lnkbt;
14552ed133dbSHong Zhang   IS             iperm;
1456f76d2b81SHong Zhang 
1457f76d2b81SHong Zhang   PetscFunctionBegin;
14582ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1459f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14602ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14612ed133dbSHong Zhang 
1462f76d2b81SHong Zhang   if (!perm_identity){
14632ed133dbSHong Zhang     /* check if perm is symmetric! */
14642ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
14652ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
14662ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
14672ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
14682ed133dbSHong Zhang     }
14692ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
14702ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1471f76d2b81SHong Zhang   }
147210c27e3dSHong Zhang 
147310c27e3dSHong Zhang   /* initialization */
14742ed133dbSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
147510c27e3dSHong Zhang   ui[0] = 0;
147610c27e3dSHong Zhang 
147710c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14782ed133dbSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
14792ed133dbSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
14802ed133dbSHong Zhang   il     = jl + mbs;
14812ed133dbSHong Zhang   cols   = il + mbs;
14822ed133dbSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
14832ed133dbSHong Zhang   for (i=0; i<mbs; i++){
14842ed133dbSHong Zhang     jl[i] = mbs; il[i] = 0;
1485f76d2b81SHong Zhang   }
148610c27e3dSHong Zhang 
148710c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14882ed133dbSHong Zhang   nlnk = mbs + 1;
14892ed133dbSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
149010c27e3dSHong Zhang 
14912ed133dbSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
14922ed133dbSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
149310c27e3dSHong Zhang   current_space = free_space;
149410c27e3dSHong Zhang 
14952ed133dbSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
149610c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
149710c27e3dSHong Zhang     nzk   = 0;
14982ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
14992ed133dbSHong Zhang     ncols_upper = 0;
15002ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1501622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
15022ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
15032ed133dbSHong Zhang         cols[ncols_upper] = i;
15042ed133dbSHong Zhang         ncols_upper++;
15052ed133dbSHong Zhang       }
15062ed133dbSHong Zhang     }
15072ed133dbSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
150810c27e3dSHong Zhang     nzk += nlnk;
150910c27e3dSHong Zhang 
151010c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
151110c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
151210c27e3dSHong Zhang 
151310c27e3dSHong Zhang     while (prow < k){
151410c27e3dSHong Zhang       nextprow = jl[prow];
151510c27e3dSHong Zhang       /* merge prow into k-th row */
15162ed133dbSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
151710c27e3dSHong Zhang       jmax = ui[prow+1];
151810c27e3dSHong Zhang       ncols = jmax-jmin;
15192ed133dbSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
15202ed133dbSHong Zhang       ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
152110c27e3dSHong Zhang       nzk += nlnk;
152210c27e3dSHong Zhang 
152310c27e3dSHong Zhang       /* update il and jl for prow */
152410c27e3dSHong Zhang       if (jmin < jmax){
152510c27e3dSHong Zhang         il[prow] = jmin;
15262ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
152710c27e3dSHong Zhang       }
152810c27e3dSHong Zhang       prow = nextprow;
152910c27e3dSHong Zhang     }
153010c27e3dSHong Zhang 
153110c27e3dSHong Zhang     /* if free space is not available, make more free space */
153210c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
15332ed133dbSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
153410c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
153510c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
153610c27e3dSHong Zhang       reallocs++;
153710c27e3dSHong Zhang     }
153810c27e3dSHong Zhang 
153910c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
15402ed133dbSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
154110c27e3dSHong Zhang 
154210c27e3dSHong Zhang     /* add the k-th row into il and jl */
154310c27e3dSHong Zhang     if (nzk-1 > 0){
15442ed133dbSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
154510c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
154610c27e3dSHong Zhang       il[k] = ui[k] + 1;
154710c27e3dSHong Zhang     }
15482ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
154910c27e3dSHong Zhang     current_space->array           += nzk;
155010c27e3dSHong Zhang     current_space->local_used      += nzk;
155110c27e3dSHong Zhang     current_space->local_remaining -= nzk;
155210c27e3dSHong Zhang 
155310c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
155410c27e3dSHong Zhang   }
155510c27e3dSHong Zhang 
15562ed133dbSHong Zhang   if (ai[mbs] != 0) {
155778910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[mbs]-mbs))/((PetscReal)ai[mbs]);
155878910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
155978910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
156078910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
156110c27e3dSHong Zhang   } else {
156278910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
156310c27e3dSHong Zhang   }
156410c27e3dSHong Zhang 
156510c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
156610c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
156710c27e3dSHong Zhang 
156810c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
15692ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
157010c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
157110c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
157210c27e3dSHong Zhang 
157310c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15742ed133dbSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
15752ed133dbSHong Zhang   B    = *fact;
15762ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
15772ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
157810c27e3dSHong Zhang 
15792ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
158010c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
158110c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
158210c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
158310c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
158410c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15852ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
158610c27e3dSHong Zhang   b->j    = uj;
158710c27e3dSHong Zhang   b->i    = ui;
158810c27e3dSHong Zhang   b->diag = 0;
158910c27e3dSHong Zhang   b->ilen = 0;
159010c27e3dSHong Zhang   b->imax = 0;
159110c27e3dSHong Zhang   b->row  = perm;
159210c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
159310c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
159410c27e3dSHong Zhang   b->icol = perm;
159510c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15962ed133dbSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
15972ed133dbSHong Zhang   PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
15982ed133dbSHong Zhang   b->maxnz = b->nz = ui[mbs];
159910c27e3dSHong Zhang 
16002ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
16012ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
16022ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
16032ed133dbSHong Zhang   if (ai[mbs] != 0) {
16042ed133dbSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
160510c27e3dSHong Zhang   } else {
16062ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
160710c27e3dSHong Zhang   }
160839e3d630SHong Zhang   (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
16092ed133dbSHong Zhang   if (perm_identity){
161010c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
161110c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
16122ed133dbSHong Zhang   }
1613f76d2b81SHong Zhang   PetscFunctionReturn(0);
1614f76d2b81SHong Zhang }
1615