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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)NULL,help)); 51*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL)); 52*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL)); 53*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 54*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL)); 55*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL)); 56*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL)); 57*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 58c7a4214aSPierre Jolivet M = size*m; 59*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 60*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m*dim,&coords)); 61*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 62*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValuesReal(rdm,m*dim,coords)); 63*9566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(M*dim,&gcoords)); 64*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B)); 65*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 66*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 67*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B,rdm)); 68*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&begin,NULL)); 69*9566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(gcoords+begin*dim,coords,m*dim)); 70*9566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 71*9566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A)); 72*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,sym)); 73*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 74*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 75*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 76*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-A_view")); 77*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A,is,NULL)); 78*9566063dSJacob Faibussowitsch PetscCall(ISDuplicate(is[0],is+1)); 79*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,1,is,2)); 80*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A,2)); 81*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,1,is+1,1)); 82*9566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is[1],&bs)); 832c71b3e2SJacob Faibussowitsch PetscCheckFalse(bs != 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs); 84*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A,1)); 85*9566063dSJacob Faibussowitsch PetscCall(ISEqual(is[0],is[1],&flg)); 8628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets"); 87*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(is)); 88*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(is+1)); 89*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&right,&left)); 90*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right,rdm)); 91*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,right,left)); 92*9566063dSJacob Faibussowitsch PetscCall(MatHtoolGetPermutationSource(A,&iss)); 93*9566063dSJacob Faibussowitsch PetscCall(MatHtoolGetPermutationTarget(A,&ist)); 94*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&perm)); 95*9566063dSJacob Faibussowitsch PetscCall(VecCopy(left,perm)); 96*9566063dSJacob Faibussowitsch PetscCall(VecPermute(perm,ist,PETSC_FALSE)); 97*9566063dSJacob Faibussowitsch PetscCall(VecPermute(right,iss,PETSC_FALSE)); 98*9566063dSJacob Faibussowitsch PetscCall(MatHtoolUsePermutation(A,PETSC_FALSE)); 99*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,right,left)); 100*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(left,-1.0,perm)); 101*9566063dSJacob Faibussowitsch PetscCall(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); 103*9566063dSJacob Faibussowitsch PetscCall(MatHtoolUsePermutation(A,PETSC_TRUE)); 104*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&perm)); 105*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 106*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 107*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ist)); 108*9566063dSJacob 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; 111*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 112*9566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); 113*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym)); 114*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 115*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 116*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 117*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(B,NULL,"-B_view")); 118*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P)); 119*9566063dSJacob Faibussowitsch PetscCall(MatNorm(P,NORM_FROBENIUS,&relative)); 120*9566063dSJacob Faibussowitsch PetscCall(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R)); 121*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN)); 122*9566063dSJacob Faibussowitsch PetscCall(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); 124*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 125*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 126*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 127c7a4214aSPierre Jolivet } else { 128*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D)); 129*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D,NULL,"-D_view")); 130*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,D,10,&flg)); 13128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 132*9566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A,D,10,&flg)); 13328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 134*9566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A,D,10,&flg)); 13528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 136*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&begin,NULL)); 137*9566063dSJacob 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); 140*9566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(D,&ptr)); 141*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,D,10,&flg)); 14228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 143*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 144*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT)); 145*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(AT,D,10,&flg)); 14628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 147*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_REUSE_MATRIX,&AT)); 148*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(AT,D,10,&flg)); 14928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 150*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN)); 151*9566063dSJacob Faibussowitsch PetscCall(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); 153*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 154*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 155*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 156*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 157*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 158*9566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A,B,P,10,&flg)); 15928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 160*9566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(A,B,P,10,&flg)); 16128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx"); 162*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 163*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 164c7a4214aSPierre Jolivet if (n) { 165*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n*dim,&scoords)); 166*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValuesReal(rdm,n*dim,scoords)); 167c7a4214aSPierre Jolivet N = n; 168*9566063dSJacob Faibussowitsch PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 169*9566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(N*dim,&gscoords)); 170*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 171*9566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(gscoords+begin*dim,scoords,n*dim)); 172*9566063dSJacob Faibussowitsch PetscCallMPI(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; 176*9566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R)); 177*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(R)); 178*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 179*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 180*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(R,NULL,"-R_view")); 181*9566063dSJacob Faibussowitsch PetscCall(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D)); 182*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(D,NULL,"-D_view")); 183*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(R,D,10,&flg)); 18428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 185*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 186*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(R,MAT_INITIAL_MATRIX,&RT)); 187*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(RT,D,10,&flg)); 18828b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 189*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(R,MAT_REUSE_MATRIX,&RT)); 190*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(RT,D,10,&flg)); 19128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 192*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RT)); 193*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 194*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B)); 195*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 196*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 197*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(B,rdm)); 198*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 199*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 200*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 201*9566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(R,B,P,10,&flg)); 20228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 203*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 204*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 205*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(R,&right,&left)); 206*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(right,rdm)); 207*9566063dSJacob Faibussowitsch PetscCall(MatMult(R,right,left)); 208*9566063dSJacob Faibussowitsch PetscCall(MatHtoolGetPermutationSource(R,&iss)); 209*9566063dSJacob Faibussowitsch PetscCall(MatHtoolGetPermutationTarget(R,&ist)); 210*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&perm)); 211*9566063dSJacob Faibussowitsch PetscCall(VecCopy(left,perm)); 212*9566063dSJacob Faibussowitsch PetscCall(VecPermute(perm,ist,PETSC_FALSE)); 213*9566063dSJacob Faibussowitsch PetscCall(VecPermute(right,iss,PETSC_FALSE)); 214*9566063dSJacob Faibussowitsch PetscCall(MatHtoolUsePermutation(R,PETSC_FALSE)); 215*9566063dSJacob Faibussowitsch PetscCall(MatMult(R,right,left)); 216*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(left,-1.0,perm)); 217*9566063dSJacob Faibussowitsch PetscCall(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); 219*9566063dSJacob Faibussowitsch PetscCall(MatHtoolUsePermutation(R,PETSC_TRUE)); 220*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&perm)); 221*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 222*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 223*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ist)); 224*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iss)); 225*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 226*9566063dSJacob Faibussowitsch PetscCall(PetscFree(gscoords)); 227*9566063dSJacob Faibussowitsch PetscCall(PetscFree(scoords)); 228c7a4214aSPierre Jolivet } 229c7a4214aSPierre Jolivet } 230*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 231*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 232*9566063dSJacob Faibussowitsch PetscCall(PetscFree(gcoords)); 233*9566063dSJacob Faibussowitsch PetscCall(PetscFree(coords)); 234*9566063dSJacob 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