1b2f1ab58SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h> 5b2f1ab58SBarry Smith 6b2f1ab58SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 7b2f1ab58SBarry Smith { 8b2f1ab58SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 9b2f1ab58SBarry Smith Mat_SeqSBAIJ *b; 10b2f1ab58SBarry Smith PetscErrorCode ierr; 11ace3abfcSBarry Smith PetscBool perm_identity,missing; 12b2f1ab58SBarry Smith PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui; 13b2f1ab58SBarry Smith const PetscInt *rip,*riip; 14b2f1ab58SBarry Smith PetscInt j; 15b2f1ab58SBarry Smith PetscInt d; 16b2f1ab58SBarry Smith PetscInt ncols,*cols,*uj; 17b2f1ab58SBarry Smith PetscReal fill=info->fill,levels=info->levels; 18b2f1ab58SBarry Smith IS iperm; 19b2f1ab58SBarry Smith spbas_matrix Pattern_0, Pattern_P; 20b2f1ab58SBarry Smith 21b2f1ab58SBarry Smith PetscFunctionBegin; 22e32f2f54SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 23b2f1ab58SBarry Smith ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 24e32f2f54SBarry Smith if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 25b2f1ab58SBarry Smith ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 26b2f1ab58SBarry Smith ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 27b2f1ab58SBarry Smith 28b2f1ab58SBarry Smith /* ICC(0) without matrix ordering: simply copies fill pattern */ 29b2f1ab58SBarry Smith if (!levels && perm_identity) { 30854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&ui);CHKERRQ(ierr); 31b2f1ab58SBarry Smith ui[0] = 0; 32b2f1ab58SBarry Smith 33b2f1ab58SBarry Smith for (i=0; i<am; i++) { 34b2f1ab58SBarry Smith ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 35b2f1ab58SBarry Smith } 36854ce69bSBarry Smith ierr = PetscMalloc1(ui[am]+1,&uj);CHKERRQ(ierr); 37b2f1ab58SBarry Smith cols = uj; 38b2f1ab58SBarry Smith for (i=0; i<am; i++) { 39b2f1ab58SBarry Smith aj = a->j + a->diag[i]; 40b2f1ab58SBarry Smith ncols = ui[i+1] - ui[i]; 41b2f1ab58SBarry Smith for (j=0; j<ncols; j++) *cols++ = *aj++; 42b2f1ab58SBarry Smith } 43b2f1ab58SBarry Smith } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 44b2f1ab58SBarry Smith ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 45b2f1ab58SBarry Smith ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 46b2f1ab58SBarry Smith 4771d55d18SBarry Smith /* Create spbas_matrix for pattern */ 48b2f1ab58SBarry Smith ierr = spbas_pattern_only(am, am, ai, aj, &Pattern_0);CHKERRQ(ierr); 49b2f1ab58SBarry Smith 5071d55d18SBarry Smith /* Apply the permutation */ 51b2f1ab58SBarry Smith ierr = spbas_apply_reordering(&Pattern_0, rip, riip);CHKERRQ(ierr); 52b2f1ab58SBarry Smith 5371d55d18SBarry Smith /* Raise the power */ 5471d55d18SBarry Smith ierr = spbas_power(Pattern_0, (int) levels+1, &Pattern_P);CHKERRQ(ierr); 55b2f1ab58SBarry Smith ierr = spbas_delete(Pattern_0);CHKERRQ(ierr); 56b2f1ab58SBarry Smith 5771d55d18SBarry Smith /* Keep only upper triangle of pattern */ 5822d28d08SBarry Smith ierr = spbas_keep_upper(&Pattern_P);CHKERRQ(ierr); 59b2f1ab58SBarry Smith 6071d55d18SBarry Smith /* Convert to Sparse Row Storage */ 610298fd71SBarry Smith ierr = spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj);CHKERRQ(ierr); 62b2f1ab58SBarry Smith ierr = spbas_delete(Pattern_P);CHKERRQ(ierr); 63b2f1ab58SBarry Smith } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 64b2f1ab58SBarry Smith 65b2f1ab58SBarry Smith /* put together the new matrix in MATSEQSBAIJ format */ 66b2f1ab58SBarry Smith 67b2f1ab58SBarry Smith b = (Mat_SeqSBAIJ*)(fact)->data; 68b2f1ab58SBarry Smith b->singlemalloc = PETSC_FALSE; 692205254eSKarl Rupp 70854ce69bSBarry Smith ierr = PetscMalloc1(ui[am]+1,&b->a);CHKERRQ(ierr); 712205254eSKarl Rupp 72b2f1ab58SBarry Smith b->j = uj; 73b2f1ab58SBarry Smith b->i = ui; 74f4259b30SLisandro Dalcin b->diag = NULL; 75f4259b30SLisandro Dalcin b->ilen = NULL; 76f4259b30SLisandro Dalcin b->imax = NULL; 77b2f1ab58SBarry Smith b->row = perm; 78b2f1ab58SBarry Smith b->col = perm; 792205254eSKarl Rupp 80b2f1ab58SBarry Smith ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 81b2f1ab58SBarry Smith ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 822205254eSKarl Rupp 83b2f1ab58SBarry Smith b->icol = iperm; 84b2f1ab58SBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 85854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&b->solve_work);CHKERRQ(ierr); 863bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)(fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 87b2f1ab58SBarry Smith b->maxnz = b->nz = ui[am]; 88b2f1ab58SBarry Smith b->free_a = PETSC_TRUE; 89b2f1ab58SBarry Smith b->free_ij = PETSC_TRUE; 90b2f1ab58SBarry Smith 91b2f1ab58SBarry Smith (fact)->info.factor_mallocs = reallocs; 92b2f1ab58SBarry Smith (fact)->info.fill_ratio_given = fill; 93b2f1ab58SBarry Smith if (ai[am] != 0) { 94b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 95b2f1ab58SBarry Smith } else { 96b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 97b2f1ab58SBarry Smith } 98b2f1ab58SBarry Smith /* (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */ 99b2f1ab58SBarry Smith PetscFunctionReturn(0); 100b2f1ab58SBarry Smith } 101b2f1ab58SBarry Smith 102b2f1ab58SBarry Smith 103b2f1ab58SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info) 104b2f1ab58SBarry Smith { 105b2f1ab58SBarry Smith Mat C = B; 106b2f1ab58SBarry Smith Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 107b2f1ab58SBarry Smith IS ip=b->row,iip = b->icol; 108b2f1ab58SBarry Smith PetscErrorCode ierr; 109b2f1ab58SBarry Smith const PetscInt *rip,*riip; 110b2f1ab58SBarry Smith PetscInt mbs=A->rmap->n,*bi=b->i,*bj=b->j; 111b2f1ab58SBarry Smith 112b2f1ab58SBarry Smith MatScalar *ba = b->a; 113c88783f4SHong Zhang PetscReal shiftnz = info->shiftamount; 114d36aa34eSBarry Smith PetscReal droptol = -1; 115ace3abfcSBarry Smith PetscBool perm_identity; 116b2f1ab58SBarry Smith spbas_matrix Pattern, matrix_L,matrix_LT; 1179767453cSBarry Smith PetscReal mem_reduction; 118b2f1ab58SBarry Smith 119b2f1ab58SBarry Smith PetscFunctionBegin; 12071d55d18SBarry Smith /* Reduce memory requirements: erase values of B-matrix */ 121b2f1ab58SBarry Smith ierr = PetscFree(ba);CHKERRQ(ierr); 12271d55d18SBarry Smith /* Compress (maximum) sparseness pattern of B-matrix */ 1239767453cSBarry Smith ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr); 124b2f1ab58SBarry Smith ierr = PetscFree(bi);CHKERRQ(ierr); 125b2f1ab58SBarry Smith ierr = PetscFree(bj);CHKERRQ(ierr); 126b2f1ab58SBarry Smith 12757622a8eSBarry Smith ierr = PetscInfo1(NULL," compression rate for spbas_compress_pattern %g \n",(double)mem_reduction);CHKERRQ(ierr); 128b2f1ab58SBarry Smith 1294e1ff37aSBarry Smith /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */ 130b2f1ab58SBarry Smith ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 131b2f1ab58SBarry Smith ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 132b2f1ab58SBarry Smith 133b2f1ab58SBarry Smith if (info->usedt) { 134b2f1ab58SBarry Smith droptol = info->dt; 135b2f1ab58SBarry Smith } 136b2f1ab58SBarry Smith for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) 137b2f1ab58SBarry Smith { 13871d55d18SBarry Smith ierr = spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr); 1392205254eSKarl Rupp if (ierr == NEGATIVE_DIAGONAL) { 140b2f1ab58SBarry Smith shiftnz *= 1.5; 1419767453cSBarry Smith if (shiftnz < 1e-5) shiftnz=1e-5; 14257622a8eSBarry Smith ierr = PetscInfo1(NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n",(double)shiftnz);CHKERRQ(ierr); 143b2f1ab58SBarry Smith } 144b2f1ab58SBarry Smith } 145b2f1ab58SBarry Smith ierr = spbas_delete(Pattern);CHKERRQ(ierr); 146b2f1ab58SBarry Smith 1476712e2f1SBarry Smith ierr = PetscInfo1(NULL," memory_usage for spbas_incomplete_cholesky %g bytes per row\n", (double)(PetscReal) (spbas_memory_requirement(matrix_LT)/ (PetscReal) mbs));CHKERRQ(ierr); 148b2f1ab58SBarry Smith 149b2f1ab58SBarry Smith ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 150b2f1ab58SBarry Smith ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 151b2f1ab58SBarry Smith 15271d55d18SBarry Smith /* Convert spbas_matrix to compressed row storage */ 153b2f1ab58SBarry Smith ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr); 154b2f1ab58SBarry Smith ierr = spbas_delete(matrix_LT);CHKERRQ(ierr); 155b2f1ab58SBarry Smith ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr); 156b2f1ab58SBarry Smith b->i =bi; b->j=bj; b->a=ba; 157b2f1ab58SBarry Smith ierr = spbas_delete(matrix_L);CHKERRQ(ierr); 158b2f1ab58SBarry Smith 15971d55d18SBarry Smith /* Set the appropriate solution functions */ 160b2f1ab58SBarry Smith ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 161b2f1ab58SBarry Smith if (perm_identity) { 162b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 163b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 164b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 165b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 166b2f1ab58SBarry Smith } else { 167b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_inplace; 168b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; 169b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; 170b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; 171b2f1ab58SBarry Smith } 172b2f1ab58SBarry Smith 173b2f1ab58SBarry Smith C->assembled = PETSC_TRUE; 174b2f1ab58SBarry Smith C->preallocated = PETSC_TRUE; 1752205254eSKarl Rupp 176b2f1ab58SBarry Smith ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 177b2f1ab58SBarry Smith PetscFunctionReturn(0); 178b2f1ab58SBarry Smith } 179b2f1ab58SBarry Smith 180ea799195SBarry Smith PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A,MatSolverType *type) 1813dfa2556SBarry Smith { 1823dfa2556SBarry Smith PetscFunctionBegin; 1833dfa2556SBarry Smith *type = MATSOLVERBAS; 1843dfa2556SBarry Smith PetscFunctionReturn(0); 1853dfa2556SBarry Smith } 1863dfa2556SBarry Smith 187cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat *B) 188b2f1ab58SBarry Smith { 189b2f1ab58SBarry Smith PetscInt n = A->rmap->n; 190b2f1ab58SBarry Smith PetscErrorCode ierr; 191b2f1ab58SBarry Smith 192b2f1ab58SBarry Smith PetscFunctionBegin; 193ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 194b2f1ab58SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 195b2f1ab58SBarry Smith if (ftype == MAT_FACTOR_ICC) { 196b2f1ab58SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 1970298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 1982205254eSKarl Rupp 199b2f1ab58SBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas; 200b2f1ab58SBarry Smith (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas; 2013ca39a21SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)*B,"MatFactorGetSolverType_C",MatFactorGetSolverType_seqaij_bas);CHKERRQ(ierr); 202e32f2f54SBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 203d5f3da31SBarry Smith (*B)->factortype = ftype; 20400c67f3bSHong Zhang 20500c67f3bSHong Zhang ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr); 20600c67f3bSHong Zhang ierr = PetscStrallocpy(MATSOLVERBAS,&(*B)->solvertype);CHKERRQ(ierr); 207*f73b0415SBarry Smith (*B)->canuseordering = PETSC_TRUE; 208b2f1ab58SBarry Smith PetscFunctionReturn(0); 209b2f1ab58SBarry Smith } 210