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 50*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)NULL,help)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL)); 575f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 58c7a4214aSPierre Jolivet M = size*m; 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*dim,&coords)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValuesReal(rdm,m*dim,coords)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(M*dim,&gcoords)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&begin,NULL)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(gcoords+begin*dim,coords,m*dim)); 705f80ce2aSJacob Faibussowitsch CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,sym)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-A_view")); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(A,is,NULL)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(ISDuplicate(is[0],is+1)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,1,is,2)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,2)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,1,is+1,1)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetBlockSize(is[1],&bs)); 832c71b3e2SJacob Faibussowitsch PetscCheckFalse(bs != 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,1)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(ISEqual(is[0],is[1],&flg)); 8628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets"); 875f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(is)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(is+1)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&right,&left)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(right,rdm)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,right,left)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationSource(A,&iss)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationTarget(A,&ist)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&perm)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,perm)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(perm,ist,PETSC_FALSE)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(right,iss,PETSC_FALSE)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(A,PETSC_FALSE)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,right,left)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(left,-1.0,perm)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(left,NORM_INFINITY,&norm)); 1022c71b3e2SJacob 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); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(A,PETSC_TRUE)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&perm)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ist)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(B,NULL,"-B_view")); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(P,NORM_FROBENIUS,&relative)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(R,NORM_INFINITY,&norm)); 1232c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(norm/relative) > epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 127c7a4214aSPierre Jolivet } else { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(D,NULL,"-D_view")); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,D,10,&flg)); 13128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(A,D,10,&flg)); 13328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,D,10,&flg)); 13528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&begin,NULL)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(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); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(D,&ptr)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,D,10,&flg)); 14228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&AT)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(AT,D,10,&flg)); 14628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_REUSE_MATRIX,&AT)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(AT,D,10,&flg)); 14928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(D,NORM_INFINITY,&norm)); 1522c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(norm) > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AT)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,P,10,&flg)); 15928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMultEqual(A,B,P,10,&flg)); 16128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx"); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 164c7a4214aSPierre Jolivet if (n) { 1655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n*dim,&scoords)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValuesReal(rdm,n*dim,scoords)); 167c7a4214aSPierre Jolivet N = n; 1685f80ce2aSJacob Faibussowitsch CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(N*dim,&gscoords)); 1705f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(gscoords+begin*dim,scoords,n*dim)); 1725f80ce2aSJacob Faibussowitsch CHKERRMPI(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; 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(R)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(R,NULL,"-R_view")); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(D,NULL,"-D_view")); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(R,D,10,&flg)); 18428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&RT)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(RT,D,10,&flg)); 18828b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(R,MAT_REUSE_MATRIX,&RT)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(RT,D,10,&flg)); 19128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RT)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(R,B,P,10,&flg)); 20228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(R,&right,&left)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(right,rdm)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(R,right,left)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationSource(R,&iss)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationTarget(R,&ist)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&perm)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,perm)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(perm,ist,PETSC_FALSE)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(right,iss,PETSC_FALSE)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(R,PETSC_FALSE)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(R,right,left)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(left,-1.0,perm)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(left,NORM_INFINITY,&norm)); 2182c71b3e2SJacob 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); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(R,PETSC_TRUE)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&perm)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 2225f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ist)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iss)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gscoords)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(scoords)); 228c7a4214aSPierre Jolivet } 229c7a4214aSPierre Jolivet } 2305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gcoords)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(coords)); 234*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 235*b122ec5aSJacob 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