xref: /petsc/src/mat/tests/ex34.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(PetscInitialize(&argc, &argv, (char*)0,help));
15   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16   PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must use 2 processors");
17 
18   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
19   CHKERRQ(MatSetType(A, MATMPIAIJ));
20   CHKERRQ(MatSetSizes(A, 1, 1, 2, 2));
21   CHKERRQ(MatMPIAIJSetPreallocationCSR(A, ia, ja, a));
22 
23   CHKERRQ(PetscCalloc1(4 * lda,&data));
24   for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0;
25 
26   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B));
27   CHKERRQ(MatSetOptionsPrefix(B,"b_"));
28   CHKERRQ(MatSetFromOptions(B));
29   CHKERRQ(MatDenseSetLDA(B, lda));
30 
31   /* Test MatMatMult() */
32   CHKERRQ(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
33   CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
34 
35   CHKERRQ(MatMatMultEqual(A,B,C,10,&flg));
36   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C");
37 
38   /* Test user-provided mpidense matrix product */
39   CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1));
40   CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1));
41   CHKERRQ(MatMatMultEqual(A,B,C1,10,&flg));
42   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1");
43 
44   CHKERRQ(MatDestroy(&C1));
45   CHKERRQ(MatDestroy(&C));
46 
47   /* Test MatTransposeMatMult() */
48   CHKERRQ(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
49   CHKERRQ(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
50 
51   CHKERRQ(MatTransposeMatMultEqual(A,B,C,10,&flg));
52   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()");
53   CHKERRQ(MatDestroy(&C));
54 
55   CHKERRQ(MatDestroy(&B));
56   CHKERRQ(MatDestroy(&A));
57   CHKERRQ(PetscFree(data));
58 
59   CHKERRQ(PetscFinalize());
60   return 0;
61 }
62 
63 /*TEST
64 
65    test:
66       suffix: 1
67       nsize: 2
68       output_file: output/ex34.out
69 
70    test:
71       suffix: 1_cuda
72       requires: cuda
73       nsize: 2
74       args: -b_mat_type mpidensecuda
75       output_file: output/ex34.out
76 TEST*/
77