xref: /petsc/src/mat/tests/ex34.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() and MatTransposeMatMult() for MPIAIJ and MPIDENSE matrices. \n\
2ad16ce7aSStefano Zampini                       Sequential part of mpidense matrix allows changes made by MatDenseSetLDA(). \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petsc.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc, char ** argv)
7c4762a1bSJed Brown {
8ad16ce7aSStefano Zampini   Mat            A, B, C, C1;
9c4762a1bSJed Brown   PetscMPIInt    size;
10c4762a1bSJed Brown   PetscInt       i,ia[2] = { 0, 2 }, ja[2] = { 0, 1 }, lda = 4;
11c4762a1bSJed Brown   PetscScalar    a[2] = { 1.0, 1.0 }, *data;
12c4762a1bSJed Brown   PetscBool      flg;
13c4762a1bSJed Brown 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, (char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
162c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must use 2 processors");
17c4762a1bSJed Brown 
185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A, MATMPIAIJ));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, 1, 1, 2, 2));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocationCSR(A, ia, ja, a));
22c4762a1bSJed Brown 
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(4 * lda,&data));
24c4762a1bSJed Brown   for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0;
25c4762a1bSJed Brown 
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(B,"b_"));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseSetLDA(B, lda));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Test MatMatMult() */
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
34c4762a1bSJed Brown 
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,B,C,10,&flg));
3628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C");
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Test user-provided mpidense matrix product */
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(A,B,C1,10,&flg));
4228b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1");
43c4762a1bSJed Brown 
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C1));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Test MatTransposeMatMult() */
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMultEqual(A,B,C,10,&flg));
5228b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()");
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
54c4762a1bSJed Brown 
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(data));
58c4762a1bSJed Brown 
59*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
60*b122ec5aSJacob Faibussowitsch   return 0;
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63c4762a1bSJed Brown /*TEST
64c4762a1bSJed Brown 
65c4762a1bSJed Brown    test:
66c4762a1bSJed Brown       suffix: 1
67c4762a1bSJed Brown       nsize: 2
68c4762a1bSJed Brown       output_file: output/ex34.out
69a5225ed3SStefano Zampini 
70a5225ed3SStefano Zampini    test:
71a5225ed3SStefano Zampini       suffix: 1_cuda
72a5225ed3SStefano Zampini       requires: cuda
73a5225ed3SStefano Zampini       nsize: 2
74a5225ed3SStefano Zampini       args: -b_mat_type mpidensecuda
75a5225ed3SStefano Zampini       output_file: output/ex34.out
76c4762a1bSJed Brown TEST*/
77