xref: /petsc/src/mat/tests/ex238.c (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 static char help[] = "Creates MatSeqBAIJ matrix of given BS for timing tests of MatMult().\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **args)
6 {
7   Mat         A;
8   Vec         x, y;
9   PetscInt    m = 50000, bs = 12, i, j, k, l, row, col, M, its = 25;
10   PetscScalar rval, *vals;
11   PetscRandom rdm;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &args, NULL, help));
15   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
16   PetscCall(PetscOptionsGetInt(NULL, NULL, "-its", &its, NULL));
17   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
18   M = m * bs;
19   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 27, NULL, &A));
20   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
21 
22   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
23   PetscCall(PetscRandomSetFromOptions(rdm));
24   PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &x));
25   PetscCall(VecDuplicate(x, &y));
26 
27   /* For each block row insert at most 27 blocks */
28   PetscCall(PetscMalloc1(bs * bs, &vals));
29   for (i = 0; i < m; i++) {
30     row = i;
31     for (j = 0; j < 27; j++) {
32       PetscCall(PetscRandomGetValue(rdm, &rval));
33       col = (PetscInt)(PetscRealPart(rval) * m);
34       for (k = 0; k < bs; k++) {
35         for (l = 0; l < bs; l++) {
36           PetscCall(PetscRandomGetValue(rdm, &rval));
37           vals[k * bs + l] = rval;
38         }
39       }
40       PetscCall(MatSetValuesBlocked(A, 1, &row, 1, &col, vals, INSERT_VALUES));
41     }
42   }
43   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
44   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
45   PetscCall(PetscFree(vals));
46 
47   /* Time MatMult(), MatMultAdd() */
48   for (i = 0; i < its; i++) {
49     PetscCall(VecSetRandom(x, rdm));
50     PetscCall(MatMult(A, x, y));
51     PetscCall(VecSetRandom(x, rdm));
52     PetscCall(VecSetRandom(y, rdm));
53     PetscCall(MatMultAdd(A, x, y, y));
54   }
55 
56   PetscCall(MatDestroy(&A));
57   PetscCall(VecDestroy(&x));
58   PetscCall(VecDestroy(&y));
59   PetscCall(PetscRandomDestroy(&rdm));
60   PetscCall(PetscFinalize());
61   return 0;
62 }
63 
64 /*TEST
65 
66    testset:
67      requires: defined(PETSC_USING_64BIT_PTR)
68      output_file: output/empty.out
69      test:
70        suffix: 1
71        args: -mat_block_size 1 -mat_size 1000 -its 2
72      test:
73        suffix: 2
74        args: -mat_block_size 2 -mat_size 1000 -its 2
75      test:
76        suffix: 4
77        args: -mat_block_size 4 -mat_size 1000 -its 2
78      test:
79        suffix: 5
80        args: -mat_block_size 5 -mat_size 1000 -its 2
81      test:
82        suffix: 6
83        args: -mat_block_size 6 -mat_size 1000 -its 2
84      test:
85        suffix: 8
86        args: -mat_block_size 8 -mat_size 1000 -its 2
87      test:
88        suffix: 12
89        args: -mat_block_size 12 -mat_size 1000 -its 2
90      test:
91        suffix: 15
92        args: -mat_block_size 15 -mat_size 1000 -its 2
93 
94 TEST*/
95