1 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2 3 /*@ 4 MatComputeBandwidth - Calculate the full bandwidth of the matrix, meaning the width 2k+1 where k diagonals on either side are sufficient to contain all the matrix nonzeros. 5 6 Collective on Mat 7 8 Input Parameters: 9 + A - The Mat 10 - fraction - An optional percentage of the Frobenius norm of the matrix that the bandwidth should enclose 11 12 Output Parameter: 13 . bw - The matrix bandwidth 14 15 Level: beginner 16 17 .seealso: `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()` 18 @*/ 19 PetscErrorCode MatComputeBandwidth(Mat A, PetscReal fraction, PetscInt *bw) { 20 PetscInt lbw[2] = {0, 0}, gbw[2]; 21 PetscInt rStart, rEnd, r; 22 23 PetscFunctionBegin; 24 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 25 PetscValidLogicalCollectiveReal(A, fraction, 2); 26 PetscValidIntPointer(bw, 3); 27 PetscCheck(!(fraction > 0.0) || !(fraction < 1.0), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "We do not yet support a fractional bandwidth"); 28 PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 29 for (r = rStart; r < rEnd; ++r) { 30 const PetscInt *cols; 31 PetscInt ncols; 32 33 PetscCall(MatGetRow(A, r, &ncols, &cols, NULL)); 34 if (ncols) { 35 lbw[0] = PetscMax(lbw[0], r - cols[0]); 36 lbw[1] = PetscMax(lbw[1], cols[ncols - 1] - r); 37 } 38 PetscCall(MatRestoreRow(A, r, &ncols, &cols, NULL)); 39 } 40 PetscCall(MPIU_Allreduce(lbw, gbw, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A))); 41 *bw = 2 * PetscMax(gbw[0], gbw[1]) + 1; 42 PetscFunctionReturn(0); 43 } 44