xref: /petsc/src/mat/impls/htool/htool.hpp (revision cab00e0df7744c4524b5b7245bc76cce28bdc10a)
1c7a4214aSPierre Jolivet #include <petsc/private/matimpl.h>
2c7a4214aSPierre Jolivet #include <htool/misc/petsc.hpp>
3c7a4214aSPierre Jolivet 
4*cab00e0dSPierre Jolivet class WrapperHtool : public htool::VirtualGenerator<PetscScalar> {
5c7a4214aSPierre Jolivet   PetscInt        dim;
6c7a4214aSPierre Jolivet   MatHtoolKernel& kernel;
7c7a4214aSPierre Jolivet   void*           ctx;
8c7a4214aSPierre Jolivet 
9c7a4214aSPierre Jolivet   public:
10*cab00e0dSPierre Jolivet   WrapperHtool(PetscInt M,PetscInt N,PetscInt sdim,MatHtoolKernel& g,void* kernelctx) : VirtualGenerator(M,N), dim(sdim), kernel(g), ctx(kernelctx) { }
1198e73e17SPierre Jolivet   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr) const {
1298e73e17SPierre Jolivet     PetscErrorCode ierr;
1398e73e17SPierre Jolivet 
1498e73e17SPierre Jolivet     PetscFunctionBegin;
1598e73e17SPierre Jolivet     ierr = kernel(dim,M,N,rows,cols,ptr,ctx);CHKERRABORT(PETSC_COMM_SELF,ierr);
1698e73e17SPierre Jolivet     PetscFunctionReturnVoid();
1798e73e17SPierre Jolivet   }
18c7a4214aSPierre Jolivet };
19c7a4214aSPierre Jolivet 
20c7a4214aSPierre Jolivet struct Mat_Htool {
21c7a4214aSPierre Jolivet   PetscInt               dim;
22c7a4214aSPierre Jolivet   PetscReal              *gcoords_target;
23c7a4214aSPierre Jolivet   PetscReal              *gcoords_source;
24c7a4214aSPierre Jolivet   PetscScalar            *work_target;
25c7a4214aSPierre Jolivet   PetscScalar            *work_source;
26c7a4214aSPierre Jolivet   PetscScalar            s;
27c7a4214aSPierre Jolivet   PetscInt               bs[2];
28c7a4214aSPierre Jolivet   PetscReal              epsilon;
29c7a4214aSPierre Jolivet   PetscReal              eta;
30c7a4214aSPierre Jolivet   PetscInt               depth[2];
31c7a4214aSPierre Jolivet   MatHtoolCompressorType compressor;
3298e73e17SPierre Jolivet   MatHtoolClusteringType clustering;
33c7a4214aSPierre Jolivet   MatHtoolKernel         kernel;
34c7a4214aSPierre Jolivet   void*                  kernelctx;
35c7a4214aSPierre Jolivet   WrapperHtool           *wrapper;
3698e73e17SPierre Jolivet   htool::VirtualHMatrix<PetscScalar> *hmatrix;
37c7a4214aSPierre Jolivet };
38c7a4214aSPierre Jolivet 
39c7a4214aSPierre Jolivet struct MatHtoolKernelTranspose {
40c7a4214aSPierre Jolivet   Mat            A;
41c7a4214aSPierre Jolivet   MatHtoolKernel kernel;
42c7a4214aSPierre Jolivet   void*          kernelctx;
43c7a4214aSPierre Jolivet };
44