xref: /petsc/src/mat/impls/htool/htool.hpp (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c7a4214aSPierre Jolivet #include <petsc/private/matimpl.h>
2c7a4214aSPierre Jolivet #include <htool/misc/petsc.hpp>
3c7a4214aSPierre Jolivet 
4cab00e0dSPierre Jolivet class WrapperHtool : public htool::VirtualGenerator<PetscScalar> {
5c7a4214aSPierre Jolivet   PetscInt        dim;
6c7a4214aSPierre Jolivet   MatHtoolKernel &kernel;
7c7a4214aSPierre Jolivet   void           *ctx;
8c7a4214aSPierre Jolivet 
9c7a4214aSPierre Jolivet public:
10cab00e0dSPierre Jolivet   WrapperHtool(PetscInt M, PetscInt N, PetscInt sdim, MatHtoolKernel &g, void *kernelctx) : VirtualGenerator(M, N), dim(sdim), kernel(g), ctx(kernelctx) { }
11*d71ae5a4SJacob Faibussowitsch   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr) const
12*d71ae5a4SJacob Faibussowitsch   {
13f782551bSPierre Jolivet #if !PetscDefined(HAVE_OPENMP)
1498e73e17SPierre Jolivet     PetscFunctionBegin;
15f782551bSPierre Jolivet #endif
169566063dSJacob Faibussowitsch     PetscCallAbort(PETSC_COMM_SELF, kernel(dim, M, N, rows, cols, ptr, ctx));
17f782551bSPierre Jolivet #if !PetscDefined(HAVE_OPENMP)
1898e73e17SPierre Jolivet     PetscFunctionReturnVoid();
19f782551bSPierre Jolivet #endif
2098e73e17SPierre Jolivet   }
21c7a4214aSPierre Jolivet };
22c7a4214aSPierre Jolivet 
23c7a4214aSPierre Jolivet struct Mat_Htool {
24c7a4214aSPierre Jolivet   PetscInt                            dim;
25c7a4214aSPierre Jolivet   PetscReal                          *gcoords_target;
26c7a4214aSPierre Jolivet   PetscReal                          *gcoords_source;
27c7a4214aSPierre Jolivet   PetscScalar                        *work_target;
28c7a4214aSPierre Jolivet   PetscScalar                        *work_source;
29c7a4214aSPierre Jolivet   PetscScalar                         s;
30c7a4214aSPierre Jolivet   PetscInt                            bs[2];
31c7a4214aSPierre Jolivet   PetscReal                           epsilon;
32c7a4214aSPierre Jolivet   PetscReal                           eta;
33c7a4214aSPierre Jolivet   PetscInt                            depth[2];
34c7a4214aSPierre Jolivet   MatHtoolCompressorType              compressor;
3598e73e17SPierre Jolivet   MatHtoolClusteringType              clustering;
36c7a4214aSPierre Jolivet   MatHtoolKernel                      kernel;
37c7a4214aSPierre Jolivet   void                               *kernelctx;
38c7a4214aSPierre Jolivet   WrapperHtool                       *wrapper;
3998e73e17SPierre Jolivet   htool::VirtualHMatrix<PetscScalar> *hmatrix;
40c7a4214aSPierre Jolivet };
41c7a4214aSPierre Jolivet 
42c7a4214aSPierre Jolivet struct MatHtoolKernelTranspose {
43c7a4214aSPierre Jolivet   Mat            A;
44c7a4214aSPierre Jolivet   MatHtoolKernel kernel;
45c7a4214aSPierre Jolivet   void          *kernelctx;
46c7a4214aSPierre Jolivet };
47