1b2f1ab58SBarry Smith 29767453cSBarry Smith 3b2f1ab58SBarry Smith #include <petsc.h> 4b2f1ab58SBarry Smith #include <petscksp.h> 5b2f1ab58SBarry Smith #include "private/kspimpl.h" 6b2f1ab58SBarry Smith #include "petscpc.h" 7b2f1ab58SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 8b2f1ab58SBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 9b2f1ab58SBarry Smith #include "../src/mat/impls/aij/seq/bas/spbas.h" 10b2f1ab58SBarry Smith 11b2f1ab58SBarry Smith #undef __FUNCT__ 12b2f1ab58SBarry Smith #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ_Bas" 13b2f1ab58SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 14b2f1ab58SBarry Smith { 15b2f1ab58SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 16b2f1ab58SBarry Smith Mat_SeqSBAIJ *b; 17b2f1ab58SBarry Smith PetscErrorCode ierr; 18b2f1ab58SBarry Smith PetscTruth perm_identity,missing; 19b2f1ab58SBarry Smith PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui; 20b2f1ab58SBarry Smith const PetscInt *rip,*riip; 21b2f1ab58SBarry Smith PetscInt j; 22b2f1ab58SBarry Smith PetscInt d; 23b2f1ab58SBarry Smith PetscInt ncols,*cols,*uj; 24b2f1ab58SBarry Smith PetscReal fill=info->fill,levels=info->levels; 25b2f1ab58SBarry Smith IS iperm; 26b2f1ab58SBarry Smith spbas_matrix Pattern_0, Pattern_P; 27b2f1ab58SBarry Smith 28b2f1ab58SBarry Smith PetscFunctionBegin; 29b2f1ab58SBarry 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); 30b2f1ab58SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 31b2f1ab58SBarry Smith if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 32b2f1ab58SBarry Smith ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 33b2f1ab58SBarry Smith ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 34b2f1ab58SBarry Smith 35b2f1ab58SBarry Smith 36b2f1ab58SBarry Smith /* ICC(0) without matrix ordering: simply copies fill pattern */ 37b2f1ab58SBarry Smith if (!levels && perm_identity) { 38b2f1ab58SBarry Smith ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 39b2f1ab58SBarry Smith ui[0] = 0; 40b2f1ab58SBarry Smith 41b2f1ab58SBarry Smith for (i=0; i<am; i++) { 42b2f1ab58SBarry Smith ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 43b2f1ab58SBarry Smith } 44b2f1ab58SBarry Smith ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 45b2f1ab58SBarry Smith cols = uj; 46b2f1ab58SBarry Smith for (i=0; i<am; i++) { 47b2f1ab58SBarry Smith aj = a->j + a->diag[i]; 48b2f1ab58SBarry Smith ncols = ui[i+1] - ui[i]; 49b2f1ab58SBarry Smith for (j=0; j<ncols; j++) *cols++ = *aj++; 50b2f1ab58SBarry Smith } 51b2f1ab58SBarry Smith } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 52b2f1ab58SBarry Smith ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 53b2f1ab58SBarry Smith ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 54b2f1ab58SBarry Smith 55*71d55d18SBarry Smith /* Create spbas_matrix for pattern */ 56b2f1ab58SBarry Smith ierr = spbas_pattern_only(am, am, ai, aj, &Pattern_0);CHKERRQ(ierr); 57b2f1ab58SBarry Smith 58*71d55d18SBarry Smith /* Apply the permutation */ 59b2f1ab58SBarry Smith ierr = spbas_apply_reordering( &Pattern_0, rip, riip);CHKERRQ(ierr); 60b2f1ab58SBarry Smith 61*71d55d18SBarry Smith /* Raise the power */ 62*71d55d18SBarry Smith ierr = spbas_power( Pattern_0, (int) levels+1, &Pattern_P);CHKERRQ(ierr); 63b2f1ab58SBarry Smith ierr = spbas_delete( Pattern_0 );CHKERRQ(ierr); 64b2f1ab58SBarry Smith 65*71d55d18SBarry Smith /* Keep only upper triangle of pattern */ 66b2f1ab58SBarry Smith ierr = spbas_keep_upper( &Pattern_P ); 67b2f1ab58SBarry Smith 68*71d55d18SBarry Smith /* Convert to Sparse Row Storage */ 69c328eeadSBarry Smith ierr = spbas_matrix_to_crs(Pattern_P, PETSC_NULL, &ui, &uj);CHKERRQ(ierr); 70b2f1ab58SBarry Smith ierr = spbas_delete(Pattern_P);CHKERRQ(ierr); 71b2f1ab58SBarry Smith } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 72b2f1ab58SBarry Smith 73b2f1ab58SBarry Smith /* put together the new matrix in MATSEQSBAIJ format */ 74b2f1ab58SBarry Smith 75b2f1ab58SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 76b2f1ab58SBarry Smith b->singlemalloc = PETSC_FALSE; 77b2f1ab58SBarry Smith ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 78b2f1ab58SBarry Smith b->j = uj; 79b2f1ab58SBarry Smith b->i = ui; 80b2f1ab58SBarry Smith b->diag = 0; 81b2f1ab58SBarry Smith b->ilen = 0; 82b2f1ab58SBarry Smith b->imax = 0; 83b2f1ab58SBarry Smith b->row = perm; 84b2f1ab58SBarry Smith b->col = perm; 85b2f1ab58SBarry Smith ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 86b2f1ab58SBarry Smith ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 87b2f1ab58SBarry Smith b->icol = iperm; 88b2f1ab58SBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 89b2f1ab58SBarry Smith ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 90b2f1ab58SBarry Smith ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 91b2f1ab58SBarry Smith b->maxnz = b->nz = ui[am]; 92b2f1ab58SBarry Smith b->free_a = PETSC_TRUE; 93b2f1ab58SBarry Smith b->free_ij = PETSC_TRUE; 94b2f1ab58SBarry Smith 95b2f1ab58SBarry Smith (fact)->info.factor_mallocs = reallocs; 96b2f1ab58SBarry Smith (fact)->info.fill_ratio_given = fill; 97b2f1ab58SBarry Smith if (ai[am] != 0) { 98b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 99b2f1ab58SBarry Smith } else { 100b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 101b2f1ab58SBarry Smith } 102b2f1ab58SBarry Smith /* (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */ 103b2f1ab58SBarry Smith PetscFunctionReturn(0); 104b2f1ab58SBarry Smith } 105b2f1ab58SBarry Smith 106b2f1ab58SBarry Smith 107b2f1ab58SBarry Smith #undef __FUNCT__ 108b2f1ab58SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_Bas" 109b2f1ab58SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info) 110b2f1ab58SBarry Smith { 111b2f1ab58SBarry Smith Mat C = B; 112b2f1ab58SBarry Smith Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 113b2f1ab58SBarry Smith IS ip=b->row,iip = b->icol; 114b2f1ab58SBarry Smith PetscErrorCode ierr; 115b2f1ab58SBarry Smith const PetscInt *rip,*riip; 116b2f1ab58SBarry Smith PetscInt mbs=A->rmap->n,*bi=b->i,*bj=b->j; 117b2f1ab58SBarry Smith 118b2f1ab58SBarry Smith MatScalar *ba=b->a; 119b2f1ab58SBarry Smith PetscReal shiftnz = info->shiftnz; 120d36aa34eSBarry Smith PetscReal droptol = -1; 121b2f1ab58SBarry Smith PetscTruth perm_identity; 122b2f1ab58SBarry Smith spbas_matrix Pattern, matrix_L,matrix_LT; 1239767453cSBarry Smith PetscReal mem_reduction; 124b2f1ab58SBarry Smith 125b2f1ab58SBarry Smith PetscFunctionBegin; 126*71d55d18SBarry Smith /* Reduce memory requirements: erase values of B-matrix */ 127b2f1ab58SBarry Smith ierr = PetscFree(ba);CHKERRQ(ierr); 128*71d55d18SBarry Smith /* Compress (maximum) sparseness pattern of B-matrix */ 1299767453cSBarry Smith ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr); 130b2f1ab58SBarry Smith ierr = PetscFree(bi);CHKERRQ(ierr); 131b2f1ab58SBarry Smith ierr = PetscFree(bj);CHKERRQ(ierr); 132b2f1ab58SBarry Smith 1339767453cSBarry Smith ierr = PetscInfo1(PETSC_NULL," compression rate for spbas_compress_pattern %G \n",mem_reduction);CHKERRQ(ierr); 134b2f1ab58SBarry Smith 135*71d55d18SBarry Smith /* Make Cholesky decompositions with larger Manteuffel shifts until no more 136*71d55d18SBarry Smith negative diagonals are found. */ 137b2f1ab58SBarry Smith ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 138b2f1ab58SBarry Smith ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 139b2f1ab58SBarry Smith 140b2f1ab58SBarry Smith if (info->usedt) { 141b2f1ab58SBarry Smith droptol = info->dt; 142b2f1ab58SBarry Smith } 143b2f1ab58SBarry Smith for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL; ) 144b2f1ab58SBarry Smith { 145*71d55d18SBarry Smith ierr = spbas_incomplete_cholesky( A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr); 146b2f1ab58SBarry Smith if (ierr == NEGATIVE_DIAGONAL) 147b2f1ab58SBarry Smith { 148b2f1ab58SBarry Smith shiftnz *= 1.5; 1499767453cSBarry Smith if (shiftnz < 1e-5) shiftnz=1e-5; 1509767453cSBarry Smith ierr = PetscInfo1(PETSC_NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%G\n",shiftnz);CHKERRQ(ierr); 151b2f1ab58SBarry Smith } 152b2f1ab58SBarry Smith } 153b2f1ab58SBarry Smith ierr = spbas_delete(Pattern);CHKERRQ(ierr); 154b2f1ab58SBarry Smith 1559767453cSBarry Smith ierr = PetscInfo1(PETSC_NULL," memory_usage for spbas_incomplete_cholesky %G bytes per row\n", 1569767453cSBarry Smith (PetscReal) spbas_memory_requirement( matrix_LT)/ (PetscReal) mbs);CHKERRQ(ierr); 157b2f1ab58SBarry Smith 158b2f1ab58SBarry Smith ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 159b2f1ab58SBarry Smith ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 160b2f1ab58SBarry Smith 161*71d55d18SBarry Smith /* Convert spbas_matrix to compressed row storage */ 162b2f1ab58SBarry Smith ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr); 163b2f1ab58SBarry Smith ierr = spbas_delete(matrix_LT);CHKERRQ(ierr); 164b2f1ab58SBarry Smith ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr); 165b2f1ab58SBarry Smith b->i=bi; b->j=bj; b->a=ba; 166b2f1ab58SBarry Smith ierr = spbas_delete(matrix_L);CHKERRQ(ierr); 167b2f1ab58SBarry Smith 168*71d55d18SBarry Smith /* Set the appropriate solution functions */ 169b2f1ab58SBarry Smith ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 170b2f1ab58SBarry Smith if (perm_identity){ 171b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 172b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 173b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 174b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 175b2f1ab58SBarry Smith } else { 176b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_inplace; 177b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; 178b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; 179b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; 180b2f1ab58SBarry Smith } 181b2f1ab58SBarry Smith 182b2f1ab58SBarry Smith C->assembled = PETSC_TRUE; 183b2f1ab58SBarry Smith C->preallocated = PETSC_TRUE; 184b2f1ab58SBarry Smith ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 185b2f1ab58SBarry Smith PetscFunctionReturn(0); 186b2f1ab58SBarry Smith } 187b2f1ab58SBarry Smith 1881e364ec1SMatthew G Knepley EXTERN_C_BEGIN 189b2f1ab58SBarry Smith #undef __FUNCT__ 190b2f1ab58SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_bas" 191b2f1ab58SBarry Smith PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B) 192b2f1ab58SBarry Smith { 193b2f1ab58SBarry Smith PetscInt n = A->rmap->n; 194b2f1ab58SBarry Smith PetscErrorCode ierr; 195b2f1ab58SBarry Smith 196b2f1ab58SBarry Smith PetscFunctionBegin; 197b2f1ab58SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 198b2f1ab58SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 199b2f1ab58SBarry Smith if (ftype == MAT_FACTOR_ICC) { 200b2f1ab58SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 201b2f1ab58SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 202b2f1ab58SBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas; 203b2f1ab58SBarry Smith (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas; 204b2f1ab58SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 205b2f1ab58SBarry Smith (*B)->factor = ftype; 206b2f1ab58SBarry Smith PetscFunctionReturn(0); 207b2f1ab58SBarry Smith } 2081e364ec1SMatthew G Knepley EXTERN_C_END 209