xref: /petsc/src/mat/tests/ex109.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
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,AT;
8   PetscInt       i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an;
9   PetscScalar    v;
10   PetscRandom    r;
11   PetscBool      equal=PETSC_FALSE,flg;
12   PetscReal      fill = 1.0,norm;
13   PetscMPIInt    size;
14 
15   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
16   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
17   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
18   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
19 
20   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
21   PetscCall(PetscRandomSetFromOptions(r));
22 
23   /* Create a aij matrix A */
24   M    = N = m*n;
25   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
26   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
27   PetscCall(MatSetType(A,MATAIJ));
28   PetscCall(MatSetFromOptions(A));
29   PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
30   PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
31 
32   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
33   am   = Iend - Istart;
34   for (Ii=Istart; Ii<Iend; Ii++) {
35     v = -1.0; i = Ii/n; j = Ii - i*n;
36     if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
37     if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
38     if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
39     if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
40     v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
41   }
42   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
43   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
44 
45   /* Create a dense matrix B */
46   PetscCall(MatGetLocalSize(A,&am,&an));
47   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
48   PetscCall(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M));
49   PetscCall(MatSetType(B,MATDENSE));
50   PetscCall(MatSeqDenseSetPreallocation(B,NULL));
51   PetscCall(MatMPIDenseSetPreallocation(B,NULL));
52   PetscCall(MatSetFromOptions(B));
53   PetscCall(MatSetRandom(B,r));
54   PetscCall(PetscRandomDestroy(&r));
55 
56   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
57   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
58   PetscCall(MatSetType(C,MATDENSE));
59   PetscCall(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N));
60   PetscCall(MatSetFromOptions(C));
61   PetscCall(MatSetUp(C));
62   PetscCall(MatZeroEntries(C));
63   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
64   PetscCall(MatNorm(C,NORM_INFINITY,&norm));
65   PetscCall(MatDestroy(&C));
66 
67   /* Test C = A*B (aij*dense) */
68   PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
69   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
70 
71   /* Test developer API */
72   PetscCall(MatProductCreate(A,B,NULL,&D));
73   PetscCall(MatProductSetType(D,MATPRODUCT_AB));
74   PetscCall(MatProductSetAlgorithm(D,"default"));
75   PetscCall(MatProductSetFill(D,fill));
76   PetscCall(MatProductSetFromOptions(D));
77   PetscCall(MatProductSymbolic(D));
78   for (i=0; i<2; i++) {
79     PetscCall(MatProductNumeric(D));
80   }
81   PetscCall(MatEqual(C,D,&equal));
82   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
83   PetscCall(MatDestroy(&D));
84 
85   /* Test D = AT*B (transpose(aij)*dense) */
86   PetscCall(MatCreateTranspose(A,&AT));
87   PetscCall(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D));
88   PetscCall(MatMatMultEqual(AT,B,D,10,&equal));
89   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)");
90   PetscCall(MatDestroy(&D));
91   PetscCall(MatDestroy(&AT));
92 
93   /* Test D = C*A (dense*aij) */
94   PetscCall(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D));
95   PetscCall(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
96   PetscCall(MatMatMultEqual(C,A,D,10,&equal));
97   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)");
98   PetscCall(MatDestroy(&D));
99 
100   /* Test D = A*C (aij*dense) */
101   PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D));
102   PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D));
103   PetscCall(MatMatMultEqual(A,C,D,10,&equal));
104   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)");
105   PetscCall(MatDestroy(&D));
106 
107   /* Test D = B*C (dense*dense) */
108   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
109   if (size == 1) {
110     PetscCall(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
111     PetscCall(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D));
112     PetscCall(MatMatMultEqual(B,C,D,10,&equal));
113     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)");
114     PetscCall(MatDestroy(&D));
115   }
116 
117   /* Test D = B*C^T (dense*dense) */
118   PetscCall(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
119   PetscCall(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D));
120   PetscCall(MatMatTransposeMultEqual(B,C,D,10,&equal));
121   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)");
122   PetscCall(MatDestroy(&D));
123 
124   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
125   flg = PETSC_FALSE;
126   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg));
127   if (flg) {
128     /* user driver */
129     PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B));
130   } else {
131     /* clear internal data structures related with previous products to avoid circular references */
132     PetscCall(MatProductClear(A));
133     PetscCall(MatProductClear(B));
134     PetscCall(MatProductClear(C));
135     PetscCall(MatProductCreateWithMat(A,C,NULL,B));
136     PetscCall(MatProductSetType(B,MATPRODUCT_AB));
137     PetscCall(MatProductSetFromOptions(B));
138     PetscCall(MatProductSymbolic(B));
139     PetscCall(MatProductNumeric(B));
140   }
141 
142   PetscCall(MatDestroy(&C));
143   PetscCall(MatDestroy(&B));
144   PetscCall(MatDestroy(&A));
145   PetscCall(PetscFinalize());
146   return 0;
147 }
148 
149 /*TEST
150 
151    test:
152       args: -M 10 -N 10
153       output_file: output/ex109.out
154 
155    test:
156       suffix: 2
157       nsize: 3
158       output_file: output/ex109.out
159 
160    test:
161       suffix: 3
162       nsize: 2
163       args: -matmattransmult_mpidense_mpidense_via cyclic
164       output_file: output/ex109.out
165 
166    test:
167       suffix: 4
168       args: -test_userAPI
169       output_file: output/ex109.out
170 
171    test:
172       suffix: 5
173       nsize: 3
174       args: -test_userAPI
175       output_file: output/ex109.out
176 
177 TEST*/
178