1be1d678aSKris Buschelman #define PETSCMAT_DLL 2289bc588SBarry Smith 37c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 47c4f633dSBarry Smith #include "../src/inline/dot.h" 5*e09e08ccSBarry Smith #define PETSC_USE_WHILE_KERNELS 67c4f633dSBarry Smith #include "../src/inline/spops.h" 71393bc97SHong Zhang #include "petscbt.h" 87c4f633dSBarry Smith #include "../src/mat/utils/freespace.h" 93b2fbd54SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 128fa24674SBarry Smith /* 138fa24674SBarry Smith Computes an ordering to get most of the large numerical values in the lower triangular part of the matrix 148fa24674SBarry Smith */ 15dfbe8321SBarry Smith PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 163b2fbd54SBarry Smith { 178fa24674SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)mat->data; 188fa24674SBarry Smith PetscErrorCode ierr; 198fa24674SBarry Smith PetscInt i,j,jj,k, kk,n = mat->rmap->n, current, newcurrent,*order; 208fa24674SBarry Smith const PetscInt *ai = a->i, *aj = a->j; 218fa24674SBarry Smith const PetscScalar *aa = a->a; 228fa24674SBarry Smith PetscTruth *done; 238fa24674SBarry Smith PetscReal best,past,future; 243a40ed3dSBarry Smith 258fa24674SBarry Smith PetscFunctionBegin; 268fa24674SBarry Smith /* pick initial row */ 278fa24674SBarry Smith best = -1; 288fa24674SBarry Smith for (i=0; i<n; i++) { 298fa24674SBarry Smith future = 0; 308fa24674SBarry Smith for (j=ai[i]; j<ai[i+1]; j++) { 318fa24674SBarry Smith if (aj[j] != i) future += PetscAbsScalar(aa[j]); else past = PetscAbsScalar(aa[j]); 328fa24674SBarry Smith } 338fa24674SBarry Smith if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 348fa24674SBarry Smith if (past/future > best) { 358fa24674SBarry Smith best = past/future; 368fa24674SBarry Smith current = i; 378fa24674SBarry Smith } 388fa24674SBarry Smith } 398fa24674SBarry Smith 408fa24674SBarry Smith ierr = PetscMalloc(n*sizeof(PetscTruth),&done);CHKERRQ(ierr); 418fa24674SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&order);CHKERRQ(ierr); 428fa24674SBarry Smith ierr = PetscMemzero(done,n*sizeof(PetscTruth));CHKERRQ(ierr); 438fa24674SBarry Smith order[0] = current; 448fa24674SBarry Smith for (i=0; i<n-1; i++) { 458fa24674SBarry Smith done[current] = PETSC_TRUE; 468fa24674SBarry Smith best = -1; 478fa24674SBarry Smith /* loop over all neighbors of current pivot */ 488fa24674SBarry Smith for (j=ai[current]; j<ai[current+1]; j++) { 498fa24674SBarry Smith jj = aj[j]; 508fa24674SBarry Smith if (done[jj]) continue; 518fa24674SBarry Smith /* loop over columns of potential next row computing weights for below and above diagonal */ 528fa24674SBarry Smith past = future = 0.0; 538fa24674SBarry Smith for (k=ai[jj]; k<ai[jj+1]; k++) { 548fa24674SBarry Smith kk = aj[k]; 558fa24674SBarry Smith if (done[kk]) past += PetscAbsScalar(aa[k]); 568fa24674SBarry Smith else if (kk != jj) future += PetscAbsScalar(aa[k]); 578fa24674SBarry Smith } 588fa24674SBarry Smith if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 598fa24674SBarry Smith if (past/future > best) { 608fa24674SBarry Smith best = past/future; 618fa24674SBarry Smith newcurrent = jj; 628fa24674SBarry Smith } 638fa24674SBarry Smith } 648fa24674SBarry Smith if (best == -1) { /* no neighbors to select from so select best of all that remain */ 658fa24674SBarry Smith best = -1; 668fa24674SBarry Smith for (k=0; k<n; k++) { 678fa24674SBarry Smith if (done[k]) continue; 688fa24674SBarry Smith future = 0; 698fa24674SBarry Smith past = 0; 708fa24674SBarry Smith for (j=ai[k]; j<ai[k+1]; j++) { 718fa24674SBarry Smith kk = aj[j]; 728fa24674SBarry Smith if (done[kk]) past += PetscAbsScalar(aa[j]); 738fa24674SBarry Smith else if (kk != k) future += PetscAbsScalar(aa[j]); 748fa24674SBarry Smith } 758fa24674SBarry Smith if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 768fa24674SBarry Smith if (past/future > best) { 778fa24674SBarry Smith best = past/future; 788fa24674SBarry Smith newcurrent = k; 798fa24674SBarry Smith } 808fa24674SBarry Smith } 818fa24674SBarry Smith } 828fa24674SBarry Smith if (current == newcurrent) SETERRQ(PETSC_ERR_PLIB,"newcurrent cannot be current"); 838fa24674SBarry Smith current = newcurrent; 848fa24674SBarry Smith order[i+1] = current; 858fa24674SBarry Smith } 868fa24674SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,n,order,irow);CHKERRQ(ierr); 878fa24674SBarry Smith *icol = *irow; 888fa24674SBarry Smith ierr = PetscObjectReference((PetscObject)*irow);CHKERRQ(ierr); 898fa24674SBarry Smith ierr = PetscFree(done);CHKERRQ(ierr); 908fa24674SBarry Smith ierr = PetscFree(order);CHKERRQ(ierr); 91d64ed03dSBarry Smith PetscFunctionReturn(0); 923b2fbd54SBarry Smith } 933b2fbd54SBarry Smith 94e631078cSBarry Smith EXTERN_C_BEGIN 954a2ae208SSatish Balay #undef __FUNCT__ 96db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc" 97db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 98db4efbfdSBarry Smith { 99db4efbfdSBarry Smith PetscFunctionBegin; 100db4efbfdSBarry Smith *flg = PETSC_TRUE; 101db4efbfdSBarry Smith PetscFunctionReturn(0); 102db4efbfdSBarry Smith } 103db4efbfdSBarry Smith EXTERN_C_END 104db4efbfdSBarry Smith 105db4efbfdSBarry Smith EXTERN_C_BEGIN 106db4efbfdSBarry Smith #undef __FUNCT__ 107b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_petsc" 108b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B) 109b24902e0SBarry Smith { 110d0f46423SBarry Smith PetscInt n = A->rmap->n; 111b24902e0SBarry Smith PetscErrorCode ierr; 112b24902e0SBarry Smith 113b24902e0SBarry Smith PetscFunctionBegin; 114b24902e0SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 115b24902e0SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 116599ef60dSHong Zhang if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT){ 117b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr); 118db4efbfdSBarry Smith (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJ; 11935bd34faSBarry Smith (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 1201da05e5aSHong Zhang (*B)->ops->iludtfactor = MatILUDTFactor_SeqAIJ; 121b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 122b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 123b24902e0SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1245c9eb25fSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ; 1255c9eb25fSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ; 126b24902e0SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 127db4efbfdSBarry Smith (*B)->factor = ftype; 128b24902e0SBarry Smith PetscFunctionReturn(0); 129b24902e0SBarry Smith } 130e631078cSBarry Smith EXTERN_C_END 131b24902e0SBarry Smith 132b24902e0SBarry Smith #undef __FUNCT__ 133b9617806SBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 1340481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 135289bc588SBarry Smith { 136416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 137289bc588SBarry Smith IS isicol; 1386849ba73SBarry Smith PetscErrorCode ierr; 1395d0c19d7SBarry Smith const PetscInt *r,*ic; 1405d0c19d7SBarry Smith PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j; 1411393bc97SHong Zhang PetscInt *bi,*bj,*ajtmp; 1421393bc97SHong Zhang PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 1439e046878SBarry Smith PetscReal f; 1444875ba38SHong Zhang PetscInt nlnk,*lnk,k,**bi_ptr; 145a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1461393bc97SHong Zhang PetscBT lnkbt; 147289bc588SBarry Smith 1483a40ed3dSBarry Smith PetscFunctionBegin; 149d0f46423SBarry Smith if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 1504c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 151ac121b3dSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 152ac121b3dSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 153289bc588SBarry Smith 154289bc588SBarry Smith /* get new row pointers */ 1551393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1561393bc97SHong Zhang bi[0] = 0; 1571393bc97SHong Zhang 1581393bc97SHong Zhang /* bdiag is location of diagonal in factor */ 1591393bc97SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1601393bc97SHong Zhang bdiag[0] = 0; 1611393bc97SHong Zhang 1621393bc97SHong Zhang /* linked list for storing column indices of the active row */ 1631393bc97SHong Zhang nlnk = n + 1; 1641393bc97SHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1651393bc97SHong Zhang 16635baf27eSSatish Balay ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 1671393bc97SHong Zhang 1681393bc97SHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 169d3d32019SBarry Smith f = info->fill; 170a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 1711393bc97SHong Zhang current_space = free_space; 172289bc588SBarry Smith 173289bc588SBarry Smith for (i=0; i<n; i++) { 1741393bc97SHong Zhang /* copy previous fill into linked list */ 175289bc588SBarry Smith nzi = 0; 1761393bc97SHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 1773b4a8b6dSBarry Smith if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1781393bc97SHong Zhang ajtmp = aj + ai[r[i]]; 179afcc9904SHong Zhang ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1801393bc97SHong Zhang nzi += nlnk; 1811393bc97SHong Zhang 1821393bc97SHong Zhang /* add pivot rows into linked list */ 1831393bc97SHong Zhang row = lnk[n]; 1841393bc97SHong Zhang while (row < i) { 185d1170560SHong Zhang nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 186d1170560SHong Zhang ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 187d1170560SHong Zhang ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 1881393bc97SHong Zhang nzi += nlnk; 1891393bc97SHong Zhang row = lnk[row]; 190289bc588SBarry Smith } 1911393bc97SHong Zhang bi[i+1] = bi[i] + nzi; 1921393bc97SHong Zhang im[i] = nzi; 1931393bc97SHong Zhang 1941393bc97SHong Zhang /* mark bdiag */ 1951393bc97SHong Zhang nzbd = 0; 1961393bc97SHong Zhang nnz = nzi; 1971393bc97SHong Zhang k = lnk[n]; 1981393bc97SHong Zhang while (nnz-- && k < i){ 1991393bc97SHong Zhang nzbd++; 2001393bc97SHong Zhang k = lnk[k]; 2011393bc97SHong Zhang } 2021393bc97SHong Zhang bdiag[i] = bi[i] + nzbd; 2031393bc97SHong Zhang 2041393bc97SHong Zhang /* if free space is not available, make more free space */ 2051393bc97SHong Zhang if (current_space->local_remaining<nzi) { 2061393bc97SHong Zhang nnz = (n - i)*nzi; /* estimated and max additional space needed */ 207a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 2081393bc97SHong Zhang reallocs++; 2091393bc97SHong Zhang } 2101393bc97SHong Zhang 2111393bc97SHong Zhang /* copy data into free space, then initialize lnk */ 2121393bc97SHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 2131393bc97SHong Zhang bi_ptr[i] = current_space->array; 2141393bc97SHong Zhang current_space->array += nzi; 2151393bc97SHong Zhang current_space->local_used += nzi; 2161393bc97SHong Zhang current_space->local_remaining -= nzi; 217289bc588SBarry Smith } 2186cf91177SBarry Smith #if defined(PETSC_USE_INFO) 219464e8d44SSatish Balay if (ai[n] != 0) { 2201393bc97SHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 221ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 222ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 223ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 224ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 22505bf559cSSatish Balay } else { 226ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 22705bf559cSSatish Balay } 22863ba0a88SBarry Smith #endif 2292fb3ed76SBarry Smith 230898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 231898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2321987afe7SBarry Smith 2331393bc97SHong Zhang /* destroy list of free space and other temporary array(s) */ 2341393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 235a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 2361393bc97SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 23735baf27eSSatish Balay ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 238289bc588SBarry Smith 239289bc588SBarry Smith /* put together the new matrix */ 240719d5645SBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 241719d5645SBarry Smith ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 242719d5645SBarry Smith b = (Mat_SeqAIJ*)(B)->data; 243e6b907acSBarry Smith b->free_a = PETSC_TRUE; 244e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 2457c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 2461393bc97SHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 2471393bc97SHong Zhang b->j = bj; 2481393bc97SHong Zhang b->i = bi; 2491393bc97SHong Zhang b->diag = bdiag; 250416022c9SBarry Smith b->ilen = 0; 251416022c9SBarry Smith b->imax = 0; 252416022c9SBarry Smith b->row = isrow; 253416022c9SBarry Smith b->col = iscol; 254c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 255c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 25682bf6240SBarry Smith b->icol = isicol; 25787828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2581393bc97SHong Zhang 2591393bc97SHong Zhang /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 260719d5645SBarry Smith ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 2611393bc97SHong Zhang b->maxnz = b->nz = bi[n] ; 2627fd98bd6SLois Curfman McInnes 263719d5645SBarry Smith (B)->factor = MAT_FACTOR_LU; 264719d5645SBarry Smith (B)->info.factor_mallocs = reallocs; 265719d5645SBarry Smith (B)->info.fill_ratio_given = f; 266703deb49SSatish Balay 2679f7d0b68SHong Zhang if (ai[n]) { 268719d5645SBarry Smith (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 269eea2913aSSatish Balay } else { 270719d5645SBarry Smith (B)->info.fill_ratio_needed = 0.0; 271eea2913aSSatish Balay } 272719d5645SBarry Smith (B)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 273719d5645SBarry Smith (B)->ops->solve = MatSolve_SeqAIJ; 274719d5645SBarry Smith (B)->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 275db4efbfdSBarry Smith /* switch to inodes if appropriate */ 276719d5645SBarry Smith ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr); 2773a40ed3dSBarry Smith PetscFunctionReturn(0); 278289bc588SBarry Smith } 2791393bc97SHong Zhang 2805b5bc046SBarry Smith /* 2815b5bc046SBarry Smith Trouble in factorization, should we dump the original matrix? 2825b5bc046SBarry Smith */ 2835b5bc046SBarry Smith #undef __FUNCT__ 2845b5bc046SBarry Smith #define __FUNCT__ "MatFactorDumpMatrix" 2855b5bc046SBarry Smith PetscErrorCode MatFactorDumpMatrix(Mat A) 2865b5bc046SBarry Smith { 2875b5bc046SBarry Smith PetscErrorCode ierr; 28890d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE; 2895b5bc046SBarry Smith 2905b5bc046SBarry Smith PetscFunctionBegin; 29190d69ab7SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_factor_dump_on_error",&flg,PETSC_NULL);CHKERRQ(ierr); 2925b5bc046SBarry Smith if (flg) { 2935b5bc046SBarry Smith PetscViewer viewer; 2945b5bc046SBarry Smith char filename[PETSC_MAX_PATH_LEN]; 2955b5bc046SBarry Smith 2965b5bc046SBarry Smith ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr); 2977adad957SLisandro Dalcin ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 2985b5bc046SBarry Smith ierr = MatView(A,viewer);CHKERRQ(ierr); 2995b5bc046SBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 3005b5bc046SBarry Smith } 3015b5bc046SBarry Smith PetscFunctionReturn(0); 3025b5bc046SBarry Smith } 3035b5bc046SBarry Smith 304db4efbfdSBarry Smith extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec); 305db4efbfdSBarry Smith 3060a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 3074a2ae208SSatish Balay #undef __FUNCT__ 3084a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 3090481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info) 310289bc588SBarry Smith { 311719d5645SBarry Smith Mat C=B; 312416022c9SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 31382bf6240SBarry Smith IS isrow = b->row,isicol = b->icol; 3146849ba73SBarry Smith PetscErrorCode ierr; 3155d0c19d7SBarry Smith const PetscInt *r,*ic,*ics; 3169f7d0b68SHong Zhang PetscInt i,j,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 3179f7d0b68SHong Zhang PetscInt *ajtmp,*bjtmp,nz,row,*diag_offset = b->diag,diag,*pj; 3189f7d0b68SHong Zhang MatScalar *rtmp,*pc,multiplier,*v,*pv,d,*aa=a->a; 31933f9893dSHong Zhang PetscReal rs=0.0; 320b3bf805bSHong Zhang LUShift_Ctx sctx; 32142898933SHong Zhang PetscInt newshift,*ddiag; 322289bc588SBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 32478b31e54SBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 32578b31e54SBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 326b78a477dSHong Zhang ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3279f7d0b68SHong Zhang ics = ic; 328289bc588SBarry Smith 329ada7143aSSatish Balay sctx.shift_top = 0; 330ada7143aSSatish Balay sctx.nshift_max = 0; 331ada7143aSSatish Balay sctx.shift_lo = 0; 332ada7143aSSatish Balay sctx.shift_hi = 0; 3330481f469SBarry Smith sctx.shift_fraction = 0; 334ada7143aSSatish Balay 335118f5789SBarry Smith /* if both shift schemes are chosen by user, only use info->shiftpd */ 3361a890ab1SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 33742898933SHong Zhang ddiag = a->diag; 3389f7d0b68SHong Zhang sctx.shift_top = info->zeropivot; 3396cc28720Svictorle for (i=0; i<n; i++) { 3409f95998fSHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 3419f7d0b68SHong Zhang d = (aa)[ddiag[i]]; 3421a890ab1SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 3439f7d0b68SHong Zhang v = aa+ai[i]; 3449f7d0b68SHong Zhang nz = ai[i+1] - ai[i]; 345e24cfe64SHong Zhang for (j=0; j<nz; j++) 3461a890ab1SHong Zhang rs += PetscAbsScalar(v[j]); 3471a890ab1SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 3486cc28720Svictorle } 349118f5789SBarry Smith sctx.shift_top *= 1.1; 350d2276718SHong Zhang sctx.nshift_max = 5; 351d2276718SHong Zhang sctx.shift_lo = 0.; 352d2276718SHong Zhang sctx.shift_hi = 1.; 353a0a356efSHong Zhang } 354a0a356efSHong Zhang 3559f7d0b68SHong Zhang sctx.shift_amount = 0.0; 356a0a356efSHong Zhang sctx.nshift = 0; 357085256b3SBarry Smith do { 358d2276718SHong Zhang sctx.lushift = PETSC_FALSE; 359289bc588SBarry Smith for (i=0; i<n; i++){ 360b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 361b3bf805bSHong Zhang bjtmp = bj + bi[i]; 3629f7d0b68SHong Zhang for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 363289bc588SBarry Smith 364289bc588SBarry Smith /* load in initial (unfactored row) */ 3659f7d0b68SHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 3669f7d0b68SHong Zhang ajtmp = aj + ai[r[i]]; 3679f7d0b68SHong Zhang v = aa + ai[r[i]]; 368085256b3SBarry Smith for (j=0; j<nz; j++) { 369b3bf805bSHong Zhang rtmp[ics[ajtmp[j]]] = v[j]; 370085256b3SBarry Smith } 371d2276718SHong Zhang rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 372289bc588SBarry Smith 373b3bf805bSHong Zhang row = *bjtmp++; 374f2582111SSatish Balay while (row < i) { 3758c37ef55SBarry Smith pc = rtmp + row; 3768c37ef55SBarry Smith if (*pc != 0.0) { 377010a20caSHong Zhang pv = b->a + diag_offset[row]; 378010a20caSHong Zhang pj = b->j + diag_offset[row] + 1; 37935aab85fSBarry Smith multiplier = *pc / *pv++; 3808c37ef55SBarry Smith *pc = multiplier; 381b3bf805bSHong Zhang nz = bi[row+1] - diag_offset[row] - 1; 3829f7d0b68SHong Zhang for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 383dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 384289bc588SBarry Smith } 385b3bf805bSHong Zhang row = *bjtmp++; 386289bc588SBarry Smith } 387416022c9SBarry Smith /* finished row so stick it into b->a */ 388b3bf805bSHong Zhang pv = b->a + bi[i] ; 389b3bf805bSHong Zhang pj = b->j + bi[i] ; 390b3bf805bSHong Zhang nz = bi[i+1] - bi[i]; 391b3bf805bSHong Zhang diag = diag_offset[i] - bi[i]; 3929f7d0b68SHong Zhang rs = -PetscAbsScalar(pv[diag]); 393d3d32019SBarry Smith for (j=0; j<nz; j++) { 3949f7d0b68SHong Zhang pv[j] = rtmp[pj[j]]; 3959f7d0b68SHong Zhang rs += PetscAbsScalar(pv[j]); 396d3d32019SBarry Smith } 397d2276718SHong Zhang 398b3bf805bSHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 399d2276718SHong Zhang sctx.rs = rs; 400d2276718SHong Zhang sctx.pv = pv[diag]; 401c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 4025b5bc046SBarry Smith if (newshift == 1) break; 40335aab85fSBarry Smith } 404d2276718SHong Zhang 4050481f469SBarry Smith if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 4066cc28720Svictorle /* 4076c037d1bSvictorle * if no shift in this attempt & shifting & started shifting & can refine, 4086cc28720Svictorle * then try lower shift 4096cc28720Svictorle */ 4100481f469SBarry Smith sctx.shift_hi = sctx.shift_fraction; 4110481f469SBarry Smith sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 4120481f469SBarry Smith sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 413d2276718SHong Zhang sctx.lushift = PETSC_TRUE; 414d2276718SHong Zhang sctx.nshift++; 4156cc28720Svictorle } 416d2276718SHong Zhang } while (sctx.lushift); 4170f11f9aeSSatish Balay 41835aab85fSBarry Smith /* invert diagonal entries for simplier triangular solves */ 41935aab85fSBarry Smith for (i=0; i<n; i++) { 420010a20caSHong Zhang b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 42135aab85fSBarry Smith } 422606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 42378b31e54SBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 42478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 425db4efbfdSBarry Smith if (b->inode.use) { 426db4efbfdSBarry Smith C->ops->solve = MatSolve_Inode; 427db4efbfdSBarry Smith } else { 4288d8c7c43SBarry Smith PetscTruth row_identity, col_identity; 4298d8c7c43SBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 4308d8c7c43SBarry Smith ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 4318d8c7c43SBarry Smith if (row_identity && col_identity) { 4328d8c7c43SBarry Smith C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 4338d8c7c43SBarry Smith } else { 434db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqAIJ; 435db4efbfdSBarry Smith } 4368d8c7c43SBarry Smith } 437db4efbfdSBarry Smith C->ops->solveadd = MatSolveAdd_SeqAIJ; 438db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 439db4efbfdSBarry Smith C->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ; 440de26497bSBarry Smith C->ops->matsolve = MatMatSolve_SeqAIJ; 441c456f294SBarry Smith C->assembled = PETSC_TRUE; 4425c9eb25fSBarry Smith C->preallocated = PETSC_TRUE; 443d0f46423SBarry Smith ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 444d2276718SHong Zhang if (sctx.nshift){ 445e738e422SBarry Smith if (info->shiftpd) { 4460481f469SBarry Smith ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 447e738e422SBarry Smith } else if (info->shiftnz) { 448e738e422SBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 4496cc28720Svictorle } 4500697b6caSHong Zhang } 4513a40ed3dSBarry Smith PetscFunctionReturn(0); 452289bc588SBarry Smith } 4530889b5dcSvictorle 454137fb511SHong Zhang /* 455137fb511SHong Zhang This routine implements inplace ILU(0) with row or/and column permutations. 456137fb511SHong Zhang Input: 457137fb511SHong Zhang A - original matrix 458137fb511SHong Zhang Output; 459137fb511SHong Zhang A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i] 460137fb511SHong Zhang a->j (col index) is permuted by the inverse of colperm, then sorted 461137fb511SHong Zhang a->a reordered accordingly with a->j 462137fb511SHong Zhang a->diag (ptr to diagonal elements) is updated. 463137fb511SHong Zhang */ 464137fb511SHong Zhang #undef __FUNCT__ 465137fb511SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm" 4660481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info) 467137fb511SHong Zhang { 468b051339dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 469b051339dSHong Zhang IS isrow = a->row,isicol = a->icol; 470137fb511SHong Zhang PetscErrorCode ierr; 4715d0c19d7SBarry Smith const PetscInt *r,*ic,*ics; 4725d0c19d7SBarry Smith PetscInt i,j,n=A->rmap->n,*ai=a->i,*aj=a->j; 4735d0c19d7SBarry Smith PetscInt *ajtmp,nz,row; 474b051339dSHong Zhang PetscInt *diag = a->diag,nbdiag,*pj; 475a77337e4SBarry Smith PetscScalar *rtmp,*pc,multiplier,d; 476a77337e4SBarry Smith MatScalar *v,*pv; 477137fb511SHong Zhang PetscReal rs; 478137fb511SHong Zhang LUShift_Ctx sctx; 479b051339dSHong Zhang PetscInt newshift; 480137fb511SHong Zhang 481137fb511SHong Zhang PetscFunctionBegin; 482719d5645SBarry Smith if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address"); 483137fb511SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 484137fb511SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 485137fb511SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 486137fb511SHong Zhang ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 487b051339dSHong Zhang ics = ic; 488137fb511SHong Zhang 489137fb511SHong Zhang sctx.shift_top = 0; 490137fb511SHong Zhang sctx.nshift_max = 0; 491137fb511SHong Zhang sctx.shift_lo = 0; 492137fb511SHong Zhang sctx.shift_hi = 0; 4930481f469SBarry Smith sctx.shift_fraction = 0; 494137fb511SHong Zhang 495137fb511SHong Zhang /* if both shift schemes are chosen by user, only use info->shiftpd */ 496137fb511SHong Zhang if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 497137fb511SHong Zhang sctx.shift_top = 0; 498137fb511SHong Zhang for (i=0; i<n; i++) { 499137fb511SHong Zhang /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 500b051339dSHong Zhang d = (a->a)[diag[i]]; 501137fb511SHong Zhang rs = -PetscAbsScalar(d) - PetscRealPart(d); 502b051339dSHong Zhang v = a->a+ai[i]; 503b051339dSHong Zhang nz = ai[i+1] - ai[i]; 504137fb511SHong Zhang for (j=0; j<nz; j++) 505137fb511SHong Zhang rs += PetscAbsScalar(v[j]); 506137fb511SHong Zhang if (rs>sctx.shift_top) sctx.shift_top = rs; 507137fb511SHong Zhang } 508137fb511SHong Zhang if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot; 509137fb511SHong Zhang sctx.shift_top *= 1.1; 510137fb511SHong Zhang sctx.nshift_max = 5; 511137fb511SHong Zhang sctx.shift_lo = 0.; 512137fb511SHong Zhang sctx.shift_hi = 1.; 513137fb511SHong Zhang } 514137fb511SHong Zhang 515137fb511SHong Zhang sctx.shift_amount = 0; 516137fb511SHong Zhang sctx.nshift = 0; 517137fb511SHong Zhang do { 518137fb511SHong Zhang sctx.lushift = PETSC_FALSE; 519137fb511SHong Zhang for (i=0; i<n; i++){ 520b051339dSHong Zhang /* load in initial unfactored row */ 521b051339dSHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 522b051339dSHong Zhang ajtmp = aj + ai[r[i]]; 523b051339dSHong Zhang v = a->a + ai[r[i]]; 524137fb511SHong Zhang /* sort permuted ajtmp and values v accordingly */ 525b051339dSHong Zhang for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]]; 526137fb511SHong Zhang ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr); 527137fb511SHong Zhang 528b051339dSHong Zhang diag[r[i]] = ai[r[i]]; 529137fb511SHong Zhang for (j=0; j<nz; j++) { 530137fb511SHong Zhang rtmp[ajtmp[j]] = v[j]; 531b051339dSHong Zhang if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */ 532137fb511SHong Zhang } 533137fb511SHong Zhang rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 534137fb511SHong Zhang 535b051339dSHong Zhang row = *ajtmp++; 536137fb511SHong Zhang while (row < i) { 537137fb511SHong Zhang pc = rtmp + row; 538137fb511SHong Zhang if (*pc != 0.0) { 539b051339dSHong Zhang pv = a->a + diag[r[row]]; 540b051339dSHong Zhang pj = aj + diag[r[row]] + 1; 541137fb511SHong Zhang 542137fb511SHong Zhang multiplier = *pc / *pv++; 543137fb511SHong Zhang *pc = multiplier; 544b051339dSHong Zhang nz = ai[r[row]+1] - diag[r[row]] - 1; 545b051339dSHong Zhang for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 546dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 547137fb511SHong Zhang } 548b051339dSHong Zhang row = *ajtmp++; 549137fb511SHong Zhang } 550b051339dSHong Zhang /* finished row so overwrite it onto a->a */ 551b051339dSHong Zhang pv = a->a + ai[r[i]] ; 552b051339dSHong Zhang pj = aj + ai[r[i]] ; 553b051339dSHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 554b051339dSHong Zhang nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */ 555137fb511SHong Zhang 556137fb511SHong Zhang rs = 0.0; 557137fb511SHong Zhang for (j=0; j<nz; j++) { 558b051339dSHong Zhang pv[j] = rtmp[pj[j]]; 559b051339dSHong Zhang if (j != nbdiag) rs += PetscAbsScalar(pv[j]); 560137fb511SHong Zhang } 561137fb511SHong Zhang 562137fb511SHong Zhang /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 563137fb511SHong Zhang sctx.rs = rs; 564b051339dSHong Zhang sctx.pv = pv[nbdiag]; 565c6b1f410SHong Zhang ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 566137fb511SHong Zhang if (newshift == 1) break; 567137fb511SHong Zhang } 568137fb511SHong Zhang 5690481f469SBarry Smith if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 570137fb511SHong Zhang /* 571137fb511SHong Zhang * if no shift in this attempt & shifting & started shifting & can refine, 572137fb511SHong Zhang * then try lower shift 573137fb511SHong Zhang */ 5740481f469SBarry Smith sctx.shift_hi = sctx.shift_fraction; 5750481f469SBarry Smith sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 5760481f469SBarry Smith sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 577137fb511SHong Zhang sctx.lushift = PETSC_TRUE; 578137fb511SHong Zhang sctx.nshift++; 579137fb511SHong Zhang } 580137fb511SHong Zhang } while (sctx.lushift); 581137fb511SHong Zhang 582137fb511SHong Zhang /* invert diagonal entries for simplier triangular solves */ 583137fb511SHong Zhang for (i=0; i<n; i++) { 584b051339dSHong Zhang a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]]; 585137fb511SHong Zhang } 586137fb511SHong Zhang 587137fb511SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 588137fb511SHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 589137fb511SHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 590b051339dSHong Zhang A->ops->solve = MatSolve_SeqAIJ_InplaceWithPerm; 591db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqAIJ; 592db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 593db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ; 594b051339dSHong Zhang A->assembled = PETSC_TRUE; 5955c9eb25fSBarry Smith A->preallocated = PETSC_TRUE; 596d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr); 597137fb511SHong Zhang if (sctx.nshift){ 598e738e422SBarry Smith if (info->shiftpd) { 5990481f469SBarry Smith ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 600e738e422SBarry Smith } else if (info->shiftnz) { 601e738e422SBarry Smith ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 602137fb511SHong Zhang } 603137fb511SHong Zhang } 604137fb511SHong Zhang PetscFunctionReturn(0); 605137fb511SHong Zhang } 606137fb511SHong Zhang 6070a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 6084a2ae208SSatish Balay #undef __FUNCT__ 6094a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqAIJ" 6100481f469SBarry Smith PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 611da3a660dSBarry Smith { 612dfbe8321SBarry Smith PetscErrorCode ierr; 613416022c9SBarry Smith Mat C; 614416022c9SBarry Smith 6153a40ed3dSBarry Smith PetscFunctionBegin; 61643244d56SBarry Smith ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 617719d5645SBarry Smith ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 618719d5645SBarry Smith ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 619db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 620db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 621273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 62252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 624da3a660dSBarry Smith } 6250a6dbcc7SLois Curfman McInnes /* ----------------------------------------------------------- */ 6261d098868SHong Zhang 6274a2ae208SSatish Balay #undef __FUNCT__ 6284a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ" 629dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 6308c37ef55SBarry Smith { 631416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 632416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 6336849ba73SBarry Smith PetscErrorCode ierr; 6345d0c19d7SBarry Smith PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 6355d0c19d7SBarry Smith PetscInt nz; 6365d0c19d7SBarry Smith const PetscInt *rout,*cout,*r,*c; 637dd6ea824SBarry Smith PetscScalar *x,*tmp,*tmps,sum; 638d9fead3dSBarry Smith const PetscScalar *b; 639dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 6408c37ef55SBarry Smith 6413a40ed3dSBarry Smith PetscFunctionBegin; 6423a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 64391bf9adeSBarry Smith 644d9fead3dSBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 6451ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 646416022c9SBarry Smith tmp = a->solve_work; 6478c37ef55SBarry Smith 6483b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 6493b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 6508c37ef55SBarry Smith 6518c37ef55SBarry Smith /* forward solve the lower triangular */ 6528c37ef55SBarry Smith tmp[0] = b[*r++]; 653010a20caSHong Zhang tmps = tmp; 6548c37ef55SBarry Smith for (i=1; i<n; i++) { 655010a20caSHong Zhang v = aa + ai[i] ; 656010a20caSHong Zhang vi = aj + ai[i] ; 65753e7425aSBarry Smith nz = a->diag[i] - ai[i]; 65853e7425aSBarry Smith sum = b[*r++]; 659ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 6608c37ef55SBarry Smith tmp[i] = sum; 6618c37ef55SBarry Smith } 6628c37ef55SBarry Smith 6638c37ef55SBarry Smith /* backward solve the upper triangular */ 6648c37ef55SBarry Smith for (i=n-1; i>=0; i--){ 665010a20caSHong Zhang v = aa + a->diag[i] + 1; 666010a20caSHong Zhang vi = aj + a->diag[i] + 1; 667416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 6688c37ef55SBarry Smith sum = tmp[i]; 669ed480e8bSBarry Smith SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 670010a20caSHong Zhang x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 6718c37ef55SBarry Smith } 6728c37ef55SBarry Smith 6733b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 6743b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 675d9fead3dSBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 6761ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 677dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 6783a40ed3dSBarry Smith PetscFunctionReturn(0); 6798c37ef55SBarry Smith } 680026e39d0SSatish Balay 6812bebee5dSHong Zhang #undef __FUNCT__ 6822bebee5dSHong Zhang #define __FUNCT__ "MatMatSolve_SeqAIJ" 6832bebee5dSHong Zhang PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X) 6842bebee5dSHong Zhang { 6852bebee5dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 6862bebee5dSHong Zhang IS iscol = a->col,isrow = a->row; 6872bebee5dSHong Zhang PetscErrorCode ierr; 6885d0c19d7SBarry Smith PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 6895d0c19d7SBarry Smith PetscInt nz,neq; 6905d0c19d7SBarry Smith const PetscInt *rout,*cout,*r,*c; 691dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,*tmps,sum; 692dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 693db4efbfdSBarry Smith PetscTruth bisdense,xisdense; 6942bebee5dSHong Zhang 6952bebee5dSHong Zhang PetscFunctionBegin; 6962bebee5dSHong Zhang if (!n) PetscFunctionReturn(0); 6972bebee5dSHong Zhang 698db4efbfdSBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr); 699db4efbfdSBarry Smith if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix"); 700db4efbfdSBarry Smith ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr); 701db4efbfdSBarry Smith if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix"); 702db4efbfdSBarry Smith 7032bebee5dSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 7042bebee5dSHong Zhang ierr = MatGetArray(X,&x);CHKERRQ(ierr); 7052bebee5dSHong Zhang 7062bebee5dSHong Zhang tmp = a->solve_work; 7072bebee5dSHong Zhang ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 7082bebee5dSHong Zhang ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 7092bebee5dSHong Zhang 710a36b22ccSSatish Balay for (neq=0; neq<B->cmap->n; neq++){ 7112bebee5dSHong Zhang /* forward solve the lower triangular */ 7122bebee5dSHong Zhang tmp[0] = b[r[0]]; 7132bebee5dSHong Zhang tmps = tmp; 7142bebee5dSHong Zhang for (i=1; i<n; i++) { 7152bebee5dSHong Zhang v = aa + ai[i] ; 7162bebee5dSHong Zhang vi = aj + ai[i] ; 7172bebee5dSHong Zhang nz = a->diag[i] - ai[i]; 7182bebee5dSHong Zhang sum = b[r[i]]; 7192bebee5dSHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 7202bebee5dSHong Zhang tmp[i] = sum; 7212bebee5dSHong Zhang } 7222bebee5dSHong Zhang /* backward solve the upper triangular */ 7232bebee5dSHong Zhang for (i=n-1; i>=0; i--){ 7242bebee5dSHong Zhang v = aa + a->diag[i] + 1; 7252bebee5dSHong Zhang vi = aj + a->diag[i] + 1; 7262bebee5dSHong Zhang nz = ai[i+1] - a->diag[i] - 1; 7272bebee5dSHong Zhang sum = tmp[i]; 7282bebee5dSHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 7292bebee5dSHong Zhang x[c[i]] = tmp[i] = sum*aa[a->diag[i]]; 7302bebee5dSHong Zhang } 7312bebee5dSHong Zhang 7322bebee5dSHong Zhang b += n; 7332bebee5dSHong Zhang x += n; 7342bebee5dSHong Zhang } 7352bebee5dSHong Zhang ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 7362bebee5dSHong Zhang ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 7372bebee5dSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 7382bebee5dSHong Zhang ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 739dc0b31edSSatish Balay ierr = PetscLogFlops(B->cmap->n*(2.0*a->nz - n));CHKERRQ(ierr); 7402bebee5dSHong Zhang PetscFunctionReturn(0); 7412bebee5dSHong Zhang } 7422bebee5dSHong Zhang 743137fb511SHong Zhang #undef __FUNCT__ 744137fb511SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm" 745137fb511SHong Zhang PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx) 746137fb511SHong Zhang { 747137fb511SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 748137fb511SHong Zhang IS iscol = a->col,isrow = a->row; 749137fb511SHong Zhang PetscErrorCode ierr; 7505d0c19d7SBarry Smith const PetscInt *r,*c,*rout,*cout; 7515d0c19d7SBarry Smith PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 7525d0c19d7SBarry Smith PetscInt nz,row; 753dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,*tmps,sum; 754dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 755137fb511SHong Zhang 756137fb511SHong Zhang PetscFunctionBegin; 757137fb511SHong Zhang if (!n) PetscFunctionReturn(0); 758137fb511SHong Zhang 759137fb511SHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 760137fb511SHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 761137fb511SHong Zhang tmp = a->solve_work; 762137fb511SHong Zhang 763137fb511SHong Zhang ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 764137fb511SHong Zhang ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 765137fb511SHong Zhang 766137fb511SHong Zhang /* forward solve the lower triangular */ 767137fb511SHong Zhang tmp[0] = b[*r++]; 768137fb511SHong Zhang tmps = tmp; 769137fb511SHong Zhang for (row=1; row<n; row++) { 770137fb511SHong Zhang i = rout[row]; /* permuted row */ 771137fb511SHong Zhang v = aa + ai[i] ; 772137fb511SHong Zhang vi = aj + ai[i] ; 773137fb511SHong Zhang nz = a->diag[i] - ai[i]; 774137fb511SHong Zhang sum = b[*r++]; 775137fb511SHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 776137fb511SHong Zhang tmp[row] = sum; 777137fb511SHong Zhang } 778137fb511SHong Zhang 779137fb511SHong Zhang /* backward solve the upper triangular */ 780137fb511SHong Zhang for (row=n-1; row>=0; row--){ 781137fb511SHong Zhang i = rout[row]; /* permuted row */ 782137fb511SHong Zhang v = aa + a->diag[i] + 1; 783137fb511SHong Zhang vi = aj + a->diag[i] + 1; 784137fb511SHong Zhang nz = ai[i+1] - a->diag[i] - 1; 785137fb511SHong Zhang sum = tmp[row]; 786137fb511SHong Zhang SPARSEDENSEMDOT(sum,tmps,v,vi,nz); 787137fb511SHong Zhang x[*c--] = tmp[row] = sum*aa[a->diag[i]]; 788137fb511SHong Zhang } 789137fb511SHong Zhang 790137fb511SHong Zhang ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 791137fb511SHong Zhang ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 792137fb511SHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 793137fb511SHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 794dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 795137fb511SHong Zhang PetscFunctionReturn(0); 796137fb511SHong Zhang } 797137fb511SHong Zhang 798930ae53cSBarry Smith /* ----------------------------------------------------------- */ 7994a2ae208SSatish Balay #undef __FUNCT__ 8004a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 801dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 802930ae53cSBarry Smith { 803930ae53cSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 8046849ba73SBarry Smith PetscErrorCode ierr; 805d0f46423SBarry Smith PetscInt n = A->rmap->n; 806356650c2SBarry Smith const PetscInt *ai = a->i,*aj = a->j,*adiag = a->diag,*vi; 807356650c2SBarry Smith PetscScalar *x; 808356650c2SBarry Smith const PetscScalar *b; 809dd6ea824SBarry Smith const MatScalar *aa = a->a; 810aa482453SBarry Smith #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 811356650c2SBarry Smith PetscInt adiag_i,i,nz,ai_i; 812dd6ea824SBarry Smith const MatScalar *v; 813dd6ea824SBarry Smith PetscScalar sum; 814d85afc42SSatish Balay #endif 815930ae53cSBarry Smith 8163a40ed3dSBarry Smith PetscFunctionBegin; 8173a40ed3dSBarry Smith if (!n) PetscFunctionReturn(0); 818930ae53cSBarry Smith 819356650c2SBarry Smith ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 8201ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 821930ae53cSBarry Smith 822aa482453SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 8231c4feecaSSatish Balay fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 8241c4feecaSSatish Balay #else 825930ae53cSBarry Smith /* forward solve the lower triangular */ 826930ae53cSBarry Smith x[0] = b[0]; 827930ae53cSBarry Smith for (i=1; i<n; i++) { 828930ae53cSBarry Smith ai_i = ai[i]; 829930ae53cSBarry Smith v = aa + ai_i; 830930ae53cSBarry Smith vi = aj + ai_i; 831930ae53cSBarry Smith nz = adiag[i] - ai_i; 832930ae53cSBarry Smith sum = b[i]; 833930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 834930ae53cSBarry Smith x[i] = sum; 835930ae53cSBarry Smith } 836930ae53cSBarry Smith 837930ae53cSBarry Smith /* backward solve the upper triangular */ 838930ae53cSBarry Smith for (i=n-1; i>=0; i--){ 839930ae53cSBarry Smith adiag_i = adiag[i]; 840d708749eSBarry Smith v = aa + adiag_i + 1; 841d708749eSBarry Smith vi = aj + adiag_i + 1; 842930ae53cSBarry Smith nz = ai[i+1] - adiag_i - 1; 843930ae53cSBarry Smith sum = x[i]; 844930ae53cSBarry Smith while (nz--) sum -= *v++ * x[*vi++]; 845930ae53cSBarry Smith x[i] = sum*aa[adiag_i]; 846930ae53cSBarry Smith } 8471c4feecaSSatish Balay #endif 848dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 849356650c2SBarry Smith ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 8501ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 8513a40ed3dSBarry Smith PetscFunctionReturn(0); 852930ae53cSBarry Smith } 853930ae53cSBarry Smith 8544a2ae208SSatish Balay #undef __FUNCT__ 8554a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqAIJ" 856dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 857da3a660dSBarry Smith { 858416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 859416022c9SBarry Smith IS iscol = a->col,isrow = a->row; 8606849ba73SBarry Smith PetscErrorCode ierr; 8615d0c19d7SBarry Smith PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 8625d0c19d7SBarry Smith PetscInt nz; 8635d0c19d7SBarry Smith const PetscInt *rout,*cout,*r,*c; 864dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,sum; 865dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 866da3a660dSBarry Smith 8673a40ed3dSBarry Smith PetscFunctionBegin; 86878b31e54SBarry Smith if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 869da3a660dSBarry Smith 8701ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 8711ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 872416022c9SBarry Smith tmp = a->solve_work; 873da3a660dSBarry Smith 8743b2fbd54SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 8753b2fbd54SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 876da3a660dSBarry Smith 877da3a660dSBarry Smith /* forward solve the lower triangular */ 878da3a660dSBarry Smith tmp[0] = b[*r++]; 879da3a660dSBarry Smith for (i=1; i<n; i++) { 880010a20caSHong Zhang v = aa + ai[i] ; 881010a20caSHong Zhang vi = aj + ai[i] ; 882416022c9SBarry Smith nz = a->diag[i] - ai[i]; 883da3a660dSBarry Smith sum = b[*r++]; 884010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 885da3a660dSBarry Smith tmp[i] = sum; 886da3a660dSBarry Smith } 887da3a660dSBarry Smith 888da3a660dSBarry Smith /* backward solve the upper triangular */ 889da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 890010a20caSHong Zhang v = aa + a->diag[i] + 1; 891010a20caSHong Zhang vi = aj + a->diag[i] + 1; 892416022c9SBarry Smith nz = ai[i+1] - a->diag[i] - 1; 893da3a660dSBarry Smith sum = tmp[i]; 894010a20caSHong Zhang while (nz--) sum -= *v++ * tmp[*vi++ ]; 895010a20caSHong Zhang tmp[i] = sum*aa[a->diag[i]]; 896da3a660dSBarry Smith x[*c--] += tmp[i]; 897da3a660dSBarry Smith } 898da3a660dSBarry Smith 8993b2fbd54SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 9003b2fbd54SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 9011ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 9021ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 903dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 904898183ecSLois Curfman McInnes 9053a40ed3dSBarry Smith PetscFunctionReturn(0); 906da3a660dSBarry Smith } 907da3a660dSBarry Smith /* -------------------------------------------------------------------*/ 9084a2ae208SSatish Balay #undef __FUNCT__ 9094a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 910dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 911da3a660dSBarry Smith { 912416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 9132235e667SBarry Smith IS iscol = a->col,isrow = a->row; 9146849ba73SBarry Smith PetscErrorCode ierr; 9155d0c19d7SBarry Smith const PetscInt *rout,*cout,*r,*c; 9165d0c19d7SBarry Smith PetscInt i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 9175d0c19d7SBarry Smith PetscInt nz,*diag = a->diag; 918dd6ea824SBarry Smith PetscScalar *x,*b,*tmp,s1; 919dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 920da3a660dSBarry Smith 9213a40ed3dSBarry Smith PetscFunctionBegin; 9221ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 9231ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 924416022c9SBarry Smith tmp = a->solve_work; 925da3a660dSBarry Smith 9262235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 9272235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 928da3a660dSBarry Smith 929da3a660dSBarry Smith /* copy the b into temp work space according to permutation */ 9302235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 931da3a660dSBarry Smith 932da3a660dSBarry Smith /* forward solve the U^T */ 933da3a660dSBarry Smith for (i=0; i<n; i++) { 934010a20caSHong Zhang v = aa + diag[i] ; 935010a20caSHong Zhang vi = aj + diag[i] + 1; 936f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 937c38d4ed2SBarry Smith s1 = tmp[i]; 938ef66eb69SBarry Smith s1 *= (*v++); /* multiply by inverse of diagonal entry */ 939da3a660dSBarry Smith while (nz--) { 940010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*s1; 941da3a660dSBarry Smith } 942c38d4ed2SBarry Smith tmp[i] = s1; 943da3a660dSBarry Smith } 944da3a660dSBarry Smith 945da3a660dSBarry Smith /* backward solve the L^T */ 946da3a660dSBarry Smith for (i=n-1; i>=0; i--){ 947010a20caSHong Zhang v = aa + diag[i] - 1 ; 948010a20caSHong Zhang vi = aj + diag[i] - 1 ; 949f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 950f1af5d2fSBarry Smith s1 = tmp[i]; 951da3a660dSBarry Smith while (nz--) { 952010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*s1; 953da3a660dSBarry Smith } 954da3a660dSBarry Smith } 955da3a660dSBarry Smith 956da3a660dSBarry Smith /* copy tmp into x according to permutation */ 957da3a660dSBarry Smith for (i=0; i<n; i++) x[r[i]] = tmp[i]; 958da3a660dSBarry Smith 9592235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 9602235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 9611ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 9621ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 963da3a660dSBarry Smith 964dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz-A->cmap->n);CHKERRQ(ierr); 9653a40ed3dSBarry Smith PetscFunctionReturn(0); 966da3a660dSBarry Smith } 967da3a660dSBarry Smith 9684a2ae208SSatish Balay #undef __FUNCT__ 9694a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 970dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 971da3a660dSBarry Smith { 972416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 9732235e667SBarry Smith IS iscol = a->col,isrow = a->row; 9746849ba73SBarry Smith PetscErrorCode ierr; 9755d0c19d7SBarry Smith const PetscInt *r,*c,*rout,*cout; 9765d0c19d7SBarry Smith PetscInt i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 9775d0c19d7SBarry Smith PetscInt nz,*diag = a->diag; 978dd6ea824SBarry Smith PetscScalar *x,*b,*tmp; 979dd6ea824SBarry Smith const MatScalar *aa = a->a,*v; 9806abc6512SBarry Smith 9813a40ed3dSBarry Smith PetscFunctionBegin; 98223bc6035SMatthew Knepley if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 9836abc6512SBarry Smith 9841ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 9851ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 986416022c9SBarry Smith tmp = a->solve_work; 9876abc6512SBarry Smith 9882235e667SBarry Smith ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 9892235e667SBarry Smith ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 9906abc6512SBarry Smith 9916abc6512SBarry Smith /* copy the b into temp work space according to permutation */ 9922235e667SBarry Smith for (i=0; i<n; i++) tmp[i] = b[c[i]]; 9936abc6512SBarry Smith 9946abc6512SBarry Smith /* forward solve the U^T */ 9956abc6512SBarry Smith for (i=0; i<n; i++) { 996010a20caSHong Zhang v = aa + diag[i] ; 997010a20caSHong Zhang vi = aj + diag[i] + 1; 998f1af5d2fSBarry Smith nz = ai[i+1] - diag[i] - 1; 9996abc6512SBarry Smith tmp[i] *= *v++; 10006abc6512SBarry Smith while (nz--) { 1001010a20caSHong Zhang tmp[*vi++ ] -= (*v++)*tmp[i]; 10026abc6512SBarry Smith } 10036abc6512SBarry Smith } 10046abc6512SBarry Smith 10056abc6512SBarry Smith /* backward solve the L^T */ 10066abc6512SBarry Smith for (i=n-1; i>=0; i--){ 1007010a20caSHong Zhang v = aa + diag[i] - 1 ; 1008010a20caSHong Zhang vi = aj + diag[i] - 1 ; 1009f1af5d2fSBarry Smith nz = diag[i] - ai[i]; 10106abc6512SBarry Smith while (nz--) { 1011010a20caSHong Zhang tmp[*vi-- ] -= (*v--)*tmp[i]; 10126abc6512SBarry Smith } 10136abc6512SBarry Smith } 10146abc6512SBarry Smith 10156abc6512SBarry Smith /* copy tmp into x according to permutation */ 10166abc6512SBarry Smith for (i=0; i<n; i++) x[r[i]] += tmp[i]; 10176abc6512SBarry Smith 10182235e667SBarry Smith ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 10192235e667SBarry Smith ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 10201ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 10211ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 10226abc6512SBarry Smith 1023dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 10243a40ed3dSBarry Smith PetscFunctionReturn(0); 1025da3a660dSBarry Smith } 10269e25ed09SBarry Smith /* ----------------------------------------------------------------*/ 10273c5fc038SBarry Smith EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth); 1028719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption); 102908480c60SBarry Smith 10304a2ae208SSatish Balay #undef __FUNCT__ 10314a2ae208SSatish Balay #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 10320481f469SBarry Smith PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 10339e25ed09SBarry Smith { 1034416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 10359e25ed09SBarry Smith IS isicol; 10366849ba73SBarry Smith PetscErrorCode ierr; 10375d0c19d7SBarry Smith const PetscInt *r,*ic; 10385d0c19d7SBarry Smith PetscInt n=A->rmap->n,*ai=a->i,*aj=a->j,d; 10395a8e39fbSHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 1040052f0c41SBarry Smith PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 10415a8e39fbSHong Zhang PetscInt i,levels,diagonal_fill; 104277c4ece6SBarry Smith PetscTruth col_identity,row_identity; 1043329f5518SBarry Smith PetscReal f; 10445a8e39fbSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 10455a8e39fbSHong Zhang PetscBT lnkbt; 10465a8e39fbSHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 1047a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1048a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 104909f38230SBarry Smith PetscTruth missing; 10509e25ed09SBarry Smith 10513a40ed3dSBarry Smith PetscFunctionBegin; 1052d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 10536cf997b0SBarry Smith f = info->fill; 105438baddfdSBarry Smith levels = (PetscInt)info->levels; 105538baddfdSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 10564c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 105782bf6240SBarry Smith 105808480c60SBarry Smith /* special case that simply copies fill pattern */ 1059be0abb6dSBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1060be0abb6dSBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 106186bacbd2SBarry Smith if (!levels && row_identity && col_identity) { 1062719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr); 1063719d5645SBarry Smith fact->factor = MAT_FACTOR_ILU; 1064719d5645SBarry Smith (fact)->info.factor_mallocs = 0; 1065719d5645SBarry Smith (fact)->info.fill_ratio_given = info->fill; 1066719d5645SBarry Smith (fact)->info.fill_ratio_needed = 1.0; 1067719d5645SBarry Smith b = (Mat_SeqAIJ*)(fact)->data; 10688ef3462fSBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 106909f38230SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 107008480c60SBarry Smith b->row = isrow; 107108480c60SBarry Smith b->col = iscol; 107282bf6240SBarry Smith b->icol = isicol; 1073719d5645SBarry Smith ierr = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1074c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1075c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1076719d5645SBarry Smith (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1077719d5645SBarry Smith ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 10783a40ed3dSBarry Smith PetscFunctionReturn(0); 107908480c60SBarry Smith } 108008480c60SBarry Smith 1081898183ecSLois Curfman McInnes ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1082898183ecSLois Curfman McInnes ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 10839e25ed09SBarry Smith 10849e25ed09SBarry Smith /* get new row pointers */ 10855a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 10865a8e39fbSHong Zhang bi[0] = 0; 10875a8e39fbSHong Zhang /* bdiag is location of diagonal in factor */ 10885a8e39fbSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 10895a8e39fbSHong Zhang bdiag[0] = 0; 10906cf997b0SBarry Smith 10915a8e39fbSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 10925a8e39fbSHong Zhang bjlvl_ptr = (PetscInt**)(bj_ptr + n); 10935a8e39fbSHong Zhang 10945a8e39fbSHong Zhang /* create a linked list for storing column indices of the active row */ 10955a8e39fbSHong Zhang nlnk = n + 1; 10965a8e39fbSHong Zhang ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 10975a8e39fbSHong Zhang 10985a8e39fbSHong Zhang /* initial FreeSpace size is f*(ai[n]+1) */ 1099a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 11005a8e39fbSHong Zhang current_space = free_space; 1101a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 11025a8e39fbSHong Zhang current_space_lvl = free_space_lvl; 11035a8e39fbSHong Zhang 11045a8e39fbSHong Zhang for (i=0; i<n; i++) { 11055a8e39fbSHong Zhang nzi = 0; 11066cf997b0SBarry Smith /* copy current row into linked list */ 11075a8e39fbSHong Zhang nnz = ai[r[i]+1] - ai[r[i]]; 11083b4a8b6dSBarry Smith if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 11095a8e39fbSHong Zhang cols = aj + ai[r[i]]; 11105a8e39fbSHong Zhang lnk[i] = -1; /* marker to indicate if diagonal exists */ 11115a8e39fbSHong Zhang ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 11125a8e39fbSHong Zhang nzi += nlnk; 11136cf997b0SBarry Smith 11146cf997b0SBarry Smith /* make sure diagonal entry is included */ 11155a8e39fbSHong Zhang if (diagonal_fill && lnk[i] == -1) { 11166cf997b0SBarry Smith fm = n; 11175a8e39fbSHong Zhang while (lnk[fm] < i) fm = lnk[fm]; 11185a8e39fbSHong Zhang lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 11195a8e39fbSHong Zhang lnk[fm] = i; 11205a8e39fbSHong Zhang lnk_lvl[i] = 0; 11215a8e39fbSHong Zhang nzi++; dcount++; 11226cf997b0SBarry Smith } 11236cf997b0SBarry Smith 11245a8e39fbSHong Zhang /* add pivot rows into the active row */ 11255a8e39fbSHong Zhang nzbd = 0; 11265a8e39fbSHong Zhang prow = lnk[n]; 11275a8e39fbSHong Zhang while (prow < i) { 11285a8e39fbSHong Zhang nnz = bdiag[prow]; 11295a8e39fbSHong Zhang cols = bj_ptr[prow] + nnz + 1; 11305a8e39fbSHong Zhang cols_lvl = bjlvl_ptr[prow] + nnz + 1; 11315a8e39fbSHong Zhang nnz = bi[prow+1] - bi[prow] - nnz - 1; 11325a8e39fbSHong Zhang ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 11335a8e39fbSHong Zhang nzi += nlnk; 11345a8e39fbSHong Zhang prow = lnk[prow]; 11355a8e39fbSHong Zhang nzbd++; 113656beaf04SBarry Smith } 11375a8e39fbSHong Zhang bdiag[i] = nzbd; 11385a8e39fbSHong Zhang bi[i+1] = bi[i] + nzi; 1139ecf371e4SBarry Smith 11405a8e39fbSHong Zhang /* if free space is not available, make more free space */ 11415a8e39fbSHong Zhang if (current_space->local_remaining<nzi) { 11425a8e39fbSHong Zhang nnz = nzi*(n - i); /* estimated and max additional space needed */ 1143a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 1144a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 11455a8e39fbSHong Zhang reallocs++; 11465a8e39fbSHong Zhang } 1147ecf371e4SBarry Smith 11485a8e39fbSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 11495a8e39fbSHong Zhang ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 11505a8e39fbSHong Zhang bj_ptr[i] = current_space->array; 11515a8e39fbSHong Zhang bjlvl_ptr[i] = current_space_lvl->array; 11525a8e39fbSHong Zhang 11535a8e39fbSHong Zhang /* make sure the active row i has diagonal entry */ 11545a8e39fbSHong Zhang if (*(bj_ptr[i]+bdiag[i]) != i) { 115577431f27SBarry Smith SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 1156d7ee0231SBarry Smith try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 11576cf997b0SBarry Smith } 11585a8e39fbSHong Zhang 11595a8e39fbSHong Zhang current_space->array += nzi; 11605a8e39fbSHong Zhang current_space->local_used += nzi; 11615a8e39fbSHong Zhang current_space->local_remaining -= nzi; 11625a8e39fbSHong Zhang current_space_lvl->array += nzi; 11635a8e39fbSHong Zhang current_space_lvl->local_used += nzi; 11645a8e39fbSHong Zhang current_space_lvl->local_remaining -= nzi; 11659e25ed09SBarry Smith } 11665a8e39fbSHong Zhang 1167898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1168898183ecSLois Curfman McInnes ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 11695a8e39fbSHong Zhang 11705a8e39fbSHong Zhang /* destroy list of free space and other temporary arrays */ 11715a8e39fbSHong Zhang ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1172a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 11735a8e39fbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1174a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 11755a8e39fbSHong Zhang ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 11769e25ed09SBarry Smith 11776cf91177SBarry Smith #if defined(PETSC_USE_INFO) 1178f64065bbSBarry Smith { 11795a8e39fbSHong Zhang PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1180ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 1181ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1182ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr); 1183ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 1184335d9088SBarry Smith if (diagonal_fill) { 1185ae15b995SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr); 1186335d9088SBarry Smith } 1187f64065bbSBarry Smith } 118863ba0a88SBarry Smith #endif 1189bbb6d6a8SBarry Smith 11909e25ed09SBarry Smith /* put together the new matrix */ 119107dbf95fSBarry Smith ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1192719d5645SBarry Smith ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr); 1193719d5645SBarry Smith b = (Mat_SeqAIJ*)(fact)->data; 1194e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1195e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 11967c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 1197052f0c41SBarry Smith ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 11985a8e39fbSHong Zhang b->j = bj; 11995a8e39fbSHong Zhang b->i = bi; 12005a8e39fbSHong Zhang for (i=0; i<n; i++) bdiag[i] += bi[i]; 12015a8e39fbSHong Zhang b->diag = bdiag; 1202416022c9SBarry Smith b->ilen = 0; 1203416022c9SBarry Smith b->imax = 0; 1204416022c9SBarry Smith b->row = isrow; 1205416022c9SBarry Smith b->col = iscol; 1206c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1207c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 120882bf6240SBarry Smith b->icol = isicol; 120987828ca2SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1210416022c9SBarry Smith /* In b structure: Free imax, ilen, old a, old j. 12115a8e39fbSHong Zhang Allocate bdiag, solve_work, new a, new j */ 1212719d5645SBarry Smith ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 12135a8e39fbSHong Zhang b->maxnz = b->nz = bi[n] ; 1214719d5645SBarry Smith (fact)->info.factor_mallocs = reallocs; 1215719d5645SBarry Smith (fact)->info.fill_ratio_given = f; 1216719d5645SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1217719d5645SBarry Smith (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1218719d5645SBarry Smith ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 12193a40ed3dSBarry Smith PetscFunctionReturn(0); 12209e25ed09SBarry Smith } 1221d5d45c9bSBarry Smith 12227c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 1223a6175056SHong Zhang #undef __FUNCT__ 1224a6175056SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 12250481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info) 1226a6175056SHong Zhang { 1227719d5645SBarry Smith Mat C = B; 12280968510aSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 12292ed133dbSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 12309bfd6278SHong Zhang IS ip=b->row,iip = b->icol; 12312ed133dbSHong Zhang PetscErrorCode ierr; 12325d0c19d7SBarry Smith const PetscInt *rip,*riip; 12335d0c19d7SBarry Smith PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol; 12342ed133dbSHong Zhang PetscInt *ai=a->i,*aj=a->j; 12352ed133dbSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1236622e2028SHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1237540703acSHong Zhang PetscReal zeropivot,rs,shiftnz; 1238fbf22428SSatish Balay PetscReal shiftpd; 1239540703acSHong Zhang ChShift_Ctx sctx; 1240540703acSHong Zhang PetscInt newshift; 1241db4efbfdSBarry Smith PetscTruth perm_identity; 1242d5d45c9bSBarry Smith 1243a6175056SHong Zhang PetscFunctionBegin; 1244117f1891Stmunson 1245540703acSHong Zhang shiftnz = info->shiftnz; 1246540703acSHong Zhang shiftpd = info->shiftpd; 1247ee45ca4aSHong Zhang zeropivot = info->zeropivot; 1248ee45ca4aSHong Zhang 12492ed133dbSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 12509bfd6278SHong Zhang ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 12512ed133dbSHong Zhang 12522ed133dbSHong Zhang /* initialization */ 12532ed133dbSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 12542ed133dbSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 12552ed133dbSHong Zhang jl = il + mbs; 12562ed133dbSHong Zhang rtmp = (MatScalar*)(jl + mbs); 12572ed133dbSHong Zhang 1258540703acSHong Zhang sctx.shift_amount = 0; 1259540703acSHong Zhang sctx.nshift = 0; 12602ed133dbSHong Zhang do { 1261540703acSHong Zhang sctx.chshift = PETSC_FALSE; 12622ed133dbSHong Zhang for (i=0; i<mbs; i++) { 12632ed133dbSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1264a6175056SHong Zhang } 12652ed133dbSHong Zhang 12662ed133dbSHong Zhang for (k = 0; k<mbs; k++){ 1267c05c3958SHong Zhang bval = ba + bi[k]; 12682ed133dbSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 12692ed133dbSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 12702ed133dbSHong Zhang for (j = jmin; j < jmax; j++){ 12719bfd6278SHong Zhang col = riip[aj[j]]; 12722ed133dbSHong Zhang if (col >= k){ /* only take upper triangular entry */ 12732ed133dbSHong Zhang rtmp[col] = aa[j]; 1274c05c3958SHong Zhang *bval++ = 0.0; /* for in-place factorization */ 12752ed133dbSHong Zhang } 12762ed133dbSHong Zhang } 127739e3d630SHong Zhang /* shift the diagonal of the matrix */ 1278540703acSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 12792ed133dbSHong Zhang 12802ed133dbSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 12812ed133dbSHong Zhang dk = rtmp[k]; 12822ed133dbSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 12832ed133dbSHong Zhang 12842ed133dbSHong Zhang while (i < k){ 12852ed133dbSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 12862ed133dbSHong Zhang 12872ed133dbSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 12882ed133dbSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 12892ed133dbSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 12902ed133dbSHong Zhang dk += uikdi*ba[ili]; 12912ed133dbSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 12922ed133dbSHong Zhang 12932ed133dbSHong Zhang /* add multiple of row i to k-th row */ 12942ed133dbSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 12952ed133dbSHong Zhang if (jmin < jmax){ 12962ed133dbSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 12972ed133dbSHong Zhang /* update il and jl for row i */ 12982ed133dbSHong Zhang il[i] = jmin; 12992ed133dbSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 13002ed133dbSHong Zhang } 13012ed133dbSHong Zhang i = nexti; 13022ed133dbSHong Zhang } 13032ed133dbSHong Zhang 13040697b6caSHong Zhang /* shift the diagonals when zero pivot is detected */ 13050697b6caSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 13060697b6caSHong Zhang rs = 0.0; 13073c9e8598SHong Zhang jmin = bi[k]+1; 13083c9e8598SHong Zhang nz = bi[k+1] - jmin; 13093c9e8598SHong Zhang bcol = bj + jmin; 13103c9e8598SHong Zhang while (nz--){ 131189f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 131289f845c8SHong Zhang bcol++; 13133c9e8598SHong Zhang } 1314540703acSHong Zhang 1315540703acSHong Zhang sctx.rs = rs; 1316540703acSHong Zhang sctx.pv = dk; 13175b5bc046SBarry Smith ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 1318117f1891Stmunson 1319117f1891Stmunson if (newshift == 1) { 1320117f1891Stmunson if (!sctx.shift_amount) { 1321117f1891Stmunson sctx.shift_amount = 1e-5; 1322117f1891Stmunson } 1323117f1891Stmunson break; 1324117f1891Stmunson } 13253c9e8598SHong Zhang 13263c9e8598SHong Zhang /* copy data into U(k,:) */ 132739e3d630SHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 132839e3d630SHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 132939e3d630SHong Zhang if (jmin < jmax) { 133039e3d630SHong Zhang for (j=jmin; j<jmax; j++){ 133139e3d630SHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 13323c9e8598SHong Zhang } 133339e3d630SHong Zhang /* add the k-th row into il and jl */ 13343c9e8598SHong Zhang il[k] = jmin; 13353c9e8598SHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 13363c9e8598SHong Zhang } 13373c9e8598SHong Zhang } 1338540703acSHong Zhang } while (sctx.chshift); 13393c9e8598SHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 13403c9e8598SHong Zhang 134139e3d630SHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 13429bfd6278SHong Zhang ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 1343db4efbfdSBarry Smith 1344db4efbfdSBarry Smith ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 1345db4efbfdSBarry Smith if (perm_identity){ 1346719d5645SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1347719d5645SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1348719d5645SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1349719d5645SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1350db4efbfdSBarry Smith } else { 1351719d5645SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1; 1352719d5645SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 1353719d5645SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 1354719d5645SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 1355db4efbfdSBarry Smith } 1356db4efbfdSBarry Smith 13573c9e8598SHong Zhang C->assembled = PETSC_TRUE; 13583c9e8598SHong Zhang C->preallocated = PETSC_TRUE; 1359d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 1360540703acSHong Zhang if (sctx.nshift){ 1361540703acSHong Zhang if (shiftnz) { 13621e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1363540703acSHong Zhang } else if (shiftpd) { 13641e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 13650697b6caSHong Zhang } 13663c9e8598SHong Zhang } 13673c9e8598SHong Zhang PetscFunctionReturn(0); 13683c9e8598SHong Zhang } 1369a6175056SHong Zhang 1370a6175056SHong Zhang #undef __FUNCT__ 1371a6175056SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 13720481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1373a6175056SHong Zhang { 13740968510aSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1375ed59401aSHong Zhang Mat_SeqSBAIJ *b; 1376dfbe8321SBarry Smith PetscErrorCode ierr; 137758ebbce7SBarry Smith PetscTruth perm_identity,missing; 13785d0c19d7SBarry Smith PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui; 13795d0c19d7SBarry Smith const PetscInt *rip,*riip; 1380ed59401aSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 138158ebbce7SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,d; 13825a8e39fbSHong Zhang PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1383ed59401aSHong Zhang PetscReal fill=info->fill,levels=info->levels; 1384a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1385a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 13867a48dd6fSHong Zhang PetscBT lnkbt; 1387b635d36bSHong Zhang IS iperm; 1388a6175056SHong Zhang 1389a6175056SHong Zhang PetscFunctionBegin; 1390d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 139158ebbce7SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 139258ebbce7SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 1393653a6975SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1394b635d36bSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 13957a48dd6fSHong Zhang 139639e3d630SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 139739e3d630SHong Zhang ui[0] = 0; 139839e3d630SHong Zhang 1399b635d36bSHong Zhang /* ICC(0) without matrix ordering: simply copies fill pattern */ 1400622e2028SHong Zhang if (!levels && perm_identity) { 140158ebbce7SBarry Smith 1402ed59401aSHong Zhang for (i=0; i<am; i++) { 140339e3d630SHong Zhang ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1404ed59401aSHong Zhang } 140539e3d630SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 140639e3d630SHong Zhang cols = uj; 1407ed59401aSHong Zhang for (i=0; i<am; i++) { 1408ed59401aSHong Zhang aj = a->j + a->diag[i]; 140939e3d630SHong Zhang ncols = ui[i+1] - ui[i]; 141039e3d630SHong Zhang for (j=0; j<ncols; j++) *cols++ = *aj++; 1411ed59401aSHong Zhang } 141239e3d630SHong Zhang } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 1413b635d36bSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 1414b635d36bSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1415b635d36bSHong Zhang 14167a48dd6fSHong Zhang /* initialization */ 14175a8e39fbSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 14187a48dd6fSHong Zhang 14197a48dd6fSHong Zhang /* jl: linked list for storing indices of the pivot rows 14207a48dd6fSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 14211393bc97SHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 14227a48dd6fSHong Zhang il = jl + am; 14237a48dd6fSHong Zhang uj_ptr = (PetscInt**)(il + am); 14247a48dd6fSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 14257a48dd6fSHong Zhang for (i=0; i<am; i++){ 14267a48dd6fSHong Zhang jl[i] = am; il[i] = 0; 14277a48dd6fSHong Zhang } 14287a48dd6fSHong Zhang 14297a48dd6fSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 14307a48dd6fSHong Zhang nlnk = am + 1; 14312ed133dbSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 14327a48dd6fSHong Zhang 14337a48dd6fSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1434a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 14357a48dd6fSHong Zhang current_space = free_space; 1436a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 14377a48dd6fSHong Zhang current_space_lvl = free_space_lvl; 14387a48dd6fSHong Zhang 14397a48dd6fSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 14407a48dd6fSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 14417a48dd6fSHong Zhang nzk = 0; 1442622e2028SHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 14433b4a8b6dSBarry Smith if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 1444622e2028SHong Zhang ncols_upper = 0; 1445622e2028SHong Zhang for (j=0; j<ncols; j++){ 1446b635d36bSHong Zhang i = *(aj + ai[rip[k]] + j); /* unpermuted column index */ 1447b635d36bSHong Zhang if (riip[i] >= k){ /* only take upper triangular entry */ 14485a8e39fbSHong Zhang ajtmp[ncols_upper] = i; 1449622e2028SHong Zhang ncols_upper++; 1450622e2028SHong Zhang } 1451622e2028SHong Zhang } 1452b635d36bSHong Zhang ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 14537a48dd6fSHong Zhang nzk += nlnk; 14547a48dd6fSHong Zhang 14557a48dd6fSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 14567a48dd6fSHong Zhang prow = jl[k]; /* 1st pivot row */ 14577a48dd6fSHong Zhang 14587a48dd6fSHong Zhang while (prow < k){ 14597a48dd6fSHong Zhang nextprow = jl[prow]; 14607a48dd6fSHong Zhang 14617a48dd6fSHong Zhang /* merge prow into k-th row */ 14627a48dd6fSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 14637a48dd6fSHong Zhang jmax = ui[prow+1]; 14647a48dd6fSHong Zhang ncols = jmax-jmin; 1465ed59401aSHong Zhang i = jmin - ui[prow]; 1466ed59401aSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 14675a8e39fbSHong Zhang uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 14685a8e39fbSHong Zhang j = *(uj - 1); 14695a8e39fbSHong Zhang ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 14707a48dd6fSHong Zhang nzk += nlnk; 14717a48dd6fSHong Zhang 14727a48dd6fSHong Zhang /* update il and jl for prow */ 14737a48dd6fSHong Zhang if (jmin < jmax){ 14747a48dd6fSHong Zhang il[prow] = jmin; 14757a48dd6fSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 14767a48dd6fSHong Zhang } 14777a48dd6fSHong Zhang prow = nextprow; 14787a48dd6fSHong Zhang } 14797a48dd6fSHong Zhang 14807a48dd6fSHong Zhang /* if free space is not available, make more free space */ 14817a48dd6fSHong Zhang if (current_space->local_remaining<nzk) { 14827a48dd6fSHong Zhang i = am - k + 1; /* num of unfactored rows */ 14837a48dd6fSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1484a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1485a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 14867a48dd6fSHong Zhang reallocs++; 14877a48dd6fSHong Zhang } 14887a48dd6fSHong Zhang 14897a48dd6fSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 1490b635d36bSHong Zhang if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k); 14912ed133dbSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 14927a48dd6fSHong Zhang 14937a48dd6fSHong Zhang /* add the k-th row into il and jl */ 149439e3d630SHong Zhang if (nzk > 1){ 14957a48dd6fSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 14967a48dd6fSHong Zhang jl[k] = jl[i]; jl[i] = k; 14977a48dd6fSHong Zhang il[k] = ui[k] + 1; 14987a48dd6fSHong Zhang } 14997a48dd6fSHong Zhang uj_ptr[k] = current_space->array; 15007a48dd6fSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 15017a48dd6fSHong Zhang 15027a48dd6fSHong Zhang current_space->array += nzk; 15037a48dd6fSHong Zhang current_space->local_used += nzk; 15047a48dd6fSHong Zhang current_space->local_remaining -= nzk; 15057a48dd6fSHong Zhang 15067a48dd6fSHong Zhang current_space_lvl->array += nzk; 15077a48dd6fSHong Zhang current_space_lvl->local_used += nzk; 15087a48dd6fSHong Zhang current_space_lvl->local_remaining -= nzk; 15097a48dd6fSHong Zhang 15107a48dd6fSHong Zhang ui[k+1] = ui[k] + nzk; 15117a48dd6fSHong Zhang } 15127a48dd6fSHong Zhang 15136cf91177SBarry Smith #if defined(PETSC_USE_INFO) 15147a48dd6fSHong Zhang if (ai[am] != 0) { 151539e3d630SHong Zhang PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 1516ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1517ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1518ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 15197a48dd6fSHong Zhang } else { 1520ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 15217a48dd6fSHong Zhang } 152263ba0a88SBarry Smith #endif 15237a48dd6fSHong Zhang 15247a48dd6fSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1525b635d36bSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 15267a48dd6fSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 15275a8e39fbSHong Zhang ierr = PetscFree(ajtmp);CHKERRQ(ierr); 15287a48dd6fSHong Zhang 15297a48dd6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 15307a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1531a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 15322ed133dbSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1533a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 15347a48dd6fSHong Zhang 153539e3d630SHong Zhang } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 153639e3d630SHong Zhang 15377a48dd6fSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 15387a48dd6fSHong Zhang 1539719d5645SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 15407a48dd6fSHong Zhang b->singlemalloc = PETSC_FALSE; 15417a48dd6fSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 15427a48dd6fSHong Zhang b->j = uj; 15437a48dd6fSHong Zhang b->i = ui; 15447a48dd6fSHong Zhang b->diag = 0; 15457a48dd6fSHong Zhang b->ilen = 0; 15467a48dd6fSHong Zhang b->imax = 0; 15477a48dd6fSHong Zhang b->row = perm; 1548b635d36bSHong Zhang b->col = perm; 1549b635d36bSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1550b635d36bSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1551b635d36bSHong Zhang b->icol = iperm; 15527a48dd6fSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 15537a48dd6fSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1554719d5645SBarry Smith ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 15557a48dd6fSHong Zhang b->maxnz = b->nz = ui[am]; 1556e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1557e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 15587a48dd6fSHong Zhang 1559719d5645SBarry Smith (fact)->info.factor_mallocs = reallocs; 1560719d5645SBarry Smith (fact)->info.fill_ratio_given = fill; 15617a48dd6fSHong Zhang if (ai[am] != 0) { 1562719d5645SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 15637a48dd6fSHong Zhang } else { 1564719d5645SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 15657a48dd6fSHong Zhang } 1566719d5645SBarry Smith (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1567a6175056SHong Zhang PetscFunctionReturn(0); 1568a6175056SHong Zhang } 1569d5d45c9bSBarry Smith 1570f76d2b81SHong Zhang #undef __FUNCT__ 1571f76d2b81SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 15720481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1573f76d2b81SHong Zhang { 1574f76d2b81SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 157510c27e3dSHong Zhang Mat_SeqSBAIJ *b; 1576dfbe8321SBarry Smith PetscErrorCode ierr; 1577f76d2b81SHong Zhang PetscTruth perm_identity; 157810c27e3dSHong Zhang PetscReal fill = info->fill; 15795d0c19d7SBarry Smith const PetscInt *rip,*riip; 15805d0c19d7SBarry Smith PetscInt i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow; 158110c27e3dSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 15822ed133dbSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1583a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 158410c27e3dSHong Zhang PetscBT lnkbt; 15852ed133dbSHong Zhang IS iperm; 1586f76d2b81SHong Zhang 1587f76d2b81SHong Zhang PetscFunctionBegin; 1588d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 15892ed133dbSHong Zhang /* check whether perm is the identity mapping */ 1590f76d2b81SHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 15912ed133dbSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 15922ed133dbSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 15939bfd6278SHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 159410c27e3dSHong Zhang 159510c27e3dSHong Zhang /* initialization */ 15961393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 159710c27e3dSHong Zhang ui[0] = 0; 159810c27e3dSHong Zhang 159910c27e3dSHong Zhang /* jl: linked list for storing indices of the pivot rows 16001393bc97SHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 16011393bc97SHong Zhang ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 16021393bc97SHong Zhang il = jl + am; 16031393bc97SHong Zhang cols = il + am; 16041393bc97SHong Zhang ui_ptr = (PetscInt**)(cols + am); 16051393bc97SHong Zhang for (i=0; i<am; i++){ 16061393bc97SHong Zhang jl[i] = am; il[i] = 0; 1607f76d2b81SHong Zhang } 160810c27e3dSHong Zhang 160910c27e3dSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 16101393bc97SHong Zhang nlnk = am + 1; 16111393bc97SHong Zhang ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 161210c27e3dSHong Zhang 16131393bc97SHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 1614a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 161510c27e3dSHong Zhang current_space = free_space; 161610c27e3dSHong Zhang 16171393bc97SHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 161810c27e3dSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 161910c27e3dSHong Zhang nzk = 0; 16202ed133dbSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 16213b4a8b6dSBarry Smith if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 16222ed133dbSHong Zhang ncols_upper = 0; 16232ed133dbSHong Zhang for (j=0; j<ncols; j++){ 16249bfd6278SHong Zhang i = riip[*(aj + ai[rip[k]] + j)]; 16252ed133dbSHong Zhang if (i >= k){ /* only take upper triangular entry */ 16262ed133dbSHong Zhang cols[ncols_upper] = i; 16272ed133dbSHong Zhang ncols_upper++; 16282ed133dbSHong Zhang } 16292ed133dbSHong Zhang } 16301393bc97SHong Zhang ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 163110c27e3dSHong Zhang nzk += nlnk; 163210c27e3dSHong Zhang 163310c27e3dSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 163410c27e3dSHong Zhang prow = jl[k]; /* 1st pivot row */ 163510c27e3dSHong Zhang 163610c27e3dSHong Zhang while (prow < k){ 163710c27e3dSHong Zhang nextprow = jl[prow]; 163810c27e3dSHong Zhang /* merge prow into k-th row */ 16391393bc97SHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 164010c27e3dSHong Zhang jmax = ui[prow+1]; 164110c27e3dSHong Zhang ncols = jmax-jmin; 16421393bc97SHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 16435a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 164410c27e3dSHong Zhang nzk += nlnk; 164510c27e3dSHong Zhang 164610c27e3dSHong Zhang /* update il and jl for prow */ 164710c27e3dSHong Zhang if (jmin < jmax){ 164810c27e3dSHong Zhang il[prow] = jmin; 16492ed133dbSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 165010c27e3dSHong Zhang } 165110c27e3dSHong Zhang prow = nextprow; 165210c27e3dSHong Zhang } 165310c27e3dSHong Zhang 165410c27e3dSHong Zhang /* if free space is not available, make more free space */ 165510c27e3dSHong Zhang if (current_space->local_remaining<nzk) { 16561393bc97SHong Zhang i = am - k + 1; /* num of unfactored rows */ 165710c27e3dSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1658a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 165910c27e3dSHong Zhang reallocs++; 166010c27e3dSHong Zhang } 166110c27e3dSHong Zhang 166210c27e3dSHong Zhang /* copy data into free space, then initialize lnk */ 16631393bc97SHong Zhang ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 166410c27e3dSHong Zhang 166510c27e3dSHong Zhang /* add the k-th row into il and jl */ 166610c27e3dSHong Zhang if (nzk-1 > 0){ 16671393bc97SHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 166810c27e3dSHong Zhang jl[k] = jl[i]; jl[i] = k; 166910c27e3dSHong Zhang il[k] = ui[k] + 1; 167010c27e3dSHong Zhang } 16712ed133dbSHong Zhang ui_ptr[k] = current_space->array; 167210c27e3dSHong Zhang current_space->array += nzk; 167310c27e3dSHong Zhang current_space->local_used += nzk; 167410c27e3dSHong Zhang current_space->local_remaining -= nzk; 167510c27e3dSHong Zhang 167610c27e3dSHong Zhang ui[k+1] = ui[k] + nzk; 167710c27e3dSHong Zhang } 167810c27e3dSHong Zhang 16796cf91177SBarry Smith #if defined(PETSC_USE_INFO) 16801393bc97SHong Zhang if (ai[am] != 0) { 16811393bc97SHong Zhang PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 1682ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1683ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1684ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 168510c27e3dSHong Zhang } else { 1686ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 168710c27e3dSHong Zhang } 168863ba0a88SBarry Smith #endif 168910c27e3dSHong Zhang 169010c27e3dSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 16919bfd6278SHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 169210c27e3dSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 169310c27e3dSHong Zhang 169410c27e3dSHong Zhang /* destroy list of free space and other temporary array(s) */ 16951393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1696a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 169710c27e3dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 169810c27e3dSHong Zhang 169910c27e3dSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 170010c27e3dSHong Zhang 1701719d5645SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 170210c27e3dSHong Zhang b->singlemalloc = PETSC_FALSE; 1703e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1704e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 17051393bc97SHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 170610c27e3dSHong Zhang b->j = uj; 170710c27e3dSHong Zhang b->i = ui; 170810c27e3dSHong Zhang b->diag = 0; 170910c27e3dSHong Zhang b->ilen = 0; 171010c27e3dSHong Zhang b->imax = 0; 171110c27e3dSHong Zhang b->row = perm; 17129bfd6278SHong Zhang b->col = perm; 17139bfd6278SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 17149bfd6278SHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 17159bfd6278SHong Zhang b->icol = iperm; 171610c27e3dSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 17171393bc97SHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1718719d5645SBarry Smith ierr = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 17191393bc97SHong Zhang b->maxnz = b->nz = ui[am]; 172010c27e3dSHong Zhang 1721719d5645SBarry Smith (fact)->info.factor_mallocs = reallocs; 1722719d5645SBarry Smith (fact)->info.fill_ratio_given = fill; 17231393bc97SHong Zhang if (ai[am] != 0) { 1724719d5645SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 172510c27e3dSHong Zhang } else { 1726719d5645SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 172710c27e3dSHong Zhang } 1728719d5645SBarry Smith (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1729f76d2b81SHong Zhang PetscFunctionReturn(0); 1730f76d2b81SHong Zhang } 1731599ef60dSHong Zhang 1732599ef60dSHong Zhang #undef __FUNCT__ 17331d098868SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_iludt" 17341d098868SHong Zhang PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat A,Vec bb,Vec xx) 1735599ef60dSHong Zhang { 17361d098868SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 17371d098868SHong Zhang PetscErrorCode ierr; 17381d098868SHong Zhang PetscInt n = A->rmap->n; 17391d098868SHong Zhang const PetscInt *ai = a->i,*aj = a->j,*adiag = a->diag,*vi; 17401d098868SHong Zhang PetscScalar *x,sum; 17411d098868SHong Zhang const PetscScalar *b; 17421d098868SHong Zhang const MatScalar *aa = a->a,*v; 17431d098868SHong Zhang PetscInt i,nz; 17441d098868SHong Zhang 17451d098868SHong Zhang PetscFunctionBegin; 17461d098868SHong Zhang if (!n) PetscFunctionReturn(0); 17471d098868SHong Zhang 17481d098868SHong Zhang ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 17491d098868SHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 17501d098868SHong Zhang 17511d098868SHong Zhang /* forward solve the lower triangular */ 17521d098868SHong Zhang x[0] = b[0]; 17531d098868SHong Zhang v = aa; 17541d098868SHong Zhang vi = aj; 17551d098868SHong Zhang for (i=1; i<n; i++) { 17561d098868SHong Zhang nz = ai[i+1] - ai[i]; 17571d098868SHong Zhang sum = b[i]; 1758*e09e08ccSBarry Smith SPARSEDENSEMDOT(sum,x,v,vi,nz); 1759*e09e08ccSBarry Smith /* while (nz--) sum -= *v++ * x[*vi++];*/ 17601d098868SHong Zhang x[i] = sum; 17611d098868SHong Zhang } 17621d098868SHong Zhang 17631d098868SHong Zhang /* backward solve the upper triangular */ 17641d098868SHong Zhang v = aa + adiag[n] + 1; 17651d098868SHong Zhang vi = aj + adiag[n] + 1; 17661d098868SHong Zhang for (i=n-1; i>=0; i--){ 17671d098868SHong Zhang nz = adiag[i] - adiag[i+1] - 1; 17681d098868SHong Zhang sum = x[i]; 1769*e09e08ccSBarry Smith SPARSEDENSEMDOT(sum,x,v,vi,nz); 1770*e09e08ccSBarry Smith /* while (nz--) sum -= *v++ * x[*vi++]; */ 17711d098868SHong Zhang x[i] = sum*aa[adiag[i]]; 17721d098868SHong Zhang v++; vi++; 17731d098868SHong Zhang } 17741d098868SHong Zhang 17751d098868SHong Zhang ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 17761d098868SHong Zhang ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 17771d098868SHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 17781d098868SHong Zhang PetscFunctionReturn(0); 17791d098868SHong Zhang } 17801d098868SHong Zhang 17811d098868SHong Zhang #undef __FUNCT__ 17821d098868SHong Zhang #define __FUNCT__ "MatSolve_SeqAIJ_iludt" 17831d098868SHong Zhang PetscErrorCode MatSolve_SeqAIJ_iludt(Mat A,Vec bb,Vec xx) 17841d098868SHong Zhang { 17851d098868SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 17861d098868SHong Zhang IS iscol = a->col,isrow = a->row; 17871d098868SHong Zhang PetscErrorCode ierr; 17881d098868SHong Zhang PetscInt i,n=A->rmap->n,*vi,*ai = a->i,*aj = a->j,*adiag=a->diag; 17891d098868SHong Zhang PetscInt nz; 17901d098868SHong Zhang const PetscInt *rout,*cout,*r,*c; 17911d098868SHong Zhang PetscScalar *x,*tmp,*tmps; 17921d098868SHong Zhang const PetscScalar *b; 17931d098868SHong Zhang const MatScalar *aa = a->a,*v; 17941d098868SHong Zhang 17951d098868SHong Zhang PetscFunctionBegin; 17961d098868SHong Zhang if (!n) PetscFunctionReturn(0); 17971d098868SHong Zhang 17981d098868SHong Zhang ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 17991d098868SHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 18001d098868SHong Zhang tmp = a->solve_work; 18011d098868SHong Zhang 18021d098868SHong Zhang ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 18031d098868SHong Zhang ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 18041d098868SHong Zhang 18051d098868SHong Zhang /* forward solve the lower triangular */ 18061d098868SHong Zhang tmp[0] = b[*r++]; 18071d098868SHong Zhang tmps = tmp; 18081d098868SHong Zhang v = aa; 18091d098868SHong Zhang vi = aj; 18101d098868SHong Zhang for (i=1; i<n; i++) { 18111d098868SHong Zhang nz = ai[i+1] - ai[i]; 18121d098868SHong Zhang tmp[i] = b[*r++]; 18131d098868SHong Zhang SPARSEDENSEMDOT(tmp[i],tmps,v,vi,nz); 18141d098868SHong Zhang v += nz; vi += nz; 18151d098868SHong Zhang } 18161d098868SHong Zhang 18171d098868SHong Zhang /* backward solve the upper triangular */ 18181d098868SHong Zhang v = aa + adiag[n] + 1; 18191d098868SHong Zhang vi = aj + adiag[n] + 1; 18201d098868SHong Zhang for (i=n-1; i>=0; i--){ 18211d098868SHong Zhang nz = adiag[i] - adiag[i+1] - 1; 18221d098868SHong Zhang SPARSEDENSEMDOT(tmp[i],tmps,v,vi,nz); 18231d098868SHong Zhang x[*c--] = tmp[i] = tmp[i]*aa[adiag[i]]; 18241d098868SHong Zhang v += nz+1; vi += nz+1; 18251d098868SHong Zhang } 18261d098868SHong Zhang 18271d098868SHong Zhang ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 18281d098868SHong Zhang ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 18291d098868SHong Zhang ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 18301d098868SHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 18311d098868SHong Zhang ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr); 18321d098868SHong Zhang PetscFunctionReturn(0); 18331d098868SHong Zhang } 18341d098868SHong Zhang 18351d098868SHong Zhang #undef __FUNCT__ 18361d098868SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 18371d098868SHong Zhang PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 18381d098868SHong Zhang { 18391d098868SHong Zhang Mat B = *fact; 1840599ef60dSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b; 1841599ef60dSHong Zhang IS isicol; 1842599ef60dSHong Zhang PetscErrorCode ierr; 18431d098868SHong Zhang const PetscInt *r,*ic; 18441d098868SHong Zhang PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 18451d098868SHong Zhang PetscInt *bi,*bj,*bdiag; 18461d098868SHong Zhang PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 18471d098868SHong Zhang PetscInt nlnk,*lnk; 18481d098868SHong Zhang PetscBT lnkbt; 18491d098868SHong Zhang PetscTruth row_identity,icol_identity,both_identity; 18501d098868SHong Zhang MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp; 18511d098868SHong Zhang const PetscInt *ics; 18521d098868SHong Zhang PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 18531d098868SHong Zhang PetscReal dt=info->dt,shift=info->shiftinblocks; 18541d098868SHong Zhang PetscInt nnz_max; 1855599ef60dSHong Zhang PetscTruth missing; 1856599ef60dSHong Zhang 1857599ef60dSHong Zhang PetscFunctionBegin; 18581d098868SHong Zhang /* ------- symbolic factorization, can be reused ---------*/ 18591d098868SHong Zhang ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 18601d098868SHong Zhang if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 18611d098868SHong Zhang adiag=a->diag; 1862599ef60dSHong Zhang 1863599ef60dSHong Zhang ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1864599ef60dSHong Zhang 1865599ef60dSHong Zhang /* bdiag is location of diagonal in factor */ 1866599ef60dSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1867599ef60dSHong Zhang 18681d098868SHong Zhang /* allocate row pointers bi */ 18691d098868SHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 18701d098868SHong Zhang 18711d098868SHong Zhang /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 1872599ef60dSHong Zhang dtcount = (PetscInt)info->dtcount; 18733c2a7987SHong Zhang if (dtcount > n/2) dtcount = n/2; 18741d098868SHong Zhang nnz_max = ai[n]+2*n*dtcount+2; 18751d098868SHong Zhang ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 18761d098868SHong Zhang ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1877599ef60dSHong Zhang 1878599ef60dSHong Zhang /* put together the new matrix */ 1879599ef60dSHong Zhang ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1880599ef60dSHong Zhang ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 1881599ef60dSHong Zhang b = (Mat_SeqAIJ*)(B)->data; 1882599ef60dSHong Zhang b->free_a = PETSC_TRUE; 1883599ef60dSHong Zhang b->free_ij = PETSC_TRUE; 1884599ef60dSHong Zhang b->singlemalloc = PETSC_FALSE; 1885599ef60dSHong Zhang b->a = ba; 1886599ef60dSHong Zhang b->j = bj; 1887599ef60dSHong Zhang b->i = bi; 1888599ef60dSHong Zhang b->diag = bdiag; 1889599ef60dSHong Zhang b->ilen = 0; 1890599ef60dSHong Zhang b->imax = 0; 1891599ef60dSHong Zhang b->row = isrow; 1892599ef60dSHong Zhang b->col = iscol; 1893599ef60dSHong Zhang ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1894599ef60dSHong Zhang ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1895599ef60dSHong Zhang b->icol = isicol; 1896599ef60dSHong Zhang ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1897599ef60dSHong Zhang 18981d098868SHong Zhang ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 18991d098868SHong Zhang b->maxnz = nnz_max; 1900599ef60dSHong Zhang 1901599ef60dSHong Zhang (B)->factor = MAT_FACTOR_ILUDT; 1902599ef60dSHong Zhang (B)->info.factor_mallocs = 0; 19031d098868SHong Zhang (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]); 19041d098868SHong Zhang CHKMEMQ; 19051d098868SHong Zhang /* ------- end of symbolic factorization ---------*/ 1906599ef60dSHong Zhang 1907599ef60dSHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1908599ef60dSHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1909599ef60dSHong Zhang ics = ic; 1910599ef60dSHong Zhang 1911599ef60dSHong Zhang /* linked list for storing column indices of the active row */ 1912599ef60dSHong Zhang nlnk = n + 1; 1913599ef60dSHong Zhang ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 19141d098868SHong Zhang 19151d098868SHong Zhang /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 19161d098868SHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr); 19171d098868SHong Zhang jtmp = im + n; 19181d098868SHong Zhang /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 19198369b28dSHong Zhang ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 19208369b28dSHong Zhang vtmp = rtmp + n; 1921599ef60dSHong Zhang 1922599ef60dSHong Zhang bi[0] = 0; 19231d098868SHong Zhang bdiag[0] = nnz_max-1; /* location of diagonal in factor B */ 1924599ef60dSHong Zhang for (i=0; i<n; i++) { 1925599ef60dSHong Zhang /* copy initial fill into linked list */ 1926599ef60dSHong Zhang nzi = 0; /* nonzeros for active row i */ 19278369b28dSHong Zhang nzi = ai[r[i]+1] - ai[r[i]]; 19288369b28dSHong Zhang if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 19298369b28dSHong Zhang nzi_al = adiag[r[i]] - ai[r[i]]; 19308369b28dSHong Zhang nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 1931599ef60dSHong Zhang ajtmp = aj + ai[r[i]]; 19328369b28dSHong Zhang ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1933599ef60dSHong Zhang 1934599ef60dSHong Zhang /* load in initial (unfactored row) */ 1935599ef60dSHong Zhang ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 19361d098868SHong Zhang aatmp = a->a + ai[r[i]]; 19378369b28dSHong Zhang for (j=0; j<nzi; j++) { 19381d098868SHong Zhang rtmp[ics[*ajtmp++]] = *aatmp++; 1939599ef60dSHong Zhang } 1940599ef60dSHong Zhang 1941599ef60dSHong Zhang /* add pivot rows into linked list */ 1942599ef60dSHong Zhang row = lnk[n]; 1943599ef60dSHong Zhang while (row < i) { 19441d098868SHong Zhang nzi_bl = bi[row+1] - bi[row] + 1; 19451d098868SHong Zhang bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 19461d098868SHong Zhang ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 1947599ef60dSHong Zhang nzi += nlnk; 1948599ef60dSHong Zhang row = lnk[row]; 1949599ef60dSHong Zhang } 1950599ef60dSHong Zhang 19518369b28dSHong Zhang /* copy data from lnk into jtmp, then initialize lnk */ 19528369b28dSHong Zhang ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 1953599ef60dSHong Zhang 1954599ef60dSHong Zhang /* numerical factorization */ 19558369b28dSHong Zhang bjtmp = jtmp; 1956599ef60dSHong Zhang row = *bjtmp++; /* 1st pivot row */ 1957599ef60dSHong Zhang while (row < i) { 1958599ef60dSHong Zhang pc = rtmp + row; 19593c2a7987SHong Zhang pv = ba + bdiag[row]; /* 1./(diag of the pivot row) */ 19603c2a7987SHong Zhang multiplier = (*pc) * (*pv); 1961599ef60dSHong Zhang *pc = multiplier; 19623c2a7987SHong Zhang if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */ 19631d098868SHong Zhang pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 19641d098868SHong Zhang pv = ba + bdiag[row+1] + 1; 19651d098868SHong Zhang /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */ 19661d098868SHong Zhang nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 19671d098868SHong Zhang for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 1968dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 1969599ef60dSHong Zhang } 1970599ef60dSHong Zhang row = *bjtmp++; 1971599ef60dSHong Zhang } 1972599ef60dSHong Zhang 19738369b28dSHong Zhang /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 19748369b28dSHong Zhang nzi_bl = 0; j = 0; 19758369b28dSHong Zhang while (jtmp[j] < i){ /* Note: jtmp is sorted */ 19768369b28dSHong Zhang vtmp[j] = rtmp[jtmp[j]]; 19778369b28dSHong Zhang nzi_bl++; j++; 19788369b28dSHong Zhang } 19798369b28dSHong Zhang nzi_bu = nzi - nzi_bl -1; 19808369b28dSHong Zhang while (j < nzi){ 19818369b28dSHong Zhang vtmp[j] = rtmp[jtmp[j]]; 19828369b28dSHong Zhang j++; 19838369b28dSHong Zhang } 19848369b28dSHong Zhang 1985599ef60dSHong Zhang bjtmp = bj + bi[i]; 1986599ef60dSHong Zhang batmp = ba + bi[i]; 19878369b28dSHong Zhang /* apply level dropping rule to L part */ 19888369b28dSHong Zhang ncut = nzi_al + dtcount; 19898369b28dSHong Zhang if (ncut < nzi_bl){ 19908369b28dSHong Zhang ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr); 19918369b28dSHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 1992599ef60dSHong Zhang } else { 19938369b28dSHong Zhang ncut = nzi_bl; 19948369b28dSHong Zhang } 19958369b28dSHong Zhang for (j=0; j<ncut; j++){ 19968369b28dSHong Zhang bjtmp[j] = jtmp[j]; 19978369b28dSHong Zhang batmp[j] = vtmp[j]; 19986da515a1SHong Zhang /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */ 19998369b28dSHong Zhang } 20001d098868SHong Zhang bi[i+1] = bi[i] + ncut; 20018369b28dSHong Zhang nzi = ncut + 1; 20028369b28dSHong Zhang 20038369b28dSHong Zhang /* apply level dropping rule to U part */ 20048369b28dSHong Zhang ncut = nzi_au + dtcount; 20058369b28dSHong Zhang if (ncut < nzi_bu){ 20068369b28dSHong Zhang ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 20078369b28dSHong Zhang ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 20088369b28dSHong Zhang } else { 20098369b28dSHong Zhang ncut = nzi_bu; 20108369b28dSHong Zhang } 20118369b28dSHong Zhang nzi += ncut; 20121d098868SHong Zhang 20131d098868SHong Zhang /* mark bdiagonal */ 20141d098868SHong Zhang bdiag[i+1] = bdiag[i] - (ncut + 1); 20151d098868SHong Zhang bjtmp = bj + bdiag[i]; 20161d098868SHong Zhang batmp = ba + bdiag[i]; 20171d098868SHong Zhang *bjtmp = i; 20181d098868SHong Zhang *batmp = rtmp[i]; 2019b78a477dSHong Zhang if (*batmp == 0.0) *batmp = dt+shift; 2020b78a477dSHong Zhang *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */ 20216da515a1SHong Zhang /* printf(" (%d,%g),",*bjtmp,*batmp); */ 20221d098868SHong Zhang 20231d098868SHong Zhang bjtmp = bj + bdiag[i+1]+1; 20241d098868SHong Zhang batmp = ba + bdiag[i+1]+1; 20258369b28dSHong Zhang for (k=0; k<ncut; k++){ 20261d098868SHong Zhang bjtmp[k] = jtmp[nzi_bl+1+k]; 20271d098868SHong Zhang batmp[k] = vtmp[nzi_bl+1+k]; 20286da515a1SHong Zhang /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */ 2029599ef60dSHong Zhang } 20306da515a1SHong Zhang /* printf("\n"); */ 2031599ef60dSHong Zhang 2032599ef60dSHong Zhang im[i] = nzi; /* used by PetscLLAddSortedLU() */ 2033599ef60dSHong Zhang /* 20341d098868SHong Zhang printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]); 20358369b28dSHong Zhang printf(" ----------------------------\n"); 2036599ef60dSHong Zhang */ 20378369b28dSHong Zhang } /* for (i=0; i<n; i++) */ 20381d098868SHong Zhang /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */ 20391d098868SHong Zhang if (bi[n] >= bdiag[n]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[n],bdiag[n]); 2040599ef60dSHong Zhang 2041599ef60dSHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2042599ef60dSHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2043599ef60dSHong Zhang 2044599ef60dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2045599ef60dSHong Zhang ierr = PetscFree(im);CHKERRQ(ierr); 2046599ef60dSHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 2047599ef60dSHong Zhang 2048599ef60dSHong Zhang ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr); 20491d098868SHong Zhang b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n]; 2050599ef60dSHong Zhang 2051599ef60dSHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 2052599ef60dSHong Zhang ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 2053599ef60dSHong Zhang both_identity = (PetscTruth) (row_identity && icol_identity); 2054599ef60dSHong Zhang if (row_identity && icol_identity) { 20551d098868SHong Zhang B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt; 2056599ef60dSHong Zhang } else { 20571d098868SHong Zhang B->ops->solve = MatSolve_SeqAIJ_iludt; 2058599ef60dSHong Zhang } 2059599ef60dSHong Zhang 2060b78a477dSHong Zhang B->ops->lufactorsymbolic = MatILUDTFactorSymbolic_SeqAIJ; 2061b78a477dSHong Zhang B->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 20621d098868SHong Zhang B->ops->solveadd = 0; 20631d098868SHong Zhang B->ops->solvetranspose = 0; 20641d098868SHong Zhang B->ops->solvetransposeadd = 0; 20651d098868SHong Zhang B->ops->matsolve = 0; 2066599ef60dSHong Zhang B->assembled = PETSC_TRUE; 2067599ef60dSHong Zhang B->preallocated = PETSC_TRUE; 2068599ef60dSHong Zhang PetscFunctionReturn(0); 2069599ef60dSHong Zhang } 2070599ef60dSHong Zhang 20713c2a7987SHong Zhang /* a wraper of MatILUDTFactor_SeqAIJ() */ 20723c2a7987SHong Zhang #undef __FUNCT__ 20733c2a7987SHong Zhang #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ" 20743c2a7987SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 20753c2a7987SHong Zhang { 20763c2a7987SHong Zhang PetscErrorCode ierr; 20773c2a7987SHong Zhang 20783c2a7987SHong Zhang PetscFunctionBegin; 20796da515a1SHong Zhang /* printf("MatILUDTFactorSymbolic_SeqAIJ is called...\n"); */ 20803c2a7987SHong Zhang ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr); 20813c2a7987SHong Zhang 20823c2a7987SHong Zhang fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 20833c2a7987SHong Zhang PetscFunctionReturn(0); 20843c2a7987SHong Zhang } 20853c2a7987SHong Zhang 20863c2a7987SHong Zhang /* 20873c2a7987SHong Zhang same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors 20883c2a7987SHong Zhang - intend to replace existing MatLUFactorNumeric_SeqAIJ() 20893c2a7987SHong Zhang */ 20903c2a7987SHong Zhang #undef __FUNCT__ 20913c2a7987SHong Zhang #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ" 20923c2a7987SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info) 20933c2a7987SHong Zhang { 20943c2a7987SHong Zhang Mat C=fact; 20953c2a7987SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 20963c2a7987SHong Zhang IS isrow = b->row,isicol = b->icol; 20973c2a7987SHong Zhang PetscErrorCode ierr; 20983c2a7987SHong Zhang const PetscInt *r,*ic,*ics; 2099b78a477dSHong Zhang PetscInt i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 2100b78a477dSHong Zhang PetscInt *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj; 21013c2a7987SHong Zhang MatScalar *rtmp,*pc,multiplier,*v,*pv,*aa=a->a; 2102b78a477dSHong Zhang PetscReal dt=info->dt,shift=info->shiftinblocks; 2103b78a477dSHong Zhang PetscTruth row_identity, col_identity; 21043c2a7987SHong Zhang 21053c2a7987SHong Zhang PetscFunctionBegin; 21066da515a1SHong Zhang /* printf("MatILUDTFactorNumeric_SeqAIJ is called...\n"); */ 21073c2a7987SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 21083c2a7987SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 2109b78a477dSHong Zhang ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 21103c2a7987SHong Zhang ics = ic; 21113c2a7987SHong Zhang 21123c2a7987SHong Zhang for (i=0; i<n; i++){ 2113b78a477dSHong Zhang /* initialize rtmp array */ 2114b78a477dSHong Zhang nzl = bi[i+1] - bi[i]; /* num of nozeros in L(i,:) */ 21153c2a7987SHong Zhang bjtmp = bj + bi[i]; 2116b78a477dSHong Zhang for (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0; 2117b78a477dSHong Zhang rtmp[i] = 0.0; 2118b78a477dSHong Zhang nzu = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */ 2119b78a477dSHong Zhang bjtmp = bj + bdiag[i+1] + 1; 2120b78a477dSHong Zhang for (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0; 21213c2a7987SHong Zhang 2122b78a477dSHong Zhang /* load in initial unfactored row of A */ 21236da515a1SHong Zhang /* printf("row %d\n",i); */ 21243c2a7987SHong Zhang nz = ai[r[i]+1] - ai[r[i]]; 21253c2a7987SHong Zhang ajtmp = aj + ai[r[i]]; 21263c2a7987SHong Zhang v = aa + ai[r[i]]; 21273c2a7987SHong Zhang for (j=0; j<nz; j++) { 2128b78a477dSHong Zhang rtmp[ics[*ajtmp++]] = v[j]; 21296da515a1SHong Zhang /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */ 21303c2a7987SHong Zhang } 21316da515a1SHong Zhang /* printf("\n"); */ 21323c2a7987SHong Zhang 2133b78a477dSHong Zhang /* numerical factorization */ 2134b78a477dSHong Zhang bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */ 2135b78a477dSHong Zhang nzl = bi[i+1] - bi[i]; /* num of entries in L(i,:) */ 2136b78a477dSHong Zhang k = 0; 2137b78a477dSHong Zhang while (k < nzl){ 21383c2a7987SHong Zhang row = *bjtmp++; 21396da515a1SHong Zhang /* printf(" prow %d\n",row); */ 21403c2a7987SHong Zhang pc = rtmp + row; 2141b78a477dSHong Zhang pv = b->a + bdiag[row]; /* 1./(diag of the pivot row) */ 2142b78a477dSHong Zhang multiplier = (*pc) * (*pv); 21433c2a7987SHong Zhang *pc = multiplier; 2144b78a477dSHong Zhang if (PetscAbsScalar(multiplier) > dt){ 2145b78a477dSHong Zhang pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 2146b78a477dSHong Zhang pv = b->a + bdiag[row+1] + 1; 2147b78a477dSHong Zhang nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 2148b78a477dSHong Zhang for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 21496da515a1SHong Zhang /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */ 21503c2a7987SHong Zhang } 2151b78a477dSHong Zhang k++; 21523c2a7987SHong Zhang } 21533c2a7987SHong Zhang 2154b78a477dSHong Zhang /* finished row so stick it into b->a */ 2155b78a477dSHong Zhang /* L-part */ 2156b78a477dSHong Zhang pv = b->a + bi[i] ; 2157b78a477dSHong Zhang pj = bj + bi[i] ; 2158b78a477dSHong Zhang nzl = bi[i+1] - bi[i]; 2159b78a477dSHong Zhang for (j=0; j<nzl; j++) { 2160b78a477dSHong Zhang pv[j] = rtmp[pj[j]]; 21616da515a1SHong Zhang /* printf(" (%d,%g),",pj[j],pv[j]); */ 21623c2a7987SHong Zhang } 2163b78a477dSHong Zhang 2164b78a477dSHong Zhang /* diagonal: invert diagonal entries for simplier triangular solves */ 2165b78a477dSHong Zhang if (rtmp[i] == 0.0) rtmp[i] = dt+shift; 2166b78a477dSHong Zhang b->a[bdiag[i]] = 1.0/rtmp[i]; 21676da515a1SHong Zhang /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */ 2168b78a477dSHong Zhang 2169b78a477dSHong Zhang /* U-part */ 2170b78a477dSHong Zhang pv = b->a + bdiag[i+1] + 1; 2171b78a477dSHong Zhang pj = bj + bdiag[i+1] + 1; 2172b78a477dSHong Zhang nzu = bdiag[i] - bdiag[i+1] - 1; 2173b78a477dSHong Zhang for (j=0; j<nzu; j++) { 2174b78a477dSHong Zhang pv[j] = rtmp[pj[j]]; 21756da515a1SHong Zhang /* printf(" (%d,%g),",pj[j],pv[j]); */ 2176b78a477dSHong Zhang } 21776da515a1SHong Zhang /* printf("\n"); */ 2178b78a477dSHong Zhang } 2179b78a477dSHong Zhang 21803c2a7987SHong Zhang ierr = PetscFree(rtmp);CHKERRQ(ierr); 21813c2a7987SHong Zhang ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 21823c2a7987SHong Zhang ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2183b78a477dSHong Zhang 21843c2a7987SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 21853c2a7987SHong Zhang ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 21863c2a7987SHong Zhang if (row_identity && col_identity) { 2187b78a477dSHong Zhang C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt; 21883c2a7987SHong Zhang } else { 2189b78a477dSHong Zhang C->ops->solve = MatSolve_SeqAIJ_iludt; 21903c2a7987SHong Zhang } 2191b78a477dSHong Zhang C->ops->solveadd = 0; 2192b78a477dSHong Zhang C->ops->solvetranspose = 0; 2193b78a477dSHong Zhang C->ops->solvetransposeadd = 0; 2194b78a477dSHong Zhang C->ops->matsolve = 0; 21953c2a7987SHong Zhang C->assembled = PETSC_TRUE; 21963c2a7987SHong Zhang C->preallocated = PETSC_TRUE; 21973c2a7987SHong Zhang ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 21983c2a7987SHong Zhang PetscFunctionReturn(0); 21993c2a7987SHong Zhang } 2200