xref: /petsc/src/mat/impls/htool/htool.hpp (revision b5a865d84738d7e14adf4bbe5280b4a418a03baa)
1 #pragma once
2 
3 #include <petsc/private/matimpl.h>
4 #include <petscmathtool.h>
5 
6 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wsign-compare")
7 #include <htool/clustering/clustering.hpp>
8 #include <htool/hmatrix/lrmat/SVD.hpp>
9 #include <htool/hmatrix/lrmat/fullACA.hpp>
10 #include <htool/hmatrix/hmatrix_distributed_output.hpp>
11 #include <htool/hmatrix/linalg/factorization.hpp>
12 #include <htool/distributed_operator/utility.hpp>
13 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
14 
15 class WrapperHtool : public htool::VirtualGenerator<PetscScalar> {
16   MatHtoolKernelFn *&kernel;
17   PetscInt           sdim;
18   void              *ctx;
19 
20 public:
21   WrapperHtool(PetscInt dim, MatHtoolKernelFn *&g, void *kernelctx) : VirtualGenerator<PetscScalar>(), kernel(g), sdim(dim), ctx(kernelctx) { }
22   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr) const
23   {
24 #if !PetscDefined(HAVE_OPENMP)
25     PetscFunctionBegin;
26 #endif
27     PetscCallAbort(PETSC_COMM_SELF, kernel(sdim, M, N, rows, cols, ptr, ctx));
28 #if !PetscDefined(HAVE_OPENMP)
29     PetscFunctionReturnVoid();
30 #endif
31   }
32 };
33 
34 struct Mat_Htool {
35   PetscInt                                                            dim;
36   PetscReal                                                          *gcoords_target;
37   PetscReal                                                          *gcoords_source;
38   PetscInt                                                            min_cluster_size;
39   PetscReal                                                           epsilon;
40   PetscReal                                                           eta;
41   PetscInt                                                            depth[2];
42   PetscBool                                                           block_tree_consistency;
43   MatHtoolCompressorType                                              compressor;
44   MatHtoolClusteringType                                              clustering;
45   MatHtoolKernelFn                                                   *kernel;
46   void                                                               *kernelctx;
47   WrapperHtool                                                       *wrapper;
48   std::unique_ptr<htool::Cluster<PetscReal>>                          target_cluster;
49   std::unique_ptr<htool::Cluster<PetscReal>>                          source_cluster;
50   std::unique_ptr<htool::DistributedOperatorFromHMatrix<PetscScalar>> distributed_operator_holder;
51 };
52 
53 struct MatHtoolKernelTranspose {
54   Mat               A;
55   MatHtoolKernelFn *kernel;
56   void             *kernelctx;
57 };
58