xref: /petsc/src/mat/impls/htool/htool.cxx (revision 462c564db5c15796f1f07c2b71ff29812cdd1aad)
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));
30c7a4214aSPierre Jolivet   a->hmatrix->copy_local_diagonal(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;
3966b4df27SPierre Jolivet   PetscScalar *ptr, scale, shift;
40c7a4214aSPierre Jolivet   PetscBool    flg;
41c7a4214aSPierre Jolivet 
42c7a4214aSPierre Jolivet   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
445f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
45c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
469566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
47c7a4214aSPierre Jolivet   if (!B) {
4866b4df27SPierre 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));
49f22e26b7SPierre Jolivet     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B));
509566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B, &ptr));
51c7a4214aSPierre Jolivet     a->hmatrix->copy_local_diagonal_block(ptr);
529566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B, &ptr));
539566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A, B));
549566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B));
55c7a4214aSPierre Jolivet     *b = B;
569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
5766b4df27SPierre Jolivet     PetscCall(MatShift(*b, shift));
5866b4df27SPierre Jolivet     PetscCall(MatScale(*b, scale));
59c9a4fd69SAudic XU   } else {
6066b4df27SPierre 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));
61c9a4fd69SAudic XU     *b = B;
62c9a4fd69SAudic XU   }
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c7a4214aSPierre Jolivet }
65c7a4214aSPierre Jolivet 
66d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y)
67d71ae5a4SJacob Faibussowitsch {
68c9a4fd69SAudic XU   Mat_Htool         *a;
69c7a4214aSPierre Jolivet   const PetscScalar *in;
70c7a4214aSPierre Jolivet   PetscScalar       *out;
71c7a4214aSPierre Jolivet 
72c7a4214aSPierre Jolivet   PetscFunctionBegin;
73c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
76c7a4214aSPierre Jolivet   a->hmatrix->mvprod_local_to_local(in, out);
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80c7a4214aSPierre Jolivet }
81c7a4214aSPierre Jolivet 
82d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y)
83d71ae5a4SJacob Faibussowitsch {
84c9a4fd69SAudic XU   Mat_Htool         *a;
858b8ddd11SPierre Jolivet   const PetscScalar *in;
868b8ddd11SPierre Jolivet   PetscScalar       *out;
878b8ddd11SPierre Jolivet 
888b8ddd11SPierre Jolivet   PetscFunctionBegin;
89c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
928b8ddd11SPierre Jolivet   a->hmatrix->mvprod_transp_local_to_local(in, out);
939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
968b8ddd11SPierre Jolivet }
978b8ddd11SPierre Jolivet 
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov)
99d71ae5a4SJacob Faibussowitsch {
100c7a4214aSPierre Jolivet   std::set<PetscInt> set;
101c7a4214aSPierre Jolivet   const PetscInt    *idx;
1028f308287SPierre Jolivet   PetscInt          *oidx, size, bs[2];
103c7a4214aSPierre Jolivet   PetscMPIInt        csize;
104c7a4214aSPierre Jolivet 
105c7a4214aSPierre Jolivet   PetscFunctionBegin;
1069566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, bs, bs + 1));
1078f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
108c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < is_max; ++i) {
109c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
110c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1119566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize));
1125f80ce2aSJacob Faibussowitsch     PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported parallel IS");
1139566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i], &size));
1149566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx));
115c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < size; ++j) {
116c7a4214aSPierre Jolivet       set.insert(idx[j]);
117c7a4214aSPierre Jolivet       for (PetscInt k = 1; k <= ov; ++k) {                                              /* for each layer of overlap      */
118c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k);                                    /* do not insert negative indices */
119c7a4214aSPierre 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 */
120c7a4214aSPierre Jolivet       }
121c7a4214aSPierre Jolivet     }
1229566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i], &idx));
1239566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is + i));
1248f308287SPierre Jolivet     if (bs[0] > 1) {
1258f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) {
1268f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1278f308287SPierre Jolivet         std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]);
1288f308287SPierre Jolivet         set.insert(block.cbegin(), block.cend());
1298f308287SPierre Jolivet       }
1308f308287SPierre Jolivet     }
131c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &oidx));
133c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
134c7a4214aSPierre Jolivet     oidx -= size;
1359566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i));
136c7a4214aSPierre Jolivet   }
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138c7a4214aSPierre Jolivet }
139c7a4214aSPierre Jolivet 
140d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
141d71ae5a4SJacob Faibussowitsch {
142c9a4fd69SAudic XU   Mat_Htool         *a;
143c7a4214aSPierre Jolivet   Mat                D, B, BT;
144c7a4214aSPierre Jolivet   const PetscScalar *copy;
14566b4df27SPierre Jolivet   PetscScalar       *ptr, shift, scale;
146c7a4214aSPierre Jolivet   const PetscInt    *idxr, *idxc, *it;
147c7a4214aSPierre Jolivet   PetscInt           nrow, m, i;
148c7a4214aSPierre Jolivet   PetscBool          flg;
149c7a4214aSPierre Jolivet 
150c7a4214aSPierre Jolivet   PetscFunctionBegin;
15166b4df27SPierre 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));
152c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
15348a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
154c7a4214aSPierre Jolivet   for (i = 0; i < n; ++i) {
1559566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &nrow));
1569566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &m));
1579566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &idxr));
1589566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &idxc));
159f22e26b7SPierre Jolivet     if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i));
1609566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr));
161c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1629566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A, &flg));
163c7a4214aSPierre Jolivet       if (flg) {
1649566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i], &flg));
165c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
166c7a4214aSPierre Jolivet           it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart);
167c7a4214aSPierre Jolivet           if (it != idxr + nrow && *it == A->rmap->rstart) {    /* rmap->rstart in IS? */
168c7a4214aSPierre Jolivet             if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
1699371c9d4SSatish Balay               for (PetscInt j = 0; j < A->rmap->n && flg; ++j)
1709371c9d4SSatish Balay                 if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE;
171c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
172c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
173c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
174c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
175c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
176c7a4214aSPierre Jolivet                  */
177c7a4214aSPierre Jolivet                 m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */
178c9a4fd69SAudic XU                 PetscCall(MatGetDiagonalBlock(A, &D));
1799566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D, &copy));
1809371c9d4SSatish 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 */ }
1819566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D, &copy));
182c7a4214aSPierre Jolivet                 if (m) {
183c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */
184c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
185b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
1869566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B));
1879566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
1889566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT));
1899566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
190b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
1919566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
192c7a4214aSPierre Jolivet                     } else {
1937fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
1949566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
195c7a4214aSPierre Jolivet                     }
1969566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
1979566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
198c7a4214aSPierre Jolivet                   } else {
199c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */
200c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow);
201c7a4214aSPierre Jolivet                     }
202c7a4214aSPierre Jolivet                   }
203c7a4214aSPierre Jolivet                 }
204c7a4214aSPierre Jolivet                 if (m + A->rmap->n != nrow) {
205c7a4214aSPierre 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 */
206c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
207b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
2089566063dSJacob 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));
2099566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
2109566063dSJacob 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));
2119566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
212b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
2139566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
214c7a4214aSPierre Jolivet                     } else {
2157fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
2169566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
217c7a4214aSPierre Jolivet                     }
2189566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2199566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
220c7a4214aSPierre Jolivet                   } else {
221c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */
222c7a4214aSPierre 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);
223c7a4214aSPierre Jolivet                     }
224c7a4214aSPierre Jolivet                   }
225c7a4214aSPierre Jolivet                 }
226c7a4214aSPierre Jolivet               } /* complete local diagonal block not in IS */
227c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
228c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE;   /* rmap->rstart not in IS */
229c7a4214aSPierre Jolivet         } /* unsorted IS */
230c7a4214aSPierre Jolivet       }
231c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE;                                       /* different row and column IS */
232c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */
2339566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &idxr));
2349566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &idxc));
2359566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr));
23666b4df27SPierre Jolivet     PetscCall(MatShift((*submat)[i], shift));
23766b4df27SPierre Jolivet     PetscCall(MatScale((*submat)[i], scale));
238c7a4214aSPierre Jolivet   }
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240c7a4214aSPierre Jolivet }
241c7a4214aSPierre Jolivet 
242d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A)
243d71ae5a4SJacob Faibussowitsch {
244c9a4fd69SAudic XU   Mat_Htool               *a;
245c7a4214aSPierre Jolivet   PetscContainer           container;
246c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
247c7a4214aSPierre Jolivet 
248c7a4214aSPierre Jolivet   PetscFunctionBegin;
249c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
250f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr));
251f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr));
252f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr));
253f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr));
254f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr));
255f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr));
256f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr));
257f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr));
258f22e26b7SPierre Jolivet   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr));
2599566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container));
260c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2619566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
2629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
2639566063dSJacob Faibussowitsch     PetscCall(PetscFree(kernelt));
2649566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
265f22e26b7SPierre Jolivet     PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr));
266c7a4214aSPierre Jolivet   }
26748a46eb9SPierre Jolivet   if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source));
2689566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->work_source, a->work_target));
270c7a4214aSPierre Jolivet   delete a->wrapper;
271c7a4214aSPierre Jolivet   delete a->hmatrix;
272c9a4fd69SAudic XU   PetscCall(PetscFree(a));
273c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable()
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275c7a4214aSPierre Jolivet }
276c7a4214aSPierre Jolivet 
277d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv)
278d71ae5a4SJacob Faibussowitsch {
279c9a4fd69SAudic XU   Mat_Htool  *a;
28066b4df27SPierre Jolivet   PetscScalar shift, scale;
281c7a4214aSPierre Jolivet   PetscBool   flg;
282c7a4214aSPierre Jolivet 
283c7a4214aSPierre Jolivet   PetscFunctionBegin;
284c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg));
286c7a4214aSPierre Jolivet   if (flg) {
28766b4df27SPierre 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));
2889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->hmatrix->get_symmetry_type()));
28966b4df27SPierre Jolivet     if (PetscAbsScalar(scale - 1.0) > PETSC_MACHINE_EPSILON) {
290c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
29166b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(scale), (double)PetscImaginaryPart(scale)));
292c7a4214aSPierre Jolivet #else
29366b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)scale));
294c9a4fd69SAudic XU #endif
295c9a4fd69SAudic XU     }
29666b4df27SPierre Jolivet     if (PetscAbsScalar(shift) > PETSC_MACHINE_EPSILON) {
297c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX)
29866b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(shift), (double)PetscImaginaryPart(shift)));
299c9a4fd69SAudic XU #else
30066b4df27SPierre Jolivet       PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)shift));
301c7a4214aSPierre Jolivet #endif
302c7a4214aSPierre Jolivet     }
3039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum cluster size: %" PetscInt_FMT "\n", a->bs[0]));
3049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "maximum block size: %" PetscInt_FMT "\n", a->bs[1]));
3059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon));
3069566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta));
3079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0]));
3089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1]));
3099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor]));
3109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering]));
3119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", a->hmatrix->get_infos("Compression_ratio").c_str()));
3129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", a->hmatrix->get_infos("Space_saving").c_str()));
3139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "number of dense (resp. low rank) matrices: %s (resp. %s)\n", a->hmatrix->get_infos("Number_of_dmat").c_str(), a->hmatrix->get_infos("Number_of_lrmat").c_str()));
3149371c9d4SSatish Balay     PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) dense block sizes: (%s, %s, %s)\n", a->hmatrix->get_infos("Dense_block_size_min").c_str(), a->hmatrix->get_infos("Dense_block_size_mean").c_str(),
3159371c9d4SSatish Balay                                      a->hmatrix->get_infos("Dense_block_size_max").c_str()));
3169371c9d4SSatish Balay     PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) low rank block sizes: (%s, %s, %s)\n", a->hmatrix->get_infos("Low_rank_block_size_min").c_str(), a->hmatrix->get_infos("Low_rank_block_size_mean").c_str(),
3179371c9d4SSatish Balay                                      a->hmatrix->get_infos("Low_rank_block_size_max").c_str()));
3189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) ranks: (%s, %s, %s)\n", a->hmatrix->get_infos("Rank_min").c_str(), a->hmatrix->get_infos("Rank_mean").c_str(), a->hmatrix->get_infos("Rank_max").c_str()));
319c7a4214aSPierre Jolivet   }
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
321c7a4214aSPierre Jolivet }
322c7a4214aSPierre Jolivet 
323c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
324d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
325d71ae5a4SJacob Faibussowitsch {
326c9a4fd69SAudic XU   Mat_Htool   *a;
32766b4df27SPierre Jolivet   PetscScalar  shift, scale;
328c7a4214aSPierre Jolivet   PetscInt    *idxc;
329c7a4214aSPierre Jolivet   PetscBLASInt one = 1, bn;
330c7a4214aSPierre Jolivet 
331c7a4214aSPierre Jolivet   PetscFunctionBegin;
33266b4df27SPierre 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));
333c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
334c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
335c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, &idxc));
337c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i;
338c7a4214aSPierre Jolivet   }
339c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
340c7a4214aSPierre Jolivet   if (v) {
3419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, v));
342c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
343cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
3449566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N, &bn));
34566b4df27SPierre Jolivet     PetscCallBLAS("BLASscal", BLASscal_(&bn, &scale, *v, &one));
34666b4df27SPierre Jolivet     if (row < A->cmap->N) (*v)[row] += shift;
347c7a4214aSPierre Jolivet   }
34848a46eb9SPierre Jolivet   if (!idx) PetscCall(PetscFree(idxc));
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
350c7a4214aSPierre Jolivet }
351c7a4214aSPierre Jolivet 
352f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v)
353d71ae5a4SJacob Faibussowitsch {
354c7a4214aSPierre Jolivet   PetscFunctionBegin;
35548a46eb9SPierre Jolivet   if (idx) PetscCall(PetscFree(*idx));
35648a46eb9SPierre Jolivet   if (v) PetscCall(PetscFree(*v));
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358c7a4214aSPierre Jolivet }
359c7a4214aSPierre Jolivet 
360d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject)
361d71ae5a4SJacob Faibussowitsch {
362c9a4fd69SAudic XU   Mat_Htool *a;
363c7a4214aSPierre Jolivet   PetscInt   n;
364c7a4214aSPierre Jolivet   PetscBool  flg;
365c7a4214aSPierre Jolivet 
366c7a4214aSPierre Jolivet   PetscFunctionBegin;
367c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
368d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Htool options");
369f22e26b7SPierre Jolivet   PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", nullptr, a->bs[0], a->bs, nullptr));
370f22e26b7SPierre Jolivet   PetscCall(PetscOptionsInt("-mat_htool_max_block_size", "Maximal number of coefficients in a dense block", nullptr, a->bs[1], a->bs + 1, nullptr));
371f22e26b7SPierre Jolivet   PetscCall(PetscOptionsReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr));
372f22e26b7SPierre Jolivet   PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr));
373f22e26b7SPierre Jolivet   PetscCall(PetscOptionsInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr));
374f22e26b7SPierre Jolivet   PetscCall(PetscOptionsInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", nullptr, a->depth[1], a->depth + 1, nullptr));
375c7a4214aSPierre Jolivet   n = 0;
376dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg));
377c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
37898e73e17SPierre Jolivet   n = 0;
379dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg));
38098e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
381d0609cedSBarry Smith   PetscOptionsHeadEnd();
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383c7a4214aSPierre Jolivet }
384c7a4214aSPierre Jolivet 
385cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType)
386d71ae5a4SJacob Faibussowitsch {
387c9a4fd69SAudic XU   Mat_Htool                                                   *a;
388c7a4214aSPierre Jolivet   const PetscInt                                              *ranges;
389c7a4214aSPierre Jolivet   PetscInt                                                    *offset;
390c7a4214aSPierre Jolivet   PetscMPIInt                                                  size;
391b94d7dedSBarry Smith   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U';
392cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                        *generator = nullptr;
39398e73e17SPierre Jolivet   std::shared_ptr<htool::VirtualCluster>                       t, s = nullptr;
3943b9338faSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr;
395c7a4214aSPierre Jolivet 
396c7a4214aSPierre Jolivet   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite));
398c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
399c7a4214aSPierre Jolivet   delete a->wrapper;
400c7a4214aSPierre Jolivet   delete a->hmatrix;
4019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
4029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &offset));
4039566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
404c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < size; ++i) {
405c7a4214aSPierre Jolivet     offset[2 * i]     = ranges[i];
406c7a4214aSPierre Jolivet     offset[2 * i + 1] = ranges[i + 1] - ranges[i];
407c7a4214aSPierre Jolivet   }
40898e73e17SPierre Jolivet   switch (a->clustering) {
409d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
410d71ae5a4SJacob Faibussowitsch     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
411d71ae5a4SJacob Faibussowitsch     break;
412d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
413d71ae5a4SJacob Faibussowitsch     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
414d71ae5a4SJacob Faibussowitsch     break;
415d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
416d71ae5a4SJacob Faibussowitsch     t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
417d71ae5a4SJacob Faibussowitsch     break;
418d71ae5a4SJacob Faibussowitsch   default:
419d71ae5a4SJacob Faibussowitsch     t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
42098e73e17SPierre Jolivet   }
421c7a4214aSPierre Jolivet   t->set_minclustersize(a->bs[0]);
4220429413eSPierre Jolivet   t->build(A->rmap->N, a->gcoords_target, offset, -1, PetscObjectComm((PetscObject)A));
423c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx);
424c7a4214aSPierre Jolivet   else {
425f22e26b7SPierre Jolivet     a->wrapper = nullptr;
426cab00e0dSPierre Jolivet     generator  = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx);
427c7a4214aSPierre Jolivet   }
428c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4299566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A, &ranges));
430c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < size; ++i) {
431c7a4214aSPierre Jolivet       offset[2 * i]     = ranges[i];
432c7a4214aSPierre Jolivet       offset[2 * i + 1] = ranges[i + 1] - ranges[i];
433c7a4214aSPierre Jolivet     }
43498e73e17SPierre Jolivet     switch (a->clustering) {
435d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC:
436d71ae5a4SJacob Faibussowitsch       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
437d71ae5a4SJacob Faibussowitsch       break;
438d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC:
439d71ae5a4SJacob Faibussowitsch       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim);
440d71ae5a4SJacob Faibussowitsch       break;
441d71ae5a4SJacob Faibussowitsch     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR:
442d71ae5a4SJacob Faibussowitsch       s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim);
443d71ae5a4SJacob Faibussowitsch       break;
444d71ae5a4SJacob Faibussowitsch     default:
445d71ae5a4SJacob Faibussowitsch       s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
44698e73e17SPierre Jolivet     }
447c7a4214aSPierre Jolivet     s->set_minclustersize(a->bs[0]);
4480429413eSPierre Jolivet     s->build(A->cmap->N, a->gcoords_source, offset, -1, PetscObjectComm((PetscObject)A));
449c7a4214aSPierre Jolivet     S = uplo = 'N';
450c7a4214aSPierre Jolivet   }
4519566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
452c7a4214aSPierre Jolivet   switch (a->compressor) {
453d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_FULL_ACA:
454d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::fullACA<PetscScalar>>();
455d71ae5a4SJacob Faibussowitsch     break;
456d71ae5a4SJacob Faibussowitsch   case MAT_HTOOL_COMPRESSOR_SVD:
457d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::SVD<PetscScalar>>();
458d71ae5a4SJacob Faibussowitsch     break;
459d71ae5a4SJacob Faibussowitsch   default:
460d71ae5a4SJacob Faibussowitsch     compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
461c7a4214aSPierre Jolivet   }
4620429413eSPierre Jolivet   a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar> *>(new htool::HMatrix<PetscScalar>(t, s ? s : t, a->epsilon, a->eta, S, uplo, -1, PetscObjectComm((PetscObject)A)));
4633b9338faSPierre Jolivet   a->hmatrix->set_compression(compressor);
464c7a4214aSPierre Jolivet   a->hmatrix->set_maxblocksize(a->bs[1]);
465c7a4214aSPierre Jolivet   a->hmatrix->set_mintargetdepth(a->depth[0]);
466c7a4214aSPierre Jolivet   a->hmatrix->set_minsourcedepth(a->depth[1]);
4673b9338faSPierre Jolivet   if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target, a->gcoords_source);
4683b9338faSPierre Jolivet   else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target);
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
470c7a4214aSPierre Jolivet }
471c7a4214aSPierre Jolivet 
472d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C)
473d71ae5a4SJacob Faibussowitsch {
474c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
475c9a4fd69SAudic XU   Mat_Htool         *a;
476c7a4214aSPierre Jolivet   const PetscScalar *in;
477c7a4214aSPierre Jolivet   PetscScalar       *out;
4788b8ddd11SPierre Jolivet   PetscInt           N, lda;
479c7a4214aSPierre Jolivet 
480c7a4214aSPierre Jolivet   PetscFunctionBegin;
481c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
482f22e26b7SPierre Jolivet   PetscCall(MatGetSize(C, nullptr, &N));
4839566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &lda));
4845f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
4859566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B, &in));
4869566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &out));
487c9a4fd69SAudic XU   PetscCall(MatShellGetContext(product->A, &a));
4888b8ddd11SPierre Jolivet   switch (product->type) {
489d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AB:
490d71ae5a4SJacob Faibussowitsch     a->hmatrix->mvprod_local_to_local(in, out, N);
491d71ae5a4SJacob Faibussowitsch     break;
492d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AtB:
493d71ae5a4SJacob Faibussowitsch     a->hmatrix->mvprod_transp_local_to_local(in, out, N);
494d71ae5a4SJacob Faibussowitsch     break;
495d71ae5a4SJacob Faibussowitsch   default:
496d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]);
497c7a4214aSPierre Jolivet   }
4989566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &out));
4999566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B, &in));
5003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
501c7a4214aSPierre Jolivet }
502c7a4214aSPierre Jolivet 
503d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C)
504d71ae5a4SJacob Faibussowitsch {
505c7a4214aSPierre Jolivet   Mat_Product *product = C->product;
506c7a4214aSPierre Jolivet   Mat          A, B;
507c7a4214aSPierre Jolivet   PetscBool    flg;
508c7a4214aSPierre Jolivet 
509c7a4214aSPierre Jolivet   PetscFunctionBegin;
510c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
511c7a4214aSPierre Jolivet   A = product->A;
512c7a4214aSPierre Jolivet   B = product->B;
5139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, ""));
51497713f8cSPierre 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);
51597713f8cSPierre Jolivet   if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) {
51697713f8cSPierre Jolivet     if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
51797713f8cSPierre Jolivet     else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
5188b8ddd11SPierre Jolivet   }
5199566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
5209566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
5219566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
5229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
5239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
524f22e26b7SPierre Jolivet   C->ops->productsymbolic = nullptr;
525c7a4214aSPierre Jolivet   C->ops->productnumeric  = MatProductNumeric_Htool;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527c7a4214aSPierre Jolivet }
528c7a4214aSPierre Jolivet 
529d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C)
530d71ae5a4SJacob Faibussowitsch {
531c7a4214aSPierre Jolivet   PetscFunctionBegin;
532c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
5338b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535c7a4214aSPierre Jolivet }
536c7a4214aSPierre Jolivet 
537d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix)
538d71ae5a4SJacob Faibussowitsch {
539c9a4fd69SAudic XU   Mat_Htool *a;
540c7a4214aSPierre Jolivet 
541c7a4214aSPierre Jolivet   PetscFunctionBegin;
542c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
543c7a4214aSPierre Jolivet   *hmatrix = a->hmatrix;
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
545c7a4214aSPierre Jolivet }
546c7a4214aSPierre Jolivet 
547c7a4214aSPierre Jolivet /*@C
54811a5261eSBarry Smith   MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`.
549c7a4214aSPierre Jolivet 
550cc4c1da9SBarry Smith   No Fortran Support, No C Support
551cc4c1da9SBarry Smith 
552c7a4214aSPierre Jolivet   Input Parameter:
553c7a4214aSPierre Jolivet . A - hierarchical matrix
554c7a4214aSPierre Jolivet 
555c7a4214aSPierre Jolivet   Output Parameter:
556c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix
557c7a4214aSPierre Jolivet 
558c7a4214aSPierre Jolivet   Level: advanced
559c7a4214aSPierre Jolivet 
5601cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`
561c7a4214aSPierre Jolivet @*/
562d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix)
563d71ae5a4SJacob Faibussowitsch {
564c7a4214aSPierre Jolivet   PetscFunctionBegin;
565c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
5664f572ea9SToby Isaac   PetscAssertPointer(hmatrix, 2);
567cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::VirtualHMatrix<PetscScalar> **), (A, hmatrix));
5683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
569c7a4214aSPierre Jolivet }
570c7a4214aSPierre Jolivet 
5718434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
572d71ae5a4SJacob Faibussowitsch {
573c9a4fd69SAudic XU   Mat_Htool *a;
574c7a4214aSPierre Jolivet 
575c7a4214aSPierre Jolivet   PetscFunctionBegin;
576c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
577c7a4214aSPierre Jolivet   a->kernel    = kernel;
578c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
579c7a4214aSPierre Jolivet   delete a->wrapper;
580c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx);
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582c7a4214aSPierre Jolivet }
583c7a4214aSPierre Jolivet 
584c7a4214aSPierre Jolivet /*@C
58511a5261eSBarry Smith   MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`.
586c7a4214aSPierre Jolivet 
587cc4c1da9SBarry Smith   Collective, No Fortran Support
588cc4c1da9SBarry Smith 
589c7a4214aSPierre Jolivet   Input Parameters:
590c7a4214aSPierre Jolivet + A         - hierarchical matrix
5912ef1f0ffSBarry Smith . kernel    - computational kernel (or `NULL`)
5922ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
593c7a4214aSPierre Jolivet 
594c7a4214aSPierre Jolivet   Level: advanced
595c7a4214aSPierre Jolivet 
5961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()`
597c7a4214aSPierre Jolivet @*/
598cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx)
599d71ae5a4SJacob Faibussowitsch {
600c7a4214aSPierre Jolivet   PetscFunctionBegin;
601c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
602c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 2);
6034f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 3);
6048434afd1SBarry Smith   PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx));
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606c7a4214aSPierre Jolivet }
607c7a4214aSPierre Jolivet 
608d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is)
609d71ae5a4SJacob Faibussowitsch {
610c9a4fd69SAudic XU   Mat_Htool            *a;
61198e73e17SPierre Jolivet   std::vector<PetscInt> source;
61298e73e17SPierre Jolivet 
61398e73e17SPierre Jolivet   PetscFunctionBegin;
614c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
615a9918087SPierre Jolivet   source = a->hmatrix->get_source_cluster()->get_local_perm();
6169566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), source.size(), source.data(), PETSC_COPY_VALUES, is));
6179566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61998e73e17SPierre Jolivet }
62098e73e17SPierre Jolivet 
621cc4c1da9SBarry Smith /*@
62211a5261eSBarry Smith   MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix.
62398e73e17SPierre Jolivet 
62498e73e17SPierre Jolivet   Input Parameter:
62598e73e17SPierre Jolivet . A - hierarchical matrix
62698e73e17SPierre Jolivet 
62798e73e17SPierre Jolivet   Output Parameter:
62898e73e17SPierre Jolivet . is - permutation
62998e73e17SPierre Jolivet 
63098e73e17SPierre Jolivet   Level: advanced
63198e73e17SPierre Jolivet 
6321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()`
63398e73e17SPierre Jolivet @*/
634cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is)
635d71ae5a4SJacob Faibussowitsch {
63698e73e17SPierre Jolivet   PetscFunctionBegin;
63798e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6384f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
639cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is));
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64198e73e17SPierre Jolivet }
64298e73e17SPierre Jolivet 
643d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is)
644d71ae5a4SJacob Faibussowitsch {
645c9a4fd69SAudic XU   Mat_Htool            *a;
64698e73e17SPierre Jolivet   std::vector<PetscInt> target;
64798e73e17SPierre Jolivet 
64898e73e17SPierre Jolivet   PetscFunctionBegin;
649c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
650a9918087SPierre Jolivet   target = a->hmatrix->get_target_cluster()->get_local_perm();
6519566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), target.size(), target.data(), PETSC_COPY_VALUES, is));
6529566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
6533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65498e73e17SPierre Jolivet }
65598e73e17SPierre Jolivet 
656cc4c1da9SBarry Smith /*@
65711a5261eSBarry Smith   MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix.
65898e73e17SPierre Jolivet 
65998e73e17SPierre Jolivet   Input Parameter:
66098e73e17SPierre Jolivet . A - hierarchical matrix
66198e73e17SPierre Jolivet 
66298e73e17SPierre Jolivet   Output Parameter:
66398e73e17SPierre Jolivet . is - permutation
66498e73e17SPierre Jolivet 
66598e73e17SPierre Jolivet   Level: advanced
66698e73e17SPierre Jolivet 
6671cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()`
66898e73e17SPierre Jolivet @*/
669cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is)
670d71ae5a4SJacob Faibussowitsch {
67198e73e17SPierre Jolivet   PetscFunctionBegin;
67298e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
6734f572ea9SToby Isaac   if (!is) PetscAssertPointer(is, 2);
674cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67698e73e17SPierre Jolivet }
67798e73e17SPierre Jolivet 
678d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use)
679d71ae5a4SJacob Faibussowitsch {
680c9a4fd69SAudic XU   Mat_Htool *a;
68198e73e17SPierre Jolivet 
68298e73e17SPierre Jolivet   PetscFunctionBegin;
683c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
68498e73e17SPierre Jolivet   a->hmatrix->set_use_permutation(use);
6853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68698e73e17SPierre Jolivet }
68798e73e17SPierre Jolivet 
688cc4c1da9SBarry Smith /*@
68911a5261eSBarry Smith   MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation.
69098e73e17SPierre Jolivet 
69198e73e17SPierre Jolivet   Input Parameters:
69298e73e17SPierre Jolivet + A   - hierarchical matrix
69398e73e17SPierre Jolivet - use - Boolean value
69498e73e17SPierre Jolivet 
69598e73e17SPierre Jolivet   Level: advanced
69698e73e17SPierre Jolivet 
6971cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()`
69898e73e17SPierre Jolivet @*/
699cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use)
700d71ae5a4SJacob Faibussowitsch {
70198e73e17SPierre Jolivet   PetscFunctionBegin;
70298e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
70398e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A, use, 2);
704cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use));
7053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70698e73e17SPierre Jolivet }
70798e73e17SPierre Jolivet 
708cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B)
709d71ae5a4SJacob Faibussowitsch {
710c7a4214aSPierre Jolivet   Mat          C;
711c9a4fd69SAudic XU   Mat_Htool   *a;
71266b4df27SPierre Jolivet   PetscScalar *array, shift, scale;
713c7a4214aSPierre Jolivet   PetscInt     lda;
714c7a4214aSPierre Jolivet 
715c7a4214aSPierre Jolivet   PetscFunctionBegin;
71666b4df27SPierre 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));
717c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
718c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
719c7a4214aSPierre Jolivet     C = *B;
7205f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions");
7219566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C, &lda));
7225f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
723c7a4214aSPierre Jolivet   } else {
7249566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
7259566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
7269566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, MATDENSE));
7279566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7289566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
729c7a4214aSPierre Jolivet   }
7309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
7319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
7329566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
733c7a4214aSPierre Jolivet   a->hmatrix->copy_local_dense_perm(array);
7349566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
73566b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
73666b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
737c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
7389566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &C));
739c7a4214aSPierre Jolivet   } else *B = C;
7403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
741c7a4214aSPierre Jolivet }
742c7a4214aSPierre Jolivet 
743d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx)
744d71ae5a4SJacob Faibussowitsch {
745c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx;
74698e73e17SPierre Jolivet   PetscScalar             *tmp;
74798e73e17SPierre Jolivet 
74898e73e17SPierre Jolivet   PetscFunctionBegin;
7493ba16761SJacob Faibussowitsch   PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx));
7509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M * N, &tmp));
7519566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp, ptr, M * N));
75298e73e17SPierre Jolivet   for (PetscInt i = 0; i < M; ++i) {
75398e73e17SPierre Jolivet     for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N];
75498e73e17SPierre Jolivet   }
7559566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757c7a4214aSPierre Jolivet }
758c7a4214aSPierre Jolivet 
759c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
760d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B)
761d71ae5a4SJacob Faibussowitsch {
762c7a4214aSPierre Jolivet   Mat                      C;
763c9a4fd69SAudic XU   Mat_Htool               *a, *c;
76466b4df27SPierre Jolivet   PetscScalar              shift, scale;
765c7a4214aSPierre Jolivet   PetscInt                 M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n;
766c7a4214aSPierre Jolivet   PetscContainer           container;
767c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
768c7a4214aSPierre Jolivet 
769c7a4214aSPierre Jolivet   PetscFunctionBegin;
77066b4df27SPierre 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));
771c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
7727fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
7735f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported");
774c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
7759566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
7769566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, n, m, N, M));
7779566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
7789566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7799566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C), &container));
7809566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
7819566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container, kernelt));
7829566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)C, "KernelTranspose", (PetscObject)container));
783c7a4214aSPierre Jolivet   } else {
784c7a4214aSPierre Jolivet     C = *B;
7859566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container));
7865f80ce2aSJacob Faibussowitsch     PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first");
7879566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
788c7a4214aSPierre Jolivet   }
789c9a4fd69SAudic XU   PetscCall(MatShellGetContext(C, &c));
790c7a4214aSPierre Jolivet   c->dim = a->dim;
79166b4df27SPierre Jolivet   PetscCall(MatShift(C, shift));
79266b4df27SPierre Jolivet   PetscCall(MatScale(C, scale));
79398e73e17SPierre Jolivet   c->kernel = GenEntriesTranspose;
794c7a4214aSPierre Jolivet   if (kernelt->A != A) {
7959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
796c7a4214aSPierre Jolivet     kernelt->A = A;
7979566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
798c7a4214aSPierre Jolivet   }
799c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
800c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
801c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
802c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
8039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target));
8049566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim));
805c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
8069566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source));
8079566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim));
808c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
8099566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target));
810c7a4214aSPierre Jolivet   }
8119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
8129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
813c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
815c7a4214aSPierre Jolivet }
816c7a4214aSPierre Jolivet 
817c7a4214aSPierre Jolivet /*@C
81811a5261eSBarry Smith   MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel.
819c7a4214aSPierre Jolivet 
820cc4c1da9SBarry Smith   Collective, No Fortran Support
821cc4c1da9SBarry Smith 
822c7a4214aSPierre Jolivet   Input Parameters:
823c7a4214aSPierre Jolivet + comm          - MPI communicator
8242ef1f0ffSBarry Smith . m             - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
8252ef1f0ffSBarry Smith . n             - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
8262ef1f0ffSBarry Smith . M             - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
8272ef1f0ffSBarry Smith . N             - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
828c7a4214aSPierre Jolivet . spacedim      - dimension of the space coordinates
829c7a4214aSPierre Jolivet . coords_target - coordinates of the target
830c7a4214aSPierre Jolivet . coords_source - coordinates of the source
8312ef1f0ffSBarry Smith . kernel        - computational kernel (or `NULL`)
8322ef1f0ffSBarry Smith - kernelctx     - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
833c7a4214aSPierre Jolivet 
834c7a4214aSPierre Jolivet   Output Parameter:
835c7a4214aSPierre Jolivet . B - matrix
836c7a4214aSPierre Jolivet 
837c7a4214aSPierre Jolivet   Options Database Keys:
83811a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`>                                                     - minimal leaf size in cluster tree
83911a5261eSBarry Smith . -mat_htool_max_block_size <`PetscInt`>                                                       - maximal number of coefficients in a dense block
84011a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`>                                                             - relative error in Frobenius norm when approximating a block
84111a5261eSBarry Smith . -mat_htool_eta <`PetscReal`>                                                                 - admissibility condition tolerance
84211a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the rows
84311a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`>                                                     - minimal cluster tree depth associated with the columns
84498e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD>                                          - type of compression
84598e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
846c7a4214aSPierre Jolivet 
847c7a4214aSPierre Jolivet   Level: intermediate
848c7a4214aSPierre Jolivet 
8491cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
850c7a4214aSPierre Jolivet @*/
8518434afd1SBarry 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)
852d71ae5a4SJacob Faibussowitsch {
853c7a4214aSPierre Jolivet   Mat        A;
854c7a4214aSPierre Jolivet   Mat_Htool *a;
855c7a4214aSPierre Jolivet 
856c7a4214aSPierre Jolivet   PetscFunctionBegin;
8579566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
858c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, spacedim, 6);
8594f572ea9SToby Isaac   PetscAssertPointer(coords_target, 7);
8604f572ea9SToby Isaac   PetscAssertPointer(coords_source, 8);
861c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 9);
8624f572ea9SToby Isaac   if (!kernel) PetscAssertPointer(kernelctx, 10);
8639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, n, M, N));
8649566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATHTOOL));
8659566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
866c9a4fd69SAudic XU   PetscCall(MatShellGetContext(A, &a));
867c7a4214aSPierre Jolivet   a->dim       = spacedim;
868c7a4214aSPierre Jolivet   a->kernel    = kernel;
869c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
8709566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target));
8719566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim));
872*462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */
873c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
8749566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source));
8759566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim));
876*462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */
877c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
8789566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target));
879c7a4214aSPierre Jolivet   *B = A;
8803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
881c7a4214aSPierre Jolivet }
882c7a4214aSPierre Jolivet 
883c7a4214aSPierre Jolivet /*MC
884c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
885c7a4214aSPierre Jolivet 
8862ef1f0ffSBarry Smith   Use `./configure --download-htool` to install PETSc to use Htool.
887c7a4214aSPierre Jolivet 
8882ef1f0ffSBarry Smith    Options Database Key:
8892ef1f0ffSBarry Smith .     -mat_type htool - matrix type to `MATHTOOL`
890c7a4214aSPierre Jolivet 
891c7a4214aSPierre Jolivet    Level: beginner
892c7a4214aSPierre Jolivet 
8931cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()`
894c7a4214aSPierre Jolivet M*/
895d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A)
896d71ae5a4SJacob Faibussowitsch {
897c7a4214aSPierre Jolivet   Mat_Htool *a;
898c7a4214aSPierre Jolivet 
899c7a4214aSPierre Jolivet   PetscFunctionBegin;
900c9a4fd69SAudic XU   PetscCall(MatSetType(A, MATSHELL));
9014dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
902c9a4fd69SAudic XU   PetscCall(MatShellSetContext(A, a));
903c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Htool));
904c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Htool));
905c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Htool));
906c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool));
907c9a4fd69SAudic XU   if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool));
908c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
909c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
910c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_Htool));
911c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_Htool));
912c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (void (*)(void))MatGetRow_Htool));
913c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (void (*)(void))MatRestoreRow_Htool));
914c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_Htool));
915c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_Htool));
916c9a4fd69SAudic XU   PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_Htool));
917c7a4214aSPierre Jolivet   a->dim            = 0;
918f22e26b7SPierre Jolivet   a->gcoords_target = nullptr;
919f22e26b7SPierre Jolivet   a->gcoords_source = nullptr;
920c7a4214aSPierre Jolivet   a->bs[0]          = 10;
921c7a4214aSPierre Jolivet   a->bs[1]          = 1000000;
922c7a4214aSPierre Jolivet   a->epsilon        = PetscSqrtReal(PETSC_SMALL);
923c7a4214aSPierre Jolivet   a->eta            = 10.0;
924c7a4214aSPierre Jolivet   a->depth[0]       = 0;
925c7a4214aSPierre Jolivet   a->depth[1]       = 0;
926c7a4214aSPierre Jolivet   a->compressor     = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
9279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool));
9289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool));
9299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense));
9309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense));
9319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool));
9329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool));
9339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool));
9349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool));
9359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool));
936c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable));
937c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable));
938c9a4fd69SAudic XU   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable));
939c9a4fd69SAudic XU   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL));
9403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
941c7a4214aSPierre Jolivet }
942