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