1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2c13ced23SMatthew G. Knepley
3ddc19b50SMatthew G. Knepley /*@
4ddc19b50SMatthew G. Knepley 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.
5c13ced23SMatthew G. Knepley
6c3339decSBarry Smith Collective
7c13ced23SMatthew G. Knepley
8c13ced23SMatthew G. Knepley Input Parameters:
92ef1f0ffSBarry Smith + A - The `Mat`
10c13ced23SMatthew G. Knepley - fraction - An optional percentage of the Frobenius norm of the matrix that the bandwidth should enclose
11c13ced23SMatthew G. Knepley
12c13ced23SMatthew G. Knepley Output Parameter:
13c13ced23SMatthew G. Knepley . bw - The matrix bandwidth
14c13ced23SMatthew G. Knepley
15c13ced23SMatthew G. Knepley Level: beginner
16c13ced23SMatthew G. Knepley
17db781477SPatrick Sanan .seealso: `DMPlexCreate()`, `DMPlexSetConeSize()`, `DMPlexSetChart()`
18c13ced23SMatthew G. Knepley @*/
MatComputeBandwidth(Mat A,PetscReal fraction,PetscInt * bw)19d71ae5a4SJacob Faibussowitsch PetscErrorCode MatComputeBandwidth(Mat A, PetscReal fraction, PetscInt *bw)
20d71ae5a4SJacob Faibussowitsch {
211ff9db0eSJed Brown PetscInt lbw[2] = {0, 0}, gbw[2];
22c13ced23SMatthew G. Knepley PetscInt rStart, rEnd, r;
23c13ced23SMatthew G. Knepley
24c13ced23SMatthew G. Knepley PetscFunctionBegin;
25c13ced23SMatthew G. Knepley PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
26c13ced23SMatthew G. Knepley PetscValidLogicalCollectiveReal(A, fraction, 2);
274f572ea9SToby Isaac PetscAssertPointer(bw, 3);
2808401ef6SPierre Jolivet PetscCheck(!(fraction > 0.0) || !(fraction < 1.0), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "We do not yet support a fractional bandwidth");
299566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
30c13ced23SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) {
31c13ced23SMatthew G. Knepley const PetscInt *cols;
32c13ced23SMatthew G. Knepley PetscInt ncols;
33c13ced23SMatthew G. Knepley
349566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &cols, NULL));
35c13ced23SMatthew G. Knepley if (ncols) {
36c13ced23SMatthew G. Knepley lbw[0] = PetscMax(lbw[0], r - cols[0]);
37c13ced23SMatthew G. Knepley lbw[1] = PetscMax(lbw[1], cols[ncols - 1] - r);
38c13ced23SMatthew G. Knepley }
399566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &cols, NULL));
40c13ced23SMatthew G. Knepley }
41*462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(lbw, gbw, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A)));
42c13ced23SMatthew G. Knepley *bw = 2 * PetscMax(gbw[0], gbw[1]) + 1;
433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
44c13ced23SMatthew G. Knepley }
45