xref: /petsc/src/mat/impls/htool/htool.cxx (revision 49abdd8a111d9c2ef7fc48ade253ef64e07f9b37) !
1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/
2c7a4214aSPierre Jolivet #include <petscblaslapack.h>
3c7a4214aSPierre Jolivet #include <set>
4c7a4214aSPierre Jolivet 
5c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = {"sympartialACA", "fullACA", "SVD"};
698e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = {"PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric"};
7c7a4214aSPierre Jolivet const char        HtoolCitation[]           = "@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"
16c7a4214aSPierre Jolivet                                               "}\n";
17c7a4214aSPierre Jolivet static PetscBool  HtoolCite                 = PETSC_FALSE;
18c7a4214aSPierre Jolivet 
19d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Htool(Mat A, Vec v)
20d71ae5a4SJacob Faibussowitsch {
21c9a4fd69SAudic XU   Mat_Htool   *a;
22c7a4214aSPierre Jolivet   PetscScalar *x;
23c7a4214aSPierre Jolivet   PetscBool    flg;
24c7a4214aSPierre Jolivet 
25c7a4214aSPierre Jolivet   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
275f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
28c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(v, &x));
301dd4f53aSPierre Jolivet   PetscStackCallExternalVoid("copy_diagonal_in_user_numbering", htool::copy_diagonal_in_user_numbering(a->distributed_operator_holder->hmatrix, x));
319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(v, &x));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33c7a4214aSPierre Jolivet }
34c7a4214aSPierre Jolivet 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b)
36d71ae5a4SJacob Faibussowitsch {
37c9a4fd69SAudic XU   Mat_Htool                 *a;
38c7a4214aSPierre Jolivet   Mat                        B;
391dd4f53aSPierre Jolivet   PetscScalar               *ptr, shift, scale;
40c7a4214aSPierre Jolivet   PetscBool                  flg;
411dd4f53aSPierre Jolivet   PetscMPIInt                rank;
421dd4f53aSPierre Jolivet   htool::Cluster<PetscReal> *source_cluster = nullptr;
43c7a4214aSPierre Jolivet 
44c7a4214aSPierre Jolivet   PetscFunctionBegin;
459566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
465f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
47c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
49c7a4214aSPierre Jolivet   if (!B) {
5066b4df27SPierre 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));
51f22e26b7SPierre Jolivet     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B));
529566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B, &ptr));
531dd4f53aSPierre Jolivet     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
541dd4f53aSPierre Jolivet     source_cluster = a->source_cluster ? a->source_cluster.get() : a->target_cluster.get();
551dd4f53aSPierre 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));
569566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B, &ptr));
579566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A, B));
589566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B));
59c7a4214aSPierre Jolivet     *b = B;
609566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
6166b4df27SPierre Jolivet     PetscCall(MatShift(*b, shift));
6266b4df27SPierre Jolivet     PetscCall(MatScale(*b, scale));
63c9a4fd69SAudic XU   } else {
6466b4df27SPierre 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));
65c9a4fd69SAudic XU     *b = B;
66c9a4fd69SAudic XU   }
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68c7a4214aSPierre Jolivet }
69c7a4214aSPierre Jolivet 
70d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y)
71d71ae5a4SJacob Faibussowitsch {
72c9a4fd69SAudic XU   Mat_Htool         *a;
73c7a4214aSPierre Jolivet   const PetscScalar *in;
74c7a4214aSPierre Jolivet   PetscScalar       *out;
75c7a4214aSPierre Jolivet 
76c7a4214aSPierre Jolivet   PetscFunctionBegin;
77c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
799566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
801dd4f53aSPierre Jolivet   a->distributed_operator_holder->distributed_operator.vector_product_local_to_local(in, out, nullptr);
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84c7a4214aSPierre Jolivet }
85c7a4214aSPierre Jolivet 
86d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y)
87d71ae5a4SJacob Faibussowitsch {
88c9a4fd69SAudic XU   Mat_Htool         *a;
898b8ddd11SPierre Jolivet   const PetscScalar *in;
908b8ddd11SPierre Jolivet   PetscScalar       *out;
918b8ddd11SPierre Jolivet 
928b8ddd11SPierre Jolivet   PetscFunctionBegin;
93c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
961dd4f53aSPierre Jolivet   a->distributed_operator_holder->distributed_operator.vector_product_transp_local_to_local(in, out, nullptr);
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1008b8ddd11SPierre Jolivet }
1018b8ddd11SPierre Jolivet 
102d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov)
103d71ae5a4SJacob Faibussowitsch {
104c7a4214aSPierre Jolivet   std::set<PetscInt> set;
105c7a4214aSPierre Jolivet   const PetscInt    *idx;
1068f308287SPierre Jolivet   PetscInt          *oidx, size, bs[2];
107c7a4214aSPierre Jolivet   PetscMPIInt        csize;
108c7a4214aSPierre Jolivet 
109c7a4214aSPierre Jolivet   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, bs, bs + 1));
1118f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
112c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < is_max; ++i) {
113c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
114c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1159566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize));
1161dd4f53aSPierre Jolivet     PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel IS");
1179566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i], &size));
1189566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx));
119c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < size; ++j) {
120c7a4214aSPierre Jolivet       set.insert(idx[j]);
121c7a4214aSPierre Jolivet       for (PetscInt k = 1; k <= ov; ++k) {                                              /* for each layer of overlap      */
122c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k);                                    /* do not insert negative indices */
123c7a4214aSPierre 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 */
124c7a4214aSPierre Jolivet       }
125c7a4214aSPierre Jolivet     }
1269566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i], &idx));
1279566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is + i));
1288f308287SPierre Jolivet     if (bs[0] > 1) {
1298f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) {
1308f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1318f308287SPierre Jolivet         std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]);
1328f308287SPierre Jolivet         set.insert(block.cbegin(), block.cend());
1338f308287SPierre Jolivet       }
1348f308287SPierre Jolivet     }
135c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &oidx));
137c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
138c7a4214aSPierre Jolivet     oidx -= size;
1399566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i));
140c7a4214aSPierre Jolivet   }
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142c7a4214aSPierre Jolivet }
143c7a4214aSPierre Jolivet 
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
145d71ae5a4SJacob Faibussowitsch {
146c9a4fd69SAudic XU   Mat_Htool         *a;
147c7a4214aSPierre Jolivet   Mat                D, B, BT;
148c7a4214aSPierre Jolivet   const PetscScalar *copy;
14966b4df27SPierre Jolivet   PetscScalar       *ptr, shift, scale;
150c7a4214aSPierre Jolivet   const PetscInt    *idxr, *idxc, *it;
151c7a4214aSPierre Jolivet   PetscInt           nrow, m, i;
152c7a4214aSPierre Jolivet   PetscBool          flg;
153c7a4214aSPierre Jolivet 
154c7a4214aSPierre Jolivet   PetscFunctionBegin;
15566b4df27SPierre 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));
156c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
15748a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
158c7a4214aSPierre Jolivet   for (i = 0; i < n; ++i) {
1599566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &nrow));
1609566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &m));
1619566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &idxr));
1629566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &idxc));
163f22e26b7SPierre Jolivet     if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i));
1649566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr));
165c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1669566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A, &flg));
167c7a4214aSPierre Jolivet       if (flg) {
1689566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i], &flg));
169c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
170c7a4214aSPierre Jolivet           it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart);
171c7a4214aSPierre Jolivet           if (it != idxr + nrow && *it == A->rmap->rstart) {    /* rmap->rstart in IS? */
172c7a4214aSPierre Jolivet             if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
1739371c9d4SSatish Balay               for (PetscInt j = 0; j < A->rmap->n && flg; ++j)
1749371c9d4SSatish Balay                 if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE;
175c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
176c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
177c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
178c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
179c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
180c7a4214aSPierre Jolivet                  */
181c7a4214aSPierre Jolivet                 m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */
182c9a4fd69SAudic XU                 PetscCall(MatGetDiagonalBlock(A, &D));
1839566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D, &copy));
1849371c9d4SSatish Balay                 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 */ }
1859566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D, &copy));
186c7a4214aSPierre Jolivet                 if (m) {
187c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */
188c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
189b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
1909566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B));
1919566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
1929566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT));
1939566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
194b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
1959566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
196c7a4214aSPierre Jolivet                     } else {
1977fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
1989566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
199c7a4214aSPierre Jolivet                     }
2009566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2019566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
202c7a4214aSPierre Jolivet                   } else {
203c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */
204c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow);
205c7a4214aSPierre Jolivet                     }
206c7a4214aSPierre Jolivet                   }
207c7a4214aSPierre Jolivet                 }
208c7a4214aSPierre Jolivet                 if (m + A->rmap->n != nrow) {
209c7a4214aSPierre 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 */
210c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
211b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
2129566063dSJacob 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));
2139566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
2149566063dSJacob 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));
2159566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
216b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
2179566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
218c7a4214aSPierre Jolivet                     } else {
2197fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
2209566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
221c7a4214aSPierre Jolivet                     }
2229566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2239566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
224c7a4214aSPierre Jolivet                   } else {
225c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */
226c7a4214aSPierre 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);
227c7a4214aSPierre Jolivet                     }
228c7a4214aSPierre Jolivet                   }
229c7a4214aSPierre Jolivet                 }
230c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
231c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
232c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE;   /* rmap->rstart not in IS */
233c7a4214aSPierre Jolivet         } /* unsorted IS */
234c7a4214aSPierre Jolivet       }
235c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE;                                       /* different row and column IS */
236c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */
2379566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &idxr));
2389566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &idxc));
2399566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr));
24066b4df27SPierre Jolivet     PetscCall(MatShift((*submat)[i], shift));
24166b4df27SPierre Jolivet     PetscCall(MatScale((*submat)[i], scale));
242c7a4214aSPierre Jolivet   }
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244c7a4214aSPierre Jolivet }
245c7a4214aSPierre Jolivet 
246d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A)
247d71ae5a4SJacob Faibussowitsch {
248c9a4fd69SAudic XU   Mat_Htool               *a;
249c7a4214aSPierre Jolivet   PetscContainer           container;
250c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
251c7a4214aSPierre Jolivet 
252c7a4214aSPierre Jolivet   PetscFunctionBegin;
253c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
254f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr));
255f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr));
256f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr));
257f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr));
258f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr));
259f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr));
260f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr));
261f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr));
262f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr));
2639566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container));
264c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2659566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
2669566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
267f22e26b7SPierre Jolivet     PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr));
268c7a4214aSPierre Jolivet   }
26948a46eb9SPierre Jolivet   if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
2719566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->work_source, a->work_target));
272c7a4214aSPierre Jolivet   delete a->wrapper;
2731dd4f53aSPierre Jolivet   a->target_cluster.reset();
2741dd4f53aSPierre Jolivet   a->source_cluster.reset();
2751dd4f53aSPierre Jolivet   a->distributed_operator_holder.reset();
276c9a4fd69SAudic XU   PetscCall(PetscFree(a));
277c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable()
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279c7a4214aSPierre Jolivet }
280c7a4214aSPierre Jolivet 
281d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv)
282d71ae5a4SJacob Faibussowitsch {
283c9a4fd69SAudic XU   Mat_Htool                         *a;
28466b4df27SPierre Jolivet   PetscScalar                        shift, scale;
285c7a4214aSPierre Jolivet   PetscBool                          flg;
2861dd4f53aSPierre Jolivet   std::map<std::string, std::string> hmatrix_information;
287c7a4214aSPierre Jolivet 
288c7a4214aSPierre Jolivet   PetscFunctionBegin;
289c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
2901dd4f53aSPierre Jolivet   hmatrix_information = htool::get_distributed_hmatrix_information(a->distributed_operator_holder->hmatrix, PetscObjectComm((PetscObject)A));
2919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg));
292c7a4214aSPierre Jolivet   if (flg) {
29366b4df27SPierre 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));
2941dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->distributed_operator_holder->distributed_operator.get_symmetry_type()));
29566b4df27SPierre Jolivet     if (PetscAbsScalar(scale - 1.0) > PETSC_MACHINE_EPSILON) {
296c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
29766b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(scale), (double)PetscImaginaryPart(scale)));
298c7a4214aSPierre Jolivet #else
29966b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)scale));
300c9a4fd69SAudic XU #endif
301c9a4fd69SAudic XU     }
30266b4df27SPierre Jolivet     if (PetscAbsScalar(shift) > PETSC_MACHINE_EPSILON) {
303c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX)
30466b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(shift), (double)PetscImaginaryPart(shift)));
305c9a4fd69SAudic XU #else
30666b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)shift));
307c7a4214aSPierre Jolivet #endif
308c7a4214aSPierre Jolivet     }
3091dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "minimum cluster size: %" PetscInt_FMT "\n", a->min_cluster_size));
3109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon));
3119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta));
3129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0]));
3139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1]));
3149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor]));
3159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering]));
3161dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", hmatrix_information["Compression_ratio"].c_str()));
3171dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", hmatrix_information["Space_saving"].c_str()));
3181dd4f53aSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(pv, "block tree consistency: %s\n", PetscBools[a->distributed_operator_holder->hmatrix.is_block_tree_consistent()]));
3191dd4f53aSPierre 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()));
3201dd4f53aSPierre Jolivet     PetscCall(
3211dd4f53aSPierre 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()));
3221dd4f53aSPierre 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(),
3231dd4f53aSPierre Jolivet                                      hmatrix_information["Low_rank_block_size_max"].c_str()));
3241dd4f53aSPierre 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()));
325c7a4214aSPierre Jolivet   }
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327c7a4214aSPierre Jolivet }
328c7a4214aSPierre Jolivet 
329c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
330d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
331d71ae5a4SJacob Faibussowitsch {
332c9a4fd69SAudic XU   Mat_Htool   *a;
33366b4df27SPierre Jolivet   PetscScalar  shift, scale;
334c7a4214aSPierre Jolivet   PetscInt    *idxc;
335c7a4214aSPierre Jolivet   PetscBLASInt one = 1, bn;
336c7a4214aSPierre Jolivet 
337c7a4214aSPierre Jolivet   PetscFunctionBegin;
33866b4df27SPierre 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));
339c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
340c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
341c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, &idxc));
343c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i;
344c7a4214aSPierre Jolivet   }
345c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
346c7a4214aSPierre Jolivet   if (v) {
3479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, v));
348c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
349cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
3509566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N, &bn));
35166b4df27SPierre Jolivet     PetscCallBLAS("BLASscal", BLASscal_(&bn, &scale, *v, &one));
35266b4df27SPierre Jolivet     if (row < A->cmap->N) (*v)[row] += shift;
353c7a4214aSPierre Jolivet   }
35448a46eb9SPierre Jolivet   if (!idx) PetscCall(PetscFree(idxc));
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356c7a4214aSPierre Jolivet }
357c7a4214aSPierre Jolivet 
358f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v)
359d71ae5a4SJacob Faibussowitsch {
360c7a4214aSPierre Jolivet   PetscFunctionBegin;
36148a46eb9SPierre Jolivet   if (idx) PetscCall(PetscFree(*idx));
36248a46eb9SPierre Jolivet   if (v) PetscCall(PetscFree(*v));
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
364c7a4214aSPierre Jolivet }
365c7a4214aSPierre Jolivet 
366d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject)
367d71ae5a4SJacob Faibussowitsch {
368c9a4fd69SAudic XU   Mat_Htool *a;
369c7a4214aSPierre Jolivet   PetscInt   n;
370c7a4214aSPierre Jolivet   PetscBool  flg;
371c7a4214aSPierre Jolivet 
372c7a4214aSPierre Jolivet   PetscFunctionBegin;
373c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
374d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Htool options");
3751dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", nullptr, a->min_cluster_size, &a->min_cluster_size, nullptr, 0));
3761dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr, 0.0));
377f22e26b7SPierre Jolivet   PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr));
3781dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBoundedInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr, 0));
3791dd4f53aSPierre 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));
3801dd4f53aSPierre Jolivet   PetscCall(PetscOptionsBool("-mat_htool_block_tree_consistency", "Block tree consistency", nullptr, a->block_tree_consistency, &a->block_tree_consistency, nullptr));
3811dd4f53aSPierre Jolivet 
382c7a4214aSPierre Jolivet   n = 0;
383dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg));
384c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
38598e73e17SPierre Jolivet   n = 0;
386dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg));
38798e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
388d0609cedSBarry Smith   PetscOptionsHeadEnd();
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
390c7a4214aSPierre Jolivet }
391c7a4214aSPierre Jolivet 
392cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType)
393d71ae5a4SJacob Faibussowitsch {
394c9a4fd69SAudic XU   Mat_Htool                                                   *a;
395c7a4214aSPierre Jolivet   const PetscInt                                              *ranges;
396c7a4214aSPierre Jolivet   PetscInt                                                    *offset;
3971dd4f53aSPierre Jolivet   PetscMPIInt                                                  size, rank;
398b94d7dedSBarry Smith   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U';
399cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                        *generator = nullptr;
4001dd4f53aSPierre Jolivet   htool::ClusterTreeBuilder<PetscReal>                         recursive_build_strategy;
4011dd4f53aSPierre Jolivet   htool::Cluster<PetscReal>                                   *source_cluster;
4021dd4f53aSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor;
403c7a4214aSPierre Jolivet 
404c7a4214aSPierre Jolivet   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite));
406c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
407c7a4214aSPierre Jolivet   delete a->wrapper;
4081dd4f53aSPierre Jolivet   a->target_cluster.reset();
4091dd4f53aSPierre Jolivet   a->source_cluster.reset();
4101dd4f53aSPierre Jolivet   a->distributed_operator_holder.reset();
4111dd4f53aSPierre Jolivet   // clustering
4129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &offset));
4149566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
415c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < size; ++i) {
416c7a4214aSPierre Jolivet     offset[2 * i]     = ranges[i];
417c7a4214aSPierre Jolivet     offset[2 * i + 1] = ranges[i + 1] - ranges[i];
418c7a4214aSPierre Jolivet   }
41998e73e17SPierre Jolivet   switch (a->clustering) {
420d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
4211dd4f53aSPierre Jolivet     recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>());
4221dd4f53aSPierre Jolivet     recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>());
423d71ae5a4SJacob Faibussowitsch     break;
424d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
4251dd4f53aSPierre Jolivet     recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>());
4261dd4f53aSPierre Jolivet     recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>());
427d71ae5a4SJacob Faibussowitsch     break;
428d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
4291dd4f53aSPierre Jolivet     recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>());
4301dd4f53aSPierre Jolivet     recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>());
431d71ae5a4SJacob Faibussowitsch     break;
432d71ae5a4SJacob Faibussowitsch   default:
4331dd4f53aSPierre Jolivet     recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>());
4341dd4f53aSPierre Jolivet     recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>());
43598e73e17SPierre Jolivet   }
4361dd4f53aSPierre Jolivet   recursive_build_strategy.set_minclustersize(a->min_cluster_size);
4371dd4f53aSPierre Jolivet   a->target_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree(A->rmap->N, a->dim, a->gcoords_target, 2, size, offset));
438c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4399566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A, &ranges));
440c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < size; ++i) {
441c7a4214aSPierre Jolivet       offset[2 * i]     = ranges[i];
442c7a4214aSPierre Jolivet       offset[2 * i + 1] = ranges[i + 1] - ranges[i];
443c7a4214aSPierre Jolivet     }
44498e73e17SPierre Jolivet     switch (a->clustering) {
445d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
4461dd4f53aSPierre Jolivet       recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>());
4471dd4f53aSPierre Jolivet       recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>());
448d71ae5a4SJacob Faibussowitsch       break;
449d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
4501dd4f53aSPierre Jolivet       recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>());
4511dd4f53aSPierre Jolivet       recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>());
452d71ae5a4SJacob Faibussowitsch       break;
453d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
4541dd4f53aSPierre Jolivet       recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>());
4551dd4f53aSPierre Jolivet       recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>());
456d71ae5a4SJacob Faibussowitsch       break;
457d71ae5a4SJacob Faibussowitsch     default:
4581dd4f53aSPierre Jolivet       recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>());
4591dd4f53aSPierre Jolivet       recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>());
46098e73e17SPierre Jolivet     }
4611dd4f53aSPierre Jolivet     recursive_build_strategy.set_minclustersize(a->min_cluster_size);
4621dd4f53aSPierre Jolivet     a->source_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree(A->cmap->N, a->dim, a->gcoords_source, 2, size, offset));
463c7a4214aSPierre Jolivet     S = uplo       = 'N';
4641dd4f53aSPierre Jolivet     source_cluster = a->source_cluster.get();
4651dd4f53aSPierre Jolivet   } else source_cluster = a->target_cluster.get();
4669566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
4671dd4f53aSPierre Jolivet   // generator
4681dd4f53aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx);
4691dd4f53aSPierre Jolivet   else {
4701dd4f53aSPierre Jolivet     a->wrapper = nullptr;
4711dd4f53aSPierre Jolivet     generator  = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx);
4721dd4f53aSPierre Jolivet   }
4731dd4f53aSPierre Jolivet   // compressor
474c7a4214aSPierre Jolivet   switch (a->compressor) {
475d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
476d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::fullACA<PetscScalar>>();
477d71ae5a4SJacob Faibussowitsch     break;
478d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_SVD:
479d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::SVD<PetscScalar>>();
480d71ae5a4SJacob Faibussowitsch     break;
481d71ae5a4SJacob Faibussowitsch   default:
482d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
483c7a4214aSPierre Jolivet   }
4841dd4f53aSPierre Jolivet   // local hierarchical matrix
4851dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
4861dd4f53aSPierre Jolivet   auto hmatrix_builder = htool::HMatrixTreeBuilder<PetscScalar>(*a->target_cluster, *source_cluster, a->epsilon, a->eta, S, uplo, -1, rank, rank);
4871dd4f53aSPierre Jolivet   hmatrix_builder.set_low_rank_generator(compressor);
4881dd4f53aSPierre Jolivet   hmatrix_builder.set_minimal_target_depth(a->depth[0]);
4891dd4f53aSPierre Jolivet   hmatrix_builder.set_minimal_source_depth(a->depth[1]);
4901dd4f53aSPierre 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");
4911dd4f53aSPierre Jolivet   hmatrix_builder.set_block_tree_consistency(a->block_tree_consistency);
4921dd4f53aSPierre Jolivet   a->distributed_operator_holder = std::make_unique<htool::DistributedOperatorFromHMatrix<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, *a->target_cluster, *source_cluster, hmatrix_builder, PetscObjectComm((PetscObject)A));
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
494c7a4214aSPierre Jolivet }
495c7a4214aSPierre Jolivet 
496d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C)
497d71ae5a4SJacob Faibussowitsch {
498c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
499c9a4fd69SAudic XU   Mat_Htool         *a;
500c7a4214aSPierre Jolivet   const PetscScalar *in;
501c7a4214aSPierre Jolivet   PetscScalar       *out;
5028b8ddd11SPierre Jolivet   PetscInt           N, lda;
503c7a4214aSPierre Jolivet 
504c7a4214aSPierre Jolivet   PetscFunctionBegin;
505c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
506f22e26b7SPierre Jolivet   PetscCall(MatGetSize(C, nullptr, &N));
5079566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &lda));
5085f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
5099566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B, &in));
5109566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &out));
511c9a4fd69SAudic XU   PetscCall(MatShellGetContext(product->A, &a));
5128b8ddd11SPierre Jolivet   switch (product->type) {
513d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AB:
5141dd4f53aSPierre Jolivet     a->distributed_operator_holder->distributed_operator.matrix_product_local_to_local(in, out, N, nullptr);
515d71ae5a4SJacob Faibussowitsch     break;
516d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AtB:
5171dd4f53aSPierre Jolivet     a->distributed_operator_holder->distributed_operator.matrix_product_transp_local_to_local(in, out, N, nullptr);
518d71ae5a4SJacob Faibussowitsch     break;
519d71ae5a4SJacob Faibussowitsch   default:
520d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]);
521c7a4214aSPierre Jolivet   }
5229566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &out));
5239566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B, &in));
5243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525c7a4214aSPierre Jolivet }
526c7a4214aSPierre Jolivet 
527d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C)
528d71ae5a4SJacob Faibussowitsch {
529c7a4214aSPierre Jolivet   Mat_Product *product = C->product;
530c7a4214aSPierre Jolivet   Mat          A, B;
531c7a4214aSPierre Jolivet   PetscBool    flg;
532c7a4214aSPierre Jolivet 
533c7a4214aSPierre Jolivet   PetscFunctionBegin;
534c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
535c7a4214aSPierre Jolivet   A = product->A;
536c7a4214aSPierre Jolivet   B = product->B;
5379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, ""));
53897713f8cSPierre 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);
53997713f8cSPierre Jolivet   if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
54097713f8cSPierre Jolivet     if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
54197713f8cSPierre Jolivet     else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
5428b8ddd11SPierre Jolivet   }
5439566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
5449566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
5459566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
5469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
5479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
548f22e26b7SPierre Jolivet   C->ops->productsymbolic = nullptr;
549c7a4214aSPierre Jolivet   C->ops->productnumeric  = MatProductNumeric_Htool;
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551c7a4214aSPierre Jolivet }
552c7a4214aSPierre Jolivet 
553d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
554d71ae5a4SJacob Faibussowitsch {
555c7a4214aSPierre Jolivet   PetscFunctionBegin;
556c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
5578b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
559c7a4214aSPierre Jolivet }
560c7a4214aSPierre Jolivet 
5611dd4f53aSPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator)
562d71ae5a4SJacob Faibussowitsch {
563c9a4fd69SAudic XU   Mat_Htool *a;
564c7a4214aSPierre Jolivet 
565c7a4214aSPierre Jolivet   PetscFunctionBegin;
566c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
5671dd4f53aSPierre Jolivet   *distributed_operator = &a->distributed_operator_holder->distributed_operator;
5683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
569c7a4214aSPierre Jolivet }
570c7a4214aSPierre Jolivet 
571c7a4214aSPierre Jolivet /*@C
57211a5261eSBarry Smith   MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`.
573c7a4214aSPierre Jolivet 
574cc4c1da9SBarry Smith   No Fortran Support, No C Support
575cc4c1da9SBarry Smith 
576c7a4214aSPierre Jolivet   Input Parameter:
577c7a4214aSPierre Jolivet . A - hierarchical matrix
578c7a4214aSPierre Jolivet 
579c7a4214aSPierre Jolivet   Output Parameter:
5801dd4f53aSPierre Jolivet . distributed_operator - opaque pointer to a Htool virtual matrix
581c7a4214aSPierre Jolivet 
582c7a4214aSPierre Jolivet   Level: advanced
583c7a4214aSPierre Jolivet 
5841cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`
585c7a4214aSPierre Jolivet @*/
5861dd4f53aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator)
587d71ae5a4SJacob Faibussowitsch {
588c7a4214aSPierre Jolivet   PetscFunctionBegin;
589c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
5901dd4f53aSPierre Jolivet   PetscAssertPointer(distributed_operator, 2);
5911dd4f53aSPierre Jolivet   PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::DistributedOperator<PetscScalar> **), (A, distributed_operator));
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593c7a4214aSPierre Jolivet }
594c7a4214aSPierre Jolivet 
5958434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
596d71ae5a4SJacob Faibussowitsch {
597c9a4fd69SAudic XU   Mat_Htool *a;
598c7a4214aSPierre Jolivet 
599c7a4214aSPierre Jolivet   PetscFunctionBegin;
600c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
601c7a4214aSPierre Jolivet   a->kernel    = kernel;
602c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
603c7a4214aSPierre Jolivet   delete a->wrapper;
6041dd4f53aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx);
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606c7a4214aSPierre Jolivet }
607c7a4214aSPierre Jolivet 
608c7a4214aSPierre Jolivet /*@C
60911a5261eSBarry Smith   MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`.
610c7a4214aSPierre Jolivet 
611cc4c1da9SBarry Smith   Collective, No Fortran Support
612cc4c1da9SBarry Smith 
613c7a4214aSPierre Jolivet   Input Parameters:
614c7a4214aSPierre Jolivet + A         - hierarchical matrix
6152ef1f0ffSBarry Smith . kernel    - computational kernel (or `NULL`)
6162ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
617c7a4214aSPierre Jolivet 
618c7a4214aSPierre Jolivet   Level: advanced
619c7a4214aSPierre Jolivet 
6201cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()`
621c7a4214aSPierre Jolivet @*/
622cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
623d71ae5a4SJacob Faibussowitsch {
624c7a4214aSPierre Jolivet   PetscFunctionBegin;
625c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
626c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 2);
6274f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 3);
6288434afd1SBarry Smith   PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx));
6293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
630c7a4214aSPierre Jolivet }
631c7a4214aSPierre Jolivet 
632d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is)
633d71ae5a4SJacob Faibussowitsch {
634c9a4fd69SAudic XU   Mat_Htool                       *a;
6351dd4f53aSPierre Jolivet   PetscMPIInt                      rank;
6361dd4f53aSPierre Jolivet   const std::vector<PetscInt>     *source;
6371dd4f53aSPierre Jolivet   const htool::Cluster<PetscReal> *local_source_cluster;
63898e73e17SPierre Jolivet 
63998e73e17SPierre Jolivet   PetscFunctionBegin;
640c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
6411dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
6421dd4f53aSPierre Jolivet   local_source_cluster = a->source_cluster ? &a->source_cluster->get_cluster_on_partition(rank) : &a->target_cluster->get_cluster_on_partition(rank);
6431dd4f53aSPierre Jolivet   source               = &local_source_cluster->get_permutation();
6441dd4f53aSPierre Jolivet   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), local_source_cluster->get_size(), source->data() + local_source_cluster->get_offset(), PETSC_COPY_VALUES, is));
6459566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64798e73e17SPierre Jolivet }
64898e73e17SPierre Jolivet 
649cc4c1da9SBarry Smith /*@
65011a5261eSBarry Smith   MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix.
65198e73e17SPierre Jolivet 
65298e73e17SPierre Jolivet   Input Parameter:
65398e73e17SPierre Jolivet . A - hierarchical matrix
65498e73e17SPierre Jolivet 
65598e73e17SPierre Jolivet   Output Parameter:
65698e73e17SPierre Jolivet . is - permutation
65798e73e17SPierre Jolivet 
65898e73e17SPierre Jolivet   Level: advanced
65998e73e17SPierre Jolivet 
6601cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()`
66198e73e17SPierre Jolivet @*/
662cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is)
663d71ae5a4SJacob Faibussowitsch {
66498e73e17SPierre Jolivet   PetscFunctionBegin;
66598e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6664f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
667cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is));
6683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66998e73e17SPierre Jolivet }
67098e73e17SPierre Jolivet 
671d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is)
672d71ae5a4SJacob Faibussowitsch {
673c9a4fd69SAudic XU   Mat_Htool                   *a;
6741dd4f53aSPierre Jolivet   const std::vector<PetscInt> *target;
6751dd4f53aSPierre Jolivet   PetscMPIInt                  rank;
67698e73e17SPierre Jolivet 
67798e73e17SPierre Jolivet   PetscFunctionBegin;
678c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
6791dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
6801dd4f53aSPierre Jolivet   target = &a->target_cluster->get_permutation();
6811dd4f53aSPierre 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));
6829566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68498e73e17SPierre Jolivet }
68598e73e17SPierre Jolivet 
686cc4c1da9SBarry Smith /*@
68711a5261eSBarry Smith   MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix.
68898e73e17SPierre Jolivet 
68998e73e17SPierre Jolivet   Input Parameter:
69098e73e17SPierre Jolivet . A - hierarchical matrix
69198e73e17SPierre Jolivet 
69298e73e17SPierre Jolivet   Output Parameter:
69398e73e17SPierre Jolivet . is - permutation
69498e73e17SPierre Jolivet 
69598e73e17SPierre Jolivet   Level: advanced
69698e73e17SPierre Jolivet 
6971cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()`
69898e73e17SPierre Jolivet @*/
699cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is)
700d71ae5a4SJacob Faibussowitsch {
70198e73e17SPierre Jolivet   PetscFunctionBegin;
70298e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
7034f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
704cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is));
7053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70698e73e17SPierre Jolivet }
70798e73e17SPierre Jolivet 
708d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use)
709d71ae5a4SJacob Faibussowitsch {
710c9a4fd69SAudic XU   Mat_Htool *a;
71198e73e17SPierre Jolivet 
71298e73e17SPierre Jolivet   PetscFunctionBegin;
713c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
7141dd4f53aSPierre Jolivet   a->distributed_operator_holder->distributed_operator.use_permutation() = use;
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71698e73e17SPierre Jolivet }
71798e73e17SPierre Jolivet 
718cc4c1da9SBarry Smith /*@
71911a5261eSBarry Smith   MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation.
72098e73e17SPierre Jolivet 
72198e73e17SPierre Jolivet   Input Parameters:
72298e73e17SPierre Jolivet + A   - hierarchical matrix
72398e73e17SPierre Jolivet - use - Boolean value
72498e73e17SPierre Jolivet 
72598e73e17SPierre Jolivet   Level: advanced
72698e73e17SPierre Jolivet 
7271cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()`
72898e73e17SPierre Jolivet @*/
729cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use)
730d71ae5a4SJacob Faibussowitsch {
73198e73e17SPierre Jolivet   PetscFunctionBegin;
73298e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
73398e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A, use, 2);
734cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use));
7353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73698e73e17SPierre Jolivet }
73798e73e17SPierre Jolivet 
738cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B)
739d71ae5a4SJacob Faibussowitsch {
740c7a4214aSPierre Jolivet   Mat          C;
741c9a4fd69SAudic XU   Mat_Htool   *a;
74266b4df27SPierre Jolivet   PetscScalar *array, shift, scale;
743c7a4214aSPierre Jolivet   PetscInt     lda;
744c7a4214aSPierre Jolivet 
745c7a4214aSPierre Jolivet   PetscFunctionBegin;
74666b4df27SPierre 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));
747c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
748c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
749c7a4214aSPierre Jolivet     C = *B;
7505f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions");
7519566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C, &lda));
7525f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
753c7a4214aSPierre Jolivet   } else {
7549566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
7559566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
7569566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, MATDENSE));
7579566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7589566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
759c7a4214aSPierre Jolivet   }
7609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
7619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
7629566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
7631dd4f53aSPierre Jolivet   htool::copy_to_dense_in_user_numbering(a->distributed_operator_holder->hmatrix, array);
7649566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
76566b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
76666b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
767c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
7689566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &C));
769c7a4214aSPierre Jolivet   } else *B = C;
7703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
771c7a4214aSPierre Jolivet }
772c7a4214aSPierre Jolivet 
773d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx)
774d71ae5a4SJacob Faibussowitsch {
775c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx;
77698e73e17SPierre Jolivet   PetscScalar             *tmp;
77798e73e17SPierre Jolivet 
77898e73e17SPierre Jolivet   PetscFunctionBegin;
7793ba16761SJacob Faibussowitsch   PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx));
7809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M * N, &tmp));
7819566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp, ptr, M * N));
78298e73e17SPierre Jolivet   for (PetscInt i = 0; i < M; ++i) {
78398e73e17SPierre Jolivet     for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N];
78498e73e17SPierre Jolivet   }
7859566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
7863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
787c7a4214aSPierre Jolivet }
788c7a4214aSPierre Jolivet 
789c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
790d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B)
791d71ae5a4SJacob Faibussowitsch {
792c7a4214aSPierre Jolivet   Mat                      C;
793c9a4fd69SAudic XU   Mat_Htool               *a, *c;
79466b4df27SPierre Jolivet   PetscScalar              shift, scale;
795c7a4214aSPierre Jolivet   PetscInt                 M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n;
796c7a4214aSPierre Jolivet   PetscContainer           container;
797c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
798c7a4214aSPierre Jolivet 
799c7a4214aSPierre Jolivet   PetscFunctionBegin;
80066b4df27SPierre 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));
801c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
8027fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
8035f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported");
804c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8059566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
8069566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, n, m, N, M));
8079566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
8089566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
8099566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
810*49abdd8aSBarry Smith     PetscCall(PetscObjectContainerCompose((PetscObject)C, "KernelTranspose", kernelt, PetscCtxDestroyDefault));
811c7a4214aSPierre Jolivet   } else {
812c7a4214aSPierre Jolivet     C = *B;
8139566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container));
8145f80ce2aSJacob Faibussowitsch     PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first");
8159566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
816c7a4214aSPierre Jolivet   }
817c9a4fd69SAudic XU   PetscCall(MatShellGetContext(C, &c));
818c7a4214aSPierre Jolivet   c->dim = a->dim;
81966b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
82066b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
82198e73e17SPierre Jolivet   c->kernel = GenEntriesTranspose;
822c7a4214aSPierre Jolivet   if (kernelt->A != A) {
8239566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
824c7a4214aSPierre Jolivet     kernelt->A = A;
8259566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
826c7a4214aSPierre Jolivet   }
827c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
828c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
829c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
830c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8319566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target));
8329566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim));
833c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
8349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source));
8359566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim));
836c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
8379566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target));
838c7a4214aSPierre Jolivet   }
8399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
8409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
841c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
8423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
843c7a4214aSPierre Jolivet }
844c7a4214aSPierre Jolivet 
8451dd4f53aSPierre Jolivet static PetscErrorCode MatDestroy_Factor(Mat F)
8461dd4f53aSPierre Jolivet {
8471dd4f53aSPierre Jolivet   PetscContainer               container;
8481dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *A;
8491dd4f53aSPierre Jolivet 
8501dd4f53aSPierre Jolivet   PetscFunctionBegin;
8511dd4f53aSPierre Jolivet   PetscCall(PetscObjectQuery((PetscObject)F, "HMatrix", (PetscObject *)&container));
8521dd4f53aSPierre Jolivet   if (container) {
8531dd4f53aSPierre Jolivet     PetscCall(PetscContainerGetPointer(container, (void **)&A));
8541dd4f53aSPierre Jolivet     delete A;
8551dd4f53aSPierre Jolivet     PetscCall(PetscObjectCompose((PetscObject)F, "HMatrix", nullptr));
8561dd4f53aSPierre Jolivet   }
8571dd4f53aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", nullptr));
8581dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
8591dd4f53aSPierre Jolivet }
8601dd4f53aSPierre Jolivet 
8611dd4f53aSPierre Jolivet static PetscErrorCode MatFactorGetSolverType_Htool(Mat, MatSolverType *type)
8621dd4f53aSPierre Jolivet {
8631dd4f53aSPierre Jolivet   PetscFunctionBegin;
8641dd4f53aSPierre Jolivet   *type = MATSOLVERHTOOL;
8651dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
8661dd4f53aSPierre Jolivet }
8671dd4f53aSPierre Jolivet 
8681dd4f53aSPierre Jolivet template <char trans>
8691dd4f53aSPierre Jolivet static inline PetscErrorCode MatSolve_Private(Mat A, htool::Matrix<PetscScalar> &X)
8701dd4f53aSPierre Jolivet {
8711dd4f53aSPierre Jolivet   PetscContainer               container;
8721dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *B;
8731dd4f53aSPierre Jolivet 
8741dd4f53aSPierre Jolivet   PetscFunctionBegin;
8751dd4f53aSPierre 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");
8761dd4f53aSPierre Jolivet   PetscCall(PetscObjectQuery((PetscObject)A, "HMatrix", (PetscObject *)&container));
8771dd4f53aSPierre 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");
8781dd4f53aSPierre Jolivet   PetscCall(PetscContainerGetPointer(container, (void **)&B));
8791dd4f53aSPierre Jolivet   if (A->factortype == MAT_FACTOR_LU) htool::lu_solve(trans, *B, X);
8801dd4f53aSPierre Jolivet   else htool::cholesky_solve('L', *B, X);
8811dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
8821dd4f53aSPierre Jolivet }
8831dd4f53aSPierre Jolivet 
8841dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type * = nullptr>
8851dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type b, Type x)
8861dd4f53aSPierre Jolivet {
8871dd4f53aSPierre Jolivet   PetscInt                   n;
8881dd4f53aSPierre Jolivet   htool::Matrix<PetscScalar> v;
8891dd4f53aSPierre Jolivet   PetscScalar               *array;
8901dd4f53aSPierre Jolivet 
8911dd4f53aSPierre Jolivet   PetscFunctionBegin;
8921dd4f53aSPierre Jolivet   PetscCall(VecGetLocalSize(b, &n));
8931dd4f53aSPierre Jolivet   PetscCall(VecCopy(b, x));
8941dd4f53aSPierre Jolivet   PetscCall(VecGetArrayWrite(x, &array));
8951dd4f53aSPierre Jolivet   v.assign(n, 1, array, false);
8961dd4f53aSPierre Jolivet   PetscCall(VecRestoreArrayWrite(x, &array));
8971dd4f53aSPierre Jolivet   PetscCall(MatSolve_Private<trans>(A, v));
8981dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
8991dd4f53aSPierre Jolivet }
9001dd4f53aSPierre Jolivet 
9011dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type * = nullptr>
9021dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type B, Type X)
9031dd4f53aSPierre Jolivet {
9041dd4f53aSPierre Jolivet   PetscInt                   m, N;
9051dd4f53aSPierre Jolivet   htool::Matrix<PetscScalar> v;
9061dd4f53aSPierre Jolivet   PetscScalar               *array;
9071dd4f53aSPierre Jolivet 
9081dd4f53aSPierre Jolivet   PetscFunctionBegin;
9091dd4f53aSPierre Jolivet   PetscCall(MatGetLocalSize(B, &m, nullptr));
9101dd4f53aSPierre Jolivet   PetscCall(MatGetLocalSize(B, nullptr, &N));
9111dd4f53aSPierre Jolivet   PetscCall(MatCopy(B, X, SAME_NONZERO_PATTERN));
9121dd4f53aSPierre Jolivet   PetscCall(MatDenseGetArrayWrite(X, &array));
9131dd4f53aSPierre Jolivet   v.assign(m, N, array, false);
9141dd4f53aSPierre Jolivet   PetscCall(MatDenseRestoreArrayWrite(X, &array));
9151dd4f53aSPierre Jolivet   PetscCall(MatSolve_Private<trans>(A, v));
9161dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9171dd4f53aSPierre Jolivet }
9181dd4f53aSPierre Jolivet 
9191dd4f53aSPierre Jolivet template <MatFactorType ftype>
9201dd4f53aSPierre Jolivet static PetscErrorCode MatFactorNumeric_Htool(Mat F, Mat A, const MatFactorInfo *)
9211dd4f53aSPierre Jolivet {
9221dd4f53aSPierre Jolivet   Mat_Htool                   *a;
9231dd4f53aSPierre Jolivet   htool::HMatrix<PetscScalar> *B;
9241dd4f53aSPierre Jolivet 
9251dd4f53aSPierre Jolivet   PetscFunctionBegin;
9261dd4f53aSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
9271dd4f53aSPierre Jolivet   B = new htool::HMatrix<PetscScalar>(a->distributed_operator_holder->hmatrix);
9281dd4f53aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) htool::lu_factorization(*B);
9291dd4f53aSPierre Jolivet   else htool::cholesky_factorization('L', *B);
9301dd4f53aSPierre Jolivet   PetscCall(PetscObjectContainerCompose((PetscObject)F, "HMatrix", B, nullptr));
9311dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9321dd4f53aSPierre Jolivet }
9331dd4f53aSPierre Jolivet 
9341dd4f53aSPierre Jolivet template <MatFactorType ftype>
9351dd4f53aSPierre Jolivet PetscErrorCode MatFactorSymbolic_Htool(Mat F, Mat)
9361dd4f53aSPierre Jolivet {
9371dd4f53aSPierre Jolivet   PetscFunctionBegin;
9381dd4f53aSPierre Jolivet   F->preallocated  = PETSC_TRUE;
9391dd4f53aSPierre Jolivet   F->assembled     = PETSC_TRUE;
9401dd4f53aSPierre Jolivet   F->ops->solve    = MatSolve_Htool<'N', Vec>;
9411dd4f53aSPierre Jolivet   F->ops->matsolve = MatSolve_Htool<'N', Mat>;
9421dd4f53aSPierre Jolivet   if (!PetscDefined(USE_COMPLEX) || ftype == MAT_FACTOR_LU) {
9431dd4f53aSPierre Jolivet     F->ops->solvetranspose    = MatSolve_Htool<'T', Vec>;
9441dd4f53aSPierre Jolivet     F->ops->matsolvetranspose = MatSolve_Htool<'T', Mat>;
9451dd4f53aSPierre Jolivet   }
9461dd4f53aSPierre Jolivet   F->ops->destroy = MatDestroy_Factor;
9471dd4f53aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_LU>;
9481dd4f53aSPierre Jolivet   else F->ops->choleskyfactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_CHOLESKY>;
9491dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9501dd4f53aSPierre Jolivet }
9511dd4f53aSPierre Jolivet 
9521dd4f53aSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Htool(Mat F, Mat A, IS, IS, const MatFactorInfo *)
9531dd4f53aSPierre Jolivet {
9541dd4f53aSPierre Jolivet   PetscFunctionBegin;
9551dd4f53aSPierre Jolivet   PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_LU>(F, A));
9561dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9571dd4f53aSPierre Jolivet }
9581dd4f53aSPierre Jolivet 
9591dd4f53aSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Htool(Mat F, Mat A, IS, const MatFactorInfo *)
9601dd4f53aSPierre Jolivet {
9611dd4f53aSPierre Jolivet   PetscFunctionBegin;
9621dd4f53aSPierre Jolivet   PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_CHOLESKY>(F, A));
9631dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9641dd4f53aSPierre Jolivet }
9651dd4f53aSPierre Jolivet 
9661dd4f53aSPierre Jolivet static PetscErrorCode MatGetFactor_htool_htool(Mat A, MatFactorType ftype, Mat *F)
9671dd4f53aSPierre Jolivet {
9681dd4f53aSPierre Jolivet   Mat         B;
9691dd4f53aSPierre Jolivet   Mat_Htool  *a;
9701dd4f53aSPierre Jolivet   PetscMPIInt size;
9711dd4f53aSPierre Jolivet 
9721dd4f53aSPierre Jolivet   PetscFunctionBegin;
9731dd4f53aSPierre Jolivet   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
9741dd4f53aSPierre Jolivet   PetscCall(MatShellGetContext(A, &a));
9751dd4f53aSPierre Jolivet   PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel MatGetFactor()");
9761dd4f53aSPierre Jolivet   PetscCheck(a->block_tree_consistency, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot factor a MatHtool with inconsistent block tree");
9771dd4f53aSPierre Jolivet   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
9781dd4f53aSPierre Jolivet   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
9791dd4f53aSPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &((PetscObject)B)->type_name));
9801dd4f53aSPierre Jolivet   PetscCall(MatSetUp(B));
9811dd4f53aSPierre Jolivet 
9821dd4f53aSPierre Jolivet   B->ops->getinfo    = MatGetInfo_External;
9831dd4f53aSPierre Jolivet   B->factortype      = ftype;
9841dd4f53aSPierre Jolivet   B->trivialsymbolic = PETSC_TRUE;
9851dd4f53aSPierre Jolivet 
9861dd4f53aSPierre Jolivet   if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_Htool;
9871dd4f53aSPierre Jolivet   else if (ftype == MAT_FACTOR_CHOLESKY) B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_Htool;
9881dd4f53aSPierre Jolivet 
9891dd4f53aSPierre Jolivet   PetscCall(PetscFree(B->solvertype));
9901dd4f53aSPierre Jolivet   PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &B->solvertype));
9911dd4f53aSPierre Jolivet 
9921dd4f53aSPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_Htool));
9931dd4f53aSPierre Jolivet   *F = B;
9941dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
9951dd4f53aSPierre Jolivet }
9961dd4f53aSPierre Jolivet 
9971dd4f53aSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Htool(void)
9981dd4f53aSPierre Jolivet {
9991dd4f53aSPierre Jolivet   PetscFunctionBegin;
10001dd4f53aSPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_LU, MatGetFactor_htool_htool));
10011dd4f53aSPierre Jolivet   PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_CHOLESKY, MatGetFactor_htool_htool));
10021dd4f53aSPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
10031dd4f53aSPierre Jolivet }
10041dd4f53aSPierre Jolivet 
1005c7a4214aSPierre Jolivet /*@C
100611a5261eSBarry Smith   MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel.
1007c7a4214aSPierre Jolivet 
1008cc4c1da9SBarry Smith   Collective, No Fortran Support
1009cc4c1da9SBarry Smith 
1010c7a4214aSPierre Jolivet   Input Parameters:
1011c7a4214aSPierre Jolivet + comm          - MPI communicator
10122ef1f0ffSBarry Smith . m             - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
10132ef1f0ffSBarry Smith . n             - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
10142ef1f0ffSBarry Smith . M             - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
10152ef1f0ffSBarry Smith . N             - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1016c7a4214aSPierre Jolivet . spacedim      - dimension of the space coordinates
1017c7a4214aSPierre Jolivet . coords_target - coordinates of the target
1018c7a4214aSPierre Jolivet . coords_source - coordinates of the source
10192ef1f0ffSBarry Smith . kernel        - computational kernel (or `NULL`)
10202ef1f0ffSBarry Smith - kernelctx     - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
1021c7a4214aSPierre Jolivet 
1022c7a4214aSPierre Jolivet   Output Parameter:
1023c7a4214aSPierre Jolivet . B - matrix
1024c7a4214aSPierre Jolivet 
1025c7a4214aSPierre Jolivet   Options Database Keys:
102611a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`>                                                     - minimal leaf size in cluster tree
102711a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`>                                                             - relative error in Frobenius norm when approximating a block
102811a5261eSBarry Smith . -mat_htool_eta <`PetscReal`>                                                                 - admissibility condition tolerance
102911a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the rows
103011a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the columns
10311dd4f53aSPierre Jolivet . -mat_htool_block_tree_consistency <`PetscBool`>                                              - block tree consistency
103298e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD>                                          - type of compression
103398e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
1034c7a4214aSPierre Jolivet 
1035c7a4214aSPierre Jolivet   Level: intermediate
1036c7a4214aSPierre Jolivet 
10371cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
1038c7a4214aSPierre Jolivet @*/
10398434afd1SBarry 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)
1040d71ae5a4SJacob Faibussowitsch {
1041c7a4214aSPierre Jolivet   Mat        A;
1042c7a4214aSPierre Jolivet   Mat_Htool *a;
1043c7a4214aSPierre Jolivet 
1044c7a4214aSPierre Jolivet   PetscFunctionBegin;
10459566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
1046c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, spacedim, 6);
10474f572ea9SToby Isaac   PetscAssertPointer(coords_target, 7);
10484f572ea9SToby Isaac   PetscAssertPointer(coords_source, 8);
1049c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 9);
10504f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 10);
10519566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, n, M, N));
10529566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATHTOOL));
10539566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1054c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
1055c7a4214aSPierre Jolivet   a->dim       = spacedim;
1056c7a4214aSPierre Jolivet   a->kernel    = kernel;
1057c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
10589566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target));
10599566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim));
1060462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */
1061c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
10629566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source));
10639566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim));
1064462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */
1065c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
10669566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target));
1067c7a4214aSPierre Jolivet   *B = A;
10683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1069c7a4214aSPierre Jolivet }
1070c7a4214aSPierre Jolivet 
1071c7a4214aSPierre Jolivet /*MC
1072c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
1073c7a4214aSPierre Jolivet 
10742ef1f0ffSBarry Smith   Use `./configure --download-htool` to install PETSc to use Htool.
1075c7a4214aSPierre Jolivet 
10762ef1f0ffSBarry Smith    Options Database Key:
10772ef1f0ffSBarry Smith .     -mat_type htool - matrix type to `MATHTOOL`
1078c7a4214aSPierre Jolivet 
1079c7a4214aSPierre Jolivet    Level: beginner
1080c7a4214aSPierre Jolivet 
10811cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()`
1082c7a4214aSPierre Jolivet M*/
1083d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
1084d71ae5a4SJacob Faibussowitsch {
1085c7a4214aSPierre Jolivet   Mat_Htool *a;
1086c7a4214aSPierre Jolivet 
1087c7a4214aSPierre Jolivet   PetscFunctionBegin;
1088c9a4fd69SAudic XU   PetscCall(MatSetType(A, MATSHELL));
10894dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
1090c9a4fd69SAudic XU   PetscCall(MatShellSetContext(A, a));
1091c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Htool));
1092c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Htool));
1093c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Htool));
1094c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool));
1095c9a4fd69SAudic XU   if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool));
1096c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
1097c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
1098c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_Htool));
1099c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_Htool));
1100c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (void (*)(void))MatGetRow_Htool));
1101c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (void (*)(void))MatRestoreRow_Htool));
1102c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_Htool));
1103c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_Htool));
1104c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_Htool));
1105c7a4214aSPierre Jolivet   a->dim                    = 0;
1106f22e26b7SPierre Jolivet   a->gcoords_target         = nullptr;
1107f22e26b7SPierre Jolivet   a->gcoords_source         = nullptr;
11081dd4f53aSPierre Jolivet   a->min_cluster_size       = 10;
1109c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
1110c7a4214aSPierre Jolivet   a->eta                    = 10.0;
1111c7a4214aSPierre Jolivet   a->depth[0]               = 0;
1112c7a4214aSPierre Jolivet   a->depth[1]               = 0;
11131dd4f53aSPierre Jolivet   a->block_tree_consistency = PETSC_TRUE;
1114c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
11159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool));
11169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool));
11179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense));
11189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense));
11199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool));
11209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool));
11219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool));
11229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool));
11239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool));
1124c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
1125c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
1126c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
1127c9a4fd69SAudic XU   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL));
11283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1129c7a4214aSPierre Jolivet }
1130