xref: /petsc/src/mat/tests/ex241.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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   PetscErrorCode ierr;
50c7a4214aSPierre Jolivet 
51c7a4214aSPierre Jolivet   ierr = PetscInitialize(&argc,&argv,(char*)NULL,help);if (ierr) return ierr;
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL));
585f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
59c7a4214aSPierre Jolivet   M = size*m;
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(m*dim,&coords));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomGetValuesReal(rdm,m*dim,coords));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(M*dim,&gcoords));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(B,rdm));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(B,&begin,NULL));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(gcoords+begin*dim,coords,m*dim));
715f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,sym));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-A_view"));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(A,is,NULL));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDuplicate(is[0],is+1));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,1,is,2));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,2));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,1,is+1,1));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetBlockSize(is[1],&bs));
842c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bs != 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs);
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,1));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISEqual(is[0],is[1],&flg));
87*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets");
885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(is));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(is+1));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&right,&left));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(right,rdm));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,right,left));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHtoolGetPermutationSource(A,&iss));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHtoolGetPermutationTarget(A,&ist));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(left,&perm));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(left,perm));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPermute(perm,ist,PETSC_FALSE));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecPermute(right,iss,PETSC_FALSE));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHtoolUsePermutation(A,PETSC_FALSE));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,right,left));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(left,-1.0,perm));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(left,NORM_INFINITY,&norm));
1032c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHtoolUsePermutation(A,PETSC_TRUE));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&perm));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&left));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&right));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ist));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iss));
110c7a4214aSPierre 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 */
111c7a4214aSPierre Jolivet     PetscReal relative;
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(B));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(B,NULL,"-B_view"));
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P));
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(P,NORM_FROBENIUS,&relative));
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R));
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN));
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(R,NORM_INFINITY,&norm));
1242c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal(norm/relative) > epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon);
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&R));
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&P));
128c7a4214aSPierre Jolivet   } else {
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D));
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatViewFromOptions(D,NULL,"-D_view"));
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,D,10,&flg));
132*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTransposeEqual(A,D,10,&flg));
134*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultAddEqual(A,D,10,&flg));
136*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx");
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOwnershipRange(B,&begin,NULL));
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseGetArrayWrite(D,&ptr));
13998e73e17SPierre Jolivet     for (PetscInt i = begin; i < m+begin; ++i)
14098e73e17SPierre Jolivet       for (PetscInt j = 0; j < M; ++j) GenEntries(dim,1,1,&i,&j,ptr+i-begin+j*m,gcoords);
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseRestoreArrayWrite(D,&ptr));
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,D,10,&flg));
143*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx");
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
1465f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(AT,D,10,&flg));
147*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTranspose(A,MAT_REUSE_MATRIX,&AT));
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(AT,D,10,&flg));
150*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx");
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN));
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(D,NORM_INFINITY,&norm));
1532c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL);
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AT));
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P));
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
1595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatMultEqual(A,B,P,10,&flg));
160*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px");
1615f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMultEqual(A,B,P,10,&flg));
162*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx");
1635f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&B));
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&P));
165c7a4214aSPierre Jolivet     if (n) {
1665f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(n*dim,&scoords));
1675f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValuesReal(rdm,n*dim,scoords));
168c7a4214aSPierre Jolivet       N = n;
1695f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
1705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc1(N*dim,&gscoords));
1715f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD));
1725f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(gscoords+begin*dim,scoords,n*dim));
1735f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD));
17498e73e17SPierre Jolivet       kernel = GenEntriesRectangular;
175c7a4214aSPierre Jolivet       ctx[0] = gcoords;
176c7a4214aSPierre Jolivet       ctx[1] = gscoords;
1775f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R));
1785f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetFromOptions(R));
1795f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY));
1805f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY));
1815f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(R,NULL,"-R_view"));
1825f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D));
1835f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(D,NULL,"-D_view"));
1845f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultEqual(R,D,10,&flg));
185*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx");
1865f80ce2aSJacob Faibussowitsch       CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D));
1875f80ce2aSJacob Faibussowitsch       CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&RT));
1885f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultEqual(RT,D,10,&flg));
189*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
1905f80ce2aSJacob Faibussowitsch       CHKERRQ(MatTranspose(R,MAT_REUSE_MATRIX,&RT));
1915f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMultEqual(RT,D,10,&flg));
192*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx");
1935f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&RT));
1945f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&D));
1955f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B));
1965f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1975f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1985f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetRandom(B,rdm));
1995f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P));
2005f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
2015f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
2025f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMatMultEqual(R,B,P,10,&flg));
203*28b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px");
2045f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B));
2055f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&P));
2065f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(R,&right,&left));
2075f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSetRandom(right,rdm));
2085f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(R,right,left));
2095f80ce2aSJacob Faibussowitsch       CHKERRQ(MatHtoolGetPermutationSource(R,&iss));
2105f80ce2aSJacob Faibussowitsch       CHKERRQ(MatHtoolGetPermutationTarget(R,&ist));
2115f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(left,&perm));
2125f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(left,perm));
2135f80ce2aSJacob Faibussowitsch       CHKERRQ(VecPermute(perm,ist,PETSC_FALSE));
2145f80ce2aSJacob Faibussowitsch       CHKERRQ(VecPermute(right,iss,PETSC_FALSE));
2155f80ce2aSJacob Faibussowitsch       CHKERRQ(MatHtoolUsePermutation(R,PETSC_FALSE));
2165f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(R,right,left));
2175f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(left,-1.0,perm));
2185f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(left,NORM_INFINITY,&norm));
2192c71b3e2SJacob Faibussowitsch       PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL);
2205f80ce2aSJacob Faibussowitsch       CHKERRQ(MatHtoolUsePermutation(R,PETSC_TRUE));
2215f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&perm));
2225f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&left));
2235f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&right));
2245f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&ist));
2255f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&iss));
2265f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&R));
2275f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(gscoords));
2285f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(scoords));
229c7a4214aSPierre Jolivet     }
230c7a4214aSPierre Jolivet   }
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(gcoords));
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(coords));
235c7a4214aSPierre Jolivet   ierr = PetscFinalize();
236c7a4214aSPierre Jolivet   return ierr;
237c7a4214aSPierre Jolivet }
238c7a4214aSPierre Jolivet 
239c7a4214aSPierre Jolivet /*TEST
240c7a4214aSPierre Jolivet 
241c7a4214aSPierre Jolivet    build:
242c7a4214aSPierre Jolivet       requires: htool
243c7a4214aSPierre Jolivet 
244c7a4214aSPierre Jolivet    test:
245c7a4214aSPierre Jolivet       requires: htool
246c7a4214aSPierre Jolivet       suffix: 1
247c7a4214aSPierre Jolivet       nsize: 4
248c7a4214aSPierre Jolivet       args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output}
249c7a4214aSPierre Jolivet       output_file: output/ex101.out
250c7a4214aSPierre Jolivet 
251c7a4214aSPierre Jolivet    test:
252c7a4214aSPierre Jolivet       requires: htool
253c7a4214aSPierre Jolivet       suffix: 2
254c7a4214aSPierre Jolivet       nsize: 4
25598e73e17SPierre 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}
256c7a4214aSPierre Jolivet       output_file: output/ex101.out
257c7a4214aSPierre Jolivet 
258c7a4214aSPierre Jolivet TEST*/
259