Actual source code: bandwidth.c
petsc-3.7.3 2016-08-01
1: #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
5: /*@
6: 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.
8: Collective on Mat
10: Input Parameters:
11: + A - The Mat
12: - fraction - An optional percentage of the Frobenius norm of the matrix that the bandwidth should enclose
14: Output Parameter:
15: . bw - The matrix bandwidth
17: Level: beginner
19: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
20: @*/
21: PetscErrorCode MatComputeBandwidth(Mat A, PetscReal fraction, PetscInt *bw)
22: {
23: PetscInt lbw[2] = {0, 0}, gbw[2];
24: PetscInt rStart, rEnd, r;
31: if ((fraction > 0.0) && (fraction < 1.0)) SETERRQ(PetscObjectComm((PetscObject) A), PETSC_ERR_SUP, "We do not yet support a fractional bandwidth");
32: MatGetOwnershipRange(A, &rStart, &rEnd);
33: for (r = rStart; r < rEnd; ++r) {
34: const PetscInt *cols;
35: PetscInt ncols;
37: MatGetRow(A, r, &ncols, &cols, NULL);
38: if (ncols) {
39: lbw[0] = PetscMax(lbw[0], r - cols[0]);
40: lbw[1] = PetscMax(lbw[1], cols[ncols-1] - r);
41: }
42: MatRestoreRow(A, r, &ncols, &cols, NULL);
43: }
44: MPIU_Allreduce(lbw, gbw, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) A));
45: *bw = 2*PetscMax(gbw[0], gbw[1]) + 1;
46: return(0);
47: }