1 static char help[] = "Tests MATHTOOL\n\n"; 2 3 #include <petscmat.h> 4 5 static PetscErrorCode GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx) 6 { 7 PetscInt d,j,k; 8 PetscReal diff = 0.0,*coords = (PetscReal*)(ctx); 9 10 PetscFunctionBeginUser; 11 for (j = 0; j < M; j++) { 12 for (k = 0; k < N; k++) { 13 diff = 0.0; 14 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]); 15 ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff)); 16 } 17 } 18 PetscFunctionReturn(0); 19 } 20 21 static PetscErrorCode GenEntriesRectangular(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt *J,const PetscInt *K,PetscScalar *ptr,void *ctx) 22 { 23 PetscInt d,j,k; 24 PetscReal diff = 0.0,**coords = (PetscReal**)(ctx); 25 26 PetscFunctionBeginUser; 27 for (j = 0; j < M; j++) { 28 for (k = 0; k < N; k++) { 29 diff = 0.0; 30 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]); 31 ptr[j+M*k] = 1.0/(1.0e-2 + PetscSqrtReal(diff)); 32 } 33 } 34 PetscFunctionReturn(0); 35 } 36 37 int main(int argc,char **argv) 38 { 39 Mat A,AT,D,B,P,R,RT; 40 PetscInt m = 100,dim = 3,M,K = 10,begin,n = 0,N,bs; 41 PetscMPIInt size; 42 PetscScalar *ptr; 43 PetscReal *coords,*gcoords,*scoords,*gscoords,*(ctx[2]),norm,epsilon; 44 MatHtoolKernel kernel = GenEntries; 45 PetscBool flg,sym = PETSC_FALSE; 46 PetscRandom rdm; 47 IS iss,ist,is[2]; 48 Vec right,left,perm; 49 50 PetscCall(PetscInitialize(&argc,&argv,(char*)NULL,help)); 51 PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_local",&m,NULL)); 52 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_local",&n,NULL)); 53 PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 54 PetscCall(PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL)); 55 PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&sym,NULL)); 56 PetscCall(PetscOptionsGetReal(NULL,NULL,"-mat_htool_epsilon",&epsilon,NULL)); 57 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 58 M = size*m; 59 PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 60 PetscCall(PetscMalloc1(m*dim,&coords)); 61 PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 62 PetscCall(PetscRandomGetValuesReal(rdm,m*dim,coords)); 63 PetscCall(PetscCalloc1(M*dim,&gcoords)); 64 PetscCall(MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,NULL,&B)); 65 PetscCall(MatSetRandom(B,rdm)); 66 PetscCall(MatGetOwnershipRange(B,&begin,NULL)); 67 PetscCall(PetscArraycpy(gcoords+begin*dim,coords,m*dim)); 68 PetscCall(MPIU_Allreduce(MPI_IN_PLACE,gcoords,M*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 69 PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&A)); 70 PetscCall(MatSetOption(A,MAT_SYMMETRIC,sym)); 71 PetscCall(MatSetFromOptions(A)); 72 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 73 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 74 PetscCall(MatViewFromOptions(A,NULL,"-A_view")); 75 PetscCall(MatGetOwnershipIS(A,is,NULL)); 76 PetscCall(ISDuplicate(is[0],is+1)); 77 PetscCall(MatIncreaseOverlap(A,1,is,2)); 78 PetscCall(MatSetBlockSize(A,2)); 79 PetscCall(MatIncreaseOverlap(A,1,is+1,1)); 80 PetscCall(ISGetBlockSize(is[1],&bs)); 81 PetscCheck(bs == 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect block size %" PetscInt_FMT " != 2",bs); 82 PetscCall(MatSetBlockSize(A,1)); 83 PetscCall(ISEqual(is[0],is[1],&flg)); 84 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unequal index sets"); 85 PetscCall(ISDestroy(is)); 86 PetscCall(ISDestroy(is+1)); 87 PetscCall(MatCreateVecs(A,&right,&left)); 88 PetscCall(VecSetRandom(right,rdm)); 89 PetscCall(MatMult(A,right,left)); 90 PetscCall(MatHtoolGetPermutationSource(A,&iss)); 91 PetscCall(MatHtoolGetPermutationTarget(A,&ist)); 92 PetscCall(VecDuplicate(left,&perm)); 93 PetscCall(VecCopy(left,perm)); 94 PetscCall(VecPermute(perm,ist,PETSC_FALSE)); 95 PetscCall(VecPermute(right,iss,PETSC_FALSE)); 96 PetscCall(MatHtoolUsePermutation(A,PETSC_FALSE)); 97 PetscCall(MatMult(A,right,left)); 98 PetscCall(VecAXPY(left,-1.0,perm)); 99 PetscCall(VecNorm(left,NORM_INFINITY,&norm)); 100 PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL); 101 PetscCall(MatHtoolUsePermutation(A,PETSC_TRUE)); 102 PetscCall(VecDestroy(&perm)); 103 PetscCall(VecDestroy(&left)); 104 PetscCall(VecDestroy(&right)); 105 PetscCall(ISDestroy(&ist)); 106 PetscCall(ISDestroy(&iss)); 107 if (PetscAbsReal(epsilon) >= PETSC_SMALL) { /* when there is compression, it is more difficult to check against MATDENSE, so just compare symmetric and nonsymmetric assemblies */ 108 PetscReal relative; 109 PetscCall(MatDestroy(&B)); 110 PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,m,M,M,dim,coords,coords,kernel,gcoords,&B)); 111 PetscCall(MatSetOption(B,MAT_SYMMETRIC,(PetscBool)!sym)); 112 PetscCall(MatSetFromOptions(B)); 113 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 114 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 115 PetscCall(MatViewFromOptions(B,NULL,"-B_view")); 116 PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&P)); 117 PetscCall(MatNorm(P,NORM_FROBENIUS,&relative)); 118 PetscCall(MatConvert(B,MATDENSE,MAT_INITIAL_MATRIX,&R)); 119 PetscCall(MatAXPY(R,-1.0,P,SAME_NONZERO_PATTERN)); 120 PetscCall(MatNorm(R,NORM_INFINITY,&norm)); 121 PetscCheck(PetscAbsReal(norm/relative) <= epsilon,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A(!symmetric)-A(symmetric)|| = %g (> %g)",(double)PetscAbsReal(norm/relative),(double)epsilon); 122 PetscCall(MatDestroy(&B)); 123 PetscCall(MatDestroy(&R)); 124 PetscCall(MatDestroy(&P)); 125 } else { 126 PetscCall(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&D)); 127 PetscCall(MatViewFromOptions(D,NULL,"-D_view")); 128 PetscCall(MatMultEqual(A,D,10,&flg)); 129 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 130 PetscCall(MatMultTransposeEqual(A,D,10,&flg)); 131 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 132 PetscCall(MatMultAddEqual(A,D,10,&flg)); 133 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"y+Ax != y+Dx"); 134 PetscCall(MatGetOwnershipRange(B,&begin,NULL)); 135 PetscCall(MatDenseGetArrayWrite(D,&ptr)); 136 for (PetscInt i = begin; i < m+begin; ++i) 137 for (PetscInt j = 0; j < M; ++j) GenEntries(dim,1,1,&i,&j,ptr+i-begin+j*m,gcoords); 138 PetscCall(MatDenseRestoreArrayWrite(D,&ptr)); 139 PetscCall(MatMultEqual(A,D,10,&flg)); 140 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Ax != Dx"); 141 PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 142 PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT)); 143 PetscCall(MatMultEqual(AT,D,10,&flg)); 144 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 145 PetscCall(MatTranspose(A,MAT_REUSE_MATRIX,&AT)); 146 PetscCall(MatMultEqual(AT,D,10,&flg)); 147 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^Tx != D^Tx"); 148 PetscCall(MatAXPY(D,-1.0,AT,SAME_NONZERO_PATTERN)); 149 PetscCall(MatNorm(D,NORM_INFINITY,&norm)); 150 PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||A-D|| = %g (> %g)",(double)norm,(double)PETSC_SMALL); 151 PetscCall(MatDestroy(&AT)); 152 PetscCall(MatDestroy(&D)); 153 PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 154 PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 155 PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 156 PetscCall(MatMatMultEqual(A,B,P,10,&flg)); 157 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ABx != Px"); 158 PetscCall(MatTransposeMatMultEqual(A,B,P,10,&flg)); 159 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"A^TBx != P^Tx"); 160 PetscCall(MatDestroy(&B)); 161 PetscCall(MatDestroy(&P)); 162 if (n) { 163 PetscCall(PetscMalloc1(n*dim,&scoords)); 164 PetscCall(PetscRandomGetValuesReal(rdm,n*dim,scoords)); 165 N = n; 166 PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&N,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 167 PetscCall(PetscCalloc1(N*dim,&gscoords)); 168 PetscCallMPI(MPI_Exscan(&n,&begin,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD)); 169 PetscCall(PetscArraycpy(gscoords+begin*dim,scoords,n*dim)); 170 PetscCall(MPIU_Allreduce(MPI_IN_PLACE,gscoords,N*dim,MPIU_REAL,MPI_SUM,PETSC_COMM_WORLD)); 171 kernel = GenEntriesRectangular; 172 ctx[0] = gcoords; 173 ctx[1] = gscoords; 174 PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD,m,n,M,N,dim,coords,scoords,kernel,ctx,&R)); 175 PetscCall(MatSetFromOptions(R)); 176 PetscCall(MatAssemblyBegin(R,MAT_FINAL_ASSEMBLY)); 177 PetscCall(MatAssemblyEnd(R,MAT_FINAL_ASSEMBLY)); 178 PetscCall(MatViewFromOptions(R,NULL,"-R_view")); 179 PetscCall(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&D)); 180 PetscCall(MatViewFromOptions(D,NULL,"-D_view")); 181 PetscCall(MatMultEqual(R,D,10,&flg)); 182 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Rx != Dx"); 183 PetscCall(MatTranspose(D,MAT_INPLACE_MATRIX,&D)); 184 PetscCall(MatTranspose(R,MAT_INITIAL_MATRIX,&RT)); 185 PetscCall(MatMultEqual(RT,D,10,&flg)); 186 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 187 PetscCall(MatTranspose(R,MAT_REUSE_MATRIX,&RT)); 188 PetscCall(MatMultEqual(RT,D,10,&flg)); 189 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"R^Tx != D^Tx"); 190 PetscCall(MatDestroy(&RT)); 191 PetscCall(MatDestroy(&D)); 192 PetscCall(MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_DETERMINE,K,NULL,&B)); 193 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 194 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 195 PetscCall(MatSetRandom(B,rdm)); 196 PetscCall(MatMatMult(R,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&P)); 197 PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 198 PetscCall(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 199 PetscCall(MatMatMultEqual(R,B,P,10,&flg)); 200 PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"RBx != Px"); 201 PetscCall(MatDestroy(&B)); 202 PetscCall(MatDestroy(&P)); 203 PetscCall(MatCreateVecs(R,&right,&left)); 204 PetscCall(VecSetRandom(right,rdm)); 205 PetscCall(MatMult(R,right,left)); 206 PetscCall(MatHtoolGetPermutationSource(R,&iss)); 207 PetscCall(MatHtoolGetPermutationTarget(R,&ist)); 208 PetscCall(VecDuplicate(left,&perm)); 209 PetscCall(VecCopy(left,perm)); 210 PetscCall(VecPermute(perm,ist,PETSC_FALSE)); 211 PetscCall(VecPermute(right,iss,PETSC_FALSE)); 212 PetscCall(MatHtoolUsePermutation(R,PETSC_FALSE)); 213 PetscCall(MatMult(R,right,left)); 214 PetscCall(VecAXPY(left,-1.0,perm)); 215 PetscCall(VecNorm(left,NORM_INFINITY,&norm)); 216 PetscCheck(PetscAbsReal(norm) <= PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"||y(with permutation)-y(without permutation)|| = %g (> %g)",(double)PetscAbsReal(norm),(double)PETSC_SMALL); 217 PetscCall(MatHtoolUsePermutation(R,PETSC_TRUE)); 218 PetscCall(VecDestroy(&perm)); 219 PetscCall(VecDestroy(&left)); 220 PetscCall(VecDestroy(&right)); 221 PetscCall(ISDestroy(&ist)); 222 PetscCall(ISDestroy(&iss)); 223 PetscCall(MatDestroy(&R)); 224 PetscCall(PetscFree(gscoords)); 225 PetscCall(PetscFree(scoords)); 226 } 227 } 228 PetscCall(PetscRandomDestroy(&rdm)); 229 PetscCall(MatDestroy(&A)); 230 PetscCall(PetscFree(gcoords)); 231 PetscCall(PetscFree(coords)); 232 PetscCall(PetscFinalize()); 233 return 0; 234 } 235 236 /*TEST 237 238 build: 239 requires: htool 240 241 test: 242 requires: htool 243 suffix: 1 244 nsize: 4 245 args: -m_local 80 -n_local 25 -mat_htool_epsilon 1.0e-11 -symmetric {{false true}shared output} 246 output_file: output/ex101.out 247 248 test: 249 requires: htool 250 suffix: 2 251 nsize: 4 252 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} 253 output_file: output/ex101.out 254 255 TEST*/ 256