xref: /petsc/src/mat/impls/aij/seq/aijfact.c (revision 2ed133db5fd632f78dd160c512dfe86a6e065f3e)
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;
43838baddfdSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj,ndamp = 0, nshift=0;
43987828ca2SBarry Smith   PetscScalar    *rtmp,*v,*pc,multiplier,*pv,*rtmps;
4400cf777b8SBarry Smith   PetscReal      damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs,d;
4410cf777b8SBarry Smith   PetscReal      row_shift,shift_fraction,shift_amount,shift_lo=0., shift_hi=1., shift_top=0.;
4420cf777b8SBarry Smith   PetscTruth     damp,lushift;
443289bc588SBarry Smith 
4443a40ed3dSBarry Smith   PetscFunctionBegin;
44578b31e54SBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
44678b31e54SBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
44787828ca2SBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr);
44887828ca2SBarry Smith   ierr  = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr);
449010a20caSHong Zhang   rtmps = rtmp; ics = ic;
450289bc588SBarry Smith 
4516cc28720Svictorle   if (!a->diag) {
4520cf777b8SBarry Smith     ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
4530cf777b8SBarry Smith   }
4540cf777b8SBarry Smith 
4550cf777b8SBarry Smith   if (b->lu_shift) { /* set max shift */
45638baddfdSBarry Smith     PetscInt *aai = a->i,*ddiag = a->diag;
4570cf777b8SBarry Smith     shift_top = 0;
4586cc28720Svictorle     for (i=0; i<n; i++) {
459010a20caSHong Zhang       d = PetscAbsScalar((a->a)[ddiag[i]]);
4606cc28720Svictorle       /* calculate amt of shift needed for this row */
4616c037d1bSvictorle       if (d<=0) {
4623aeef889SHong Zhang         row_shift = 0;
4633aeef889SHong Zhang       } else {
4643aeef889SHong Zhang         row_shift = -2*d;
4653aeef889SHong Zhang       }
466010a20caSHong Zhang       v = a->a+aai[i];
4678a2e2d9cSvictorle       for (j=0; j<aai[i+1]-aai[i]; j++)
4686cc28720Svictorle 	row_shift += PetscAbsScalar(v[j]);
4696cc28720Svictorle       if (row_shift>shift_top) shift_top = row_shift;
4706cc28720Svictorle     }
4716cc28720Svictorle   }
4726cc28720Svictorle 
4736cc28720Svictorle   shift_fraction = 0; shift_amount = 0;
474085256b3SBarry Smith   do {
4758a5eb818SBarry Smith     damp    = PETSC_FALSE;
4766cc28720Svictorle     lushift = PETSC_FALSE;
477289bc588SBarry Smith     for (i=0; i<n; i++) {
478289bc588SBarry Smith       nz    = ai[i+1] - ai[i];
479010a20caSHong Zhang       ajtmp = aj + ai[i];
48048e94559SBarry Smith       for  (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0;
481289bc588SBarry Smith 
482289bc588SBarry Smith       /* load in initial (unfactored row) */
483416022c9SBarry Smith       nz       = a->i[r[i]+1] - a->i[r[i]];
484010a20caSHong Zhang       ajtmpold = a->j + a->i[r[i]];
485010a20caSHong Zhang       v        = a->a + a->i[r[i]];
486085256b3SBarry Smith       for (j=0; j<nz; j++) {
487085256b3SBarry Smith         rtmp[ics[ajtmpold[j]]] = v[j];
488085256b3SBarry Smith       }
489e2dd4fc4Svictorle       rtmp[ics[r[i]]] += damping + shift_amount; /* damp the diagonal of the matrix */
490289bc588SBarry Smith 
491010a20caSHong Zhang       row = *ajtmp++ ;
492f2582111SSatish Balay       while  (row < i) {
4938c37ef55SBarry Smith         pc = rtmp + row;
4948c37ef55SBarry Smith         if (*pc != 0.0) {
495010a20caSHong Zhang           pv         = b->a + diag_offset[row];
496010a20caSHong Zhang           pj         = b->j + diag_offset[row] + 1;
49735aab85fSBarry Smith           multiplier = *pc / *pv++;
4988c37ef55SBarry Smith           *pc        = multiplier;
499f2582111SSatish Balay           nz         = ai[row+1] - diag_offset[row] - 1;
500f2582111SSatish Balay           for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j];
501b0a32e0cSBarry Smith           PetscLogFlops(2*nz);
502289bc588SBarry Smith         }
503010a20caSHong Zhang         row = *ajtmp++;
504289bc588SBarry Smith       }
505416022c9SBarry Smith       /* finished row so stick it into b->a */
506010a20caSHong Zhang       pv   = b->a + ai[i] ;
507010a20caSHong Zhang       pj   = b->j + ai[i] ;
5088c37ef55SBarry Smith       nz   = ai[i+1] - ai[i];
50935aab85fSBarry Smith       diag = diag_offset[i] - ai[i];
510d3d32019SBarry Smith       /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */
511d3d32019SBarry Smith       rs   = 0.0;
512d3d32019SBarry Smith       for (j=0; j<nz; j++) {
513d3d32019SBarry Smith         pv[j] = rtmps[pj[j]];
514d3d32019SBarry Smith         if (j != diag) rs += PetscAbsScalar(pv[j]);
515d3d32019SBarry Smith       }
5166cc28720Svictorle #define MAX_NSHIFT 5
517b2da817dSvictorle       if (PetscRealPart(pv[diag]) <= zeropivot*rs && b->lu_shift) {
5186cc28720Svictorle 	if (nshift>MAX_NSHIFT) {
519e005ede5SBarry Smith 	  SETERRQ(PETSC_ERR_CONV_FAILED,"Unable to determine shift to enforce positive definite preconditioner");
5206cc28720Svictorle 	} else if (nshift==MAX_NSHIFT) {
5216cc28720Svictorle 	  shift_fraction = shift_hi;
5226c037d1bSvictorle 	  lushift      = PETSC_FALSE;
5236cc28720Svictorle 	} else {
5246cc28720Svictorle 	  shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.;
5256c037d1bSvictorle 	  lushift      = PETSC_TRUE;
5266cc28720Svictorle 	}
5276cc28720Svictorle 	shift_amount = shift_fraction * shift_top;
5280cf777b8SBarry Smith         nshift++;
5290cf777b8SBarry Smith         break;
5306cc28720Svictorle       }
5319e29f15eSvictorle       if (PetscAbsScalar(pv[diag]) <= zeropivot*rs) {
532d3d32019SBarry Smith         if (damping) {
5338a5eb818SBarry Smith           if (ndamp) damping *= 2.0;
534085256b3SBarry Smith           damp = PETSC_TRUE;
535085256b3SBarry Smith           ndamp++;
536085256b3SBarry Smith           break;
537085256b3SBarry Smith         } else {
53877431f27SBarry Smith           SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs);
539d708749eSBarry Smith         }
54035aab85fSBarry Smith       }
54135aab85fSBarry Smith     }
5426cc28720Svictorle     if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) {
5436cc28720Svictorle       /*
5446c037d1bSvictorle        * if no shift in this attempt & shifting & started shifting & can refine,
5456cc28720Svictorle        * then try lower shift
5466cc28720Svictorle        */
5470cf777b8SBarry Smith       shift_hi       = shift_fraction;
5480cf777b8SBarry Smith       shift_fraction = (shift_hi+shift_lo)/2.;
5496cc28720Svictorle       shift_amount   = shift_fraction * shift_top;
5500cf777b8SBarry Smith       lushift        = PETSC_TRUE;
5510cf777b8SBarry Smith       nshift++;
5526cc28720Svictorle     }
5536cc28720Svictorle   } while (damp || lushift);
5540f11f9aeSSatish Balay 
55535aab85fSBarry Smith   /* invert diagonal entries for simplier triangular solves */
55635aab85fSBarry Smith   for (i=0; i<n; i++) {
557010a20caSHong Zhang     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
55835aab85fSBarry Smith   }
55935aab85fSBarry Smith 
560606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
56178b31e54SBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
56278b31e54SBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
563416022c9SBarry Smith   C->factor = FACTOR_LU;
5647dda7e10SBarry Smith   (*B)->ops->lufactornumeric   =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
565c456f294SBarry Smith   C->assembled = PETSC_TRUE;
566b0a32e0cSBarry Smith   PetscLogFlops(C->n);
567d3d32019SBarry Smith   if (ndamp) {
56877431f27SBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %D damping value %g\n",ndamp,damping);
569085256b3SBarry Smith   }
5706cc28720Svictorle   if (nshift) {
5716cc28720Svictorle     b->lu_shift_fraction = shift_fraction;
57277431f27SBarry Smith     PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %D\n",shift_fraction,shift_top,nshift);
5736cc28720Svictorle   }
5743a40ed3dSBarry Smith   PetscFunctionReturn(0);
575289bc588SBarry Smith }
5760889b5dcSvictorle 
5770889b5dcSvictorle #undef __FUNCT__
5780889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ"
579dfbe8321SBarry Smith PetscErrorCode MatUsePETSc_SeqAIJ(Mat A)
5800889b5dcSvictorle {
5810889b5dcSvictorle   PetscFunctionBegin;
5820889b5dcSvictorle   A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ;
5830889b5dcSvictorle   A->ops->lufactornumeric  = MatLUFactorNumeric_SeqAIJ;
5840889b5dcSvictorle   PetscFunctionReturn(0);
5850889b5dcSvictorle }
5860889b5dcSvictorle 
5870889b5dcSvictorle 
5880a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
5894a2ae208SSatish Balay #undef __FUNCT__
5904a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ"
591dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
592da3a660dSBarry Smith {
593dfbe8321SBarry Smith   PetscErrorCode ierr;
594416022c9SBarry Smith   Mat            C;
595416022c9SBarry Smith 
5963a40ed3dSBarry Smith   PetscFunctionBegin;
5979e046878SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
598f2582111SSatish Balay   ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr);
599273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
600440f4c53SSatish Balay   PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);
6013a40ed3dSBarry Smith   PetscFunctionReturn(0);
602da3a660dSBarry Smith }
6030a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */
6044a2ae208SSatish Balay #undef __FUNCT__
6054a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ"
606dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
6078c37ef55SBarry Smith {
608416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
609416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
6106849ba73SBarry Smith   PetscErrorCode ierr;
61138baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
61238baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
61387828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*tmps,*aa = a->a,sum,*v;
6148c37ef55SBarry Smith 
6153a40ed3dSBarry Smith   PetscFunctionBegin;
6163a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
61791bf9adeSBarry Smith 
6181ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6191ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
620416022c9SBarry Smith   tmp  = a->solve_work;
6218c37ef55SBarry Smith 
6223b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
6233b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
6248c37ef55SBarry Smith 
6258c37ef55SBarry Smith   /* forward solve the lower triangular */
6268c37ef55SBarry Smith   tmp[0] = b[*r++];
627010a20caSHong Zhang   tmps   = tmp;
6288c37ef55SBarry Smith   for (i=1; i<n; i++) {
629010a20caSHong Zhang     v   = aa + ai[i] ;
630010a20caSHong Zhang     vi  = aj + ai[i] ;
63153e7425aSBarry Smith     nz  = a->diag[i] - ai[i];
63253e7425aSBarry Smith     sum = b[*r++];
633ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
6348c37ef55SBarry Smith     tmp[i] = sum;
6358c37ef55SBarry Smith   }
6368c37ef55SBarry Smith 
6378c37ef55SBarry Smith   /* backward solve the upper triangular */
6388c37ef55SBarry Smith   for (i=n-1; i>=0; i--){
639010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
640010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
641416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
6428c37ef55SBarry Smith     sum = tmp[i];
643ed480e8bSBarry Smith     SPARSEDENSEMDOT(sum,tmps,v,vi,nz);
644010a20caSHong Zhang     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
6458c37ef55SBarry Smith   }
6468c37ef55SBarry Smith 
6473b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
6483b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
6491ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
6501ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
651b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
6523a40ed3dSBarry Smith   PetscFunctionReturn(0);
6538c37ef55SBarry Smith }
654026e39d0SSatish Balay 
655930ae53cSBarry Smith /* ----------------------------------------------------------- */
6564a2ae208SSatish Balay #undef __FUNCT__
6574a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering"
658dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
659930ae53cSBarry Smith {
660930ae53cSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
6616849ba73SBarry Smith   PetscErrorCode ierr;
66238baddfdSBarry Smith   PetscInt       n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag;
663362ced78SSatish Balay   PetscScalar    *x,*b,*aa = a->a;
664aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
66538baddfdSBarry Smith   PetscInt       adiag_i,i,*vi,nz,ai_i;
666362ced78SSatish Balay   PetscScalar    *v,sum;
667d85afc42SSatish Balay #endif
668930ae53cSBarry Smith 
6693a40ed3dSBarry Smith   PetscFunctionBegin;
6703a40ed3dSBarry Smith   if (!n) PetscFunctionReturn(0);
671930ae53cSBarry Smith 
6721ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
6731ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
674930ae53cSBarry Smith 
675aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
6761c4feecaSSatish Balay   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
6771c4feecaSSatish Balay #else
678930ae53cSBarry Smith   /* forward solve the lower triangular */
679930ae53cSBarry Smith   x[0] = b[0];
680930ae53cSBarry Smith   for (i=1; i<n; i++) {
681930ae53cSBarry Smith     ai_i = ai[i];
682930ae53cSBarry Smith     v    = aa + ai_i;
683930ae53cSBarry Smith     vi   = aj + ai_i;
684930ae53cSBarry Smith     nz   = adiag[i] - ai_i;
685930ae53cSBarry Smith     sum  = b[i];
686930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
687930ae53cSBarry Smith     x[i] = sum;
688930ae53cSBarry Smith   }
689930ae53cSBarry Smith 
690930ae53cSBarry Smith   /* backward solve the upper triangular */
691930ae53cSBarry Smith   for (i=n-1; i>=0; i--){
692930ae53cSBarry Smith     adiag_i = adiag[i];
693d708749eSBarry Smith     v       = aa + adiag_i + 1;
694d708749eSBarry Smith     vi      = aj + adiag_i + 1;
695930ae53cSBarry Smith     nz      = ai[i+1] - adiag_i - 1;
696930ae53cSBarry Smith     sum     = x[i];
697930ae53cSBarry Smith     while (nz--) sum -= *v++ * x[*vi++];
698930ae53cSBarry Smith     x[i]    = sum*aa[adiag_i];
699930ae53cSBarry Smith   }
7001c4feecaSSatish Balay #endif
701b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->n);
7021ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7031ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7043a40ed3dSBarry Smith   PetscFunctionReturn(0);
705930ae53cSBarry Smith }
706930ae53cSBarry Smith 
7074a2ae208SSatish Balay #undef __FUNCT__
7084a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ"
709dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
710da3a660dSBarry Smith {
711416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
712416022c9SBarry Smith   IS             iscol = a->col,isrow = a->row;
7136849ba73SBarry Smith   PetscErrorCode ierr;
71438baddfdSBarry Smith   PetscInt       *r,*c,i, n = A->m,*vi,*ai = a->i,*aj = a->j;
71538baddfdSBarry Smith   PetscInt       nz,*rout,*cout;
71687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,sum,*v;
717da3a660dSBarry Smith 
7183a40ed3dSBarry Smith   PetscFunctionBegin;
71978b31e54SBarry Smith   if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);}
720da3a660dSBarry Smith 
7211ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
723416022c9SBarry Smith   tmp  = a->solve_work;
724da3a660dSBarry Smith 
7253b2fbd54SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7263b2fbd54SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
727da3a660dSBarry Smith 
728da3a660dSBarry Smith   /* forward solve the lower triangular */
729da3a660dSBarry Smith   tmp[0] = b[*r++];
730da3a660dSBarry Smith   for (i=1; i<n; i++) {
731010a20caSHong Zhang     v   = aa + ai[i] ;
732010a20caSHong Zhang     vi  = aj + ai[i] ;
733416022c9SBarry Smith     nz  = a->diag[i] - ai[i];
734da3a660dSBarry Smith     sum = b[*r++];
735010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
736da3a660dSBarry Smith     tmp[i] = sum;
737da3a660dSBarry Smith   }
738da3a660dSBarry Smith 
739da3a660dSBarry Smith   /* backward solve the upper triangular */
740da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
741010a20caSHong Zhang     v   = aa + a->diag[i] + 1;
742010a20caSHong Zhang     vi  = aj + a->diag[i] + 1;
743416022c9SBarry Smith     nz  = ai[i+1] - a->diag[i] - 1;
744da3a660dSBarry Smith     sum = tmp[i];
745010a20caSHong Zhang     while (nz--) sum -= *v++ * tmp[*vi++ ];
746010a20caSHong Zhang     tmp[i] = sum*aa[a->diag[i]];
747da3a660dSBarry Smith     x[*c--] += tmp[i];
748da3a660dSBarry Smith   }
749da3a660dSBarry Smith 
7503b2fbd54SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
7513b2fbd54SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
7521ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
7531ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
754b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
755898183ecSLois Curfman McInnes 
7563a40ed3dSBarry Smith   PetscFunctionReturn(0);
757da3a660dSBarry Smith }
758da3a660dSBarry Smith /* -------------------------------------------------------------------*/
7594a2ae208SSatish Balay #undef __FUNCT__
7604a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ"
761dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
762da3a660dSBarry Smith {
763416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
7642235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
7656849ba73SBarry Smith   PetscErrorCode ierr;
76638baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
76738baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
76887828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v,s1;
769da3a660dSBarry Smith 
7703a40ed3dSBarry Smith   PetscFunctionBegin;
7711ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
7721ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
773416022c9SBarry Smith   tmp  = a->solve_work;
774da3a660dSBarry Smith 
7752235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
7762235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
777da3a660dSBarry Smith 
778da3a660dSBarry Smith   /* copy the b into temp work space according to permutation */
7792235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
780da3a660dSBarry Smith 
781da3a660dSBarry Smith   /* forward solve the U^T */
782da3a660dSBarry Smith   for (i=0; i<n; i++) {
783010a20caSHong Zhang     v   = aa + diag[i] ;
784010a20caSHong Zhang     vi  = aj + diag[i] + 1;
785f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
786c38d4ed2SBarry Smith     s1  = tmp[i];
787ef66eb69SBarry Smith     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
788da3a660dSBarry Smith     while (nz--) {
789010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*s1;
790da3a660dSBarry Smith     }
791c38d4ed2SBarry Smith     tmp[i] = s1;
792da3a660dSBarry Smith   }
793da3a660dSBarry Smith 
794da3a660dSBarry Smith   /* backward solve the L^T */
795da3a660dSBarry Smith   for (i=n-1; i>=0; i--){
796010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
797010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
798f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
799f1af5d2fSBarry Smith     s1  = tmp[i];
800da3a660dSBarry Smith     while (nz--) {
801010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*s1;
802da3a660dSBarry Smith     }
803da3a660dSBarry Smith   }
804da3a660dSBarry Smith 
805da3a660dSBarry Smith   /* copy tmp into x according to permutation */
806da3a660dSBarry Smith   for (i=0; i<n; i++) x[r[i]] = tmp[i];
807da3a660dSBarry Smith 
8082235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8092235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8101ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8111ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
812da3a660dSBarry Smith 
813b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz-A->n);
8143a40ed3dSBarry Smith   PetscFunctionReturn(0);
815da3a660dSBarry Smith }
816da3a660dSBarry Smith 
8174a2ae208SSatish Balay #undef __FUNCT__
8184a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ"
819dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
820da3a660dSBarry Smith {
821416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
8222235e667SBarry Smith   IS             iscol = a->col,isrow = a->row;
8236849ba73SBarry Smith   PetscErrorCode ierr;
82438baddfdSBarry Smith   PetscInt       *r,*c,i,n = A->m,*vi,*ai = a->i,*aj = a->j;
82538baddfdSBarry Smith   PetscInt       nz,*rout,*cout,*diag = a->diag;
82687828ca2SBarry Smith   PetscScalar    *x,*b,*tmp,*aa = a->a,*v;
8276abc6512SBarry Smith 
8283a40ed3dSBarry Smith   PetscFunctionBegin;
82923bc6035SMatthew Knepley   if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);}
8306abc6512SBarry Smith 
8311ebc52fbSHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
8321ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
833416022c9SBarry Smith   tmp = a->solve_work;
8346abc6512SBarry Smith 
8352235e667SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
8362235e667SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
8376abc6512SBarry Smith 
8386abc6512SBarry Smith   /* copy the b into temp work space according to permutation */
8392235e667SBarry Smith   for (i=0; i<n; i++) tmp[i] = b[c[i]];
8406abc6512SBarry Smith 
8416abc6512SBarry Smith   /* forward solve the U^T */
8426abc6512SBarry Smith   for (i=0; i<n; i++) {
843010a20caSHong Zhang     v   = aa + diag[i] ;
844010a20caSHong Zhang     vi  = aj + diag[i] + 1;
845f1af5d2fSBarry Smith     nz  = ai[i+1] - diag[i] - 1;
8466abc6512SBarry Smith     tmp[i] *= *v++;
8476abc6512SBarry Smith     while (nz--) {
848010a20caSHong Zhang       tmp[*vi++ ] -= (*v++)*tmp[i];
8496abc6512SBarry Smith     }
8506abc6512SBarry Smith   }
8516abc6512SBarry Smith 
8526abc6512SBarry Smith   /* backward solve the L^T */
8536abc6512SBarry Smith   for (i=n-1; i>=0; i--){
854010a20caSHong Zhang     v   = aa + diag[i] - 1 ;
855010a20caSHong Zhang     vi  = aj + diag[i] - 1 ;
856f1af5d2fSBarry Smith     nz  = diag[i] - ai[i];
8576abc6512SBarry Smith     while (nz--) {
858010a20caSHong Zhang       tmp[*vi-- ] -= (*v--)*tmp[i];
8596abc6512SBarry Smith     }
8606abc6512SBarry Smith   }
8616abc6512SBarry Smith 
8626abc6512SBarry Smith   /* copy tmp into x according to permutation */
8636abc6512SBarry Smith   for (i=0; i<n; i++) x[r[i]] += tmp[i];
8646abc6512SBarry Smith 
8652235e667SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
8662235e667SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
8671ebc52fbSHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
8681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8696abc6512SBarry Smith 
870b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
872da3a660dSBarry Smith }
8739e25ed09SBarry Smith /* ----------------------------------------------------------------*/
874dfbe8321SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat);
87508480c60SBarry Smith 
8764a2ae208SSatish Balay #undef __FUNCT__
8774a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ"
878dfbe8321SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
8799e25ed09SBarry Smith {
880416022c9SBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
8819e25ed09SBarry Smith   IS             isicol;
8826849ba73SBarry Smith   PetscErrorCode ierr;
88338baddfdSBarry Smith   PetscInt       *r,*ic,prow,n = A->m,*ai = a->i,*aj = a->j;
88438baddfdSBarry Smith   PetscInt       *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev;
885418422e8SSatish Balay   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,len, reallocs = 0,dcount = 0;
88638baddfdSBarry Smith   PetscInt       incrlev,nnz,i,levels,diagonal_fill;
88777c4ece6SBarry Smith   PetscTruth     col_identity,row_identity;
888329f5518SBarry Smith   PetscReal      f;
8899e25ed09SBarry Smith 
8903a40ed3dSBarry Smith   PetscFunctionBegin;
8916cf997b0SBarry Smith   f             = info->fill;
89238baddfdSBarry Smith   levels        = (PetscInt)info->levels;
89338baddfdSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
8944c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
89582bf6240SBarry Smith 
89608480c60SBarry Smith   /* special case that simply copies fill pattern */
897be0abb6dSBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
898be0abb6dSBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
89986bacbd2SBarry Smith   if (!levels && row_identity && col_identity) {
9002e8a6d31SBarry Smith     ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr);
90108480c60SBarry Smith     (*fact)->factor = FACTOR_LU;
90208480c60SBarry Smith     b               = (Mat_SeqAIJ*)(*fact)->data;
90308480c60SBarry Smith     if (!b->diag) {
9047c922b88SBarry Smith       ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90508480c60SBarry Smith     }
9067c922b88SBarry Smith     ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr);
90708480c60SBarry Smith     b->row              = isrow;
90808480c60SBarry Smith     b->col              = iscol;
90982bf6240SBarry Smith     b->icol             = isicol;
910a2e34c3dSBarry Smith     b->lu_damping       = info->damping;
91187828ca2SBarry Smith     b->lu_zeropivot     = info->zeropivot;
9122cea7109SBarry Smith     b->lu_shift         = info->shift;
9132cea7109SBarry Smith     b->lu_shift_fraction= info->shift_fraction;
91487828ca2SBarry Smith     ierr                = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
915f830108cSBarry Smith     (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
916c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
917c38d4ed2SBarry Smith     ierr                = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
9183a40ed3dSBarry Smith     PetscFunctionReturn(0);
91908480c60SBarry Smith   }
92008480c60SBarry Smith 
921898183ecSLois Curfman McInnes   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
922898183ecSLois Curfman McInnes   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
9239e25ed09SBarry Smith 
9249e25ed09SBarry Smith   /* get new row pointers */
92538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&ainew);CHKERRQ(ierr);
926010a20caSHong Zhang   ainew[0] = 0;
9279e25ed09SBarry Smith   /* don't know how many column pointers are needed so estimate */
92838baddfdSBarry Smith   jmax = (PetscInt)(f*(ai[n]+1));
92938baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajnew);CHKERRQ(ierr);
9309e25ed09SBarry Smith   /* ajfill is level of fill for each fill entry */
93138baddfdSBarry Smith   ierr = PetscMalloc((jmax)*sizeof(PetscInt),&ajfill);CHKERRQ(ierr);
9329e25ed09SBarry Smith   /* fill is a linked list of nonzeros in active row */
93338baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
93456beaf04SBarry Smith   /* im is level for each filled value */
93538baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
93656beaf04SBarry Smith   /* dloc is location of diagonal in factor */
93738baddfdSBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&dloc);CHKERRQ(ierr);
93856beaf04SBarry Smith   dloc[0]  = 0;
93956beaf04SBarry Smith   for (prow=0; prow<n; prow++) {
9406cf997b0SBarry Smith 
9416cf997b0SBarry Smith     /* copy current row into linked list */
94256beaf04SBarry Smith     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
94329bbc08cSBarry Smith     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix");
944010a20caSHong Zhang     xi      = aj + ai[r[prow]] ;
9459e25ed09SBarry Smith     fill[n]    = n;
946435faa5fSBarry Smith     fill[prow] = -1; /* marker to indicate if diagonal exists */
9479e25ed09SBarry Smith     while (nz--) {
9489e25ed09SBarry Smith       fm  = n;
949010a20caSHong Zhang       idx = ic[*xi++ ];
9509e25ed09SBarry Smith       do {
9519e25ed09SBarry Smith         m  = fm;
9529e25ed09SBarry Smith         fm = fill[m];
9539e25ed09SBarry Smith       } while (fm < idx);
9549e25ed09SBarry Smith       fill[m]   = idx;
9559e25ed09SBarry Smith       fill[idx] = fm;
95656beaf04SBarry Smith       im[idx]   = 0;
9579e25ed09SBarry Smith     }
9586cf997b0SBarry Smith 
9596cf997b0SBarry Smith     /* make sure diagonal entry is included */
960435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
9616cf997b0SBarry Smith       fm = n;
962435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
963435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
964435faa5fSBarry Smith       fill[fm]   = prow;
9656cf997b0SBarry Smith       im[prow]   = 0;
9666cf997b0SBarry Smith       nzf++;
967335d9088SBarry Smith       dcount++;
9686cf997b0SBarry Smith     }
9696cf997b0SBarry Smith 
97056beaf04SBarry Smith     nzi = 0;
9719e25ed09SBarry Smith     row = fill[n];
97256beaf04SBarry Smith     while (row < prow) {
97356beaf04SBarry Smith       incrlev = im[row] + 1;
97456beaf04SBarry Smith       nz      = dloc[row];
975010a20caSHong Zhang       xi      = ajnew  + ainew[row]  + nz + 1;
976010a20caSHong Zhang       flev    = ajfill + ainew[row]  + nz + 1;
977416022c9SBarry Smith       nnz     = ainew[row+1] - ainew[row] - nz - 1;
97856beaf04SBarry Smith       fm      = row;
97956beaf04SBarry Smith       while (nnz-- > 0) {
980010a20caSHong Zhang         idx = *xi++ ;
98156beaf04SBarry Smith         if (*flev + incrlev > levels) {
98256beaf04SBarry Smith           flev++;
98356beaf04SBarry Smith           continue;
98456beaf04SBarry Smith         }
98556beaf04SBarry Smith         do {
9869e25ed09SBarry Smith           m  = fm;
9879e25ed09SBarry Smith           fm = fill[m];
98856beaf04SBarry Smith         } while (fm < idx);
9899e25ed09SBarry Smith         if (fm != idx) {
99056beaf04SBarry Smith           im[idx]   = *flev + incrlev;
9919e25ed09SBarry Smith           fill[m]   = idx;
9929e25ed09SBarry Smith           fill[idx] = fm;
9939e25ed09SBarry Smith           fm        = idx;
99456beaf04SBarry Smith           nzf++;
995ecf371e4SBarry Smith         } else {
99656beaf04SBarry Smith           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
9979e25ed09SBarry Smith         }
99856beaf04SBarry Smith         flev++;
9999e25ed09SBarry Smith       }
10009e25ed09SBarry Smith       row = fill[row];
100156beaf04SBarry Smith       nzi++;
10029e25ed09SBarry Smith     }
10039e25ed09SBarry Smith     /* copy new filled row into permanent storage */
100456beaf04SBarry Smith     ainew[prow+1] = ainew[prow] + nzf;
1005010a20caSHong Zhang     if (ainew[prow+1] > jmax) {
1006ecf371e4SBarry Smith 
1007ecf371e4SBarry Smith       /* estimate how much additional space we will need */
1008ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
1009ecf371e4SBarry Smith       /* just double the memory each time */
101038baddfdSBarry Smith       /*  maxadd = (PetscInt)((f*(ai[n]+!shift)*(n-prow+5))/n); */
101138baddfdSBarry Smith       PetscInt maxadd = jmax;
101256beaf04SBarry Smith       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1013bbb6d6a8SBarry Smith       jmax += maxadd;
1014ecf371e4SBarry Smith 
1015ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
101638baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
101738baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1018606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
101956beaf04SBarry Smith       ajnew  = xi;
102038baddfdSBarry Smith       ierr   = PetscMalloc(jmax*sizeof(PetscInt),&xi);CHKERRQ(ierr);
102138baddfdSBarry Smith       ierr   = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(PetscInt));CHKERRQ(ierr);
1022606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
102356beaf04SBarry Smith       ajfill = xi;
1024418422e8SSatish Balay       reallocs++; /* count how many times we realloc */
10259e25ed09SBarry Smith     }
1026010a20caSHong Zhang     xi          = ajnew + ainew[prow] ;
1027010a20caSHong Zhang     flev        = ajfill + ainew[prow] ;
102856beaf04SBarry Smith     dloc[prow]  = nzi;
10299e25ed09SBarry Smith     fm          = fill[n];
103056beaf04SBarry Smith     while (nzf--) {
1031010a20caSHong Zhang       *xi++   = fm ;
103256beaf04SBarry Smith       *flev++ = im[fm];
10339e25ed09SBarry Smith       fm      = fill[fm];
10349e25ed09SBarry Smith     }
10356cf997b0SBarry Smith     /* make sure row has diagonal entry */
1036010a20caSHong Zhang     if (ajnew[ainew[prow]+dloc[prow]] != prow) {
103777431f27SBarry Smith       SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
10386cf997b0SBarry Smith     try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow);
10396cf997b0SBarry Smith     }
10409e25ed09SBarry Smith   }
1041606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
1042898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1043898183ecSLois Curfman McInnes   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1044606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
1045606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
10469e25ed09SBarry Smith 
1047f64065bbSBarry Smith   {
1048329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
1049418422e8SSatish Balay     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af);
1050c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af);
1051c0d6ac57SBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af);
1052b0a32e0cSBarry Smith     PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n");
1053335d9088SBarry Smith     if (diagonal_fill) {
105477431f27SBarry Smith       PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %D missing diagonals",dcount);
1055335d9088SBarry Smith     }
1056f64065bbSBarry Smith   }
1057bbb6d6a8SBarry Smith 
10589e25ed09SBarry Smith   /* put together the new matrix */
1059f204ca49SKris Buschelman   ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr);
1060f204ca49SKris Buschelman   ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr);
1061f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr);
1062b0a32e0cSBarry Smith   PetscLogObjectParent(*fact,isicol);
1063416022c9SBarry Smith   b = (Mat_SeqAIJ*)(*fact)->data;
1064606d414cSSatish Balay   ierr = PetscFree(b->imax);CHKERRQ(ierr);
10657c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
1066010a20caSHong Zhang   len = (ainew[n] )*sizeof(PetscScalar);
10679e25ed09SBarry Smith   /* the next line frees the default space generated by the Create() */
1068606d414cSSatish Balay   ierr = PetscFree(b->a);CHKERRQ(ierr);
1069606d414cSSatish Balay   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1070b0a32e0cSBarry Smith   ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr);
1071416022c9SBarry Smith   b->j          = ajnew;
1072416022c9SBarry Smith   b->i          = ainew;
107356beaf04SBarry Smith   for (i=0; i<n; i++) dloc[i] += ainew[i];
1074416022c9SBarry Smith   b->diag       = dloc;
1075416022c9SBarry Smith   b->ilen       = 0;
1076416022c9SBarry Smith   b->imax       = 0;
1077416022c9SBarry Smith   b->row        = isrow;
1078416022c9SBarry Smith   b->col        = iscol;
1079c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1080c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
108182bf6240SBarry Smith   b->icol       = isicol;
108287828ca2SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1083416022c9SBarry Smith   /* In b structure:  Free imax, ilen, old a, old j.
108456beaf04SBarry Smith      Allocate dloc, solve_work, new a, new j */
108538baddfdSBarry Smith   PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
1086010a20caSHong Zhang   b->maxnz             = b->nz = ainew[n] ;
1087a2e34c3dSBarry Smith   b->lu_damping        = info->damping;
10882cea7109SBarry Smith   b->lu_shift          = info->shift;
10892cea7109SBarry Smith   b->lu_shift_fraction = info->shift_fraction;
109087828ca2SBarry Smith   b->lu_zeropivot = info->zeropivot;
10919e25ed09SBarry Smith   (*fact)->factor = FACTOR_LU;
10923a7fca6bSBarry Smith   ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr);
10933a7fca6bSBarry Smith   (*fact)->ops->lufactornumeric =  A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */
1094ae068f58SLois Curfman McInnes 
1095418422e8SSatish Balay   (*fact)->info.factor_mallocs    = reallocs;
1096ae068f58SLois Curfman McInnes   (*fact)->info.fill_ratio_given  = f;
1097329f5518SBarry Smith   (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
10983a40ed3dSBarry Smith   PetscFunctionReturn(0);
10999e25ed09SBarry Smith }
1100d5d45c9bSBarry Smith 
11013c9e8598SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
1102a6175056SHong Zhang #undef __FUNCT__
1103a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ"
1104*2ed133dbSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *B)
1105a6175056SHong Zhang {
1106*2ed133dbSHong Zhang   Mat            C = *B;
11070968510aSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1108*2ed133dbSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
1109*2ed133dbSHong Zhang   IS             ip=b->row;
1110*2ed133dbSHong Zhang   PetscErrorCode ierr;
1111*2ed133dbSHong Zhang   PetscInt       *rip,i,j,mbs=A->m,*bi=b->i,*bj=b->j,*bcol;
1112*2ed133dbSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
1113*2ed133dbSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1114*2ed133dbSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1115*2ed133dbSHong Zhang   PetscReal      damping=b->factor_damping,zeropivot=b->factor_zeropivot,shift_amount;
1116*2ed133dbSHong Zhang   PetscTruth     damp,chshift;
1117*2ed133dbSHong Zhang   PetscInt       nshift=0,ndamp=0;
1118d5d45c9bSBarry Smith 
1119a6175056SHong Zhang   PetscFunctionBegin;
1120*2ed133dbSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1121*2ed133dbSHong Zhang 
1122*2ed133dbSHong Zhang   /* initialization */
1123*2ed133dbSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1124*2ed133dbSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1125*2ed133dbSHong Zhang   jl   = il + mbs;
1126*2ed133dbSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
1127*2ed133dbSHong Zhang 
1128*2ed133dbSHong Zhang   shift_amount = 0;
1129*2ed133dbSHong Zhang   do {
1130*2ed133dbSHong Zhang     damp = PETSC_FALSE;
1131*2ed133dbSHong Zhang     chshift = PETSC_FALSE;
1132*2ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
1133*2ed133dbSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1134a6175056SHong Zhang     }
1135*2ed133dbSHong Zhang 
1136*2ed133dbSHong Zhang     for (k = 0; k<mbs; k++){
1137*2ed133dbSHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
1138*2ed133dbSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1139*2ed133dbSHong Zhang       PetscScalar *baval = ba + bi[k];
1140*2ed133dbSHong Zhang       for (j = jmin; j < jmax; j++){
1141*2ed133dbSHong Zhang         col = rip[aj[j]];
1142*2ed133dbSHong Zhang         if (col >= k){ /* only take upper triangular entry */
1143*2ed133dbSHong Zhang           rtmp[col] = aa[j];
1144*2ed133dbSHong Zhang           *baval++  = 0.0; /* for in-place factorization */
1145*2ed133dbSHong Zhang         }
1146*2ed133dbSHong Zhang       }
1147*2ed133dbSHong Zhang       /* damp the diagonal of the matrix */
1148*2ed133dbSHong Zhang       if (ndamp||nshift) rtmp[k] += damping+shift_amount;
1149*2ed133dbSHong Zhang 
1150*2ed133dbSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1151*2ed133dbSHong Zhang       dk = rtmp[k];
1152*2ed133dbSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
1153*2ed133dbSHong Zhang 
1154*2ed133dbSHong Zhang       while (i < k){
1155*2ed133dbSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1156*2ed133dbSHong Zhang 
1157*2ed133dbSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1158*2ed133dbSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1159*2ed133dbSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1160*2ed133dbSHong Zhang         dk += uikdi*ba[ili];
1161*2ed133dbSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1162*2ed133dbSHong Zhang 
1163*2ed133dbSHong Zhang         /* add multiple of row i to k-th row */
1164*2ed133dbSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1165*2ed133dbSHong Zhang         if (jmin < jmax){
1166*2ed133dbSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1167*2ed133dbSHong Zhang           /* update il and jl for row i */
1168*2ed133dbSHong Zhang           il[i] = jmin;
1169*2ed133dbSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1170*2ed133dbSHong Zhang         }
1171*2ed133dbSHong Zhang         i = nexti;
1172*2ed133dbSHong Zhang       }
1173*2ed133dbSHong Zhang 
1174*2ed133dbSHong Zhang       if (PetscRealPart(dk) < zeropivot && b->factor_shift){
1175*2ed133dbSHong Zhang 	/* calculate a shift that would make this row diagonally dominant */
1176*2ed133dbSHong Zhang 	PetscReal rs = PetscAbs(PetscRealPart(dk));
1177*2ed133dbSHong Zhang 	jmin      = bi[k]+1;
1178*2ed133dbSHong Zhang 	nz        = bi[k+1] - jmin;
1179*2ed133dbSHong Zhang 	if (nz){
1180*2ed133dbSHong Zhang 	  bcol = bj + jmin;
1181*2ed133dbSHong Zhang 	  bval = ba + jmin;
1182*2ed133dbSHong Zhang 	  while (nz--){
1183*2ed133dbSHong Zhang 	    rs += PetscAbsScalar(rtmp[*bcol++]);
1184*2ed133dbSHong Zhang 	  }
1185*2ed133dbSHong Zhang 	}
1186*2ed133dbSHong Zhang 	/* if this shift is less than the previous, just up the previous
1187*2ed133dbSHong Zhang 	   one by a bit */
1188*2ed133dbSHong Zhang 	shift_amount = PetscMax(rs,1.1*shift_amount);
1189*2ed133dbSHong Zhang 	chshift  = PETSC_TRUE;
1190*2ed133dbSHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
1191*2ed133dbSHong Zhang 	   we increase the shift until it converges. There is no guarantee that
1192*2ed133dbSHong Zhang 	   this algorithm converges faster or slower, or is better or worse
1193*2ed133dbSHong Zhang 	   than the ILU algorithm. */
1194*2ed133dbSHong Zhang 	nshift++;
1195*2ed133dbSHong Zhang 	break;
1196*2ed133dbSHong Zhang       }
1197*2ed133dbSHong Zhang       if (PetscRealPart(dk) < zeropivot){
1198*2ed133dbSHong Zhang         if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1);
1199*2ed133dbSHong Zhang         if (damping > 0.0) {
1200*2ed133dbSHong Zhang           if (ndamp) damping *= 2.0;
1201*2ed133dbSHong Zhang           damp = PETSC_TRUE;
1202*2ed133dbSHong Zhang           ndamp++;
1203*2ed133dbSHong Zhang           break;
1204*2ed133dbSHong Zhang         } else if (PetscAbsScalar(dk) < zeropivot){
1205*2ed133dbSHong Zhang           SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g",k,PetscRealPart(dk),zeropivot);
1206*2ed133dbSHong Zhang         } else {
1207*2ed133dbSHong Zhang           PetscLogInfo((PetscObject)A,"Negative pivot %g in row %D of Cholesky factorization\n",PetscRealPart(dk),k);
1208*2ed133dbSHong Zhang         }
1209*2ed133dbSHong Zhang       }
1210*2ed133dbSHong Zhang 
1211*2ed133dbSHong Zhang       /* copy data into U(k,:) */
1212*2ed133dbSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1213*2ed133dbSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
1214*2ed133dbSHong Zhang       if (jmin < jmax) {
1215*2ed133dbSHong Zhang         for (j=jmin; j<jmax; j++){
1216*2ed133dbSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1217*2ed133dbSHong Zhang         }
1218*2ed133dbSHong Zhang         /* add the k-th row into il and jl */
1219*2ed133dbSHong Zhang         il[k] = jmin;
1220*2ed133dbSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1221*2ed133dbSHong Zhang       }
1222*2ed133dbSHong Zhang     }
1223*2ed133dbSHong Zhang   } while (damp||chshift);
1224*2ed133dbSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
1225*2ed133dbSHong Zhang 
1226*2ed133dbSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1227*2ed133dbSHong Zhang   C->factor       = FACTOR_CHOLESKY;
1228*2ed133dbSHong Zhang   C->assembled    = PETSC_TRUE;
1229*2ed133dbSHong Zhang   C->preallocated = PETSC_TRUE;
1230*2ed133dbSHong Zhang   PetscLogFlops(C->m);
1231*2ed133dbSHong Zhang   if (ndamp) {
1232*2ed133dbSHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqAIJ: number of damping tries %D damping value %g\n",ndamp,damping);
1233*2ed133dbSHong Zhang   }
1234*2ed133dbSHong Zhang   if (nshift) {
1235*2ed133dbSHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ diagonal shifted %D shifts\n",nshift);
1236*2ed133dbSHong Zhang   }
1237a6175056SHong Zhang   PetscFunctionReturn(0);
1238a6175056SHong Zhang }
12393c9e8598SHong Zhang #undef __FUNCT__
12403c9e8598SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering"
12413c9e8598SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering(Mat A,Mat *fact)
12423c9e8598SHong Zhang {
12433c9e8598SHong Zhang   Mat            C = *fact;
12443c9e8598SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
12453c9e8598SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
12463c9e8598SHong Zhang   PetscErrorCode ierr;
12473c9e8598SHong Zhang   PetscInt       i,j,am = A->m;
12483c9e8598SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
12493c9e8598SHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz,ndamp = 0;
12503c9e8598SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
12513c9e8598SHong Zhang   PetscReal      damping=b->factor_damping, zeropivot=b->factor_zeropivot,shift_amount;
12523c9e8598SHong Zhang   PetscTruth     damp,chshift;
12533c9e8598SHong Zhang   PetscInt       nshift=0;
12543c9e8598SHong Zhang 
12553c9e8598SHong Zhang   PetscFunctionBegin;
12563c9e8598SHong Zhang   /* initialization */
12573c9e8598SHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
12583c9e8598SHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
12593c9e8598SHong Zhang   jl   = il + am;
12603c9e8598SHong Zhang   rtmp = (MatScalar*)(jl + am);
12613c9e8598SHong Zhang 
12623c9e8598SHong Zhang   shift_amount = 0;
12633c9e8598SHong Zhang   do {
12643c9e8598SHong Zhang     damp = PETSC_FALSE;
12653c9e8598SHong Zhang     chshift = PETSC_FALSE;
12663c9e8598SHong Zhang     for (i=0; i<am; i++) {
12673c9e8598SHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
12683c9e8598SHong Zhang     }
12693c9e8598SHong Zhang 
12703c9e8598SHong Zhang     for (k = 0; k<am; k++){
12713c9e8598SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
12723c9e8598SHong Zhang       nz   = ai[k+1] - ai[k];
12733c9e8598SHong Zhang       acol = aj + ai[k];
12743c9e8598SHong Zhang       aval = aa + ai[k];
12753c9e8598SHong Zhang       bval = ba + bi[k];
12763c9e8598SHong Zhang       while (nz -- ){
12773c9e8598SHong Zhang         if (*acol < k) { /* skip lower triangular entries */
12783c9e8598SHong Zhang           acol++; aval++;
12793c9e8598SHong Zhang         } else {
12803c9e8598SHong Zhang           rtmp[*acol++] = *aval++;
12813c9e8598SHong Zhang           *bval++       = 0.0; /* for in-place factorization */
12823c9e8598SHong Zhang         }
12833c9e8598SHong Zhang       }
12843c9e8598SHong Zhang       /* damp the diagonal of the matrix */
12853c9e8598SHong Zhang       if (ndamp||nshift) rtmp[k] += damping+shift_amount;
12863c9e8598SHong Zhang 
12873c9e8598SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
12883c9e8598SHong Zhang       dk = rtmp[k];
12893c9e8598SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
12903c9e8598SHong Zhang 
12913c9e8598SHong Zhang       while (i < k){
12923c9e8598SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
12933c9e8598SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
12943c9e8598SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
12953c9e8598SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
12963c9e8598SHong Zhang         dk   += uikdi*ba[ili];
12973c9e8598SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
12983c9e8598SHong Zhang 
12993c9e8598SHong Zhang         /* add multiple of row i to k-th row ... */
13003c9e8598SHong Zhang         jmin = ili + 1;
13013c9e8598SHong Zhang         nz   = bi[i+1] - jmin;
13023c9e8598SHong Zhang         if (nz > 0){
13033c9e8598SHong Zhang           bcol = bj + jmin;
13043c9e8598SHong Zhang           bval = ba + jmin;
13053c9e8598SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
13063c9e8598SHong Zhang           /* update il and jl for i-th row */
13073c9e8598SHong Zhang           il[i] = jmin;
13083c9e8598SHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
13093c9e8598SHong Zhang         }
13103c9e8598SHong Zhang         i = nexti;
13113c9e8598SHong Zhang       }
13123c9e8598SHong Zhang 
13133c9e8598SHong Zhang       if (PetscRealPart(dk) < zeropivot && b->factor_shift){
13143c9e8598SHong Zhang 	/* calculate a shift that would make this row diagonally dominant */
13153c9e8598SHong Zhang 	PetscReal rs = PetscAbs(PetscRealPart(dk));
13163c9e8598SHong Zhang 	jmin      = bi[k]+1;
13173c9e8598SHong Zhang 	nz        = bi[k+1] - jmin;
13183c9e8598SHong Zhang 	if (nz){
13193c9e8598SHong Zhang 	  bcol = bj + jmin;
13203c9e8598SHong Zhang 	  bval = ba + jmin;
13213c9e8598SHong Zhang 	  while (nz--){
13223c9e8598SHong Zhang 	    rs += PetscAbsScalar(rtmp[*bcol++]);
13233c9e8598SHong Zhang 	  }
13243c9e8598SHong Zhang 	}
13253c9e8598SHong Zhang 	/* if this shift is less than the previous, just up the previous
13263c9e8598SHong Zhang 	   one by a bit */
13273c9e8598SHong Zhang 	shift_amount = PetscMax(rs,1.1*shift_amount);
13283c9e8598SHong Zhang 	chshift  = PETSC_TRUE;
13293c9e8598SHong Zhang 	/* Unlike in the ILU case there is no exit condition on nshift:
13303c9e8598SHong Zhang 	   we increase the shift until it converges. There is no guarantee that
13313c9e8598SHong Zhang 	   this algorithm converges faster or slower, or is better or worse
13323c9e8598SHong Zhang 	   than the ILU algorithm. */
13333c9e8598SHong Zhang 	nshift++;
13343c9e8598SHong Zhang 	break;
13353c9e8598SHong Zhang       }
13363c9e8598SHong Zhang       if (PetscRealPart(dk) < zeropivot){
13373c9e8598SHong Zhang         if (damping == (PetscReal) PETSC_DECIDE) damping = -PetscRealPart(dk)/(k+1);
13383c9e8598SHong Zhang         if (damping > 0.0) {
13393c9e8598SHong Zhang           if (ndamp) damping *= 2.0;
13403c9e8598SHong Zhang           damp = PETSC_TRUE;
13413c9e8598SHong Zhang           ndamp++;
13423c9e8598SHong Zhang           break;
13433c9e8598SHong Zhang         } else if (PetscAbsScalar(dk) < zeropivot){
13443c9e8598SHong Zhang           SETERRQ3(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g",k,PetscRealPart(dk),zeropivot);
13453c9e8598SHong Zhang         } else {
13463c9e8598SHong Zhang           PetscLogInfo((PetscObject)A,"Negative pivot %g in row %D of Cholesky factorization\n",PetscRealPart(dk),k);
13473c9e8598SHong Zhang         }
13483c9e8598SHong Zhang       }
13493c9e8598SHong Zhang 
13503c9e8598SHong Zhang       /* copy data into U(k,:) */
13513c9e8598SHong Zhang       ba[bi[k]] = 1.0/dk;
13523c9e8598SHong Zhang       jmin      = bi[k]+1;
13533c9e8598SHong Zhang       nz        = bi[k+1] - jmin;
13543c9e8598SHong Zhang       if (nz){
13553c9e8598SHong Zhang         bcol = bj + jmin;
13563c9e8598SHong Zhang         bval = ba + jmin;
13573c9e8598SHong Zhang         while (nz--){
13583c9e8598SHong Zhang           *bval++       = rtmp[*bcol];
13593c9e8598SHong Zhang           rtmp[*bcol++] = 0.0;
13603c9e8598SHong Zhang         }
13613c9e8598SHong Zhang         /* add k-th row into il and jl */
13623c9e8598SHong Zhang         il[k] = jmin;
13633c9e8598SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
13643c9e8598SHong Zhang       }
13653c9e8598SHong Zhang     }
13663c9e8598SHong Zhang   } while (damp||chshift);
13673c9e8598SHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
13683c9e8598SHong Zhang 
13693c9e8598SHong Zhang   C->factor       = FACTOR_CHOLESKY;
13703c9e8598SHong Zhang   C->assembled    = PETSC_TRUE;
13713c9e8598SHong Zhang   C->preallocated = PETSC_TRUE;
13723c9e8598SHong Zhang   PetscLogFlops(C->m);
13733c9e8598SHong Zhang   if (ndamp) {
13743c9e8598SHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumerical_SeqAIJ_NaturalOrdering: number of damping tries %D damping value %g\n",ndamp,damping);
13753c9e8598SHong Zhang   }
13763c9e8598SHong Zhang   if (nshift) {
13773c9e8598SHong Zhang     PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering diagonal shifted %D shifts\n",nshift);
13783c9e8598SHong Zhang   }
13793c9e8598SHong Zhang   PetscFunctionReturn(0);
13803c9e8598SHong Zhang }
1381a6175056SHong Zhang 
13827a48dd6fSHong Zhang #include "petscbt.h"
13837a48dd6fSHong Zhang #include "src/mat/utils/freespace.h"
1384a6175056SHong Zhang #undef __FUNCT__
1385a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ"
1386dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1387a6175056SHong Zhang {
13880968510aSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1389ed59401aSHong Zhang   Mat_SeqSBAIJ   *b;
1390ed59401aSHong Zhang   Mat            B;
1391dfbe8321SBarry Smith   PetscErrorCode ierr;
1392653a6975SHong Zhang   PetscTruth     perm_identity;
1393ed59401aSHong Zhang   PetscInt       reallocs=0,*rip,i,*ai,*aj,am=A->m,*ui;
1394ed59401aSHong Zhang   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
13957a48dd6fSHong Zhang   PetscInt       nlnk,*lnk,*lnk_lvl,ncols,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
1396ed59401aSHong Zhang   PetscReal      fill=info->fill,levels=info->levels;
13977a48dd6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
13987a48dd6fSHong Zhang   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
13997a48dd6fSHong Zhang   PetscBT        lnkbt;
1400a6175056SHong Zhang 
1401a6175056SHong Zhang   PetscFunctionBegin;
1402653a6975SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1403653a6975SHong Zhang   if (!perm_identity){
1404e005ede5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Non-identity permutation is not supported yet");
1405653a6975SHong Zhang   }
14067a48dd6fSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
14077a48dd6fSHong Zhang   ai = a->i; aj = a->j;
14087a48dd6fSHong Zhang   ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr);
14097a48dd6fSHong Zhang 
1410ed59401aSHong Zhang   /* levels = 0: simply copies fill pattern */
1411ed59401aSHong Zhang   if (!levels) {
1412ed59401aSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1413ed59401aSHong Zhang     for (i=0; i<am; i++) {
1414ed59401aSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1415ed59401aSHong Zhang     }
1416ed59401aSHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1417ed59401aSHong Zhang     B = *fact;
1418ed59401aSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1419ed59401aSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
1420ed59401aSHong Zhang 
1421ed59401aSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
1422ed59401aSHong Zhang     uj = b->j;
1423ed59401aSHong Zhang     for (i=0; i<am; i++) {
1424ed59401aSHong Zhang       aj = a->j + a->diag[i];
1425ed59401aSHong Zhang       for (j=0; j<ui[i]; j++){
1426ed59401aSHong Zhang         *uj++ = *aj++;
1427ed59401aSHong Zhang       }
1428ed59401aSHong Zhang       b->ilen[i] = ui[i];
1429ed59401aSHong Zhang     }
1430ed59401aSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
1431ed59401aSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1432ed59401aSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1433ed59401aSHong Zhang 
1434ed59401aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14353c9e8598SHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)B,"MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering","dummyname",(FCNVOID)B->ops->choleskyfactornumeric);CHKERRQ(ierr);
14363c9e8598SHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
1437ed59401aSHong Zhang     PetscFunctionReturn(0);
1438ed59401aSHong Zhang   }
1439ed59401aSHong Zhang 
14407a48dd6fSHong Zhang   /* initialization */
14417a48dd6fSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
14427a48dd6fSHong Zhang   ui[0] = 0;
14437a48dd6fSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
14447a48dd6fSHong Zhang 
14457a48dd6fSHong Zhang   /* jl: linked list for storing indices of the pivot rows
14467a48dd6fSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
14477a48dd6fSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
14487a48dd6fSHong Zhang   il         = jl + am;
14497a48dd6fSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
14507a48dd6fSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
14517a48dd6fSHong Zhang   for (i=0; i<am; i++){
14527a48dd6fSHong Zhang     jl[i] = am; il[i] = 0;
14537a48dd6fSHong Zhang   }
14547a48dd6fSHong Zhang 
14557a48dd6fSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
14567a48dd6fSHong Zhang   nlnk = am + 1;
1457*2ed133dbSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
14587a48dd6fSHong Zhang 
14597a48dd6fSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
14607a48dd6fSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
14617a48dd6fSHong Zhang   current_space = free_space;
14627a48dd6fSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
14637a48dd6fSHong Zhang   current_space_lvl = free_space_lvl;
14647a48dd6fSHong Zhang 
14657a48dd6fSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
14667a48dd6fSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
14677a48dd6fSHong Zhang     nzk   = 0;
1468ed59401aSHong Zhang     ncols = ai[rip[k]+1] - a->diag[rip[k]]; /* changes for !perm_identity */
14697a48dd6fSHong Zhang     cols  = aj + a->diag[rip[k]];
14707a48dd6fSHong Zhang     for (j=0; j<ncols; j++) cols_lvl[j] = -1;  /* initialize level for nonzero entries */
1471*2ed133dbSHong Zhang     ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
14727a48dd6fSHong Zhang     nzk += nlnk;
14737a48dd6fSHong Zhang 
14747a48dd6fSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
14757a48dd6fSHong Zhang     prow = jl[k]; /* 1st pivot row */
14767a48dd6fSHong Zhang 
14777a48dd6fSHong Zhang     while (prow < k){
14787a48dd6fSHong Zhang       nextprow = jl[prow];
14797a48dd6fSHong Zhang 
14807a48dd6fSHong Zhang       /* merge prow into k-th row */
14817a48dd6fSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
14827a48dd6fSHong Zhang       jmax = ui[prow+1];
14837a48dd6fSHong Zhang       ncols = jmax-jmin;
1484ed59401aSHong Zhang       i     = jmin - ui[prow];
1485ed59401aSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1486ed59401aSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
1487*2ed133dbSHong Zhang       ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
14887a48dd6fSHong Zhang       nzk += nlnk;
14897a48dd6fSHong Zhang 
14907a48dd6fSHong Zhang       /* update il and jl for prow */
14917a48dd6fSHong Zhang       if (jmin < jmax){
14927a48dd6fSHong Zhang         il[prow] = jmin;
14937a48dd6fSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
14947a48dd6fSHong Zhang       }
14957a48dd6fSHong Zhang       prow = nextprow;
14967a48dd6fSHong Zhang     }
14977a48dd6fSHong Zhang 
14987a48dd6fSHong Zhang     /* if free space is not available, make more free space */
14997a48dd6fSHong Zhang     if (current_space->local_remaining<nzk) {
15007a48dd6fSHong Zhang       i = am - k + 1; /* num of unfactored rows */
15017a48dd6fSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
15027a48dd6fSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
15037a48dd6fSHong Zhang       ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
15047a48dd6fSHong Zhang       reallocs++;
15057a48dd6fSHong Zhang     }
15067a48dd6fSHong Zhang 
15077a48dd6fSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
1508*2ed133dbSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
15097a48dd6fSHong Zhang 
15107a48dd6fSHong Zhang     /* add the k-th row into il and jl */
15117a48dd6fSHong Zhang     if (nzk-1 > 0){
15127a48dd6fSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
15137a48dd6fSHong Zhang       jl[k] = jl[i]; jl[i] = k;
15147a48dd6fSHong Zhang       il[k] = ui[k] + 1;
15157a48dd6fSHong Zhang     }
15167a48dd6fSHong Zhang     uj_ptr[k]     = current_space->array;
15177a48dd6fSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
15187a48dd6fSHong Zhang 
15197a48dd6fSHong Zhang     current_space->array           += nzk;
15207a48dd6fSHong Zhang     current_space->local_used      += nzk;
15217a48dd6fSHong Zhang     current_space->local_remaining -= nzk;
15227a48dd6fSHong Zhang 
15237a48dd6fSHong Zhang     current_space_lvl->array           += nzk;
15247a48dd6fSHong Zhang     current_space_lvl->local_used      += nzk;
15257a48dd6fSHong Zhang     current_space_lvl->local_remaining -= nzk;
15267a48dd6fSHong Zhang 
15277a48dd6fSHong Zhang     ui[k+1] = ui[k] + nzk;
15287a48dd6fSHong Zhang   }
15297a48dd6fSHong Zhang 
15307a48dd6fSHong Zhang   if (ai[am] != 0) {
15317a48dd6fSHong Zhang     PetscReal af = ((PetscReal)ui[am])/((PetscReal)ai[am]);
15327a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
15337a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af);
15347a48dd6fSHong Zhang     PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
15357a48dd6fSHong Zhang   } else {
1536ed59401aSHong Zhang      PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Empty matrix.\n");
15377a48dd6fSHong Zhang   }
15387a48dd6fSHong Zhang 
15397a48dd6fSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
15407a48dd6fSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
15417a48dd6fSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
15427a48dd6fSHong Zhang 
15437a48dd6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
15447a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
15457a48dd6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1546*2ed133dbSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15477a48dd6fSHong Zhang   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
15487a48dd6fSHong Zhang 
15497a48dd6fSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
15507a48dd6fSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
1551ed59401aSHong Zhang   B = *fact;
1552ed59401aSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1553ed59401aSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
15547a48dd6fSHong Zhang 
1555ed59401aSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
15567a48dd6fSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
15577a48dd6fSHong Zhang   b->singlemalloc = PETSC_FALSE;
15587a48dd6fSHong Zhang   /* the next line frees the default space generated by the Create() */
15597a48dd6fSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
15607a48dd6fSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
15617a48dd6fSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
15627a48dd6fSHong Zhang   b->j    = uj;
15637a48dd6fSHong Zhang   b->i    = ui;
15647a48dd6fSHong Zhang   b->diag = 0;
15657a48dd6fSHong Zhang   b->ilen = 0;
15667a48dd6fSHong Zhang   b->imax = 0;
15677a48dd6fSHong Zhang   b->row  = perm;
15687a48dd6fSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
15697a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15707a48dd6fSHong Zhang   b->icol = perm;
15717a48dd6fSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
15727a48dd6fSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1573ed59401aSHong Zhang   PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
15747a48dd6fSHong Zhang   b->maxnz = b->nz = ui[am];
15757a48dd6fSHong Zhang 
1576ed59401aSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1577ed59401aSHong Zhang   B->info.factor_mallocs    = reallocs;
1578ed59401aSHong Zhang   B->info.fill_ratio_given  = fill;
15797a48dd6fSHong Zhang   if (ai[am] != 0) {
1580ed59401aSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
15817a48dd6fSHong Zhang   } else {
1582ed59401aSHong Zhang     B->info.fill_ratio_needed = 0.0;
15837a48dd6fSHong Zhang   }
15843c9e8598SHong Zhang 
1585ed59401aSHong Zhang   B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1586ed59401aSHong Zhang   B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
15873c9e8598SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject) B,"MatCholeskyFactorNumeric_NaturalOrdering","dummyname",(FCNVOID)B->ops->choleskyfactornumeric);CHKERRQ(ierr);
15883c9e8598SHong Zhang   B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
1589a6175056SHong Zhang   PetscFunctionReturn(0);
1590a6175056SHong Zhang }
1591d5d45c9bSBarry Smith 
1592f76d2b81SHong Zhang #undef __FUNCT__
1593f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ"
1594dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
1595f76d2b81SHong Zhang {
1596f76d2b81SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
159710c27e3dSHong Zhang   Mat_SeqSBAIJ   *b;
1598*2ed133dbSHong Zhang   Mat            B;
1599dfbe8321SBarry Smith   PetscErrorCode ierr;
1600f76d2b81SHong Zhang   PetscTruth     perm_identity;
160110c27e3dSHong Zhang   PetscReal      fill = info->fill;
1602*2ed133dbSHong Zhang   PetscInt       *rip,*riip,i,mbs=A->m,*ai=a->i,*aj=a->j,reallocs=0,prow;
160310c27e3dSHong Zhang   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1604*2ed133dbSHong Zhang   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
160510c27e3dSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
160610c27e3dSHong Zhang   PetscBT        lnkbt;
1607*2ed133dbSHong Zhang   IS             iperm;
1608f76d2b81SHong Zhang 
1609f76d2b81SHong Zhang   PetscFunctionBegin;
1610*2ed133dbSHong Zhang   /* check whether perm is the identity mapping */
1611f76d2b81SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1612*2ed133dbSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1613*2ed133dbSHong Zhang 
1614f76d2b81SHong Zhang   if (!perm_identity){
1615*2ed133dbSHong Zhang     /* check if perm is symmetric! */
1616*2ed133dbSHong Zhang     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
1617*2ed133dbSHong Zhang     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
1618*2ed133dbSHong Zhang     for (i=0; i<mbs; i++) {
1619*2ed133dbSHong Zhang       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
1620*2ed133dbSHong Zhang     }
1621*2ed133dbSHong Zhang     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
1622*2ed133dbSHong Zhang     ierr = ISDestroy(iperm);CHKERRQ(ierr);
1623f76d2b81SHong Zhang   }
162410c27e3dSHong Zhang 
162510c27e3dSHong Zhang   /* initialization */
1626*2ed133dbSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
162710c27e3dSHong Zhang   ui[0] = 0;
162810c27e3dSHong Zhang 
162910c27e3dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
1630*2ed133dbSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1631*2ed133dbSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
1632*2ed133dbSHong Zhang   il     = jl + mbs;
1633*2ed133dbSHong Zhang   cols   = il + mbs;
1634*2ed133dbSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
1635*2ed133dbSHong Zhang   for (i=0; i<mbs; i++){
1636*2ed133dbSHong Zhang     jl[i] = mbs; il[i] = 0;
1637f76d2b81SHong Zhang   }
163810c27e3dSHong Zhang 
163910c27e3dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
1640*2ed133dbSHong Zhang   nlnk = mbs + 1;
1641*2ed133dbSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
164210c27e3dSHong Zhang 
1643*2ed133dbSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1644*2ed133dbSHong Zhang   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
164510c27e3dSHong Zhang   current_space = free_space;
164610c27e3dSHong Zhang 
1647*2ed133dbSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
164810c27e3dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
164910c27e3dSHong Zhang     nzk   = 0;
1650*2ed133dbSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
1651*2ed133dbSHong Zhang     ncols_upper = 0;
1652*2ed133dbSHong Zhang     for (j=0; j<ncols; j++){
1653*2ed133dbSHong Zhang       i = *(aj + ai[rip[k]] + j);
1654*2ed133dbSHong Zhang       i = rip[i];
1655*2ed133dbSHong Zhang       if (i >= k){ /* only take upper triangular entry */
1656*2ed133dbSHong Zhang         cols[ncols_upper] = i;
1657*2ed133dbSHong Zhang         ncols_upper++;
1658*2ed133dbSHong Zhang       }
1659*2ed133dbSHong Zhang     }
1660*2ed133dbSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
166110c27e3dSHong Zhang     nzk += nlnk;
166210c27e3dSHong Zhang 
166310c27e3dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
166410c27e3dSHong Zhang     prow = jl[k]; /* 1st pivot row */
166510c27e3dSHong Zhang 
166610c27e3dSHong Zhang     while (prow < k){
166710c27e3dSHong Zhang       nextprow = jl[prow];
166810c27e3dSHong Zhang       /* merge prow into k-th row */
1669*2ed133dbSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
167010c27e3dSHong Zhang       jmax = ui[prow+1];
167110c27e3dSHong Zhang       ncols = jmax-jmin;
1672*2ed133dbSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
1673*2ed133dbSHong Zhang       ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
167410c27e3dSHong Zhang       nzk += nlnk;
167510c27e3dSHong Zhang 
167610c27e3dSHong Zhang       /* update il and jl for prow */
167710c27e3dSHong Zhang       if (jmin < jmax){
167810c27e3dSHong Zhang         il[prow] = jmin;
1679*2ed133dbSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
168010c27e3dSHong Zhang       }
168110c27e3dSHong Zhang       prow = nextprow;
168210c27e3dSHong Zhang     }
168310c27e3dSHong Zhang 
168410c27e3dSHong Zhang     /* if free space is not available, make more free space */
168510c27e3dSHong Zhang     if (current_space->local_remaining<nzk) {
1686*2ed133dbSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
168710c27e3dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
168810c27e3dSHong Zhang       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
168910c27e3dSHong Zhang       reallocs++;
169010c27e3dSHong Zhang     }
169110c27e3dSHong Zhang 
169210c27e3dSHong Zhang     /* copy data into free space, then initialize lnk */
1693*2ed133dbSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
169410c27e3dSHong Zhang 
169510c27e3dSHong Zhang     /* add the k-th row into il and jl */
169610c27e3dSHong Zhang     if (nzk-1 > 0){
1697*2ed133dbSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
169810c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
169910c27e3dSHong Zhang       il[k] = ui[k] + 1;
170010c27e3dSHong Zhang     }
1701*2ed133dbSHong Zhang     ui_ptr[k] = current_space->array;
170210c27e3dSHong Zhang     current_space->array           += nzk;
170310c27e3dSHong Zhang     current_space->local_used      += nzk;
170410c27e3dSHong Zhang     current_space->local_remaining -= nzk;
170510c27e3dSHong Zhang 
170610c27e3dSHong Zhang     ui[k+1] = ui[k] + nzk;
170710c27e3dSHong Zhang   }
170810c27e3dSHong Zhang 
1709*2ed133dbSHong Zhang   if (ai[mbs] != 0) {
1710*2ed133dbSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1711*2ed133dbSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af);
1712*2ed133dbSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af);
1713*2ed133dbSHong Zhang     PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af);
171410c27e3dSHong Zhang   } else {
1715*2ed133dbSHong Zhang      PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n");
171610c27e3dSHong Zhang   }
171710c27e3dSHong Zhang 
171810c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
171910c27e3dSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
172010c27e3dSHong Zhang 
172110c27e3dSHong Zhang   /* destroy list of free space and other temporary array(s) */
1722*2ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
172310c27e3dSHong Zhang   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
172410c27e3dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
172510c27e3dSHong Zhang 
172610c27e3dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1727*2ed133dbSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
1728*2ed133dbSHong Zhang   B    = *fact;
1729*2ed133dbSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
1730*2ed133dbSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
173110c27e3dSHong Zhang 
1732*2ed133dbSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
173310c27e3dSHong Zhang   ierr = PetscFree(b->imax);CHKERRQ(ierr);
173410c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
173510c27e3dSHong Zhang   /* the next line frees the default space generated by the Create() */
173610c27e3dSHong Zhang   ierr = PetscFree(b->a);CHKERRQ(ierr);
173710c27e3dSHong Zhang   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
1738*2ed133dbSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
173910c27e3dSHong Zhang   b->j    = uj;
174010c27e3dSHong Zhang   b->i    = ui;
174110c27e3dSHong Zhang   b->diag = 0;
174210c27e3dSHong Zhang   b->ilen = 0;
174310c27e3dSHong Zhang   b->imax = 0;
174410c27e3dSHong Zhang   b->row  = perm;
174510c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
174610c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
174710c27e3dSHong Zhang   b->icol = perm;
174810c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1749*2ed133dbSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1750*2ed133dbSHong Zhang   PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
1751*2ed133dbSHong Zhang   b->maxnz = b->nz = ui[mbs];
175210c27e3dSHong Zhang 
1753*2ed133dbSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
1754*2ed133dbSHong Zhang   B->info.factor_mallocs    = reallocs;
1755*2ed133dbSHong Zhang   B->info.fill_ratio_given  = fill;
1756*2ed133dbSHong Zhang   if (ai[mbs] != 0) {
1757*2ed133dbSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
175810c27e3dSHong Zhang   } else {
1759*2ed133dbSHong Zhang     B->info.fill_ratio_needed = 0.0;
176010c27e3dSHong Zhang   }
1761*2ed133dbSHong Zhang   if (perm_identity){
176210c27e3dSHong Zhang     (*fact)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
176310c27e3dSHong Zhang     (*fact)->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1764*2ed133dbSHong Zhang     ierr = PetscObjectComposeFunction((PetscObject)(*fact),"MatCholeskyFactorNumeric_NaturalOrdering","dummyname",(FCNVOID)(*fact)->ops->choleskyfactornumeric);CHKERRQ(ierr);
17653c9e8598SHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_NaturalOrdering;
1766*2ed133dbSHong Zhang   } else {
1767*2ed133dbSHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
1768*2ed133dbSHong Zhang   }
1769f76d2b81SHong Zhang   PetscFunctionReturn(0);
1770f76d2b81SHong Zhang }
1771