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; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL)); 58*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 59c7a4214aSPierre Jolivet M = size*m; 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*dim,&coords)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValuesReal(rdm,m*dim,coords)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(M*dim,&gcoords)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&begin,NULL)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(gcoords+begin*dim,coords,m*dim)); 71*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,sym)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-A_view")); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(A,is,NULL)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDuplicate(is[0],is+1)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,1,is,2)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,2)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,1,is+1,1)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetBlockSize(is[1],&bs)); 842c71b3e2SJacob Faibussowitsch PetscCheckFalse(bs != 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,1)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISEqual(is[0],is[1],&flg)); 872c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets"); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(is)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(is+1)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&right,&left)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(right,rdm)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,right,left)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationSource(A,&iss)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationTarget(A,&ist)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&perm)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,perm)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(perm,ist,PETSC_FALSE)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(right,iss,PETSC_FALSE)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(A,PETSC_FALSE)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,right,left)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(left,-1.0,perm)); 102*5f80ce2aSJacob 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); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(A,PETSC_TRUE)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&perm)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ist)); 109*5f80ce2aSJacob 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; 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(B,NULL,"-B_view")); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(P,NORM_FROBENIUS,&relative)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN)); 123*5f80ce2aSJacob 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); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 128c7a4214aSPierre Jolivet } else { 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(D,NULL,"-D_view")); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,D,10,&flg)); 1322c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(A,D,10,&flg)); 1342c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,D,10,&flg)); 1362c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(B,&begin,NULL)); 138*5f80ce2aSJacob 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); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(D,&ptr)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,D,10,&flg)); 1432c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&AT)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(AT,D,10,&flg)); 1472c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_REUSE_MATRIX,&AT)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(AT,D,10,&flg)); 1502c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN)); 152*5f80ce2aSJacob 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); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AT)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,P,10,&flg)); 1602c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMultEqual(A,B,P,10,&flg)); 1622c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx"); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 165c7a4214aSPierre Jolivet if (n) { 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n*dim,&scoords)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValuesReal(rdm,n*dim,scoords)); 168c7a4214aSPierre Jolivet N = n; 169*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(N*dim,&gscoords)); 171*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(gscoords+begin*dim,scoords,n*dim)); 173*5f80ce2aSJacob 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; 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(R)); 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(R,NULL,"-R_view")); 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D)); 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(D,NULL,"-D_view")); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(R,D,10,&flg)); 1852c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(R,MAT_INITIAL_MATRIX,&RT)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(RT,D,10,&flg)); 1892c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(R,MAT_REUSE_MATRIX,&RT)); 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(RT,D,10,&flg)); 1922c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RT)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(B,rdm)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(R,B,P,10,&flg)); 2032c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 204*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(R,&right,&left)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(right,rdm)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(R,right,left)); 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationSource(R,&iss)); 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolGetPermutationTarget(R,&ist)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&perm)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,perm)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(perm,ist,PETSC_FALSE)); 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPermute(right,iss,PETSC_FALSE)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(R,PETSC_FALSE)); 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(R,right,left)); 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(left,-1.0,perm)); 218*5f80ce2aSJacob 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); 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHtoolUsePermutation(R,PETSC_TRUE)); 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&perm)); 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 223*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 224*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ist)); 225*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iss)); 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gscoords)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(scoords)); 229c7a4214aSPierre Jolivet } 230c7a4214aSPierre Jolivet } 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(gcoords)); 234*5f80ce2aSJacob 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