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 */ 5386bacbd2SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,MatILUInfo *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" 2759e046878SBarry Smith int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatLUInfo *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; 407c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 408c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 40982bf6240SBarry Smith b->icol = isicol; 41087828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 411416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 4127fd98bd6SLois Curfman McInnes Allocate idnew, solve_work, new a, new j */ 41387828ca2SBarry Smith PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(PetscScalar))); 414416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 4157fd98bd6SLois Curfman McInnes 416329f5518SBarry Smith (*B)->factor = FACTOR_LU; 417ae068f58SLois Curfman McInnes (*B)->info.factor_mallocs = realloc; 418ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 4193a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr); 4207dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 421703deb49SSatish Balay 422e93ccc4dSBarry Smith if (ai[n] != 0) { 423329f5518SBarry Smith (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 424eea2913aSSatish Balay } else { 425eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 426eea2913aSSatish Balay } 4273a40ed3dSBarry Smith PetscFunctionReturn(0); 428289bc588SBarry Smith } 4290a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4303a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth); 43141c01911SSatish Balay 4324a2ae208SSatish Balay #undef __FUNCT__ 4334a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 434416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B) 435289bc588SBarry Smith { 436416022c9SBarry Smith Mat C = *B; 437416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data; 43882bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 439273d9f13SBarry Smith int *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j; 4401987afe7SBarry Smith int *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift; 441085256b3SBarry Smith int *diag_offset = b->diag,diag; 442085256b3SBarry Smith int *pj,ndamp = 0; 44387828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 444d3d32019SBarry Smith PetscReal damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs; 4458a5eb818SBarry Smith PetscTruth damp; 446289bc588SBarry Smith 4473a40ed3dSBarry Smith PetscFunctionBegin; 44878b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 44978b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 45087828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 45187828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 4520cf60383SBarry Smith rtmps = rtmp + shift; ics = ic + shift; 453289bc588SBarry Smith 454085256b3SBarry Smith do { 4558a5eb818SBarry Smith damp = PETSC_FALSE; 456289bc588SBarry Smith for (i=0; i<n; i++) { 457289bc588SBarry Smith nz = ai[i+1] - ai[i]; 458dbb450caSBarry Smith ajtmp = aj + ai[i] + shift; 45948e94559SBarry Smith for (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0; 460289bc588SBarry Smith 461289bc588SBarry Smith /* load in initial (unfactored row) */ 462416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 463416022c9SBarry Smith ajtmpold = a->j + a->i[r[i]] + shift; 464416022c9SBarry Smith v = a->a + a->i[r[i]] + shift; 465085256b3SBarry Smith for (j=0; j<nz; j++) { 466085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] = v[j]; 4678a5eb818SBarry Smith if (ndamp && ajtmpold[j] == r[i]) { /* damp the diagonal of the matrix */ 468085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] += damping; 469085256b3SBarry Smith } 470085256b3SBarry Smith } 471289bc588SBarry Smith 472dbb450caSBarry Smith row = *ajtmp++ + shift; 473f2582111SSatish Balay while (row < i) { 4748c37ef55SBarry Smith pc = rtmp + row; 4758c37ef55SBarry Smith if (*pc != 0.0) { 4761987afe7SBarry Smith pv = b->a + diag_offset[row] + shift; 4771987afe7SBarry Smith pj = b->j + diag_offset[row] + (!shift); 47835aab85fSBarry Smith multiplier = *pc / *pv++; 4798c37ef55SBarry Smith *pc = multiplier; 480f2582111SSatish Balay nz = ai[row+1] - diag_offset[row] - 1; 481f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 482b0a32e0cSBarry Smith PetscLogFlops(2*nz); 483289bc588SBarry Smith } 484dbb450caSBarry Smith row = *ajtmp++ + shift; 485289bc588SBarry Smith } 486416022c9SBarry Smith /* finished row so stick it into b->a */ 487416022c9SBarry Smith pv = b->a + ai[i] + shift; 488416022c9SBarry Smith pj = b->j + ai[i] + shift; 4898c37ef55SBarry Smith nz = ai[i+1] - ai[i]; 49035aab85fSBarry Smith diag = diag_offset[i] - ai[i]; 491d3d32019SBarry Smith /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 492d3d32019SBarry Smith rs = 0.0; 493d3d32019SBarry Smith for (j=0; j<nz; j++) { 494d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 495d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 496d3d32019SBarry Smith } 497d3d32019SBarry Smith if (PetscAbsScalar(pv[diag]) < zeropivot*rs) { 498d3d32019SBarry Smith if (damping) { 4998a5eb818SBarry Smith if (ndamp) damping *= 2.0; 500085256b3SBarry Smith damp = PETSC_TRUE; 501085256b3SBarry Smith ndamp++; 502085256b3SBarry Smith break; 503085256b3SBarry Smith } else { 50463bb2aa6SBarry Smith SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs); 505d708749eSBarry Smith } 50635aab85fSBarry Smith } 50735aab85fSBarry Smith } 508085256b3SBarry Smith } while (damp); 5090f11f9aeSSatish Balay 51035aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 51135aab85fSBarry Smith for (i=0; i<n; i++) { 51235aab85fSBarry Smith b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift]; 51335aab85fSBarry Smith } 51435aab85fSBarry Smith 515606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 51678b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 51778b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 518416022c9SBarry Smith C->factor = FACTOR_LU; 5197dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 520c456f294SBarry Smith C->assembled = PETSC_TRUE; 521b0a32e0cSBarry Smith PetscLogFlops(C->n); 522d3d32019SBarry Smith if (ndamp) { 523b0a32e0cSBarry Smith PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping); 524085256b3SBarry Smith } 5253a40ed3dSBarry Smith PetscFunctionReturn(0); 526289bc588SBarry Smith } 5270a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5284a2ae208SSatish Balay #undef __FUNCT__ 5294a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 5309e046878SBarry Smith int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatLUInfo *info) 531da3a660dSBarry Smith { 532273d9f13SBarry Smith int ierr; 533416022c9SBarry Smith Mat C; 534416022c9SBarry Smith 5353a40ed3dSBarry Smith PetscFunctionBegin; 5369e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 537f2582111SSatish Balay ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr); 538273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 539440f4c53SSatish Balay PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol); 5403a40ed3dSBarry Smith PetscFunctionReturn(0); 541da3a660dSBarry Smith } 5420a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5434a2ae208SSatish Balay #undef __FUNCT__ 5444a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 545416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5468c37ef55SBarry Smith { 547416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 548416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 549273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 5503b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 55187828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 5528c37ef55SBarry Smith 5533a40ed3dSBarry Smith PetscFunctionBegin; 5543a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 55591bf9adeSBarry Smith 556e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 557e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 558416022c9SBarry Smith tmp = a->solve_work; 5598c37ef55SBarry Smith 5603b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 5613b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 5628c37ef55SBarry Smith 5638c37ef55SBarry Smith /* forward solve the lower triangular */ 5648c37ef55SBarry Smith tmp[0] = b[*r++]; 5650cf60383SBarry Smith tmps = tmp + shift; 5668c37ef55SBarry Smith for (i=1; i<n; i++) { 567dbb450caSBarry Smith v = aa + ai[i] + shift; 568dbb450caSBarry Smith vi = aj + ai[i] + shift; 56953e7425aSBarry Smith nz = a->diag[i] - ai[i]; 57053e7425aSBarry Smith sum = b[*r++]; 571ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 5728c37ef55SBarry Smith tmp[i] = sum; 5738c37ef55SBarry Smith } 5748c37ef55SBarry Smith 5758c37ef55SBarry Smith /* backward solve the upper triangular */ 5768c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 577416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 578416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 579416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 5808c37ef55SBarry Smith sum = tmp[i]; 581ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 582416022c9SBarry Smith x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift]; 5838c37ef55SBarry Smith } 5848c37ef55SBarry Smith 5853b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 5863b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 587cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 588cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 589b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 5903a40ed3dSBarry Smith PetscFunctionReturn(0); 5918c37ef55SBarry Smith } 592026e39d0SSatish Balay 593930ae53cSBarry Smith /* ----------------------------------------------------------- */ 5944a2ae208SSatish Balay #undef __FUNCT__ 5954a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 596930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 597930ae53cSBarry Smith { 598930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 599273d9f13SBarry Smith int n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr; 600362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 601aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 602d85afc42SSatish Balay int adiag_i,i,*vi,nz,ai_i; 603362ced78SSatish Balay PetscScalar *v,sum; 604d85afc42SSatish Balay #endif 605930ae53cSBarry Smith 6063a40ed3dSBarry Smith PetscFunctionBegin; 6073a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 608930ae53cSBarry Smith if (a->indexshift) { 6093a40ed3dSBarry Smith ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr); 6103a40ed3dSBarry Smith PetscFunctionReturn(0); 611930ae53cSBarry Smith } 612930ae53cSBarry Smith 613e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 614e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 615930ae53cSBarry Smith 616aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6171c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6181c4feecaSSatish Balay #else 619930ae53cSBarry Smith /* forward solve the lower triangular */ 620930ae53cSBarry Smith x[0] = b[0]; 621930ae53cSBarry Smith for (i=1; i<n; i++) { 622930ae53cSBarry Smith ai_i = ai[i]; 623930ae53cSBarry Smith v = aa + ai_i; 624930ae53cSBarry Smith vi = aj + ai_i; 625930ae53cSBarry Smith nz = adiag[i] - ai_i; 626930ae53cSBarry Smith sum = b[i]; 627930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 628930ae53cSBarry Smith x[i] = sum; 629930ae53cSBarry Smith } 630930ae53cSBarry Smith 631930ae53cSBarry Smith /* backward solve the upper triangular */ 632930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 633930ae53cSBarry Smith adiag_i = adiag[i]; 634d708749eSBarry Smith v = aa + adiag_i + 1; 635d708749eSBarry Smith vi = aj + adiag_i + 1; 636930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 637930ae53cSBarry Smith sum = x[i]; 638930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 639930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 640930ae53cSBarry Smith } 6411c4feecaSSatish Balay #endif 642b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 643cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 644cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 646930ae53cSBarry Smith } 647930ae53cSBarry Smith 6484a2ae208SSatish Balay #undef __FUNCT__ 6494a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 650416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 651da3a660dSBarry Smith { 652416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 653416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 654273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 6553b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 65687828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 657da3a660dSBarry Smith 6583a40ed3dSBarry Smith PetscFunctionBegin; 65978b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 660da3a660dSBarry Smith 661e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 662e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 663416022c9SBarry Smith tmp = a->solve_work; 664da3a660dSBarry Smith 6653b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6663b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 667da3a660dSBarry Smith 668da3a660dSBarry Smith /* forward solve the lower triangular */ 669da3a660dSBarry Smith tmp[0] = b[*r++]; 670da3a660dSBarry Smith for (i=1; i<n; i++) { 671dbb450caSBarry Smith v = aa + ai[i] + shift; 672dbb450caSBarry Smith vi = aj + ai[i] + shift; 673416022c9SBarry Smith nz = a->diag[i] - ai[i]; 674da3a660dSBarry Smith sum = b[*r++]; 675dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 676da3a660dSBarry Smith tmp[i] = sum; 677da3a660dSBarry Smith } 678da3a660dSBarry Smith 679da3a660dSBarry Smith /* backward solve the upper triangular */ 680da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 681416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 682416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 683416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 684da3a660dSBarry Smith sum = tmp[i]; 685dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 686416022c9SBarry Smith tmp[i] = sum*aa[a->diag[i]+shift]; 687da3a660dSBarry Smith x[*c--] += tmp[i]; 688da3a660dSBarry Smith } 689da3a660dSBarry Smith 6903b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6913b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 692cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 693cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 694b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 695898183ecSLois Curfman McInnes 6963a40ed3dSBarry Smith PetscFunctionReturn(0); 697da3a660dSBarry Smith } 698da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 6994a2ae208SSatish Balay #undef __FUNCT__ 7004a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 7017c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 702da3a660dSBarry Smith { 703416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7042235e667SBarry 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; 706f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 70787828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 708da3a660dSBarry Smith 7093a40ed3dSBarry Smith PetscFunctionBegin; 710e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 711e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 712416022c9SBarry Smith tmp = a->solve_work; 713da3a660dSBarry Smith 7142235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7152235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 716da3a660dSBarry Smith 717da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7182235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 719da3a660dSBarry Smith 720da3a660dSBarry Smith /* forward solve the U^T */ 721da3a660dSBarry Smith for (i=0; i<n; i++) { 722f1af5d2fSBarry Smith v = aa + diag[i] + shift; 723f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 724f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 725c38d4ed2SBarry Smith s1 = tmp[i]; 726ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 727da3a660dSBarry Smith while (nz--) { 728f1af5d2fSBarry Smith tmp[*vi++ + shift] -= (*v++)*s1; 729da3a660dSBarry Smith } 730c38d4ed2SBarry Smith tmp[i] = s1; 731da3a660dSBarry Smith } 732da3a660dSBarry Smith 733da3a660dSBarry Smith /* backward solve the L^T */ 734da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 735f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 736f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 737f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 738f1af5d2fSBarry Smith s1 = tmp[i]; 739da3a660dSBarry Smith while (nz--) { 740f1af5d2fSBarry Smith tmp[*vi-- + shift] -= (*v--)*s1; 741da3a660dSBarry Smith } 742da3a660dSBarry Smith } 743da3a660dSBarry Smith 744da3a660dSBarry Smith /* copy tmp into x according to permutation */ 745da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 746da3a660dSBarry Smith 7472235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7482235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 749cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 750cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 751da3a660dSBarry Smith 752b0a32e0cSBarry Smith PetscLogFlops(2*a->nz-A->n); 7533a40ed3dSBarry Smith PetscFunctionReturn(0); 754da3a660dSBarry Smith } 755da3a660dSBarry Smith 7564a2ae208SSatish Balay #undef __FUNCT__ 7574a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 7587c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 759da3a660dSBarry Smith { 760416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7612235e667SBarry Smith IS iscol = a->col,isrow = a->row; 762273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 763f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 76487828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 7656abc6512SBarry Smith 7663a40ed3dSBarry Smith PetscFunctionBegin; 7676abc6512SBarry Smith if (zz != xx) VecCopy(zz,xx); 7686abc6512SBarry Smith 769e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 770e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 771416022c9SBarry Smith tmp = a->solve_work; 7726abc6512SBarry Smith 7732235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7742235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 7756abc6512SBarry Smith 7766abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 7772235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 7786abc6512SBarry Smith 7796abc6512SBarry Smith /* forward solve the U^T */ 7806abc6512SBarry Smith for (i=0; i<n; i++) { 781f1af5d2fSBarry Smith v = aa + diag[i] + shift; 782f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 783f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 7846abc6512SBarry Smith tmp[i] *= *v++; 7856abc6512SBarry Smith while (nz--) { 786dbb450caSBarry Smith tmp[*vi++ + shift] -= (*v++)*tmp[i]; 7876abc6512SBarry Smith } 7886abc6512SBarry Smith } 7896abc6512SBarry Smith 7906abc6512SBarry Smith /* backward solve the L^T */ 7916abc6512SBarry Smith for (i=n-1; i>=0; i--){ 792f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 793f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 794f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 7956abc6512SBarry Smith while (nz--) { 796dbb450caSBarry Smith tmp[*vi-- + shift] -= (*v--)*tmp[i]; 7976abc6512SBarry Smith } 7986abc6512SBarry Smith } 7996abc6512SBarry Smith 8006abc6512SBarry Smith /* copy tmp into x according to permutation */ 8016abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8026abc6512SBarry Smith 8032235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8042235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 805cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 806cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8076abc6512SBarry Smith 808b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 8093a40ed3dSBarry Smith PetscFunctionReturn(0); 810da3a660dSBarry Smith } 8119e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 812ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat); 81308480c60SBarry Smith 8144a2ae208SSatish Balay #undef __FUNCT__ 8154a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 8166cf997b0SBarry Smith int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatILUInfo *info,Mat *fact) 8179e25ed09SBarry Smith { 818416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8199e25ed09SBarry Smith IS isicol; 820273d9f13SBarry Smith int *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j; 82156beaf04SBarry Smith int *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev; 822335d9088SBarry Smith int *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0; 8236cf997b0SBarry Smith int incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill; 82477c4ece6SBarry Smith PetscTruth col_identity,row_identity; 825329f5518SBarry Smith PetscReal f; 8269e25ed09SBarry Smith 8273a40ed3dSBarry Smith PetscFunctionBegin; 8286cf997b0SBarry Smith f = info->fill; 8290cd17293SBarry Smith levels = (int)info->levels; 8300cd17293SBarry Smith diagonal_fill = (int)info->diagonal_fill; 8314c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 83282bf6240SBarry Smith 83308480c60SBarry Smith /* special case that simply copies fill pattern */ 834be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 835be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 83686bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8372e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 83808480c60SBarry Smith (*fact)->factor = FACTOR_LU; 83908480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 84008480c60SBarry Smith if (!b->diag) { 8417c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 84208480c60SBarry Smith } 8437c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 84408480c60SBarry Smith b->row = isrow; 84508480c60SBarry Smith b->col = iscol; 84682bf6240SBarry Smith b->icol = isicol; 847a2e34c3dSBarry Smith b->lu_damping = info->damping; 84887828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 84987828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 850f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 851c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 852c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 8533a40ed3dSBarry Smith PetscFunctionReturn(0); 85408480c60SBarry Smith } 85508480c60SBarry Smith 856898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 857898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 8589e25ed09SBarry Smith 8599e25ed09SBarry Smith /* get new row pointers */ 860b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 861dbb450caSBarry Smith ainew[0] = -shift; 8629e25ed09SBarry Smith /* don't know how many column pointers are needed so estimate */ 863dbb450caSBarry Smith jmax = (int)(f*(ai[n]+!shift)); 864b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 8659e25ed09SBarry Smith /* ajfill is level of fill for each fill entry */ 866b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr); 8679e25ed09SBarry Smith /* fill is a linked list of nonzeros in active row */ 868b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr); 86956beaf04SBarry Smith /* im is level for each filled value */ 870b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr); 87156beaf04SBarry Smith /* dloc is location of diagonal in factor */ 872b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr); 87356beaf04SBarry Smith dloc[0] = 0; 87456beaf04SBarry Smith for (prow=0; prow<n; prow++) { 8756cf997b0SBarry Smith 8766cf997b0SBarry Smith /* copy current row into linked list */ 87756beaf04SBarry Smith nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 87829bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 879dbb450caSBarry Smith xi = aj + ai[r[prow]] + shift; 8809e25ed09SBarry Smith fill[n] = n; 881435faa5fSBarry Smith fill[prow] = -1; /* marker to indicate if diagonal exists */ 8829e25ed09SBarry Smith while (nz--) { 8839e25ed09SBarry Smith fm = n; 884dbb450caSBarry Smith idx = ic[*xi++ + shift]; 8859e25ed09SBarry Smith do { 8869e25ed09SBarry Smith m = fm; 8879e25ed09SBarry Smith fm = fill[m]; 8889e25ed09SBarry Smith } while (fm < idx); 8899e25ed09SBarry Smith fill[m] = idx; 8909e25ed09SBarry Smith fill[idx] = fm; 89156beaf04SBarry Smith im[idx] = 0; 8929e25ed09SBarry Smith } 8936cf997b0SBarry Smith 8946cf997b0SBarry Smith /* make sure diagonal entry is included */ 895435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 8966cf997b0SBarry Smith fm = n; 897435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 898435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 899435faa5fSBarry Smith fill[fm] = prow; 9006cf997b0SBarry Smith im[prow] = 0; 9016cf997b0SBarry Smith nzf++; 902335d9088SBarry Smith dcount++; 9036cf997b0SBarry Smith } 9046cf997b0SBarry Smith 90556beaf04SBarry Smith nzi = 0; 9069e25ed09SBarry Smith row = fill[n]; 90756beaf04SBarry Smith while (row < prow) { 90856beaf04SBarry Smith incrlev = im[row] + 1; 90956beaf04SBarry Smith nz = dloc[row]; 9106cf997b0SBarry Smith xi = ajnew + ainew[row] + shift + nz + 1; 911dbb450caSBarry Smith flev = ajfill + ainew[row] + shift + nz + 1; 912416022c9SBarry Smith nnz = ainew[row+1] - ainew[row] - nz - 1; 91356beaf04SBarry Smith fm = row; 91456beaf04SBarry Smith while (nnz-- > 0) { 915dbb450caSBarry Smith idx = *xi++ + shift; 91656beaf04SBarry Smith if (*flev + incrlev > levels) { 91756beaf04SBarry Smith flev++; 91856beaf04SBarry Smith continue; 91956beaf04SBarry Smith } 92056beaf04SBarry Smith do { 9219e25ed09SBarry Smith m = fm; 9229e25ed09SBarry Smith fm = fill[m]; 92356beaf04SBarry Smith } while (fm < idx); 9249e25ed09SBarry Smith if (fm != idx) { 92556beaf04SBarry Smith im[idx] = *flev + incrlev; 9269e25ed09SBarry Smith fill[m] = idx; 9279e25ed09SBarry Smith fill[idx] = fm; 9289e25ed09SBarry Smith fm = idx; 92956beaf04SBarry Smith nzf++; 930ecf371e4SBarry Smith } else { 93156beaf04SBarry Smith if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 9329e25ed09SBarry Smith } 93356beaf04SBarry Smith flev++; 9349e25ed09SBarry Smith } 9359e25ed09SBarry Smith row = fill[row]; 93656beaf04SBarry Smith nzi++; 9379e25ed09SBarry Smith } 9389e25ed09SBarry Smith /* copy new filled row into permanent storage */ 93956beaf04SBarry Smith ainew[prow+1] = ainew[prow] + nzf; 940d7e8b826SBarry Smith if (ainew[prow+1] > jmax-shift) { 941ecf371e4SBarry Smith 942ecf371e4SBarry Smith /* estimate how much additional space we will need */ 943ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 944ecf371e4SBarry Smith /* just double the memory each time */ 945ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */ 946ecf371e4SBarry Smith int maxadd = jmax; 94756beaf04SBarry Smith if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 948bbb6d6a8SBarry Smith jmax += maxadd; 949ecf371e4SBarry Smith 950ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 951b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 952549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 953606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 95456beaf04SBarry Smith ajnew = xi; 955b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 956549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 957606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 95856beaf04SBarry Smith ajfill = xi; 959ecf371e4SBarry Smith realloc++; /* count how many times we realloc */ 9609e25ed09SBarry Smith } 961dbb450caSBarry Smith xi = ajnew + ainew[prow] + shift; 962dbb450caSBarry Smith flev = ajfill + ainew[prow] + shift; 96356beaf04SBarry Smith dloc[prow] = nzi; 9649e25ed09SBarry Smith fm = fill[n]; 96556beaf04SBarry Smith while (nzf--) { 966dbb450caSBarry Smith *xi++ = fm - shift; 96756beaf04SBarry Smith *flev++ = im[fm]; 9689e25ed09SBarry Smith fm = fill[fm]; 9699e25ed09SBarry Smith } 9706cf997b0SBarry Smith /* make sure row has diagonal entry */ 9716cf997b0SBarry Smith if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) { 97229bbc08cSBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\ 9736cf997b0SBarry Smith try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow); 9746cf997b0SBarry Smith } 9759e25ed09SBarry Smith } 976606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 977898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 978898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 979606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 980606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 9819e25ed09SBarry Smith 982f64065bbSBarry Smith { 983329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 984b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 985c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af); 986c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af); 987b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"); 988335d9088SBarry Smith if (diagonal_fill) { 989b1bcba4aSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount); 990335d9088SBarry Smith } 991f64065bbSBarry Smith } 992bbb6d6a8SBarry Smith 9939e25ed09SBarry Smith /* put together the new matrix */ 994b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr); 995b0a32e0cSBarry Smith PetscLogObjectParent(*fact,isicol); 996416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 997606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 9987c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 99987828ca2SBarry Smith len = (ainew[n] + shift)*sizeof(PetscScalar); 10009e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1001606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1002606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1003b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 1004416022c9SBarry Smith b->j = ajnew; 1005416022c9SBarry Smith b->i = ainew; 100656beaf04SBarry Smith for (i=0; i<n; i++) dloc[i] += ainew[i]; 1007416022c9SBarry Smith b->diag = dloc; 1008416022c9SBarry Smith b->ilen = 0; 1009416022c9SBarry Smith b->imax = 0; 1010416022c9SBarry Smith b->row = isrow; 1011416022c9SBarry Smith b->col = iscol; 1012c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1013c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 101482bf6240SBarry Smith b->icol = isicol; 101587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1016416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 101756beaf04SBarry Smith Allocate dloc, solve_work, new a, new j */ 101887828ca2SBarry Smith PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(PetscScalar))); 1019416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 1020a2e34c3dSBarry Smith b->lu_damping = info->damping; 102187828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 10229e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 10233a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 10243a7fca6bSBarry Smith (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 1025ae068f58SLois Curfman McInnes 1026ae068f58SLois Curfman McInnes (*fact)->info.factor_mallocs = realloc; 1027ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 1028329f5518SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 10293a40ed3dSBarry Smith PetscFunctionReturn(0); 10309e25ed09SBarry Smith } 1031d5d45c9bSBarry Smith 1032a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1033653a6975SHong Zhang typedef struct { 1034653a6975SHong Zhang int levels; 1035653a6975SHong Zhang } Mat_SeqAIJ_SeqSBAIJ; 1036653a6975SHong Zhang 1037a6175056SHong Zhang #undef __FUNCT__ 1038a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1039a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact) 1040a6175056SHong Zhang { 10410968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1042a6175056SHong Zhang int ierr; 1043653a6975SHong Zhang Mat_SeqAIJ_SeqSBAIJ *ptr = (Mat_SeqAIJ_SeqSBAIJ*)(*fact)->spptr; 1044d5d45c9bSBarry Smith 1045a6175056SHong Zhang PetscFunctionBegin; 10460968510aSHong Zhang if (!a->sbaijMat){ 10470968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1048a6175056SHong Zhang } 1049064503c5SHong Zhang ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(a->sbaijMat,fact);CHKERRQ(ierr); 1050064503c5SHong Zhang if (ptr->levels > 0){ 10510968510aSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 10520968510aSHong Zhang } 1053064503c5SHong Zhang a->sbaijMat = PETSC_NULL; 1054653a6975SHong Zhang ierr = PetscFree(ptr);CHKERRQ(ierr); 1055653a6975SHong Zhang 1056a6175056SHong Zhang PetscFunctionReturn(0); 1057a6175056SHong Zhang } 1058a6175056SHong Zhang 1059a6175056SHong Zhang #undef __FUNCT__ 1060a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1061*15e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1062a6175056SHong Zhang { 10630968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1064653a6975SHong Zhang Mat_SeqSBAIJ *b; 1065a6175056SHong Zhang int ierr; 1066653a6975SHong Zhang PetscTruth perm_identity; 1067653a6975SHong Zhang Mat_SeqAIJ_SeqSBAIJ *ptr; 1068b00f7748SHong Zhang int levels = info->levels; 1069b00f7748SHong Zhang PetscReal fill = info->fill; 1070a6175056SHong Zhang 1071a6175056SHong Zhang PetscFunctionBegin; 1072653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1073653a6975SHong Zhang if (!perm_identity){ 1074653a6975SHong Zhang SETERRQ(1,"Non-identity permutation is not supported yet"); 1075653a6975SHong Zhang } 10760968510aSHong Zhang if (!a->sbaijMat){ 10770968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 10780968510aSHong Zhang } 1079a6175056SHong Zhang 1080653a6975SHong Zhang if (levels > 0){ 1081b00f7748SHong Zhang ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 1082653a6975SHong Zhang 1083653a6975SHong Zhang } else { /* in-place icc(0) */ 1084653a6975SHong Zhang (*fact) = a->sbaijMat; 1085653a6975SHong Zhang (*fact)->factor = FACTOR_CHOLESKY; 1086653a6975SHong Zhang b = (Mat_SeqSBAIJ*)(*fact)->data; 1087653a6975SHong Zhang b->row = perm; 1088653a6975SHong Zhang b->icol = perm; 1089653a6975SHong Zhang ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1090653a6975SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1091653a6975SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1092653a6975SHong Zhang } 1093653a6975SHong Zhang 1094653a6975SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1095064503c5SHong Zhang (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1096653a6975SHong Zhang 1097653a6975SHong Zhang ierr = PetscNew(Mat_SeqAIJ_SeqSBAIJ,&ptr);CHKERRQ(ierr); 1098653a6975SHong Zhang (*fact)->spptr = (void*)ptr; 1099653a6975SHong Zhang ptr->levels = levels; /* to be used by CholeskyFactorNumeric_SeqAIJ() */ 1100653a6975SHong Zhang 1101a6175056SHong Zhang PetscFunctionReturn(0); 1102a6175056SHong Zhang } 1103d5d45c9bSBarry Smith 1104