xref: /petsc/src/mat/tests/ex184.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
6c4762a1bSJed Brown int main(int argc, char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown     Mat            A,A_inv;
9c4762a1bSJed Brown     PetscMPIInt    rank,size;
10c4762a1bSJed Brown     PetscInt       M,m,bs,rstart,rend,j,x,y;
11c4762a1bSJed Brown     PetscInt*      dnnz;
12c4762a1bSJed Brown     PetscScalar    *v;
13c4762a1bSJed Brown     Vec            X, Y;
14c4762a1bSJed Brown     PetscReal      norm;
15c4762a1bSJed Brown 
16*327415f7SBarry Smith     PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
20c4762a1bSJed Brown 
21d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex184","Mat");
22c4762a1bSJed Brown     M=8;
239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&M,NULL));
24c4762a1bSJed Brown     bs=3;
259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
26d0609cedSBarry Smith     PetscOptionsEnd();
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
299566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
309566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M*bs,M*bs));
319566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(A,bs));
329566063dSJacob Faibussowitsch     PetscCall(MatSetUp(A));
339566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A,&m,NULL));
349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m/bs,&dnnz));
35c4762a1bSJed Brown     for (j = 0; j < m/bs; j++) {
36c4762a1bSJed Brown         dnnz[j] = 1;
37c4762a1bSJed Brown     }
389566063dSJacob Faibussowitsch     PetscCall(MatXAIJSetPreallocation(A,bs,dnnz,NULL,NULL,NULL));
399566063dSJacob Faibussowitsch     PetscCall(PetscFree(dnnz));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs*bs,&v));
429566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
43c4762a1bSJed Brown     for (j = rstart/bs; j < rend/bs; j++) {
44c4762a1bSJed Brown         for (x = 0; x < bs; x++) {
45c4762a1bSJed Brown             for (y = 0; y < bs; y++) {
46c4762a1bSJed Brown                 if (x == y) {
47c4762a1bSJed Brown                     v[y+bs*x] = 2*bs;
48c4762a1bSJed Brown                 } else {
49c4762a1bSJed Brown                     v[y+bs*x] = -1 * (x < y) - 2 * (x > y);
50c4762a1bSJed Brown                 }
51c4762a1bSJed Brown             }
52c4762a1bSJed Brown         }
539566063dSJacob Faibussowitsch         PetscCall(MatSetValuesBlocked(A,1,&j,1,&j,v,INSERT_VALUES));
54c4762a1bSJed Brown     }
559566063dSJacob Faibussowitsch     PetscCall(PetscFree(v));
569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
579566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     /* check that A  = inv(inv(A)) */
609566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD,&A_inv));
619566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A_inv));
629566063dSJacob Faibussowitsch     PetscCall(MatInvertBlockDiagonalMat(A,A_inv));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown     /* Test A_inv * A on a random vector */
659566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(A, &X, &Y));
669566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(X, NULL));
679566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, Y));
689566063dSJacob Faibussowitsch     PetscCall(VecScale(X, -1));
699566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A_inv, Y, X, X));
709566063dSJacob Faibussowitsch     PetscCall(VecNorm(X, NORM_MAX, &norm));
71c4762a1bSJed Brown     if (norm > PETSC_SMALL) {
729566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error exceeds tolerance.\nInverse of block diagonal A\n"));
739566063dSJacob Faibussowitsch         PetscCall(MatView(A_inv,PETSC_VIEWER_STDOUT_WORLD));
74c4762a1bSJed Brown     }
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
779566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_inv));
789566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X));
799566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Y));
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch     PetscCall(PetscFinalize());
82b122ec5aSJacob Faibussowitsch     return 0;
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /*TEST
86c4762a1bSJed Brown   test:
87c4762a1bSJed Brown     suffix: seqaij
88c4762a1bSJed Brown     args: -mat_type seqaij -mat_size 12 -mat_block_size 3
89c4762a1bSJed Brown     nsize: 1
90c4762a1bSJed Brown   test:
91c4762a1bSJed Brown     suffix: seqbaij
92c4762a1bSJed Brown     args: -mat_type seqbaij -mat_size 12 -mat_block_size 3
93c4762a1bSJed Brown     nsize: 1
94c4762a1bSJed Brown   test:
95c4762a1bSJed Brown     suffix: mpiaij
96c4762a1bSJed Brown     args: -mat_type mpiaij -mat_size 12 -mat_block_size 3
97c4762a1bSJed Brown     nsize: 2
98c4762a1bSJed Brown   test:
99c4762a1bSJed Brown     suffix: mpibaij
100c4762a1bSJed Brown     args: -mat_type mpibaij -mat_size 12 -mat_block_size 3
101c4762a1bSJed Brown     nsize: 2
102c4762a1bSJed Brown TEST*/
103