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