xref: /petsc/src/mat/tests/ex241.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1c7a4214aSPierre Jolivet static char help[] = "Tests MATHTOOL\n\n";
2c7a4214aSPierre Jolivet 
3c7a4214aSPierre Jolivet #include <petscmat.h>
4c7a4214aSPierre Jolivet 
598e73e17SPierre Jolivet static PetscErrorCode GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx)
6c7a4214aSPierre Jolivet {
798e73e17SPierre Jolivet   PetscInt  d,j,k;
8c7a4214aSPierre Jolivet   PetscReal diff = 0.0,*coords = (PetscReal*)(ctx);
9c7a4214aSPierre Jolivet 
10c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
1198e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
1298e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
1398e73e17SPierre Jolivet       diff = 0.0;
1498e73e17SPierre Jolivet       for (d = 0; d < sdim; d++) diff += (coords[J[j]*sdim+d] - coords[K[k]*sdim+d]) * (coords[J[j]*sdim+d] - coords[K[k]*sdim+d]);
1598e73e17SPierre Jolivet       ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
1698e73e17SPierre Jolivet     }
1798e73e17SPierre Jolivet   }
1898e73e17SPierre Jolivet   PetscFunctionReturn(0);
19c7a4214aSPierre Jolivet }
20c7a4214aSPierre Jolivet 
2198e73e17SPierre Jolivet static PetscErrorCode GenEntriesRectangular(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx)
22c7a4214aSPierre Jolivet {
2398e73e17SPierre Jolivet   PetscInt  d,j,k;
24c7a4214aSPierre Jolivet   PetscReal diff = 0.0,**coords = (PetscReal**)(ctx);
25c7a4214aSPierre Jolivet 
26c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
2798e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
2898e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
2998e73e17SPierre Jolivet       diff = 0.0;
3098e73e17SPierre Jolivet       for (d = 0; d < sdim; d++) diff += (coords[0][J[j]*sdim+d] - coords[1][K[k]*sdim+d]) * (coords[0][J[j]*sdim+d] - coords[1][K[k]*sdim+d]);
3198e73e17SPierre Jolivet       ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff));
3298e73e17SPierre Jolivet     }
3398e73e17SPierre Jolivet   }
3498e73e17SPierre Jolivet   PetscFunctionReturn(0);
35c7a4214aSPierre Jolivet }
36c7a4214aSPierre Jolivet 
37c7a4214aSPierre Jolivet int main(int argc,char **argv)
38c7a4214aSPierre Jolivet {
39c7a4214aSPierre Jolivet   Mat            A,AT,D,B,P,R,RT;
408f308287SPierre Jolivet   PetscInt       m = 100,dim = 3,M,K = 10,begin,n = 0,N,bs;
41c7a4214aSPierre Jolivet   PetscMPIInt    size;
42c7a4214aSPierre Jolivet   PetscScalar    *ptr;
43c7a4214aSPierre Jolivet   PetscReal      *coords,*gcoords,*scoords,*gscoords,*(ctx[2]),norm,epsilon;
4498e73e17SPierre Jolivet   MatHtoolKernel kernel = GenEntries;
45c7a4214aSPierre Jolivet   PetscBool      flg,sym = PETSC_FALSE;
46c7a4214aSPierre Jolivet   PetscRandom    rdm;
478f308287SPierre Jolivet   IS             iss,ist,is[2];
4898e73e17SPierre Jolivet   Vec            right,left,perm;
49c7a4214aSPierre Jolivet 
509566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)NULL,help));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL));
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL));
579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
58c7a4214aSPierre Jolivet   M = size*m;
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m*dim,&coords));
619566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
629566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal(rdm,m*dim,coords));
639566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(M*dim,&gcoords));
649566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B));
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B,rdm));
689566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B,&begin,NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(gcoords+begin*dim,coords,m*dim));
701c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD));
719566063dSJacob Faibussowitsch   PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A));
729566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_SYMMETRIC,sym));
739566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
769566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A,NULL,"-A_view"));
779566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A,is,NULL));
789566063dSJacob Faibussowitsch   PetscCall(ISDuplicate(is[0],is+1));
799566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A,1,is,2));
809566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A,2));
819566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A,1,is+1,1));
829566063dSJacob Faibussowitsch   PetscCall(ISGetBlockSize(is[1],&bs));
83*08401ef6SPierre Jolivet   PetscCheck(bs == 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs);
849566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A,1));
859566063dSJacob Faibussowitsch   PetscCall(ISEqual(is[0],is[1],&flg));
8628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets");
879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(is));
889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(is+1));
899566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&right,&left));
909566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(right,rdm));
919566063dSJacob Faibussowitsch   PetscCall(MatMult(A,right,left));
929566063dSJacob Faibussowitsch   PetscCall(MatHtoolGetPermutationSource(A,&iss));
939566063dSJacob Faibussowitsch   PetscCall(MatHtoolGetPermutationTarget(A,&ist));
949566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(left,&perm));
959566063dSJacob Faibussowitsch   PetscCall(VecCopy(left,perm));
969566063dSJacob Faibussowitsch   PetscCall(VecPermute(perm,ist,PETSC_FALSE));
979566063dSJacob Faibussowitsch   PetscCall(VecPermute(right,iss,PETSC_FALSE));
989566063dSJacob Faibussowitsch   PetscCall(MatHtoolUsePermutation(A,PETSC_FALSE));
999566063dSJacob Faibussowitsch   PetscCall(MatMult(A,right,left));
1009566063dSJacob Faibussowitsch   PetscCall(VecAXPY(left,-1.0,perm));
1019566063dSJacob Faibussowitsch   PetscCall(VecNorm(left,NORM_INFINITY,&norm));
102*08401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
1039566063dSJacob Faibussowitsch   PetscCall(MatHtoolUsePermutation(A,PETSC_TRUE));
1049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&perm));
1059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
1069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
1079566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ist));
1089566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iss));
109c7a4214aSPierre Jolivet   if (PetscAbsReal(epsilon) >= PETSC_SMALL) { /* when there is compression, it is more difficult to check against MATDENSE, so just compare symmetric and nonsymmetric assemblies */
110c7a4214aSPierre Jolivet     PetscReal relative;
1119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
1129566063dSJacob Faibussowitsch     PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B));
1139566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym));
1149566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(B));
1159566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1169566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1179566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(B,NULL,"-B_view"));
1189566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P));
1199566063dSJacob Faibussowitsch     PetscCall(MatNorm(P,NORM_FROBENIUS,&relative));
1209566063dSJacob Faibussowitsch     PetscCall(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R));
1219566063dSJacob Faibussowitsch     PetscCall(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN));
1229566063dSJacob Faibussowitsch     PetscCall(MatNorm(R,NORM_INFINITY,&norm));
123*08401ef6SPierre Jolivet     PetscCheck(PetscAbsReal(norm/relative) <= epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
1249566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
1259566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&R));
1269566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&P));
127c7a4214aSPierre Jolivet   } else {
1289566063dSJacob Faibussowitsch     PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D));
1299566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(D,NULL,"-D_view"));
1309566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,D,10,&flg));
13128b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
1329566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeEqual(A,D,10,&flg));
13328b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1349566063dSJacob Faibussowitsch     PetscCall(MatMultAddEqual(A,D,10,&flg));
13528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx");
1369566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(B,&begin,NULL));
1379566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(D,&ptr));
13898e73e17SPierre Jolivet     for (PetscInt i = begin; i < m+begin; ++i)
13998e73e17SPierre Jolivet       for (PetscInt j = 0; j < M; ++j) GenEntries(dim,1,1,&i,&j,ptr+i-begin+j*m,gcoords);
1409566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(D,&ptr));
1419566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A,D,10,&flg));
14228b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
1439566063dSJacob Faibussowitsch     PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D));
1449566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
1459566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(AT,D,10,&flg));
14628b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1479566063dSJacob Faibussowitsch     PetscCall(MatTranspose(A,MAT_REUSE_MATRIX,&AT));
1489566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(AT,D,10,&flg));
14928b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1509566063dSJacob Faibussowitsch     PetscCall(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN));
1519566063dSJacob Faibussowitsch     PetscCall(MatNorm(D,NORM_INFINITY,&norm));
152*08401ef6SPierre Jolivet     PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL);
1539566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
1549566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1559566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P));
1569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
1579566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
1589566063dSJacob Faibussowitsch     PetscCall(MatMatMultEqual(A,B,P,10,&flg));
15928b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px");
1609566063dSJacob Faibussowitsch     PetscCall(MatTransposeMatMultEqual(A,B,P,10,&flg));
16128b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx");
1629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
1639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&P));
164c7a4214aSPierre Jolivet     if (n) {
1659566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n*dim,&scoords));
1669566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValuesReal(rdm,n*dim,scoords));
167c7a4214aSPierre Jolivet       N = n;
1681c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
1699566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(N*dim,&gscoords));
1709566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
1719566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(gscoords+begin*dim,scoords,n*dim));
1721c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD));
17398e73e17SPierre Jolivet       kernel = GenEntriesRectangular;
174c7a4214aSPierre Jolivet       ctx[0] = gcoords;
175c7a4214aSPierre Jolivet       ctx[1] = gscoords;
1769566063dSJacob Faibussowitsch       PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R));
1779566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(R));
1789566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY));
1799566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY));
1809566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(R,NULL,"-R_view"));
1819566063dSJacob Faibussowitsch       PetscCall(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D));
1829566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(D,NULL,"-D_view"));
1839566063dSJacob Faibussowitsch       PetscCall(MatMultEqual(R,D,10,&flg));
18428b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx");
1859566063dSJacob Faibussowitsch       PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D));
1869566063dSJacob Faibussowitsch       PetscCall(MatTranspose(R,MAT_INITIAL_MATRIX,&RT));
1879566063dSJacob Faibussowitsch       PetscCall(MatMultEqual(RT,D,10,&flg));
18828b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
1899566063dSJacob Faibussowitsch       PetscCall(MatTranspose(R,MAT_REUSE_MATRIX,&RT));
1909566063dSJacob Faibussowitsch       PetscCall(MatMultEqual(RT,D,10,&flg));
19128b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
1929566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RT));
1939566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&D));
1949566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B));
1959566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1969566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1979566063dSJacob Faibussowitsch       PetscCall(MatSetRandom(B,rdm));
1989566063dSJacob Faibussowitsch       PetscCall(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P));
1999566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
2009566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
2019566063dSJacob Faibussowitsch       PetscCall(MatMatMultEqual(R,B,P,10,&flg));
20228b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px");
2039566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
2049566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&P));
2059566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(R,&right,&left));
2069566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(right,rdm));
2079566063dSJacob Faibussowitsch       PetscCall(MatMult(R,right,left));
2089566063dSJacob Faibussowitsch       PetscCall(MatHtoolGetPermutationSource(R,&iss));
2099566063dSJacob Faibussowitsch       PetscCall(MatHtoolGetPermutationTarget(R,&ist));
2109566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(left,&perm));
2119566063dSJacob Faibussowitsch       PetscCall(VecCopy(left,perm));
2129566063dSJacob Faibussowitsch       PetscCall(VecPermute(perm,ist,PETSC_FALSE));
2139566063dSJacob Faibussowitsch       PetscCall(VecPermute(right,iss,PETSC_FALSE));
2149566063dSJacob Faibussowitsch       PetscCall(MatHtoolUsePermutation(R,PETSC_FALSE));
2159566063dSJacob Faibussowitsch       PetscCall(MatMult(R,right,left));
2169566063dSJacob Faibussowitsch       PetscCall(VecAXPY(left,-1.0,perm));
2179566063dSJacob Faibussowitsch       PetscCall(VecNorm(left,NORM_INFINITY,&norm));
218*08401ef6SPierre Jolivet       PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
2199566063dSJacob Faibussowitsch       PetscCall(MatHtoolUsePermutation(R,PETSC_TRUE));
2209566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&perm));
2219566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&left));
2229566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&right));
2239566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&ist));
2249566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&iss));
2259566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&R));
2269566063dSJacob Faibussowitsch       PetscCall(PetscFree(gscoords));
2279566063dSJacob Faibussowitsch       PetscCall(PetscFree(scoords));
228c7a4214aSPierre Jolivet     }
229c7a4214aSPierre Jolivet   }
2309566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
2319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree(gcoords));
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree(coords));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
235b122ec5aSJacob Faibussowitsch   return 0;
236c7a4214aSPierre Jolivet }
237c7a4214aSPierre Jolivet 
238c7a4214aSPierre Jolivet /*TEST
239c7a4214aSPierre Jolivet 
240c7a4214aSPierre Jolivet    build:
241c7a4214aSPierre Jolivet       requires: htool
242c7a4214aSPierre Jolivet 
243c7a4214aSPierre Jolivet    test:
244c7a4214aSPierre Jolivet       requires: htool
245c7a4214aSPierre Jolivet       suffix: 1
246c7a4214aSPierre Jolivet       nsize: 4
247c7a4214aSPierre Jolivet       args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output}
248c7a4214aSPierre Jolivet       output_file: output/ex101.out
249c7a4214aSPierre Jolivet 
250c7a4214aSPierre Jolivet    test:
251c7a4214aSPierre Jolivet       requires: htool
252c7a4214aSPierre Jolivet       suffix: 2
253c7a4214aSPierre Jolivet       nsize: 4
25498e73e17SPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2 -mat_htool_compressor {{sympartialACA fullACA SVD}shared output} -mat_htool_clustering {{PCARegular PCAGeometric BoundingBox1Regular BoundingBox1Geometric}shared output}
255c7a4214aSPierre Jolivet       output_file: output/ex101.out
256c7a4214aSPierre Jolivet 
257c7a4214aSPierre Jolivet TEST*/
258