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