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 1101a6175056SHong Zhang #undef __FUNCT__ 1102a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1103dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact) 1104a6175056SHong Zhang { 11050968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 110679416396SBarry Smith PetscErrorCode ierr,(*f)(Mat,Mat*); 1107d5d45c9bSBarry Smith 1108a6175056SHong Zhang PetscFunctionBegin; 11090968510aSHong Zhang if (!a->sbaijMat){ 11100968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1111a6175056SHong Zhang } 1112d07a9264SSatish Balay ierr = PetscObjectQueryFunction((PetscObject) *fact,"MatCholeskyFactorNumeric",(void (**)(void))&f);CHKERRQ(ierr); 111379416396SBarry Smith ierr = (*f)(a->sbaijMat,fact);CHKERRQ(ierr); 11140968510aSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 1115064503c5SHong Zhang a->sbaijMat = PETSC_NULL; 1116a6175056SHong Zhang PetscFunctionReturn(0); 1117a6175056SHong Zhang } 1118a6175056SHong Zhang 1119*7a48dd6fSHong Zhang #include "petscbt.h" 1120*7a48dd6fSHong Zhang #include "src/mat/utils/freespace.h" 1121*7a48dd6fSHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1122a6175056SHong Zhang #undef __FUNCT__ 1123a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1124dfbe8321SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1125a6175056SHong Zhang { 11260968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1127dfbe8321SBarry Smith PetscErrorCode ierr; 1128653a6975SHong Zhang PetscTruth perm_identity; 1129*7a48dd6fSHong Zhang Mat_SeqSBAIJ *b; 1130*7a48dd6fSHong Zhang PetscInt *rip,i,*ai = a->i,*aj = a->j; 1131*7a48dd6fSHong Zhang PetscInt reallocs=0; 1132*7a48dd6fSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl; 1133*7a48dd6fSHong Zhang PetscInt prow; 1134*7a48dd6fSHong Zhang PetscInt *il,nextprow; 1135*7a48dd6fSHong Zhang PetscReal fill = info->fill,levels = info->levels; 1136*7a48dd6fSHong Zhang PetscInt am=A->m,*ui; 1137*7a48dd6fSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl,ncols,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 1138*7a48dd6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1139*7a48dd6fSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 1140*7a48dd6fSHong Zhang PetscBT lnkbt; 1141a6175056SHong Zhang 1142a6175056SHong Zhang PetscFunctionBegin; 1143653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1144653a6975SHong Zhang if (!perm_identity){ 1145e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Non-identity permutation is not supported yet"); 1146653a6975SHong Zhang } 1147*7a48dd6fSHong Zhang 1148*7a48dd6fSHong Zhang /* levels = 0: simply copies fill pattern */ 1149*7a48dd6fSHong Zhang if (!levels) { 1150*7a48dd6fSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,fact);CHKERRQ(ierr); /* don't need copy numerical values! */ 1151*7a48dd6fSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 1152*7a48dd6fSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1153*7a48dd6fSHong Zhang 1154*7a48dd6fSHong Zhang ierr = PetscObjectComposeFunction((PetscObject) *fact,"MatCholeskyFactorNumeric","dummyname",(FCNVOID)(*fact)->ops->choleskyfactornumeric);CHKERRQ(ierr); 1155*7a48dd6fSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1156*7a48dd6fSHong Zhang PetscFunctionReturn(0); 11570968510aSHong Zhang } 1158a6175056SHong Zhang 1159*7a48dd6fSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1160*7a48dd6fSHong Zhang ai = a->i; aj = a->j; 1161*7a48dd6fSHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 1162*7a48dd6fSHong Zhang 1163*7a48dd6fSHong Zhang /* initialization */ 1164*7a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1165*7a48dd6fSHong Zhang ui[0] = 0; 1166*7a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 1167*7a48dd6fSHong Zhang 1168*7a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 1169*7a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1170*7a48dd6fSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 1171*7a48dd6fSHong Zhang il = jl + am; 1172*7a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 1173*7a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 1174*7a48dd6fSHong Zhang for (i=0; i<am; i++){ 1175*7a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 1176*7a48dd6fSHong Zhang } 1177*7a48dd6fSHong Zhang 1178*7a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 1179*7a48dd6fSHong Zhang nlnk = am + 1; 1180*7a48dd6fSHong Zhang ierr = PetscLLCreate_PermutedLeveled(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1181*7a48dd6fSHong Zhang 1182*7a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1183*7a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 1184*7a48dd6fSHong Zhang current_space = free_space; 1185*7a48dd6fSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 1186*7a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 1187*7a48dd6fSHong Zhang 1188*7a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 1189*7a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 1190*7a48dd6fSHong Zhang nzk = 0; 1191*7a48dd6fSHong Zhang ncols = ai[rip[k]+1] - a->diag[rip[k]]; /* diag would change for !perm_identity */ 1192*7a48dd6fSHong Zhang cols = aj + a->diag[rip[k]]; 1193*7a48dd6fSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = -1; /* initialize level for nonzero entries */ 1194*7a48dd6fSHong Zhang ierr = PetscLLAdd_PermutedLeveled(ncols,cols,levels,cols_lvl,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1195*7a48dd6fSHong Zhang nzk += nlnk; 1196*7a48dd6fSHong Zhang 1197*7a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 1198*7a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 1199*7a48dd6fSHong Zhang 1200*7a48dd6fSHong Zhang while (prow < k){ 1201*7a48dd6fSHong Zhang nextprow = jl[prow]; 1202*7a48dd6fSHong Zhang 1203*7a48dd6fSHong Zhang /* merge prow into k-th row */ 1204*7a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 1205*7a48dd6fSHong Zhang jmax = ui[prow+1]; 1206*7a48dd6fSHong Zhang ncols = jmax-jmin; 1207*7a48dd6fSHong Zhang cols = uj_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 1208*7a48dd6fSHong Zhang /* ierr = PetscLLAdd(ncols,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); */ 1209*7a48dd6fSHong Zhang /* cols_lvl = uj_lvl_ptr[prow] + jmin - ui[prow]; */ 1210*7a48dd6fSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + jmin - ui[prow] + j); 1211*7a48dd6fSHong Zhang ierr = PetscLLAdd_PermutedLeveled(ncols,cols,levels,cols_lvl,am,rip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1212*7a48dd6fSHong Zhang nzk += nlnk; 1213*7a48dd6fSHong Zhang 1214*7a48dd6fSHong Zhang /* update il and jl for prow */ 1215*7a48dd6fSHong Zhang if (jmin < jmax){ 1216*7a48dd6fSHong Zhang il[prow] = jmin; 1217*7a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 1218*7a48dd6fSHong Zhang } 1219*7a48dd6fSHong Zhang prow = nextprow; 1220*7a48dd6fSHong Zhang } 1221*7a48dd6fSHong Zhang 1222*7a48dd6fSHong Zhang /* if free space is not available, make more free space */ 1223*7a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 1224*7a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 1225*7a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1226*7a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 1227*7a48dd6fSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 1228*7a48dd6fSHong Zhang reallocs++; 1229*7a48dd6fSHong Zhang } 1230*7a48dd6fSHong Zhang 1231*7a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 1232*7a48dd6fSHong Zhang ierr = PetscLLClean_PermutedLeveled(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1233*7a48dd6fSHong Zhang 1234*7a48dd6fSHong Zhang /* add the k-th row into il and jl */ 1235*7a48dd6fSHong Zhang if (nzk-1 > 0){ 1236*7a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 1237*7a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 1238*7a48dd6fSHong Zhang il[k] = ui[k] + 1; 1239*7a48dd6fSHong Zhang } 1240*7a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 1241*7a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 1242*7a48dd6fSHong Zhang 1243*7a48dd6fSHong Zhang current_space->array += nzk; 1244*7a48dd6fSHong Zhang current_space->local_used += nzk; 1245*7a48dd6fSHong Zhang current_space->local_remaining -= nzk; 1246*7a48dd6fSHong Zhang 1247*7a48dd6fSHong Zhang current_space_lvl->array += nzk; 1248*7a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 1249*7a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 1250*7a48dd6fSHong Zhang 1251*7a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 1252*7a48dd6fSHong Zhang } 1253*7a48dd6fSHong Zhang 1254*7a48dd6fSHong Zhang if (ai[am] != 0) { 1255*7a48dd6fSHong Zhang PetscReal af = ((PetscReal)ui[am])/((PetscReal)ai[am]); 1256*7a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 1257*7a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af); 1258*7a48dd6fSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 1259*7a48dd6fSHong Zhang } else { 1260*7a48dd6fSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"); 1261*7a48dd6fSHong Zhang } 1262*7a48dd6fSHong Zhang 1263*7a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1264*7a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 1265*7a48dd6fSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 1266*7a48dd6fSHong Zhang 1267*7a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 1268*7a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1269*7a48dd6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1270*7a48dd6fSHong Zhang ierr = PetscLLDestroy_PermutedLeveled(lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1271*7a48dd6fSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 1272*7a48dd6fSHong Zhang 1273*7a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1274*7a48dd6fSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 1275*7a48dd6fSHong Zhang ierr = MatSetType(*fact,MATSEQSBAIJ);CHKERRQ(ierr); 1276*7a48dd6fSHong Zhang ierr = MatSeqSBAIJSetPreallocation(*fact,1,0,PETSC_NULL);CHKERRQ(ierr); 1277*7a48dd6fSHong Zhang 1278*7a48dd6fSHong Zhang b = (Mat_SeqSBAIJ*)(*fact)->data; 1279*7a48dd6fSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 1280*7a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 1281*7a48dd6fSHong Zhang /* the next line frees the default space generated by the Create() */ 1282*7a48dd6fSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 1283*7a48dd6fSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1284*7a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1285*7a48dd6fSHong Zhang b->j = uj; 1286*7a48dd6fSHong Zhang b->i = ui; 1287*7a48dd6fSHong Zhang b->diag = 0; 1288*7a48dd6fSHong Zhang b->ilen = 0; 1289*7a48dd6fSHong Zhang b->imax = 0; 1290*7a48dd6fSHong Zhang b->row = perm; 1291*7a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1292*7a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1293*7a48dd6fSHong Zhang b->icol = perm; 1294*7a48dd6fSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1295*7a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1296*7a48dd6fSHong Zhang /* In b structure: Free imax, ilen, old a, old j. 1297*7a48dd6fSHong Zhang Allocate idnew, solve_work, new a, new j */ 1298*7a48dd6fSHong Zhang PetscLogObjectMemory(*fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))); 1299*7a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 1300*7a48dd6fSHong Zhang 1301*7a48dd6fSHong Zhang (*fact)->factor = FACTOR_CHOLESKY; 1302*7a48dd6fSHong Zhang (*fact)->info.factor_mallocs = reallocs; 1303*7a48dd6fSHong Zhang (*fact)->info.fill_ratio_given = fill; 1304*7a48dd6fSHong Zhang if (ai[am] != 0) { 1305*7a48dd6fSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 1306*7a48dd6fSHong Zhang } else { 1307*7a48dd6fSHong Zhang (*fact)->info.fill_ratio_needed = 0.0; 1308*7a48dd6fSHong Zhang } 1309*7a48dd6fSHong Zhang /* Should the followings be moved to MatCholeskyFactorNumeric_SeqAIJ()? */ 1310*7a48dd6fSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 1311*7a48dd6fSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1312*7a48dd6fSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1313*7a48dd6fSHong Zhang 1314*7a48dd6fSHong Zhang /* -------- end of new ---------------------------*/ 1315d07a9264SSatish Balay ierr = PetscObjectComposeFunction((PetscObject) *fact,"MatCholeskyFactorNumeric","dummyname",(FCNVOID)(*fact)->ops->choleskyfactornumeric);CHKERRQ(ierr); 1316653a6975SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1317653a6975SHong Zhang 1318a6175056SHong Zhang PetscFunctionReturn(0); 1319a6175056SHong Zhang } 1320d5d45c9bSBarry Smith 1321f76d2b81SHong Zhang #undef __FUNCT__ 1322f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1323dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1324f76d2b81SHong Zhang { 1325f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 132610c27e3dSHong Zhang Mat_SeqSBAIJ *b; 1327dfbe8321SBarry Smith PetscErrorCode ierr; 1328f76d2b81SHong Zhang PetscTruth perm_identity; 132910c27e3dSHong Zhang PetscReal fill = info->fill; 133010c27e3dSHong Zhang PetscInt *rip,i,am=A->m,*ai,*aj,reallocs=0,prow; 133110c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 133210c27e3dSHong Zhang PetscInt nlnk,*lnk,ncols,*cols,*uj,**uj_ptr; 133310c27e3dSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 133410c27e3dSHong Zhang PetscBT lnkbt; 1335f76d2b81SHong Zhang 1336f76d2b81SHong Zhang PetscFunctionBegin; 1337f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1338f76d2b81SHong Zhang if (!perm_identity){ 1339e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Non-identity permutation is not supported yet"); 1340f76d2b81SHong Zhang } 134110c27e3dSHong Zhang 134210c27e3dSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 134310c27e3dSHong Zhang ai = a->i; aj = a->j; 134410c27e3dSHong Zhang ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 134510c27e3dSHong Zhang 134610c27e3dSHong Zhang /* initialization */ 134710c27e3dSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 134810c27e3dSHong Zhang ui[0] = 0; 134910c27e3dSHong Zhang 135010c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 135110c27e3dSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 135210c27e3dSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 135310c27e3dSHong Zhang il = jl + am; 135410c27e3dSHong Zhang uj_ptr = (PetscInt**)(il + am); 135510c27e3dSHong Zhang for (i=0; i<am; i++){ 135610c27e3dSHong Zhang jl[i] = am; il[i] = 0; 1357f76d2b81SHong Zhang } 135810c27e3dSHong Zhang 135910c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 136010c27e3dSHong Zhang nlnk = am + 1; 136110c27e3dSHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 136210c27e3dSHong Zhang 136310c27e3dSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 136410c27e3dSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 136510c27e3dSHong Zhang current_space = free_space; 136610c27e3dSHong Zhang 136710c27e3dSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 136810c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 136910c27e3dSHong Zhang nzk = 0; 1370*7a48dd6fSHong Zhang ncols = ai[rip[k]+1] - a->diag[rip[k]]; /* diag would change for !perm_identity */ 137110c27e3dSHong Zhang cols = aj + a->diag[rip[k]]; 137210c27e3dSHong Zhang ierr = PetscLLAdd(ncols,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 137310c27e3dSHong Zhang nzk += nlnk; 137410c27e3dSHong Zhang 137510c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 137610c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 137710c27e3dSHong Zhang 137810c27e3dSHong Zhang while (prow < k){ 137910c27e3dSHong Zhang nextprow = jl[prow]; 138010c27e3dSHong Zhang 138110c27e3dSHong Zhang /* merge prow into k-th row */ 138210c27e3dSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 138310c27e3dSHong Zhang jmax = ui[prow+1]; 138410c27e3dSHong Zhang ncols = jmax-jmin; 138510c27e3dSHong Zhang cols = uj_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 138610c27e3dSHong Zhang ierr = PetscLLAdd(ncols,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 138710c27e3dSHong Zhang nzk += nlnk; 138810c27e3dSHong Zhang 138910c27e3dSHong Zhang /* update il and jl for prow */ 139010c27e3dSHong Zhang if (jmin < jmax){ 139110c27e3dSHong Zhang il[prow] = jmin; 139210c27e3dSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 139310c27e3dSHong Zhang } 139410c27e3dSHong Zhang prow = nextprow; 139510c27e3dSHong Zhang } 139610c27e3dSHong Zhang 139710c27e3dSHong Zhang /* if free space is not available, make more free space */ 139810c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 139910c27e3dSHong Zhang i = am - k + 1; /* num of unfactored rows */ 140010c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 140110c27e3dSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 140210c27e3dSHong Zhang reallocs++; 140310c27e3dSHong Zhang } 140410c27e3dSHong Zhang 140510c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 140610c27e3dSHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 140710c27e3dSHong Zhang 140810c27e3dSHong Zhang /* add the k-th row into il and jl */ 140910c27e3dSHong Zhang if (nzk-1 > 0){ 141010c27e3dSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 141110c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 141210c27e3dSHong Zhang il[k] = ui[k] + 1; 141310c27e3dSHong Zhang } 141410c27e3dSHong Zhang uj_ptr[k] = current_space->array; 141510c27e3dSHong Zhang current_space->array += nzk; 141610c27e3dSHong Zhang current_space->local_used += nzk; 141710c27e3dSHong Zhang current_space->local_remaining -= nzk; 141810c27e3dSHong Zhang 141910c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 142010c27e3dSHong Zhang } 142110c27e3dSHong Zhang 142210c27e3dSHong Zhang if (ai[am] != 0) { 142310c27e3dSHong Zhang PetscReal af = ((PetscReal)ui[am])/((PetscReal)ai[am]); 142410c27e3dSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 142510c27e3dSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Run with -pc_cholesky_fill %g or use \n",af); 142610c27e3dSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 142710c27e3dSHong Zhang } else { 142810c27e3dSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqAIJ:Empty matrix.\n"); 142910c27e3dSHong Zhang } 143010c27e3dSHong Zhang 143110c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 143210c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 143310c27e3dSHong Zhang 143410c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 143510c27e3dSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 143610c27e3dSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 143710c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 143810c27e3dSHong Zhang 143910c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 1440*7a48dd6fSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 144110c27e3dSHong Zhang ierr = MatSetType(*fact,MATSEQSBAIJ);CHKERRQ(ierr); 144210c27e3dSHong Zhang ierr = MatSeqSBAIJSetPreallocation(*fact,1,0,PETSC_NULL);CHKERRQ(ierr); 144310c27e3dSHong Zhang 144410c27e3dSHong Zhang b = (Mat_SeqSBAIJ*)(*fact)->data; 144510c27e3dSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 144610c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 144710c27e3dSHong Zhang /* the next line frees the default space generated by the Create() */ 144810c27e3dSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 144910c27e3dSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 145010c27e3dSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 145110c27e3dSHong Zhang b->j = uj; 145210c27e3dSHong Zhang b->i = ui; 145310c27e3dSHong Zhang b->diag = 0; 145410c27e3dSHong Zhang b->ilen = 0; 145510c27e3dSHong Zhang b->imax = 0; 145610c27e3dSHong Zhang b->row = perm; 145710c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 145810c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 145910c27e3dSHong Zhang b->icol = perm; 146010c27e3dSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 146110c27e3dSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 146210c27e3dSHong Zhang /* In b structure: Free imax, ilen, old a, old j. 146310c27e3dSHong Zhang Allocate idnew, solve_work, new a, new j */ 146410c27e3dSHong Zhang PetscLogObjectMemory(*fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))); 146510c27e3dSHong Zhang b->maxnz = b->nz = ui[am]; 146610c27e3dSHong Zhang 146710c27e3dSHong Zhang (*fact)->factor = FACTOR_CHOLESKY; 146810c27e3dSHong Zhang (*fact)->info.factor_mallocs = reallocs; 146910c27e3dSHong Zhang (*fact)->info.fill_ratio_given = fill; 147010c27e3dSHong Zhang if (ai[am] != 0) { 147110c27e3dSHong Zhang (*fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 147210c27e3dSHong Zhang } else { 147310c27e3dSHong Zhang (*fact)->info.fill_ratio_needed = 0.0; 147410c27e3dSHong Zhang } 1475*7a48dd6fSHong Zhang /* Should the followings be moved to MatCholeskyFactorNumeric_SeqAIJ()? */ 147610c27e3dSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 147710c27e3dSHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 147810c27e3dSHong Zhang (*fact)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 147910c27e3dSHong Zhang 148010c27e3dSHong Zhang /* -------- end of new ---------------------------*/ 1481d07a9264SSatish Balay ierr = PetscObjectComposeFunction((PetscObject) *fact,"MatCholeskyFactorNumeric","dummyname",(FCNVOID)(*fact)->ops->choleskyfactornumeric);CHKERRQ(ierr); 1482f76d2b81SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 148310c27e3dSHong Zhang 1484f76d2b81SHong Zhang PetscFunctionReturn(0); 1485f76d2b81SHong Zhang } 1486