xref: /petsc/src/mat/tests/ex109.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c20d7725SJed Brown   Mat            A,B,C,D,AT;
8c4762a1bSJed Brown   PetscInt       i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an;
9c4762a1bSJed Brown   PetscScalar    v;
10c4762a1bSJed Brown   PetscRandom    r;
113fbca975SStefano Zampini   PetscBool      equal=PETSC_FALSE,flg;
12c20d7725SJed Brown   PetscReal      fill = 1.0,norm;
13c4762a1bSJed Brown   PetscMPIInt    size;
14c4762a1bSJed Brown 
15*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
16*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
17*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
19c4762a1bSJed Brown 
20*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
21*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* Create a aij matrix A */
24c4762a1bSJed Brown   M    = N = m*n;
25*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
26*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
27*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
28*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
29*9566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
30*9566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
31c4762a1bSJed Brown 
32*9566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
33c4762a1bSJed Brown   am   = Iend - Istart;
34c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
35c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
36*9566063dSJacob Faibussowitsch     if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
37*9566063dSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
38*9566063dSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
39*9566063dSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
40*9566063dSJacob Faibussowitsch     v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
41c4762a1bSJed Brown   }
42*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
43*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Create a dense matrix B */
46*9566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&am,&an));
47*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
48*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M));
49*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATDENSE));
50*9566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(B,NULL));
51*9566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(B,NULL));
52*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
53*9566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B,r));
54*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
55*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
56*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
57c4762a1bSJed Brown 
58c20d7725SJed Brown   /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */
59*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
60*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATDENSE));
61*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N));
62*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
63*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
64*9566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
65*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
66*9566063dSJacob Faibussowitsch   PetscCall(MatNorm(C,NORM_INFINITY,&norm));
67*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
68c20d7725SJed Brown 
69c4762a1bSJed Brown   /* Test C = A*B (aij*dense) */
70*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C));
71*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C));
72c4762a1bSJed Brown 
73c20d7725SJed Brown   /* Test developer API */
74*9566063dSJacob Faibussowitsch   PetscCall(MatProductCreate(A,B,NULL,&D));
75*9566063dSJacob Faibussowitsch   PetscCall(MatProductSetType(D,MATPRODUCT_AB));
76*9566063dSJacob Faibussowitsch   PetscCall(MatProductSetAlgorithm(D,"default"));
77*9566063dSJacob Faibussowitsch   PetscCall(MatProductSetFill(D,fill));
78*9566063dSJacob Faibussowitsch   PetscCall(MatProductSetFromOptions(D));
79*9566063dSJacob Faibussowitsch   PetscCall(MatProductSymbolic(D));
80c4762a1bSJed Brown   for (i=0; i<2; i++) {
81*9566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(D));
82c4762a1bSJed Brown   }
83*9566063dSJacob Faibussowitsch   PetscCall(MatEqual(C,D,&equal));
8428b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D");
85*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
86c4762a1bSJed Brown 
87c20d7725SJed Brown   /* Test D = AT*B (transpose(aij)*dense) */
88*9566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(A,&AT));
89*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D));
90*9566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(AT,B,D,10,&equal));
9128b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)");
92*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
93*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AT));
94c20d7725SJed Brown 
95c4762a1bSJed Brown   /* Test D = C*A (dense*aij) */
96*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D));
97*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
98*9566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(C,A,D,10,&equal));
9928b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)");
100*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Test D = A*C (aij*dense) */
103*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D));
104*9566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D));
105*9566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A,C,D,10,&equal));
10628b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)");
107*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Test D = B*C (dense*dense) */
110*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
111c4762a1bSJed Brown   if (size == 1) {
112*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
113*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D));
114*9566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(B,C,D,10,&equal));
11528b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)");
116*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c20d7725SJed Brown   /* Test D = B*C^T (dense*dense) */
120*9566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D));
121*9566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D));
122*9566063dSJacob Faibussowitsch   PetscCall(MatMatTransposeMultEqual(B,C,D,10,&equal));
12328b400f6SJacob Faibussowitsch   PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)");
124*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
125c20d7725SJed Brown 
1263fbca975SStefano Zampini   /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */
1273fbca975SStefano Zampini   flg = PETSC_FALSE;
128*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg));
1293fbca975SStefano Zampini   if (flg) {
1303fbca975SStefano Zampini     /* user driver */
131*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B));
1323fbca975SStefano Zampini   } else {
1333fbca975SStefano Zampini     /* clear internal data structures related with previous products to avoid circular references */
134*9566063dSJacob Faibussowitsch     PetscCall(MatProductClear(A));
135*9566063dSJacob Faibussowitsch     PetscCall(MatProductClear(B));
136*9566063dSJacob Faibussowitsch     PetscCall(MatProductClear(C));
137*9566063dSJacob Faibussowitsch     PetscCall(MatProductCreateWithMat(A,C,NULL,B));
138*9566063dSJacob Faibussowitsch     PetscCall(MatProductSetType(B,MATPRODUCT_AB));
139*9566063dSJacob Faibussowitsch     PetscCall(MatProductSetFromOptions(B));
140*9566063dSJacob Faibussowitsch     PetscCall(MatProductSymbolic(B));
141*9566063dSJacob Faibussowitsch     PetscCall(MatProductNumeric(B));
1423fbca975SStefano Zampini   }
1433fbca975SStefano Zampini 
144*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
145*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
146*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
147*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
148b122ec5aSJacob Faibussowitsch   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /*TEST
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       args: -M 10 -N 10
155c4762a1bSJed Brown       output_file: output/ex109.out
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       suffix: 2
159c4762a1bSJed Brown       nsize: 3
160c4762a1bSJed Brown       output_file: output/ex109.out
161c4762a1bSJed Brown 
162c20d7725SJed Brown    test:
163c20d7725SJed Brown       suffix: 3
164c20d7725SJed Brown       nsize: 2
165c20d7725SJed Brown       args: -matmattransmult_mpidense_mpidense_via cyclic
166c20d7725SJed Brown       output_file: output/ex109.out
167c20d7725SJed Brown 
1683fbca975SStefano Zampini    test:
1693fbca975SStefano Zampini       suffix: 4
1703fbca975SStefano Zampini       args: -test_userAPI
1713fbca975SStefano Zampini       output_file: output/ex109.out
1723fbca975SStefano Zampini 
1733fbca975SStefano Zampini    test:
1743fbca975SStefano Zampini       suffix: 5
1753fbca975SStefano Zampini       nsize: 3
1763fbca975SStefano Zampini       args: -test_userAPI
1773fbca975SStefano Zampini       output_file: output/ex109.out
1783fbca975SStefano Zampini 
179c4762a1bSJed Brown TEST*/
180