xref: /petsc/src/mat/tests/ex246.cxx (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c7a4214aSPierre Jolivet static char help[] = "Tests MATHTOOL with a derived htool::IMatrix<PetscScalar> class\n\n";
2c7a4214aSPierre Jolivet 
3c7a4214aSPierre Jolivet #include <petscmat.h>
41dd4f53aSPierre Jolivet #include <htool/hmatrix/interfaces/virtual_generator.hpp>
5c7a4214aSPierre Jolivet 
GenEntries(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt * J,const PetscInt * K,PetscScalar * ptr,PetscCtx ctx)6*2a8381b2SBarry Smith static PetscErrorCode GenEntries(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *J, const PetscInt *K, PetscScalar *ptr, PetscCtx 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   }
193ba16761SJacob 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:
MyIMatrix(PetscInt spacedim,PetscReal * gcoords)271dd4f53aSPierre Jolivet   MyIMatrix(PetscInt spacedim, PetscReal *gcoords) : htool::VirtualGenerator<PetscScalar>(), coords(gcoords), sdim(spacedim) { }
28c7a4214aSPierre Jolivet 
copy_submatrix(PetscInt M,PetscInt N,const PetscInt * J,const PetscInt * K,PetscScalar * ptr) const291dd4f53aSPierre Jolivet   virtual 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 
main(int argc,char ** argv)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;
501dd4f53aSPierre Jolivet   PetscBool         flg, sym = PETSC_FALSE;
51c7a4214aSPierre Jolivet   PetscRandom       rdm;
528434afd1SBarry Smith   MatHtoolKernelFn *kernel = GenEntries;
53c7a4214aSPierre Jolivet   MyIMatrix        *imatrix;
54c7a4214aSPierre Jolivet 
55327415f7SBarry Smith   PetscFunctionBeginUser;
561dd4f53aSPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, 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));
70462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, gcoords, M * dim, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD));
711dd4f53aSPierre Jolivet   imatrix = new MyIMatrix(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"));
841dd4f53aSPierre Jolivet   PetscCall(MatMultEqual(A, B, 10, &flg));
851dd4f53aSPierre Jolivet   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Ax != Bx");
869566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &P));
879566063dSJacob Faibussowitsch   PetscCall(MatNorm(P, NORM_FROBENIUS, &relative));
889566063dSJacob Faibussowitsch   PetscCall(MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &R));
899566063dSJacob Faibussowitsch   PetscCall(MatAXPY(R, -1.0, P, SAME_NONZERO_PATTERN));
909566063dSJacob Faibussowitsch   PetscCall(MatNorm(R, NORM_INFINITY, &norm));
9108401ef6SPierre Jolivet   PetscCheck(PetscAbsReal(norm / relative) <= epsilon, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "||A(!symmetric)-A(symmetric)|| = %g (> %g)", (double)PetscAbsReal(norm / relative), (double)epsilon);
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&R));
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
959566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
979566063dSJacob Faibussowitsch   PetscCall(PetscFree(gcoords));
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(coords));
99c7a4214aSPierre Jolivet   delete imatrix;
1009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
101b122ec5aSJacob Faibussowitsch   return 0;
102c7a4214aSPierre Jolivet }
103c7a4214aSPierre Jolivet 
104c7a4214aSPierre Jolivet /*TEST
105c7a4214aSPierre Jolivet 
106c7a4214aSPierre Jolivet    build:
107c7a4214aSPierre Jolivet       requires: htool
1081dd4f53aSPierre Jolivet    testset:
109c7a4214aSPierre Jolivet       requires: htool
1103886731fSPierre Jolivet       output_file: output/empty.out
111c7a4214aSPierre Jolivet       nsize: 4
1121dd4f53aSPierre Jolivet       args: -m_local 120 -mat_htool_epsilon 1.0e-2
1131dd4f53aSPierre Jolivet       test:
1141dd4f53aSPierre Jolivet         args: -symmetric
1151dd4f53aSPierre Jolivet       test:
1161dd4f53aSPierre Jolivet         args: -mat_htool_block_tree_consistency false
117c7a4214aSPierre Jolivet 
118c7a4214aSPierre Jolivet TEST*/
119