xref: /petsc/src/mat/impls/htool/htool.hpp (revision 2bd2c0250c3b101ad40b9219faf1f20554191d26)
1a4963045SJacob Faibussowitsch #pragma once
23ea99036SJacob Faibussowitsch 
3c7a4214aSPierre Jolivet #include <petsc/private/matimpl.h>
41dd4f53aSPierre Jolivet #include <petscmathtool.h>
51dd4f53aSPierre Jolivet 
61dd4f53aSPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wsign-compare")
7*c095613aSPierre Jolivet #include <htool/clustering/tree_builder/tree_builder.hpp>
81dd4f53aSPierre Jolivet #include <htool/hmatrix/lrmat/SVD.hpp>
91dd4f53aSPierre Jolivet #include <htool/hmatrix/lrmat/fullACA.hpp>
10*c095613aSPierre Jolivet #include <htool/hmatrix/lrmat/recompressed_low_rank_generator.hpp>
111dd4f53aSPierre Jolivet #include <htool/hmatrix/hmatrix_distributed_output.hpp>
121dd4f53aSPierre Jolivet #include <htool/hmatrix/linalg/factorization.hpp>
131dd4f53aSPierre Jolivet #include <htool/distributed_operator/utility.hpp>
14*c095613aSPierre Jolivet #include <htool/distributed_operator/linalg/add_distributed_operator_vector_product_local_to_local.hpp>
15*c095613aSPierre Jolivet #include <htool/distributed_operator/linalg/add_distributed_operator_matrix_product_local_to_local.hpp>
161dd4f53aSPierre Jolivet PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
17c7a4214aSPierre Jolivet 
18cab00e0dSPierre Jolivet class WrapperHtool : public htool::VirtualGenerator<PetscScalar> {
198434afd1SBarry Smith   MatHtoolKernelFn *&kernel;
201dd4f53aSPierre Jolivet   PetscInt           sdim;
21c7a4214aSPierre Jolivet   void              *ctx;
22c7a4214aSPierre Jolivet 
23c7a4214aSPierre Jolivet public:
WrapperHtool(PetscInt dim,MatHtoolKernelFn * & g,void * kernelctx)241dd4f53aSPierre Jolivet   WrapperHtool(PetscInt dim, MatHtoolKernelFn *&g, void *kernelctx) : VirtualGenerator<PetscScalar>(), kernel(g), sdim(dim), ctx(kernelctx) { }
copy_submatrix(PetscInt M,PetscInt N,const PetscInt * rows,const PetscInt * cols,PetscScalar * ptr) const25d71ae5a4SJacob Faibussowitsch   void copy_submatrix(PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr) const
26d71ae5a4SJacob Faibussowitsch   {
27f782551bSPierre Jolivet #if !PetscDefined(HAVE_OPENMP)
2898e73e17SPierre Jolivet     PetscFunctionBegin;
29f782551bSPierre Jolivet #endif
301dd4f53aSPierre Jolivet     PetscCallAbort(PETSC_COMM_SELF, kernel(sdim, M, N, rows, cols, ptr, ctx));
31f782551bSPierre Jolivet #if !PetscDefined(HAVE_OPENMP)
3298e73e17SPierre Jolivet     PetscFunctionReturnVoid();
33f782551bSPierre Jolivet #endif
3498e73e17SPierre Jolivet   }
35c7a4214aSPierre Jolivet };
36c7a4214aSPierre Jolivet 
37c7a4214aSPierre Jolivet struct Mat_Htool {
38c7a4214aSPierre Jolivet   PetscInt                                                         dim;
39c7a4214aSPierre Jolivet   PetscReal                                                       *gcoords_target;
40c7a4214aSPierre Jolivet   PetscReal                                                       *gcoords_source;
41*c095613aSPierre Jolivet   PetscInt                                                         max_cluster_leaf_size;
42c7a4214aSPierre Jolivet   PetscReal                                                        epsilon;
43c7a4214aSPierre Jolivet   PetscReal                                                        eta;
44c7a4214aSPierre Jolivet   PetscInt                                                         depth[2];
451dd4f53aSPierre Jolivet   PetscBool                                                        block_tree_consistency;
46*c095613aSPierre Jolivet   PetscBool                                                        permutation;
47*c095613aSPierre Jolivet   PetscBool                                                        recompression;
48c7a4214aSPierre Jolivet   MatHtoolCompressorType                                           compressor;
4998e73e17SPierre Jolivet   MatHtoolClusteringType                                           clustering;
508434afd1SBarry Smith   MatHtoolKernelFn                                                *kernel;
51c7a4214aSPierre Jolivet   void                                                            *kernelctx;
52c7a4214aSPierre Jolivet   WrapperHtool                                                    *wrapper;
531dd4f53aSPierre Jolivet   std::unique_ptr<htool::Cluster<PetscReal>>                       target_cluster;
541dd4f53aSPierre Jolivet   std::unique_ptr<htool::Cluster<PetscReal>>                       source_cluster;
55*c095613aSPierre Jolivet   std::unique_ptr<htool::DefaultApproximationBuilder<PetscScalar>> distributed_operator_holder;
56c7a4214aSPierre Jolivet };
57c7a4214aSPierre Jolivet 
58c7a4214aSPierre Jolivet struct MatHtoolKernelTranspose {
59c7a4214aSPierre Jolivet   Mat               A;
608434afd1SBarry Smith   MatHtoolKernelFn *kernel;
61c7a4214aSPierre Jolivet   void             *kernelctx;
62c7a4214aSPierre Jolivet };
63