xref: /petsc/src/mat/tests/ex122.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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   PetscRandom    r;
10   PetscBool      equal=PETSC_FALSE;
11   PetscReal      fill = 1.0;
12   PetscInt       nza,am,an;
13 
14   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
15   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
16   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
17   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
18 
19   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
20   PetscCall(PetscRandomSetFromOptions(r));
21 
22   /* create a aij matrix A */
23   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
24   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M));
25   PetscCall(MatSetType(A,MATAIJ));
26   nza  = (PetscInt)(.3*M); /* num of nozeros in each row of A */
27   PetscCall(MatSeqAIJSetPreallocation(A,nza,NULL));
28   PetscCall(MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL));
29   PetscCall(MatSetRandom(A,r));
30 
31   /* create a dense matrix B */
32   PetscCall(MatGetLocalSize(A,&am,&an));
33   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
34   PetscCall(MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE));
35   PetscCall(MatSetType(B,MATDENSE));
36   PetscCall(MatSeqDenseSetPreallocation(B,NULL));
37   PetscCall(MatMPIDenseSetPreallocation(B,NULL));
38   PetscCall(MatSetRandom(B,r));
39   PetscCall(PetscRandomDestroy(&r));
40   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
41   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
42 
43   /* Test MatMatMult() */
44   PetscCall(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C));
45   PetscCall(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
46   PetscCall(MatMatMultEqual(B,A,C,10,&equal));
47   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A");
48 
49   PetscCall(MatDestroy(&C));
50   PetscCall(MatDestroy(&B));
51   PetscCall(MatDestroy(&A));
52   PetscCall(PetscFinalize());
53   return 0;
54 }
55 
56 /*TEST
57 
58    test:
59       output_file: output/ex122.out
60 
61 TEST*/
62