xref: /petsc/src/mat/tests/ex122.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **argv) {
6   Mat         A, B, C;
7   PetscInt    M = 10, N = 5;
8   PetscRandom r;
9   PetscBool   equal = PETSC_FALSE;
10   PetscReal   fill  = 1.0;
11   PetscInt    nza, am, an;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
15   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
16   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
17   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL));
18 
19   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r));
20   PetscCall(PetscRandomSetFromOptions(r));
21 
22   /* create a aij matrix A */
23   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
24   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, M));
25   PetscCall(MatSetType(A, MATAIJ));
26   nza = (PetscInt)(.3 * M); /* num of nozeros in each row of A */
27   PetscCall(MatSeqAIJSetPreallocation(A, nza, NULL));
28   PetscCall(MatMPIAIJSetPreallocation(A, nza, NULL, nza, NULL));
29   PetscCall(MatSetRandom(A, r));
30 
31   /* create a dense matrix B */
32   PetscCall(MatGetLocalSize(A, &am, &an));
33   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
34   PetscCall(MatSetSizes(B, PETSC_DECIDE, am, N, PETSC_DECIDE));
35   PetscCall(MatSetType(B, MATDENSE));
36   PetscCall(MatSeqDenseSetPreallocation(B, NULL));
37   PetscCall(MatMPIDenseSetPreallocation(B, NULL));
38   PetscCall(MatSetRandom(B, r));
39   PetscCall(PetscRandomDestroy(&r));
40 
41   /* Test MatMatMult() */
42   PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C));
43   PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C));
44   PetscCall(MatMatMultEqual(B, A, C, 10, &equal));
45   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "C != B*A");
46 
47   PetscCall(MatDestroy(&C));
48   PetscCall(MatDestroy(&B));
49   PetscCall(MatDestroy(&A));
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55 
56    test:
57       output_file: output/ex122.out
58 
59 TEST*/
60