xref: /petsc/src/mat/tests/ex70.c (revision 89669be4d29968dc8d4c19ce1b69194a6a561ea4)
1 #include <petscmat.h>
2 
3 static char help[] = "Tests MatMat operations with MAT_REUSE_MATRIX and already allocated dense result.\n\n";
4 
5 static PetscScalar MAGIC_NUMBER = 12345;
6 
7 static PetscErrorCode CheckLocal(Mat A, Mat B, PetscScalar *a, PetscScalar *b)
8 {
9   PetscBool      wA = PETSC_FALSE, wB = PETSC_FALSE;
10   PetscBool      wAv = PETSC_FALSE, wBv = PETSC_FALSE;
11   PetscInt       lda,i,j,m,n;
12 
13   PetscFunctionBegin;
14   if (a) {
15     const PetscScalar *Aa;
16     PetscCall(MatDenseGetArrayRead(A,&Aa));
17     wA   = (PetscBool)(a != Aa);
18     PetscCall(MatDenseGetLDA(A,&lda));
19     PetscCall(MatGetLocalSize(A,&m,&n));
20     for (j=0;j<n;j++) {
21       for (i=m;i<lda;i++) {
22         if (Aa[j*lda +i] != MAGIC_NUMBER) wAv = PETSC_TRUE;
23       }
24     }
25     PetscCall(MatDenseRestoreArrayRead(A,&Aa));
26   }
27   if (b) {
28     const PetscScalar *Bb;
29     PetscCall(MatDenseGetArrayRead(B,&Bb));
30     wB   = (PetscBool)(b != Bb);
31     PetscCall(MatDenseGetLDA(B,&lda));
32     PetscCall(MatGetLocalSize(B,&m,&n));
33     for (j=0;j<n;j++) {
34       for (i=m;i<lda;i++) {
35         if (Bb[j*lda +i] != MAGIC_NUMBER) wBv = PETSC_TRUE;
36       }
37     }
38     PetscCall(MatDenseRestoreArrayRead(B,&Bb));
39   }
40   PetscCheck(!wA && !wB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array in first Mat? %d, Wrong array in second Mat? %d",wA,wB);
41   PetscCheck(!wAv && !wBv,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong data in first Mat? %d, Wrong data in second Mat? %d",wAv,wBv);
42   PetscFunctionReturn(0);
43 }
44 
45 typedef struct {
46   Mat A;
47   Mat P;
48   Mat R;
49 } proj_data;
50 
51 PetscErrorCode proj_destroy(void *ctx)
52 {
53   proj_data      *userdata = (proj_data*)ctx;
54 
55   PetscFunctionBegin;
56   PetscCheck(userdata,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing userdata");
57   PetscCall(MatDestroy(&userdata->A));
58   PetscCall(MatDestroy(&userdata->P));
59   PetscCall(MatDestroy(&userdata->R));
60   PetscCall(PetscFree(userdata));
61   PetscFunctionReturn(0);
62 }
63 
64 PetscErrorCode proj_mult(Mat S, Vec X, Vec Y)
65 {
66   Mat            A,R,P;
67   Vec            Ax,Ay;
68   Vec            Px,Py;
69   proj_data      *userdata;
70 
71   PetscFunctionBegin;
72   PetscCall(MatShellGetContext(S,&userdata));
73   PetscCheck(userdata,PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing userdata");
74   A = userdata->A;
75   R = userdata->R;
76   P = userdata->P;
77   PetscCheck(A,PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing matrix");
78   PetscCheck(R || P,PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Missing projectors");
79   PetscCheck(!R || !P,PetscObjectComm((PetscObject)S),PETSC_ERR_PLIB,"Both projectors");
80   PetscCall(MatCreateVecs(A,&Ax,&Ay));
81   if (R) {
82     PetscCall(MatCreateVecs(R,&Py,&Px));
83   } else {
84     PetscCall(MatCreateVecs(P,&Px,&Py));
85   }
86   PetscCall(VecCopy(X,Px));
87   if (P) {
88     PetscCall(MatMult(P,Px,Py));
89   } else {
90     PetscCall(MatMultTranspose(R,Px,Py));
91   }
92   PetscCall(VecCopy(Py,Ax));
93   PetscCall(MatMult(A,Ax,Ay));
94   PetscCall(VecCopy(Ay,Py));
95   if (P) {
96     PetscCall(MatMultTranspose(P,Py,Px));
97   } else {
98     PetscCall(MatMult(R,Py,Px));
99   }
100   PetscCall(VecCopy(Px,Y));
101   PetscCall(VecDestroy(&Px));
102   PetscCall(VecDestroy(&Py));
103   PetscCall(VecDestroy(&Ax));
104   PetscCall(VecDestroy(&Ay));
105   PetscFunctionReturn(0);
106 }
107 
108 PetscErrorCode MyPtShellPMultSymbolic(Mat S, Mat P, Mat PtAP, void** ctx)
109 {
110   proj_data      *userdata;
111 
112   PetscFunctionBegin;
113   PetscCall(PetscNew(&userdata));
114   PetscCall(MatShellSetContext(PtAP,userdata));
115   *ctx = (void *)userdata;
116   PetscFunctionReturn(0);
117 }
118 
119 PetscErrorCode MyPtShellPMultNumeric(Mat S, Mat P, Mat PtAP, void *ctx)
120 {
121   Mat            A;
122   proj_data      *userdata = (proj_data*)ctx;
123 
124   PetscFunctionBegin;
125   PetscCall(MatShellGetContext(S,&A));
126   PetscCall(PetscObjectReference((PetscObject)A));
127   PetscCall(PetscObjectReference((PetscObject)P));
128   PetscCall(MatDestroy(&userdata->A));
129   PetscCall(MatDestroy(&userdata->P));
130   PetscCall(MatDestroy(&userdata->R));
131   userdata->A = A;
132   userdata->P = P;
133   PetscCall(MatShellSetOperation(PtAP,MATOP_MULT,(void (*)(void))proj_mult));
134   PetscCall(MatSetUp(PtAP));
135   PetscCall(MatAssemblyBegin(PtAP,MAT_FINAL_ASSEMBLY));
136   PetscCall(MatAssemblyEnd(PtAP,MAT_FINAL_ASSEMBLY));
137   PetscFunctionReturn(0);
138 }
139 
140 PetscErrorCode MyRShellRtMultSymbolic(Mat S, Mat R, Mat RARt, void **ctx)
141 {
142   proj_data      *userdata;
143 
144   PetscFunctionBegin;
145   PetscCall(PetscNew(&userdata));
146   PetscCall(MatShellSetContext(RARt,userdata));
147   *ctx = (void *)userdata;
148   PetscFunctionReturn(0);
149 }
150 
151 PetscErrorCode MyRShellRtMultNumeric(Mat S, Mat R, Mat RARt, void *ctx)
152 {
153   Mat            A;
154   proj_data      *userdata = (proj_data*)ctx;
155 
156   PetscFunctionBegin;
157   PetscCall(MatShellGetContext(S,&A));
158   PetscCall(PetscObjectReference((PetscObject)A));
159   PetscCall(PetscObjectReference((PetscObject)R));
160   PetscCall(MatDestroy(&userdata->A));
161   PetscCall(MatDestroy(&userdata->P));
162   PetscCall(MatDestroy(&userdata->R));
163   userdata->A = A;
164   userdata->R = R;
165   PetscCall(MatShellSetOperation(RARt,MATOP_MULT,(void (*)(void))proj_mult));
166   PetscCall(MatSetUp(RARt));
167   PetscCall(MatAssemblyBegin(RARt,MAT_FINAL_ASSEMBLY));
168   PetscCall(MatAssemblyEnd(RARt,MAT_FINAL_ASSEMBLY));
169   PetscFunctionReturn(0);
170 }
171 
172 PetscErrorCode MyMatShellMatMultNumeric(Mat S, Mat B, Mat C, void *ctx)
173 {
174   Mat            A;
175 
176   PetscFunctionBegin;
177   PetscCall(MatShellGetContext(S,&A));
178   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C));
179   PetscFunctionReturn(0);
180 }
181 
182 PetscErrorCode MyMatTransposeShellMatMultNumeric(Mat S, Mat B, Mat C, void *ctx)
183 {
184   Mat            A;
185 
186   PetscFunctionBegin;
187   PetscCall(MatShellGetContext(S,&A));
188   PetscCall(MatTransposeMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C));
189   PetscFunctionReturn(0);
190 }
191 
192 PetscErrorCode MyMatShellMatTransposeMultNumeric(Mat S, Mat B, Mat C, void *ctx)
193 {
194   Mat            A;
195 
196   PetscFunctionBegin;
197   PetscCall(MatShellGetContext(S,&A));
198   PetscCall(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C));
199   PetscFunctionReturn(0);
200 }
201 
202 int main(int argc,char **args)
203 {
204   Mat            X,B,A,Bt,T,T2,PtAP = NULL,RARt = NULL, R = NULL;
205   Vec            r,l,rs,ls;
206   PetscInt       m,n,k,M = 10,N = 10,K = 5, ldx = 3, ldb = 5, ldr = 4;
207   const char     *deft = MATAIJ;
208   char           mattype[256];
209   PetscBool      flg,symm = PETSC_FALSE,testtt = PETSC_TRUE, testnest = PETSC_TRUE, testtranspose = PETSC_TRUE, testcircular = PETSC_FALSE, local = PETSC_TRUE;
210   PetscBool      testhtranspose = PETSC_TRUE;
211   PetscBool      xgpu = PETSC_FALSE, bgpu = PETSC_FALSE, testshellops = PETSC_FALSE, testproj = PETSC_TRUE, testrart = PETSC_TRUE, testmatmatt = PETSC_TRUE, testmattmat = PETSC_TRUE;
212   PetscScalar    *dataX = NULL,*dataB = NULL, *dataR = NULL, *dataBt = NULL;
213   PetscScalar    *aX,*aB,*aBt;
214   PetscReal      err;
215 
216   PetscCall(PetscInitialize(&argc,&args,NULL,help));
217   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
218   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
219   PetscCall(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL));
220   PetscCall(PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL));
221   PetscCall(PetscOptionsGetBool(NULL,NULL,"-local",&local,NULL));
222   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ldx",&ldx,NULL));
223   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ldb",&ldb,NULL));
224   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ldr",&ldr,NULL));
225   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testtranspose",&testtranspose,NULL));
226   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testnest",&testnest,NULL));
227   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testtt",&testtt,NULL));
228   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testcircular",&testcircular,NULL));
229   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testshellops",&testshellops,NULL));
230   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testproj",&testproj,NULL));
231   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testrart",&testrart,NULL));
232   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testmatmatt",&testmatmatt,NULL));
233   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testmattmat",&testmattmat,NULL));
234   PetscCall(PetscOptionsGetBool(NULL,NULL,"-xgpu",&xgpu,NULL));
235   PetscCall(PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL));
236   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-magic_number",&MAGIC_NUMBER,NULL));
237   if (M != N) testproj = PETSC_FALSE;
238 
239   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
240   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
241   PetscCall(MatSetType(A,MATAIJ));
242   PetscCall(MatSetUp(A));
243   PetscCall(MatSetRandom(A,NULL));
244   if (M==N && symm) {
245     Mat AT;
246 
247     PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
248     PetscCall(MatAXPY(A,1.0,AT,DIFFERENT_NONZERO_PATTERN));
249     PetscCall(MatDestroy(&AT));
250     PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
251   }
252   PetscCall(MatViewFromOptions(A,NULL,"-A_init_view"));
253   PetscOptionsBegin(PETSC_COMM_WORLD,"","","");
254   PetscCall(PetscOptionsFList("-A_mat_type","Matrix type","MatSetType",MatList,deft,mattype,256,&flg));
255   PetscOptionsEnd();
256   if (flg) {
257     Mat A2;
258 
259     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2));
260     PetscCall(MatConvert(A,mattype,MAT_INPLACE_MATRIX,&A));
261     PetscCall(MatMultEqual(A,A2,10,&flg));
262     if (!flg) {
263       Mat AE,A2E;
264 
265       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with convert\n"));
266       PetscCall(MatComputeOperator(A,MATDENSE,&AE));
267       PetscCall(MatComputeOperator(A2,MATDENSE,&A2E));
268       PetscCall(MatView(AE,NULL));
269       PetscCall(MatView(A2E,NULL));
270       PetscCall(MatAXPY(A2E,-1.0,A,SAME_NONZERO_PATTERN));
271       PetscCall(MatView(A2E,NULL));
272       PetscCall(MatDestroy(&A2E));
273       PetscCall(MatDestroy(&AE));
274     }
275     PetscCall(MatDestroy(&A2));
276   }
277   PetscCall(MatViewFromOptions(A,NULL,"-A_view"));
278 
279   PetscCall(MatGetLocalSize(A,&m,&n));
280   if (local) {
281     PetscInt i;
282 
283     PetscCall(PetscMalloc1((m+ldx)*K,&dataX));
284     PetscCall(PetscMalloc1((n+ldb)*K,&dataB));
285     for (i=0;i<(m+ldx)*K;i++) dataX[i] = MAGIC_NUMBER;
286     for (i=0;i<(n+ldb)*K;i++) dataB[i] = MAGIC_NUMBER;
287   }
288   PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,K,dataB,&B));
289   PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,dataX,&X));
290   if (local) {
291     PetscCall(MatDenseSetLDA(X,m+ldx));
292     PetscCall(MatDenseSetLDA(B,n+ldb));
293   }
294   PetscCall(MatGetLocalSize(B,NULL,&k));
295   if (local) {
296     PetscInt i;
297 
298     PetscCall(PetscMalloc1((k+ldr)*N,&dataBt));
299     for (i=0;i<(k+ldr)*N;i++) dataBt[i] = MAGIC_NUMBER;
300   }
301   PetscCall(MatCreateDense(PETSC_COMM_WORLD,k,n,K,N,dataBt,&Bt));
302   if (local) {
303     PetscCall(MatDenseSetLDA(Bt,k+ldr));
304   }
305 
306   /* store pointer to dense data for testing */
307   PetscCall(MatDenseGetArrayRead(B,(const PetscScalar**)&dataB));
308   PetscCall(MatDenseGetArrayRead(X,(const PetscScalar**)&dataX));
309   PetscCall(MatDenseGetArrayRead(Bt,(const PetscScalar**)&dataBt));
310   aX   = dataX;
311   aB   = dataB;
312   aBt  = dataBt;
313   PetscCall(MatDenseRestoreArrayRead(Bt,(const PetscScalar**)&dataBt));
314   PetscCall(MatDenseRestoreArrayRead(B,(const PetscScalar**)&dataB));
315   PetscCall(MatDenseRestoreArrayRead(X,(const PetscScalar**)&dataX));
316   if (local) {
317     dataX  = aX;
318     dataB  = aB;
319     dataBt = aBt;
320   }
321 
322   PetscCall(MatSetRandom(X,NULL));
323   PetscCall(MatSetRandom(B,NULL));
324   PetscCall(MatSetRandom(Bt,NULL));
325   PetscCall(CheckLocal(X,NULL,aX,NULL));
326   PetscCall(CheckLocal(Bt,B,aBt,aB));
327 
328   /* convert to CUDA if needed */
329   if (bgpu) {
330     PetscCall(MatConvert(B,MATDENSECUDA,MAT_INPLACE_MATRIX,&B));
331     PetscCall(MatConvert(Bt,MATDENSECUDA,MAT_INPLACE_MATRIX,&Bt));
332   }
333   if (xgpu) {
334     PetscCall(MatConvert(X,MATDENSECUDA,MAT_INPLACE_MATRIX,&X));
335   }
336   PetscCall(CheckLocal(B,X,aB,aX));
337 
338   /* Test MatDenseGetSubMatrix */
339   {
340     Mat B2,T3,T4;
341 
342     PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2));
343     PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&T4));
344     PetscCall(MatSetRandom(T4,NULL));
345     PetscCall(MatAXPY(B2,1.0,T4,SAME_NONZERO_PATTERN));
346     PetscCall(MatDenseGetSubMatrix(B,PetscMin(1,K-1),PetscMin(2,K),&T));
347     PetscCall(MatDenseGetSubMatrix(T4,PetscMin(1,K-1),PetscMin(2,K),&T2));
348     PetscCall(MatDenseGetSubMatrix(B2,PetscMin(1,K-1),PetscMin(2,K),&T3));
349     PetscCall(MatAXPY(T,1.0,T2,SAME_NONZERO_PATTERN));
350     PetscCall(MatAXPY(T3,-1.0,T,SAME_NONZERO_PATTERN));
351     PetscCall(MatNorm(T3,NORM_FROBENIUS,&err));
352     if (err > PETSC_SMALL) {
353       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetSubMatrix\n"));
354       PetscCall(MatView(T3,NULL));
355     }
356     PetscCall(MatDenseRestoreSubMatrix(B,&T));
357     PetscCall(MatDenseRestoreSubMatrix(T4,&T2));
358     PetscCall(MatDenseRestoreSubMatrix(B2,&T3));
359     PetscCall(CheckLocal(B,NULL,aB,NULL));
360     PetscCall(MatDestroy(&B2));
361     PetscCall(MatDestroy(&T4));
362   }
363 
364   /* Test reusing a previously allocated dense buffer */
365   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
366   PetscCall(CheckLocal(B,X,aB,aX));
367   PetscCall(MatMatMultEqual(A,B,X,10,&flg));
368   if (!flg) {
369     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage\n"));
370     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
371     PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
372     PetscCall(MatView(T,NULL));
373     PetscCall(MatDestroy(&T));
374   }
375 
376   /* Test MatTransposeMat and MatMatTranspose */
377   if (testmattmat) {
378     PetscCall(MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
379     PetscCall(CheckLocal(B,X,aB,aX));
380     PetscCall(MatTransposeMatMultEqual(A,X,B,10,&flg));
381     if (!flg) {
382       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatTransposeMat)\n"));
383       PetscCall(MatTransposeMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B));
384       PetscCall(MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN));
385       PetscCall(MatView(T,NULL));
386       PetscCall(MatDestroy(&T));
387     }
388   }
389   if (testmatmatt) {
390     PetscCall(MatMatTransposeMult(A,Bt,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
391     PetscCall(CheckLocal(Bt,X,aBt,aX));
392     PetscCall(MatMatTransposeMultEqual(A,Bt,X,10,&flg));
393     if (!flg) {
394       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatMatTranspose)\n"));
395       PetscCall(MatMatTransposeMult(A,Bt,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
396       PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
397       PetscCall(MatView(T,NULL));
398       PetscCall(MatDestroy(&T));
399     }
400   }
401 
402   /* Test projection operations (PtAP and RARt) */
403   if (testproj) {
404     PetscCall(MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP));
405     PetscCall(MatPtAPMultEqual(A,B,PtAP,10,&flg));
406     if (!flg) {
407       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP\n"));
408       PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
409       PetscCall(MatTransposeMatMult(B,T,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T2));
410       PetscCall(MatAXPY(T2,-1.0,PtAP,SAME_NONZERO_PATTERN));
411       PetscCall(MatView(T2,NULL));
412       PetscCall(MatDestroy(&T2));
413       PetscCall(MatDestroy(&T));
414     }
415     PetscCall(PetscMalloc1((k+ldr)*M,&dataR));
416     PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,m,K,M,dataR,&R));
417     PetscCall(MatDenseSetLDA(R,k+ldr));
418     PetscCall(MatSetRandom(R,NULL));
419     if (testrart) { /* fails for AIJCUSPARSE because RA operation is not defined */
420       PetscCall(MatRARt(A,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RARt));
421       PetscCall(MatRARtMultEqual(A,R,RARt,10,&flg));
422       if (!flg) {
423         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with RARt\n"));
424         PetscCall(MatMatTransposeMult(A,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
425         PetscCall(MatMatMult(R,T,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T2));
426         PetscCall(MatAXPY(T2,-1.0,RARt,SAME_NONZERO_PATTERN));
427         PetscCall(MatView(T2,NULL));
428         PetscCall(MatDestroy(&T2));
429         PetscCall(MatDestroy(&T));
430       }
431     }
432   }
433 
434   /* Test MatDenseGetColumnVec and friends */
435   PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
436   PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
437   PetscCall(MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&T2));
438   for (k=0;k<K;k++) {
439     Vec Xv,Tv,T2v;
440 
441     PetscCall(MatDenseGetColumnVecRead(X,k,&Xv));
442     PetscCall(MatDenseGetColumnVec(T,k,&Tv));
443     PetscCall(MatDenseGetColumnVecWrite(T2,k,&T2v));
444     PetscCall(VecCopy(Xv,T2v));
445     PetscCall(VecAXPY(Tv,-1.,Xv));
446     PetscCall(MatDenseRestoreColumnVecRead(X,k,&Xv));
447     PetscCall(MatDenseRestoreColumnVec(T,k,&Tv));
448     PetscCall(MatDenseRestoreColumnVecWrite(T2,k,&T2v));
449   }
450   PetscCall(MatNorm(T,NORM_FROBENIUS,&err));
451   if (err > PETSC_SMALL) {
452     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVec\n"));
453     PetscCall(MatView(T,NULL));
454   }
455   PetscCall(MatAXPY(T2,-1.,X,SAME_NONZERO_PATTERN));
456   PetscCall(MatNorm(T2,NORM_FROBENIUS,&err));
457   if (err > PETSC_SMALL) {
458     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVecWrite\n"));
459     PetscCall(MatView(T2,NULL));
460   }
461   PetscCall(MatDestroy(&T));
462   PetscCall(MatDestroy(&T2));
463 
464   /* Test with MatShell */
465   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&T));
466   PetscCall(MatConvert(T,MATSHELL,MAT_INITIAL_MATRIX,&T2));
467   PetscCall(MatDestroy(&T));
468 
469   /* scale matrix */
470   PetscCall(MatScale(A,2.0));
471   PetscCall(MatScale(T2,2.0));
472   PetscCall(MatCreateVecs(A,&r,&l));
473   PetscCall(VecSetRandom(r,NULL));
474   PetscCall(VecSetRandom(l,NULL));
475   PetscCall(MatCreateVecs(T2,&rs,&ls));
476   PetscCall(VecCopy(r,rs));
477   PetscCall(VecCopy(l,ls));
478   if (testproj) {
479     PetscCall(MatDiagonalScale(A,r,r));
480     PetscCall(MatDiagonalScale(T2,rs,rs));
481   } else {
482     PetscCall(MatDiagonalScale(A,l,r));
483     PetscCall(MatDiagonalScale(T2,ls,rs));
484   }
485   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&T));
486   PetscCall(MatAXPY(A,4.5,T,SAME_NONZERO_PATTERN));
487   PetscCall(MatAXPY(T2,4.5,T,DIFFERENT_NONZERO_PATTERN));
488   PetscCall(MatMultEqual(T2,A,10,&flg));
489   if (!flg) {
490     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with MATSHELL (MatMult)\n"));
491   }
492   PetscCall(MatMultTransposeEqual(T2,A,10,&flg));
493   if (!flg) {
494     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with MATSHELL (MatMultTranspose)\n"));
495   }
496   PetscCall(MatDestroy(&T));
497   PetscCall(VecDestroy(&ls));
498   PetscCall(VecDestroy(&rs));
499   PetscCall(VecDestroy(&l));
500   PetscCall(VecDestroy(&r));
501 
502   /* recompute projections, test reusage */
503   if (PtAP) PetscCall(MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&PtAP));
504   if (RARt) PetscCall(MatRARt(A,R,MAT_REUSE_MATRIX,PETSC_DEFAULT,&RARt));
505   if (testshellops) { /* test callbacks for user defined MatProducts */
506     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_AB,NULL,MyMatShellMatMultNumeric,NULL,MATDENSE,MATDENSE));
507     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_AB,NULL,MyMatShellMatMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA));
508     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_AtB,NULL,MyMatTransposeShellMatMultNumeric,NULL,MATDENSE,MATDENSE));
509     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_AtB,NULL,MyMatTransposeShellMatMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA));
510     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_ABt,NULL,MyMatShellMatTransposeMultNumeric,NULL,MATDENSE,MATDENSE));
511     PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_ABt,NULL,MyMatShellMatTransposeMultNumeric,NULL,MATDENSECUDA,MATDENSECUDA));
512     if (testproj) {
513       PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_PtAP,MyPtShellPMultSymbolic,MyPtShellPMultNumeric,proj_destroy,MATDENSE,MATSHELL));
514       PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_PtAP,MyPtShellPMultSymbolic,MyPtShellPMultNumeric,proj_destroy,MATDENSECUDA,MATSHELL));
515       PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_RARt,MyRShellRtMultSymbolic,MyRShellRtMultNumeric,proj_destroy,MATDENSE,MATSHELL));
516       PetscCall(MatShellSetMatProductOperation(T2,MATPRODUCT_RARt,MyRShellRtMultSymbolic,MyRShellRtMultNumeric,proj_destroy,MATDENSECUDA,MATSHELL));
517     }
518   }
519   PetscCall(CheckLocal(B,X,aB,aX));
520   /* we either use the shell operations or the loop over columns code, applying the operator */
521   PetscCall(MatMatMult(T2,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
522   PetscCall(CheckLocal(B,X,aB,aX));
523   PetscCall(MatMatMultEqual(T2,B,X,10,&flg));
524   if (!flg) {
525     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MATSHELL)\n"));
526     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
527     PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
528     PetscCall(MatView(T,NULL));
529     PetscCall(MatDestroy(&T));
530   }
531   if (testproj) {
532     PetscCall(MatPtAPMultEqual(T2,B,PtAP,10,&flg));
533     if (!flg) {
534       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP (MATSHELL)\n"));
535     }
536     if (testshellops) { /* projections fail if the product operations are not specified */
537       PetscCall(MatPtAP(T2,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
538       PetscCall(MatPtAP(T2,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&T));
539       PetscCall(MatPtAPMultEqual(T2,B,T,10,&flg));
540       if (!flg) {
541         Mat TE;
542 
543         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with PtAP (MATSHELL user defined)\n"));
544         PetscCall(MatComputeOperator(T,MATDENSE,&TE));
545         PetscCall(MatView(TE,NULL));
546         PetscCall(MatView(PtAP,NULL));
547         PetscCall(MatAXPY(TE,-1.0,PtAP,SAME_NONZERO_PATTERN));
548         PetscCall(MatView(TE,NULL));
549         PetscCall(MatDestroy(&TE));
550       }
551       PetscCall(MatDestroy(&T));
552     }
553     if (RARt) {
554       PetscCall(MatRARtMultEqual(T2,R,RARt,10,&flg));
555       if (!flg) {
556         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with RARt (MATSHELL)\n"));
557       }
558     }
559     if (testshellops) {
560       PetscCall(MatRARt(T2,R,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
561       PetscCall(MatRARt(T2,R,MAT_REUSE_MATRIX,PETSC_DEFAULT,&T));
562       PetscCall(MatRARtMultEqual(T2,R,T,10,&flg));
563       if (!flg) {
564         Mat TE;
565 
566         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with RARt (MATSHELL user defined)\n"));
567         PetscCall(MatComputeOperator(T,MATDENSE,&TE));
568         PetscCall(MatView(TE,NULL));
569         if (RARt) {
570           PetscCall(MatView(RARt,NULL));
571           PetscCall(MatAXPY(TE,-1.0,RARt,SAME_NONZERO_PATTERN));
572           PetscCall(MatView(TE,NULL));
573         }
574         PetscCall(MatDestroy(&TE));
575       }
576       PetscCall(MatDestroy(&T));
577     }
578   }
579 
580   if (testmattmat) { /* we either use the shell operations or the loop over columns code applying the transposed operator */
581     PetscCall(MatTransposeMatMult(T2,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
582     PetscCall(CheckLocal(B,X,aB,aX));
583     PetscCall(MatTransposeMatMultEqual(T2,X,B,10,&flg));
584     if (!flg) {
585       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatTranspose, MATSHELL)\n"));
586       PetscCall(MatTransposeMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
587       PetscCall(MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN));
588       PetscCall(MatView(T,NULL));
589       PetscCall(MatDestroy(&T));
590     }
591   }
592   if (testmatmatt && testshellops) { /* only when shell operations are set */
593     PetscCall(MatMatTransposeMult(T2,Bt,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
594     PetscCall(CheckLocal(Bt,X,aBt,aX));
595     PetscCall(MatMatTransposeMultEqual(T2,Bt,X,10,&flg));
596     if (!flg) {
597       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatMatTranspose, MATSHELL)\n"));
598       PetscCall(MatMatTransposeMult(A,Bt,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
599       PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
600       PetscCall(MatView(T,NULL));
601       PetscCall(MatDestroy(&T));
602     }
603   }
604   PetscCall(MatDestroy(&T2));
605 
606   if (testnest) { /* test with MatNest */
607     Mat NA;
608 
609     PetscCall(MatCreateNest(PETSC_COMM_WORLD,1,NULL,1,NULL,&A,&NA));
610     PetscCall(MatViewFromOptions(NA,NULL,"-NA_view"));
611     PetscCall(MatMatMult(NA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
612     PetscCall(CheckLocal(B,X,aB,aX));
613     PetscCall(MatMatMultEqual(NA,B,X,10,&flg));
614     if (!flg) {
615       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with Nest\n"));
616       PetscCall(MatMatMult(NA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
617       PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
618       PetscCall(MatView(T,NULL));
619       PetscCall(MatDestroy(&T));
620     }
621     PetscCall(MatDestroy(&NA));
622   }
623 
624   if (testtranspose) { /* test with Transpose */
625     Mat TA;
626 
627     PetscCall(MatCreateTranspose(A,&TA));
628     PetscCall(MatMatMult(TA,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
629     PetscCall(CheckLocal(B,X,aB,aX));
630     PetscCall(MatMatMultEqual(TA,X,B,10,&flg));
631     if (!flg) {
632       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose\n"));
633       PetscCall(MatMatMult(TA,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
634       PetscCall(MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN));
635       PetscCall(MatView(T,NULL));
636       PetscCall(MatDestroy(&T));
637     }
638     PetscCall(MatDestroy(&TA));
639   }
640 
641   if (testhtranspose) { /* test with Hermitian Transpose */
642     Mat TA;
643 
644     PetscCall(MatCreateHermitianTranspose(A,&TA));
645     PetscCall(MatMatMult(TA,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
646     PetscCall(CheckLocal(B,X,aB,aX));
647     PetscCall(MatMatMultEqual(TA,X,B,10,&flg));
648     if (!flg) {
649       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose\n"));
650       PetscCall(MatMatMult(TA,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
651       PetscCall(MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN));
652       PetscCall(MatView(T,NULL));
653       PetscCall(MatDestroy(&T));
654     }
655     PetscCall(MatDestroy(&TA));
656   }
657 
658   if (testtt) { /* test with Transpose(Transpose) */
659     Mat TA, TTA;
660 
661     PetscCall(MatCreateTranspose(A,&TA));
662     PetscCall(MatCreateTranspose(TA,&TTA));
663     PetscCall(MatMatMult(TTA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
664     PetscCall(CheckLocal(B,X,aB,aX));
665     PetscCall(MatMatMultEqual(TTA,B,X,10,&flg));
666     if (!flg) {
667       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose(Transpose)\n"));
668       PetscCall(MatMatMult(TTA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
669       PetscCall(MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN));
670       PetscCall(MatView(T,NULL));
671       PetscCall(MatDestroy(&T));
672     }
673     PetscCall(MatDestroy(&TA));
674     PetscCall(MatDestroy(&TTA));
675   }
676 
677   if (testcircular) { /* test circular */
678     Mat AB;
679 
680     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AB));
681     PetscCall(MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X));
682     PetscCall(CheckLocal(B,X,aB,aX));
683     if (M == N && N == K) {
684       PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
685     } else {
686       PetscCall(MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B));
687     }
688     PetscCall(CheckLocal(B,X,aB,aX));
689     PetscCall(MatDestroy(&AB));
690   }
691 
692   /* Test by Pierre Jolivet */
693   {
694     Mat C,D,D2,AtA;
695     PetscCall(MatCreateNormal(A,&AtA));
696     PetscCall(MatDuplicate(X,MAT_DO_NOT_COPY_VALUES,&C));
697     PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&D));
698     PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&D2));
699     PetscCall(MatSetRandom(B,NULL));
700     PetscCall(MatSetRandom(C,NULL));
701     PetscCall(MatSetRandom(D,NULL));
702     PetscCall(MatSetRandom(D2,NULL));
703     PetscCall(MatProductCreateWithMat(A,B,NULL,C));
704     PetscCall(MatProductSetType(C,MATPRODUCT_AB));
705     PetscCall(MatProductSetFromOptions(C));
706     PetscCall(MatProductSymbolic(C));
707     PetscCall(MatProductCreateWithMat(A,C,NULL,D));
708     PetscCall(MatProductSetType(D, MATPRODUCT_AtB));
709     PetscCall(MatProductSetFromOptions(D));
710     PetscCall(MatProductSymbolic(D));
711     PetscCall(MatProductNumeric(C));
712     PetscCall(MatProductNumeric(D));
713     PetscCall(MatMatMultEqual(AtA,B,D,10,&flg));
714     if (!flg) {
715       PetscCall(MatMatMult(AtA,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T));
716       PetscCall(MatAXPY(T,-1.0,D,SAME_NONZERO_PATTERN));
717       PetscCall(MatView(T,NULL));
718       PetscCall(MatDestroy(&T));
719     }
720     PetscCall(MatDestroy(&C));
721     PetscCall(MatDestroy(&D));
722     PetscCall(MatDestroy(&D2));
723     PetscCall(MatDestroy(&AtA));
724   }
725 
726   PetscCall(MatDestroy(&X));
727   PetscCall(MatDestroy(&Bt));
728   PetscCall(MatDestroy(&B));
729   PetscCall(MatDestroy(&A));
730   PetscCall(MatDestroy(&R));
731   PetscCall(MatDestroy(&PtAP));
732   PetscCall(MatDestroy(&RARt));
733   PetscCall(PetscFree(dataX));
734   PetscCall(PetscFree(dataB));
735   PetscCall(PetscFree(dataR));
736   PetscCall(PetscFree(dataBt));
737   PetscCall(PetscFinalize());
738   return 0;
739 }
740 
741 /*TEST
742 
743   test:
744     suffix: 1
745     args: -local {{0 1}} -testshellops
746 
747   test:
748     output_file: output/ex70_1.out
749     requires: cuda
750     suffix: 1_cuda
751     args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{seqaijcusparse seqaij}} -testshellops {{0 1}}
752 
753   test:
754     output_file: output/ex70_1.out
755     nsize: 2
756     suffix: 1_par
757     args: -local {{0 1}} -testmatmatt 0
758 
759   test:
760     output_file: output/ex70_1.out
761     requires: cuda
762     nsize: 2
763     suffix: 1_par_cuda
764     args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{mpiaijcusparse mpiaij}} -testnest 0 -testmatmatt 0 -matmatmult_Bbn 3
765 
766   test:
767     output_file: output/ex70_1.out
768     suffix: 2
769     nsize: 1
770     args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}}
771 
772   testset:
773     requires: cuda
774     output_file: output/ex70_1.out
775     nsize: 1
776     args: -M 7 -N 9 -K 2 -local {{0 1}} -testnest 0 -A_mat_type {{seqdensecuda seqdense}} -xgpu {{0 1}} -bgpu {{0 1}}
777     test:
778       requires: !complex
779       suffix: 2_cuda_real
780     test:
781       # complex+single gives a little bigger error in the MatDenseGetColumnVec test
782       requires: complex !single
783       suffix: 2_cuda_complex
784 
785   test:
786     output_file: output/ex70_1.out
787     suffix: 2_par
788     nsize: 2
789     args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular -testmatmatt 0
790 
791   test:
792     requires: cuda
793     output_file: output/ex70_1.out
794     suffix: 2_par_cuda
795     nsize: 2
796     args: -M 11 -N 9 -K 1 -local {{0 1}} -testcircular 0 -A_mat_type mpiaijcusparse -xgpu -bgpu -testnest 0 -testmatmatt 0
797 
798   test:
799     output_file: output/ex70_1.out
800     suffix: 3
801     nsize: {{1 3}}
802     args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type sbaij -symm -testproj 0 -testmatmatt 0
803 
804   test:
805     output_file: output/ex70_1.out
806     suffix: 4
807     nsize: 1
808     args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular
809 
810   test:
811     output_file: output/ex70_1.out
812     suffix: 5
813     nsize: {{2 4}}
814     args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular -testmatmatt 0
815 
816   test:
817     output_file: output/ex70_1.out
818     suffix: 6
819     nsize: 1
820     args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular
821 
822   test:
823     output_file: output/ex70_1.out
824     suffix: 7
825     nsize: 1
826     args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type dense -testnest -testcircular
827 
828 TEST*/
829