xref: /petsc/src/mat/impls/htool/htool.cxx (revision 7530f0189ecda25d023fdcb28178244332f3efb5)
1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/
2c7a4214aSPierre Jolivet #include <set>
3c7a4214aSPierre Jolivet 
4c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = {"sympartialACA", "fullACA", "SVD"};
598e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = {"PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric"};
6*c895e7bbSPierre Jolivet // clang-format off
7*c895e7bbSPierre Jolivet const char       *HtoolCitations[2]         = {"@article{marchand2020two,\n"
8c7a4214aSPierre Jolivet                                                "  Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n"
9c7a4214aSPierre Jolivet                                                "  Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n"
10c7a4214aSPierre Jolivet                                                "  Year = {2020},\n"
11c7a4214aSPierre Jolivet                                                "  Publisher = {Elsevier},\n"
12c7a4214aSPierre Jolivet                                                "  Journal = {Numerische Mathematik},\n"
13c7a4214aSPierre Jolivet                                                "  Volume = {146},\n"
14c7a4214aSPierre Jolivet                                                "  Pages = {597--628},\n"
15c7a4214aSPierre Jolivet                                                "  Url = {https://github.com/htool-ddm/htool}\n"
16*c895e7bbSPierre Jolivet                                                "}\n",
17*c895e7bbSPierre Jolivet                                                "@article{Marchand2026,\n"
18*c895e7bbSPierre Jolivet                                                "  Author = {Marchand, Pierre and Tournier, Pierre-Henri and Jolivet, Pierre},\n"
19*c895e7bbSPierre Jolivet                                                "  Title = {{Htool-DDM}: A {C++} library for parallel solvers and compressed linear systems},\n"
20*c895e7bbSPierre Jolivet                                                "  Year = {2026},\n"
21*c895e7bbSPierre Jolivet                                                "  Publisher = {The Open Journal},\n"
22*c895e7bbSPierre Jolivet                                                "  Journal = {Journal of Open Source Software},\n"
23*c895e7bbSPierre Jolivet                                                "  Volume = {11},\n"
24*c895e7bbSPierre Jolivet                                                "  Number = {118},\n"
25*c895e7bbSPierre Jolivet                                                "  Pages = {9279},\n"
26*c895e7bbSPierre Jolivet                                                "  Url = {https://doi.org/10.21105/joss.09279}\n"
27*c895e7bbSPierre Jolivet                                                "}\n"};
28*c895e7bbSPierre Jolivet static PetscBool  HtoolCite[2]              = {PETSC_FALSE, PETSC_FALSE};
29*c895e7bbSPierre Jolivet // clang-format on
30c7a4214aSPierre Jolivet 
MatGetDiagonal_Htool(Mat A,Vec v)31d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Htool(Mat A, Vec v)
32d71ae5a4SJacob Faibussowitsch {
33c9a4fd69SAudic XU   Mat_Htool   *a;
34c7a4214aSPierre Jolivet   PetscScalar *x;
35c7a4214aSPierre Jolivet   PetscBool    flg;
36c7a4214aSPierre Jolivet 
37c7a4214aSPierre Jolivet   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
395f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
40c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(v, &x));
421dd4f53aSPierre Jolivet   PetscStackCallExternalVoid("copy_diagonal_in_user_numbering", htool::copy_diagonal_in_user_numbering(a->distributed_operator_holder->hmatrix, x));
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(v, &x));
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45c7a4214aSPierre Jolivet }
46c7a4214aSPierre Jolivet 
MatGetDiagonalBlock_Htool(Mat A,Mat * b)47d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b)
48d71ae5a4SJacob Faibussowitsch {
49c9a4fd69SAudic XU   Mat_Htool                 *a;
50c7a4214aSPierre Jolivet   Mat                        B;
511dd4f53aSPierre Jolivet   PetscScalar               *ptr, shift, scale;
52c7a4214aSPierre Jolivet   PetscBool                  flg;
531dd4f53aSPierre Jolivet   PetscMPIInt                rank;
541dd4f53aSPierre Jolivet   htool::Cluster<PetscReal> *source_cluster = nullptr;
55c7a4214aSPierre Jolivet 
56c7a4214aSPierre Jolivet   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
585f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
59c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
61c7a4214aSPierre Jolivet   if (!B) {
6266b4df27SPierre Jolivet     PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
63f22e26b7SPierre Jolivet     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B));
649566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B, &ptr));
651dd4f53aSPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
661dd4f53aSPierre Jolivet     source_cluster = a->source_cluster ? a->source_cluster.get() : a->target_cluster.get();
671dd4f53aSPierre Jolivet     PetscStackCallExternalVoid("copy_to_dense_in_user_numbering", htool::copy_to_dense_in_user_numbering(*a->distributed_operator_holder->hmatrix.get_sub_hmatrix(a->target_cluster->get_cluster_on_partition(rank), source_cluster->get_cluster_on_partition(rank)), ptr));
689566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B, &ptr));
699566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A, B));
709566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B));
71c7a4214aSPierre Jolivet     *b = B;
729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
7366b4df27SPierre Jolivet     PetscCall(MatShift(*b, shift));
7466b4df27SPierre Jolivet     PetscCall(MatScale(*b, scale));
75c9a4fd69SAudic XU   } else {
7666b4df27SPierre Jolivet     PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
77c9a4fd69SAudic XU     *b = B;
78c9a4fd69SAudic XU   }
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80c7a4214aSPierre Jolivet }
81c7a4214aSPierre Jolivet 
MatMult_Htool(Mat A,Vec x,Vec y)82d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y)
83d71ae5a4SJacob Faibussowitsch {
84c9a4fd69SAudic XU   Mat_Htool         *a;
85c7a4214aSPierre Jolivet   const PetscScalar *in;
86c7a4214aSPierre Jolivet   PetscScalar       *out;
87c7a4214aSPierre Jolivet 
88c7a4214aSPierre Jolivet   PetscFunctionBegin;
89c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
92c095613aSPierre Jolivet   if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_vector_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr);
93c095613aSPierre Jolivet   else htool::internal_add_distributed_operator_vector_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr);
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97c7a4214aSPierre Jolivet }
98c7a4214aSPierre Jolivet 
MatMultTranspose_Htool(Mat A,Vec x,Vec y)99d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y)
100d71ae5a4SJacob Faibussowitsch {
101c9a4fd69SAudic XU   Mat_Htool         *a;
1028b8ddd11SPierre Jolivet   const PetscScalar *in;
1038b8ddd11SPierre Jolivet   PetscScalar       *out;
1048b8ddd11SPierre Jolivet 
1058b8ddd11SPierre Jolivet   PetscFunctionBegin;
106c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
109c095613aSPierre Jolivet   if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_vector_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr);
110c095613aSPierre Jolivet   else htool::internal_add_distributed_operator_vector_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr);
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1148b8ddd11SPierre Jolivet }
1158b8ddd11SPierre Jolivet 
MatIncreaseOverlap_Htool(Mat A,PetscInt is_max,IS is[],PetscInt ov)116d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov)
117d71ae5a4SJacob Faibussowitsch {
118c7a4214aSPierre Jolivet   std::set<PetscInt> set;
119c7a4214aSPierre Jolivet   const PetscInt    *idx;
1208f308287SPierre Jolivet   PetscInt          *oidx, size, bs[2];
121c7a4214aSPierre Jolivet   PetscMPIInt        csize;
122c7a4214aSPierre Jolivet 
123c7a4214aSPierre Jolivet   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, bs, bs + 1));
1258f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
126c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < is_max; ++i) {
127c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
128c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1299566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize));
1301dd4f53aSPierre Jolivet     PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel IS");
1319566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i], &size));
1329566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx));
133c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < size; ++j) {
134c7a4214aSPierre Jolivet       set.insert(idx[j]);
135c7a4214aSPierre Jolivet       for (PetscInt k = 1; k <= ov; ++k) {                                              /* for each layer of overlap      */
136c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k);                                    /* do not insert negative indices */
137c7a4214aSPierre Jolivet         if (idx[j] + k < A->rmap->N && idx[j] + k < A->cmap->N) set.insert(idx[j] + k); /* do not insert indices greater than the dimension of A */
138c7a4214aSPierre Jolivet       }
139c7a4214aSPierre Jolivet     }
1409566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i], &idx));
1419566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is + i));
1428f308287SPierre Jolivet     if (bs[0] > 1) {
1438f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) {
1448f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1458f308287SPierre Jolivet         std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]);
1468f308287SPierre Jolivet         set.insert(block.cbegin(), block.cend());
1478f308287SPierre Jolivet       }
1488f308287SPierre Jolivet     }
149c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &oidx));
151c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
152c7a4214aSPierre Jolivet     oidx -= size;
1539566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i));
154c7a4214aSPierre Jolivet   }
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156c7a4214aSPierre Jolivet }
157c7a4214aSPierre Jolivet 
MatCreateSubMatrices_Htool(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat * submat[])158d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
159d71ae5a4SJacob Faibussowitsch {
160c9a4fd69SAudic XU   Mat_Htool         *a;
161c7a4214aSPierre Jolivet   Mat                D, B, BT;
162c7a4214aSPierre Jolivet   const PetscScalar *copy;
16366b4df27SPierre Jolivet   PetscScalar       *ptr, shift, scale;
164c7a4214aSPierre Jolivet   const PetscInt    *idxr, *idxc, *it;
165c7a4214aSPierre Jolivet   PetscInt           nrow, m, i;
166c7a4214aSPierre Jolivet   PetscBool          flg;
167c7a4214aSPierre Jolivet 
168c7a4214aSPierre Jolivet   PetscFunctionBegin;
16966b4df27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
170c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
17148a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
172c7a4214aSPierre Jolivet   for (i = 0; i < n; ++i) {
1739566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &nrow));
1749566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &m));
1759566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &idxr));
1769566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &idxc));
177f22e26b7SPierre Jolivet     if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i));
1789566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr));
179c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1809566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A, &flg));
181c7a4214aSPierre Jolivet       if (flg) {
1829566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i], &flg));
183c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
184c7a4214aSPierre Jolivet           it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart);
185c7a4214aSPierre Jolivet           if (it != idxr + nrow && *it == A->rmap->rstart) {    /* rmap->rstart in IS? */
186c7a4214aSPierre Jolivet             if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
1879371c9d4SSatish Balay               for (PetscInt j = 0; j < A->rmap->n && flg; ++j)
1889371c9d4SSatish Balay                 if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE;
189c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
190c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
191c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
192c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
193c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
194c7a4214aSPierre Jolivet                  */
195c7a4214aSPierre Jolivet                 m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */
196c9a4fd69SAudic XU                 PetscCall(MatGetDiagonalBlock(A, &D));
1979566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D, &copy));
198ac530a7eSPierre Jolivet                 for (PetscInt k = 0; k < A->rmap->n; ++k) PetscCall(PetscArraycpy(ptr + (m + k) * nrow + m, copy + k * A->rmap->n, A->rmap->n)); /* block D from above */
1999566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D, &copy));
200c7a4214aSPierre Jolivet                 if (m) {
201c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */
202c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
203b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
2049566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B));
2059566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
2069566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT));
2079566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
208b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
2099566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
210c7a4214aSPierre Jolivet                     } else {
2117fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
2129566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
213c7a4214aSPierre Jolivet                     }
2149566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2159566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
216c7a4214aSPierre Jolivet                   } else {
217c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */
218c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow);
219c7a4214aSPierre Jolivet                     }
220c7a4214aSPierre Jolivet                   }
221c7a4214aSPierre Jolivet                 }
222c7a4214aSPierre Jolivet                 if (m + A->rmap->n != nrow) {
223c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow, std::distance(it + A->rmap->n, idxr + nrow), idxr, idxc + m + A->rmap->n, ptr + (m + A->rmap->n) * nrow); /* vertical block E from above */
224c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
225b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
2269566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), ptr + (m + A->rmap->n) * nrow + m, &B));
2279566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
2289566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, ptr + m * nrow + m + A->rmap->n, &BT));
2299566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
230b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
2319566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
232c7a4214aSPierre Jolivet                     } else {
2337fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
2349566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
235c7a4214aSPierre Jolivet                     }
2369566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2379566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
238c7a4214aSPierre Jolivet                   } else {
239c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */
240c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(std::distance(it + A->rmap->n, idxr + nrow), 1, it + A->rmap->n, idxc + m + k, ptr + (m + k) * nrow + m + A->rmap->n);
241c7a4214aSPierre Jolivet                     }
242c7a4214aSPierre Jolivet                   }
243c7a4214aSPierre Jolivet                 }
244c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
245c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
246c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE;   /* rmap->rstart not in IS */
247c7a4214aSPierre Jolivet         } /* unsorted IS */
248c7a4214aSPierre Jolivet       }
249c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE;                                       /* different row and column IS */
250c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */
2519566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &idxr));
2529566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &idxc));
2539566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr));
25466b4df27SPierre Jolivet     PetscCall(MatShift((*submat)[i], shift));
25566b4df27SPierre Jolivet     PetscCall(MatScale((*submat)[i], scale));
256c7a4214aSPierre Jolivet   }
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258c7a4214aSPierre Jolivet }
259c7a4214aSPierre Jolivet 
MatDestroy_Htool(Mat A)260d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A)
261d71ae5a4SJacob Faibussowitsch {
262c9a4fd69SAudic XU   Mat_Htool               *a;
263c7a4214aSPierre Jolivet   PetscContainer           container;
264c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
265c7a4214aSPierre Jolivet 
266c7a4214aSPierre Jolivet   PetscFunctionBegin;
267c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
268f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr));
269f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr));
270f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr));
271f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr));
272f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr));
273f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr));
274f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr));
275f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr));
276f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr));
277c095613aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUseRecompression_C", nullptr));
2789566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container));
279c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2802a8381b2SBarry Smith     PetscCall(PetscContainerGetPointer(container, &kernelt));
2819566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
282f22e26b7SPierre Jolivet     PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr));
283c7a4214aSPierre Jolivet   }
28448a46eb9SPierre Jolivet   if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
286c7a4214aSPierre Jolivet   delete a->wrapper;
2871dd4f53aSPierre Jolivet   a->target_cluster.reset();
2881dd4f53aSPierre Jolivet   a->source_cluster.reset();
2891dd4f53aSPierre Jolivet   a->distributed_operator_holder.reset();
290c9a4fd69SAudic XU   PetscCall(PetscFree(a));
291c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable()
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293c7a4214aSPierre Jolivet }
294c7a4214aSPierre Jolivet 
MatView_Htool(Mat A,PetscViewer pv)295d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv)
296d71ae5a4SJacob Faibussowitsch {
297c9a4fd69SAudic XU   Mat_Htool                         *a;
29866b4df27SPierre Jolivet   PetscScalar                        shift, scale;
299c7a4214aSPierre Jolivet   PetscBool                          flg;
3001dd4f53aSPierre Jolivet   std::map<std::string, std::string> hmatrix_information;
301c7a4214aSPierre Jolivet 
302c7a4214aSPierre Jolivet   PetscFunctionBegin;
303c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
3041dd4f53aSPierre Jolivet   hmatrix_information = htool::get_distributed_hmatrix_information(a->distributed_operator_holder->hmatrix, PetscObjectComm((PetscObject)A));
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg));
306c7a4214aSPierre Jolivet   if (flg) {
30766b4df27SPierre Jolivet     PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
308c095613aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->distributed_operator_holder->block_diagonal_hmatrix->get_symmetry()));
30966b4df27SPierre Jolivet     if (PetscAbsScalar(scale - 1.0) > PETSC_MACHINE_EPSILON) {
310c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
31166b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(scale), (double)PetscImaginaryPart(scale)));
312c7a4214aSPierre Jolivet #else
31366b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)scale));
314c9a4fd69SAudic XU #endif
315c9a4fd69SAudic XU     }
31666b4df27SPierre Jolivet     if (PetscAbsScalar(shift) > PETSC_MACHINE_EPSILON) {
317c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX)
31866b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(shift), (double)PetscImaginaryPart(shift)));
319c9a4fd69SAudic XU #else
32066b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)shift));
321c7a4214aSPierre Jolivet #endif
322c7a4214aSPierre Jolivet     }
323c095613aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "maximal cluster leaf size: %" PetscInt_FMT "\n", a->max_cluster_leaf_size));
3249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon));
3259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta));
3269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0]));
3279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1]));
3289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor]));
3299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering]));
3301dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", hmatrix_information["Compression_ratio"].c_str()));
3311dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", hmatrix_information["Space_saving"].c_str()));
3321dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "block tree consistency: %s\n", PetscBools[a->distributed_operator_holder->hmatrix.is_block_tree_consistent()]));
333c095613aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "recompression: %s\n", PetscBools[a->recompression]));
3341dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "number of dense (resp. low rank) matrices: %s (resp. %s)\n", hmatrix_information["Number_of_dense_blocks"].c_str(), hmatrix_information["Number_of_low_rank_blocks"].c_str()));
3351dd4f53aSPierre Jolivet     PetscCall(
3361dd4f53aSPierre Jolivet       PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) dense block sizes: (%s, %s, %s)\n", hmatrix_information["Dense_block_size_min"].c_str(), hmatrix_information["Dense_block_size_mean"].c_str(), hmatrix_information["Dense_block_size_max"].c_str()));
3371dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) low rank block sizes: (%s, %s, %s)\n", hmatrix_information["Low_rank_block_size_min"].c_str(), hmatrix_information["Low_rank_block_size_mean"].c_str(),
3381dd4f53aSPierre Jolivet                                      hmatrix_information["Low_rank_block_size_max"].c_str()));
3391dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) ranks: (%s, %s, %s)\n", hmatrix_information["Rank_min"].c_str(), hmatrix_information["Rank_mean"].c_str(), hmatrix_information["Rank_max"].c_str()));
340c7a4214aSPierre Jolivet   }
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
342c7a4214aSPierre Jolivet }
343c7a4214aSPierre Jolivet 
344c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
MatGetRow_Htool(Mat A,PetscInt row,PetscInt * nz,PetscInt ** idx,PetscScalar ** v)345d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
346d71ae5a4SJacob Faibussowitsch {
347c9a4fd69SAudic XU   Mat_Htool   *a;
34866b4df27SPierre Jolivet   PetscScalar  shift, scale;
349c7a4214aSPierre Jolivet   PetscInt    *idxc;
350c7a4214aSPierre Jolivet   PetscBLASInt one = 1, bn;
351c7a4214aSPierre Jolivet 
352c7a4214aSPierre Jolivet   PetscFunctionBegin;
35366b4df27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
354c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
355c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
356c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, &idxc));
358c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i;
359c7a4214aSPierre Jolivet   }
360c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
361c7a4214aSPierre Jolivet   if (v) {
3629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, v));
363c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
364cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
3659566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N, &bn));
366ebd42cecSPierre Jolivet     PetscCallCXX(htool::Blas<PetscScalar>::scal(&bn, &scale, *v, &one));
36766b4df27SPierre Jolivet     if (row < A->cmap->N) (*v)[row] += shift;
368c7a4214aSPierre Jolivet   }
36948a46eb9SPierre Jolivet   if (!idx) PetscCall(PetscFree(idxc));
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
371c7a4214aSPierre Jolivet }
372c7a4214aSPierre Jolivet 
MatRestoreRow_Htool(Mat,PetscInt,PetscInt *,PetscInt ** idx,PetscScalar ** v)373f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v)
374d71ae5a4SJacob Faibussowitsch {
375c7a4214aSPierre Jolivet   PetscFunctionBegin;
37648a46eb9SPierre Jolivet   if (idx) PetscCall(PetscFree(*idx));
37748a46eb9SPierre Jolivet   if (v) PetscCall(PetscFree(*v));
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379c7a4214aSPierre Jolivet }
380c7a4214aSPierre Jolivet 
MatSetFromOptions_Htool(Mat A,PetscOptionItems PetscOptionsObject)381ce78bad3SBarry Smith static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems PetscOptionsObject)
382d71ae5a4SJacob Faibussowitsch {
383c9a4fd69SAudic XU   Mat_Htool *a;
384c7a4214aSPierre Jolivet   PetscInt   n;
385c7a4214aSPierre Jolivet   PetscBool  flg;
386c7a4214aSPierre Jolivet 
387c7a4214aSPierre Jolivet   PetscFunctionBegin;
388c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
389d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Htool options");
390c095613aSPierre Jolivet   PetscCall(PetscOptionsBoundedInt("-mat_htool_max_cluster_leaf_size", "Maximal leaf size in cluster tree", nullptr, a->max_cluster_leaf_size, &a->max_cluster_leaf_size, nullptr, 0));
3911dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr, 0.0));
392f22e26b7SPierre Jolivet   PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr));
3931dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr, 0));
3941dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", nullptr, a->depth[1], a->depth + 1, nullptr, 0));
3951dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBool("-mat_htool_block_tree_consistency", "Block tree consistency", nullptr, a->block_tree_consistency, &a->block_tree_consistency, nullptr));
396c095613aSPierre Jolivet   PetscCall(PetscOptionsBool("-mat_htool_recompression", "Use recompression", nullptr, a->recompression, &a->recompression, nullptr));
3971dd4f53aSPierre Jolivet 
398c7a4214aSPierre Jolivet   n = 0;
399dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg));
400c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
40198e73e17SPierre Jolivet   n = 0;
402dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg));
40398e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
404d0609cedSBarry Smith   PetscOptionsHeadEnd();
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406c7a4214aSPierre Jolivet }
407c7a4214aSPierre Jolivet 
MatAssemblyEnd_Htool(Mat A,MatAssemblyType)408cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType)
409d71ae5a4SJacob Faibussowitsch {
410c9a4fd69SAudic XU   Mat_Htool                                                           *a;
411c7a4214aSPierre Jolivet   const PetscInt                                                      *ranges;
412c7a4214aSPierre Jolivet   PetscInt                                                            *offset;
4131dd4f53aSPierre Jolivet   PetscMPIInt                                                          size, rank;
414b94d7dedSBarry Smith   char                                                                 S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U';
415cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                                *generator = nullptr;
4161dd4f53aSPierre Jolivet   htool::ClusterTreeBuilder<PetscReal>                                 recursive_build_strategy;
4171dd4f53aSPierre Jolivet   htool::Cluster<PetscReal>                                           *source_cluster;
418c095613aSPierre Jolivet   std::shared_ptr<htool::VirtualInternalLowRankGenerator<PetscScalar>> compressor;
419c7a4214aSPierre Jolivet 
420c7a4214aSPierre Jolivet   PetscFunctionBegin;
421*c895e7bbSPierre Jolivet   for (size_t i = 0; i < PETSC_STATIC_ARRAY_LENGTH(HtoolCite); ++i) PetscCall(PetscCitationsRegister(HtoolCitations[i], HtoolCite + i));
422c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
423c7a4214aSPierre Jolivet   delete a->wrapper;
4241dd4f53aSPierre Jolivet   a->target_cluster.reset();
4251dd4f53aSPierre Jolivet   a->source_cluster.reset();
4261dd4f53aSPierre Jolivet   a->distributed_operator_holder.reset();
4271dd4f53aSPierre Jolivet   // clustering
4289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &offset));
4309566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
431c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < size; ++i) {
432c7a4214aSPierre Jolivet     offset[2 * i]     = ranges[i];
433c7a4214aSPierre Jolivet     offset[2 * i + 1] = ranges[i + 1] - ranges[i];
434c7a4214aSPierre Jolivet   }
43598e73e17SPierre Jolivet   switch (a->clustering) {
436d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
437c095613aSPierre Jolivet     recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::GeometricSplitting<PetscReal>>>());
438d71ae5a4SJacob Faibussowitsch     break;
439d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
440c095613aSPierre Jolivet     recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::GeometricSplitting<PetscReal>>>());
441d71ae5a4SJacob Faibussowitsch     break;
442d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
443c095613aSPierre Jolivet     recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::RegularSplitting<PetscReal>>>());
444d71ae5a4SJacob Faibussowitsch     break;
445d71ae5a4SJacob Faibussowitsch   default:
446c095613aSPierre Jolivet     recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::RegularSplitting<PetscReal>>>());
44798e73e17SPierre Jolivet   }
448c095613aSPierre Jolivet   recursive_build_strategy.set_maximal_leaf_size(a->max_cluster_leaf_size);
449c095613aSPierre Jolivet   a->target_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree_from_local_partition(A->rmap->N, a->dim, a->gcoords_target, 2, size, offset));
450c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4519566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A, &ranges));
452c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < size; ++i) {
453c7a4214aSPierre Jolivet       offset[2 * i]     = ranges[i];
454c7a4214aSPierre Jolivet       offset[2 * i + 1] = ranges[i + 1] - ranges[i];
455c7a4214aSPierre Jolivet     }
45698e73e17SPierre Jolivet     switch (a->clustering) {
457d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
458c095613aSPierre Jolivet       recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::GeometricSplitting<PetscReal>>>());
459d71ae5a4SJacob Faibussowitsch       break;
460d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
461c095613aSPierre Jolivet       recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::GeometricSplitting<PetscReal>>>());
462d71ae5a4SJacob Faibussowitsch       break;
463d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
464c095613aSPierre Jolivet       recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::RegularSplitting<PetscReal>>>());
465d71ae5a4SJacob Faibussowitsch       break;
466d71ae5a4SJacob Faibussowitsch     default:
467c095613aSPierre Jolivet       recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::RegularSplitting<PetscReal>>>());
46898e73e17SPierre Jolivet     }
469c095613aSPierre Jolivet     recursive_build_strategy.set_maximal_leaf_size(a->max_cluster_leaf_size);
470c095613aSPierre Jolivet     a->source_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree_from_local_partition(A->cmap->N, a->dim, a->gcoords_source, 2, size, offset));
471c7a4214aSPierre Jolivet     S = uplo       = 'N';
4721dd4f53aSPierre Jolivet     source_cluster = a->source_cluster.get();
4731dd4f53aSPierre Jolivet   } else source_cluster = a->target_cluster.get();
4749566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
4751dd4f53aSPierre Jolivet   // generator
4761dd4f53aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx);
4771dd4f53aSPierre Jolivet   else {
4781dd4f53aSPierre Jolivet     a->wrapper = nullptr;
4791dd4f53aSPierre Jolivet     generator  = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx);
4801dd4f53aSPierre Jolivet   }
4811dd4f53aSPierre Jolivet   // compressor
482c7a4214aSPierre Jolivet   switch (a->compressor) {
483d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
484c095613aSPierre Jolivet     compressor = std::make_shared<htool::fullACA<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data());
485d71ae5a4SJacob Faibussowitsch     break;
486d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_SVD:
487c095613aSPierre Jolivet     compressor = std::make_shared<htool::SVD<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data());
488d71ae5a4SJacob Faibussowitsch     break;
489d71ae5a4SJacob Faibussowitsch   default:
490c095613aSPierre Jolivet     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data());
491c7a4214aSPierre Jolivet   }
4921dd4f53aSPierre Jolivet   // local hierarchical matrix
4931dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
494c095613aSPierre Jolivet   auto hmatrix_builder = htool::HMatrixTreeBuilder<PetscScalar>(a->epsilon, a->eta, S, uplo);
495c095613aSPierre Jolivet   if (a->recompression) {
496c095613aSPierre Jolivet     std::shared_ptr<htool::VirtualInternalLowRankGenerator<PetscScalar>> RecompressedLowRankGenerator = std::make_shared<htool::RecompressedLowRankGenerator<PetscScalar>>(*compressor, std::function<void(htool::LowRankMatrix<PetscScalar> &)>(htool::SVD_recompression<PetscScalar>));
497c095613aSPierre Jolivet     hmatrix_builder.set_low_rank_generator(RecompressedLowRankGenerator);
498c095613aSPierre Jolivet   } else {
4991dd4f53aSPierre Jolivet     hmatrix_builder.set_low_rank_generator(compressor);
500c095613aSPierre Jolivet   }
5011dd4f53aSPierre Jolivet   hmatrix_builder.set_minimal_target_depth(a->depth[0]);
5021dd4f53aSPierre Jolivet   hmatrix_builder.set_minimal_source_depth(a->depth[1]);
5031dd4f53aSPierre Jolivet   PetscCheck(a->block_tree_consistency || (!a->block_tree_consistency && !(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE)), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot have a MatHtool with inconsistent block tree which is either symmetric or Hermitian");
5041dd4f53aSPierre Jolivet   hmatrix_builder.set_block_tree_consistency(a->block_tree_consistency);
505c095613aSPierre Jolivet   a->distributed_operator_holder = std::make_unique<htool::DefaultApproximationBuilder<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, *a->target_cluster, *source_cluster, hmatrix_builder, PetscObjectComm((PetscObject)A));
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
507c7a4214aSPierre Jolivet }
508c7a4214aSPierre Jolivet 
MatProductNumeric_Htool(Mat C)509d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C)
510d71ae5a4SJacob Faibussowitsch {
511c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
512c9a4fd69SAudic XU   Mat_Htool         *a;
513c7a4214aSPierre Jolivet   const PetscScalar *in;
514c7a4214aSPierre Jolivet   PetscScalar       *out;
5158b8ddd11SPierre Jolivet   PetscInt           N, lda;
516c7a4214aSPierre Jolivet 
517c7a4214aSPierre Jolivet   PetscFunctionBegin;
518c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
519f22e26b7SPierre Jolivet   PetscCall(MatGetSize(C, nullptr, &N));
5209566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &lda));
5215f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
5229566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B, &in));
5239566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &out));
524c9a4fd69SAudic XU   PetscCall(MatShellGetContext(product->A, &a));
5258b8ddd11SPierre Jolivet   switch (product->type) {
526d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AB:
527c095613aSPierre Jolivet     if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_matrix_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr);
528c095613aSPierre Jolivet     else htool::internal_add_distributed_operator_matrix_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr);
529d71ae5a4SJacob Faibussowitsch     break;
530d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AtB:
531c095613aSPierre Jolivet     if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_matrix_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr);
532c095613aSPierre Jolivet     else htool::internal_add_distributed_operator_matrix_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr);
533d71ae5a4SJacob Faibussowitsch     break;
534d71ae5a4SJacob Faibussowitsch   default:
535d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]);
536c7a4214aSPierre Jolivet   }
5379566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &out));
5389566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B, &in));
5393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
540c7a4214aSPierre Jolivet }
541c7a4214aSPierre Jolivet 
MatProductSymbolic_Htool(Mat C)542d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C)
543d71ae5a4SJacob Faibussowitsch {
544c7a4214aSPierre Jolivet   Mat_Product *product = C->product;
545c7a4214aSPierre Jolivet   Mat          A, B;
546c7a4214aSPierre Jolivet   PetscBool    flg;
547c7a4214aSPierre Jolivet 
548c7a4214aSPierre Jolivet   PetscFunctionBegin;
549c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
550c7a4214aSPierre Jolivet   A = product->A;
551c7a4214aSPierre Jolivet   B = product->B;
5529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, ""));
55397713f8cSPierre Jolivet   PetscCheck(flg && (product->type == MATPRODUCT_AB || product->type == MATPRODUCT_AtB), PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "ProductType %s not supported for %s", MatProductTypes[product->type], ((PetscObject)product->B)->type_name);
55497713f8cSPierre Jolivet   if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
55597713f8cSPierre Jolivet     if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
55697713f8cSPierre Jolivet     else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
5578b8ddd11SPierre Jolivet   }
5589566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
5599566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
5609566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
5619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
5629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
563f22e26b7SPierre Jolivet   C->ops->productsymbolic = nullptr;
564c7a4214aSPierre Jolivet   C->ops->productnumeric  = MatProductNumeric_Htool;
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
566c7a4214aSPierre Jolivet }
567c7a4214aSPierre Jolivet 
MatProductSetFromOptions_Htool(Mat C)568d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
569d71ae5a4SJacob Faibussowitsch {
570c7a4214aSPierre Jolivet   PetscFunctionBegin;
571c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
5728b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
574c7a4214aSPierre Jolivet }
575c7a4214aSPierre Jolivet 
MatHtoolGetHierarchicalMat_Htool(Mat A,const htool::DistributedOperator<PetscScalar> ** distributed_operator)5761dd4f53aSPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator)
577d71ae5a4SJacob Faibussowitsch {
578c9a4fd69SAudic XU   Mat_Htool *a;
579c7a4214aSPierre Jolivet 
580c7a4214aSPierre Jolivet   PetscFunctionBegin;
581c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
5821dd4f53aSPierre Jolivet   *distributed_operator = &a->distributed_operator_holder->distributed_operator;
5833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
584c7a4214aSPierre Jolivet }
585c7a4214aSPierre Jolivet 
586c7a4214aSPierre Jolivet /*@C
58711a5261eSBarry Smith   MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`.
588c7a4214aSPierre Jolivet 
589cc4c1da9SBarry Smith   No Fortran Support, No C Support
590cc4c1da9SBarry Smith 
591c7a4214aSPierre Jolivet   Input Parameter:
592c7a4214aSPierre Jolivet . A - hierarchical matrix
593c7a4214aSPierre Jolivet 
594c7a4214aSPierre Jolivet   Output Parameter:
5951dd4f53aSPierre Jolivet . distributed_operator - opaque pointer to a Htool virtual matrix
596c7a4214aSPierre Jolivet 
597c7a4214aSPierre Jolivet   Level: advanced
598c7a4214aSPierre Jolivet 
5991cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`
600c7a4214aSPierre Jolivet @*/
MatHtoolGetHierarchicalMat(Mat A,const htool::DistributedOperator<PetscScalar> ** distributed_operator)6011dd4f53aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator)
602d71ae5a4SJacob Faibussowitsch {
603c7a4214aSPierre Jolivet   PetscFunctionBegin;
604c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6051dd4f53aSPierre Jolivet   PetscAssertPointer(distributed_operator, 2);
6061dd4f53aSPierre Jolivet   PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::DistributedOperator<PetscScalar> **), (A, distributed_operator));
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
608c7a4214aSPierre Jolivet }
609c7a4214aSPierre Jolivet 
MatHtoolSetKernel_Htool(Mat A,MatHtoolKernelFn * kernel,void * kernelctx)6108434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
611d71ae5a4SJacob Faibussowitsch {
612c9a4fd69SAudic XU   Mat_Htool *a;
613c7a4214aSPierre Jolivet 
614c7a4214aSPierre Jolivet   PetscFunctionBegin;
615c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
616c7a4214aSPierre Jolivet   a->kernel    = kernel;
617c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
618c7a4214aSPierre Jolivet   delete a->wrapper;
6191dd4f53aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx);
6203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
621c7a4214aSPierre Jolivet }
622c7a4214aSPierre Jolivet 
623c7a4214aSPierre Jolivet /*@C
62411a5261eSBarry Smith   MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`.
625c7a4214aSPierre Jolivet 
626cc4c1da9SBarry Smith   Collective, No Fortran Support
627cc4c1da9SBarry Smith 
628c7a4214aSPierre Jolivet   Input Parameters:
629c7a4214aSPierre Jolivet + A         - hierarchical matrix
6302ef1f0ffSBarry Smith . kernel    - computational kernel (or `NULL`)
6312ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
632c7a4214aSPierre Jolivet 
633c7a4214aSPierre Jolivet   Level: advanced
634c7a4214aSPierre Jolivet 
6351cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()`
636c7a4214aSPierre Jolivet @*/
MatHtoolSetKernel(Mat A,MatHtoolKernelFn * kernel,void * kernelctx)637cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
638d71ae5a4SJacob Faibussowitsch {
639c7a4214aSPierre Jolivet   PetscFunctionBegin;
640c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
641c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 2);
6424f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 3);
6438434afd1SBarry Smith   PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx));
6443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
645c7a4214aSPierre Jolivet }
646c7a4214aSPierre Jolivet 
MatHtoolGetPermutationSource_Htool(Mat A,IS * is)647d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is)
648d71ae5a4SJacob Faibussowitsch {
649c9a4fd69SAudic XU   Mat_Htool                       *a;
6501dd4f53aSPierre Jolivet   PetscMPIInt                      rank;
6511dd4f53aSPierre Jolivet   const std::vector<PetscInt>     *source;
6521dd4f53aSPierre Jolivet   const htool::Cluster<PetscReal> *local_source_cluster;
65398e73e17SPierre Jolivet 
65498e73e17SPierre Jolivet   PetscFunctionBegin;
655c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
6561dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
6571dd4f53aSPierre Jolivet   local_source_cluster = a->source_cluster ? &a->source_cluster->get_cluster_on_partition(rank) : &a->target_cluster->get_cluster_on_partition(rank);
6581dd4f53aSPierre Jolivet   source               = &local_source_cluster->get_permutation();
6591dd4f53aSPierre Jolivet   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), local_source_cluster->get_size(), source->data() + local_source_cluster->get_offset(), PETSC_COPY_VALUES, is));
6609566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66298e73e17SPierre Jolivet }
66398e73e17SPierre Jolivet 
664cc4c1da9SBarry Smith /*@
66511a5261eSBarry Smith   MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix.
66698e73e17SPierre Jolivet 
66798e73e17SPierre Jolivet   Input Parameter:
66898e73e17SPierre Jolivet . A - hierarchical matrix
66998e73e17SPierre Jolivet 
67098e73e17SPierre Jolivet   Output Parameter:
67198e73e17SPierre Jolivet . is - permutation
67298e73e17SPierre Jolivet 
67398e73e17SPierre Jolivet   Level: advanced
67498e73e17SPierre Jolivet 
6751cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()`
67698e73e17SPierre Jolivet @*/
MatHtoolGetPermutationSource(Mat A,IS * is)677cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is)
678d71ae5a4SJacob Faibussowitsch {
67998e73e17SPierre Jolivet   PetscFunctionBegin;
68098e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6814f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
682cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is));
6833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68498e73e17SPierre Jolivet }
68598e73e17SPierre Jolivet 
MatHtoolGetPermutationTarget_Htool(Mat A,IS * is)686d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is)
687d71ae5a4SJacob Faibussowitsch {
688c9a4fd69SAudic XU   Mat_Htool                   *a;
6891dd4f53aSPierre Jolivet   const std::vector<PetscInt> *target;
6901dd4f53aSPierre Jolivet   PetscMPIInt                  rank;
69198e73e17SPierre Jolivet 
69298e73e17SPierre Jolivet   PetscFunctionBegin;
693c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
6941dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
6951dd4f53aSPierre Jolivet   target = &a->target_cluster->get_permutation();
6961dd4f53aSPierre Jolivet   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), a->target_cluster->get_cluster_on_partition(rank).get_size(), target->data() + a->target_cluster->get_cluster_on_partition(rank).get_offset(), PETSC_COPY_VALUES, is));
6979566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69998e73e17SPierre Jolivet }
70098e73e17SPierre Jolivet 
701cc4c1da9SBarry Smith /*@
70211a5261eSBarry Smith   MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix.
70398e73e17SPierre Jolivet 
70498e73e17SPierre Jolivet   Input Parameter:
70598e73e17SPierre Jolivet . A - hierarchical matrix
70698e73e17SPierre Jolivet 
70798e73e17SPierre Jolivet   Output Parameter:
70898e73e17SPierre Jolivet . is - permutation
70998e73e17SPierre Jolivet 
71098e73e17SPierre Jolivet   Level: advanced
71198e73e17SPierre Jolivet 
7121cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()`
71398e73e17SPierre Jolivet @*/
MatHtoolGetPermutationTarget(Mat A,IS * is)714cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is)
715d71ae5a4SJacob Faibussowitsch {
71698e73e17SPierre Jolivet   PetscFunctionBegin;
71798e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
7184f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
719cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is));
7203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72198e73e17SPierre Jolivet }
72298e73e17SPierre Jolivet 
MatHtoolUsePermutation_Htool(Mat A,PetscBool use)723d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use)
724d71ae5a4SJacob Faibussowitsch {
725c9a4fd69SAudic XU   Mat_Htool *a;
72698e73e17SPierre Jolivet 
72798e73e17SPierre Jolivet   PetscFunctionBegin;
728c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
729c095613aSPierre Jolivet   a->permutation = use;
7303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73198e73e17SPierre Jolivet }
73298e73e17SPierre Jolivet 
733cc4c1da9SBarry Smith /*@
73411a5261eSBarry Smith   MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation.
73598e73e17SPierre Jolivet 
73698e73e17SPierre Jolivet   Input Parameters:
73798e73e17SPierre Jolivet + A   - hierarchical matrix
73898e73e17SPierre Jolivet - use - Boolean value
73998e73e17SPierre Jolivet 
74098e73e17SPierre Jolivet   Level: advanced
74198e73e17SPierre Jolivet 
7421cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()`
74398e73e17SPierre Jolivet @*/
MatHtoolUsePermutation(Mat A,PetscBool use)744cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use)
745d71ae5a4SJacob Faibussowitsch {
74698e73e17SPierre Jolivet   PetscFunctionBegin;
74798e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
74898e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A, use, 2);
749cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use));
7503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75198e73e17SPierre Jolivet }
75298e73e17SPierre Jolivet 
MatHtoolUseRecompression_Htool(Mat A,PetscBool use)753c095613aSPierre Jolivet static PetscErrorCode MatHtoolUseRecompression_Htool(Mat A, PetscBool use)
754c095613aSPierre Jolivet {
755c095613aSPierre Jolivet   Mat_Htool *a;
756c095613aSPierre Jolivet 
757c095613aSPierre Jolivet   PetscFunctionBegin;
758c095613aSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
759c095613aSPierre Jolivet   a->recompression = use;
760c095613aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
761c095613aSPierre Jolivet }
762c095613aSPierre Jolivet 
763c095613aSPierre Jolivet /*@
764c095613aSPierre Jolivet   MatHtoolUseRecompression - Sets whether a `MATHTOOL` matrix should use recompression.
765c095613aSPierre Jolivet 
766c095613aSPierre Jolivet   Input Parameters:
767c095613aSPierre Jolivet + A   - hierarchical matrix
768c095613aSPierre Jolivet - use - Boolean value
769c095613aSPierre Jolivet 
770c095613aSPierre Jolivet   Level: advanced
771c095613aSPierre Jolivet 
772c095613aSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MATHTOOL`
773c095613aSPierre Jolivet @*/
MatHtoolUseRecompression(Mat A,PetscBool use)774c095613aSPierre Jolivet PetscErrorCode MatHtoolUseRecompression(Mat A, PetscBool use)
775c095613aSPierre Jolivet {
776c095613aSPierre Jolivet   PetscFunctionBegin;
777c095613aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
778c095613aSPierre Jolivet   PetscValidLogicalCollectiveBool(A, use, 2);
779c095613aSPierre Jolivet   PetscTryMethod(A, "MatHtoolUseRecompression_C", (Mat, PetscBool), (A, use));
780c095613aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
781c095613aSPierre Jolivet }
782c095613aSPierre Jolivet 
MatConvert_Htool_Dense(Mat A,MatType,MatReuse reuse,Mat * B)783cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B)
784d71ae5a4SJacob Faibussowitsch {
785c7a4214aSPierre Jolivet   Mat          C;
786c9a4fd69SAudic XU   Mat_Htool   *a;
78766b4df27SPierre Jolivet   PetscScalar *array, shift, scale;
788c7a4214aSPierre Jolivet   PetscInt     lda;
789c7a4214aSPierre Jolivet 
790c7a4214aSPierre Jolivet   PetscFunctionBegin;
79166b4df27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
792c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
793c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
794c7a4214aSPierre Jolivet     C = *B;
7955f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions");
7969566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C, &lda));
7975f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
798c7a4214aSPierre Jolivet   } else {
7999566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
8009566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
8019566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, MATDENSE));
8029566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
8039566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
804c7a4214aSPierre Jolivet   }
8059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
8069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
8079566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
8081dd4f53aSPierre Jolivet   htool::copy_to_dense_in_user_numbering(a->distributed_operator_holder->hmatrix, array);
8099566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
81066b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
81166b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
812ac530a7eSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &C));
813ac530a7eSPierre Jolivet   else *B = C;
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
815c7a4214aSPierre Jolivet }
816c7a4214aSPierre Jolivet 
GenEntriesTranspose(PetscInt sdim,PetscInt M,PetscInt N,const PetscInt * rows,const PetscInt * cols,PetscScalar * ptr,PetscCtx ctx)8172a8381b2SBarry Smith static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, PetscCtx ctx)
818d71ae5a4SJacob Faibussowitsch {
819c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx;
82098e73e17SPierre Jolivet   PetscScalar             *tmp;
82198e73e17SPierre Jolivet 
82298e73e17SPierre Jolivet   PetscFunctionBegin;
8233ba16761SJacob Faibussowitsch   PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx));
8249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M * N, &tmp));
8259566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp, ptr, M * N));
82698e73e17SPierre Jolivet   for (PetscInt i = 0; i < M; ++i) {
82798e73e17SPierre Jolivet     for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N];
82898e73e17SPierre Jolivet   }
8299566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
831c7a4214aSPierre Jolivet }
832c7a4214aSPierre Jolivet 
833c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
MatTranspose_Htool(Mat A,MatReuse reuse,Mat * B)834d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B)
835d71ae5a4SJacob Faibussowitsch {
836c7a4214aSPierre Jolivet   Mat                      C;
837c9a4fd69SAudic XU   Mat_Htool               *a, *c;
83866b4df27SPierre Jolivet   PetscScalar              shift, scale;
839c7a4214aSPierre Jolivet   PetscInt                 M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n;
840c7a4214aSPierre Jolivet   PetscContainer           container;
841c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
842c7a4214aSPierre Jolivet 
843c7a4214aSPierre Jolivet   PetscFunctionBegin;
84466b4df27SPierre Jolivet   PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED));
845c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
8467fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
8475f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported");
848c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8499566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
8509566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, n, m, N, M));
8519566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
8529566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
8539566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
85449abdd8aSBarry Smith     PetscCall(PetscObjectContainerCompose((PetscObject)C, "KernelTranspose", kernelt, PetscCtxDestroyDefault));
855c7a4214aSPierre Jolivet   } else {
856c7a4214aSPierre Jolivet     C = *B;
8579566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container));
8585f80ce2aSJacob Faibussowitsch     PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first");
8592a8381b2SBarry Smith     PetscCall(PetscContainerGetPointer(container, &kernelt));
860c7a4214aSPierre Jolivet   }
861c9a4fd69SAudic XU   PetscCall(MatShellGetContext(C, &c));
862c7a4214aSPierre Jolivet   c->dim = a->dim;
86366b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
86466b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
86598e73e17SPierre Jolivet   c->kernel = GenEntriesTranspose;
866c7a4214aSPierre Jolivet   if (kernelt->A != A) {
8679566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
868c7a4214aSPierre Jolivet     kernelt->A = A;
8699566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
870c7a4214aSPierre Jolivet   }
871c7a4214aSPierre Jolivet   kernelt->kernel           = a->kernel;
872c7a4214aSPierre Jolivet   kernelt->kernelctx        = a->kernelctx;
873c7a4214aSPierre Jolivet   c->kernelctx              = kernelt;
874c095613aSPierre Jolivet   c->max_cluster_leaf_size  = a->max_cluster_leaf_size;
875de647125SPierre Marchand   c->epsilon                = a->epsilon;
876de647125SPierre Marchand   c->eta                    = a->eta;
877de647125SPierre Marchand   c->block_tree_consistency = a->block_tree_consistency;
878c095613aSPierre Jolivet   c->permutation            = a->permutation;
879c095613aSPierre Jolivet   c->recompression          = a->recompression;
880de647125SPierre Marchand   c->compressor             = a->compressor;
881de647125SPierre Marchand   c->clustering             = a->clustering;
882c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target));
8849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim));
885c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
8869566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source));
8879566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim));
888c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
889c7a4214aSPierre Jolivet   }
8909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
8919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
892c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
8933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
894c7a4214aSPierre Jolivet }
895c7a4214aSPierre Jolivet 
MatDestroy_Factor(Mat F)8961dd4f53aSPierre Jolivet static PetscErrorCode MatDestroy_Factor(Mat F)
8971dd4f53aSPierre Jolivet {
8981dd4f53aSPierre Jolivet   PetscContainer               container;
8991dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *A;
9001dd4f53aSPierre Jolivet 
9011dd4f53aSPierre Jolivet   PetscFunctionBegin;
9021dd4f53aSPierre Jolivet   PetscCall(PetscObjectQuery((PetscObject)F, "HMatrix", (PetscObject *)&container));
9031dd4f53aSPierre Jolivet   if (container) {
9042a8381b2SBarry Smith     PetscCall(PetscContainerGetPointer(container, &A));
9051dd4f53aSPierre Jolivet     delete A;
9061dd4f53aSPierre Jolivet     PetscCall(PetscObjectCompose((PetscObject)F, "HMatrix", nullptr));
9071dd4f53aSPierre Jolivet   }
9081dd4f53aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", nullptr));
9091dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9101dd4f53aSPierre Jolivet }
9111dd4f53aSPierre Jolivet 
MatFactorGetSolverType_Htool(Mat,MatSolverType * type)9121dd4f53aSPierre Jolivet static PetscErrorCode MatFactorGetSolverType_Htool(Mat, MatSolverType *type)
9131dd4f53aSPierre Jolivet {
9141dd4f53aSPierre Jolivet   PetscFunctionBegin;
9151dd4f53aSPierre Jolivet   *type = MATSOLVERHTOOL;
9161dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9171dd4f53aSPierre Jolivet }
9181dd4f53aSPierre Jolivet 
9191dd4f53aSPierre Jolivet template <char trans>
MatSolve_Private(Mat A,htool::Matrix<PetscScalar> & X)9201dd4f53aSPierre Jolivet static inline PetscErrorCode MatSolve_Private(Mat A, htool::Matrix<PetscScalar> &X)
9211dd4f53aSPierre Jolivet {
9221dd4f53aSPierre Jolivet   PetscContainer               container;
9231dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *B;
9241dd4f53aSPierre Jolivet 
9251dd4f53aSPierre Jolivet   PetscFunctionBegin;
9261dd4f53aSPierre Jolivet   PetscCheck(A->factortype == MAT_FACTOR_LU || A->factortype == MAT_FACTOR_CHOLESKY, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_UNKNOWN_TYPE, "Only MAT_LU_FACTOR and MAT_CHOLESKY_FACTOR are supported");
9271dd4f53aSPierre Jolivet   PetscCall(PetscObjectQuery((PetscObject)A, "HMatrix", (PetscObject *)&container));
9281dd4f53aSPierre Jolivet   PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call Mat%sFactorNumeric() before Mat%sSolve%s()", A->factortype == MAT_FACTOR_LU ? "LU" : "Cholesky", X.nb_cols() == 1 ? "" : "Mat", trans == 'N' ? "" : "Transpose");
9292a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(container, &B));
9301dd4f53aSPierre Jolivet   if (A->factortype == MAT_FACTOR_LU) htool::lu_solve(trans, *B, X);
9311dd4f53aSPierre Jolivet   else htool::cholesky_solve('L', *B, X);
9321dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9331dd4f53aSPierre Jolivet }
9341dd4f53aSPierre Jolivet 
9351dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type * = nullptr>
MatSolve_Htool(Mat A,Type b,Type x)9361dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type b, Type x)
9371dd4f53aSPierre Jolivet {
9381dd4f53aSPierre Jolivet   PetscInt                   n;
9391dd4f53aSPierre Jolivet   htool::Matrix<PetscScalar> v;
9401dd4f53aSPierre Jolivet   PetscScalar               *array;
9411dd4f53aSPierre Jolivet 
9421dd4f53aSPierre Jolivet   PetscFunctionBegin;
9431dd4f53aSPierre Jolivet   PetscCall(VecGetLocalSize(b, &n));
9441dd4f53aSPierre Jolivet   PetscCall(VecCopy(b, x));
9451dd4f53aSPierre Jolivet   PetscCall(VecGetArrayWrite(x, &array));
9461dd4f53aSPierre Jolivet   v.assign(n, 1, array, false);
9471dd4f53aSPierre Jolivet   PetscCall(VecRestoreArrayWrite(x, &array));
9481dd4f53aSPierre Jolivet   PetscCall(MatSolve_Private<trans>(A, v));
9491dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9501dd4f53aSPierre Jolivet }
9511dd4f53aSPierre Jolivet 
9521dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type * = nullptr>
MatSolve_Htool(Mat A,Type B,Type X)9531dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type B, Type X)
9541dd4f53aSPierre Jolivet {
9551dd4f53aSPierre Jolivet   PetscInt                   m, N;
9561dd4f53aSPierre Jolivet   htool::Matrix<PetscScalar> v;
9571dd4f53aSPierre Jolivet   PetscScalar               *array;
9581dd4f53aSPierre Jolivet 
9591dd4f53aSPierre Jolivet   PetscFunctionBegin;
9601dd4f53aSPierre Jolivet   PetscCall(MatGetLocalSize(B, &m, nullptr));
9611dd4f53aSPierre Jolivet   PetscCall(MatGetLocalSize(B, nullptr, &N));
9621dd4f53aSPierre Jolivet   PetscCall(MatCopy(B, X, SAME_NONZERO_PATTERN));
9631dd4f53aSPierre Jolivet   PetscCall(MatDenseGetArrayWrite(X, &array));
9641dd4f53aSPierre Jolivet   v.assign(m, N, array, false);
9651dd4f53aSPierre Jolivet   PetscCall(MatDenseRestoreArrayWrite(X, &array));
9661dd4f53aSPierre Jolivet   PetscCall(MatSolve_Private<trans>(A, v));
9671dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9681dd4f53aSPierre Jolivet }
9691dd4f53aSPierre Jolivet 
9701dd4f53aSPierre Jolivet template <MatFactorType ftype>
MatFactorNumeric_Htool(Mat F,Mat A,const MatFactorInfo *)9711dd4f53aSPierre Jolivet static PetscErrorCode MatFactorNumeric_Htool(Mat F, Mat A, const MatFactorInfo *)
9721dd4f53aSPierre Jolivet {
9731dd4f53aSPierre Jolivet   Mat_Htool                   *a;
9741dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *B;
9751dd4f53aSPierre Jolivet 
9761dd4f53aSPierre Jolivet   PetscFunctionBegin;
9771dd4f53aSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
9781dd4f53aSPierre Jolivet   B = new htool::HMatrix<PetscScalar>(a->distributed_operator_holder->hmatrix);
979c095613aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) htool::sequential_lu_factorization(*B);
980c095613aSPierre Jolivet   else htool::sequential_cholesky_factorization('L', *B);
9811dd4f53aSPierre Jolivet   PetscCall(PetscObjectContainerCompose((PetscObject)F, "HMatrix", B, nullptr));
9821dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9831dd4f53aSPierre Jolivet }
9841dd4f53aSPierre Jolivet 
9851dd4f53aSPierre Jolivet template <MatFactorType ftype>
MatFactorSymbolic_Htool(Mat F,Mat)9861dd4f53aSPierre Jolivet PetscErrorCode MatFactorSymbolic_Htool(Mat F, Mat)
9871dd4f53aSPierre Jolivet {
9881dd4f53aSPierre Jolivet   PetscFunctionBegin;
9891dd4f53aSPierre Jolivet   F->preallocated  = PETSC_TRUE;
9901dd4f53aSPierre Jolivet   F->assembled     = PETSC_TRUE;
9911dd4f53aSPierre Jolivet   F->ops->solve    = MatSolve_Htool<'N', Vec>;
9921dd4f53aSPierre Jolivet   F->ops->matsolve = MatSolve_Htool<'N', Mat>;
9931dd4f53aSPierre Jolivet   if (!PetscDefined(USE_COMPLEX) || ftype == MAT_FACTOR_LU) {
9941dd4f53aSPierre Jolivet     F->ops->solvetranspose    = MatSolve_Htool<'T', Vec>;
9951dd4f53aSPierre Jolivet     F->ops->matsolvetranspose = MatSolve_Htool<'T', Mat>;
9961dd4f53aSPierre Jolivet   }
9971dd4f53aSPierre Jolivet   F->ops->destroy = MatDestroy_Factor;
9981dd4f53aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_LU>;
9991dd4f53aSPierre Jolivet   else F->ops->choleskyfactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_CHOLESKY>;
10001dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10011dd4f53aSPierre Jolivet }
10021dd4f53aSPierre Jolivet 
MatLUFactorSymbolic_Htool(Mat F,Mat A,IS,IS,const MatFactorInfo *)10031dd4f53aSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Htool(Mat F, Mat A, IS, IS, const MatFactorInfo *)
10041dd4f53aSPierre Jolivet {
10051dd4f53aSPierre Jolivet   PetscFunctionBegin;
10061dd4f53aSPierre Jolivet   PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_LU>(F, A));
10071dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10081dd4f53aSPierre Jolivet }
10091dd4f53aSPierre Jolivet 
MatCholeskyFactorSymbolic_Htool(Mat F,Mat A,IS,const MatFactorInfo *)10101dd4f53aSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Htool(Mat F, Mat A, IS, const MatFactorInfo *)
10111dd4f53aSPierre Jolivet {
10121dd4f53aSPierre Jolivet   PetscFunctionBegin;
10131dd4f53aSPierre Jolivet   PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_CHOLESKY>(F, A));
10141dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10151dd4f53aSPierre Jolivet }
10161dd4f53aSPierre Jolivet 
MatGetFactor_htool_htool(Mat A,MatFactorType ftype,Mat * F)10171dd4f53aSPierre Jolivet static PetscErrorCode MatGetFactor_htool_htool(Mat A, MatFactorType ftype, Mat *F)
10181dd4f53aSPierre Jolivet {
10191dd4f53aSPierre Jolivet   Mat         B;
10201dd4f53aSPierre Jolivet   Mat_Htool  *a;
10211dd4f53aSPierre Jolivet   PetscMPIInt size;
10221dd4f53aSPierre Jolivet 
10231dd4f53aSPierre Jolivet   PetscFunctionBegin;
10241dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
10251dd4f53aSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
10261dd4f53aSPierre Jolivet   PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel MatGetFactor()");
10271dd4f53aSPierre Jolivet   PetscCheck(a->block_tree_consistency, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot factor a MatHtool with inconsistent block tree");
10281dd4f53aSPierre Jolivet   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
10291dd4f53aSPierre Jolivet   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
10301dd4f53aSPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &((PetscObject)B)->type_name));
10311dd4f53aSPierre Jolivet   PetscCall(MatSetUp(B));
10321dd4f53aSPierre Jolivet 
10331dd4f53aSPierre Jolivet   B->ops->getinfo    = MatGetInfo_External;
10341dd4f53aSPierre Jolivet   B->factortype      = ftype;
10351dd4f53aSPierre Jolivet   B->trivialsymbolic = PETSC_TRUE;
10361dd4f53aSPierre Jolivet 
10371dd4f53aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_Htool;
10381dd4f53aSPierre Jolivet   else if (ftype == MAT_FACTOR_CHOLESKY) B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_Htool;
10391dd4f53aSPierre Jolivet 
10401dd4f53aSPierre Jolivet   PetscCall(PetscFree(B->solvertype));
10411dd4f53aSPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &B->solvertype));
10421dd4f53aSPierre Jolivet 
10431dd4f53aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_Htool));
10441dd4f53aSPierre Jolivet   *F = B;
10451dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10461dd4f53aSPierre Jolivet }
10471dd4f53aSPierre Jolivet 
MatSolverTypeRegister_Htool(void)10481dd4f53aSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Htool(void)
10491dd4f53aSPierre Jolivet {
10501dd4f53aSPierre Jolivet   PetscFunctionBegin;
10511dd4f53aSPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_LU, MatGetFactor_htool_htool));
10521dd4f53aSPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_CHOLESKY, MatGetFactor_htool_htool));
10531dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10541dd4f53aSPierre Jolivet }
10551dd4f53aSPierre Jolivet 
1056c7a4214aSPierre Jolivet /*@C
105711a5261eSBarry Smith   MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel.
1058c7a4214aSPierre Jolivet 
1059cc4c1da9SBarry Smith   Collective, No Fortran Support
1060cc4c1da9SBarry Smith 
1061c7a4214aSPierre Jolivet   Input Parameters:
1062c7a4214aSPierre Jolivet + comm          - MPI communicator
10632ef1f0ffSBarry Smith . m             - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
10642ef1f0ffSBarry Smith . n             - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
10652ef1f0ffSBarry Smith . M             - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
10662ef1f0ffSBarry Smith . N             - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1067c7a4214aSPierre Jolivet . spacedim      - dimension of the space coordinates
1068c7a4214aSPierre Jolivet . coords_target - coordinates of the target
1069c7a4214aSPierre Jolivet . coords_source - coordinates of the source
10702ef1f0ffSBarry Smith . kernel        - computational kernel (or `NULL`)
10712ef1f0ffSBarry Smith - kernelctx     - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
1072c7a4214aSPierre Jolivet 
1073c7a4214aSPierre Jolivet   Output Parameter:
1074c7a4214aSPierre Jolivet . B - matrix
1075c7a4214aSPierre Jolivet 
1076c7a4214aSPierre Jolivet   Options Database Keys:
1077c095613aSPierre Jolivet + -mat_htool_max_cluster_leaf_size <`PetscInt`>                                                - maximal leaf size in cluster tree
107811a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`>                                                             - relative error in Frobenius norm when approximating a block
107911a5261eSBarry Smith . -mat_htool_eta <`PetscReal`>                                                                 - admissibility condition tolerance
108011a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the rows
108111a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the columns
10821dd4f53aSPierre Jolivet . -mat_htool_block_tree_consistency <`PetscBool`>                                              - block tree consistency
1083c095613aSPierre Jolivet . -mat_htool_recompression <`PetscBool`>                                                       - use recompression
108498e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD>                                          - type of compression
108598e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
1086c7a4214aSPierre Jolivet 
1087c7a4214aSPierre Jolivet   Level: intermediate
1088c7a4214aSPierre Jolivet 
10891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
1090c7a4214aSPierre Jolivet @*/
MatCreateHtoolFromKernel(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt spacedim,const PetscReal coords_target[],const PetscReal coords_source[],MatHtoolKernelFn * kernel,void * kernelctx,Mat * B)10918434afd1SBarry Smith PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords_target[], const PetscReal coords_source[], MatHtoolKernelFn *kernel, void *kernelctx, Mat *B)
1092d71ae5a4SJacob Faibussowitsch {
1093c7a4214aSPierre Jolivet   Mat        A;
1094c7a4214aSPierre Jolivet   Mat_Htool *a;
1095c7a4214aSPierre Jolivet 
1096c7a4214aSPierre Jolivet   PetscFunctionBegin;
10979566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
1098c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, spacedim, 6);
10994f572ea9SToby Isaac   PetscAssertPointer(coords_target, 7);
11004f572ea9SToby Isaac   PetscAssertPointer(coords_source, 8);
1101c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 9);
11024f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 10);
11039566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, n, M, N));
11049566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATHTOOL));
11059566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1106c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
1107c7a4214aSPierre Jolivet   a->dim       = spacedim;
1108c7a4214aSPierre Jolivet   a->kernel    = kernel;
1109c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
11109566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target));
11119566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim));
1112462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */
1113c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
11149566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source));
11159566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim));
1116462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */
1117c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
1118c7a4214aSPierre Jolivet   *B = A;
11193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1120c7a4214aSPierre Jolivet }
1121c7a4214aSPierre Jolivet 
1122c7a4214aSPierre Jolivet /*MC
1123c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
1124c7a4214aSPierre Jolivet 
11252ef1f0ffSBarry Smith   Use `./configure --download-htool` to install PETSc to use Htool.
1126c7a4214aSPierre Jolivet 
11272ef1f0ffSBarry Smith    Options Database Key:
11282ef1f0ffSBarry Smith .     -mat_type htool - matrix type to `MATHTOOL`
1129c7a4214aSPierre Jolivet 
1130c7a4214aSPierre Jolivet    Level: beginner
1131c7a4214aSPierre Jolivet 
11321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()`
1133c7a4214aSPierre Jolivet M*/
MatCreate_Htool(Mat A)1134d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
1135d71ae5a4SJacob Faibussowitsch {
1136c7a4214aSPierre Jolivet   Mat_Htool *a;
1137c7a4214aSPierre Jolivet 
1138c7a4214aSPierre Jolivet   PetscFunctionBegin;
1139c9a4fd69SAudic XU   PetscCall(MatSetType(A, MATSHELL));
11404dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
1141c9a4fd69SAudic XU   PetscCall(MatShellSetContext(A, a));
114257d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_Htool));
114357d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (PetscErrorCodeFn *)MatGetDiagonalBlock_Htool));
114457d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Htool));
114557d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Htool));
114657d50842SBarry Smith   if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Htool));
1147c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
1148c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
114957d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (PetscErrorCodeFn *)MatView_Htool));
115057d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (PetscErrorCodeFn *)MatSetFromOptions_Htool));
115157d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (PetscErrorCodeFn *)MatGetRow_Htool));
115257d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (PetscErrorCodeFn *)MatRestoreRow_Htool));
115357d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (PetscErrorCodeFn *)MatAssemblyEnd_Htool));
115457d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (PetscErrorCodeFn *)MatTranspose_Htool));
115557d50842SBarry Smith   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_Htool));
1156c7a4214aSPierre Jolivet   a->dim                    = 0;
1157f22e26b7SPierre Jolivet   a->gcoords_target         = nullptr;
1158f22e26b7SPierre Jolivet   a->gcoords_source         = nullptr;
1159c095613aSPierre Jolivet   a->max_cluster_leaf_size  = 10;
1160c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
1161c7a4214aSPierre Jolivet   a->eta                    = 10.0;
1162c7a4214aSPierre Jolivet   a->depth[0]               = 0;
1163c7a4214aSPierre Jolivet   a->depth[1]               = 0;
11641dd4f53aSPierre Jolivet   a->block_tree_consistency = PETSC_TRUE;
1165c095613aSPierre Jolivet   a->permutation            = PETSC_TRUE;
1166c095613aSPierre Jolivet   a->recompression          = PETSC_FALSE;
1167c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
11689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool));
11699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool));
11709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense));
11719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense));
11729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool));
11739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool));
11749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool));
11759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool));
11769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool));
1177c095613aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUseRecompression_C", MatHtoolUseRecompression_Htool));
1178c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
1179c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
1180c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
1181c9a4fd69SAudic XU   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL));
11823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1183c7a4214aSPierre Jolivet }
1184