Actual source code: bandwidth.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1:  #include <petsc/private/matimpl.h>

  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.

  6:   Collective on Mat

  8:   Input Parameters:
  9: + A - The Mat
 10: - fraction - An optional percentage of the Frobenius norm of the matrix that the bandwidth should enclose

 12:   Output Parameter:
 13: . bw - The matrix bandwidth

 15:   Level: beginner

 17: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
 18: @*/
 19: PetscErrorCode MatComputeBandwidth(Mat A, PetscReal fraction, PetscInt *bw)
 20: {
 21:   PetscInt       lbw[2] = {0, 0}, gbw[2];
 22:   PetscInt       rStart, rEnd, r;

 29:   if ((fraction > 0.0) && (fraction < 1.0)) SETERRQ(PetscObjectComm((PetscObject) A), PETSC_ERR_SUP, "We do not yet support a fractional bandwidth");
 30:   MatGetOwnershipRange(A, &rStart, &rEnd);
 31:   for (r = rStart; r < rEnd; ++r) {
 32:     const PetscInt *cols;
 33:     PetscInt        ncols;

 35:     MatGetRow(A, r, &ncols, &cols, NULL);
 36:     if (ncols) {
 37:       lbw[0] = PetscMax(lbw[0], r - cols[0]);
 38:       lbw[1] = PetscMax(lbw[1], cols[ncols-1] - r);
 39:     }
 40:     MatRestoreRow(A, r, &ncols, &cols, NULL);
 41:   }
 42:   MPIU_Allreduce(lbw, gbw, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) A));
 43:   *bw = 2*PetscMax(gbw[0], gbw[1]) + 1;
 44:   return(0);
 45: }