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" 41c4feecaSSatish Balay #include "src/inline/dot.h" 5ed480e8bSBarry Smith #include "src/inline/spops.h" 63b2fbd54SBarry Smith 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 9234db7c0SBarry Smith int MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 103b2fbd54SBarry Smith { 113a40ed3dSBarry Smith PetscFunctionBegin; 123a40ed3dSBarry Smith 1329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Code not written"); 14aa482453SBarry Smith #if !defined(PETSC_USE_DEBUG) 15d64ed03dSBarry Smith PetscFunctionReturn(0); 16d64ed03dSBarry Smith #endif 173b2fbd54SBarry Smith } 183b2fbd54SBarry Smith 1986bacbd2SBarry Smith 20ca44d042SBarry Smith EXTERN int MatMarkDiagonal_SeqAIJ(Mat); 213a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth); 2286bacbd2SBarry Smith 2387828ca2SBarry Smith EXTERN int SPARSEKIT2dperm(int*,PetscScalar*,int*,int*,PetscScalar*,int*,int*,int*,int*,int*); 2487828ca2SBarry Smith EXTERN int SPARSEKIT2ilutp(int*,PetscScalar*,int*,int*,int*,PetscReal,PetscReal*,int*,PetscScalar*,int*,int*,int*,PetscScalar*,int*,int*,int*); 2587828ca2SBarry Smith EXTERN int SPARSEKIT2msrcsr(int*,PetscScalar*,int*,PetscScalar*,int*,int*,PetscScalar*,int*); 2686bacbd2SBarry Smith 274a2ae208SSatish Balay #undef __FUNCT__ 284a2ae208SSatish Balay #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 2986bacbd2SBarry Smith /* ------------------------------------------------------------ 3086bacbd2SBarry Smith 3186bacbd2SBarry Smith This interface was contribed by Tony Caola 3286bacbd2SBarry Smith 3386bacbd2SBarry Smith This routine is an interface to the pivoting drop-tolerance 345bd2ddc7SBarry Smith ILU routine written by Yousef Saad (saad@cs.umn.edu) as part of 355bd2ddc7SBarry Smith SPARSEKIT2. 365bd2ddc7SBarry Smith 375bd2ddc7SBarry Smith The SPARSEKIT2 routines used here are covered by the GNU 385bd2ddc7SBarry Smith copyright; see the file gnu in this directory. 3986bacbd2SBarry Smith 4086bacbd2SBarry Smith Thanks to Prof. Saad, Dr. Hysom, and Dr. Smith for their 4186bacbd2SBarry Smith help in getting this routine ironed out. 4286bacbd2SBarry Smith 435bd2ddc7SBarry Smith The major drawback to this routine is that if info->fill is 445bd2ddc7SBarry Smith not large enough it fails rather than allocating more space; 455bd2ddc7SBarry Smith this can be fixed by hacking/improving the f2c version of 465bd2ddc7SBarry Smith Yousef Saad's code. 4786bacbd2SBarry Smith 4886bacbd2SBarry Smith ------------------------------------------------------------ 4986bacbd2SBarry Smith */ 50b380c88cSHong Zhang int MatILUDTFactor_SeqAIJ(Mat A,MatFactorInfo *info,IS isrow,IS iscol,Mat *fact) 5186bacbd2SBarry Smith { 5260ec86bdSBarry Smith #if defined(PETSC_AVOID_GNUCOPYRIGHT_CODE) 5360ec86bdSBarry Smith PetscFunctionBegin; 5460ec86bdSBarry Smith SETERRQ(1,"This distribution does not include GNU Copyright code\n\ 5560ec86bdSBarry Smith You can obtain the drop tolerance routines by installing PETSc from\n\ 5660ec86bdSBarry Smith www.mcs.anl.gov/petsc\n"); 5760ec86bdSBarry Smith #else 5886bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 5907d2056aSBarry Smith IS iscolf,isicol,isirow; 6086bacbd2SBarry Smith PetscTruth reorder; 61273d9f13SBarry Smith int *c,*r,*ic,ierr,i,n = A->m; 62329f5518SBarry Smith int *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 63313ae024SBarry Smith int *ordcol,*iwk,*iperm,*jw; 64b19713cbSBarry Smith int jmax,lfill,job,*o_i,*o_j; 6587828ca2SBarry Smith PetscScalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 66329f5518SBarry Smith PetscReal permtol,af; 6786bacbd2SBarry Smith 6886bacbd2SBarry Smith PetscFunctionBegin; 6986bacbd2SBarry Smith 7086bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 7186bacbd2SBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax); 7286bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7333259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 746faa4630SBarry Smith lfill = (int)(info->dtcount/2.0); 756faa4630SBarry Smith jmax = (int)(info->fill*a->nz); 7686bacbd2SBarry Smith permtol = info->dtcol; 7786bacbd2SBarry Smith 7886bacbd2SBarry Smith 7986bacbd2SBarry Smith /* ------------------------------------------------------------ 8086bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 8186bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8286bacbd2SBarry Smith then de-reordered so it is in it's original format. 8386bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8486bacbd2SBarry Smith the original matrix and allocate more storage. . . 8586bacbd2SBarry Smith ------------------------------------------------------------ 8686bacbd2SBarry Smith */ 8786bacbd2SBarry Smith 8886bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 8986bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 9086bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 9186bacbd2SBarry Smith reorder = PetscNot(reorder); 9286bacbd2SBarry Smith 9386bacbd2SBarry Smith 9486bacbd2SBarry Smith /* storage for ilu factor */ 95b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr); 96b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr); 9787828ca2SBarry Smith ierr = PetscMalloc(jmax*sizeof(PetscScalar),&new_a);CHKERRQ(ierr); 98b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr); 9986bacbd2SBarry Smith 10086bacbd2SBarry Smith /* ------------------------------------------------------------ 10186bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10286bacbd2SBarry Smith ------------------------------------------------------------ 10386bacbd2SBarry Smith */ 104b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 105b19713cbSBarry Smith old_j[i]++; 10686bacbd2SBarry Smith } 107b19713cbSBarry Smith for(i=0;i<n+1;i++) { 108b19713cbSBarry Smith old_i[i]++; 109b19713cbSBarry Smith }; 110010a20caSHong Zhang 11186bacbd2SBarry Smith 11286bacbd2SBarry Smith if (reorder) { 113c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c); CHKERRQ(ierr); 114c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r); CHKERRQ(ierr); 115c0c2c81eSBarry Smith for(i=0;i<n;i++) { 116c0c2c81eSBarry Smith r[i] = r[i]+1; 117c0c2c81eSBarry Smith c[i] = c[i]+1; 118c0c2c81eSBarry Smith } 119b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr); 120b0a32e0cSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr); 12187828ca2SBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(PetscScalar),&old_a2);CHKERRQ(ierr); 1225bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 123c0c2c81eSBarry Smith for (i=0;i<n;i++) { 124c0c2c81eSBarry Smith r[i] = r[i]-1; 125c0c2c81eSBarry Smith c[i] = c[i]-1; 126c0c2c81eSBarry Smith } 127c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 128c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 129b19713cbSBarry Smith o_a = old_a2; 130b19713cbSBarry Smith o_j = old_j2; 131b19713cbSBarry Smith o_i = old_i2; 132b19713cbSBarry Smith } else { 133b19713cbSBarry Smith o_a = old_a; 134b19713cbSBarry Smith o_j = old_j; 135b19713cbSBarry Smith o_i = old_i; 13686bacbd2SBarry Smith } 13786bacbd2SBarry Smith 13886bacbd2SBarry Smith /* ------------------------------------------------------------ 13986bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 14086bacbd2SBarry Smith ------------------------------------------------------------ 14186bacbd2SBarry Smith */ 14286bacbd2SBarry Smith 143b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr); 144b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr); 14587828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&w);CHKERRQ(ierr); 14686bacbd2SBarry Smith 14787828ca2SBarry 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); 14886bacbd2SBarry Smith if (ierr) { 14986bacbd2SBarry Smith switch (ierr) { 15029bbc08cSBarry Smith case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15129bbc08cSBarry Smith case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15229bbc08cSBarry Smith case -5: SETERRQ(1,"ilutp(), zero row encountered"); 15329bbc08cSBarry Smith case -1: SETERRQ(1,"ilutp(), input matrix may be wrong"); 15429bbc08cSBarry Smith case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax); 15529bbc08cSBarry Smith default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr); 15686bacbd2SBarry Smith } 15786bacbd2SBarry Smith } 15886bacbd2SBarry Smith 15986bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 16086bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16186bacbd2SBarry Smith 16286bacbd2SBarry Smith /* ------------------------------------------------------------ 16386bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16486bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16586bacbd2SBarry Smith ------------------------------------------------------------ 16686bacbd2SBarry Smith */ 16786bacbd2SBarry Smith 16887828ca2SBarry Smith ierr = PetscMalloc(n*sizeof(PetscScalar),&wk);CHKERRQ(ierr); 169b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr); 17086bacbd2SBarry Smith 1715bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17286bacbd2SBarry Smith 17386bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17486bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17586bacbd2SBarry Smith 17686bacbd2SBarry Smith if (reorder) { 17786bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 17886bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 17986bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 180313ae024SBarry Smith } else { 181b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 182f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 183b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 184b19713cbSBarry Smith } 185b19713cbSBarry Smith } 186b19713cbSBarry Smith 187b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 18886bacbd2SBarry Smith for (i=0; i<n+1; i++) { 189b19713cbSBarry Smith old_i[i]--; 190b19713cbSBarry Smith } 191b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 192b19713cbSBarry Smith old_j[i]--; 193b19713cbSBarry Smith } 19486bacbd2SBarry Smith 195b19713cbSBarry Smith /* Make the factored matrix 0-based */ 19686bacbd2SBarry Smith for (i=0; i<n+1; i++) { 197b19713cbSBarry Smith new_i[i]--; 198b19713cbSBarry Smith } 199b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 200b19713cbSBarry Smith new_j[i]--; 201b19713cbSBarry Smith } 20286bacbd2SBarry Smith 20386bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20486bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2054c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2064c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 207c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 20886bacbd2SBarry Smith for(i=0; i<n; i++) { 20986bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21086bacbd2SBarry Smith }; 21186bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21207d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21386bacbd2SBarry Smith 214c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 215c0c2c81eSBarry Smith 216329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 21707d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 21886bacbd2SBarry Smith 21986bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22086bacbd2SBarry Smith 221*f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 222*f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 223*f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 22486bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22586bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 22686bacbd2SBarry Smith 22786bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 22886bacbd2SBarry Smith ierr = PetscFree(b->imax);CHKERRQ(ierr); 22986bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23007d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 231b19713cbSBarry Smith /* the next line frees the default space generated by the MatCreate() */ 23286bacbd2SBarry Smith ierr = PetscFree(b->a);CHKERRQ(ierr); 23386bacbd2SBarry Smith ierr = PetscFree(b->ilen);CHKERRQ(ierr); 23486bacbd2SBarry Smith b->a = new_a; 23586bacbd2SBarry Smith b->j = new_j; 23686bacbd2SBarry Smith b->i = new_i; 23786bacbd2SBarry Smith b->ilen = 0; 23886bacbd2SBarry Smith b->imax = 0; 2391f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 240313ae024SBarry Smith b->row = isirow; 24186bacbd2SBarry Smith b->col = iscolf; 24287828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 24386bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 24486bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 24586bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 24686bacbd2SBarry Smith 24786bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 24886bacbd2SBarry Smith 24986bacbd2SBarry Smith /* check out for identical nodes. If found, use inode functions */ 2503a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 25186bacbd2SBarry Smith 252431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 253b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af); 254b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af); 255b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af); 256b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n"); 257431e710aSBarry Smith 25886bacbd2SBarry Smith PetscFunctionReturn(0); 25960ec86bdSBarry Smith #endif 26086bacbd2SBarry Smith } 26186bacbd2SBarry Smith 262289bc588SBarry Smith /* 263289bc588SBarry Smith Factorization code for AIJ format. 264289bc588SBarry Smith */ 2654a2ae208SSatish Balay #undef __FUNCT__ 266b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 267b380c88cSHong Zhang int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *B) 268289bc588SBarry Smith { 269416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 270289bc588SBarry Smith IS isicol; 271273d9f13SBarry Smith int *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j; 272010a20caSHong Zhang int *ainew,*ajnew,jmax,*fill,*ajtmp,nz; 27391bf9adeSBarry Smith int *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im; 2749e046878SBarry Smith PetscReal f; 275289bc588SBarry Smith 2763a40ed3dSBarry Smith PetscFunctionBegin; 27729bbc08cSBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2784c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 279ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 280ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 281289bc588SBarry Smith 282289bc588SBarry Smith /* get new row pointers */ 283b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 284010a20caSHong Zhang ainew[0] = 0; 285289bc588SBarry Smith /* don't know how many column pointers are needed so estimate */ 286d3d32019SBarry Smith f = info->fill; 287010a20caSHong Zhang jmax = (int)(f*ai[n]+1); 288b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 289289bc588SBarry Smith /* fill is a linked list of nonzeros in active row */ 290b0a32e0cSBarry Smith ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr); 2912fb3ed76SBarry Smith im = fill + n + 1; 292289bc588SBarry Smith /* idnew is location of diagonal in factor */ 293b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr); 294010a20caSHong Zhang idnew[0] = 0; 295289bc588SBarry Smith 296289bc588SBarry Smith for (i=0; i<n; i++) { 297289bc588SBarry Smith /* first copy previous fill into linked list */ 298289bc588SBarry Smith nnz = nz = ai[r[i]+1] - ai[r[i]]; 29929bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 300010a20caSHong Zhang ajtmp = aj + ai[r[i]]; 301289bc588SBarry Smith fill[n] = n; 302289bc588SBarry Smith while (nz--) { 303289bc588SBarry Smith fm = n; 304010a20caSHong Zhang idx = ic[*ajtmp++]; 305289bc588SBarry Smith do { 306289bc588SBarry Smith m = fm; 307289bc588SBarry Smith fm = fill[m]; 308289bc588SBarry Smith } while (fm < idx); 309289bc588SBarry Smith fill[m] = idx; 310289bc588SBarry Smith fill[idx] = fm; 311289bc588SBarry Smith } 312289bc588SBarry Smith row = fill[n]; 313289bc588SBarry Smith while (row < i) { 314010a20caSHong Zhang ajtmp = ajnew + idnew[row] + 1; 3152fb3ed76SBarry Smith nzbd = 1 + idnew[row] - ainew[row]; 3162fb3ed76SBarry Smith nz = im[row] - nzbd; 317289bc588SBarry Smith fm = row; 3182fb3ed76SBarry Smith while (nz-- > 0) { 319010a20caSHong Zhang idx = *ajtmp++ ; 3202fb3ed76SBarry Smith nzbd++; 3212fb3ed76SBarry Smith if (idx == i) im[row] = nzbd; 322289bc588SBarry Smith do { 323289bc588SBarry Smith m = fm; 324289bc588SBarry Smith fm = fill[m]; 325289bc588SBarry Smith } while (fm < idx); 326289bc588SBarry Smith if (fm != idx) { 327289bc588SBarry Smith fill[m] = idx; 328289bc588SBarry Smith fill[idx] = fm; 329289bc588SBarry Smith fm = idx; 330289bc588SBarry Smith nnz++; 331289bc588SBarry Smith } 332289bc588SBarry Smith } 333289bc588SBarry Smith row = fill[row]; 334289bc588SBarry Smith } 335289bc588SBarry Smith /* copy new filled row into permanent storage */ 336289bc588SBarry Smith ainew[i+1] = ainew[i] + nnz; 337496e697eSBarry Smith if (ainew[i+1] > jmax) { 338ecf371e4SBarry Smith 339ecf371e4SBarry Smith /* estimate how much additional space we will need */ 340ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 341ecf371e4SBarry Smith /* just double the memory each time */ 342ecf371e4SBarry Smith int maxadd = jmax; 343ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */ 344bbb6d6a8SBarry Smith if (maxadd < nnz) maxadd = (n-i)*(nnz+1); 3452fb3ed76SBarry Smith jmax += maxadd; 346ecf371e4SBarry Smith 347ecf371e4SBarry Smith /* allocate a longer ajnew */ 348b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr); 349010a20caSHong Zhang ierr = PetscMemcpy(ajtmp,ajnew,(ainew[i])*sizeof(int));CHKERRQ(ierr); 350606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 351289bc588SBarry Smith ajnew = ajtmp; 3522fb3ed76SBarry Smith realloc++; /* count how many times we realloc */ 353289bc588SBarry Smith } 354010a20caSHong Zhang ajtmp = ajnew + ainew[i]; 355289bc588SBarry Smith fm = fill[n]; 356289bc588SBarry Smith nzi = 0; 3572fb3ed76SBarry Smith im[i] = nnz; 358289bc588SBarry Smith while (nnz--) { 359289bc588SBarry Smith if (fm < i) nzi++; 360010a20caSHong Zhang *ajtmp++ = fm ; 361289bc588SBarry Smith fm = fill[fm]; 362289bc588SBarry Smith } 363289bc588SBarry Smith idnew[i] = ainew[i] + nzi; 364289bc588SBarry Smith } 365464e8d44SSatish Balay if (ai[n] != 0) { 366329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 367b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 368b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af); 369b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af); 370b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"); 37105bf559cSSatish Balay } else { 372b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"); 37305bf559cSSatish Balay } 3742fb3ed76SBarry Smith 375898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 376898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3771987afe7SBarry Smith 378606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 379289bc588SBarry Smith 380289bc588SBarry Smith /* put together the new matrix */ 381*f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,B);CHKERRQ(ierr); 382*f204ca49SKris Buschelman ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); 383*f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*B,0,PETSC_NULL);CHKERRQ(ierr); 384b0a32e0cSBarry Smith PetscLogObjectParent(*B,isicol); 385416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 386606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3877c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 388e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 389606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 390606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 391010a20caSHong Zhang ierr = PetscMalloc((ainew[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 392416022c9SBarry Smith b->j = ajnew; 393416022c9SBarry Smith b->i = ainew; 394416022c9SBarry Smith b->diag = idnew; 395416022c9SBarry Smith b->ilen = 0; 396416022c9SBarry Smith b->imax = 0; 397416022c9SBarry Smith b->row = isrow; 398416022c9SBarry Smith b->col = iscol; 399085256b3SBarry Smith b->lu_damping = info->damping; 40087828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 4012cea7109SBarry Smith b->lu_shift = info->shift; 4022cea7109SBarry Smith b->lu_shift_fraction = info->shift_fraction; 403c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 404c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 40582bf6240SBarry Smith b->icol = isicol; 40687828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 407416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 4087fd98bd6SLois Curfman McInnes Allocate idnew, solve_work, new a, new j */ 409010a20caSHong Zhang PetscLogObjectMemory(*B,(ainew[n]-n)*(sizeof(int)+sizeof(PetscScalar))); 410010a20caSHong Zhang b->maxnz = b->nz = ainew[n] ; 4117fd98bd6SLois Curfman McInnes 412329f5518SBarry Smith (*B)->factor = FACTOR_LU; 413ae068f58SLois Curfman McInnes (*B)->info.factor_mallocs = realloc; 414ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 4153a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*B,PETSC_FALSE);CHKERRQ(ierr); 4167dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 417703deb49SSatish Balay 418e93ccc4dSBarry Smith if (ai[n] != 0) { 419329f5518SBarry Smith (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 420eea2913aSSatish Balay } else { 421eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 422eea2913aSSatish Balay } 4233a40ed3dSBarry Smith PetscFunctionReturn(0); 424289bc588SBarry Smith } 4250a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 4263a7fca6bSBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth); 42741c01911SSatish Balay 4284a2ae208SSatish Balay #undef __FUNCT__ 4294a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 430416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B) 431289bc588SBarry Smith { 432416022c9SBarry Smith Mat C = *B; 433416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data; 43482bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 435273d9f13SBarry Smith int *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j; 436010a20caSHong Zhang int *ajtmpold,*ajtmp,nz,row,*ics; 4370cf777b8SBarry Smith int *diag_offset = b->diag,diag,*pj,ndamp = 0, nshift=0; 43887828ca2SBarry Smith PetscScalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 4390cf777b8SBarry Smith PetscReal damping = b->lu_damping, zeropivot = b->lu_zeropivot,rs,d; 4400cf777b8SBarry Smith PetscReal row_shift,shift_fraction,shift_amount,shift_lo=0., shift_hi=1., shift_top=0.; 4410cf777b8SBarry Smith PetscTruth damp,lushift; 442289bc588SBarry Smith 4433a40ed3dSBarry Smith PetscFunctionBegin; 44478b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 44578b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 44687828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 44787828ca2SBarry Smith ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 448010a20caSHong Zhang rtmps = rtmp; ics = ic; 449289bc588SBarry Smith 4506cc28720Svictorle if (!a->diag) { 4510cf777b8SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A); CHKERRQ(ierr); 4520cf777b8SBarry Smith } 4530cf777b8SBarry Smith 4540cf777b8SBarry Smith if (b->lu_shift) { /* set max shift */ 4550cf777b8SBarry Smith int *aai = a->i,*ddiag = a->diag; 4560cf777b8SBarry Smith shift_top = 0; 4576cc28720Svictorle for (i=0; i<n; i++) { 458010a20caSHong Zhang d = PetscAbsScalar((a->a)[ddiag[i]]); 4596cc28720Svictorle /* calculate amt of shift needed for this row */ 4606c037d1bSvictorle if (d<=0) { 4613aeef889SHong Zhang row_shift = 0; 4623aeef889SHong Zhang } else { 4633aeef889SHong Zhang row_shift = -2*d; 4643aeef889SHong Zhang } 465010a20caSHong Zhang v = a->a+aai[i]; 4668a2e2d9cSvictorle for (j=0; j<aai[i+1]-aai[i]; j++) 4676cc28720Svictorle row_shift += PetscAbsScalar(v[j]); 4686cc28720Svictorle if (row_shift>shift_top) shift_top = row_shift; 4696cc28720Svictorle } 4706cc28720Svictorle } 4716cc28720Svictorle 4726cc28720Svictorle shift_fraction = 0; shift_amount = 0; 473085256b3SBarry Smith do { 4748a5eb818SBarry Smith damp = PETSC_FALSE; 4756cc28720Svictorle lushift = PETSC_FALSE; 476289bc588SBarry Smith for (i=0; i<n; i++) { 477289bc588SBarry Smith nz = ai[i+1] - ai[i]; 478010a20caSHong Zhang ajtmp = aj + ai[i]; 47948e94559SBarry Smith for (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0; 480289bc588SBarry Smith 481289bc588SBarry Smith /* load in initial (unfactored row) */ 482416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 483010a20caSHong Zhang ajtmpold = a->j + a->i[r[i]]; 484010a20caSHong Zhang v = a->a + a->i[r[i]]; 485085256b3SBarry Smith for (j=0; j<nz; j++) { 486085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] = v[j]; 487085256b3SBarry Smith } 488e2dd4fc4Svictorle rtmp[ics[r[i]]] += damping + shift_amount; /* damp the diagonal of the matrix */ 489289bc588SBarry Smith 490010a20caSHong Zhang row = *ajtmp++ ; 491f2582111SSatish Balay while (row < i) { 4928c37ef55SBarry Smith pc = rtmp + row; 4938c37ef55SBarry Smith if (*pc != 0.0) { 494010a20caSHong Zhang pv = b->a + diag_offset[row] ; 495010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 49635aab85fSBarry Smith multiplier = *pc / *pv++; 4978c37ef55SBarry Smith *pc = multiplier; 498f2582111SSatish Balay nz = ai[row+1] - diag_offset[row] - 1; 499f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 500b0a32e0cSBarry Smith PetscLogFlops(2*nz); 501289bc588SBarry Smith } 502010a20caSHong Zhang row = *ajtmp++; 503289bc588SBarry Smith } 504416022c9SBarry Smith /* finished row so stick it into b->a */ 505010a20caSHong Zhang pv = b->a + ai[i] ; 506010a20caSHong Zhang pj = b->j + ai[i] ; 5078c37ef55SBarry Smith nz = ai[i+1] - ai[i]; 50835aab85fSBarry Smith diag = diag_offset[i] - ai[i]; 509d3d32019SBarry Smith /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 510d3d32019SBarry Smith rs = 0.0; 511d3d32019SBarry Smith for (j=0; j<nz; j++) { 512d3d32019SBarry Smith pv[j] = rtmps[pj[j]]; 513d3d32019SBarry Smith if (j != diag) rs += PetscAbsScalar(pv[j]); 514d3d32019SBarry Smith } 5156cc28720Svictorle #define MAX_NSHIFT 5 516b2da817dSvictorle if (PetscRealPart(pv[diag]) <= zeropivot*rs && b->lu_shift) { 5176cc28720Svictorle if (nshift>MAX_NSHIFT) { 5180cf777b8SBarry Smith SETERRQ(1,"Unable to determine shift to enforce positive definite preconditioner"); 5196cc28720Svictorle } else if (nshift==MAX_NSHIFT) { 5206cc28720Svictorle shift_fraction = shift_hi; 5216c037d1bSvictorle lushift = PETSC_FALSE; 5226cc28720Svictorle } else { 5236cc28720Svictorle shift_lo = shift_fraction; shift_fraction = (shift_hi+shift_lo)/2.; 5246c037d1bSvictorle lushift = PETSC_TRUE; 5256cc28720Svictorle } 5266cc28720Svictorle shift_amount = shift_fraction * shift_top; 5270cf777b8SBarry Smith nshift++; 5280cf777b8SBarry Smith break; 5296cc28720Svictorle } 5309e29f15eSvictorle if (PetscAbsScalar(pv[diag]) <= zeropivot*rs) { 531d3d32019SBarry Smith if (damping) { 5328a5eb818SBarry Smith if (ndamp) damping *= 2.0; 533085256b3SBarry Smith damp = PETSC_TRUE; 534085256b3SBarry Smith ndamp++; 535085256b3SBarry Smith break; 536085256b3SBarry Smith } else { 53763bb2aa6SBarry Smith SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d value %g tolerance %g * rs %g",i,PetscAbsScalar(pv[diag]),zeropivot,rs); 538d708749eSBarry Smith } 53935aab85fSBarry Smith } 54035aab85fSBarry Smith } 5416cc28720Svictorle if (!lushift && b->lu_shift && shift_fraction>0 && nshift<MAX_NSHIFT) { 5426cc28720Svictorle /* 5436c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 5446cc28720Svictorle * then try lower shift 5456cc28720Svictorle */ 5460cf777b8SBarry Smith shift_hi = shift_fraction; 5470cf777b8SBarry Smith shift_fraction = (shift_hi+shift_lo)/2.; 5486cc28720Svictorle shift_amount = shift_fraction * shift_top; 5490cf777b8SBarry Smith lushift = PETSC_TRUE; 5500cf777b8SBarry Smith nshift++; 5516cc28720Svictorle } 5526cc28720Svictorle } while (damp || lushift); 5530f11f9aeSSatish Balay 55435aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 55535aab85fSBarry Smith for (i=0; i<n; i++) { 556010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 55735aab85fSBarry Smith } 55835aab85fSBarry Smith 559606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 56078b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 56178b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 562416022c9SBarry Smith C->factor = FACTOR_LU; 5637dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 564c456f294SBarry Smith C->assembled = PETSC_TRUE; 565b0a32e0cSBarry Smith PetscLogFlops(C->n); 566d3d32019SBarry Smith if (ndamp) { 567b0a32e0cSBarry Smith PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping); 568085256b3SBarry Smith } 5696cc28720Svictorle if (nshift) { 5706cc28720Svictorle b->lu_shift_fraction = shift_fraction; 571fb10cecfSBarry Smith PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: diagonal shifted up by %e fraction top_value %e number shifts %d\n",shift_fraction,shift_top,nshift); 5726cc28720Svictorle } 5733a40ed3dSBarry Smith PetscFunctionReturn(0); 574289bc588SBarry Smith } 5750889b5dcSvictorle 5760889b5dcSvictorle #undef __FUNCT__ 5770889b5dcSvictorle #define __FUNCT__ "MatUsePETSc_SeqAIJ" 5780889b5dcSvictorle int MatUsePETSc_SeqAIJ(Mat A) 5790889b5dcSvictorle { 5800889b5dcSvictorle PetscFunctionBegin; 5810889b5dcSvictorle A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 5820889b5dcSvictorle A->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 5830889b5dcSvictorle PetscFunctionReturn(0); 5840889b5dcSvictorle } 5850889b5dcSvictorle 5860889b5dcSvictorle 5870a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5884a2ae208SSatish Balay #undef __FUNCT__ 5894a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 590b380c88cSHong Zhang int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 591da3a660dSBarry Smith { 592273d9f13SBarry Smith int ierr; 593416022c9SBarry Smith Mat C; 594416022c9SBarry Smith 5953a40ed3dSBarry Smith PetscFunctionBegin; 5969e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 597f2582111SSatish Balay ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr); 598273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 599440f4c53SSatish Balay PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol); 6003a40ed3dSBarry Smith PetscFunctionReturn(0); 601da3a660dSBarry Smith } 6020a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 6034a2ae208SSatish Balay #undef __FUNCT__ 6044a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 605416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 6068c37ef55SBarry Smith { 607416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 608416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 609273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 610010a20caSHong Zhang int nz,*rout,*cout; 61187828ca2SBarry Smith PetscScalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 6128c37ef55SBarry Smith 6133a40ed3dSBarry Smith PetscFunctionBegin; 6143a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 61591bf9adeSBarry Smith 616b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 617b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 618416022c9SBarry Smith tmp = a->solve_work; 6198c37ef55SBarry Smith 6203b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6213b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6228c37ef55SBarry Smith 6238c37ef55SBarry Smith /* forward solve the lower triangular */ 6248c37ef55SBarry Smith tmp[0] = b[*r++]; 625010a20caSHong Zhang tmps = tmp; 6268c37ef55SBarry Smith for (i=1; i<n; i++) { 627010a20caSHong Zhang v = aa + ai[i] ; 628010a20caSHong Zhang vi = aj + ai[i] ; 62953e7425aSBarry Smith nz = a->diag[i] - ai[i]; 63053e7425aSBarry Smith sum = b[*r++]; 631ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6328c37ef55SBarry Smith tmp[i] = sum; 6338c37ef55SBarry Smith } 6348c37ef55SBarry Smith 6358c37ef55SBarry Smith /* backward solve the upper triangular */ 6368c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 637010a20caSHong Zhang v = aa + a->diag[i] + 1; 638010a20caSHong Zhang vi = aj + a->diag[i] + 1; 639416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6408c37ef55SBarry Smith sum = tmp[i]; 641ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 642010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6438c37ef55SBarry Smith } 6448c37ef55SBarry Smith 6453b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6463b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 647b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 648b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 649b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 6503a40ed3dSBarry Smith PetscFunctionReturn(0); 6518c37ef55SBarry Smith } 652026e39d0SSatish Balay 653930ae53cSBarry Smith /* ----------------------------------------------------------- */ 6544a2ae208SSatish Balay #undef __FUNCT__ 6554a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 656930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 657930ae53cSBarry Smith { 658930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 659273d9f13SBarry Smith int n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr; 660362ced78SSatish Balay PetscScalar *x,*b,*aa = a->a; 661aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 662d85afc42SSatish Balay int adiag_i,i,*vi,nz,ai_i; 663362ced78SSatish Balay PetscScalar *v,sum; 664d85afc42SSatish Balay #endif 665930ae53cSBarry Smith 6663a40ed3dSBarry Smith PetscFunctionBegin; 6673a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 668930ae53cSBarry Smith 669b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 670b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 671930ae53cSBarry Smith 672aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6731c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6741c4feecaSSatish Balay #else 675930ae53cSBarry Smith /* forward solve the lower triangular */ 676930ae53cSBarry Smith x[0] = b[0]; 677930ae53cSBarry Smith for (i=1; i<n; i++) { 678930ae53cSBarry Smith ai_i = ai[i]; 679930ae53cSBarry Smith v = aa + ai_i; 680930ae53cSBarry Smith vi = aj + ai_i; 681930ae53cSBarry Smith nz = adiag[i] - ai_i; 682930ae53cSBarry Smith sum = b[i]; 683930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 684930ae53cSBarry Smith x[i] = sum; 685930ae53cSBarry Smith } 686930ae53cSBarry Smith 687930ae53cSBarry Smith /* backward solve the upper triangular */ 688930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 689930ae53cSBarry Smith adiag_i = adiag[i]; 690d708749eSBarry Smith v = aa + adiag_i + 1; 691d708749eSBarry Smith vi = aj + adiag_i + 1; 692930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 693930ae53cSBarry Smith sum = x[i]; 694930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 695930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 696930ae53cSBarry Smith } 6971c4feecaSSatish Balay #endif 698b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 699b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 700b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 7013a40ed3dSBarry Smith PetscFunctionReturn(0); 702930ae53cSBarry Smith } 703930ae53cSBarry Smith 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 706416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 707da3a660dSBarry Smith { 708416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 709416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 710273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 711010a20caSHong Zhang int nz,*rout,*cout; 71287828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,sum,*v; 713da3a660dSBarry Smith 7143a40ed3dSBarry Smith PetscFunctionBegin; 71578b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 716da3a660dSBarry Smith 717b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 718b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 719416022c9SBarry Smith tmp = a->solve_work; 720da3a660dSBarry Smith 7213b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7223b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 723da3a660dSBarry Smith 724da3a660dSBarry Smith /* forward solve the lower triangular */ 725da3a660dSBarry Smith tmp[0] = b[*r++]; 726da3a660dSBarry Smith for (i=1; i<n; i++) { 727010a20caSHong Zhang v = aa + ai[i] ; 728010a20caSHong Zhang vi = aj + ai[i] ; 729416022c9SBarry Smith nz = a->diag[i] - ai[i]; 730da3a660dSBarry Smith sum = b[*r++]; 731010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 732da3a660dSBarry Smith tmp[i] = sum; 733da3a660dSBarry Smith } 734da3a660dSBarry Smith 735da3a660dSBarry Smith /* backward solve the upper triangular */ 736da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 737010a20caSHong Zhang v = aa + a->diag[i] + 1; 738010a20caSHong Zhang vi = aj + a->diag[i] + 1; 739416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 740da3a660dSBarry Smith sum = tmp[i]; 741010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 742010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 743da3a660dSBarry Smith x[*c--] += tmp[i]; 744da3a660dSBarry Smith } 745da3a660dSBarry Smith 7463b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7473b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 748b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 749b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 750b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 751898183ecSLois Curfman McInnes 7523a40ed3dSBarry Smith PetscFunctionReturn(0); 753da3a660dSBarry Smith } 754da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7554a2ae208SSatish Balay #undef __FUNCT__ 7564a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 7577c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 758da3a660dSBarry Smith { 759416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7602235e667SBarry Smith IS iscol = a->col,isrow = a->row; 761273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 762010a20caSHong Zhang int nz,*rout,*cout,*diag = a->diag; 76387828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v,s1; 764da3a660dSBarry Smith 7653a40ed3dSBarry Smith PetscFunctionBegin; 766b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 767b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 768416022c9SBarry Smith tmp = a->solve_work; 769da3a660dSBarry Smith 7702235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7712235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 772da3a660dSBarry Smith 773da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7742235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 775da3a660dSBarry Smith 776da3a660dSBarry Smith /* forward solve the U^T */ 777da3a660dSBarry Smith for (i=0; i<n; i++) { 778010a20caSHong Zhang v = aa + diag[i] ; 779010a20caSHong Zhang vi = aj + diag[i] + 1; 780f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 781c38d4ed2SBarry Smith s1 = tmp[i]; 782ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 783da3a660dSBarry Smith while (nz--) { 784010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 785da3a660dSBarry Smith } 786c38d4ed2SBarry Smith tmp[i] = s1; 787da3a660dSBarry Smith } 788da3a660dSBarry Smith 789da3a660dSBarry Smith /* backward solve the L^T */ 790da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 791010a20caSHong Zhang v = aa + diag[i] - 1 ; 792010a20caSHong Zhang vi = aj + diag[i] - 1 ; 793f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 794f1af5d2fSBarry Smith s1 = tmp[i]; 795da3a660dSBarry Smith while (nz--) { 796010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 797da3a660dSBarry Smith } 798da3a660dSBarry Smith } 799da3a660dSBarry Smith 800da3a660dSBarry Smith /* copy tmp into x according to permutation */ 801da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 802da3a660dSBarry Smith 8032235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8042235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 805b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 806b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 807da3a660dSBarry Smith 808b0a32e0cSBarry Smith PetscLogFlops(2*a->nz-A->n); 8093a40ed3dSBarry Smith PetscFunctionReturn(0); 810da3a660dSBarry Smith } 811da3a660dSBarry Smith 8124a2ae208SSatish Balay #undef __FUNCT__ 8134a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 8147c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 815da3a660dSBarry Smith { 816416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8172235e667SBarry Smith IS iscol = a->col,isrow = a->row; 818273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 819010a20caSHong Zhang int nz,*rout,*cout,*diag = a->diag; 82087828ca2SBarry Smith PetscScalar *x,*b,*tmp,*aa = a->a,*v; 8216abc6512SBarry Smith 8223a40ed3dSBarry Smith PetscFunctionBegin; 8236abc6512SBarry Smith if (zz != xx) VecCopy(zz,xx); 8246abc6512SBarry Smith 825b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 826b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 827416022c9SBarry Smith tmp = a->solve_work; 8286abc6512SBarry Smith 8292235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8302235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 8316abc6512SBarry Smith 8326abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 8332235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 8346abc6512SBarry Smith 8356abc6512SBarry Smith /* forward solve the U^T */ 8366abc6512SBarry Smith for (i=0; i<n; i++) { 837010a20caSHong Zhang v = aa + diag[i] ; 838010a20caSHong Zhang vi = aj + diag[i] + 1; 839f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 8406abc6512SBarry Smith tmp[i] *= *v++; 8416abc6512SBarry Smith while (nz--) { 842010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 8436abc6512SBarry Smith } 8446abc6512SBarry Smith } 8456abc6512SBarry Smith 8466abc6512SBarry Smith /* backward solve the L^T */ 8476abc6512SBarry Smith for (i=n-1; i>=0; i--){ 848010a20caSHong Zhang v = aa + diag[i] - 1 ; 849010a20caSHong Zhang vi = aj + diag[i] - 1 ; 850f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 8516abc6512SBarry Smith while (nz--) { 852010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 8536abc6512SBarry Smith } 8546abc6512SBarry Smith } 8556abc6512SBarry Smith 8566abc6512SBarry Smith /* copy tmp into x according to permutation */ 8576abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8586abc6512SBarry Smith 8592235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8602235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 861b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 862b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 8636abc6512SBarry Smith 864b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 8653a40ed3dSBarry Smith PetscFunctionReturn(0); 866da3a660dSBarry Smith } 8679e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 868ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat); 86908480c60SBarry Smith 8704a2ae208SSatish Balay #undef __FUNCT__ 8714a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 872b380c88cSHong Zhang int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 8739e25ed09SBarry Smith { 874416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8759e25ed09SBarry Smith IS isicol; 876273d9f13SBarry Smith int *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j; 87756beaf04SBarry Smith int *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev; 878335d9088SBarry Smith int *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0; 879010a20caSHong Zhang int incrlev,nnz,i,levels,diagonal_fill; 88077c4ece6SBarry Smith PetscTruth col_identity,row_identity; 881329f5518SBarry Smith PetscReal f; 8829e25ed09SBarry Smith 8833a40ed3dSBarry Smith PetscFunctionBegin; 8846cf997b0SBarry Smith f = info->fill; 8850cd17293SBarry Smith levels = (int)info->levels; 8860cd17293SBarry Smith diagonal_fill = (int)info->diagonal_fill; 8874c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 88882bf6240SBarry Smith 88908480c60SBarry Smith /* special case that simply copies fill pattern */ 890be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 891be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 89286bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8932e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 89408480c60SBarry Smith (*fact)->factor = FACTOR_LU; 89508480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 89608480c60SBarry Smith if (!b->diag) { 8977c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 89808480c60SBarry Smith } 8997c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 90008480c60SBarry Smith b->row = isrow; 90108480c60SBarry Smith b->col = iscol; 90282bf6240SBarry Smith b->icol = isicol; 903a2e34c3dSBarry Smith b->lu_damping = info->damping; 90487828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 9052cea7109SBarry Smith b->lu_shift = info->shift; 9062cea7109SBarry Smith b->lu_shift_fraction= info->shift_fraction; 90787828ca2SBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 908f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 909c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 910c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 9113a40ed3dSBarry Smith PetscFunctionReturn(0); 91208480c60SBarry Smith } 91308480c60SBarry Smith 914898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 915898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 9169e25ed09SBarry Smith 9179e25ed09SBarry Smith /* get new row pointers */ 918b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 919010a20caSHong Zhang ainew[0] = 0; 9209e25ed09SBarry Smith /* don't know how many column pointers are needed so estimate */ 921010a20caSHong Zhang jmax = (int)(f*(ai[n]+1)); 922b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 9239e25ed09SBarry Smith /* ajfill is level of fill for each fill entry */ 924b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr); 9259e25ed09SBarry Smith /* fill is a linked list of nonzeros in active row */ 926b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr); 92756beaf04SBarry Smith /* im is level for each filled value */ 928b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr); 92956beaf04SBarry Smith /* dloc is location of diagonal in factor */ 930b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr); 93156beaf04SBarry Smith dloc[0] = 0; 93256beaf04SBarry Smith for (prow=0; prow<n; prow++) { 9336cf997b0SBarry Smith 9346cf997b0SBarry Smith /* copy current row into linked list */ 93556beaf04SBarry Smith nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 93629bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 937010a20caSHong Zhang xi = aj + ai[r[prow]] ; 9389e25ed09SBarry Smith fill[n] = n; 939435faa5fSBarry Smith fill[prow] = -1; /* marker to indicate if diagonal exists */ 9409e25ed09SBarry Smith while (nz--) { 9419e25ed09SBarry Smith fm = n; 942010a20caSHong Zhang idx = ic[*xi++ ]; 9439e25ed09SBarry Smith do { 9449e25ed09SBarry Smith m = fm; 9459e25ed09SBarry Smith fm = fill[m]; 9469e25ed09SBarry Smith } while (fm < idx); 9479e25ed09SBarry Smith fill[m] = idx; 9489e25ed09SBarry Smith fill[idx] = fm; 94956beaf04SBarry Smith im[idx] = 0; 9509e25ed09SBarry Smith } 9516cf997b0SBarry Smith 9526cf997b0SBarry Smith /* make sure diagonal entry is included */ 953435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 9546cf997b0SBarry Smith fm = n; 955435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 956435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 957435faa5fSBarry Smith fill[fm] = prow; 9586cf997b0SBarry Smith im[prow] = 0; 9596cf997b0SBarry Smith nzf++; 960335d9088SBarry Smith dcount++; 9616cf997b0SBarry Smith } 9626cf997b0SBarry Smith 96356beaf04SBarry Smith nzi = 0; 9649e25ed09SBarry Smith row = fill[n]; 96556beaf04SBarry Smith while (row < prow) { 96656beaf04SBarry Smith incrlev = im[row] + 1; 96756beaf04SBarry Smith nz = dloc[row]; 968010a20caSHong Zhang xi = ajnew + ainew[row] + nz + 1; 969010a20caSHong Zhang flev = ajfill + ainew[row] + nz + 1; 970416022c9SBarry Smith nnz = ainew[row+1] - ainew[row] - nz - 1; 97156beaf04SBarry Smith fm = row; 97256beaf04SBarry Smith while (nnz-- > 0) { 973010a20caSHong Zhang idx = *xi++ ; 97456beaf04SBarry Smith if (*flev + incrlev > levels) { 97556beaf04SBarry Smith flev++; 97656beaf04SBarry Smith continue; 97756beaf04SBarry Smith } 97856beaf04SBarry Smith do { 9799e25ed09SBarry Smith m = fm; 9809e25ed09SBarry Smith fm = fill[m]; 98156beaf04SBarry Smith } while (fm < idx); 9829e25ed09SBarry Smith if (fm != idx) { 98356beaf04SBarry Smith im[idx] = *flev + incrlev; 9849e25ed09SBarry Smith fill[m] = idx; 9859e25ed09SBarry Smith fill[idx] = fm; 9869e25ed09SBarry Smith fm = idx; 98756beaf04SBarry Smith nzf++; 988ecf371e4SBarry Smith } else { 98956beaf04SBarry Smith if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 9909e25ed09SBarry Smith } 99156beaf04SBarry Smith flev++; 9929e25ed09SBarry Smith } 9939e25ed09SBarry Smith row = fill[row]; 99456beaf04SBarry Smith nzi++; 9959e25ed09SBarry Smith } 9969e25ed09SBarry Smith /* copy new filled row into permanent storage */ 99756beaf04SBarry Smith ainew[prow+1] = ainew[prow] + nzf; 998010a20caSHong Zhang if (ainew[prow+1] > jmax) { 999ecf371e4SBarry Smith 1000ecf371e4SBarry Smith /* estimate how much additional space we will need */ 1001ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 1002ecf371e4SBarry Smith /* just double the memory each time */ 1003ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */ 1004ecf371e4SBarry Smith int maxadd = jmax; 100556beaf04SBarry Smith if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 1006bbb6d6a8SBarry Smith jmax += maxadd; 1007ecf371e4SBarry Smith 1008ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 1009b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 1010010a20caSHong Zhang ierr = PetscMemcpy(xi,ajnew,(ainew[prow])*sizeof(int));CHKERRQ(ierr); 1011606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 101256beaf04SBarry Smith ajnew = xi; 1013b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 1014010a20caSHong Zhang ierr = PetscMemcpy(xi,ajfill,(ainew[prow])*sizeof(int));CHKERRQ(ierr); 1015606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 101656beaf04SBarry Smith ajfill = xi; 1017ecf371e4SBarry Smith realloc++; /* count how many times we realloc */ 10189e25ed09SBarry Smith } 1019010a20caSHong Zhang xi = ajnew + ainew[prow] ; 1020010a20caSHong Zhang flev = ajfill + ainew[prow] ; 102156beaf04SBarry Smith dloc[prow] = nzi; 10229e25ed09SBarry Smith fm = fill[n]; 102356beaf04SBarry Smith while (nzf--) { 1024010a20caSHong Zhang *xi++ = fm ; 102556beaf04SBarry Smith *flev++ = im[fm]; 10269e25ed09SBarry Smith fm = fill[fm]; 10279e25ed09SBarry Smith } 10286cf997b0SBarry Smith /* make sure row has diagonal entry */ 1029010a20caSHong Zhang if (ajnew[ainew[prow]+dloc[prow]] != prow) { 103029bbc08cSBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\ 10316cf997b0SBarry Smith try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow); 10326cf997b0SBarry Smith } 10339e25ed09SBarry Smith } 1034606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 1035898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1036898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1037606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 1038606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 10399e25ed09SBarry Smith 1040f64065bbSBarry Smith { 1041329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 1042b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 1043c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -[sub_]pc_ilu_fill %g or use \n",af); 1044c0d6ac57SBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill([sub]pc,%g);\n",af); 1045b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"); 1046335d9088SBarry Smith if (diagonal_fill) { 1047b1bcba4aSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replaced %d missing diagonals",dcount); 1048335d9088SBarry Smith } 1049f64065bbSBarry Smith } 1050bbb6d6a8SBarry Smith 10519e25ed09SBarry Smith /* put together the new matrix */ 1052*f204ca49SKris Buschelman ierr = MatCreate(A->comm,n,n,n,n,fact);CHKERRQ(ierr); 1053*f204ca49SKris Buschelman ierr = MatSetType(*fact,A->type_name);CHKERRQ(ierr); 1054*f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(*fact,0,PETSC_NULL);CHKERRQ(ierr); 1055b0a32e0cSBarry Smith PetscLogObjectParent(*fact,isicol); 1056416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1057606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10587c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 1059010a20caSHong Zhang len = (ainew[n] )*sizeof(PetscScalar); 10609e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1061606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1062606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1063b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 1064416022c9SBarry Smith b->j = ajnew; 1065416022c9SBarry Smith b->i = ainew; 106656beaf04SBarry Smith for (i=0; i<n; i++) dloc[i] += ainew[i]; 1067416022c9SBarry Smith b->diag = dloc; 1068416022c9SBarry Smith b->ilen = 0; 1069416022c9SBarry Smith b->imax = 0; 1070416022c9SBarry Smith b->row = isrow; 1071416022c9SBarry Smith b->col = iscol; 1072c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1073c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 107482bf6240SBarry Smith b->icol = isicol; 107587828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1076416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 107756beaf04SBarry Smith Allocate dloc, solve_work, new a, new j */ 1078010a20caSHong Zhang PetscLogObjectMemory(*fact,(ainew[n]-n) * (sizeof(int)+sizeof(PetscScalar))); 1079010a20caSHong Zhang b->maxnz = b->nz = ainew[n] ; 1080a2e34c3dSBarry Smith b->lu_damping = info->damping; 10812cea7109SBarry Smith b->lu_shift = info->shift; 10822cea7109SBarry Smith b->lu_shift_fraction = info->shift_fraction; 108387828ca2SBarry Smith b->lu_zeropivot = info->zeropivot; 10849e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 10853a7fca6bSBarry Smith ierr = Mat_AIJ_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 10863a7fca6bSBarry Smith (*fact)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 1087ae068f58SLois Curfman McInnes 1088ae068f58SLois Curfman McInnes (*fact)->info.factor_mallocs = realloc; 1089ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 1090329f5518SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 10913a40ed3dSBarry Smith PetscFunctionReturn(0); 10929e25ed09SBarry Smith } 1093d5d45c9bSBarry Smith 1094a6175056SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 1095a6175056SHong Zhang #undef __FUNCT__ 1096a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1097a6175056SHong Zhang int MatCholeskyFactorNumeric_SeqAIJ(Mat A,Mat *fact) 1098a6175056SHong Zhang { 10990968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1100a6175056SHong Zhang int ierr; 1101d5d45c9bSBarry Smith 1102a6175056SHong Zhang PetscFunctionBegin; 11030968510aSHong Zhang if (!a->sbaijMat){ 11040968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1105a6175056SHong Zhang } 110603aac1b8SHong Zhang 1107b45a75daSHong Zhang ierr = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(a->sbaijMat,fact);CHKERRQ(ierr); 11080968510aSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 1109064503c5SHong Zhang a->sbaijMat = PETSC_NULL; 1110653a6975SHong Zhang 1111a6175056SHong Zhang PetscFunctionReturn(0); 1112a6175056SHong Zhang } 1113a6175056SHong Zhang 1114a6175056SHong Zhang #undef __FUNCT__ 1115a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 111615e8a5b3SHong Zhang int MatICCFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1117a6175056SHong Zhang { 11180968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 111903aac1b8SHong Zhang int ierr; 1120653a6975SHong Zhang PetscTruth perm_identity; 1121a6175056SHong Zhang 1122a6175056SHong Zhang PetscFunctionBegin; 1123653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1124653a6975SHong Zhang if (!perm_identity){ 1125653a6975SHong Zhang SETERRQ(1,"Non-identity permutation is not supported yet"); 1126653a6975SHong Zhang } 11270968510aSHong Zhang if (!a->sbaijMat){ 11280968510aSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 11290968510aSHong Zhang } 1130a6175056SHong Zhang 1131b00f7748SHong Zhang ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 1132653a6975SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1133653a6975SHong Zhang 1134a6175056SHong Zhang PetscFunctionReturn(0); 1135a6175056SHong Zhang } 1136d5d45c9bSBarry Smith 1137f76d2b81SHong Zhang #undef __FUNCT__ 1138f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1139f76d2b81SHong Zhang int MatCholeskyFactorSymbolic_SeqAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 1140f76d2b81SHong Zhang { 1141f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 114203aac1b8SHong Zhang int ierr; 1143f76d2b81SHong Zhang PetscTruth perm_identity; 1144f76d2b81SHong Zhang 1145f76d2b81SHong Zhang PetscFunctionBegin; 1146f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1147f76d2b81SHong Zhang if (!perm_identity){ 1148f76d2b81SHong Zhang SETERRQ(1,"Non-identity permutation is not supported yet"); 1149f76d2b81SHong Zhang } 1150f76d2b81SHong Zhang if (!a->sbaijMat){ 1151f76d2b81SHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 1152f76d2b81SHong Zhang } 1153f76d2b81SHong Zhang 1154f76d2b81SHong Zhang ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 1155f76d2b81SHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1156f76d2b81SHong Zhang 1157f76d2b81SHong Zhang PetscFunctionReturn(0); 1158f76d2b81SHong Zhang } 1159