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