xref: /petsc/src/mat/tests/ex122.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **argv)
6 {
7   Mat            A,B,C;
8   PetscInt       M=10,N=5;
9   PetscErrorCode ierr;
10   PetscRandom    r;
11   PetscBool      equal=PETSC_FALSE;
12   PetscReal      fill = 1.0;
13   PetscInt       nza,am,an;
14 
15   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
17   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
18   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
19 
20   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
21   CHKERRQ(PetscRandomSetFromOptions(r));
22 
23   /* create a aij matrix A */
24   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
25   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M));
26   CHKERRQ(MatSetType(A,MATAIJ));
27   nza  = (PetscInt)(.3*M); /* num of nozeros in each row of A */
28   CHKERRQ(MatSeqAIJSetPreallocation(A,nza,NULL));
29   CHKERRQ(MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL));
30   CHKERRQ(MatSetRandom(A,r));
31 
32   /* create a dense matrix B */
33   CHKERRQ(MatGetLocalSize(A,&am,&an));
34   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
35   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE));
36   CHKERRQ(MatSetType(B,MATDENSE));
37   CHKERRQ(MatSeqDenseSetPreallocation(B,NULL));
38   CHKERRQ(MatMPIDenseSetPreallocation(B,NULL));
39   CHKERRQ(MatSetRandom(B,r));
40   CHKERRQ(PetscRandomDestroy(&r));
41   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
42   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
43 
44   /* Test MatMatMult() */
45   CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C));
46   CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
47   CHKERRQ(MatMatMultEqual(B,A,C,10,&equal));
48   PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A");
49 
50   CHKERRQ(MatDestroy(&C));
51   CHKERRQ(MatDestroy(&B));
52   CHKERRQ(MatDestroy(&A));
53   ierr = PetscFinalize();
54   return ierr;
55 }
56 
57 /*TEST
58 
59    test:
60       output_file: output/ex122.out
61 
62 TEST*/
63