1c7a4214aSPierre Jolivet static char help[] = "Tests MATHTOOL with a derived htool::IMatrix<PetscScalar> class\n\n"; 2c7a4214aSPierre Jolivet 3c7a4214aSPierre Jolivet #include <petscmat.h> 4c7a4214aSPierre Jolivet #include <htool/misc/petsc.hpp> 5c7a4214aSPierre Jolivet 6*9371c9d4SSatish Balay static PetscErrorCode GenEntries(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr, void *ctx) { 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)); 16c7a4214aSPierre Jolivet } 1798e73e17SPierre Jolivet } 1898e73e17SPierre Jolivet PetscFunctionReturn(0); 1998e73e17SPierre Jolivet } 20cab00e0dSPierre Jolivet class MyIMatrix : public htool::VirtualGenerator<PetscScalar> { 21c7a4214aSPierre Jolivet private: 22c7a4214aSPierre Jolivet PetscReal *coords; 23c7a4214aSPierre Jolivet PetscInt sdim; 24*9371c9d4SSatish Balay 25c7a4214aSPierre Jolivet public: 26cab00e0dSPierre Jolivet MyIMatrix(PetscInt M, PetscInt N, PetscInt spacedim, PetscReal *gcoords) : htool::VirtualGenerator<PetscScalar>(M, N), coords(gcoords), sdim(spacedim) { } 27c7a4214aSPierre Jolivet 28*9371c9d4SSatish Balay void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr) const override { 29c7a4214aSPierre Jolivet PetscReal diff = 0.0; 30c7a4214aSPierre Jolivet 31c7a4214aSPierre Jolivet PetscFunctionBeginUser; 32c7a4214aSPierre Jolivet for (PetscInt j = 0; j < M; j++) /* could be optimized by the user how they see fit, e.g., vectorization */ 3398e73e17SPierre Jolivet for (PetscInt k = 0; k < N; k++) { 3498e73e17SPierre Jolivet diff = 0.0; 3598e73e17SPierre Jolivet for (PetscInt 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]); 3698e73e17SPierre Jolivet ptr[j + M * k] = 1.0 / (1.0e-2 + PetscSqrtReal(diff)); 3798e73e17SPierre Jolivet } 38c7a4214aSPierre Jolivet PetscFunctionReturnVoid(); 39c7a4214aSPierre Jolivet } 40c7a4214aSPierre Jolivet }; 41c7a4214aSPierre Jolivet 42*9371c9d4SSatish Balay int main(int argc, char **argv) { 43c7a4214aSPierre Jolivet Mat A, B, P, R; 44c7a4214aSPierre Jolivet PetscInt m = 100, dim = 3, M, begin = 0; 45c7a4214aSPierre Jolivet PetscMPIInt size; 46c7a4214aSPierre Jolivet PetscReal *coords, *gcoords, norm, epsilon, relative; 47c7a4214aSPierre Jolivet PetscBool sym = PETSC_FALSE; 48c7a4214aSPierre Jolivet PetscRandom rdm; 4998e73e17SPierre Jolivet MatHtoolKernel kernel = GenEntries; 50c7a4214aSPierre Jolivet MyIMatrix *imatrix; 51c7a4214aSPierre Jolivet 52327415f7SBarry Smith PetscFunctionBeginUser; 539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)NULL, help)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_local", &m, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric", &sym, NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mat_htool_epsilon", &epsilon, NULL)); 589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 59c7a4214aSPierre Jolivet M = size * m; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * dim, &coords)); 629566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm)); 639566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValuesReal(rdm, m * dim, coords)); 649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(M * dim, &gcoords)); 659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&m, &begin, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD)); 669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(gcoords + begin * dim, coords, m * dim)); 671c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, gcoords, M * dim, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD)); 68c7a4214aSPierre Jolivet imatrix = new MyIMatrix(M, M, dim, gcoords); 699566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD, m, m, M, M, dim, coords, coords, NULL, imatrix, &A)); /* block-wise assembly using htool::IMatrix<PetscScalar>::copy_submatrix() */ 709566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, sym)); 719566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 759566063dSJacob Faibussowitsch PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD, m, m, M, M, dim, coords, coords, kernel, gcoords, &B)); /* entry-wise assembly using GenEntries() */ 769566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SYMMETRIC, sym)); 779566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 809566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(B, NULL, "-B_view")); 819566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &P)); 829566063dSJacob Faibussowitsch PetscCall(MatNorm(P, NORM_FROBENIUS, &relative)); 839566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &R)); 849566063dSJacob Faibussowitsch PetscCall(MatAXPY(R, -1.0, P, SAME_NONZERO_PATTERN)); 859566063dSJacob Faibussowitsch PetscCall(MatNorm(R, NORM_INFINITY, &norm)); 8608401ef6SPierre Jolivet PetscCheck(PetscAbsReal(norm / relative) <= epsilon, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "||A(!symmetric)-A(symmetric)|| = %g (> %g)", (double)PetscAbsReal(norm / relative), (double)epsilon); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&R)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 909566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(gcoords)); 939566063dSJacob Faibussowitsch PetscCall(PetscFree(coords)); 94c7a4214aSPierre Jolivet delete imatrix; 959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 96b122ec5aSJacob Faibussowitsch return 0; 97c7a4214aSPierre Jolivet } 98c7a4214aSPierre Jolivet 99c7a4214aSPierre Jolivet /*TEST 100c7a4214aSPierre Jolivet 101c7a4214aSPierre Jolivet build: 102c7a4214aSPierre Jolivet requires: htool 103c7a4214aSPierre Jolivet test: 104c7a4214aSPierre Jolivet requires: htool 105c7a4214aSPierre Jolivet suffix: 2 106c7a4214aSPierre Jolivet nsize: 4 107c7a4214aSPierre Jolivet args: -m_local 120 -mat_htool_epsilon 1.0e-2 -symmetric {{false true}shared output} 108c7a4214aSPierre Jolivet output_file: output/ex101.out 109c7a4214aSPierre Jolivet 110c7a4214aSPierre Jolivet TEST*/ 111