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