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