173f4d377SMatthew Knepley /*$Id: aijfact.c,v 1.167 2001/09/11 16:32:26 bsmith Exp $*/ 2289bc588SBarry Smith 370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 490f02eecSBarry Smith #include "src/vec/vecimpl.h" 51c4feecaSSatish Balay #include "src/inline/dot.h" 6ed480e8bSBarry Smith #include "src/inline/spops.h" 73b2fbd54SBarry Smith 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 1091e9ee9fSBarry Smith int MatOrdering_Flow_SeqAIJ(Mat mat,MatOrderingType type,IS *irow,IS *icol) 113b2fbd54SBarry Smith { 123a40ed3dSBarry Smith PetscFunctionBegin; 133a40ed3dSBarry Smith 1429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 15aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 16d64ed03dSBarry Smith PetscFunctionReturn(0); 17d64ed03dSBarry Smith #endif 183b2fbd54SBarry Smith } 193b2fbd54SBarry Smith 2086bacbd2SBarry Smith 21ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqAIJ(Mat); 223a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth); 2386bacbd2SBarry Smith 2487828ca2SBarry Smith EXTERN int SPARSEKIT2dperm(int*,PetscScalar*,int*,int*,PetscScalar*,int*,int*,int*,int*,int*); 2587828ca2SBarry Smith EXTERN int SPARSEKIT2ilutp(int*,PetscScalar*,int*,int*,int*,PetscReal,PetscReal*,int*,PetscScalar*,int*,int*,int*,PetscScalar*,int*,int*,int*); 2687828ca2SBarry Smith EXTERN int SPARSEKIT2msrcsr(int*,PetscScalar*,int*,PetscScalar*,int*,int*,PetscScalar*,int*); 2786bacbd2SBarry Smith 284a2ae208SSatish Balay #undef __FUNCT__ 294a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 3086bacbd2SBarry Smith /* ------------------------------------------------------------ 3186bacbd2SBarry Smith 3286bacbd2SBarry Smith This interface was contribed by Tony Caola 3386bacbd2SBarry Smith 3486bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 355bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 365bd2ddc7SBarry Smith SPARSEKIT2. 375bd2ddc7SBarry Smith 385bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 395bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 4086bacbd2SBarry Smith 4186bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4286bacbd2SBarry Smith help in getting this routine ironed out. 4386bacbd2SBarry Smith 445bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 455bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 465bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 475bd2ddc7SBarry Smith Yousef Saad's code. 4886bacbd2SBarry Smith 495bd2ddc7SBarry Smith ishift = 0, for indices start at 1 505bd2ddc7SBarry Smith ishift = 1, for indices starting at 0 5186bacbd2SBarry Smith ------------------------------------------------------------ 5286bacbd2SBarry Smith */ 53b380c88cSHong Zhang int MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact) 5486bacbd2SBarry Smith { 5560ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5660ec86bdSBarry Smith PetscFunctionBegin; 5760ec86bdSBarry Smith SETERRQ(1,"This distribution does not include GNU Copyright code\n\ 5860ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5960ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 6060ec86bdSBarry Smith #else 6186bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 6207d2056aSBarry Smith IS iscolf,isicol,isirow; 6386bacbd2SBarry Smith PetscTruth reorder; 64273d9f13SBarry Smith int *c,*r,*ic,ierr,i,n = A->m; 65329f5518SBarry Smith int *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 66313ae024SBarry Smith int *ordcol,*iwk,*iperm,*jw; 675bd2ddc7SBarry Smith int ishift = !a->indexshift; 68b19713cbSBarry Smith int jmax,lfill,job,*o_i,*o_j; 6987828ca2SBarry Smith PetscScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 70329f5518SBarry Smith PetscReal permtol,af; 7186bacbd2SBarry Smith 7286bacbd2SBarry Smith PetscFunctionBegin; 7386bacbd2SBarry Smith 7486bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7586bacbd2SBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax); 7686bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7733259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 786faa4630SBarry Smith lfill = (int)(info->dtcount/2.0); 796faa4630SBarry Smith jmax = (int)(info->fill*a->nz); 8086bacbd2SBarry Smith permtol = info->dtcol; 8186bacbd2SBarry Smith 8286bacbd2SBarry Smith 8386bacbd2SBarry Smith /* ------------------------------------------------------------ 8486bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8586bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8686bacbd2SBarry Smith then de-reordered so it is in it's original format. 8786bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8886bacbd2SBarry Smith the original matrix and allocate more storage. . . 8986bacbd2SBarry Smith ------------------------------------------------------------ 9086bacbd2SBarry Smith */ 9186bacbd2SBarry Smith 9286bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 9386bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9486bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9586bacbd2SBarry Smith reorder = PetscNot(reorder); 9686bacbd2SBarry Smith 9786bacbd2SBarry Smith 9886bacbd2SBarry Smith /* storage for ilu factor */ 99b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr); 100b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr); 10187828ca2SBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr); 102b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr); 10386bacbd2SBarry Smith 10486bacbd2SBarry Smith /* ------------------------------------------------------------ 10586bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10686bacbd2SBarry Smith ------------------------------------------------------------ 10786bacbd2SBarry Smith */ 108b19713cbSBarry Smith if (ishift) { 109b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 110b19713cbSBarry Smith old_j[i]++; 11186bacbd2SBarry Smith } 112b19713cbSBarry Smith for(i=0;i<n+1;i++) { 113b19713cbSBarry Smith old_i[i]++; 114b19713cbSBarry Smith }; 11586bacbd2SBarry Smith }; 11686bacbd2SBarry Smith 11786bacbd2SBarry Smith if (reorder) { 118c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c); CHKERRQ(ierr); 119c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r); CHKERRQ(ierr); 120c0c2c81eSBarry Smith for(i=0;i<n;i++) { 121c0c2c81eSBarry Smith r[i] = r[i]+1; 122c0c2c81eSBarry Smith c[i] = c[i]+1; 123c0c2c81eSBarry Smith } 124b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr); 125b0a32e0cSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr); 12687828ca2SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr); 1275bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 128c0c2c81eSBarry Smith for (i=0;i<n;i++) { 129c0c2c81eSBarry Smith r[i] = r[i]-1; 130c0c2c81eSBarry Smith c[i] = c[i]-1; 131c0c2c81eSBarry Smith } 132c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 133c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 134b19713cbSBarry Smith o_a = old_a2; 135b19713cbSBarry Smith o_j = old_j2; 136b19713cbSBarry Smith o_i = old_i2; 137b19713cbSBarry Smith } else { 138b19713cbSBarry Smith o_a = old_a; 139b19713cbSBarry Smith o_j = old_j; 140b19713cbSBarry Smith o_i = old_i; 14186bacbd2SBarry Smith } 14286bacbd2SBarry Smith 14386bacbd2SBarry Smith /* ------------------------------------------------------------ 14486bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14586bacbd2SBarry Smith ------------------------------------------------------------ 14686bacbd2SBarry Smith */ 14786bacbd2SBarry Smith 148b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr); 149b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr); 15087828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 15186bacbd2SBarry Smith 15287828ca2SBarry Smith SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,(PetscReal)info->dt,&permtol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&ierr); 15386bacbd2SBarry Smith if (ierr) { 15486bacbd2SBarry Smith switch (ierr) { 15529bbc08cSBarry Smith case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15629bbc08cSBarry Smith case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15729bbc08cSBarry Smith case -5: SETERRQ(1,"ilutp(), zero row encountered"); 15829bbc08cSBarry Smith case -1: SETERRQ(1,"ilutp(), input matrix may be wrong"); 15929bbc08cSBarry Smith case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax); 16029bbc08cSBarry Smith default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr); 16186bacbd2SBarry Smith } 16286bacbd2SBarry Smith } 16386bacbd2SBarry Smith 16486bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16586bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16686bacbd2SBarry Smith 16786bacbd2SBarry Smith /* ------------------------------------------------------------ 16886bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16986bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 17086bacbd2SBarry Smith ------------------------------------------------------------ 17186bacbd2SBarry Smith */ 17286bacbd2SBarry Smith 17387828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 174b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr); 17586bacbd2SBarry Smith 1765bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17786bacbd2SBarry Smith 17886bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17986bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 18086bacbd2SBarry Smith 18186bacbd2SBarry Smith if (reorder) { 18286bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 18386bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 18486bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 185313ae024SBarry Smith } else { 186b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 187f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 188b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 189b19713cbSBarry Smith } 190b19713cbSBarry Smith } 191b19713cbSBarry Smith 192b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 193b19713cbSBarry Smith if (ishift) { 19486bacbd2SBarry Smith for (i=0; i<n+1; i++) { 195b19713cbSBarry Smith old_i[i]--; 196b19713cbSBarry Smith } 197b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 198b19713cbSBarry Smith old_j[i]--; 199b19713cbSBarry Smith } 20086bacbd2SBarry Smith } 20186bacbd2SBarry Smith 202b19713cbSBarry Smith /* Make the factored matrix 0-based */ 203b19713cbSBarry Smith if (ishift) { 20486bacbd2SBarry Smith for (i=0; i<n+1; i++) { 205b19713cbSBarry Smith new_i[i]--; 206b19713cbSBarry Smith } 207b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 208b19713cbSBarry Smith new_j[i]--; 209b19713cbSBarry Smith } 21086bacbd2SBarry Smith } 21186bacbd2SBarry Smith 21286bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 21386bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2144c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2154c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 216c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 21786bacbd2SBarry Smith for(i=0; i<n; i++) { 21886bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21986bacbd2SBarry Smith }; 22086bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 22107d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 22286bacbd2SBarry Smith 223c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 224c0c2c81eSBarry Smith 225329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 22607d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 22786bacbd2SBarry Smith 22886bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22986bacbd2SBarry Smith 23086bacbd2SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr); 23186bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 23286bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 23386bacbd2SBarry Smith 23486bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 23586bacbd2SBarry Smith ierr = PetscFree(b->imax);CHKERRQ(ierr); 23686bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23707d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 238b19713cbSBarry Smith /* the next line frees the default space generated by the MatCreate() */ 23986bacbd2SBarry Smith ierr = PetscFree(b->a);CHKERRQ(ierr); 24086bacbd2SBarry Smith ierr = PetscFree(b->ilen);CHKERRQ(ierr); 24186bacbd2SBarry Smith b->a = new_a; 24286bacbd2SBarry Smith b->j = new_j; 24386bacbd2SBarry Smith b->i = new_i; 24486bacbd2SBarry Smith b->ilen = 0; 24586bacbd2SBarry Smith b->imax = 0; 2461f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 247313ae024SBarry Smith b->row = isirow; 24886bacbd2SBarry Smith b->col = iscolf; 24987828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 25086bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 25186bacbd2SBarry Smith b->indexshift = a->indexshift; 25286bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 25386bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 25486bacbd2SBarry Smith 25586bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 25686bacbd2SBarry Smith 25786bacbd2SBarry Smith /* check out for identical nodes. If found, use inode functions */ 2583a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 25986bacbd2SBarry Smith 260431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 261b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af); 262b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af); 263b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af); 264b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n"); 265431e710aSBarry Smith 26686bacbd2SBarry Smith PetscFunctionReturn(0); 26760ec86bdSBarry Smith #endif 26886bacbd2SBarry Smith } 26986bacbd2SBarry Smith 270289bc588SBarry Smith /* 271289bc588SBarry Smith Factorization code for AIJ format. 272289bc588SBarry Smith */ 2734a2ae208SSatish Balay #undef __FUNCT__ 274b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 275b380c88cSHong Zhang int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B) 276289bc588SBarry Smith { 277416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 278289bc588SBarry Smith IS isicol; 279273d9f13SBarry Smith int *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j; 280416022c9SBarry Smith int *ainew,*ajnew,jmax,*fill,*ajtmp,nz,shift = a->indexshift; 28191bf9adeSBarry Smith int *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im; 2829e046878SBarry Smith PetscReal f; 283289bc588SBarry Smith 2843a40ed3dSBarry Smith PetscFunctionBegin; 28529bbc08cSBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2864c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 287ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 288ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 289289bc588SBarry Smith 290289bc588SBarry Smith /* get new row pointers */ 291b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 292dbb450caSBarry Smith ainew[0] = -shift; 293289bc588SBarry Smith /* don't know how many column pointers are needed so estimate */ 294d3d32019SBarry Smith f = info->fill; 295dbb450caSBarry Smith jmax = (int)(f*ai[n]+(!shift)); 296b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 297289bc588SBarry Smith /* fill is a linked list of nonzeros in active row */ 298b0a32e0cSBarry Smith ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr); 2992fb3ed76SBarry Smith im = fill + n + 1; 300289bc588SBarry Smith /* idnew is location of diagonal in factor */ 301b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr); 302dbb450caSBarry Smith idnew[0] = -shift; 303289bc588SBarry Smith 304289bc588SBarry Smith for (i=0; i<n; i++) { 305289bc588SBarry Smith /* first copy previous fill into linked list */ 306289bc588SBarry Smith nnz = nz = ai[r[i]+1] - ai[r[i]]; 30729bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 308dbb450caSBarry Smith ajtmp = aj + ai[r[i]] + shift; 309289bc588SBarry Smith fill[n] = n; 310289bc588SBarry Smith while (nz--) { 311289bc588SBarry Smith fm = n; 312dbb450caSBarry Smith idx = ic[*ajtmp++ + shift]; 313289bc588SBarry Smith do { 314289bc588SBarry Smith m = fm; 315289bc588SBarry Smith fm = fill[m]; 316289bc588SBarry Smith } while (fm < idx); 317289bc588SBarry Smith fill[m] = idx; 318289bc588SBarry Smith fill[idx] = fm; 319289bc588SBarry Smith } 320289bc588SBarry Smith row = fill[n]; 321289bc588SBarry Smith while (row < i) { 322dbb450caSBarry Smith ajtmp = ajnew + idnew[row] + (!shift); 3232fb3ed76SBarry Smith nzbd = 1 + idnew[row] - ainew[row]; 3242fb3ed76SBarry Smith nz = im[row] - nzbd; 325289bc588SBarry Smith fm = row; 3262fb3ed76SBarry Smith while (nz-- > 0) { 327dbb450caSBarry Smith idx = *ajtmp++ + shift; 3282fb3ed76SBarry Smith nzbd++; 3292fb3ed76SBarry Smith if (idx == i) im[row] = nzbd; 330289bc588SBarry Smith do { 331289bc588SBarry Smith m = fm; 332289bc588SBarry Smith fm = fill[m]; 333289bc588SBarry Smith } while (fm < idx); 334289bc588SBarry Smith if (fm != idx) { 335289bc588SBarry Smith fill[m] = idx; 336289bc588SBarry Smith fill[idx] = fm; 337289bc588SBarry Smith fm = idx; 338289bc588SBarry Smith nnz++; 339289bc588SBarry Smith } 340289bc588SBarry Smith } 341289bc588SBarry Smith row = fill[row]; 342289bc588SBarry Smith } 343289bc588SBarry Smith /* copy new filled row into permanent storage */ 344289bc588SBarry Smith ainew[i+1] = ainew[i] + nnz; 345496e697eSBarry Smith if (ainew[i+1] > jmax) { 346ecf371e4SBarry Smith 347ecf371e4SBarry Smith /* estimate how much additional space we will need */ 348ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 349ecf371e4SBarry Smith /* just double the memory each time */ 350ecf371e4SBarry Smith int maxadd = jmax; 351ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */ 352bbb6d6a8SBarry Smith if (maxadd < nnz) maxadd = (n-i)*(nnz+1); 3532fb3ed76SBarry Smith jmax += maxadd; 354ecf371e4SBarry Smith 355ecf371e4SBarry Smith /* allocate a longer ajnew */ 356b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr); 357549d3d68SSatish Balay ierr = PetscMemcpy(ajtmp,ajnew,(ainew[i]+shift)*sizeof(int));CHKERRQ(ierr); 358606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 359289bc588SBarry Smith ajnew = ajtmp; 3602fb3ed76SBarry Smith realloc++; /* count how many times we realloc */ 361289bc588SBarry Smith } 362dbb450caSBarry Smith ajtmp = ajnew + ainew[i] + shift; 363289bc588SBarry Smith fm = fill[n]; 364289bc588SBarry Smith nzi = 0; 3652fb3ed76SBarry Smith im[i] = nnz; 366289bc588SBarry Smith while (nnz--) { 367289bc588SBarry Smith if (fm < i) nzi++; 368dbb450caSBarry Smith *ajtmp++ = fm - shift; 369289bc588SBarry Smith fm = fill[fm]; 370289bc588SBarry Smith } 371289bc588SBarry Smith idnew[i] = ainew[i] + nzi; 372289bc588SBarry Smith } 373464e8d44SSatish Balay if (ai[n] != 0) { 374329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 375b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 376b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af); 377b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af); 378b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"); 37905bf559cSSatish Balay } else { 380b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"); 38105bf559cSSatish Balay } 3822fb3ed76SBarry Smith 383898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 384898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3851987afe7SBarry Smith 386606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 387289bc588SBarry Smith 388289bc588SBarry Smith /* put together the new matrix */ 389b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr); 390b0a32e0cSBarry Smith PetscLogObjectParent(*B,isicol); 391416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 392606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3937c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 394e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 395606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 396606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 39787828ca2SBarry Smith ierr = PetscMalloc((ainew[n]+shift+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 398416022c9SBarry Smith b->j = ajnew; 399416022c9SBarry Smith b->i = ainew; 400416022c9SBarry Smith b->diag = idnew; 401416022c9SBarry Smith b->ilen = 0; 402416022c9SBarry Smith b->imax = 0; 403416022c9SBarry Smith b->row = isrow; 404416022c9SBarry Smith b->col = iscol; 405085256b3SBarry Smith b->lu_damping = info->damping; 40687828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 4076cc28720Svictorle b->lu_shift = info->lu_shift; 4086cc28720Svictorle b->lu_shift_fraction = info->lu_shift_fraction; 409c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 410c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 41182bf6240SBarry Smith b->icol = isicol; 41287828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 413416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 4147fd98bd6SLois Curfman McInnes Allocate idnew, solve_work, new a, new j */ 41587828ca2SBarry Smith PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(PetscScalar))); 416416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 4177fd98bd6SLois Curfman McInnes 418329f5518SBarry Smith (*B)->factor = FACTOR_LU; 419ae068f58SLois Curfman McInnes (*B)->info.factor_mallocs = realloc; 420ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 4213a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr); 4227dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 423703deb49SSatish Balay 424e93ccc4dSBarry Smith if (ai[n] != 0) { 425329f5518SBarry Smith (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 426eea2913aSSatish Balay } else { 427eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 428eea2913aSSatish Balay } 4293a40ed3dSBarry Smith PetscFunctionReturn(0); 430289bc588SBarry Smith } 4310a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4323a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth); 43341c01911SSatish Balay 4344a2ae208SSatish Balay #undef __FUNCT__ 4354a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 436416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B) 437289bc588SBarry Smith { 438416022c9SBarry Smith Mat C = *B; 439416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data; 44082bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 441273d9f13SBarry Smith int *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j; 4421987afe7SBarry Smith int *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift; 4430cf777b8SBarry Smith int *diag_offset = b->diag,diag,*pj,ndamp = 0, nshift=0; 44487828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 4450cf777b8SBarry Smith PetscReal damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs,d; 4460cf777b8SBarry Smith PetscReal row_shift,shift_fraction,shift_amount,shift_lo=0., shift_hi=1., shift_top=0.; 4470cf777b8SBarry Smith PetscTruth damp,lushift; 448289bc588SBarry Smith 4493a40ed3dSBarry Smith PetscFunctionBegin; 45078b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 45178b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 45287828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 45387828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 4540cf60383SBarry Smith rtmps = rtmp + shift; ics = ic + shift; 455289bc588SBarry Smith 4566cc28720Svictorle if (!a->diag) { 4570cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A); CHKERRQ(ierr); 4580cf777b8SBarry Smith } 4590cf777b8SBarry Smith 4600cf777b8SBarry Smith if (b->lu_shift) { /* set max shift */ 4610cf777b8SBarry Smith int *aai = a->i,*ddiag = a->diag; 4620cf777b8SBarry Smith shift_top = 0; 4636cc28720Svictorle for (i=0; i<n; i++) { 4640cf777b8SBarry Smith d = PetscAbsScalar((a->a)[ddiag[i]+shift]); 4656cc28720Svictorle /* calculate amt of shift needed for this row */ 466*28bb5da7SSatish Balay if (d<0) row_shift = 0; else row_shift = -2*d; 4678a2e2d9cSvictorle v = a->a+aai[i]+shift; 4688a2e2d9cSvictorle for (j=0; j<aai[i+1]-aai[i]; j++) 4696cc28720Svictorle row_shift += PetscAbsScalar(v[j]); 4706cc28720Svictorle if (row_shift>shift_top) shift_top = row_shift; 4716cc28720Svictorle } 4726cc28720Svictorle } 4736cc28720Svictorle 4746cc28720Svictorle shift_fraction = 0; shift_amount = 0; 475085256b3SBarry Smith do { 4768a5eb818SBarry Smith damp = PETSC_FALSE; 4776cc28720Svictorle lushift = PETSC_FALSE; 478289bc588SBarry Smith for (i=0; i<n; i++) { 479289bc588SBarry Smith nz = ai[i+1] - ai[i]; 480dbb450caSBarry Smith ajtmp = aj + ai[i] + shift; 48148e94559SBarry Smith for (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0; 482289bc588SBarry Smith 483289bc588SBarry Smith /* load in initial (unfactored row) */ 484416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 485416022c9SBarry Smith ajtmpold = a->j + a->i[r[i]] + shift; 486416022c9SBarry Smith v = a->a + a->i[r[i]] + shift; 487085256b3SBarry Smith for (j=0; j<nz; j++) { 488085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] = v[j]; 4896cc28720Svictorle if ( (ndamp||nshift) && ajtmpold[j] == r[i]) { /* damp the diagonal of the matrix */ 4906cc28720Svictorle rtmp[ics[ajtmpold[j]]] += damping + shift_amount; 491085256b3SBarry Smith } 492085256b3SBarry Smith } 493289bc588SBarry Smith 494dbb450caSBarry Smith row = *ajtmp++ + shift; 495f2582111SSatish Balay while (row < i) { 4968c37ef55SBarry Smith pc = rtmp + row; 4978c37ef55SBarry Smith if (*pc != 0.0) { 4981987afe7SBarry Smith pv = b->a + diag_offset[row] + shift; 4991987afe7SBarry Smith pj = b->j + diag_offset[row] + (!shift); 50035aab85fSBarry Smith multiplier = *pc / *pv++; 5018c37ef55SBarry Smith *pc = multiplier; 502f2582111SSatish Balay nz = ai[row+1] - diag_offset[row] - 1; 503f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 504b0a32e0cSBarry Smith PetscLogFlops(2*nz); 505289bc588SBarry Smith } 506dbb450caSBarry Smith row = *ajtmp++ + shift; 507289bc588SBarry Smith } 508416022c9SBarry Smith /* finished row so stick it into b->a */ 509416022c9SBarry Smith pv = b->a + ai[i] + shift; 510416022c9SBarry Smith pj = b->j + ai[i] + shift; 5118c37ef55SBarry Smith nz = ai[i+1] - ai[i]; 51235aab85fSBarry Smith diag = diag_offset[i] - ai[i]; 513d3d32019SBarry Smith /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 514d3d32019SBarry Smith rs = 0.0; 515d3d32019SBarry Smith for (j=0; j<nz; j++) { 516d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 517d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 518d3d32019SBarry Smith } 5196cc28720Svictorle #define MAX_NSHIFT 5 5208a2e2d9cSvictorle if (PetscRealPart(pv[diag]) < zeropivot*rs && b->lu_shift) { 5216cc28720Svictorle if (nshift>MAX_NSHIFT) { 5220cf777b8SBarry Smith SETERRQ(1,"Unable to determine shift to enforce positive definite preconditioner"); 5236cc28720Svictorle } else if (nshift==MAX_NSHIFT) { 5246cc28720Svictorle shift_fraction = shift_hi; 5256cc28720Svictorle } else { 5266cc28720Svictorle shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.; 5276cc28720Svictorle } 5286cc28720Svictorle shift_amount = shift_fraction * shift_top; 5290cf777b8SBarry Smith lushift = PETSC_TRUE; 5300cf777b8SBarry Smith nshift++; 5310cf777b8SBarry Smith break; 5326cc28720Svictorle } 533d3d32019SBarry Smith if (PetscAbsScalar(pv[diag]) < zeropivot*rs) { 534d3d32019SBarry Smith if (damping) { 5358a5eb818SBarry Smith if (ndamp) damping *= 2.0; 536085256b3SBarry Smith damp = PETSC_TRUE; 537085256b3SBarry Smith ndamp++; 538085256b3SBarry Smith break; 539085256b3SBarry Smith } else { 54063bb2aa6SBarry Smith SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs); 541d708749eSBarry Smith } 54235aab85fSBarry Smith } 54335aab85fSBarry Smith } 5446cc28720Svictorle if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) { 5456cc28720Svictorle /* 5466cc28720Svictorle * if not already shifting up & shifting & started shifting & can refine, 5476cc28720Svictorle * then try lower shift 5486cc28720Svictorle */ 5490cf777b8SBarry Smith shift_hi = shift_fraction; 5500cf777b8SBarry Smith shift_fraction = (shift_hi+shift_lo)/2.; 5516cc28720Svictorle shift_amount = shift_fraction * shift_top; 5520cf777b8SBarry Smith lushift = PETSC_TRUE; 5530cf777b8SBarry Smith nshift++; 5546cc28720Svictorle } 5556cc28720Svictorle } while (damp || lushift); 5560f11f9aeSSatish Balay 55735aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 55835aab85fSBarry Smith for (i=0; i<n; i++) { 55935aab85fSBarry Smith b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift]; 56035aab85fSBarry Smith } 56135aab85fSBarry Smith 562606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 56378b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 56478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 565416022c9SBarry Smith C->factor = FACTOR_LU; 5667dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 567c456f294SBarry Smith C->assembled = PETSC_TRUE; 568b0a32e0cSBarry Smith PetscLogFlops(C->n); 569d3d32019SBarry Smith if (ndamp) { 570b0a32e0cSBarry Smith PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping); 571085256b3SBarry Smith } 5726cc28720Svictorle if (nshift) { 5736cc28720Svictorle b->lu_shift_fraction = shift_fraction; 5746cc28720Svictorle PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction\n",shift_fraction); 5756cc28720Svictorle } 5763a40ed3dSBarry Smith PetscFunctionReturn(0); 577289bc588SBarry Smith } 5780a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5794a2ae208SSatish Balay #undef __FUNCT__ 5804a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 581b380c88cSHong Zhang int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 582da3a660dSBarry Smith { 583273d9f13SBarry Smith int ierr; 584416022c9SBarry Smith Mat C; 585416022c9SBarry Smith 5863a40ed3dSBarry Smith PetscFunctionBegin; 5879e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 588f2582111SSatish Balay ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr); 589273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 590440f4c53SSatish Balay PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol); 5913a40ed3dSBarry Smith PetscFunctionReturn(0); 592da3a660dSBarry Smith } 5930a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5944a2ae208SSatish Balay #undef __FUNCT__ 5954a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 596416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5978c37ef55SBarry Smith { 598416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 599416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 600273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 6013b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 60287828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 6038c37ef55SBarry Smith 6043a40ed3dSBarry Smith PetscFunctionBegin; 6053a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 60691bf9adeSBarry Smith 607e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 608e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 609416022c9SBarry Smith tmp = a->solve_work; 6108c37ef55SBarry Smith 6113b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6123b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6138c37ef55SBarry Smith 6148c37ef55SBarry Smith /* forward solve the lower triangular */ 6158c37ef55SBarry Smith tmp[0] = b[*r++]; 6160cf60383SBarry Smith tmps = tmp + shift; 6178c37ef55SBarry Smith for (i=1; i<n; i++) { 618dbb450caSBarry Smith v = aa + ai[i] + shift; 619dbb450caSBarry Smith vi = aj + ai[i] + shift; 62053e7425aSBarry Smith nz = a->diag[i] - ai[i]; 62153e7425aSBarry Smith sum = b[*r++]; 622ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6238c37ef55SBarry Smith tmp[i] = sum; 6248c37ef55SBarry Smith } 6258c37ef55SBarry Smith 6268c37ef55SBarry Smith /* backward solve the upper triangular */ 6278c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 628416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 629416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 630416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6318c37ef55SBarry Smith sum = tmp[i]; 632ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 633416022c9SBarry Smith x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift]; 6348c37ef55SBarry Smith } 6358c37ef55SBarry Smith 6363b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6373b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 638cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 639cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 640b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 6413a40ed3dSBarry Smith PetscFunctionReturn(0); 6428c37ef55SBarry Smith } 643026e39d0SSatish Balay 644930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6454a2ae208SSatish Balay #undef __FUNCT__ 6464a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 647930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 648930ae53cSBarry Smith { 649930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 650273d9f13SBarry Smith int n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr; 651362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 652aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 653d85afc42SSatish Balay int adiag_i,i,*vi,nz,ai_i; 654362ced78SSatish Balay PetscScalar *v,sum; 655d85afc42SSatish Balay #endif 656930ae53cSBarry Smith 6573a40ed3dSBarry Smith PetscFunctionBegin; 6583a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 659930ae53cSBarry Smith if (a->indexshift) { 6603a40ed3dSBarry Smith ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr); 6613a40ed3dSBarry Smith PetscFunctionReturn(0); 662930ae53cSBarry Smith } 663930ae53cSBarry Smith 664e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 665e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 666930ae53cSBarry Smith 667aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6681c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6691c4feecaSSatish Balay #else 670930ae53cSBarry Smith /* forward solve the lower triangular */ 671930ae53cSBarry Smith x[0] = b[0]; 672930ae53cSBarry Smith for (i=1; i<n; i++) { 673930ae53cSBarry Smith ai_i = ai[i]; 674930ae53cSBarry Smith v = aa + ai_i; 675930ae53cSBarry Smith vi = aj + ai_i; 676930ae53cSBarry Smith nz = adiag[i] - ai_i; 677930ae53cSBarry Smith sum = b[i]; 678930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 679930ae53cSBarry Smith x[i] = sum; 680930ae53cSBarry Smith } 681930ae53cSBarry Smith 682930ae53cSBarry Smith /* backward solve the upper triangular */ 683930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 684930ae53cSBarry Smith adiag_i = adiag[i]; 685d708749eSBarry Smith v = aa + adiag_i + 1; 686d708749eSBarry Smith vi = aj + adiag_i + 1; 687930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 688930ae53cSBarry Smith sum = x[i]; 689930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 690930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 691930ae53cSBarry Smith } 6921c4feecaSSatish Balay #endif 693b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 694cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 695cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6963a40ed3dSBarry Smith PetscFunctionReturn(0); 697930ae53cSBarry Smith } 698930ae53cSBarry Smith 6994a2ae208SSatish Balay #undef __FUNCT__ 7004a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 701416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 702da3a660dSBarry Smith { 703416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 704416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 705273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 7063b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 70787828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 708da3a660dSBarry Smith 7093a40ed3dSBarry Smith PetscFunctionBegin; 71078b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 711da3a660dSBarry Smith 712e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 713e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 714416022c9SBarry Smith tmp = a->solve_work; 715da3a660dSBarry Smith 7163b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7173b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 718da3a660dSBarry Smith 719da3a660dSBarry Smith /* forward solve the lower triangular */ 720da3a660dSBarry Smith tmp[0] = b[*r++]; 721da3a660dSBarry Smith for (i=1; i<n; i++) { 722dbb450caSBarry Smith v = aa + ai[i] + shift; 723dbb450caSBarry Smith vi = aj + ai[i] + shift; 724416022c9SBarry Smith nz = a->diag[i] - ai[i]; 725da3a660dSBarry Smith sum = b[*r++]; 726dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 727da3a660dSBarry Smith tmp[i] = sum; 728da3a660dSBarry Smith } 729da3a660dSBarry Smith 730da3a660dSBarry Smith /* backward solve the upper triangular */ 731da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 732416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 733416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 734416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 735da3a660dSBarry Smith sum = tmp[i]; 736dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 737416022c9SBarry Smith tmp[i] = sum*aa[a->diag[i]+shift]; 738da3a660dSBarry Smith x[*c--] += tmp[i]; 739da3a660dSBarry Smith } 740da3a660dSBarry Smith 7413b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7423b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 743cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 744cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 745b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 746898183ecSLois Curfman McInnes 7473a40ed3dSBarry Smith PetscFunctionReturn(0); 748da3a660dSBarry Smith } 749da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7504a2ae208SSatish Balay #undef __FUNCT__ 7514a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 7527c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 753da3a660dSBarry Smith { 754416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7552235e667SBarry Smith IS iscol = a->col,isrow = a->row; 756273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 757f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 75887828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 759da3a660dSBarry Smith 7603a40ed3dSBarry Smith PetscFunctionBegin; 761e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 762e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 763416022c9SBarry Smith tmp = a->solve_work; 764da3a660dSBarry Smith 7652235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7662235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 767da3a660dSBarry Smith 768da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7692235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 770da3a660dSBarry Smith 771da3a660dSBarry Smith /* forward solve the U^T */ 772da3a660dSBarry Smith for (i=0; i<n; i++) { 773f1af5d2fSBarry Smith v = aa + diag[i] + shift; 774f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 775f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 776c38d4ed2SBarry Smith s1 = tmp[i]; 777ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 778da3a660dSBarry Smith while (nz--) { 779f1af5d2fSBarry Smith tmp[*vi++ + shift] -= (*v++)*s1; 780da3a660dSBarry Smith } 781c38d4ed2SBarry Smith tmp[i] = s1; 782da3a660dSBarry Smith } 783da3a660dSBarry Smith 784da3a660dSBarry Smith /* backward solve the L^T */ 785da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 786f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 787f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 788f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 789f1af5d2fSBarry Smith s1 = tmp[i]; 790da3a660dSBarry Smith while (nz--) { 791f1af5d2fSBarry Smith tmp[*vi-- + shift] -= (*v--)*s1; 792da3a660dSBarry Smith } 793da3a660dSBarry Smith } 794da3a660dSBarry Smith 795da3a660dSBarry Smith /* copy tmp into x according to permutation */ 796da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 797da3a660dSBarry Smith 7982235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7992235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 800cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 801cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 802da3a660dSBarry Smith 803b0a32e0cSBarry Smith PetscLogFlops(2*a->nz-A->n); 8043a40ed3dSBarry Smith PetscFunctionReturn(0); 805da3a660dSBarry Smith } 806da3a660dSBarry Smith 8074a2ae208SSatish Balay #undef __FUNCT__ 8084a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 8097c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 810da3a660dSBarry Smith { 811416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8122235e667SBarry Smith IS iscol = a->col,isrow = a->row; 813273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 814f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 81587828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 8166abc6512SBarry Smith 8173a40ed3dSBarry Smith PetscFunctionBegin; 8186abc6512SBarry Smith if (zz != xx) VecCopy(zz,xx); 8196abc6512SBarry Smith 820e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 821e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 822416022c9SBarry Smith tmp = a->solve_work; 8236abc6512SBarry Smith 8242235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8252235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8266abc6512SBarry Smith 8276abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8282235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8296abc6512SBarry Smith 8306abc6512SBarry Smith /* forward solve the U^T */ 8316abc6512SBarry Smith for (i=0; i<n; i++) { 832f1af5d2fSBarry Smith v = aa + diag[i] + shift; 833f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 834f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8356abc6512SBarry Smith tmp[i] *= *v++; 8366abc6512SBarry Smith while (nz--) { 837dbb450caSBarry Smith tmp[*vi++ + shift] -= (*v++)*tmp[i]; 8386abc6512SBarry Smith } 8396abc6512SBarry Smith } 8406abc6512SBarry Smith 8416abc6512SBarry Smith /* backward solve the L^T */ 8426abc6512SBarry Smith for (i=n-1; i>=0; i--){ 843f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 844f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 845f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8466abc6512SBarry Smith while (nz--) { 847dbb450caSBarry Smith tmp[*vi-- + shift] -= (*v--)*tmp[i]; 8486abc6512SBarry Smith } 8496abc6512SBarry Smith } 8506abc6512SBarry Smith 8516abc6512SBarry Smith /* copy tmp into x according to permutation */ 8526abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8536abc6512SBarry Smith 8542235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8552235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 856cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 857cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8586abc6512SBarry Smith 859b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 8603a40ed3dSBarry Smith PetscFunctionReturn(0); 861da3a660dSBarry Smith } 8629e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 863ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat); 86408480c60SBarry Smith 8654a2ae208SSatish Balay #undef __FUNCT__ 8664a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 867b380c88cSHong Zhang int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8689e25ed09SBarry Smith { 869416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8709e25ed09SBarry Smith IS isicol; 871273d9f13SBarry Smith int *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j; 87256beaf04SBarry Smith int *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev; 873335d9088SBarry Smith int *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0; 8746cf997b0SBarry Smith int incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill; 87577c4ece6SBarry Smith PetscTruth col_identity,row_identity; 876329f5518SBarry Smith PetscReal f; 8779e25ed09SBarry Smith 8783a40ed3dSBarry Smith PetscFunctionBegin; 8796cf997b0SBarry Smith f = info->fill; 8800cd17293SBarry Smith levels = (int)info->levels; 8810cd17293SBarry Smith diagonal_fill = (int)info->diagonal_fill; 8824c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 88382bf6240SBarry Smith 88408480c60SBarry Smith /* special case that simply copies fill pattern */ 885be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 886be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 88786bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8882e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 88908480c60SBarry Smith (*fact)->factor = FACTOR_LU; 89008480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 89108480c60SBarry Smith if (!b->diag) { 8927c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89308480c60SBarry Smith } 8947c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89508480c60SBarry Smith b->row = isrow; 89608480c60SBarry Smith b->col = iscol; 89782bf6240SBarry Smith b->icol = isicol; 898a2e34c3dSBarry Smith b->lu_damping = info->damping; 89987828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 9006cc28720Svictorle b->lu_shift = info->lu_shift; 9016cc28720Svictorle b->lu_shift_fraction= info->lu_shift_fraction; 90287828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 903f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 904c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 905c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 9063a40ed3dSBarry Smith PetscFunctionReturn(0); 90708480c60SBarry Smith } 90808480c60SBarry Smith 909898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 910898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 9119e25ed09SBarry Smith 9129e25ed09SBarry Smith /* get new row pointers */ 913b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 914dbb450caSBarry Smith ainew[0] = -shift; 9159e25ed09SBarry Smith /* don't know how many column pointers are needed so estimate */ 916dbb450caSBarry Smith jmax = (int)(f*(ai[n]+!shift)); 917b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 9189e25ed09SBarry Smith /* ajfill is level of fill for each fill entry */ 919b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr); 9209e25ed09SBarry Smith /* fill is a linked list of nonzeros in active row */ 921b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr); 92256beaf04SBarry Smith /* im is level for each filled value */ 923b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr); 92456beaf04SBarry Smith /* dloc is location of diagonal in factor */ 925b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr); 92656beaf04SBarry Smith dloc[0] = 0; 92756beaf04SBarry Smith for (prow=0; prow<n; prow++) { 9286cf997b0SBarry Smith 9296cf997b0SBarry Smith /* copy current row into linked list */ 93056beaf04SBarry Smith nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 93129bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 932dbb450caSBarry Smith xi = aj + ai[r[prow]] + shift; 9339e25ed09SBarry Smith fill[n] = n; 934435faa5fSBarry Smith fill[prow] = -1; /* marker to indicate if diagonal exists */ 9359e25ed09SBarry Smith while (nz--) { 9369e25ed09SBarry Smith fm = n; 937dbb450caSBarry Smith idx = ic[*xi++ + shift]; 9389e25ed09SBarry Smith do { 9399e25ed09SBarry Smith m = fm; 9409e25ed09SBarry Smith fm = fill[m]; 9419e25ed09SBarry Smith } while (fm < idx); 9429e25ed09SBarry Smith fill[m] = idx; 9439e25ed09SBarry Smith fill[idx] = fm; 94456beaf04SBarry Smith im[idx] = 0; 9459e25ed09SBarry Smith } 9466cf997b0SBarry Smith 9476cf997b0SBarry Smith /* make sure diagonal entry is included */ 948435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 9496cf997b0SBarry Smith fm = n; 950435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 951435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 952435faa5fSBarry Smith fill[fm] = prow; 9536cf997b0SBarry Smith im[prow] = 0; 9546cf997b0SBarry Smith nzf++; 955335d9088SBarry Smith dcount++; 9566cf997b0SBarry Smith } 9576cf997b0SBarry Smith 95856beaf04SBarry Smith nzi = 0; 9599e25ed09SBarry Smith row = fill[n]; 96056beaf04SBarry Smith while (row < prow) { 96156beaf04SBarry Smith incrlev = im[row] + 1; 96256beaf04SBarry Smith nz = dloc[row]; 9636cf997b0SBarry Smith xi = ajnew + ainew[row] + shift + nz + 1; 964dbb450caSBarry Smith flev = ajfill + ainew[row] + shift + nz + 1; 965416022c9SBarry Smith nnz = ainew[row+1] - ainew[row] - nz - 1; 96656beaf04SBarry Smith fm = row; 96756beaf04SBarry Smith while (nnz-- > 0) { 968dbb450caSBarry Smith idx = *xi++ + shift; 96956beaf04SBarry Smith if (*flev + incrlev > levels) { 97056beaf04SBarry Smith flev++; 97156beaf04SBarry Smith continue; 97256beaf04SBarry Smith } 97356beaf04SBarry Smith do { 9749e25ed09SBarry Smith m = fm; 9759e25ed09SBarry Smith fm = fill[m]; 97656beaf04SBarry Smith } while (fm < idx); 9779e25ed09SBarry Smith if (fm != idx) { 97856beaf04SBarry Smith im[idx] = *flev + incrlev; 9799e25ed09SBarry Smith fill[m] = idx; 9809e25ed09SBarry Smith fill[idx] = fm; 9819e25ed09SBarry Smith fm = idx; 98256beaf04SBarry Smith nzf++; 983ecf371e4SBarry Smith } else { 98456beaf04SBarry Smith if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 9859e25ed09SBarry Smith } 98656beaf04SBarry Smith flev++; 9879e25ed09SBarry Smith } 9889e25ed09SBarry Smith row = fill[row]; 98956beaf04SBarry Smith nzi++; 9909e25ed09SBarry Smith } 9919e25ed09SBarry Smith /* copy new filled row into permanent storage */ 99256beaf04SBarry Smith ainew[prow+1] = ainew[prow] + nzf; 993d7e8b826SBarry Smith if (ainew[prow+1] > jmax-shift) { 994ecf371e4SBarry Smith 995ecf371e4SBarry Smith /* estimate how much additional space we will need */ 996ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 997ecf371e4SBarry Smith /* just double the memory each time */ 998ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */ 999ecf371e4SBarry Smith int maxadd = jmax; 100056beaf04SBarry Smith if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 1001bbb6d6a8SBarry Smith jmax += maxadd; 1002ecf371e4SBarry Smith 1003ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 1004b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 1005549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 1006606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 100756beaf04SBarry Smith ajnew = xi; 1008b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 1009549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 1010606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 101156beaf04SBarry Smith ajfill = xi; 1012ecf371e4SBarry Smith realloc++; /* count how many times we realloc */ 10139e25ed09SBarry Smith } 1014dbb450caSBarry Smith xi = ajnew + ainew[prow] + shift; 1015dbb450caSBarry Smith flev = ajfill + ainew[prow] + shift; 101656beaf04SBarry Smith dloc[prow] = nzi; 10179e25ed09SBarry Smith fm = fill[n]; 101856beaf04SBarry Smith while (nzf--) { 1019dbb450caSBarry Smith *xi++ = fm - shift; 102056beaf04SBarry Smith *flev++ = im[fm]; 10219e25ed09SBarry Smith fm = fill[fm]; 10229e25ed09SBarry Smith } 10236cf997b0SBarry Smith /* make sure row has diagonal entry */ 10246cf997b0SBarry Smith if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) { 102529bbc08cSBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\ 10266cf997b0SBarry Smith try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow); 10276cf997b0SBarry Smith } 10289e25ed09SBarry Smith } 1029606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 1030898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1031898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1032606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 1033606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 10349e25ed09SBarry Smith 1035f64065bbSBarry Smith { 1036329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 1037b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 1038c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af); 1039c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af); 1040b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"); 1041335d9088SBarry Smith if (diagonal_fill) { 1042b1bcba4aSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount); 1043335d9088SBarry Smith } 1044f64065bbSBarry Smith } 1045bbb6d6a8SBarry Smith 10469e25ed09SBarry Smith /* put together the new matrix */ 1047b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr); 1048b0a32e0cSBarry Smith PetscLogObjectParent(*fact,isicol); 1049416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1050606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10517c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 105287828ca2SBarry Smith len = (ainew[n] + shift)*sizeof(PetscScalar); 10539e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1054606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1055606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1056b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 1057416022c9SBarry Smith b->j = ajnew; 1058416022c9SBarry Smith b->i = ainew; 105956beaf04SBarry Smith for (i=0; i<n; i++) dloc[i] += ainew[i]; 1060416022c9SBarry Smith b->diag = dloc; 1061416022c9SBarry Smith b->ilen = 0; 1062416022c9SBarry Smith b->imax = 0; 1063416022c9SBarry Smith b->row = isrow; 1064416022c9SBarry Smith b->col = iscol; 1065c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1066c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 106782bf6240SBarry Smith b->icol = isicol; 106887828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1069416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 107056beaf04SBarry Smith Allocate dloc, solve_work, new a, new j */ 107187828ca2SBarry Smith PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(PetscScalar))); 1072416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 1073a2e34c3dSBarry Smith b->lu_damping = info->damping; 107487828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 10759e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 10763a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 10773a7fca6bSBarry Smith (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 1078ae068f58SLois Curfman McInnes 1079ae068f58SLois Curfman McInnes (*fact)->info.factor_mallocs = realloc; 1080ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 1081329f5518SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 10823a40ed3dSBarry Smith PetscFunctionReturn(0); 10839e25ed09SBarry Smith } 1084d5d45c9bSBarry Smith 1085a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1086a6175056SHong Zhang #undef __FUNCT__ 1087a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1088a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact) 1089a6175056SHong Zhang { 10900968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1091a6175056SHong Zhang int ierr; 1092d5d45c9bSBarry Smith 1093a6175056SHong Zhang PetscFunctionBegin; 10940968510aSHong Zhang if (!a->sbaijMat){ 10950968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1096a6175056SHong Zhang } 109703aac1b8SHong Zhang 1098b45a75daSHong Zhang ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr); 10990968510aSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 1100064503c5SHong Zhang a->sbaijMat = PETSC_NULL; 1101653a6975SHong Zhang 1102a6175056SHong Zhang PetscFunctionReturn(0); 1103a6175056SHong Zhang } 1104a6175056SHong Zhang 1105a6175056SHong Zhang #undef __FUNCT__ 1106a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 110715e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1108a6175056SHong Zhang { 11090968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 111003aac1b8SHong Zhang int ierr; 1111653a6975SHong Zhang PetscTruth perm_identity; 1112a6175056SHong Zhang 1113a6175056SHong Zhang PetscFunctionBegin; 1114653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1115653a6975SHong Zhang if (!perm_identity){ 1116653a6975SHong Zhang SETERRQ(1,"Non-identity permutation is not supported yet"); 1117653a6975SHong Zhang } 11180968510aSHong Zhang if (!a->sbaijMat){ 11190968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 11200968510aSHong Zhang } 1121a6175056SHong Zhang 1122b00f7748SHong Zhang ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 1123653a6975SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1124653a6975SHong Zhang 1125a6175056SHong Zhang PetscFunctionReturn(0); 1126a6175056SHong Zhang } 1127d5d45c9bSBarry Smith 1128f76d2b81SHong Zhang #undef __FUNCT__ 1129f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1130f76d2b81SHong Zhang int MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1131f76d2b81SHong Zhang { 1132f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 113303aac1b8SHong Zhang int ierr; 1134f76d2b81SHong Zhang PetscTruth perm_identity; 1135f76d2b81SHong Zhang 1136f76d2b81SHong Zhang PetscFunctionBegin; 1137f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1138f76d2b81SHong Zhang if (!perm_identity){ 1139f76d2b81SHong Zhang SETERRQ(1,"Non-identity permutation is not supported yet"); 1140f76d2b81SHong Zhang } 1141f76d2b81SHong Zhang if (!a->sbaijMat){ 1142f76d2b81SHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1143f76d2b81SHong Zhang } 1144f76d2b81SHong Zhang 1145f76d2b81SHong Zhang ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 1146f76d2b81SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1147f76d2b81SHong Zhang 1148f76d2b81SHong Zhang PetscFunctionReturn(0); 1149f76d2b81SHong Zhang } 1150