1 static char help[] = "Example of inverting a block diagonal matrix.\n" 2 "\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc, char **args) { 7 Mat A, A_inv; 8 PetscMPIInt rank, size; 9 PetscInt M, m, bs, rstart, rend, j, x, y; 10 PetscInt *dnnz; 11 PetscScalar *v; 12 Vec X, Y; 13 PetscReal norm; 14 15 PetscFunctionBeginUser; 16 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 17 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 18 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 19 20 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex184", "Mat"); 21 M = 8; 22 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &M, NULL)); 23 bs = 3; 24 PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL)); 25 PetscOptionsEnd(); 26 27 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 28 PetscCall(MatSetFromOptions(A)); 29 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M * bs, M * bs)); 30 PetscCall(MatSetBlockSize(A, bs)); 31 PetscCall(MatSetUp(A)); 32 PetscCall(MatGetLocalSize(A, &m, NULL)); 33 PetscCall(PetscMalloc1(m / bs, &dnnz)); 34 for (j = 0; j < m / bs; j++) { dnnz[j] = 1; } 35 PetscCall(MatXAIJSetPreallocation(A, bs, dnnz, NULL, NULL, NULL)); 36 PetscCall(PetscFree(dnnz)); 37 38 PetscCall(PetscMalloc1(bs * bs, &v)); 39 PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 40 for (j = rstart / bs; j < rend / bs; j++) { 41 for (x = 0; x < bs; x++) { 42 for (y = 0; y < bs; y++) { 43 if (x == y) { 44 v[y + bs * x] = 2 * bs; 45 } else { 46 v[y + bs * x] = -1 * (x < y) - 2 * (x > y); 47 } 48 } 49 } 50 PetscCall(MatSetValuesBlocked(A, 1, &j, 1, &j, v, INSERT_VALUES)); 51 } 52 PetscCall(PetscFree(v)); 53 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 54 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 55 56 /* check that A = inv(inv(A)) */ 57 PetscCall(MatCreate(PETSC_COMM_WORLD, &A_inv)); 58 PetscCall(MatSetFromOptions(A_inv)); 59 PetscCall(MatInvertBlockDiagonalMat(A, A_inv)); 60 61 /* Test A_inv * A on a random vector */ 62 PetscCall(MatCreateVecs(A, &X, &Y)); 63 PetscCall(VecSetRandom(X, NULL)); 64 PetscCall(MatMult(A, X, Y)); 65 PetscCall(VecScale(X, -1)); 66 PetscCall(MatMultAdd(A_inv, Y, X, X)); 67 PetscCall(VecNorm(X, NORM_MAX, &norm)); 68 if (norm > PETSC_SMALL) { 69 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error exceeds tolerance.\nInverse of block diagonal A\n")); 70 PetscCall(MatView(A_inv, PETSC_VIEWER_STDOUT_WORLD)); 71 } 72 73 PetscCall(MatDestroy(&A)); 74 PetscCall(MatDestroy(&A_inv)); 75 PetscCall(VecDestroy(&X)); 76 PetscCall(VecDestroy(&Y)); 77 78 PetscCall(PetscFinalize()); 79 return 0; 80 } 81 82 /*TEST 83 test: 84 suffix: seqaij 85 args: -mat_type seqaij -mat_size 12 -mat_block_size 3 86 nsize: 1 87 test: 88 suffix: seqbaij 89 args: -mat_type seqbaij -mat_size 12 -mat_block_size 3 90 nsize: 1 91 test: 92 suffix: mpiaij 93 args: -mat_type mpiaij -mat_size 12 -mat_block_size 3 94 nsize: 2 95 test: 96 suffix: mpibaij 97 args: -mat_type mpibaij -mat_size 12 -mat_block_size 3 98 nsize: 2 99 TEST*/ 100