xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 0697b6cae9af508b3938b6fbd8a1890c7adeb4de)
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 
339ecf371e4SBarry Smith       /* estimate how much additional space we will need */
340ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
341ecf371e4SBarry Smith       /* just double the memory each time */
34238baddfdSBarry Smith       PetscInt maxadd = jmax;
343ecf371e4SBarry Smith       /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */
344bbb6d6a8SBarry Smith       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
3452fb3ed76SBarry Smith       jmax += maxadd;
346ecf371e4SBarry Smith 
347ecf371e4SBarry Smith       /* allocate a longer ajnew */
34838baddfdSBarry Smith       ierr  = PetscMalloc(jmax*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr);
34938baddfdSBarry Smith       ierr  = PetscMemcpy(ajtmp,ajnew,(ainew[i])*sizeof(PetscInt));CHKERRQ(ierr);
350606d414cSSatish Balay       ierr  = PetscFree(ajnew);CHKERRQ(ierr);
351289bc588SBarry Smith       ajnew = ajtmp;
352418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
353289bc588SBarry Smith     }
354010a20caSHong Zhang     ajtmp = ajnew + ainew[i];
355289bc588SBarry Smith     fm    = fill[n];
356289bc588SBarry Smith     nzi   = 0;
3572fb3ed76SBarry Smith     im[i] = nnz;
358289bc588SBarry Smith     while (nnz--) {
359289bc588SBarry Smith       if (fm < i) nzi++;
360010a20caSHong Zhang       *ajtmp++ = fm ;
361289bc588SBarry Smith       fm       = fill[fm];
362289bc588SBarry Smith     }
363289bc588SBarry Smith     idnew[i] = ainew[i] + nzi;
364289bc588SBarry Smith   }
365464e8d44SSatish Balay   if (ai[n] != 0) {
366329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
367418422e8SSatish Balay     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
368b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af);
369b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af);
370b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n");
37105bf559cSSatish Balay   } else {
372b0a32e0cSBarry Smith     PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n");
37305bf559cSSatish Balay   }
3742fb3ed76SBarry Smith 
375898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
376898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3771987afe7SBarry Smith 
378606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
379289bc588SBarry Smith 
380289bc588SBarry Smith   /* put together the new matrix */
381f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr);
382f204ca49SKris Buschelman   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
383f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr);
384b0a32e0cSBarry Smith   PetscLogObjectParent(*B,isicol);
385416022c9SBarry Smith   b    = (Mat_SeqAIJ*)(*B)->data;
386606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
3877c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
388e8d4e0b9SBarry Smith   /* the next line frees the default space generated by the Create() */
389606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
390606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
391010a20caSHong Zhang   ierr          = PetscMalloc((ainew[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr);
392416022c9SBarry Smith   b->j          = ajnew;
393416022c9SBarry Smith   b->i          = ainew;
394416022c9SBarry Smith   b->diag       = idnew;
395416022c9SBarry Smith   b->ilen       = 0;
396416022c9SBarry Smith   b->imax       = 0;
397416022c9SBarry Smith   b->row        = isrow;
398416022c9SBarry Smith   b->col        = iscol;
399085256b3SBarry Smith   b->lu_damping        = info->damping;
40087828ca2SBarry Smith   b->lu_zeropivot      = info->zeropivot;
4012cea7109SBarry Smith   b->lu_shift          = info->shift;
4022cea7109SBarry Smith   b->lu_shift_fraction = info->shift_fraction;
403c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
404c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
40582bf6240SBarry Smith   b->icol       = isicol;
40687828ca2SBarry Smith   ierr          = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
407416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
4087fd98bd6SLois Curfman McInnes      Allocate idnew, solve_work, new a, new j */
40938baddfdSBarry Smith   PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));
410010a20caSHong Zhang   b->maxnz = b->nz = ainew[n] ;
4117fd98bd6SLois Curfman McInnes 
412329f5518SBarry Smith   (*B)->factor                 =  FACTOR_LU;
413418422e8SSatish Balay   (*B)->info.factor_mallocs    = reallocs;
414ae068f58SLois Curfman McInnes   (*B)->info.fill_ratio_given  = f;
4153a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr);
4167dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
417703deb49SSatish Balay 
418e93ccc4dSBarry Smith   if (ai[n] != 0) {
419329f5518SBarry Smith     (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
420eea2913aSSatish Balay   } else {
421eea2913aSSatish Balay     (*B)->info.fill_ratio_needed = 0.0;
422eea2913aSSatish Balay   }
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
424289bc588SBarry Smith }
4250a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
426dfbe8321SBarry Smith EXTERN PetscErrorCode Mat_AIJ_CheckInode(Mat,PetscTruth);
42741c01911SSatish Balay 
4284a2ae208SSatish Balay #undef __FUNCT__
4294a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ"
430dfbe8321SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B)
431289bc588SBarry Smith {
432416022c9SBarry Smith   Mat            C=*B;
433416022c9SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data;
43482bf6240SBarry Smith   IS             isrow = b->row,isicol = b->icol;
4356849ba73SBarry Smith   PetscErrorCode ierr;
43638baddfdSBarry Smith   PetscInt       *r,*ic,i,j,n=A->m,*ai=b->i,*aj=b->j;
43738baddfdSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ics;
438e24cfe64SHong Zhang   PetscInt       *diag_offset = b->diag,diag,*pj,nshift=0;
43987828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
440e24cfe64SHong Zhang   PetscReal      zeropivot=b->lu_zeropivot,rs,d;
4410cf777b8SBarry Smith   PetscReal      row_shift,shift_fraction,shift_amount,shift_lo=0.,shift_hi=1.,shift_top=0.;
442e24cfe64SHong Zhang   PetscTruth     lushift;
443e24cfe64SHong Zhang   PetscReal      shift_pd=b->lu_shift,shift_nonzero=b->lu_damping;
444289bc588SBarry Smith 
4453a40ed3dSBarry Smith   PetscFunctionBegin;
44678b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44778b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44887828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44987828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
450010a20caSHong Zhang   rtmps = rtmp; ics = ic;
451289bc588SBarry Smith 
4526cc28720Svictorle   if (!a->diag) {
4530cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4540cf777b8SBarry Smith   }
455e24cfe64SHong Zhang   /* if both shift schemes are chosen by user, only use shift_pd */
456e24cfe64SHong Zhang   if (shift_pd && shift_nonzero) shift_nonzero = 0.0;
457e24cfe64SHong Zhang   if (shift_pd) { /* set shift_top=max{row_shift} */
45838baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4590cf777b8SBarry Smith     shift_top = 0;
4606cc28720Svictorle     for (i=0; i<n; i++) {
461010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4626cc28720Svictorle       /* calculate amt of shift needed for this row */
4636c037d1bSvictorle       if (d<=0) {
4643aeef889SHong Zhang         row_shift = 0;
4653aeef889SHong Zhang       } else {
4663aeef889SHong Zhang         row_shift = -2*d;
4673aeef889SHong Zhang       }
468010a20caSHong Zhang       v  = a->a+aai[i];
469e24cfe64SHong Zhang       nz = aai[i+1] - aai[i];
470e24cfe64SHong Zhang       for (j=0; j<nz; j++)
4716cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4726cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4736cc28720Svictorle     }
4746cc28720Svictorle   }
4756cc28720Svictorle 
476*0697b6caSHong Zhang   shift_fraction = b->lu_shift_fraction;
477*0697b6caSHong Zhang   shift_amount   = 0;
478085256b3SBarry Smith   do {
4796cc28720Svictorle     lushift = PETSC_FALSE;
480289bc588SBarry Smith     for (i=0; i<n; i++){
481289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
482010a20caSHong Zhang       ajtmp = aj + ai[i];
48348e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
484289bc588SBarry Smith 
485289bc588SBarry Smith       /* load in initial (unfactored row) */
486416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
487010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
488010a20caSHong Zhang       v        = a->a + a->i[r[i]];
489085256b3SBarry Smith       for (j=0; j<nz; j++) {
490085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
491085256b3SBarry Smith       }
492e24cfe64SHong Zhang       rtmp[ics[r[i]]] += shift_amount; /* shift the diagonal of the matrix */
493289bc588SBarry Smith 
494010a20caSHong Zhang       row = *ajtmp++;
495f2582111SSatish Balay       while  (row < i) {
4968c37ef55SBarry Smith         pc = rtmp + row;
4978c37ef55SBarry Smith         if (*pc != 0.0) {
498010a20caSHong Zhang           pv         = b->a + diag_offset[row];
499010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
50035aab85fSBarry Smith           multiplier = *pc / *pv++;
5018c37ef55SBarry Smith           *pc        = multiplier;
502f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
503f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
504b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
505289bc588SBarry Smith         }
506010a20caSHong Zhang         row = *ajtmp++;
507289bc588SBarry Smith       }
508416022c9SBarry Smith       /* finished row so stick it into b->a */
509010a20caSHong Zhang       pv   = b->a + ai[i] ;
510010a20caSHong Zhang       pj   = b->j + ai[i] ;
5118c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
51235aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
513d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
514d3d32019SBarry Smith       rs   = 0.0;
515d3d32019SBarry Smith       for (j=0; j<nz; j++) {
516d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
517d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
518d3d32019SBarry Smith       }
519e24cfe64SHong Zhang       /* shift the diagonals when zero pivot is detected */
5206cc28720Svictorle #define MAX_NSHIFT 5
521e24cfe64SHong Zhang       if (PetscAbsScalar(pv[diag]) <= zeropivot*rs && shift_nonzero){
522e24cfe64SHong Zhang         /* force |diag(*B)| > zeropivot*rs */
523e24cfe64SHong Zhang         if (!nshift){
524e24cfe64SHong Zhang           shift_amount = shift_nonzero;
525e24cfe64SHong Zhang         } else {
526e24cfe64SHong Zhang           shift_amount *= 2.0;
527e24cfe64SHong Zhang         }
528e24cfe64SHong Zhang         lushift = PETSC_TRUE;
529e24cfe64SHong Zhang         nshift++;
530e24cfe64SHong Zhang         break;
531e24cfe64SHong Zhang       } else if (PetscRealPart(pv[diag]) <= zeropivot*rs && shift_pd){
532e24cfe64SHong Zhang         /* force *B to be diagonally dominant */
5336cc28720Svictorle 	if (nshift>MAX_NSHIFT) {
534e005ede5SBarry Smith 	  SETERRQ(PETSC_ERR_CONV_FAILED,"Unable to determine shift to enforce positive definite preconditioner");
5356cc28720Svictorle 	} else if (nshift==MAX_NSHIFT) {
5366cc28720Svictorle 	  shift_fraction = shift_hi;
5376c037d1bSvictorle 	  lushift        = PETSC_FALSE;
5386cc28720Svictorle 	} else {
5396cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5406c037d1bSvictorle 	  lushift  = PETSC_TRUE;
5416cc28720Svictorle 	}
5426cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5430cf777b8SBarry Smith         nshift++;
5440cf777b8SBarry Smith         break;
545e24cfe64SHong Zhang       } else if (PetscAbsScalar(pv[diag]) <= zeropivot*rs){
54677431f27SBarry Smith         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
547d708749eSBarry Smith       }
54835aab85fSBarry Smith     }
549e24cfe64SHong Zhang     if (shift_pd && !lushift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5506cc28720Svictorle       /*
5516c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5526cc28720Svictorle        * then try lower shift
5536cc28720Svictorle        */
5540cf777b8SBarry Smith       shift_hi       = shift_fraction;
5550cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5566cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5570cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5580cf777b8SBarry Smith       nshift++;
5596cc28720Svictorle     }
560e24cfe64SHong Zhang   } while (lushift);
5610f11f9aeSSatish Balay 
56235aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
56335aab85fSBarry Smith   for (i=0; i<n; i++) {
564010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
56535aab85fSBarry Smith   }
56635aab85fSBarry Smith 
567606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
56878b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
56978b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
570416022c9SBarry Smith   C->factor = FACTOR_LU;
5717dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
572c456f294SBarry Smith   C->assembled = PETSC_TRUE;
573b0a32e0cSBarry Smith   PetscLogFlops(C->n);
574*0697b6caSHong Zhang   if (nshift){
575*0697b6caSHong Zhang     if (shift_nonzero) {
576*0697b6caSHong Zhang       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of shift_nonzero tries %D, shift_ amount%g\n",nshift,shift_amount);
577*0697b6caSHong Zhang     } else if (shift_pd) {
5786cc28720Svictorle       b->lu_shift_fraction = shift_fraction;
57977431f27SBarry Smith       PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %D\n",shift_fraction,shift_top,nshift);
5806cc28720Svictorle     }
581*0697b6caSHong Zhang   }
5823a40ed3dSBarry Smith   PetscFunctionReturn(0);
583289bc588SBarry Smith }
5840889b5dcSvictorle 
5850889b5dcSvictorle #undef __FUNCT__
5860889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
587dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5880889b5dcSvictorle {
5890889b5dcSvictorle   PetscFunctionBegin;
5900889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5910889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5920889b5dcSvictorle   PetscFunctionReturn(0);
5930889b5dcSvictorle }
5940889b5dcSvictorle 
5950889b5dcSvictorle 
5960a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5974a2ae208SSatish Balay #undef __FUNCT__
5984a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
599dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
600da3a660dSBarry Smith {
601dfbe8321SBarry Smith   PetscErrorCode ierr;
602416022c9SBarry Smith   Mat            C;
603416022c9SBarry Smith 
6043a40ed3dSBarry Smith   PetscFunctionBegin;
6059e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
606f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
607273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
608440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
6093a40ed3dSBarry Smith   PetscFunctionReturn(0);
610da3a660dSBarry Smith }
6110a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
6124a2ae208SSatish Balay #undef __FUNCT__
6134a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
614dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6158c37ef55SBarry Smith {
616416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
617416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6186849ba73SBarry Smith   PetscErrorCode ierr;
61938baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
62038baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
62187828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6228c37ef55SBarry Smith 
6233a40ed3dSBarry Smith   PetscFunctionBegin;
6243a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
62591bf9adeSBarry Smith 
6261ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6271ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
628416022c9SBarry Smith   tmp  = a->solve_work;
6298c37ef55SBarry Smith 
6303b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6313b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6328c37ef55SBarry Smith 
6338c37ef55SBarry Smith   /* forward solve the lower triangular */
6348c37ef55SBarry Smith   tmp[0] = b[*r++];
635010a20caSHong Zhang   tmps   = tmp;
6368c37ef55SBarry Smith   for (i=1; i<n; i++) {
637010a20caSHong Zhang     v   = aa + ai[i] ;
638010a20caSHong Zhang     vi  = aj + ai[i] ;
63953e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
64053e7425aSBarry Smith     sum = b[*r++];
641ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6428c37ef55SBarry Smith     tmp[i] = sum;
6438c37ef55SBarry Smith   }
6448c37ef55SBarry Smith 
6458c37ef55SBarry Smith   /* backward solve the upper triangular */
6468c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
647010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
648010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
649416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6508c37ef55SBarry Smith     sum = tmp[i];
651ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
652010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6538c37ef55SBarry Smith   }
6548c37ef55SBarry Smith 
6553b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6563b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6571ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
659b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
6618c37ef55SBarry Smith }
662026e39d0SSatish Balay 
663930ae53cSBarry Smith /* ----------------------------------------------------------- */
6644a2ae208SSatish Balay #undef __FUNCT__
6654a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
666dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
667930ae53cSBarry Smith {
668930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6696849ba73SBarry Smith   PetscErrorCode ierr;
67038baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
671362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
672aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
67338baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
674362ced78SSatish Balay   PetscScalar    *v,sum;
675d85afc42SSatish Balay #endif
676930ae53cSBarry Smith 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
6783a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
679930ae53cSBarry Smith 
6801ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6811ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
682930ae53cSBarry Smith 
683aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6841c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6851c4feecaSSatish Balay #else
686930ae53cSBarry Smith   /* forward solve the lower triangular */
687930ae53cSBarry Smith   x[0] = b[0];
688930ae53cSBarry Smith   for (i=1; i<n; i++) {
689930ae53cSBarry Smith     ai_i = ai[i];
690930ae53cSBarry Smith     v    = aa + ai_i;
691930ae53cSBarry Smith     vi   = aj + ai_i;
692930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
693930ae53cSBarry Smith     sum  = b[i];
694930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
695930ae53cSBarry Smith     x[i] = sum;
696930ae53cSBarry Smith   }
697930ae53cSBarry Smith 
698930ae53cSBarry Smith   /* backward solve the upper triangular */
699930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
700930ae53cSBarry Smith     adiag_i = adiag[i];
701d708749eSBarry Smith     v       = aa + adiag_i + 1;
702d708749eSBarry Smith     vi      = aj + adiag_i + 1;
703930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
704930ae53cSBarry Smith     sum     = x[i];
705930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
706930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
707930ae53cSBarry Smith   }
7081c4feecaSSatish Balay #endif
709b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
7101ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7111ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7123a40ed3dSBarry Smith   PetscFunctionReturn(0);
713930ae53cSBarry Smith }
714930ae53cSBarry Smith 
7154a2ae208SSatish Balay #undef __FUNCT__
7164a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
717dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
718da3a660dSBarry Smith {
719416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
720416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7216849ba73SBarry Smith   PetscErrorCode ierr;
72238baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
72338baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
72487828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
725da3a660dSBarry Smith 
7263a40ed3dSBarry Smith   PetscFunctionBegin;
72778b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
728da3a660dSBarry Smith 
7291ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7301ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
731416022c9SBarry Smith   tmp  = a->solve_work;
732da3a660dSBarry Smith 
7333b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7343b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
735da3a660dSBarry Smith 
736da3a660dSBarry Smith   /* forward solve the lower triangular */
737da3a660dSBarry Smith   tmp[0] = b[*r++];
738da3a660dSBarry Smith   for (i=1; i<n; i++) {
739010a20caSHong Zhang     v   = aa + ai[i] ;
740010a20caSHong Zhang     vi  = aj + ai[i] ;
741416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
742da3a660dSBarry Smith     sum = b[*r++];
743010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
744da3a660dSBarry Smith     tmp[i] = sum;
745da3a660dSBarry Smith   }
746da3a660dSBarry Smith 
747da3a660dSBarry Smith   /* backward solve the upper triangular */
748da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
749010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
750010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
751416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
752da3a660dSBarry Smith     sum = tmp[i];
753010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
754010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
755da3a660dSBarry Smith     x[*c--] += tmp[i];
756da3a660dSBarry Smith   }
757da3a660dSBarry Smith 
7583b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7593b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7601ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7611ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
762b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
763898183ecSLois Curfman McInnes 
7643a40ed3dSBarry Smith   PetscFunctionReturn(0);
765da3a660dSBarry Smith }
766da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7674a2ae208SSatish Balay #undef __FUNCT__
7684a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
769dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
770da3a660dSBarry Smith {
771416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7722235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7736849ba73SBarry Smith   PetscErrorCode ierr;
77438baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
77538baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
77687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
777da3a660dSBarry Smith 
7783a40ed3dSBarry Smith   PetscFunctionBegin;
7791ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7801ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
781416022c9SBarry Smith   tmp  = a->solve_work;
782da3a660dSBarry Smith 
7832235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7842235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
785da3a660dSBarry Smith 
786da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7872235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
788da3a660dSBarry Smith 
789da3a660dSBarry Smith   /* forward solve the U^T */
790da3a660dSBarry Smith   for (i=0; i<n; i++) {
791010a20caSHong Zhang     v   = aa + diag[i] ;
792010a20caSHong Zhang     vi  = aj + diag[i] + 1;
793f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
794c38d4ed2SBarry Smith     s1  = tmp[i];
795ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
796da3a660dSBarry Smith     while (nz--) {
797010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
798da3a660dSBarry Smith     }
799c38d4ed2SBarry Smith     tmp[i] = s1;
800da3a660dSBarry Smith   }
801da3a660dSBarry Smith 
802da3a660dSBarry Smith   /* backward solve the L^T */
803da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
804010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
805010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
806f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
807f1af5d2fSBarry Smith     s1  = tmp[i];
808da3a660dSBarry Smith     while (nz--) {
809010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
810da3a660dSBarry Smith     }
811da3a660dSBarry Smith   }
812da3a660dSBarry Smith 
813da3a660dSBarry Smith   /* copy tmp into x according to permutation */
814da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
815da3a660dSBarry Smith 
8162235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8172235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8181ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
820da3a660dSBarry Smith 
821b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8223a40ed3dSBarry Smith   PetscFunctionReturn(0);
823da3a660dSBarry Smith }
824da3a660dSBarry Smith 
8254a2ae208SSatish Balay #undef __FUNCT__
8264a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
827dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
828da3a660dSBarry Smith {
829416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8302235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8316849ba73SBarry Smith   PetscErrorCode ierr;
83238baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
83338baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
83487828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8356abc6512SBarry Smith 
8363a40ed3dSBarry Smith   PetscFunctionBegin;
83723bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8386abc6512SBarry Smith 
8391ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
841416022c9SBarry Smith   tmp = a->solve_work;
8426abc6512SBarry Smith 
8432235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8442235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8456abc6512SBarry Smith 
8466abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8472235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8486abc6512SBarry Smith 
8496abc6512SBarry Smith   /* forward solve the U^T */
8506abc6512SBarry Smith   for (i=0; i<n; i++) {
851010a20caSHong Zhang     v   = aa + diag[i] ;
852010a20caSHong Zhang     vi  = aj + diag[i] + 1;
853f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8546abc6512SBarry Smith     tmp[i] *= *v++;
8556abc6512SBarry Smith     while (nz--) {
856010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8576abc6512SBarry Smith     }
8586abc6512SBarry Smith   }
8596abc6512SBarry Smith 
8606abc6512SBarry Smith   /* backward solve the L^T */
8616abc6512SBarry Smith   for (i=n-1; i>=0; i--){
862010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
863010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
864f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8656abc6512SBarry Smith     while (nz--) {
866010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8676abc6512SBarry Smith     }
8686abc6512SBarry Smith   }
8696abc6512SBarry Smith 
8706abc6512SBarry Smith   /* copy tmp into x according to permutation */
8716abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8726abc6512SBarry Smith 
8732235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8742235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8751ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8761ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8776abc6512SBarry Smith 
878b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8793a40ed3dSBarry Smith   PetscFunctionReturn(0);
880da3a660dSBarry Smith }
8819e25ed09SBarry Smith /* ----------------------------------------------------------------*/
882dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
88308480c60SBarry Smith 
8844a2ae208SSatish Balay #undef __FUNCT__
8854a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
886dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8879e25ed09SBarry Smith {
888416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8899e25ed09SBarry Smith   IS             isicol;
8906849ba73SBarry Smith   PetscErrorCode ierr;
89138baddfdSBarry Smith   PetscInt       *r,*ic,prow,n = A->m,*ai = a->i,*aj = a->j;
89238baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
893418422e8SSatish Balay   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,len, reallocs = 0,dcount = 0;
89438baddfdSBarry Smith   PetscInt       incrlev,nnz,i,levels,diagonal_fill;
89577c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
896329f5518SBarry Smith   PetscReal      f;
8979e25ed09SBarry Smith 
8983a40ed3dSBarry Smith   PetscFunctionBegin;
8996cf997b0SBarry Smith   f             = info->fill;
90038baddfdSBarry Smith   levels        = (PetscInt)info->levels;
90138baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
9024c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
90382bf6240SBarry Smith 
90408480c60SBarry Smith   /* special case that simply copies fill pattern */
905be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
906be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
90786bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9082e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
90908480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
91008480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
91108480c60SBarry Smith     if (!b->diag) {
9127c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
91308480c60SBarry Smith     }
9147c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
91508480c60SBarry Smith     b->row              = isrow;
91608480c60SBarry Smith     b->col              = iscol;
91782bf6240SBarry Smith     b->icol             = isicol;
918a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
91987828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
9202cea7109SBarry Smith     b->lu_shift         = info->shift;
9212cea7109SBarry Smith     b->lu_shift_fraction= info->shift_fraction;
92287828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
923f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
924c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
925c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9263a40ed3dSBarry Smith     PetscFunctionReturn(0);
92708480c60SBarry Smith   }
92808480c60SBarry Smith 
929898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
930898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9319e25ed09SBarry Smith 
9329e25ed09SBarry Smith   /* get new row pointers */
93338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
934010a20caSHong Zhang   ainew[0] = 0;
9359e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
93638baddfdSBarry Smith   jmax = (PetscInt)(f*(ai[n]+1));
93738baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
9389e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
93938baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajfill);CHKERRQ(ierr);
9409e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
94138baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
94256beaf04SBarry Smith   /* im is level for each filled value */
94338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
94456beaf04SBarry Smith   /* dloc is location of diagonal in factor */
94538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&dloc);CHKERRQ(ierr);
94656beaf04SBarry Smith   dloc[0]  = 0;
94756beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9486cf997b0SBarry Smith 
9496cf997b0SBarry Smith     /* copy current row into linked list */
95056beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
95129bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
952010a20caSHong Zhang     xi      = aj + ai[r[prow]];
9539e25ed09SBarry Smith     fill[n] = n;
954435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9559e25ed09SBarry Smith     while (nz--) {
9569e25ed09SBarry Smith       fm  = n;
957010a20caSHong Zhang       idx = ic[*xi++ ];
9589e25ed09SBarry Smith       do {
9599e25ed09SBarry Smith         m  = fm;
9609e25ed09SBarry Smith         fm = fill[m];
9619e25ed09SBarry Smith       } while (fm < idx);
9629e25ed09SBarry Smith       fill[m]   = idx;
9639e25ed09SBarry Smith       fill[idx] = fm;
96456beaf04SBarry Smith       im[idx]   = 0;
9659e25ed09SBarry Smith     }
9666cf997b0SBarry Smith 
9676cf997b0SBarry Smith     /* make sure diagonal entry is included */
968435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9696cf997b0SBarry Smith       fm = n;
970435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
971435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
972435faa5fSBarry Smith       fill[fm]   = prow;
9736cf997b0SBarry Smith       im[prow]   = 0;
9746cf997b0SBarry Smith       nzf++;
975335d9088SBarry Smith       dcount++;
9766cf997b0SBarry Smith     }
9776cf997b0SBarry Smith 
97856beaf04SBarry Smith     nzi = 0;
9799e25ed09SBarry Smith     row = fill[n];
98056beaf04SBarry Smith     while (row < prow) {
98156beaf04SBarry Smith       incrlev = im[row] + 1;
98256beaf04SBarry Smith       nz      = dloc[row];
983010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
984010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
985416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
98656beaf04SBarry Smith       fm      = row;
98756beaf04SBarry Smith       while (nnz-- > 0) {
988010a20caSHong Zhang         idx = *xi++ ;
98956beaf04SBarry Smith         if (*flev + incrlev > levels) {
99056beaf04SBarry Smith           flev++;
99156beaf04SBarry Smith           continue;
99256beaf04SBarry Smith         }
99356beaf04SBarry Smith         do {
9949e25ed09SBarry Smith           m  = fm;
9959e25ed09SBarry Smith           fm = fill[m];
99656beaf04SBarry Smith         } while (fm < idx);
9979e25ed09SBarry Smith         if (fm != idx) {
99856beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9999e25ed09SBarry Smith           fill[m]   = idx;
10009e25ed09SBarry Smith           fill[idx] = fm;
10019e25ed09SBarry Smith           fm        = idx;
100256beaf04SBarry Smith           nzf++;
1003ecf371e4SBarry Smith         } else {
100456beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
10059e25ed09SBarry Smith         }
100656beaf04SBarry Smith         flev++;
10079e25ed09SBarry Smith       }
10089e25ed09SBarry Smith       row = fill[row];
100956beaf04SBarry Smith       nzi++;
10109e25ed09SBarry Smith     }
10119e25ed09SBarry Smith     /* copy new filled row into permanent storage */
101256beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
1013010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
1014ecf371e4SBarry Smith 
1015ecf371e4SBarry Smith       /* estimate how much additional space we will need */
1016ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
1017ecf371e4SBarry Smith       /* just double the memory each time */
101838baddfdSBarry Smith       /*  maxadd = (PetscInt)((f*(ai[n]+!shift)*(n-prow+5))/n); */
101938baddfdSBarry Smith       PetscInt maxadd = jmax;
102056beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1021bbb6d6a8SBarry Smith       jmax += maxadd;
1022ecf371e4SBarry Smith 
1023ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
102438baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
102538baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1026606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
102756beaf04SBarry Smith       ajnew  = xi;
102838baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
102938baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1030606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
103156beaf04SBarry Smith       ajfill = xi;
1032418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
10339e25ed09SBarry Smith     }
1034010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1035010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
103656beaf04SBarry Smith     dloc[prow]  = nzi;
10379e25ed09SBarry Smith     fm          = fill[n];
103856beaf04SBarry Smith     while (nzf--) {
1039010a20caSHong Zhang       *xi++   = fm ;
104056beaf04SBarry Smith       *flev++ = im[fm];
10419e25ed09SBarry Smith       fm      = fill[fm];
10429e25ed09SBarry Smith     }
10436cf997b0SBarry Smith     /* make sure row has diagonal entry */
1044010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
104577431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10466cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10476cf997b0SBarry Smith     }
10489e25ed09SBarry Smith   }
1049606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1050898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1051898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1052606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1053606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10549e25ed09SBarry Smith 
1055f64065bbSBarry Smith   {
1056329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1057418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1058c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1059c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1060b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1061335d9088SBarry Smith     if (diagonal_fill) {
106277431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1063335d9088SBarry Smith     }
1064f64065bbSBarry Smith   }
1065bbb6d6a8SBarry Smith 
10669e25ed09SBarry Smith   /* put together the new matrix */
1067f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1068f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1069f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
1070b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1071416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1072606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10737c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1074010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10759e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1076606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1077606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1078b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1079416022c9SBarry Smith   b->j          = ajnew;
1080416022c9SBarry Smith   b->i          = ainew;
108156beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1082416022c9SBarry Smith   b->diag       = dloc;
1083416022c9SBarry Smith   b->ilen       = 0;
1084416022c9SBarry Smith   b->imax       = 0;
1085416022c9SBarry Smith   b->row        = isrow;
1086416022c9SBarry Smith   b->col        = iscol;
1087c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1088c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
108982bf6240SBarry Smith   b->icol       = isicol;
109087828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1091416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
109256beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
109338baddfdSBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
1094010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
1095a2e34c3dSBarry Smith   b->lu_damping        = info->damping;
10962cea7109SBarry Smith   b->lu_shift          = info->shift;
10972cea7109SBarry Smith   b->lu_shift_fraction = info->shift_fraction;
109887828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10999e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
11003a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
11013a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1102ae068f58SLois Curfman McInnes 
1103418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1104ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1105329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
11063a40ed3dSBarry Smith   PetscFunctionReturn(0);
11079e25ed09SBarry Smith }
1108d5d45c9bSBarry Smith 
11093c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1110a6175056SHong Zhang #undef __FUNCT__
1111a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
11122ed133dbSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *B)
1113a6175056SHong Zhang {
11142ed133dbSHong Zhang   Mat            C = *B;
11150968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
11162ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
11172ed133dbSHong Zhang   IS             ip=b->row;
11182ed133dbSHong Zhang   PetscErrorCode ierr;
11192ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
11202ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
11212ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1122622e2028SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
11232ed133dbSHong Zhang   PetscReal      damping=b->factor_damping,zeropivot=b->factor_zeropivot,shift_amount;
11242ed133dbSHong Zhang   PetscTruth     damp,chshift;
11252ed133dbSHong Zhang   PetscInt       nshift=0,ndamp=0;
1126d5d45c9bSBarry Smith 
1127a6175056SHong Zhang   PetscFunctionBegin;
11282ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
11292ed133dbSHong Zhang 
11302ed133dbSHong Zhang   /* initialization */
11312ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
11322ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
11332ed133dbSHong Zhang   jl   = il + mbs;
11342ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
11352ed133dbSHong Zhang 
11362ed133dbSHong Zhang   shift_amount = 0;
11372ed133dbSHong Zhang   do {
11382ed133dbSHong Zhang     damp = PETSC_FALSE;
11392ed133dbSHong Zhang     chshift = PETSC_FALSE;
1140*0697b6caSHong Zhang 
11412ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
11422ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1143a6175056SHong Zhang     }
11442ed133dbSHong Zhang 
11452ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1146c05c3958SHong Zhang       bval = ba + bi[k];
11472ed133dbSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
11482ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11492ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
11502ed133dbSHong Zhang         col = rip[aj[j]];
11512ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
11522ed133dbSHong Zhang           rtmp[col] = aa[j];
1153c05c3958SHong Zhang           *bval++  = 0.0; /* for in-place factorization */
11542ed133dbSHong Zhang         }
11552ed133dbSHong Zhang       }
11562ed133dbSHong Zhang       /* damp the diagonal of the matrix */
11572ed133dbSHong Zhang       if (ndamp||nshift) rtmp[k] += damping+shift_amount;
11582ed133dbSHong Zhang 
11592ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
11602ed133dbSHong Zhang       dk = rtmp[k];
11612ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
11622ed133dbSHong Zhang 
11632ed133dbSHong Zhang       while (i < k){
11642ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
11652ed133dbSHong Zhang 
11662ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
11672ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
11682ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
11692ed133dbSHong Zhang         dk += uikdi*ba[ili];
11702ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
11712ed133dbSHong Zhang 
11722ed133dbSHong Zhang         /* add multiple of row i to k-th row */
11732ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
11742ed133dbSHong Zhang         if (jmin < jmax){
11752ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
11762ed133dbSHong Zhang           /* update il and jl for row i */
11772ed133dbSHong Zhang           il[i] = jmin;
11782ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
11792ed133dbSHong Zhang         }
11802ed133dbSHong Zhang         i = nexti;
11812ed133dbSHong Zhang       }
11822ed133dbSHong Zhang 
11832ed133dbSHong Zhang       if (PetscRealPart(dk) < zeropivot && b->factor_shift){
11842ed133dbSHong Zhang 	/* calculate a shift that would make this row diagonally dominant */
11852ed133dbSHong Zhang 	PetscReal rs = PetscAbs(PetscRealPart(dk));
11862ed133dbSHong Zhang 	jmin      = bi[k]+1;
11872ed133dbSHong Zhang 	nz        = bi[k+1] - jmin;
11882ed133dbSHong Zhang 	if (nz){
11892ed133dbSHong Zhang 	  bcol = bj + jmin;
11902ed133dbSHong Zhang 	  bval = ba + jmin;
11912ed133dbSHong Zhang 	  while (nz--){
11922ed133dbSHong Zhang 	    rs += PetscAbsScalar(rtmp[*bcol++]);
11932ed133dbSHong Zhang 	  }
11942ed133dbSHong Zhang 	}
11952ed133dbSHong Zhang 	/* if this shift is less than the previous, just up the previous
11962ed133dbSHong Zhang 	   one by a bit */
11972ed133dbSHong Zhang 	shift_amount = PetscMax(rs,1.1*shift_amount);
11982ed133dbSHong Zhang 	chshift  = PETSC_TRUE;
11992ed133dbSHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
12002ed133dbSHong Zhang 	   we increase the shift until it converges. There is no guarantee that
12012ed133dbSHong Zhang 	   this algorithm converges faster or slower, or is better or worse
12022ed133dbSHong Zhang 	   than the ILU algorithm. */
12032ed133dbSHong Zhang 	nshift++;
12042ed133dbSHong Zhang 	break;
12052ed133dbSHong Zhang       }
12062ed133dbSHong Zhang       if (PetscRealPart(dk) < zeropivot){
12072ed133dbSHong Zhang         if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1);
12082ed133dbSHong Zhang         if (damping > 0.0) {
12092ed133dbSHong Zhang           if (ndamp) damping *= 2.0;
12102ed133dbSHong Zhang           damp = PETSC_TRUE;
12112ed133dbSHong Zhang           ndamp++;
12122ed133dbSHong Zhang           break;
12132ed133dbSHong Zhang         } else if (PetscAbsScalar(dk) < zeropivot){
12142ed133dbSHong Zhang           SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g",k,PetscRealPart(dk),zeropivot);
12152ed133dbSHong Zhang         } else {
12162ed133dbSHong Zhang           PetscLogInfo((PetscObject)A,"Negative pivot %g in row %D of Cholesky factorization\n",PetscRealPart(dk),k);
12172ed133dbSHong Zhang         }
12182ed133dbSHong Zhang       }
12192ed133dbSHong Zhang 
12202ed133dbSHong Zhang       /* copy data into U(k,:) */
12212ed133dbSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
12222ed133dbSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
12232ed133dbSHong Zhang       if (jmin < jmax) {
12242ed133dbSHong Zhang         for (j=jmin; j<jmax; j++){
12252ed133dbSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
12262ed133dbSHong Zhang         }
12272ed133dbSHong Zhang         /* add the k-th row into il and jl */
12282ed133dbSHong Zhang         il[k] = jmin;
12292ed133dbSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
12302ed133dbSHong Zhang       }
12312ed133dbSHong Zhang     }
12322ed133dbSHong Zhang   } while (damp||chshift);
12332ed133dbSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
12342ed133dbSHong Zhang 
12352ed133dbSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12362ed133dbSHong Zhang   C->factor       = FACTOR_CHOLESKY;
12372ed133dbSHong Zhang   C->assembled    = PETSC_TRUE;
12382ed133dbSHong Zhang   C->preallocated = PETSC_TRUE;
12392ed133dbSHong Zhang   PetscLogFlops(C->m);
12402ed133dbSHong Zhang   if (ndamp) {
12412ed133dbSHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqAIJ: number of damping tries %D damping value %g\n",ndamp,damping);
12422ed133dbSHong Zhang   }
12432ed133dbSHong Zhang   if (nshift) {
12442ed133dbSHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ diagonal shifted %D shifts\n",nshift);
12452ed133dbSHong Zhang   }
1246a6175056SHong Zhang   PetscFunctionReturn(0);
1247a6175056SHong Zhang }
12483c9e8598SHong Zhang #undef __FUNCT__
12493c9e8598SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering"
12503c9e8598SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering(Mat A,Mat *fact)
12513c9e8598SHong Zhang {
12523c9e8598SHong Zhang   Mat            C=*fact;
12533c9e8598SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
12543c9e8598SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
12553c9e8598SHong Zhang   PetscErrorCode ierr;
12563c9e8598SHong Zhang   PetscInt       i,j,am=A->m;
12573c9e8598SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1258*0697b6caSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
12593c9e8598SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1260*0697b6caSHong Zhang   PetscReal      zeropivot=b->factor_zeropivot,shift_amount,rs;
1261*0697b6caSHong Zhang   PetscTruth     chshift;
12623c9e8598SHong Zhang   PetscInt       nshift=0;
1263*0697b6caSHong Zhang   PetscReal      shift_pd=b->factor_shift,shift_nonzero=b->factor_damping;
12643c9e8598SHong Zhang 
12653c9e8598SHong Zhang   PetscFunctionBegin;
12663c9e8598SHong Zhang   /* initialization */
12673c9e8598SHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
12683c9e8598SHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
12693c9e8598SHong Zhang   jl   = il + am;
12703c9e8598SHong Zhang   rtmp = (MatScalar*)(jl + am);
12713c9e8598SHong Zhang 
12723c9e8598SHong Zhang   shift_amount = 0;
12733c9e8598SHong Zhang   do {
12743c9e8598SHong Zhang     chshift = PETSC_FALSE;
12753c9e8598SHong Zhang     for (i=0; i<am; i++) {
12763c9e8598SHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
12773c9e8598SHong Zhang     }
12783c9e8598SHong Zhang 
12793c9e8598SHong Zhang     for (k = 0; k<am; k++){
12803c9e8598SHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
12813c9e8598SHong Zhang       nz   = ai[k+1] - ai[k];
12823c9e8598SHong Zhang       acol = aj + ai[k];
12833c9e8598SHong Zhang       aval = aa + ai[k];
12843c9e8598SHong Zhang       bval = ba + bi[k];
12853c9e8598SHong Zhang       while (nz -- ){
12863c9e8598SHong Zhang         if (*acol < k) { /* skip lower triangular entries */
12873c9e8598SHong Zhang           acol++; aval++;
12883c9e8598SHong Zhang         } else {
12893c9e8598SHong Zhang           rtmp[*acol++] = *aval++;
12903c9e8598SHong Zhang           *bval++       = 0.0; /* for in-place factorization */
12913c9e8598SHong Zhang         }
12923c9e8598SHong Zhang       }
1293*0697b6caSHong Zhang       /* shift the diagonal of the matrix */
1294*0697b6caSHong Zhang       if (nshift) rtmp[k] += shift_amount;
12953c9e8598SHong Zhang 
12963c9e8598SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
12973c9e8598SHong Zhang       dk = rtmp[k];
12983c9e8598SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
12993c9e8598SHong Zhang 
13003c9e8598SHong Zhang       while (i < k){
13013c9e8598SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
13023c9e8598SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
13033c9e8598SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
13043c9e8598SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
13053c9e8598SHong Zhang         dk   += uikdi*ba[ili];
13063c9e8598SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
13073c9e8598SHong Zhang 
13083c9e8598SHong Zhang         /* add multiple of row i to k-th row ... */
13093c9e8598SHong Zhang         jmin = ili + 1;
13103c9e8598SHong Zhang         nz   = bi[i+1] - jmin;
13113c9e8598SHong Zhang         if (nz > 0){
13123c9e8598SHong Zhang           bcol = bj + jmin;
13133c9e8598SHong Zhang           bval = ba + jmin;
13143c9e8598SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
13153c9e8598SHong Zhang           /* update il and jl for i-th row */
13163c9e8598SHong Zhang           il[i] = jmin;
13173c9e8598SHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
13183c9e8598SHong Zhang         }
13193c9e8598SHong Zhang         i = nexti;
13203c9e8598SHong Zhang       }
13213c9e8598SHong Zhang 
1322*0697b6caSHong Zhang       /* shift the diagonals when zero pivot is detected */
1323*0697b6caSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
1324*0697b6caSHong Zhang       rs   = 0.0;
13253c9e8598SHong Zhang       jmin = bi[k]+1;
13263c9e8598SHong Zhang       nz   = bi[k+1] - jmin;
13273c9e8598SHong Zhang       if (nz){
13283c9e8598SHong Zhang         bcol = bj + jmin;
13293c9e8598SHong Zhang         while (nz--){
13303c9e8598SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol++]);
13313c9e8598SHong Zhang         }
13323c9e8598SHong Zhang       }
1333*0697b6caSHong Zhang       if (PetscAbsScalar(dk) <= zeropivot*rs && shift_nonzero){
1334*0697b6caSHong Zhang         /* force |diag(*B)| > zeropivot*rs */
1335*0697b6caSHong Zhang         if (!nshift){
1336*0697b6caSHong Zhang           shift_amount = shift_nonzero;
1337*0697b6caSHong Zhang         } else {
1338*0697b6caSHong Zhang           shift_amount *= 2.0;
1339*0697b6caSHong Zhang         }
1340*0697b6caSHong Zhang         chshift = PETSC_TRUE;
1341*0697b6caSHong Zhang         nshift++;
1342*0697b6caSHong Zhang         break;
1343*0697b6caSHong Zhang       } else if (PetscRealPart(dk) <= zeropivot*rs && shift_pd){
1344*0697b6caSHong Zhang         /* calculate a shift that would make this row diagonally dominant */
1345*0697b6caSHong Zhang 	shift_amount = PetscMax(rs+PetscAbs(PetscRealPart(dk)),1.1*shift_amount);
13463c9e8598SHong Zhang 	chshift      = PETSC_TRUE;
13473c9e8598SHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
13483c9e8598SHong Zhang 	   we increase the shift until it converges. There is no guarantee that
13493c9e8598SHong Zhang 	   this algorithm converges faster or slower, or is better or worse
13503c9e8598SHong Zhang 	   than the ILU algorithm. */
13513c9e8598SHong Zhang 	nshift++;
13523c9e8598SHong Zhang 	break;
1353*0697b6caSHong Zhang       } else if (PetscAbsScalar(dk) <= zeropivot*rs){
1354*0697b6caSHong Zhang         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
13553c9e8598SHong Zhang       }
13563c9e8598SHong Zhang 
13573c9e8598SHong Zhang       /* copy data into U(k,:) */
13583c9e8598SHong Zhang       ba[bi[k]] = 1.0/dk;
13593c9e8598SHong Zhang       jmin      = bi[k]+1;
13603c9e8598SHong Zhang       nz        = bi[k+1] - jmin;
13613c9e8598SHong Zhang       if (nz){
13623c9e8598SHong Zhang         bcol = bj + jmin;
13633c9e8598SHong Zhang         bval = ba + jmin;
13643c9e8598SHong Zhang         while (nz--){
13653c9e8598SHong Zhang           *bval++       = rtmp[*bcol];
13663c9e8598SHong Zhang           rtmp[*bcol++] = 0.0;
13673c9e8598SHong Zhang         }
13683c9e8598SHong Zhang         /* add k-th row into il and jl */
13693c9e8598SHong Zhang         il[k] = jmin;
13703c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
13713c9e8598SHong Zhang       }
13723c9e8598SHong Zhang     }
1373*0697b6caSHong Zhang   } while (chshift);
13743c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
13753c9e8598SHong Zhang 
13763c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
13773c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
13783c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
13793c9e8598SHong Zhang   PetscLogFlops(C->m);
13803c9e8598SHong Zhang   if (nshift){
1381*0697b6caSHong Zhang     if (shift_nonzero) {
1382*0697b6caSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering: number of shift_nonzero tries %D, shift_amount %g\n",nshift,shift_amount);
1383*0697b6caSHong Zhang     } else if (shift_pd) {
1384*0697b6caSHong Zhang       PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering: number of shift_pd tries %D, shift_amount %g\n",nshift,shift_amount);
1385*0697b6caSHong Zhang     }
13863c9e8598SHong Zhang   }
13873c9e8598SHong Zhang   PetscFunctionReturn(0);
13883c9e8598SHong Zhang }
1389a6175056SHong Zhang 
13907a48dd6fSHong Zhang #include "petscbt.h"
13917a48dd6fSHong Zhang #include "src/mat/utils/freespace.h"
1392a6175056SHong Zhang #undef __FUNCT__
1393a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1394dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1395a6175056SHong Zhang {
13960968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1397ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1398ed59401aSHong Zhang   Mat            B;
1399dfbe8321SBarry Smith   PetscErrorCode ierr;
1400653a6975SHong Zhang   PetscTruth     perm_identity;
1401622e2028SHong Zhang   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=A->m,*ui;
1402ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
1403622e2028SHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
1404ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
14057a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
14067a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
14077a48dd6fSHong Zhang   PetscBT        lnkbt;
1408a6175056SHong Zhang 
1409a6175056SHong Zhang   PetscFunctionBegin;
1410653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
14117a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14127a48dd6fSHong Zhang 
1413622e2028SHong Zhang   /* special case that simply copies fill pattern */
1414622e2028SHong Zhang   if (!levels && perm_identity) {
1415622e2028SHong Zhang     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
1416ed59401aSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1417ed59401aSHong Zhang     for (i=0; i<am; i++) {
1418ed59401aSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1419ed59401aSHong Zhang     }
1420ed59401aSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1421ed59401aSHong Zhang     B = *fact;
1422ed59401aSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1423ed59401aSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
1424ed59401aSHong Zhang 
1425ed59401aSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
1426ed59401aSHong Zhang     uj = b->j;
1427ed59401aSHong Zhang     for (i=0; i<am; i++) {
1428ed59401aSHong Zhang       aj = a->j + a->diag[i];
1429ed59401aSHong Zhang       for (j=0; j<ui[i]; j++){
1430ed59401aSHong Zhang         *uj++ = *aj++;
1431ed59401aSHong Zhang       }
1432ed59401aSHong Zhang       b->ilen[i] = ui[i];
1433ed59401aSHong Zhang     }
1434ed59401aSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
1435ed59401aSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1436ed59401aSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1437ed59401aSHong Zhang 
1438*0697b6caSHong Zhang     b->factor_zeropivot      = info->zeropivot;
1439*0697b6caSHong Zhang     b->factor_damping        = info->damping;
1440*0697b6caSHong Zhang     b->factor_shift          = info->shift;
1441ed59401aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1442622e2028SHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14433c9e8598SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)B,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering","dummyname",(FCNVOID)B->ops->choleskyfactornumeric);CHKERRQ(ierr);
14443c9e8598SHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
1445ed59401aSHong Zhang     PetscFunctionReturn(0);
1446ed59401aSHong Zhang   }
1447ed59401aSHong Zhang 
14487a48dd6fSHong Zhang   /* initialization */
14497a48dd6fSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
14507a48dd6fSHong Zhang   ui[0] = 0;
1451622e2028SHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
14527a48dd6fSHong Zhang 
14537a48dd6fSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14547a48dd6fSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14557a48dd6fSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
14567a48dd6fSHong Zhang   il         = jl + am;
14577a48dd6fSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
14587a48dd6fSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
14597a48dd6fSHong Zhang   for (i=0; i<am; i++){
14607a48dd6fSHong Zhang     jl[i] = am; il[i] = 0;
14617a48dd6fSHong Zhang   }
14627a48dd6fSHong Zhang 
14637a48dd6fSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14647a48dd6fSHong Zhang   nlnk = am + 1;
14652ed133dbSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
14667a48dd6fSHong Zhang 
14677a48dd6fSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
14687a48dd6fSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
14697a48dd6fSHong Zhang   current_space = free_space;
14707a48dd6fSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
14717a48dd6fSHong Zhang   current_space_lvl = free_space_lvl;
14727a48dd6fSHong Zhang 
14737a48dd6fSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
14747a48dd6fSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
14757a48dd6fSHong Zhang     nzk   = 0;
1476622e2028SHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
1477622e2028SHong Zhang     ncols_upper = 0;
1478622e2028SHong Zhang     cols        = cols_lvl + am;
1479622e2028SHong Zhang     for (j=0; j<ncols; j++){
1480622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
1481622e2028SHong Zhang       if (i >= k){ /* only take upper triangular entry */
1482622e2028SHong Zhang         cols[ncols_upper] = i;
1483622e2028SHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
1484622e2028SHong Zhang         ncols_upper++;
1485622e2028SHong Zhang       }
1486622e2028SHong Zhang     }
1487622e2028SHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
14887a48dd6fSHong Zhang     nzk += nlnk;
14897a48dd6fSHong Zhang 
14907a48dd6fSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
14917a48dd6fSHong Zhang     prow = jl[k]; /* 1st pivot row */
14927a48dd6fSHong Zhang 
14937a48dd6fSHong Zhang     while (prow < k){
14947a48dd6fSHong Zhang       nextprow = jl[prow];
14957a48dd6fSHong Zhang 
14967a48dd6fSHong Zhang       /* merge prow into k-th row */
14977a48dd6fSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
14987a48dd6fSHong Zhang       jmax = ui[prow+1];
14997a48dd6fSHong Zhang       ncols = jmax-jmin;
1500ed59401aSHong Zhang       i     = jmin - ui[prow];
1501ed59401aSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1502ed59401aSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
15032ed133dbSHong Zhang       ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
15047a48dd6fSHong Zhang       nzk += nlnk;
15057a48dd6fSHong Zhang 
15067a48dd6fSHong Zhang       /* update il and jl for prow */
15077a48dd6fSHong Zhang       if (jmin < jmax){
15087a48dd6fSHong Zhang         il[prow] = jmin;
15097a48dd6fSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
15107a48dd6fSHong Zhang       }
15117a48dd6fSHong Zhang       prow = nextprow;
15127a48dd6fSHong Zhang     }
15137a48dd6fSHong Zhang 
15147a48dd6fSHong Zhang     /* if free space is not available, make more free space */
15157a48dd6fSHong Zhang     if (current_space->local_remaining<nzk) {
15167a48dd6fSHong Zhang       i = am - k + 1; /* num of unfactored rows */
15177a48dd6fSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
15187a48dd6fSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
15197a48dd6fSHong Zhang       ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
15207a48dd6fSHong Zhang       reallocs++;
15217a48dd6fSHong Zhang     }
15227a48dd6fSHong Zhang 
15237a48dd6fSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
15242ed133dbSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
15257a48dd6fSHong Zhang 
15267a48dd6fSHong Zhang     /* add the k-th row into il and jl */
15277a48dd6fSHong Zhang     if (nzk-1 > 0){
15287a48dd6fSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
15297a48dd6fSHong Zhang       jl[k] = jl[i]; jl[i] = k;
15307a48dd6fSHong Zhang       il[k] = ui[k] + 1;
15317a48dd6fSHong Zhang     }
15327a48dd6fSHong Zhang     uj_ptr[k]     = current_space->array;
15337a48dd6fSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
15347a48dd6fSHong Zhang 
15357a48dd6fSHong Zhang     current_space->array           += nzk;
15367a48dd6fSHong Zhang     current_space->local_used      += nzk;
15377a48dd6fSHong Zhang     current_space->local_remaining -= nzk;
15387a48dd6fSHong Zhang 
15397a48dd6fSHong Zhang     current_space_lvl->array           += nzk;
15407a48dd6fSHong Zhang     current_space_lvl->local_used      += nzk;
15417a48dd6fSHong Zhang     current_space_lvl->local_remaining -= nzk;
15427a48dd6fSHong Zhang 
15437a48dd6fSHong Zhang     ui[k+1] = ui[k] + nzk;
15447a48dd6fSHong Zhang   }
15457a48dd6fSHong Zhang 
15467a48dd6fSHong Zhang   if (ai[am] != 0) {
1547622e2028SHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
15487a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
15497a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
15507a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
15517a48dd6fSHong Zhang   } else {
1552ed59401aSHong Zhang      PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
15537a48dd6fSHong Zhang   }
15547a48dd6fSHong Zhang 
15557a48dd6fSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
15567a48dd6fSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
15577a48dd6fSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
15587a48dd6fSHong Zhang 
15597a48dd6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
15607a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
15617a48dd6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
15622ed133dbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15637a48dd6fSHong Zhang   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
15647a48dd6fSHong Zhang 
15657a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15667a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1567ed59401aSHong Zhang   B = *fact;
1568ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1569ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
15707a48dd6fSHong Zhang 
1571ed59401aSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
15727a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
15737a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
15747a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
15757a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
15767a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15777a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
15787a48dd6fSHong Zhang   b->j    = uj;
15797a48dd6fSHong Zhang   b->i    = ui;
15807a48dd6fSHong Zhang   b->diag = 0;
15817a48dd6fSHong Zhang   b->ilen = 0;
15827a48dd6fSHong Zhang   b->imax = 0;
15837a48dd6fSHong Zhang   b->row  = perm;
15847a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
15857a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15867a48dd6fSHong Zhang   b->icol = perm;
15877a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15887a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1589ed59401aSHong Zhang   PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
15907a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
15917a48dd6fSHong Zhang 
1592ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1593ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1594ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
15957a48dd6fSHong Zhang   if (ai[am] != 0) {
1596ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
15977a48dd6fSHong Zhang   } else {
1598ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
15997a48dd6fSHong Zhang   }
1600*0697b6caSHong Zhang 
1601*0697b6caSHong Zhang   b->factor_zeropivot      = info->zeropivot;
1602*0697b6caSHong Zhang   b->factor_damping        = info->damping;
1603*0697b6caSHong Zhang   b->factor_shift          = info->shift;
1604622e2028SHong Zhang   if (perm_identity){
1605ed59401aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1606ed59401aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
160778910aadSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)B,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering","dummyname",(FCNVOID)B->ops->choleskyfactornumeric);CHKERRQ(ierr);
16083c9e8598SHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
1609622e2028SHong Zhang   } else {
1610622e2028SHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1611622e2028SHong Zhang   }
1612a6175056SHong Zhang   PetscFunctionReturn(0);
1613a6175056SHong Zhang }
1614d5d45c9bSBarry Smith 
1615f76d2b81SHong Zhang #undef __FUNCT__
1616f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1617dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1618f76d2b81SHong Zhang {
1619f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
162010c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
16212ed133dbSHong Zhang   Mat            B;
1622dfbe8321SBarry Smith   PetscErrorCode ierr;
1623f76d2b81SHong Zhang   PetscTruth     perm_identity;
162410c27e3dSHong Zhang   PetscReal      fill = info->fill;
16252ed133dbSHong Zhang   PetscInt       *rip,*riip,i,mbs=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
162610c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
16272ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
162810c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
162910c27e3dSHong Zhang   PetscBT        lnkbt;
16302ed133dbSHong Zhang   IS             iperm;
1631f76d2b81SHong Zhang 
1632f76d2b81SHong Zhang   PetscFunctionBegin;
16332ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1634f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
16352ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
16362ed133dbSHong Zhang 
1637f76d2b81SHong Zhang   if (!perm_identity){
16382ed133dbSHong Zhang     /* check if perm is symmetric! */
16392ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
16402ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
16412ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
16422ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
16432ed133dbSHong Zhang     }
16442ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
16452ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1646f76d2b81SHong Zhang   }
164710c27e3dSHong Zhang 
164810c27e3dSHong Zhang   /* initialization */
16492ed133dbSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
165010c27e3dSHong Zhang   ui[0] = 0;
165110c27e3dSHong Zhang 
165210c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
16532ed133dbSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
16542ed133dbSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
16552ed133dbSHong Zhang   il     = jl + mbs;
16562ed133dbSHong Zhang   cols   = il + mbs;
16572ed133dbSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
16582ed133dbSHong Zhang   for (i=0; i<mbs; i++){
16592ed133dbSHong Zhang     jl[i] = mbs; il[i] = 0;
1660f76d2b81SHong Zhang   }
166110c27e3dSHong Zhang 
166210c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
16632ed133dbSHong Zhang   nlnk = mbs + 1;
16642ed133dbSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
166510c27e3dSHong Zhang 
16662ed133dbSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
16672ed133dbSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
166810c27e3dSHong Zhang   current_space = free_space;
166910c27e3dSHong Zhang 
16702ed133dbSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
167110c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
167210c27e3dSHong Zhang     nzk   = 0;
16732ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
16742ed133dbSHong Zhang     ncols_upper = 0;
16752ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1676622e2028SHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
16772ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
16782ed133dbSHong Zhang         cols[ncols_upper] = i;
16792ed133dbSHong Zhang         ncols_upper++;
16802ed133dbSHong Zhang       }
16812ed133dbSHong Zhang     }
16822ed133dbSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
168310c27e3dSHong Zhang     nzk += nlnk;
168410c27e3dSHong Zhang 
168510c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
168610c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
168710c27e3dSHong Zhang 
168810c27e3dSHong Zhang     while (prow < k){
168910c27e3dSHong Zhang       nextprow = jl[prow];
169010c27e3dSHong Zhang       /* merge prow into k-th row */
16912ed133dbSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
169210c27e3dSHong Zhang       jmax = ui[prow+1];
169310c27e3dSHong Zhang       ncols = jmax-jmin;
16942ed133dbSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
16952ed133dbSHong Zhang       ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
169610c27e3dSHong Zhang       nzk += nlnk;
169710c27e3dSHong Zhang 
169810c27e3dSHong Zhang       /* update il and jl for prow */
169910c27e3dSHong Zhang       if (jmin < jmax){
170010c27e3dSHong Zhang         il[prow] = jmin;
17012ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
170210c27e3dSHong Zhang       }
170310c27e3dSHong Zhang       prow = nextprow;
170410c27e3dSHong Zhang     }
170510c27e3dSHong Zhang 
170610c27e3dSHong Zhang     /* if free space is not available, make more free space */
170710c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
17082ed133dbSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
170910c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
171010c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
171110c27e3dSHong Zhang       reallocs++;
171210c27e3dSHong Zhang     }
171310c27e3dSHong Zhang 
171410c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
17152ed133dbSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
171610c27e3dSHong Zhang 
171710c27e3dSHong Zhang     /* add the k-th row into il and jl */
171810c27e3dSHong Zhang     if (nzk-1 > 0){
17192ed133dbSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
172010c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
172110c27e3dSHong Zhang       il[k] = ui[k] + 1;
172210c27e3dSHong Zhang     }
17232ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
172410c27e3dSHong Zhang     current_space->array           += nzk;
172510c27e3dSHong Zhang     current_space->local_used      += nzk;
172610c27e3dSHong Zhang     current_space->local_remaining -= nzk;
172710c27e3dSHong Zhang 
172810c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
172910c27e3dSHong Zhang   }
173010c27e3dSHong Zhang 
17312ed133dbSHong Zhang   if (ai[mbs] != 0) {
173278910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[mbs]-mbs))/((PetscReal)ai[mbs]);
173378910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
173478910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
173578910aadSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
173610c27e3dSHong Zhang   } else {
173778910aadSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n");
173810c27e3dSHong Zhang   }
173910c27e3dSHong Zhang 
174010c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
174110c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
174210c27e3dSHong Zhang 
174310c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
17442ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
174510c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
174610c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
174710c27e3dSHong Zhang 
174810c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
17492ed133dbSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
17502ed133dbSHong Zhang   B    = *fact;
17512ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
17522ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
175310c27e3dSHong Zhang 
17542ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
175510c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
175610c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
175710c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
175810c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
175910c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
17602ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
176110c27e3dSHong Zhang   b->j    = uj;
176210c27e3dSHong Zhang   b->i    = ui;
176310c27e3dSHong Zhang   b->diag = 0;
176410c27e3dSHong Zhang   b->ilen = 0;
176510c27e3dSHong Zhang   b->imax = 0;
176610c27e3dSHong Zhang   b->row  = perm;
176710c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
176810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
176910c27e3dSHong Zhang   b->icol = perm;
177010c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
17712ed133dbSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
17722ed133dbSHong Zhang   PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
17732ed133dbSHong Zhang   b->maxnz = b->nz = ui[mbs];
177410c27e3dSHong Zhang 
17752ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
17762ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
17772ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
17782ed133dbSHong Zhang   if (ai[mbs] != 0) {
17792ed133dbSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
178010c27e3dSHong Zhang   } else {
17812ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
178210c27e3dSHong Zhang   }
17832ed133dbSHong Zhang   if (perm_identity){
178410c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
178510c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
17862ed133dbSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)(*fact),"MatCholeskyFactorNumeric_NaturalOrdering","dummyname",(FCNVOID)(*fact)->ops->choleskyfactornumeric);CHKERRQ(ierr);
17873c9e8598SHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
17882ed133dbSHong Zhang   } else {
17892ed133dbSHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
17902ed133dbSHong Zhang   }
1791f76d2b81SHong Zhang   PetscFunctionReturn(0);
1792f76d2b81SHong Zhang }
1793