xref: /petsc/src/mat/tests/ex246.cxx (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
6d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntries(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr, void *ctx)
7d71ae5a4SJacob Faibussowitsch {
898e73e17SPierre Jolivet   PetscInt  d, j, k;
9c7a4214aSPierre Jolivet   PetscReal diff = 0.0, *coords = (PetscReal *)(ctx);
10c7a4214aSPierre Jolivet 
11c7a4214aSPierre Jolivet   PetscFunctionBeginUser;
1298e73e17SPierre Jolivet   for (j = 0; j < M; j++) {
1398e73e17SPierre Jolivet     for (k = 0; k < N; k++) {
1498e73e17SPierre Jolivet       diff = 0.0;
1598e73e17SPierre 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]);
1698e73e17SPierre Jolivet       ptr[j + M * k] = 1.0 / (1.0e-2 + PetscSqrtReal(diff));
17c7a4214aSPierre Jolivet     }
1898e73e17SPierre Jolivet   }
19*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2098e73e17SPierre Jolivet }
21cab00e0dSPierre Jolivet class MyIMatrix : public htool::VirtualGenerator<PetscScalar> {
22c7a4214aSPierre Jolivet private:
23c7a4214aSPierre Jolivet   PetscReal *coords;
24c7a4214aSPierre Jolivet   PetscInt   sdim;
259371c9d4SSatish Balay 
26c7a4214aSPierre Jolivet public:
27cab00e0dSPierre Jolivet   MyIMatrix(PetscInt M, PetscInt N, PetscInt spacedim, PetscReal *gcoords) : htool::VirtualGenerator<PetscScalar>(M, N), coords(gcoords), sdim(spacedim) { }
28c7a4214aSPierre Jolivet 
29d71ae5a4SJacob Faibussowitsch   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr) const override
30d71ae5a4SJacob Faibussowitsch   {
31c7a4214aSPierre Jolivet     PetscReal diff = 0.0;
32c7a4214aSPierre Jolivet 
33c7a4214aSPierre Jolivet     PetscFunctionBeginUser;
34c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < M; j++) /* could be optimized by the user how they see fit, e.g., vectorization */
3598e73e17SPierre Jolivet       for (PetscInt k = 0; k < N; k++) {
3698e73e17SPierre Jolivet         diff = 0.0;
3798e73e17SPierre 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]);
3898e73e17SPierre Jolivet         ptr[j + M * k] = 1.0 / (1.0e-2 + PetscSqrtReal(diff));
3998e73e17SPierre Jolivet       }
40c7a4214aSPierre Jolivet     PetscFunctionReturnVoid();
41c7a4214aSPierre Jolivet   }
42c7a4214aSPierre Jolivet };
43c7a4214aSPierre Jolivet 
44d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
45d71ae5a4SJacob Faibussowitsch {
46c7a4214aSPierre Jolivet   Mat            A, B, P, R;
47c7a4214aSPierre Jolivet   PetscInt       m = 100, dim = 3, M, begin = 0;
48c7a4214aSPierre Jolivet   PetscMPIInt    size;
49c7a4214aSPierre Jolivet   PetscReal     *coords, *gcoords, norm, epsilon, relative;
50c7a4214aSPierre Jolivet   PetscBool      sym = PETSC_FALSE;
51c7a4214aSPierre Jolivet   PetscRandom    rdm;
5298e73e17SPierre Jolivet   MatHtoolKernel kernel = GenEntries;
53c7a4214aSPierre Jolivet   MyIMatrix     *imatrix;
54c7a4214aSPierre Jolivet 
55327415f7SBarry Smith   PetscFunctionBeginUser;
569566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)NULL, help));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_local", &m, NULL));
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric", &sym, NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mat_htool_epsilon", &epsilon, NULL));
619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
62c7a4214aSPierre Jolivet   M = size * m;
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m * dim, &coords));
659566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
669566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal(rdm, m * dim, coords));
679566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(M * dim, &gcoords));
689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Exscan(&m, &begin, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD));
699566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(gcoords + begin * dim, coords, m * dim));
701c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, gcoords, M * dim, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD));
71c7a4214aSPierre Jolivet   imatrix = new MyIMatrix(M, M, dim, gcoords);
729566063dSJacob 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() */
739566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, sym));
749566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
779566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
789566063dSJacob Faibussowitsch   PetscCall(MatCreateHtoolFromKernel(PETSC_COMM_WORLD, m, m, M, M, dim, coords, coords, kernel, gcoords, &B)); /* entry-wise assembly using GenEntries() */
799566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_SYMMETRIC, sym));
809566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
839566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(B, NULL, "-B_view"));
849566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &P));
859566063dSJacob Faibussowitsch   PetscCall(MatNorm(P, NORM_FROBENIUS, &relative));
869566063dSJacob Faibussowitsch   PetscCall(MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &R));
879566063dSJacob Faibussowitsch   PetscCall(MatAXPY(R, -1.0, P, SAME_NONZERO_PATTERN));
889566063dSJacob Faibussowitsch   PetscCall(MatNorm(R, NORM_INFINITY, &norm));
8908401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(norm / relative) <= epsilon, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "||A(!symmetric)-A(symmetric)|| = %g (> %g)", (double)PetscAbsReal(norm / relative), (double)epsilon);
909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&R));
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
939566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
959566063dSJacob Faibussowitsch   PetscCall(PetscFree(gcoords));
969566063dSJacob Faibussowitsch   PetscCall(PetscFree(coords));
97c7a4214aSPierre Jolivet   delete imatrix;
989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
99b122ec5aSJacob Faibussowitsch   return 0;
100c7a4214aSPierre Jolivet }
101c7a4214aSPierre Jolivet 
102c7a4214aSPierre Jolivet /*TEST
103c7a4214aSPierre Jolivet 
104c7a4214aSPierre Jolivet    build:
105c7a4214aSPierre Jolivet       requires: htool
106c7a4214aSPierre Jolivet    test:
107c7a4214aSPierre Jolivet       requires: htool
108c7a4214aSPierre Jolivet       suffix: 2
109c7a4214aSPierre Jolivet       nsize: 4
110c7a4214aSPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2 -symmetric {{false true}shared output}
111c7a4214aSPierre Jolivet       output_file: output/ex101.out
112c7a4214aSPierre Jolivet 
113c7a4214aSPierre Jolivet TEST*/
114