xref: /petsc/src/mat/tests/ex184.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Example of inverting a block diagonal matrix.\n"
2c4762a1bSJed Brown                      "\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*9371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat          A, A_inv;
8c4762a1bSJed Brown   PetscMPIInt  rank, size;
9c4762a1bSJed Brown   PetscInt     M, m, bs, rstart, rend, j, x, y;
10c4762a1bSJed Brown   PetscInt    *dnnz;
11c4762a1bSJed Brown   PetscScalar *v;
12c4762a1bSJed Brown   Vec          X, Y;
13c4762a1bSJed Brown   PetscReal    norm;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
19c4762a1bSJed Brown 
20d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex184", "Mat");
21c4762a1bSJed Brown   M = 8;
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &M, NULL));
23c4762a1bSJed Brown   bs = 3;
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
25d0609cedSBarry Smith   PetscOptionsEnd();
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
289566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M * bs, M * bs));
309566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A, bs));
319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
329566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m / bs, &dnnz));
34*9371c9d4SSatish Balay   for (j = 0; j < m / bs; j++) { dnnz[j] = 1; }
359566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(A, bs, dnnz, NULL, NULL, NULL));
369566063dSJacob Faibussowitsch   PetscCall(PetscFree(dnnz));
37c4762a1bSJed Brown 
389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &v));
399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
40c4762a1bSJed Brown   for (j = rstart / bs; j < rend / bs; j++) {
41c4762a1bSJed Brown     for (x = 0; x < bs; x++) {
42c4762a1bSJed Brown       for (y = 0; y < bs; y++) {
43c4762a1bSJed Brown         if (x == y) {
44c4762a1bSJed Brown           v[y + bs * x] = 2 * bs;
45c4762a1bSJed Brown         } else {
46c4762a1bSJed Brown           v[y + bs * x] = -1 * (x < y) - 2 * (x > y);
47c4762a1bSJed Brown         }
48c4762a1bSJed Brown       }
49c4762a1bSJed Brown     }
509566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &j, 1, &j, v, INSERT_VALUES));
51c4762a1bSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* check that A  = inv(inv(A)) */
579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A_inv));
589566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A_inv));
599566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonalMat(A, A_inv));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Test A_inv * A on a random vector */
629566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &X, &Y));
639566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(X, NULL));
649566063dSJacob Faibussowitsch   PetscCall(MatMult(A, X, Y));
659566063dSJacob Faibussowitsch   PetscCall(VecScale(X, -1));
669566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(A_inv, Y, X, X));
679566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_MAX, &norm));
68c4762a1bSJed Brown   if (norm > PETSC_SMALL) {
699566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error exceeds tolerance.\nInverse of block diagonal A\n"));
709566063dSJacob Faibussowitsch     PetscCall(MatView(A_inv, PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A_inv));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
79b122ec5aSJacob Faibussowitsch   return 0;
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown /*TEST
83c4762a1bSJed Brown   test:
84c4762a1bSJed Brown     suffix: seqaij
85c4762a1bSJed Brown     args: -mat_type seqaij -mat_size 12 -mat_block_size 3
86c4762a1bSJed Brown     nsize: 1
87c4762a1bSJed Brown   test:
88c4762a1bSJed Brown     suffix: seqbaij
89c4762a1bSJed Brown     args: -mat_type seqbaij -mat_size 12 -mat_block_size 3
90c4762a1bSJed Brown     nsize: 1
91c4762a1bSJed Brown   test:
92c4762a1bSJed Brown     suffix: mpiaij
93c4762a1bSJed Brown     args: -mat_type mpiaij -mat_size 12 -mat_block_size 3
94c4762a1bSJed Brown     nsize: 2
95c4762a1bSJed Brown   test:
96c4762a1bSJed Brown     suffix: mpibaij
97c4762a1bSJed Brown     args: -mat_type mpibaij -mat_size 12 -mat_block_size 3
98c4762a1bSJed Brown     nsize: 2
99c4762a1bSJed Brown TEST*/
100