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