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; 40c7a4214aSPierre Jolivet PetscInt m = 100,dim = 3,M,K = 10,begin,n = 0,N; 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; 4798e73e17SPierre Jolivet IS iss,ist; 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; 52c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL);CHKERRQ(ierr); 53c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL);CHKERRQ(ierr); 54c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); 55c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr); 56c7a4214aSPierre Jolivet ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL);CHKERRQ(ierr); 57c7a4214aSPierre Jolivet ierr = PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL);CHKERRQ(ierr); 58c7a4214aSPierre Jolivet ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 59c7a4214aSPierre Jolivet M = size*m; 60c7a4214aSPierre Jolivet ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 61c7a4214aSPierre Jolivet ierr = PetscMalloc1(m*dim,&coords);CHKERRQ(ierr); 62c7a4214aSPierre Jolivet ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 63c7a4214aSPierre Jolivet ierr = PetscRandomGetValuesReal(rdm,m*dim,coords);CHKERRQ(ierr); 64c7a4214aSPierre Jolivet ierr = PetscCalloc1(M*dim,&gcoords);CHKERRQ(ierr); 65c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B);CHKERRQ(ierr); 66c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68c7a4214aSPierre Jolivet ierr = MatSetRandom(B,rdm);CHKERRQ(ierr); 69c7a4214aSPierre Jolivet ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr); 70c7a4214aSPierre Jolivet ierr = PetscArraycpy(gcoords+begin*dim,coords,m*dim);CHKERRQ(ierr); 71c7a4214aSPierre Jolivet ierr = MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 72c7a4214aSPierre Jolivet ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A);CHKERRQ(ierr); 73c7a4214aSPierre Jolivet ierr = MatSetOption(A,MAT_SYMMETRIC,sym);CHKERRQ(ierr); 74c7a4214aSPierre Jolivet ierr = MatSetFromOptions(A);CHKERRQ(ierr); 75c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77c7a4214aSPierre Jolivet ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 7898e73e17SPierre Jolivet ierr = MatCreateVecs(A,&right,&left);CHKERRQ(ierr); 7998e73e17SPierre Jolivet ierr = VecSetRandom(right,rdm);CHKERRQ(ierr); 8098e73e17SPierre Jolivet ierr = MatMult(A,right,left);CHKERRQ(ierr); 8198e73e17SPierre Jolivet ierr = MatHtoolGetPermutationSource(A,&iss);CHKERRQ(ierr); 8298e73e17SPierre Jolivet ierr = MatHtoolGetPermutationTarget(A,&ist);CHKERRQ(ierr); 8398e73e17SPierre Jolivet ierr = VecDuplicate(left,&perm);CHKERRQ(ierr); 8498e73e17SPierre Jolivet ierr = VecCopy(left,perm);CHKERRQ(ierr); 8598e73e17SPierre Jolivet ierr = VecPermute(perm,ist,PETSC_FALSE);CHKERRQ(ierr); 8698e73e17SPierre Jolivet ierr = VecPermute(right,iss,PETSC_FALSE);CHKERRQ(ierr); 8798e73e17SPierre Jolivet ierr = MatHtoolUsePermutation(A,PETSC_FALSE);CHKERRQ(ierr); 8898e73e17SPierre Jolivet ierr = MatMult(A,right,left);CHKERRQ(ierr); 8998e73e17SPierre Jolivet ierr = VecAXPY(left,-1.0,perm);CHKERRQ(ierr); 9098e73e17SPierre Jolivet ierr = VecNorm(left,NORM_INFINITY,&norm);CHKERRQ(ierr); 9198e73e17SPierre Jolivet if (PetscAbsReal(norm) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL); 9298e73e17SPierre Jolivet ierr = MatHtoolUsePermutation(A,PETSC_TRUE);CHKERRQ(ierr); 9398e73e17SPierre Jolivet ierr = VecDestroy(&perm);CHKERRQ(ierr); 9498e73e17SPierre Jolivet ierr = VecDestroy(&left);CHKERRQ(ierr); 9598e73e17SPierre Jolivet ierr = VecDestroy(&right);CHKERRQ(ierr); 9698e73e17SPierre Jolivet ierr = ISDestroy(&ist);CHKERRQ(ierr); 9798e73e17SPierre Jolivet ierr = ISDestroy(&iss);CHKERRQ(ierr); 98c7a4214aSPierre 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 */ 99c7a4214aSPierre Jolivet PetscReal relative; 100c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 101c7a4214aSPierre Jolivet ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B);CHKERRQ(ierr); 102c7a4214aSPierre Jolivet ierr = MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym);CHKERRQ(ierr); 103c7a4214aSPierre Jolivet ierr = MatSetFromOptions(B);CHKERRQ(ierr); 104c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106c7a4214aSPierre Jolivet ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 107c7a4214aSPierre Jolivet ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P);CHKERRQ(ierr); 108c7a4214aSPierre Jolivet ierr = MatNorm(P,NORM_FROBENIUS,&relative);CHKERRQ(ierr); 109c7a4214aSPierre Jolivet ierr = MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr); 110c7a4214aSPierre Jolivet ierr = MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 111c7a4214aSPierre Jolivet ierr = MatNorm(R,NORM_INFINITY,&norm);CHKERRQ(ierr); 112c7a4214aSPierre Jolivet if (PetscAbsReal(norm/relative) > epsilon) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon); 113c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 114c7a4214aSPierre Jolivet ierr = MatDestroy(&R);CHKERRQ(ierr); 115c7a4214aSPierre Jolivet ierr = MatDestroy(&P);CHKERRQ(ierr); 116c7a4214aSPierre Jolivet } else { 117c7a4214aSPierre Jolivet ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 118c7a4214aSPierre Jolivet ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr); 119c7a4214aSPierre Jolivet ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr); 120c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 121*8b8ddd11SPierre Jolivet ierr = MatMultTransposeEqual(A,D,10,&flg);CHKERRQ(ierr); 122*8b8ddd11SPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 123c7a4214aSPierre Jolivet ierr = MatMultAddEqual(A,D,10,&flg);CHKERRQ(ierr); 124c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 125c7a4214aSPierre Jolivet ierr = MatGetOwnershipRange(B,&begin,NULL);CHKERRQ(ierr); 126c7a4214aSPierre Jolivet ierr = MatDenseGetArrayWrite(D,&ptr);CHKERRQ(ierr); 12798e73e17SPierre Jolivet for (PetscInt i = begin; i < m+begin; ++i) 12898e73e17SPierre Jolivet for (PetscInt j = 0; j < M; ++j) GenEntries(dim,1,1,&i,&j,ptr+i-begin+j*m,gcoords); 129c7a4214aSPierre Jolivet ierr = MatDenseRestoreArrayWrite(D,&ptr);CHKERRQ(ierr); 130c7a4214aSPierre Jolivet ierr = MatMultEqual(A,D,10,&flg);CHKERRQ(ierr); 131c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 132c7a4214aSPierre Jolivet ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr); 133c7a4214aSPierre Jolivet ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 134c7a4214aSPierre Jolivet ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr); 135c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 136c7a4214aSPierre Jolivet ierr = MatTranspose(A,MAT_REUSE_MATRIX,&AT);CHKERRQ(ierr); 137c7a4214aSPierre Jolivet ierr = MatMultEqual(AT,D,10,&flg);CHKERRQ(ierr); 138c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 139c7a4214aSPierre Jolivet ierr = MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 140c7a4214aSPierre Jolivet ierr = MatNorm(D,NORM_INFINITY,&norm);CHKERRQ(ierr); 141c7a4214aSPierre Jolivet if (PetscAbsReal(norm) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL); 142c7a4214aSPierre Jolivet ierr = MatDestroy(&AT);CHKERRQ(ierr); 143c7a4214aSPierre Jolivet ierr = MatDestroy(&D);CHKERRQ(ierr); 144c7a4214aSPierre Jolivet ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr); 145c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147c7a4214aSPierre Jolivet ierr = MatMatMultEqual(A,B,P,10,&flg);CHKERRQ(ierr); 148c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 149*8b8ddd11SPierre Jolivet ierr = MatTransposeMatMultEqual(A,B,P,10,&flg);CHKERRQ(ierr); 150*8b8ddd11SPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx"); 151c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 152c7a4214aSPierre Jolivet ierr = MatDestroy(&P);CHKERRQ(ierr); 153c7a4214aSPierre Jolivet if (n) { 154c7a4214aSPierre Jolivet ierr = PetscMalloc1(n*dim,&scoords);CHKERRQ(ierr); 155c7a4214aSPierre Jolivet ierr = PetscRandomGetValuesReal(rdm,n*dim,scoords);CHKERRQ(ierr); 156c7a4214aSPierre Jolivet N = n; 157c7a4214aSPierre Jolivet ierr = MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 158c7a4214aSPierre Jolivet ierr = PetscCalloc1(N*dim,&gscoords);CHKERRQ(ierr); 159c7a4214aSPierre Jolivet ierr = MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 160c7a4214aSPierre Jolivet ierr = PetscArraycpy(gscoords+begin*dim,scoords,n*dim);CHKERRQ(ierr); 161c7a4214aSPierre Jolivet ierr = MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD);CHKERRMPI(ierr); 16298e73e17SPierre Jolivet kernel = GenEntriesRectangular; 163c7a4214aSPierre Jolivet ctx[0] = gcoords; 164c7a4214aSPierre Jolivet ctx[1] = gscoords; 165c7a4214aSPierre Jolivet ierr = MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R);CHKERRQ(ierr); 166c7a4214aSPierre Jolivet ierr = MatSetFromOptions(R);CHKERRQ(ierr); 167c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169c7a4214aSPierre Jolivet ierr = MatViewFromOptions(R,NULL,"-R_view");CHKERRQ(ierr); 170c7a4214aSPierre Jolivet ierr = MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D);CHKERRQ(ierr); 171c7a4214aSPierre Jolivet ierr = MatViewFromOptions(D,NULL,"-D_view");CHKERRQ(ierr); 172c7a4214aSPierre Jolivet ierr = MatMultEqual(R,D,10,&flg);CHKERRQ(ierr); 173c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 174c7a4214aSPierre Jolivet ierr = MatTranspose(D,MAT_INPLACE_MATRIX,&D);CHKERRQ(ierr); 175c7a4214aSPierre Jolivet ierr = MatTranspose(R,MAT_INITIAL_MATRIX,&RT);CHKERRQ(ierr); 176c7a4214aSPierre Jolivet ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr); 177c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 178c7a4214aSPierre Jolivet ierr = MatTranspose(R,MAT_REUSE_MATRIX,&RT);CHKERRQ(ierr); 179c7a4214aSPierre Jolivet ierr = MatMultEqual(RT,D,10,&flg);CHKERRQ(ierr); 180c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 181c7a4214aSPierre Jolivet ierr = MatDestroy(&RT);CHKERRQ(ierr); 182c7a4214aSPierre Jolivet ierr = MatDestroy(&D);CHKERRQ(ierr); 183c7a4214aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B);CHKERRQ(ierr); 184c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 185c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 186c7a4214aSPierre Jolivet ierr = MatSetRandom(B,rdm);CHKERRQ(ierr); 187c7a4214aSPierre Jolivet ierr = MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P);CHKERRQ(ierr); 188c7a4214aSPierre Jolivet ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 189c7a4214aSPierre Jolivet ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190c7a4214aSPierre Jolivet ierr = MatMatMultEqual(R,B,P,10,&flg);CHKERRQ(ierr); 191c7a4214aSPierre Jolivet if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 192c7a4214aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 193c7a4214aSPierre Jolivet ierr = MatDestroy(&P);CHKERRQ(ierr); 19498e73e17SPierre Jolivet ierr = MatCreateVecs(R,&right,&left);CHKERRQ(ierr); 19598e73e17SPierre Jolivet ierr = VecSetRandom(right,rdm);CHKERRQ(ierr); 19698e73e17SPierre Jolivet ierr = MatMult(R,right,left);CHKERRQ(ierr); 19798e73e17SPierre Jolivet ierr = MatHtoolGetPermutationSource(R,&iss);CHKERRQ(ierr); 19898e73e17SPierre Jolivet ierr = MatHtoolGetPermutationTarget(R,&ist);CHKERRQ(ierr); 19998e73e17SPierre Jolivet ierr = VecDuplicate(left,&perm);CHKERRQ(ierr); 20098e73e17SPierre Jolivet ierr = VecCopy(left,perm);CHKERRQ(ierr); 20198e73e17SPierre Jolivet ierr = VecPermute(perm,ist,PETSC_FALSE);CHKERRQ(ierr); 20298e73e17SPierre Jolivet ierr = VecPermute(right,iss,PETSC_FALSE);CHKERRQ(ierr); 20398e73e17SPierre Jolivet ierr = MatHtoolUsePermutation(R,PETSC_FALSE);CHKERRQ(ierr); 20498e73e17SPierre Jolivet ierr = MatMult(R,right,left);CHKERRQ(ierr); 20598e73e17SPierre Jolivet ierr = VecAXPY(left,-1.0,perm);CHKERRQ(ierr); 20698e73e17SPierre Jolivet ierr = VecNorm(left,NORM_INFINITY,&norm);CHKERRQ(ierr); 20798e73e17SPierre Jolivet if (PetscAbsReal(norm) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL); 20898e73e17SPierre Jolivet ierr = MatHtoolUsePermutation(R,PETSC_TRUE);CHKERRQ(ierr); 20998e73e17SPierre Jolivet ierr = VecDestroy(&perm);CHKERRQ(ierr); 21098e73e17SPierre Jolivet ierr = VecDestroy(&left);CHKERRQ(ierr); 21198e73e17SPierre Jolivet ierr = VecDestroy(&right);CHKERRQ(ierr); 21298e73e17SPierre Jolivet ierr = ISDestroy(&ist);CHKERRQ(ierr); 21398e73e17SPierre Jolivet ierr = ISDestroy(&iss);CHKERRQ(ierr); 214c7a4214aSPierre Jolivet ierr = MatDestroy(&R);CHKERRQ(ierr); 215c7a4214aSPierre Jolivet ierr = PetscFree(gscoords);CHKERRQ(ierr); 216c7a4214aSPierre Jolivet ierr = PetscFree(scoords);CHKERRQ(ierr); 217c7a4214aSPierre Jolivet } 218c7a4214aSPierre Jolivet } 219c7a4214aSPierre Jolivet ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 220c7a4214aSPierre Jolivet ierr = MatDestroy(&A);CHKERRQ(ierr); 221c7a4214aSPierre Jolivet ierr = PetscFree(gcoords);CHKERRQ(ierr); 222c7a4214aSPierre Jolivet ierr = PetscFree(coords);CHKERRQ(ierr); 223c7a4214aSPierre Jolivet ierr = PetscFinalize(); 224c7a4214aSPierre Jolivet return ierr; 225c7a4214aSPierre Jolivet } 226c7a4214aSPierre Jolivet 227c7a4214aSPierre Jolivet /*TEST 228c7a4214aSPierre Jolivet 229c7a4214aSPierre Jolivet build: 230c7a4214aSPierre Jolivet requires: htool 231c7a4214aSPierre Jolivet 232c7a4214aSPierre Jolivet test: 233c7a4214aSPierre Jolivet requires: htool 234c7a4214aSPierre Jolivet suffix: 1 235c7a4214aSPierre Jolivet nsize: 4 236c7a4214aSPierre Jolivet args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output} 237c7a4214aSPierre Jolivet output_file: output/ex101.out 238c7a4214aSPierre Jolivet 239c7a4214aSPierre Jolivet test: 240c7a4214aSPierre Jolivet requires: htool 241c7a4214aSPierre Jolivet suffix: 2 242c7a4214aSPierre Jolivet nsize: 4 24398e73e17SPierre 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} 244c7a4214aSPierre Jolivet output_file: output/ex101.out 245c7a4214aSPierre Jolivet 246c7a4214aSPierre Jolivet TEST*/ 247