1*b0a32e0cSBarry Smith /*$Id: aijfact.c,v 1.157 2000/10/24 20:25:32 bsmith Exp bsmith $*/ 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" 63b2fbd54SBarry Smith 75615d1e5SSatish Balay #undef __FUNC__ 8*b0a32e0cSBarry Smith #define __FUNC__ "MatOrdering_Flow_SeqAIJ" 991e9ee9fSBarry Smith int MatOrdering_Flow_SeqAIJ(Mat mat,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); 21ca44d042SBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat); 2286bacbd2SBarry Smith 23ca44d042SBarry Smith EXTERN int SPARSEKIT2dperm(int*,Scalar*,int*,int*,Scalar*,int*,int*,int*,int*,int*); 24ca44d042SBarry Smith EXTERN int SPARSEKIT2ilutp(int*,Scalar*,int*,int*,int*,PetscReal*,PetscReal*,int*,Scalar*,int*,int*,int*,Scalar*,int*,int*,int*); 25ca44d042SBarry Smith EXTERN int SPARSEKIT2msrcsr(int*,Scalar*,int*,Scalar*,int*,int*,Scalar*,int*); 2686bacbd2SBarry Smith 2786bacbd2SBarry Smith #undef __FUNC__ 28*b0a32e0cSBarry Smith #define __FUNC__ "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 485bd2ddc7SBarry Smith ishift = 0, for indices start at 1 495bd2ddc7SBarry Smith ishift = 1, for indices starting at 0 5086bacbd2SBarry Smith ------------------------------------------------------------ 5186bacbd2SBarry Smith */ 5286bacbd2SBarry Smith 5386bacbd2SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,MatILUInfo *info,IS isrow,IS iscol,Mat *fact) 5486bacbd2SBarry Smith { 5586bacbd2SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 5607d2056aSBarry Smith IS iscolf,isicol,isirow; 5786bacbd2SBarry Smith PetscTruth reorder; 58273d9f13SBarry Smith int *c,*r,*ic,ierr,i,n = A->m; 59329f5518SBarry Smith int *old_i = a->i,*old_j = a->j,*new_i,*old_i2 = 0,*old_j2 = 0,*new_j; 60313ae024SBarry Smith int *ordcol,*iwk,*iperm,*jw; 615bd2ddc7SBarry Smith int ishift = !a->indexshift; 62b19713cbSBarry Smith int jmax,lfill,job,*o_i,*o_j; 63329f5518SBarry Smith Scalar *old_a = a->a,*w,*new_a,*old_a2 = 0,*wk,*o_a; 64329f5518SBarry Smith PetscReal permtol,af; 6586bacbd2SBarry Smith 6686bacbd2SBarry Smith PetscFunctionBegin; 6786bacbd2SBarry Smith 6886bacbd2SBarry Smith if (info->dt == PETSC_DEFAULT) info->dt = .005; 6986bacbd2SBarry Smith if (info->dtcount == PETSC_DEFAULT) info->dtcount = (int)(1.5*a->rmax); 7086bacbd2SBarry Smith if (info->dtcol == PETSC_DEFAULT) info->dtcol = .01; 7133259db3SBarry Smith if (info->fill == PETSC_DEFAULT) info->fill = ((double)(n*(info->dtcount+1)))/a->nz; 726faa4630SBarry Smith lfill = (int)(info->dtcount/2.0); 736faa4630SBarry Smith jmax = (int)(info->fill*a->nz); 7486bacbd2SBarry Smith permtol = info->dtcol; 7586bacbd2SBarry Smith 7686bacbd2SBarry Smith 7786bacbd2SBarry Smith /* ------------------------------------------------------------ 7886bacbd2SBarry Smith If reorder=.TRUE., then the original matrix has to be 7986bacbd2SBarry Smith reordered to reflect the user selected ordering scheme, and 8086bacbd2SBarry Smith then de-reordered so it is in it's original format. 8186bacbd2SBarry Smith Because Saad's dperm() is NOT in place, we have to copy 8286bacbd2SBarry Smith the original matrix and allocate more storage. . . 8386bacbd2SBarry Smith ------------------------------------------------------------ 8486bacbd2SBarry Smith */ 8586bacbd2SBarry Smith 8686bacbd2SBarry Smith /* set reorder to true if either isrow or iscol is not identity */ 8786bacbd2SBarry Smith ierr = ISIdentity(isrow,&reorder);CHKERRQ(ierr); 8886bacbd2SBarry Smith if (reorder) {ierr = ISIdentity(iscol,&reorder);CHKERRQ(ierr);} 8986bacbd2SBarry Smith reorder = PetscNot(reorder); 9086bacbd2SBarry Smith 9186bacbd2SBarry Smith 9286bacbd2SBarry Smith /* storage for ilu factor */ 93*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&new_i);CHKERRQ(ierr); 94*b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&new_j);CHKERRQ(ierr); 95*b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(Scalar),&new_a);CHKERRQ(ierr); 96*b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(int),&ordcol);CHKERRQ(ierr); 9786bacbd2SBarry Smith 9886bacbd2SBarry Smith /* ------------------------------------------------------------ 9986bacbd2SBarry Smith Make sure that everything is Fortran formatted (1-Based) 10086bacbd2SBarry Smith ------------------------------------------------------------ 10186bacbd2SBarry Smith */ 102b19713cbSBarry Smith if (ishift) { 103b19713cbSBarry Smith for (i=old_i[0];i<old_i[n];i++) { 104b19713cbSBarry Smith old_j[i]++; 10586bacbd2SBarry Smith } 106b19713cbSBarry Smith for(i=0;i<n+1;i++) { 107b19713cbSBarry Smith old_i[i]++; 108b19713cbSBarry Smith }; 10986bacbd2SBarry Smith }; 11086bacbd2SBarry Smith 11186bacbd2SBarry Smith if (reorder) { 112c0c2c81eSBarry Smith ierr = ISGetIndices(iscol,&c); CHKERRQ(ierr); 113c0c2c81eSBarry Smith ierr = ISGetIndices(isrow,&r); CHKERRQ(ierr); 114c0c2c81eSBarry Smith for(i=0;i<n;i++) { 115c0c2c81eSBarry Smith r[i] = r[i]+1; 116c0c2c81eSBarry Smith c[i] = c[i]+1; 117c0c2c81eSBarry Smith } 118*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&old_i2);CHKERRQ(ierr); 119*b0a32e0cSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(int),&old_j2);CHKERRQ(ierr); 120*b0a32e0cSBarry Smith ierr = PetscMalloc((old_i[n]-old_i[0]+1)*sizeof(Scalar),old_a2);CHKERRQ(ierr); 1215bd2ddc7SBarry Smith job = 3; SPARSEKIT2dperm(&n,old_a,old_j,old_i,old_a2,old_j2,old_i2,r,c,&job); 122c0c2c81eSBarry Smith for (i=0;i<n;i++) { 123c0c2c81eSBarry Smith r[i] = r[i]-1; 124c0c2c81eSBarry Smith c[i] = c[i]-1; 125c0c2c81eSBarry Smith } 126c0c2c81eSBarry Smith ierr = ISRestoreIndices(iscol,&c);CHKERRQ(ierr); 127c0c2c81eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 128b19713cbSBarry Smith o_a = old_a2; 129b19713cbSBarry Smith o_j = old_j2; 130b19713cbSBarry Smith o_i = old_i2; 131b19713cbSBarry Smith } else { 132b19713cbSBarry Smith o_a = old_a; 133b19713cbSBarry Smith o_j = old_j; 134b19713cbSBarry Smith o_i = old_i; 13586bacbd2SBarry Smith } 13686bacbd2SBarry Smith 13786bacbd2SBarry Smith /* ------------------------------------------------------------ 13886bacbd2SBarry Smith Call Saad's ilutp() routine to generate the factorization 13986bacbd2SBarry Smith ------------------------------------------------------------ 14086bacbd2SBarry Smith */ 14186bacbd2SBarry Smith 142*b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&iperm);CHKERRQ(ierr); 143*b0a32e0cSBarry Smith ierr = PetscMalloc(2*n*sizeof(int),&jw);CHKERRQ(ierr); 144*b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(Scalar),&w);CHKERRQ(ierr); 14586bacbd2SBarry Smith 1465bd2ddc7SBarry Smith SPARSEKIT2ilutp(&n,o_a,o_j,o_i,&lfill,&info->dt,&permtol,&n,new_a,new_j,new_i,&jmax,w,jw,iperm,&ierr); 14786bacbd2SBarry Smith if (ierr) { 14886bacbd2SBarry Smith switch (ierr) { 14929bbc08cSBarry Smith case -3: SETERRQ2(1,"ilutp(), matrix U overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15029bbc08cSBarry Smith case -2: SETERRQ2(1,"ilutp(), matrix L overflows, need larger info->fill current fill %g space allocated %d",info->fill,jmax); 15129bbc08cSBarry Smith case -5: SETERRQ(1,"ilutp(), zero row encountered"); 15229bbc08cSBarry Smith case -1: SETERRQ(1,"ilutp(), input matrix may be wrong"); 15329bbc08cSBarry Smith case -4: SETERRQ1(1,"ilutp(), illegal info->fill value %d",jmax); 15429bbc08cSBarry Smith default: SETERRQ1(1,"ilutp(), zero pivot detected on row %d",ierr); 15586bacbd2SBarry Smith } 15686bacbd2SBarry Smith } 15786bacbd2SBarry Smith 15886bacbd2SBarry Smith ierr = PetscFree(w);CHKERRQ(ierr); 15986bacbd2SBarry Smith ierr = PetscFree(jw);CHKERRQ(ierr); 16086bacbd2SBarry Smith 16186bacbd2SBarry Smith /* ------------------------------------------------------------ 16286bacbd2SBarry Smith Saad's routine gives the result in Modified Sparse Row (msr) 16386bacbd2SBarry Smith Convert to Compressed Sparse Row format (csr) 16486bacbd2SBarry Smith ------------------------------------------------------------ 16586bacbd2SBarry Smith */ 16686bacbd2SBarry Smith 167*b0a32e0cSBarry Smith ierr = PetscMalloc(n*sizeof(Scalar),&wk);CHKERRQ(ierr); 168*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&iwk);CHKERRQ(ierr); 16986bacbd2SBarry Smith 1705bd2ddc7SBarry Smith SPARSEKIT2msrcsr(&n,new_a,new_j,new_a,new_j,new_i,wk,iwk); 17186bacbd2SBarry Smith 17286bacbd2SBarry Smith ierr = PetscFree(iwk);CHKERRQ(ierr); 17386bacbd2SBarry Smith ierr = PetscFree(wk);CHKERRQ(ierr); 17486bacbd2SBarry Smith 17586bacbd2SBarry Smith if (reorder) { 17686bacbd2SBarry Smith ierr = PetscFree(old_a2);CHKERRQ(ierr); 17786bacbd2SBarry Smith ierr = PetscFree(old_j2);CHKERRQ(ierr); 17886bacbd2SBarry Smith ierr = PetscFree(old_i2);CHKERRQ(ierr); 179313ae024SBarry Smith } else { 180b19713cbSBarry Smith /* fix permutation of old_j that the factorization introduced */ 181f170005cSBarry Smith for (i=old_i[0]; i<old_i[n]; i++) { 182b19713cbSBarry Smith old_j[i-1] = iperm[old_j[i-1]-1]; 183b19713cbSBarry Smith } 184b19713cbSBarry Smith } 185b19713cbSBarry Smith 186b801d0f9SBarry Smith /* get rid of the shift to indices starting at 1 */ 187b19713cbSBarry Smith if (ishift) { 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 } 19586bacbd2SBarry Smith 196b19713cbSBarry Smith /* Make the factored matrix 0-based */ 197b19713cbSBarry Smith if (ishift) { 19886bacbd2SBarry Smith for (i=0; i<n+1; i++) { 199b19713cbSBarry Smith new_i[i]--; 200b19713cbSBarry Smith } 201b19713cbSBarry Smith for (i=new_i[0];i<new_i[n];i++) { 202b19713cbSBarry Smith new_j[i]--; 203b19713cbSBarry Smith } 20486bacbd2SBarry Smith } 20586bacbd2SBarry Smith 20686bacbd2SBarry Smith /*-- due to the pivoting, we need to reorder iscol to correctly --*/ 20786bacbd2SBarry Smith /*-- permute the right-hand-side and solution vectors --*/ 2084c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2094c49b128SBarry Smith ierr = ISInvertPermutation(isrow,PETSC_DECIDE,&isirow);CHKERRQ(ierr); 210c0c2c81eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 21186bacbd2SBarry Smith for(i=0; i<n; i++) { 21286bacbd2SBarry Smith ordcol[i] = ic[iperm[i]-1]; 21386bacbd2SBarry Smith }; 21486bacbd2SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21507d2056aSBarry Smith ierr = ISDestroy(isicol);CHKERRQ(ierr); 21686bacbd2SBarry Smith 217c0c2c81eSBarry Smith ierr = PetscFree(iperm);CHKERRQ(ierr); 218c0c2c81eSBarry Smith 219329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ordcol,&iscolf);CHKERRQ(ierr); 22007d2056aSBarry Smith ierr = PetscFree(ordcol);CHKERRQ(ierr); 22186bacbd2SBarry Smith 22286bacbd2SBarry Smith /*----- put together the new matrix -----*/ 22386bacbd2SBarry Smith 22486bacbd2SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr); 22586bacbd2SBarry Smith (*fact)->factor = FACTOR_LU; 22686bacbd2SBarry Smith (*fact)->assembled = PETSC_TRUE; 22786bacbd2SBarry Smith 22886bacbd2SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 22986bacbd2SBarry Smith ierr = PetscFree(b->imax);CHKERRQ(ierr); 23086bacbd2SBarry Smith b->sorted = PETSC_FALSE; 23107d2056aSBarry Smith b->singlemalloc = PETSC_FALSE; 232b19713cbSBarry Smith /* the next line frees the default space generated by the MatCreate() */ 23386bacbd2SBarry Smith ierr = PetscFree(b->a);CHKERRQ(ierr); 23486bacbd2SBarry Smith ierr = PetscFree(b->ilen);CHKERRQ(ierr); 23586bacbd2SBarry Smith b->a = new_a; 23686bacbd2SBarry Smith b->j = new_j; 23786bacbd2SBarry Smith b->i = new_i; 23886bacbd2SBarry Smith b->ilen = 0; 23986bacbd2SBarry Smith b->imax = 0; 2401f9e874aSBarry Smith /* I am not sure why these are the inverses of the row and column permutations; but the other way is NO GOOD */ 241313ae024SBarry Smith b->row = isirow; 24286bacbd2SBarry Smith b->col = iscolf; 243*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Scalar),&b->solve_work);CHKERRQ(ierr); 24486bacbd2SBarry Smith b->maxnz = b->nz = new_i[n]; 24586bacbd2SBarry Smith b->indexshift = a->indexshift; 24686bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 24786bacbd2SBarry Smith (*fact)->info.factor_mallocs = 0; 24886bacbd2SBarry Smith 24986bacbd2SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); 25086bacbd2SBarry Smith 25186bacbd2SBarry Smith /* check out for identical nodes. If found, use inode functions */ 25286bacbd2SBarry Smith ierr = Mat_AIJ_CheckInode(*fact);CHKERRQ(ierr); 25386bacbd2SBarry Smith 254431e710aSBarry Smith af = ((double)b->nz)/((double)a->nz) + .001; 255*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Fill ratio:given %g needed %g\n",info->fill,af); 256*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af); 257*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:PCILUSetFill(pc,%g);\n",af); 258*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUDTFactor_SeqAIJ:for best performance.\n"); 259431e710aSBarry Smith 26086bacbd2SBarry Smith PetscFunctionReturn(0); 26186bacbd2SBarry Smith } 26286bacbd2SBarry Smith 263289bc588SBarry Smith /* 264289bc588SBarry Smith Factorization code for AIJ format. 265289bc588SBarry Smith */ 2665615d1e5SSatish Balay #undef __FUNC__ 2679e046878SBarry Smith #define __FUNC__ /*<a name="MatLUFactorSymbolic_SeqAIJ""></a>*/"MatLUFactorSymbolic_SeqAIJ" 2689e046878SBarry Smith int MatLUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatLUInfo *info,Mat *B) 269289bc588SBarry Smith { 270416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 271289bc588SBarry Smith IS isicol; 272273d9f13SBarry Smith int *r,*ic,ierr,i,n = A->m,*ai = a->i,*aj = a->j; 273416022c9SBarry Smith int *ainew,*ajnew,jmax,*fill,*ajtmp,nz,shift = a->indexshift; 27491bf9adeSBarry Smith int *idnew,idx,row,m,fm,nnz,nzi,realloc = 0,nzbd,*im; 2759e046878SBarry Smith PetscReal f; 276289bc588SBarry Smith 2773a40ed3dSBarry Smith PetscFunctionBegin; 278d3cbd50cSLois Curfman McInnes PetscValidHeaderSpecific(isrow,IS_COOKIE); 279d3cbd50cSLois Curfman McInnes PetscValidHeaderSpecific(iscol,IS_COOKIE); 28029bbc08cSBarry Smith if (A->M != A->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 2813b2fbd54SBarry Smith 2829863cb4fSBarry Smith if (!isrow) { 2839863cb4fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->M,0,1,&isrow);CHKERRQ(ierr); 2849863cb4fSBarry Smith } 2859863cb4fSBarry Smith if (!iscol) { 2869863cb4fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->M,0,1,&iscol);CHKERRQ(ierr); 2879863cb4fSBarry Smith } 2889863cb4fSBarry Smith 2894c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 290ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 291ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 292289bc588SBarry Smith 293289bc588SBarry Smith /* get new row pointers */ 294*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 295dbb450caSBarry Smith ainew[0] = -shift; 296289bc588SBarry Smith /* don't know how many column pointers are needed so estimate */ 2979e046878SBarry Smith if (info) f = info->fill; else f = 1.0; 298dbb450caSBarry Smith jmax = (int)(f*ai[n]+(!shift)); 299*b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 300289bc588SBarry Smith /* fill is a linked list of nonzeros in active row */ 301*b0a32e0cSBarry Smith ierr = PetscMalloc((2*n+1)*sizeof(int),&fill);CHKERRQ(ierr); 3022fb3ed76SBarry Smith im = fill + n + 1; 303289bc588SBarry Smith /* idnew is location of diagonal in factor */ 304*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&idnew);CHKERRQ(ierr); 305dbb450caSBarry Smith idnew[0] = -shift; 306289bc588SBarry Smith 307289bc588SBarry Smith for (i=0; i<n; i++) { 308289bc588SBarry Smith /* first copy previous fill into linked list */ 309289bc588SBarry Smith nnz = nz = ai[r[i]+1] - ai[r[i]]; 31029bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 311dbb450caSBarry Smith ajtmp = aj + ai[r[i]] + shift; 312289bc588SBarry Smith fill[n] = n; 313289bc588SBarry Smith while (nz--) { 314289bc588SBarry Smith fm = n; 315dbb450caSBarry Smith idx = ic[*ajtmp++ + shift]; 316289bc588SBarry Smith do { 317289bc588SBarry Smith m = fm; 318289bc588SBarry Smith fm = fill[m]; 319289bc588SBarry Smith } while (fm < idx); 320289bc588SBarry Smith fill[m] = idx; 321289bc588SBarry Smith fill[idx] = fm; 322289bc588SBarry Smith } 323289bc588SBarry Smith row = fill[n]; 324289bc588SBarry Smith while (row < i) { 325dbb450caSBarry Smith ajtmp = ajnew + idnew[row] + (!shift); 3262fb3ed76SBarry Smith nzbd = 1 + idnew[row] - ainew[row]; 3272fb3ed76SBarry Smith nz = im[row] - nzbd; 328289bc588SBarry Smith fm = row; 3292fb3ed76SBarry Smith while (nz-- > 0) { 330dbb450caSBarry Smith idx = *ajtmp++ + shift; 3312fb3ed76SBarry Smith nzbd++; 3322fb3ed76SBarry Smith if (idx == i) im[row] = nzbd; 333289bc588SBarry Smith do { 334289bc588SBarry Smith m = fm; 335289bc588SBarry Smith fm = fill[m]; 336289bc588SBarry Smith } while (fm < idx); 337289bc588SBarry Smith if (fm != idx) { 338289bc588SBarry Smith fill[m] = idx; 339289bc588SBarry Smith fill[idx] = fm; 340289bc588SBarry Smith fm = idx; 341289bc588SBarry Smith nnz++; 342289bc588SBarry Smith } 343289bc588SBarry Smith } 344289bc588SBarry Smith row = fill[row]; 345289bc588SBarry Smith } 346289bc588SBarry Smith /* copy new filled row into permanent storage */ 347289bc588SBarry Smith ainew[i+1] = ainew[i] + nnz; 348496e697eSBarry Smith if (ainew[i+1] > jmax) { 349ecf371e4SBarry Smith 350ecf371e4SBarry Smith /* estimate how much additional space we will need */ 351ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 352ecf371e4SBarry Smith /* just double the memory each time */ 353ecf371e4SBarry Smith int maxadd = jmax; 354ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+(!shift))*(n-i+5))/n); */ 355bbb6d6a8SBarry Smith if (maxadd < nnz) maxadd = (n-i)*(nnz+1); 3562fb3ed76SBarry Smith jmax += maxadd; 357ecf371e4SBarry Smith 358ecf371e4SBarry Smith /* allocate a longer ajnew */ 359*b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&ajtmp);CHKERRQ(ierr); 360549d3d68SSatish Balay ierr = PetscMemcpy(ajtmp,ajnew,(ainew[i]+shift)*sizeof(int));CHKERRQ(ierr); 361606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 362289bc588SBarry Smith ajnew = ajtmp; 3632fb3ed76SBarry Smith realloc++; /* count how many times we realloc */ 364289bc588SBarry Smith } 365dbb450caSBarry Smith ajtmp = ajnew + ainew[i] + shift; 366289bc588SBarry Smith fm = fill[n]; 367289bc588SBarry Smith nzi = 0; 3682fb3ed76SBarry Smith im[i] = nnz; 369289bc588SBarry Smith while (nnz--) { 370289bc588SBarry Smith if (fm < i) nzi++; 371dbb450caSBarry Smith *ajtmp++ = fm - shift; 372289bc588SBarry Smith fm = fill[fm]; 373289bc588SBarry Smith } 374289bc588SBarry Smith idnew[i] = ainew[i] + nzi; 375289bc588SBarry Smith } 376464e8d44SSatish Balay if (ai[n] != 0) { 377329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 378*b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 379*b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:Run with -pc_lu_fill %g or use \n",af); 380*b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:PCLUSetFill(pc,%g);\n",af); 381*b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ:for best performance.\n"); 38205bf559cSSatish Balay } else { 383*b0a32e0cSBarry Smith PetscLogInfo(A,"MatLUFactorSymbolic_SeqAIJ: Empty matrix\n"); 38405bf559cSSatish Balay } 3852fb3ed76SBarry Smith 386898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 387898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3881987afe7SBarry Smith 389606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 390289bc588SBarry Smith 391289bc588SBarry Smith /* put together the new matrix */ 392b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,B);CHKERRQ(ierr); 393*b0a32e0cSBarry Smith PetscLogObjectParent(*B,isicol); 394416022c9SBarry Smith b = (Mat_SeqAIJ*)(*B)->data; 395606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 3967c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 397e8d4e0b9SBarry Smith /* the next line frees the default space generated by the Create() */ 398606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 399606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 400*b0a32e0cSBarry Smith ierr = PetscMalloc((ainew[n]+shift+1)*sizeof(Scalar),&b->a);CHKERRQ(ierr); 401416022c9SBarry Smith b->j = ajnew; 402416022c9SBarry Smith b->i = ainew; 403416022c9SBarry Smith b->diag = idnew; 404416022c9SBarry Smith b->ilen = 0; 405416022c9SBarry Smith b->imax = 0; 406416022c9SBarry Smith b->row = isrow; 407416022c9SBarry Smith b->col = iscol; 408085256b3SBarry Smith if (info) { 4091d1367b7SBarry Smith b->lu_damp = (PetscTruth) ((int)info->damp); 410085256b3SBarry Smith b->lu_damping = info->damping; 411085256b3SBarry Smith } else { 412085256b3SBarry Smith b->lu_damp = PETSC_FALSE; 413085256b3SBarry Smith b->lu_damping = 0.0; 414085256b3SBarry Smith } 415c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 416c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 41782bf6240SBarry Smith b->icol = isicol; 418*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Scalar),&b->solve_work);CHKERRQ(ierr); 419416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 4207fd98bd6SLois Curfman McInnes Allocate idnew, solve_work, new a, new j */ 421*b0a32e0cSBarry Smith PetscLogObjectMemory(*B,(ainew[n]+shift-n)*(sizeof(int)+sizeof(Scalar))); 422416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 4237fd98bd6SLois Curfman McInnes 424329f5518SBarry Smith (*B)->factor = FACTOR_LU; 425ae068f58SLois Curfman McInnes (*B)->info.factor_mallocs = realloc; 426ae068f58SLois Curfman McInnes (*B)->info.fill_ratio_given = f; 4277dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 428703deb49SSatish Balay 429e93ccc4dSBarry Smith if (ai[n] != 0) { 430329f5518SBarry Smith (*B)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 431eea2913aSSatish Balay } else { 432eea2913aSSatish Balay (*B)->info.fill_ratio_needed = 0.0; 433eea2913aSSatish Balay } 4343a40ed3dSBarry Smith PetscFunctionReturn(0); 435289bc588SBarry Smith } 4360a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 437ca44d042SBarry Smith EXTERN int Mat_AIJ_CheckInode(Mat); 43841c01911SSatish Balay 4395615d1e5SSatish Balay #undef __FUNC__ 440*b0a32e0cSBarry Smith #define __FUNC__ "MatLUFactorNumeric_SeqAIJ" 441416022c9SBarry Smith int MatLUFactorNumeric_SeqAIJ(Mat A,Mat *B) 442289bc588SBarry Smith { 443416022c9SBarry Smith Mat C = *B; 444416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ *)C->data; 44582bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 446273d9f13SBarry Smith int *r,*ic,ierr,i,j,n = A->m,*ai = b->i,*aj = b->j; 4471987afe7SBarry Smith int *ajtmpold,*ajtmp,nz,row,*ics,shift = a->indexshift; 448085256b3SBarry Smith int *diag_offset = b->diag,diag; 449085256b3SBarry Smith int *pj,ndamp = 0; 450085256b3SBarry Smith Scalar *rtmp,*v,*pc,multiplier,*pv,*rtmps; 451085256b3SBarry Smith PetscReal damping = b->lu_damping; 452085256b3SBarry Smith PetscTruth damp; 453289bc588SBarry Smith 4543a40ed3dSBarry Smith PetscFunctionBegin; 45582bf6240SBarry Smith 45678b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 45778b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 458*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Scalar),&rtmp);CHKERRQ(ierr); 459549d3d68SSatish Balay ierr = PetscMemzero(rtmp,(n+1)*sizeof(Scalar));CHKERRQ(ierr); 4600cf60383SBarry Smith rtmps = rtmp + shift; ics = ic + shift; 461289bc588SBarry Smith 462085256b3SBarry Smith do { 463085256b3SBarry Smith damp = PETSC_FALSE; 464289bc588SBarry Smith for (i=0; i<n; i++) { 465289bc588SBarry Smith nz = ai[i+1] - ai[i]; 466dbb450caSBarry Smith ajtmp = aj + ai[i] + shift; 46748e94559SBarry Smith for (j=0; j<nz; j++) rtmps[ajtmp[j]] = 0.0; 468289bc588SBarry Smith 469289bc588SBarry Smith /* load in initial (unfactored row) */ 470416022c9SBarry Smith nz = a->i[r[i]+1] - a->i[r[i]]; 471416022c9SBarry Smith ajtmpold = a->j + a->i[r[i]] + shift; 472416022c9SBarry Smith v = a->a + a->i[r[i]] + shift; 473085256b3SBarry Smith for (j=0; j<nz; j++) { 474085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] = v[j]; 475085256b3SBarry Smith if (ajtmpold[j] == r[i]) { /* damp the diagonal of the matrix */ 476085256b3SBarry Smith rtmp[ics[ajtmpold[j]]] += damping; 477085256b3SBarry Smith } 478085256b3SBarry Smith } 479289bc588SBarry Smith 480dbb450caSBarry Smith row = *ajtmp++ + shift; 481f2582111SSatish Balay while (row < i) { 4828c37ef55SBarry Smith pc = rtmp + row; 4838c37ef55SBarry Smith if (*pc != 0.0) { 4841987afe7SBarry Smith pv = b->a + diag_offset[row] + shift; 4851987afe7SBarry Smith pj = b->j + diag_offset[row] + (!shift); 48635aab85fSBarry Smith multiplier = *pc / *pv++; 4878c37ef55SBarry Smith *pc = multiplier; 488f2582111SSatish Balay nz = ai[row+1] - diag_offset[row] - 1; 489f2582111SSatish Balay for (j=0; j<nz; j++) rtmps[pj[j]] -= multiplier * pv[j]; 490*b0a32e0cSBarry Smith PetscLogFlops(2*nz); 491289bc588SBarry Smith } 492dbb450caSBarry Smith row = *ajtmp++ + shift; 493289bc588SBarry Smith } 494416022c9SBarry Smith /* finished row so stick it into b->a */ 495416022c9SBarry Smith pv = b->a + ai[i] + shift; 496416022c9SBarry Smith pj = b->j + ai[i] + shift; 4978c37ef55SBarry Smith nz = ai[i+1] - ai[i]; 498416022c9SBarry Smith for (j=0; j<nz; j++) {pv[j] = rtmps[pj[j]];} 49935aab85fSBarry Smith diag = diag_offset[i] - ai[i]; 500085256b3SBarry Smith if (PetscAbsScalar(pv[diag]) < 1.e-12) { 501085256b3SBarry Smith if (b->lu_damp) { 502085256b3SBarry Smith damp = PETSC_TRUE; 503085256b3SBarry Smith if (damping) damping *= 2.0; 504085256b3SBarry Smith else damping = 1.e-12; 505085256b3SBarry Smith ndamp++; 506085256b3SBarry Smith break; 507085256b3SBarry Smith } else { 50829bbc08cSBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %d",i); 509d708749eSBarry Smith } 51035aab85fSBarry Smith } 51135aab85fSBarry Smith } 512085256b3SBarry Smith } while (damp); 5130f11f9aeSSatish Balay 51435aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 51535aab85fSBarry Smith for (i=0; i<n; i++) { 51635aab85fSBarry Smith b->a[diag_offset[i]+shift] = 1.0/b->a[diag_offset[i]+shift]; 51735aab85fSBarry Smith } 51835aab85fSBarry Smith 519606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 52078b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 52178b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 522416022c9SBarry Smith C->factor = FACTOR_LU; 52341c01911SSatish Balay ierr = Mat_AIJ_CheckInode(C);CHKERRQ(ierr); 5247dda7e10SBarry Smith (*B)->ops->lufactornumeric = A->ops->lufactornumeric; /* Use Inode variant ONLY if A has inodes */ 525c456f294SBarry Smith C->assembled = PETSC_TRUE; 526*b0a32e0cSBarry Smith PetscLogFlops(C->n); 527a2e34c3dSBarry Smith if (ndamp || b->lu_damping) { 528*b0a32e0cSBarry Smith PetscLogInfo(0,"MatLUFactorNumerical_SeqAIJ: number of damping tries %d damping value %g\n",ndamp,damping); 529085256b3SBarry Smith } 5303a40ed3dSBarry Smith PetscFunctionReturn(0); 531289bc588SBarry Smith } 5320a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5335615d1e5SSatish Balay #undef __FUNC__ 534*b0a32e0cSBarry Smith #define __FUNC__ "MatLUFactor_SeqAIJ" 5359e046878SBarry Smith int MatLUFactor_SeqAIJ(Mat A,IS row,IS col,MatLUInfo *info) 536da3a660dSBarry Smith { 537273d9f13SBarry Smith int ierr; 538416022c9SBarry Smith Mat C; 539416022c9SBarry Smith 5403a40ed3dSBarry Smith PetscFunctionBegin; 5419e046878SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 542f2582111SSatish Balay ierr = MatLUFactorNumeric(A,&C);CHKERRQ(ierr); 543273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 5443a40ed3dSBarry Smith PetscFunctionReturn(0); 545da3a660dSBarry Smith } 5460a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 5475615d1e5SSatish Balay #undef __FUNC__ 548*b0a32e0cSBarry Smith #define __FUNC__ "MatSolve_SeqAIJ" 549416022c9SBarry Smith int MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 5508c37ef55SBarry Smith { 551416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 552416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 553273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 5543b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 555416022c9SBarry Smith Scalar *x,*b,*tmp,*tmps,*aa = a->a,sum,*v; 5568c37ef55SBarry Smith 5573a40ed3dSBarry Smith PetscFunctionBegin; 5583a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 55991bf9adeSBarry Smith 560e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 561e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 562416022c9SBarry Smith tmp = a->solve_work; 5638c37ef55SBarry Smith 5643b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 5653b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 5668c37ef55SBarry Smith 5678c37ef55SBarry Smith /* forward solve the lower triangular */ 5688c37ef55SBarry Smith tmp[0] = b[*r++]; 5690cf60383SBarry Smith tmps = tmp + shift; 5708c37ef55SBarry Smith for (i=1; i<n; i++) { 571dbb450caSBarry Smith v = aa + ai[i] + shift; 572dbb450caSBarry Smith vi = aj + ai[i] + shift; 573416022c9SBarry Smith nz = a->diag[i] - ai[i]; 5748c37ef55SBarry Smith sum = b[*r++]; 5750cf60383SBarry Smith while (nz--) sum -= *v++ * tmps[*vi++]; 5768c37ef55SBarry Smith tmp[i] = sum; 5778c37ef55SBarry Smith } 5788c37ef55SBarry Smith 5798c37ef55SBarry Smith /* backward solve the upper triangular */ 5808c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 581416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 582416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 583416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 5848c37ef55SBarry Smith sum = tmp[i]; 5850cf60383SBarry Smith while (nz--) sum -= *v++ * tmps[*vi++]; 586416022c9SBarry Smith x[*c--] = tmp[i] = sum*aa[a->diag[i]+shift]; 5878c37ef55SBarry Smith } 5888c37ef55SBarry Smith 5893b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 5903b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 591cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 592cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 593*b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 5943a40ed3dSBarry Smith PetscFunctionReturn(0); 5958c37ef55SBarry Smith } 596026e39d0SSatish Balay 597930ae53cSBarry Smith /* ----------------------------------------------------------- */ 598930ae53cSBarry Smith #undef __FUNC__ 599*b0a32e0cSBarry Smith #define __FUNC__ "MatSolve_SeqAIJ_NaturalOrdering" 600930ae53cSBarry Smith int MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 601930ae53cSBarry Smith { 602930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 603273d9f13SBarry Smith int n = A->m,*ai = a->i,*aj = a->j,*adiag = a->diag,ierr; 604d85afc42SSatish Balay Scalar *x,*b,*aa = a->a,sum; 605aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 606d85afc42SSatish Balay int adiag_i,i,*vi,nz,ai_i; 607d85afc42SSatish Balay Scalar *v; 608d85afc42SSatish Balay #endif 609930ae53cSBarry Smith 6103a40ed3dSBarry Smith PetscFunctionBegin; 6113a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 612930ae53cSBarry Smith if (a->indexshift) { 6133a40ed3dSBarry Smith ierr = MatSolve_SeqAIJ(A,bb,xx);CHKERRQ(ierr); 6143a40ed3dSBarry Smith PetscFunctionReturn(0); 615930ae53cSBarry Smith } 616930ae53cSBarry Smith 617e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 618e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 619930ae53cSBarry Smith 620aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 6211c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 6221c4feecaSSatish Balay #else 623930ae53cSBarry Smith /* forward solve the lower triangular */ 624930ae53cSBarry Smith x[0] = b[0]; 625930ae53cSBarry Smith for (i=1; i<n; i++) { 626930ae53cSBarry Smith ai_i = ai[i]; 627930ae53cSBarry Smith v = aa + ai_i; 628930ae53cSBarry Smith vi = aj + ai_i; 629930ae53cSBarry Smith nz = adiag[i] - ai_i; 630930ae53cSBarry Smith sum = b[i]; 631930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 632930ae53cSBarry Smith x[i] = sum; 633930ae53cSBarry Smith } 634930ae53cSBarry Smith 635930ae53cSBarry Smith /* backward solve the upper triangular */ 636930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 637930ae53cSBarry Smith adiag_i = adiag[i]; 638d708749eSBarry Smith v = aa + adiag_i + 1; 639d708749eSBarry Smith vi = aj + adiag_i + 1; 640930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 641930ae53cSBarry Smith sum = x[i]; 642930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 643930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 644930ae53cSBarry Smith } 6451c4feecaSSatish Balay #endif 646*b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->n); 647cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 648cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6493a40ed3dSBarry Smith PetscFunctionReturn(0); 650930ae53cSBarry Smith } 651930ae53cSBarry Smith 6525615d1e5SSatish Balay #undef __FUNC__ 653*b0a32e0cSBarry Smith #define __FUNC__ "MatSolveAdd_SeqAIJ" 654416022c9SBarry Smith int MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 655da3a660dSBarry Smith { 656416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 657416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 658273d9f13SBarry Smith int *r,*c,ierr,i, n = A->m,*vi,*ai = a->i,*aj = a->j; 6593b2fbd54SBarry Smith int nz,shift = a->indexshift,*rout,*cout; 660416022c9SBarry Smith Scalar *x,*b,*tmp,*aa = a->a,sum,*v; 661da3a660dSBarry Smith 6623a40ed3dSBarry Smith PetscFunctionBegin; 66378b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 664da3a660dSBarry Smith 665e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 666e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 667416022c9SBarry Smith tmp = a->solve_work; 668da3a660dSBarry Smith 6693b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6703b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 671da3a660dSBarry Smith 672da3a660dSBarry Smith /* forward solve the lower triangular */ 673da3a660dSBarry Smith tmp[0] = b[*r++]; 674da3a660dSBarry Smith for (i=1; i<n; i++) { 675dbb450caSBarry Smith v = aa + ai[i] + shift; 676dbb450caSBarry Smith vi = aj + ai[i] + shift; 677416022c9SBarry Smith nz = a->diag[i] - ai[i]; 678da3a660dSBarry Smith sum = b[*r++]; 679dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 680da3a660dSBarry Smith tmp[i] = sum; 681da3a660dSBarry Smith } 682da3a660dSBarry Smith 683da3a660dSBarry Smith /* backward solve the upper triangular */ 684da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 685416022c9SBarry Smith v = aa + a->diag[i] + (!shift); 686416022c9SBarry Smith vi = aj + a->diag[i] + (!shift); 687416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 688da3a660dSBarry Smith sum = tmp[i]; 689dbb450caSBarry Smith while (nz--) sum -= *v++ * tmp[*vi++ + shift]; 690416022c9SBarry Smith tmp[i] = sum*aa[a->diag[i]+shift]; 691da3a660dSBarry Smith x[*c--] += tmp[i]; 692da3a660dSBarry Smith } 693da3a660dSBarry Smith 6943b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6953b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 696cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 697cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 698*b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 699898183ecSLois Curfman McInnes 7003a40ed3dSBarry Smith PetscFunctionReturn(0); 701da3a660dSBarry Smith } 702da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 7035615d1e5SSatish Balay #undef __FUNC__ 704*b0a32e0cSBarry Smith #define __FUNC__ "MatSolveTranspose_SeqAIJ" 7057c922b88SBarry Smith int MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 706da3a660dSBarry Smith { 707416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7082235e667SBarry Smith IS iscol = a->col,isrow = a->row; 709273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 710f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 711f1af5d2fSBarry Smith Scalar *x,*b,*tmp,*aa = a->a,*v,s1; 712da3a660dSBarry Smith 7133a40ed3dSBarry Smith PetscFunctionBegin; 714e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 715e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 716416022c9SBarry Smith tmp = a->solve_work; 717da3a660dSBarry Smith 7182235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7192235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 720da3a660dSBarry Smith 721da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 7222235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 723da3a660dSBarry Smith 724da3a660dSBarry Smith /* forward solve the U^T */ 725da3a660dSBarry Smith for (i=0; i<n; i++) { 726f1af5d2fSBarry Smith v = aa + diag[i] + shift; 727f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 728f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 729c38d4ed2SBarry Smith s1 = tmp[i]; 730c38d4ed2SBarry Smith s1 *= *(v++); /* multiply by inverse of diagonal entry */ 731da3a660dSBarry Smith while (nz--) { 732f1af5d2fSBarry Smith tmp[*vi++ + shift] -= (*v++)*s1; 733da3a660dSBarry Smith } 734c38d4ed2SBarry Smith tmp[i] = s1; 735da3a660dSBarry Smith } 736da3a660dSBarry Smith 737da3a660dSBarry Smith /* backward solve the L^T */ 738da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 739f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 740f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 741f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 742f1af5d2fSBarry Smith s1 = tmp[i]; 743da3a660dSBarry Smith while (nz--) { 744f1af5d2fSBarry Smith tmp[*vi-- + shift] -= (*v--)*s1; 745da3a660dSBarry Smith } 746da3a660dSBarry Smith } 747da3a660dSBarry Smith 748da3a660dSBarry Smith /* copy tmp into x according to permutation */ 749da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 750da3a660dSBarry Smith 7512235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7522235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 753cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 754cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 755da3a660dSBarry Smith 756*b0a32e0cSBarry Smith PetscLogFlops(2*a->nz-A->n); 7573a40ed3dSBarry Smith PetscFunctionReturn(0); 758da3a660dSBarry Smith } 759da3a660dSBarry Smith 7605615d1e5SSatish Balay #undef __FUNC__ 761*b0a32e0cSBarry Smith #define __FUNC__ "MatSolveTransposeAdd_SeqAIJ" 7627c922b88SBarry Smith int MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 763da3a660dSBarry Smith { 764416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 7652235e667SBarry Smith IS iscol = a->col,isrow = a->row; 766273d9f13SBarry Smith int *r,*c,ierr,i,n = A->m,*vi,*ai = a->i,*aj = a->j; 767f1af5d2fSBarry Smith int nz,shift = a->indexshift,*rout,*cout,*diag = a->diag; 768416022c9SBarry Smith Scalar *x,*b,*tmp,*aa = a->a,*v; 7696abc6512SBarry Smith 7703a40ed3dSBarry Smith PetscFunctionBegin; 7716abc6512SBarry Smith if (zz != xx) VecCopy(zz,xx); 7726abc6512SBarry Smith 773e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 774e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 775416022c9SBarry Smith tmp = a->solve_work; 7766abc6512SBarry Smith 7772235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7782235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 7796abc6512SBarry Smith 7806abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 7812235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 7826abc6512SBarry Smith 7836abc6512SBarry Smith /* forward solve the U^T */ 7846abc6512SBarry Smith for (i=0; i<n; i++) { 785f1af5d2fSBarry Smith v = aa + diag[i] + shift; 786f1af5d2fSBarry Smith vi = aj + diag[i] + (!shift); 787f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 7886abc6512SBarry Smith tmp[i] *= *v++; 7896abc6512SBarry Smith while (nz--) { 790dbb450caSBarry Smith tmp[*vi++ + shift] -= (*v++)*tmp[i]; 7916abc6512SBarry Smith } 7926abc6512SBarry Smith } 7936abc6512SBarry Smith 7946abc6512SBarry Smith /* backward solve the L^T */ 7956abc6512SBarry Smith for (i=n-1; i>=0; i--){ 796f1af5d2fSBarry Smith v = aa + diag[i] - 1 + shift; 797f1af5d2fSBarry Smith vi = aj + diag[i] - 1 + shift; 798f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 7996abc6512SBarry Smith while (nz--) { 800dbb450caSBarry Smith tmp[*vi-- + shift] -= (*v--)*tmp[i]; 8016abc6512SBarry Smith } 8026abc6512SBarry Smith } 8036abc6512SBarry Smith 8046abc6512SBarry Smith /* copy tmp into x according to permutation */ 8056abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 8066abc6512SBarry Smith 8072235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 8082235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 809cb5b572fSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 810cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8116abc6512SBarry Smith 812*b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 8133a40ed3dSBarry Smith PetscFunctionReturn(0); 814da3a660dSBarry Smith } 8159e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 816ca44d042SBarry Smith EXTERN int MatMissingDiagonal_SeqAIJ(Mat); 81708480c60SBarry Smith 8185615d1e5SSatish Balay #undef __FUNC__ 819*b0a32e0cSBarry Smith #define __FUNC__ "MatILUFactorSymbolic_SeqAIJ" 8206cf997b0SBarry Smith int MatILUFactorSymbolic_SeqAIJ(Mat A,IS isrow,IS iscol,MatILUInfo *info,Mat *fact) 8219e25ed09SBarry Smith { 822416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 8239e25ed09SBarry Smith IS isicol; 824273d9f13SBarry Smith int *r,*ic,ierr,prow,n = A->m,*ai = a->i,*aj = a->j; 82556beaf04SBarry Smith int *ainew,*ajnew,jmax,*fill,*xi,nz,*im,*ajfill,*flev; 826335d9088SBarry Smith int *dloc,idx,row,m,fm,nzf,nzi,len, realloc = 0,dcount = 0; 8276cf997b0SBarry Smith int incrlev,nnz,i,shift = a->indexshift,levels,diagonal_fill; 82877c4ece6SBarry Smith PetscTruth col_identity,row_identity; 829329f5518SBarry Smith PetscReal f; 8309e25ed09SBarry Smith 8313a40ed3dSBarry Smith PetscFunctionBegin; 8326cf997b0SBarry Smith if (info) { 8336cf997b0SBarry Smith f = info->fill; 8340cd17293SBarry Smith levels = (int)info->levels; 8350cd17293SBarry Smith diagonal_fill = (int)info->diagonal_fill; 8366cf997b0SBarry Smith } else { 8376cf997b0SBarry Smith f = 1.0; 8386cf997b0SBarry Smith levels = 0; 8396cf997b0SBarry Smith diagonal_fill = 0; 8406cf997b0SBarry Smith } 8414da77479SBarry Smith if (!isrow) { 8424da77479SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->M,0,1,&isrow);CHKERRQ(ierr); 8434da77479SBarry Smith } 8444da77479SBarry Smith if (!iscol) { 8454da77479SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->M,0,1,&iscol);CHKERRQ(ierr); 8464da77479SBarry Smith } 8474da77479SBarry Smith 8484c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 84982bf6240SBarry Smith 85008480c60SBarry Smith /* special case that simply copies fill pattern */ 851be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 852be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 85386bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 8542e8a6d31SBarry Smith ierr = MatDuplicate_SeqAIJ(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 85508480c60SBarry Smith (*fact)->factor = FACTOR_LU; 85608480c60SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 85708480c60SBarry Smith if (!b->diag) { 8587c922b88SBarry Smith ierr = MatMarkDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 85908480c60SBarry Smith } 8607c922b88SBarry Smith ierr = MatMissingDiagonal_SeqAIJ(*fact);CHKERRQ(ierr); 86108480c60SBarry Smith b->row = isrow; 86208480c60SBarry Smith b->col = iscol; 86382bf6240SBarry Smith b->icol = isicol; 864a2e34c3dSBarry Smith if (info) { 8651d1367b7SBarry Smith b->lu_damp = (PetscTruth)((int)info->damp); 866a2e34c3dSBarry Smith b->lu_damping = info->damping; 867a2e34c3dSBarry Smith } else { 868a2e34c3dSBarry Smith b->lu_damp = PETSC_FALSE; 869a2e34c3dSBarry Smith b->lu_damping = 0.0; 870a2e34c3dSBarry Smith } 871*b0a32e0cSBarry Smith ierr = PetscMalloc(((*fact)->m+1)*sizeof(Scalar),&b->solve_work);CHKERRQ(ierr); 872f830108cSBarry Smith (*fact)->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 873c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 874c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 8753a40ed3dSBarry Smith PetscFunctionReturn(0); 87608480c60SBarry Smith } 87708480c60SBarry Smith 878898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 879898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 8809e25ed09SBarry Smith 8819e25ed09SBarry Smith /* get new row pointers */ 882*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&ainew);CHKERRQ(ierr); 883dbb450caSBarry Smith ainew[0] = -shift; 8849e25ed09SBarry Smith /* don't know how many column pointers are needed so estimate */ 885dbb450caSBarry Smith jmax = (int)(f*(ai[n]+!shift)); 886*b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajnew);CHKERRQ(ierr); 8879e25ed09SBarry Smith /* ajfill is level of fill for each fill entry */ 888*b0a32e0cSBarry Smith ierr = PetscMalloc((jmax)*sizeof(int),&ajfill);CHKERRQ(ierr); 8899e25ed09SBarry Smith /* fill is a linked list of nonzeros in active row */ 890*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&fill);CHKERRQ(ierr); 89156beaf04SBarry Smith /* im is level for each filled value */ 892*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&im);CHKERRQ(ierr); 89356beaf04SBarry Smith /* dloc is location of diagonal in factor */ 894*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&dloc);CHKERRQ(ierr); 89556beaf04SBarry Smith dloc[0] = 0; 89656beaf04SBarry Smith for (prow=0; prow<n; prow++) { 8976cf997b0SBarry Smith 8986cf997b0SBarry Smith /* copy current row into linked list */ 89956beaf04SBarry Smith nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 90029bbc08cSBarry Smith if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix"); 901dbb450caSBarry Smith xi = aj + ai[r[prow]] + shift; 9029e25ed09SBarry Smith fill[n] = n; 903435faa5fSBarry Smith fill[prow] = -1; /* marker to indicate if diagonal exists */ 9049e25ed09SBarry Smith while (nz--) { 9059e25ed09SBarry Smith fm = n; 906dbb450caSBarry Smith idx = ic[*xi++ + shift]; 9079e25ed09SBarry Smith do { 9089e25ed09SBarry Smith m = fm; 9099e25ed09SBarry Smith fm = fill[m]; 9109e25ed09SBarry Smith } while (fm < idx); 9119e25ed09SBarry Smith fill[m] = idx; 9129e25ed09SBarry Smith fill[idx] = fm; 91356beaf04SBarry Smith im[idx] = 0; 9149e25ed09SBarry Smith } 9156cf997b0SBarry Smith 9166cf997b0SBarry Smith /* make sure diagonal entry is included */ 917435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 9186cf997b0SBarry Smith fm = n; 919435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 920435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 921435faa5fSBarry Smith fill[fm] = prow; 9226cf997b0SBarry Smith im[prow] = 0; 9236cf997b0SBarry Smith nzf++; 924335d9088SBarry Smith dcount++; 9256cf997b0SBarry Smith } 9266cf997b0SBarry Smith 92756beaf04SBarry Smith nzi = 0; 9289e25ed09SBarry Smith row = fill[n]; 92956beaf04SBarry Smith while (row < prow) { 93056beaf04SBarry Smith incrlev = im[row] + 1; 93156beaf04SBarry Smith nz = dloc[row]; 9326cf997b0SBarry Smith xi = ajnew + ainew[row] + shift + nz + 1; 933dbb450caSBarry Smith flev = ajfill + ainew[row] + shift + nz + 1; 934416022c9SBarry Smith nnz = ainew[row+1] - ainew[row] - nz - 1; 93556beaf04SBarry Smith fm = row; 93656beaf04SBarry Smith while (nnz-- > 0) { 937dbb450caSBarry Smith idx = *xi++ + shift; 93856beaf04SBarry Smith if (*flev + incrlev > levels) { 93956beaf04SBarry Smith flev++; 94056beaf04SBarry Smith continue; 94156beaf04SBarry Smith } 94256beaf04SBarry Smith do { 9439e25ed09SBarry Smith m = fm; 9449e25ed09SBarry Smith fm = fill[m]; 94556beaf04SBarry Smith } while (fm < idx); 9469e25ed09SBarry Smith if (fm != idx) { 94756beaf04SBarry Smith im[idx] = *flev + incrlev; 9489e25ed09SBarry Smith fill[m] = idx; 9499e25ed09SBarry Smith fill[idx] = fm; 9509e25ed09SBarry Smith fm = idx; 95156beaf04SBarry Smith nzf++; 952ecf371e4SBarry Smith } else { 95356beaf04SBarry Smith if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 9549e25ed09SBarry Smith } 95556beaf04SBarry Smith flev++; 9569e25ed09SBarry Smith } 9579e25ed09SBarry Smith row = fill[row]; 95856beaf04SBarry Smith nzi++; 9599e25ed09SBarry Smith } 9609e25ed09SBarry Smith /* copy new filled row into permanent storage */ 96156beaf04SBarry Smith ainew[prow+1] = ainew[prow] + nzf; 962d7e8b826SBarry Smith if (ainew[prow+1] > jmax-shift) { 963ecf371e4SBarry Smith 964ecf371e4SBarry Smith /* estimate how much additional space we will need */ 965ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 966ecf371e4SBarry Smith /* just double the memory each time */ 967ecf371e4SBarry Smith /* maxadd = (int)((f*(ai[n]+!shift)*(n-prow+5))/n); */ 968ecf371e4SBarry Smith int maxadd = jmax; 96956beaf04SBarry Smith if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 970bbb6d6a8SBarry Smith jmax += maxadd; 971ecf371e4SBarry Smith 972ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 973*b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 974549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajnew,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 975606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 97656beaf04SBarry Smith ajnew = xi; 977*b0a32e0cSBarry Smith ierr = PetscMalloc(jmax*sizeof(int),&xi);CHKERRQ(ierr); 978549d3d68SSatish Balay ierr = PetscMemcpy(xi,ajfill,(ainew[prow]+shift)*sizeof(int));CHKERRQ(ierr); 979606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 98056beaf04SBarry Smith ajfill = xi; 981ecf371e4SBarry Smith realloc++; /* count how many times we realloc */ 9829e25ed09SBarry Smith } 983dbb450caSBarry Smith xi = ajnew + ainew[prow] + shift; 984dbb450caSBarry Smith flev = ajfill + ainew[prow] + shift; 98556beaf04SBarry Smith dloc[prow] = nzi; 9869e25ed09SBarry Smith fm = fill[n]; 98756beaf04SBarry Smith while (nzf--) { 988dbb450caSBarry Smith *xi++ = fm - shift; 98956beaf04SBarry Smith *flev++ = im[fm]; 9909e25ed09SBarry Smith fm = fill[fm]; 9919e25ed09SBarry Smith } 9926cf997b0SBarry Smith /* make sure row has diagonal entry */ 9936cf997b0SBarry Smith if (ajnew[ainew[prow]+shift+dloc[prow]]+shift != prow) { 99429bbc08cSBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %d has missing diagonal in factored matrix\n\ 9956cf997b0SBarry Smith try running with -pc_ilu_nonzeros_along_diagonal or -pc_ilu_diagonal_fill",prow); 9966cf997b0SBarry Smith } 9979e25ed09SBarry Smith } 998606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 999898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1000898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1001606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 1002606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 10039e25ed09SBarry Smith 1004f64065bbSBarry Smith { 1005329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 1006*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Reallocs %d Fill ratio:given %g needed %g\n",realloc,f,af); 1007*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Run with -pc_ilu_fill %g or use \n",af); 1008*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:PCILUSetFill(pc,%g);\n",af); 1009*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:for best performance.\n"); 1010335d9088SBarry Smith if (diagonal_fill) { 1011*b0a32e0cSBarry Smith PetscLogInfo(A,"MatILUFactorSymbolic_SeqAIJ:Detected and replace %d missing diagonals",dcount); 1012335d9088SBarry Smith } 1013f64065bbSBarry Smith } 1014bbb6d6a8SBarry Smith 10159e25ed09SBarry Smith /* put together the new matrix */ 1016b4fd4287SBarry Smith ierr = MatCreateSeqAIJ(A->comm,n,n,0,PETSC_NULL,fact);CHKERRQ(ierr); 1017*b0a32e0cSBarry Smith PetscLogObjectParent(*fact,isicol); 1018416022c9SBarry Smith b = (Mat_SeqAIJ*)(*fact)->data; 1019606d414cSSatish Balay ierr = PetscFree(b->imax);CHKERRQ(ierr); 10207c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 1021dbb450caSBarry Smith len = (ainew[n] + shift)*sizeof(Scalar); 10229e25ed09SBarry Smith /* the next line frees the default space generated by the Create() */ 1023606d414cSSatish Balay ierr = PetscFree(b->a);CHKERRQ(ierr); 1024606d414cSSatish Balay ierr = PetscFree(b->ilen);CHKERRQ(ierr); 1025*b0a32e0cSBarry Smith ierr = PetscMalloc(len+1,&b->a);CHKERRQ(ierr); 1026416022c9SBarry Smith b->j = ajnew; 1027416022c9SBarry Smith b->i = ainew; 102856beaf04SBarry Smith for (i=0; i<n; i++) dloc[i] += ainew[i]; 1029416022c9SBarry Smith b->diag = dloc; 1030416022c9SBarry Smith b->ilen = 0; 1031416022c9SBarry Smith b->imax = 0; 1032416022c9SBarry Smith b->row = isrow; 1033416022c9SBarry Smith b->col = iscol; 1034c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1035c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 103682bf6240SBarry Smith b->icol = isicol; 1037*b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Scalar),&b->solve_work);CHKERRQ(ierr); 1038416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 103956beaf04SBarry Smith Allocate dloc, solve_work, new a, new j */ 1040*b0a32e0cSBarry Smith PetscLogObjectMemory(*fact,(ainew[n]+shift-n) * (sizeof(int)+sizeof(Scalar))); 1041416022c9SBarry Smith b->maxnz = b->nz = ainew[n] + shift; 1042a2e34c3dSBarry Smith if (info) { 10431d1367b7SBarry Smith b->lu_damp = (PetscTruth)((int)info->damp); 1044a2e34c3dSBarry Smith b->lu_damping = info->damping; 1045a2e34c3dSBarry Smith } else { 1046a2e34c3dSBarry Smith b->lu_damp = PETSC_FALSE; 1047a2e34c3dSBarry Smith b->lu_damping = 0.0; 1048a2e34c3dSBarry Smith } 10499e25ed09SBarry Smith (*fact)->factor = FACTOR_LU; 1050ae068f58SLois Curfman McInnes 1051ae068f58SLois Curfman McInnes (*fact)->info.factor_mallocs = realloc; 1052ae068f58SLois Curfman McInnes (*fact)->info.fill_ratio_given = f; 1053329f5518SBarry Smith (*fact)->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 1054329f5518SBarry Smith (*fact)->factor = FACTOR_LU; 1055ae068f58SLois Curfman McInnes 10563a40ed3dSBarry Smith PetscFunctionReturn(0); 10579e25ed09SBarry Smith } 1058d5d45c9bSBarry Smith 1059d5d45c9bSBarry Smith 1060d5d45c9bSBarry Smith 1061d5d45c9bSBarry Smith 1062