xref: /petsc/src/mat/tests/ex109.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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,D;
8   PetscInt       i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an;
9   PetscScalar    v;
10   PetscErrorCode ierr;
11   PetscRandom    r;
12   PetscBool      equal=PETSC_FALSE;
13   PetscReal      fill = 1.0;
14   PetscMPIInt    size;
15 
16   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
18   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
19   ierr = PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);CHKERRQ(ierr);
20 
21   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr);
22   ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
23 
24   /* Create a aij matrix A */
25   M    = N = m*n;
26   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
27   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
28   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
29   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
30   ierr = MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);CHKERRQ(ierr);
31   ierr = MatSeqAIJSetPreallocation(A,5,NULL);CHKERRQ(ierr);
32 
33   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
34   am   = Iend - Istart;
35   for (Ii=Istart; Ii<Iend; Ii++) {
36     v = -1.0; i = Ii/n; j = Ii - i*n;
37     if (i>0)   {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
38     if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
39     if (j>0)   {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
40     if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
41     v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
42   }
43   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
44   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45 
46   /* Create a dense matrix B */
47   ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr);
48   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
49   ierr = MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M);CHKERRQ(ierr);
50   ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr);
51   ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
52   ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr);
53   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
54   ierr = MatSetRandom(B,r);CHKERRQ(ierr);
55   ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
56   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
58 
59   /* Test C = A*B (aij*dense) */
60   ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr);
61   ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr);
62 
63   ierr = MatMatMultSymbolic(A,B,fill,&D);CHKERRQ(ierr);
64   for (i=0; i<2; i++) {
65     ierr = MatMatMultNumeric(A,B,D);CHKERRQ(ierr);
66   }
67   ierr = MatEqual(C,D,&equal);CHKERRQ(ierr);
68   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
69   ierr = MatDestroy(&D);CHKERRQ(ierr);
70 
71   /* Test D = C*A (dense*aij) */
72   ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr);
73   ierr = MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr);
74   ierr = MatDestroy(&D);CHKERRQ(ierr);
75 
76   /* Test D = A*C (aij*dense) */
77   ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr);
78   ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr);
79   ierr = MatDestroy(&D);CHKERRQ(ierr);
80 
81   /* Test D = B*C (dense*dense) */
82   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
83   if (size == 1) {
84     ierr = MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr);
85     ierr = MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr);
86     ierr = MatDestroy(&D);CHKERRQ(ierr);
87   }
88 
89   ierr = MatDestroy(&C);CHKERRQ(ierr);
90   ierr = MatDestroy(&B);CHKERRQ(ierr);
91   ierr = MatDestroy(&A);CHKERRQ(ierr);
92   ierr = PetscFinalize();
93   return ierr;
94 }
95 
96 
97 /*TEST
98 
99    test:
100       args: -M 10 -N 10
101       output_file: output/ex109.out
102 
103    test:
104       suffix: 2
105       nsize: 3
106       output_file: output/ex109.out
107 
108 TEST*/
109