xref: /petsc/src/mat/tests/ex66.c (revision 6ffe77eaecce1557e50d00ca5347a7f48e598865)
1 static char help[] = "Tests MATH2OPUS\n\n";
2 
3 #include <petscmat.h>
4 #include <petscsf.h>
5 
6 static PetscScalar GenEntry_Symm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx)
7 {
8     PetscInt  d;
9     PetscReal clength = sdim == 3 ? 0.2 : 0.1;
10     PetscReal dist, diff = 0.0;
11 
12     for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); }
13     dist = PetscSqrtReal(diff);
14     return PetscExpReal(-dist / clength);
15 }
16 
17 static PetscScalar GenEntry_Unsymm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx)
18 {
19     PetscInt  d;
20     PetscReal clength = sdim == 3 ? 0.2 : 0.1;
21     PetscReal dist, diff = 0.0, nx = 0.0, ny = 0.0;
22 
23     for (d = 0; d < sdim; d++) { nx += x[d]*x[d]; }
24     for (d = 0; d < sdim; d++) { ny += y[d]*y[d]; }
25     for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); }
26     dist = PetscSqrtReal(diff);
27     return nx > ny ? PetscExpReal(-dist / clength) : PetscExpReal(-dist / clength) + 1.;
28 }
29 
30 int main(int argc,char **argv)
31 {
32   Mat            A,B,C,D;
33   Vec            v,x,y,Ax,Ay,Bx,By;
34   PetscRandom    r;
35   PetscLayout    map;
36   PetscScalar    *Adata = NULL, *Cdata = NULL, scale = 1.0;
37   PetscReal      *coords,nA,nD,nB,err,nX,norms[3];
38   PetscInt       N, n = 64, dim = 1, i, j, nrhs = 11, lda = 0, ldc = 0, ldu = 0, nlr = 7, nt, ntrials = 2;
39   PetscMPIInt    size,rank;
40   PetscBool      testlayout = PETSC_FALSE, flg, symm = PETSC_FALSE, Asymm = PETSC_TRUE, kernel = PETSC_TRUE;
41   PetscBool      checkexpl = PETSC_FALSE, agpu = PETSC_FALSE, bgpu = PETSC_FALSE, cgpu = PETSC_FALSE, flgglob;
42   PetscBool      testtrans, testnorm, randommat = PETSC_TRUE, testorthog, testcompress, testhlru;
43   void           (*approxnormfunc)(void);
44   void           (*Anormfunc)(void);
45 
46 #if defined(PETSC_HAVE_MPI_INIT_THREAD)
47   PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE;
48 #endif
49   PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help));
50   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ng",&N,&flgglob));
51   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
52   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL));
53   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
54   PetscCall(PetscOptionsGetInt(NULL,NULL,"-lda",&lda,NULL));
55   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ldc",&ldc,NULL));
56   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nlr",&nlr,NULL));
57   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ldu",&ldu,NULL));
58   PetscCall(PetscOptionsGetInt(NULL,NULL,"-matmattrials",&ntrials,NULL));
59   PetscCall(PetscOptionsGetBool(NULL,NULL,"-randommat",&randommat,NULL));
60   if (!flgglob) PetscCall(PetscOptionsGetBool(NULL,NULL,"-testlayout",&testlayout,NULL));
61   PetscCall(PetscOptionsGetBool(NULL,NULL,"-Asymm",&Asymm,NULL));
62   PetscCall(PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL));
63   PetscCall(PetscOptionsGetBool(NULL,NULL,"-kernel",&kernel,NULL));
64   PetscCall(PetscOptionsGetBool(NULL,NULL,"-checkexpl",&checkexpl,NULL));
65   PetscCall(PetscOptionsGetBool(NULL,NULL,"-agpu",&agpu,NULL));
66   PetscCall(PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL));
67   PetscCall(PetscOptionsGetBool(NULL,NULL,"-cgpu",&cgpu,NULL));
68   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-scale",&scale,NULL));
69   if (!Asymm) symm = PETSC_FALSE;
70 
71   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
72 
73   /* Disable tests for unimplemented variants */
74   testtrans = (PetscBool)(size == 1 || symm);
75   testnorm = (PetscBool)(size == 1 || symm);
76   testorthog = (PetscBool)(size == 1 || symm);
77   testcompress = (PetscBool)(size == 1 || symm);
78   testhlru = (PetscBool)(size == 1);
79 
80   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
81   PetscCall(PetscLayoutCreate(PETSC_COMM_WORLD,&map));
82   if (testlayout) {
83     if (rank%2) n = PetscMax(2*n-5*rank,0);
84     else n = 2*n+rank;
85   }
86   if (!flgglob) {
87     PetscCall(PetscLayoutSetLocalSize(map,n));
88     PetscCall(PetscLayoutSetUp(map));
89     PetscCall(PetscLayoutGetSize(map,&N));
90   } else {
91     PetscCall(PetscLayoutSetSize(map,N));
92     PetscCall(PetscLayoutSetUp(map));
93     PetscCall(PetscLayoutGetLocalSize(map,&n));
94   }
95   PetscCall(PetscLayoutDestroy(&map));
96 
97   if (lda) {
98     PetscCall(PetscMalloc1(N*(n+lda),&Adata));
99   }
100   PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,n,N,N,Adata,&A));
101   PetscCall(MatDenseSetLDA(A,n+lda));
102 
103   /* Create random points; these are replicated in order to populate a dense matrix and to compare sequential and dense runs
104      The constructor for MATH2OPUS can take as input the distributed coordinates and replicates them internally in case
105      a kernel construction is requested */
106   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r));
107   PetscCall(PetscRandomSetFromOptions(r));
108   PetscCall(PetscRandomSetSeed(r,123456));
109   PetscCall(PetscRandomSeed(r));
110   PetscCall(PetscMalloc1(N*dim,&coords));
111   PetscCall(PetscRandomGetValuesReal(r,N*dim,coords));
112   PetscCall(PetscRandomDestroy(&r));
113 
114   if (kernel || !randommat) {
115     MatH2OpusKernel k = Asymm ? GenEntry_Symm : GenEntry_Unsymm;
116     PetscInt        ist,ien;
117 
118     PetscCall(MatGetOwnershipRange(A,&ist,&ien));
119     for (i = ist; i < ien; i++) {
120       for (j = 0; j < N; j++) {
121         PetscCall(MatSetValue(A,i,j,(*k)(dim,coords + i*dim,coords + j*dim,NULL),INSERT_VALUES));
122       }
123     }
124     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
125     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
126     if (kernel) {
127       PetscCall(MatCreateH2OpusFromKernel(PETSC_COMM_WORLD,n,n,N,N,dim,coords + ist*dim,PETSC_TRUE,k,NULL,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B));
128     } else {
129       PetscCall(MatCreateH2OpusFromMat(A,dim,coords + ist*dim,PETSC_TRUE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B));
130     }
131   } else {
132     PetscInt ist;
133 
134     PetscCall(MatGetOwnershipRange(A,&ist,NULL));
135     PetscCall(MatSetRandom(A,NULL));
136     if (Asymm) {
137       PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&B));
138       PetscCall(MatAXPY(A,1.0,B,SAME_NONZERO_PATTERN));
139       PetscCall(MatDestroy(&B));
140       PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
141     }
142     PetscCall(MatCreateH2OpusFromMat(A,dim,coords + ist*dim,PETSC_TRUE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B));
143   }
144   PetscCall(PetscFree(coords));
145   if (agpu) {
146     PetscCall(MatConvert(A,MATDENSECUDA,MAT_INPLACE_MATRIX,&A));
147   }
148   PetscCall(MatViewFromOptions(A,NULL,"-A_view"));
149 
150   PetscCall(MatSetOption(B,MAT_SYMMETRIC,symm));
151 
152   /* assemble the H-matrix */
153   PetscCall(MatBindToCPU(B,(PetscBool)!bgpu));
154   PetscCall(MatSetFromOptions(B));
155   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
156   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
157   PetscCall(MatViewFromOptions(B,NULL,"-B_view"));
158 
159   /* Test MatScale */
160   PetscCall(MatScale(A,scale));
161   PetscCall(MatScale(B,scale));
162 
163   /* Test MatMult */
164   PetscCall(MatCreateVecs(A,&Ax,&Ay));
165   PetscCall(MatCreateVecs(B,&Bx,&By));
166   PetscCall(VecSetRandom(Ax,NULL));
167   PetscCall(VecCopy(Ax,Bx));
168   PetscCall(MatMult(A,Ax,Ay));
169   PetscCall(MatMult(B,Bx,By));
170   PetscCall(VecViewFromOptions(Ay,NULL,"-mult_vec_view"));
171   PetscCall(VecViewFromOptions(By,NULL,"-mult_vec_view"));
172   PetscCall(VecNorm(Ay,NORM_INFINITY,&nX));
173   PetscCall(VecAXPY(Ay,-1.0,By));
174   PetscCall(VecViewFromOptions(Ay,NULL,"-mult_vec_view"));
175   PetscCall(VecNorm(Ay,NORM_INFINITY,&err));
176   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMult err %g\n",err/nX));
177   PetscCall(VecScale(By,-1.0));
178   PetscCall(MatMultAdd(B,Bx,By,By));
179   PetscCall(VecNorm(By,NORM_INFINITY,&err));
180   PetscCall(VecViewFromOptions(By,NULL,"-mult_vec_view"));
181   if (err > 10.*PETSC_SMALL) {
182     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMultAdd err %g\n",err));
183   }
184 
185   /* Test MatNorm */
186   PetscCall(MatNorm(A,NORM_INFINITY,&norms[0]));
187   PetscCall(MatNorm(A,NORM_1,&norms[1]));
188   norms[2] = -1.; /* NORM_2 not supported */
189   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"A Matrix norms:        infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]));
190   PetscCall(MatGetOperation(A,MATOP_NORM,&Anormfunc));
191   PetscCall(MatGetOperation(B,MATOP_NORM,&approxnormfunc));
192   PetscCall(MatSetOperation(A,MATOP_NORM,approxnormfunc));
193   PetscCall(MatNorm(A,NORM_INFINITY,&norms[0]));
194   PetscCall(MatNorm(A,NORM_1,&norms[1]));
195   PetscCall(MatNorm(A,NORM_2,&norms[2]));
196   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"A Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]));
197   if (testnorm) {
198     PetscCall(MatNorm(B,NORM_INFINITY,&norms[0]));
199     PetscCall(MatNorm(B,NORM_1,&norms[1]));
200     PetscCall(MatNorm(B,NORM_2,&norms[2]));
201   } else {
202     norms[0] = -1.;
203     norms[1] = -1.;
204     norms[2] = -1.;
205   }
206   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]));
207   PetscCall(MatSetOperation(A,MATOP_NORM,Anormfunc));
208 
209   /* Test MatDuplicate */
210   PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&D));
211   PetscCall(MatSetOption(D,MAT_SYMMETRIC,symm));
212   PetscCall(MatMultEqual(B,D,10,&flg));
213   if (!flg) {
214     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMult error after MatDuplicate\n"));
215   }
216   if (testtrans) {
217     PetscCall(MatMultTransposeEqual(B,D,10,&flg));
218     if (!flg) {
219       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMultTranspose error after MatDuplicate\n"));
220     }
221   }
222   PetscCall(MatDestroy(&D));
223 
224   if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */
225     PetscCall(VecSetRandom(Ay,NULL));
226     PetscCall(VecCopy(Ay,By));
227     PetscCall(MatMultTranspose(A,Ay,Ax));
228     PetscCall(MatMultTranspose(B,By,Bx));
229     PetscCall(VecViewFromOptions(Ax,NULL,"-multtrans_vec_view"));
230     PetscCall(VecViewFromOptions(Bx,NULL,"-multtrans_vec_view"));
231     PetscCall(VecNorm(Ax,NORM_INFINITY,&nX));
232     PetscCall(VecAXPY(Ax,-1.0,Bx));
233     PetscCall(VecViewFromOptions(Ax,NULL,"-multtrans_vec_view"));
234     PetscCall(VecNorm(Ax,NORM_INFINITY,&err));
235     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMultTranspose err %g\n",err/nX));
236     PetscCall(VecScale(Bx,-1.0));
237     PetscCall(MatMultTransposeAdd(B,By,Bx,Bx));
238     PetscCall(VecNorm(Bx,NORM_INFINITY,&err));
239     if (err > 10.*PETSC_SMALL) {
240       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMultTransposeAdd err %g\n",err));
241     }
242   }
243   PetscCall(VecDestroy(&Ax));
244   PetscCall(VecDestroy(&Ay));
245   PetscCall(VecDestroy(&Bx));
246   PetscCall(VecDestroy(&By));
247 
248   /* Test MatMatMult */
249   if (ldc) {
250     PetscCall(PetscMalloc1(nrhs*(n+ldc),&Cdata));
251   }
252   PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,nrhs,Cdata,&C));
253   PetscCall(MatDenseSetLDA(C,n+ldc));
254   PetscCall(MatSetRandom(C,NULL));
255   if (cgpu) {
256     PetscCall(MatConvert(C,MATDENSECUDA,MAT_INPLACE_MATRIX,&C));
257   }
258   for (nt = 0; nt < ntrials; nt++) {
259     PetscCall(MatMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D));
260     PetscCall(MatViewFromOptions(D,NULL,"-bc_view"));
261     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,""));
262     if (flg) {
263       PetscCall(MatCreateVecs(B,&x,&y));
264       PetscCall(MatCreateVecs(D,NULL,&v));
265       for (i = 0; i < nrhs; i++) {
266         PetscCall(MatGetColumnVector(D,v,i));
267         PetscCall(MatGetColumnVector(C,x,i));
268         PetscCall(MatMult(B,x,y));
269         PetscCall(VecAXPY(y,-1.0,v));
270         PetscCall(VecNorm(y,NORM_INFINITY,&err));
271         if (err > 10.*PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMat err %" PetscInt_FMT " %g\n",i,err));
272       }
273       PetscCall(VecDestroy(&y));
274       PetscCall(VecDestroy(&x));
275       PetscCall(VecDestroy(&v));
276     }
277   }
278   PetscCall(MatDestroy(&D));
279 
280   /* Test MatTransposeMatMult */
281   if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */
282     for (nt = 0; nt < ntrials; nt++) {
283       PetscCall(MatTransposeMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D));
284       PetscCall(MatViewFromOptions(D,NULL,"-btc_view"));
285       PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,""));
286       if (flg) {
287         PetscCall(MatCreateVecs(B,&y,&x));
288         PetscCall(MatCreateVecs(D,NULL,&v));
289         for (i = 0; i < nrhs; i++) {
290           PetscCall(MatGetColumnVector(D,v,i));
291           PetscCall(MatGetColumnVector(C,x,i));
292           PetscCall(MatMultTranspose(B,x,y));
293           PetscCall(VecAXPY(y,-1.0,v));
294           PetscCall(VecNorm(y,NORM_INFINITY,&err));
295           if (err > 10.*PETSC_SMALL) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatTransMat err %" PetscInt_FMT " %g\n",i,err));
296         }
297         PetscCall(VecDestroy(&y));
298         PetscCall(VecDestroy(&x));
299         PetscCall(VecDestroy(&v));
300       }
301     }
302     PetscCall(MatDestroy(&D));
303   }
304 
305   /* Test basis orthogonalization */
306   if (testorthog) {
307     PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&D));
308     PetscCall(MatSetOption(D,MAT_SYMMETRIC,symm));
309     PetscCall(MatH2OpusOrthogonalize(D));
310     PetscCall(MatMultEqual(B,D,10,&flg));
311     if (!flg) {
312       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMult error after basis ortogonalization\n"));
313     }
314     PetscCall(MatDestroy(&D));
315   }
316 
317   /* Test matrix compression */
318   if (testcompress) {
319     PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&D));
320     PetscCall(MatSetOption(D,MAT_SYMMETRIC,symm));
321     PetscCall(MatH2OpusCompress(D,PETSC_SMALL));
322     PetscCall(MatDestroy(&D));
323   }
324 
325   /* Test low-rank update */
326   if (testhlru) {
327     Mat         U, V;
328     PetscScalar *Udata = NULL, *Vdata = NULL;
329 
330     if (ldu) {
331       PetscCall(PetscMalloc1(nlr*(n+ldu),&Udata));
332       PetscCall(PetscMalloc1(nlr*(n+ldu+2),&Vdata));
333     }
334     PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&D));
335     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)D),n,PETSC_DECIDE,N,nlr,Udata,&U));
336     PetscCall(MatDenseSetLDA(U,n+ldu));
337     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)D),n,PETSC_DECIDE,N,nlr,Vdata,&V));
338     if (ldu) PetscCall(MatDenseSetLDA(V,n+ldu+2));
339     PetscCall(MatSetRandom(U,NULL));
340     PetscCall(MatSetRandom(V,NULL));
341     PetscCall(MatH2OpusLowRankUpdate(D,U,V,0.5));
342     PetscCall(MatH2OpusLowRankUpdate(D,U,V,-0.5));
343     PetscCall(MatMultEqual(B,D,10,&flg));
344     if (!flg) {
345       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatMult error after low-rank update\n"));
346     }
347     PetscCall(MatDestroy(&D));
348     PetscCall(MatDestroy(&U));
349     PetscCall(PetscFree(Udata));
350     PetscCall(MatDestroy(&V));
351     PetscCall(PetscFree(Vdata));
352   }
353 
354   /* check explicit operator */
355   if (checkexpl) {
356     Mat Be, Bet;
357 
358     PetscCall(MatComputeOperator(B,MATDENSE,&D));
359     PetscCall(MatDuplicate(D,MAT_COPY_VALUES,&Be));
360     PetscCall(MatNorm(D,NORM_FROBENIUS,&nB));
361     PetscCall(MatViewFromOptions(D,NULL,"-expl_view"));
362     PetscCall(MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN));
363     PetscCall(MatViewFromOptions(D,NULL,"-diff_view"));
364     PetscCall(MatNorm(D,NORM_FROBENIUS,&nD));
365     PetscCall(MatNorm(A,NORM_FROBENIUS,&nA));
366     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Approximation error %g (%g / %g, %g)\n",nD/nA,nD,nA,nB));
367     PetscCall(MatDestroy(&D));
368 
369     if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */
370       PetscCall(MatTranspose(A,MAT_INPLACE_MATRIX,&A));
371       PetscCall(MatComputeOperatorTranspose(B,MATDENSE,&D));
372       PetscCall(MatDuplicate(D,MAT_COPY_VALUES,&Bet));
373       PetscCall(MatNorm(D,NORM_FROBENIUS,&nB));
374       PetscCall(MatViewFromOptions(D,NULL,"-expl_trans_view"));
375       PetscCall(MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN));
376       PetscCall(MatViewFromOptions(D,NULL,"-diff_trans_view"));
377       PetscCall(MatNorm(D,NORM_FROBENIUS,&nD));
378       PetscCall(MatNorm(A,NORM_FROBENIUS,&nA));
379       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Approximation error transpose %g (%g / %g, %g)\n",nD/nA,nD,nA,nB));
380       PetscCall(MatDestroy(&D));
381 
382       PetscCall(MatTranspose(Bet,MAT_INPLACE_MATRIX,&Bet));
383       PetscCall(MatAXPY(Be,-1.0,Bet,SAME_NONZERO_PATTERN));
384       PetscCall(MatViewFromOptions(Be,NULL,"-diff_expl_view"));
385       PetscCall(MatNorm(Be,NORM_FROBENIUS,&nB));
386       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Approximation error B - (B^T)^T %g\n",nB));
387       PetscCall(MatDestroy(&Be));
388       PetscCall(MatDestroy(&Bet));
389     }
390   }
391   PetscCall(MatDestroy(&A));
392   PetscCall(MatDestroy(&B));
393   PetscCall(MatDestroy(&C));
394   PetscCall(PetscFree(Cdata));
395   PetscCall(PetscFree(Adata));
396   PetscCall(PetscFinalize());
397   return 0;
398 }
399 
400 /*TEST
401 
402    build:
403      requires: h2opus
404 
405 #tests from kernel
406    test:
407      requires: h2opus
408      nsize: 1
409      suffix: 1
410      args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 0
411 
412    test:
413      requires: h2opus
414      nsize: 1
415      suffix: 1_ld
416      output_file: output/ex66_1.out
417      args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -ldu 17 -symm 0 -checkexpl -bgpu 0
418 
419    test:
420      requires: h2opus cuda
421      nsize: 1
422      suffix: 1_cuda
423      output_file: output/ex66_1.out
424      args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 1
425 
426    test:
427      requires: h2opus cuda
428      nsize: 1
429      suffix: 1_cuda_ld
430      output_file: output/ex66_1.out
431      args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -ldu 17 -symm 0 -checkexpl -bgpu 1
432 
433    test:
434      requires: h2opus
435      nsize: {{2 3}}
436      suffix: 1_par
437      args: -n 64 -symm -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu 0 -cgpu 0
438 
439    test:
440      requires: h2opus cuda
441      nsize: {{2 3}}
442      suffix: 1_par_cuda
443      args: -n 64 -symm -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu {{0 1}} -cgpu {{0 1}}
444      output_file: output/ex66_1_par.out
445 
446 #tests from matrix sampling (parallel or unsymmetric not supported)
447    test:
448      requires: h2opus
449      nsize: 1
450      suffix: 2
451      args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu 0
452 
453    test:
454      requires: h2opus cuda
455      nsize: 1
456      suffix: 2_cuda
457      output_file: output/ex66_2.out
458      args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu {{0 1}} -agpu {{0 1}}
459 
460 #tests view operation
461    test:
462      requires: h2opus !cuda
463      filter: grep -v " MPI process" | grep -v "\[" | grep -v "\]"
464      nsize: {{1 2 3}}
465      suffix: view
466      args: -ng 64 -kernel 1 -dim 2 -symm 1 -checkexpl -B_view -mat_h2opus_leafsize 17 -mat_h2opus_normsamples 13 -mat_h2opus_indexmap_view ::ascii_matlab -mat_approximate_norm_samples 2 -mat_h2opus_normsamples 2
467 
468    test:
469      requires: h2opus cuda
470      filter: grep -v " MPI process" | grep -v "\[" | grep -v "\]"
471      nsize: {{1 2 3}}
472      suffix: view_cuda
473      args: -ng 64 -kernel 1 -dim 2 -symm 1 -checkexpl -bgpu -B_view -mat_h2opus_leafsize 17 -mat_h2opus_normsamples 13 -mat_h2opus_indexmap_view ::ascii_matlab -mat_approximate_norm_samples 2 -mat_h2opus_normsamples 2
474 
475 TEST*/
476