xref: /petsc/src/mat/tests/ex34.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 static char help[] = "Test MatMatMult() and MatTransposeMatMult() for MPIAIJ and MPIDENSE matrices. \n\
2                       Sequential part of mpidense matrix allows changes made by MatDenseSetLDA(). \n\n";
3 
4 #include <petsc.h>
5 
6 int main(int argc, char **argv)
7 {
8   Mat         A, B, C, C1;
9   PetscMPIInt size;
10   PetscInt    i, ia[2] = {0, 2}, ja[2] = {0, 1}, lda = 4;
11   PetscScalar a[2] = {1.0, 1.0}, *data;
12   PetscBool   flg;
13 
14   PetscFunctionBeginUser;
15   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
16   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must use 2 processors");
18 
19   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
20   PetscCall(MatSetType(A, MATMPIAIJ));
21   PetscCall(MatSetSizes(A, 1, 1, 2, 2));
22   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a));
23 
24   PetscCall(PetscCalloc1(4 * lda, &data));
25   for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0;
26 
27   PetscCall(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B));
28   PetscCall(MatSetOptionsPrefix(B, "b_"));
29   PetscCall(MatSetFromOptions(B));
30   PetscCall(MatDenseSetLDA(B, lda));
31 
32   /* Test MatMatMult() */
33   PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
34   PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
35 
36   PetscCall(MatMatMultEqual(A, B, C, 10, &flg));
37   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatMatMult() for C");
38 
39   /* Test user-provided mpidense matrix product */
40   PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &C1));
41   PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1));
42   PetscCall(MatMatMultEqual(A, B, C1, 10, &flg));
43   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatMatMult() for C1");
44 
45   PetscCall(MatDestroy(&C1));
46   PetscCall(MatDestroy(&C));
47 
48   /* Test MatTransposeMatMult() */
49   PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
50   PetscCall(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
51 
52   PetscCall(MatTransposeMatMultEqual(A, B, C, 10, &flg));
53   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatTransposeMatMult()");
54   PetscCall(MatDestroy(&C));
55 
56   PetscCall(MatDestroy(&B));
57   PetscCall(MatDestroy(&A));
58   PetscCall(PetscFree(data));
59 
60   PetscCall(PetscFinalize());
61   return 0;
62 }
63 
64 /*TEST
65 
66    test:
67       suffix: 1
68       nsize: 2
69       output_file: output/ex34.out
70 
71    test:
72       suffix: 1_cuda
73       requires: cuda
74       nsize: 2
75       args: -b_mat_type mpidensecuda
76       output_file: output/ex34.out
77 TEST*/
78