PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo *info) { Mat C = B; Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; IS ip=b->row,iip = b->icol; PetscErrorCode ierr; const PetscInt *rip,*riip; PetscInt mbs=A->rmap->n,*bi=b->i,*bj=b->j; MatScalar *ba = b->a; PetscReal shiftnz = info->shiftamount; PetscReal droptol = -1; PetscBool perm_identity; spbas_matrix Pattern, matrix_L,matrix_LT; PetscReal mem_reduction; PetscFunctionBegin; /* Reduce memory requirements: erase values of B-matrix */ ierr = PetscFree(ba);CHKERRQ(ierr); /* Compress (maximum) sparseness pattern of B-matrix */ ierr = spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS,&Pattern, &mem_reduction);CHKERRQ(ierr); ierr = PetscFree(bi);CHKERRQ(ierr); ierr = PetscFree(bj);CHKERRQ(ierr); ierr = PetscInfo1(NULL," compression rate for spbas_compress_pattern %g \n",(double)mem_reduction);CHKERRQ(ierr); /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */ ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); if (info->usedt) { droptol = info->dt; } for (ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) { ierr = spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz,&matrix_LT);CHKERRQ(ierr); if (ierr == NEGATIVE_DIAGONAL) { shiftnz *= 1.5; if (shiftnz < 1e-5) shiftnz=1e-5; ierr = PetscInfo1(NULL,"spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n",(double)shiftnz);CHKERRQ(ierr); } } ierr = spbas_delete(Pattern);CHKERRQ(ierr); 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); ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); /* Convert spbas_matrix to compressed row storage */ ierr = spbas_transpose(matrix_LT, &matrix_L);CHKERRQ(ierr); ierr = spbas_delete(matrix_LT);CHKERRQ(ierr); ierr = spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj);CHKERRQ(ierr); b->i =bi; b->j=bj; b->a=ba; ierr = spbas_delete(matrix_L);CHKERRQ(ierr); /* Set the appropriate solution functions */ ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); if (perm_identity) { (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; } else { (B)->ops->solve = MatSolve_SeqSBAIJ_1_inplace; (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; } C->assembled = PETSC_TRUE; C->preallocated = PETSC_TRUE; ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* spbas_compress_pattern: calculate a compressed sparseness pattern for a sparseness pattern given in compressed row storage. The compressed sparseness pattern may require (much) less memory. */ PetscErrorCode spbas_compress_pattern(PetscInt *irow_in, PetscInt *icol_in, PetscInt nrows, PetscInt ncols, PetscInt col_idx_type, spbas_matrix *B,PetscReal *mem_reduction) { PetscInt nnz = irow_in[nrows]; size_t mem_orig = (nrows + nnz) * sizeof(PetscInt); size_t mem_compressed; PetscErrorCode ierr; PetscInt *isort; PetscInt *icols; PetscInt row_nnz; PetscInt *ipoint; PetscBool *used; PetscInt ptr; PetscInt i,j; const PetscBool no_values = PETSC_FALSE; PetscFunctionBegin; /* Allocate the structure of the new matrix */ B->nrows = nrows; B->ncols = ncols; B->nnz = nnz; B->col_idx_type = col_idx_type; B->block_data = PETSC_TRUE; ierr = spbas_allocate_pattern(B, no_values);CHKERRQ(ierr); /* When using an offset array, set it */ if (col_idx_type==SPBAS_OFFSET_ARRAY) { for (i=0; i<nrows; i++) B->icol0[i] = icol_in[irow_in[i]]; } /* Allocate the ordering for the rows */ ierr = PetscMalloc1(nrows,&isort);CHKERRQ(ierr); ierr = PetscMalloc1(nrows,&ipoint);CHKERRQ(ierr); ierr = PetscMalloc1(nrows,&used);CHKERRQ(ierr); /* Initialize the sorting */ ierr = PetscMemzero((void*) used, nrows*sizeof(PetscBool));CHKERRQ(ierr); for (i = 0; i<nrows; i++) { B->row_nnz[i] = irow_in[i+1]-irow_in[i]; isort[i] = i; ipoint[i] = i; } /* Sort the rows so that identical columns will be next to each other */ ierr = spbas_mergesort_icols(nrows, irow_in, icol_in, col_idx_type, isort);CHKERRQ(ierr); ierr = PetscInfo(NULL,"Rows have been sorted for patterns\n");CHKERRQ(ierr); /* Replace identical rows with the first one in the list */ for (i=1; i<nrows; i++) { if (spbas_row_order_icol(isort[i-1], isort[i], irow_in, icol_in, col_idx_type) == 0) { ipoint[isort[i]] = ipoint[isort[i-1]]; } } /* Collect the rows which are used*/ for (i=0; i<nrows; i++) used[ipoint[i]] = PETSC_TRUE; /* Calculate needed memory */ B->n_alloc_icol = 0; for (i=0; i<nrows; i++) { if (used[i]) B->n_alloc_icol += B->row_nnz[i]; } ierr = PetscMalloc1(B->n_alloc_icol,&B->alloc_icol);CHKERRQ(ierr); /* Fill in the diagonal offsets for the rows which store their own data */ ptr = 0; for (i=0; i<B->nrows; i++) { if (used[i]) { B->icols[i] = &B->alloc_icol[ptr]; icols = &icol_in[irow_in[i]]; row_nnz = B->row_nnz[i]; if (col_idx_type == SPBAS_COLUMN_NUMBERS) { for (j=0; j<row_nnz; j++) { B->icols[i][j] = icols[j]; } } else if (col_idx_type == SPBAS_DIAGONAL_OFFSETS) { for (j=0; j<row_nnz; j++) { B->icols[i][j] = icols[j]-i; } } else if (col_idx_type == SPBAS_OFFSET_ARRAY) { for (j=0; j<row_nnz; j++) { B->icols[i][j] = icols[j]-icols[0]; } } ptr += B->row_nnz[i]; } } /* Point to the right places for all data */ for (i=0; i<nrows; i++) { B->icols[i] = B->icols[ipoint[i]]; } ierr = PetscInfo(NULL,"Row patterns have been compressed\n");CHKERRQ(ierr); ierr = PetscInfo1(NULL," (%g nonzeros per row)\n", (double) ((PetscReal) nnz / (PetscReal) nrows));CHKERRQ(ierr); ierr=PetscFree(isort);CHKERRQ(ierr); ierr=PetscFree(used);CHKERRQ(ierr); ierr=PetscFree(ipoint);CHKERRQ(ierr); mem_compressed = spbas_memory_requirement(*B); *mem_reduction = 100.0 * (PetscReal)(mem_orig-mem_compressed)/ (PetscReal) mem_orig; PetscFunctionReturn(0); }