xref: /petsc/src/mat/impls/htool/htool.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
199371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_Htool(Mat A, Vec v) {
20c7a4214aSPierre Jolivet   Mat_Htool   *a = (Mat_Htool *)A->data;
21c7a4214aSPierre Jolivet   PetscScalar *x;
22c7a4214aSPierre Jolivet   PetscBool    flg;
23c7a4214aSPierre Jolivet 
24c7a4214aSPierre Jolivet   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
265f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(v, &x));
28c7a4214aSPierre Jolivet   a->hmatrix->copy_local_diagonal(x);
299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(v, &x));
309566063dSJacob Faibussowitsch   PetscCall(VecScale(v, a->s));
31c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
32c7a4214aSPierre Jolivet }
33c7a4214aSPierre Jolivet 
349371c9d4SSatish Balay static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b) {
35c7a4214aSPierre Jolivet   Mat_Htool   *a = (Mat_Htool *)A->data;
36c7a4214aSPierre Jolivet   Mat          B;
37c7a4214aSPierre Jolivet   PetscScalar *ptr;
38c7a4214aSPierre Jolivet   PetscBool    flg;
39c7a4214aSPierre Jolivet 
40c7a4214aSPierre Jolivet   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &flg));
425f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported");
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */
44c7a4214aSPierre Jolivet   if (!B) {
459566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, NULL, &B));
469566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(B, &ptr));
47c7a4214aSPierre Jolivet     a->hmatrix->copy_local_diagonal_block(ptr);
489566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(B, &ptr));
499566063dSJacob Faibussowitsch     PetscCall(MatPropagateSymmetryOptions(A, B));
509566063dSJacob Faibussowitsch     PetscCall(MatScale(B, a->s));
519566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B));
52c7a4214aSPierre Jolivet     *b = B;
539566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&B));
54c7a4214aSPierre Jolivet   } else *b = B;
55c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
56c7a4214aSPierre Jolivet }
57c7a4214aSPierre Jolivet 
589371c9d4SSatish Balay static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y) {
59c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool *)A->data;
60c7a4214aSPierre Jolivet   const PetscScalar *in;
61c7a4214aSPierre Jolivet   PetscScalar       *out;
62c7a4214aSPierre Jolivet 
63c7a4214aSPierre Jolivet   PetscFunctionBegin;
649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
66c7a4214aSPierre Jolivet   a->hmatrix->mvprod_local_to_local(in, out);
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
699566063dSJacob Faibussowitsch   PetscCall(VecScale(y, a->s));
70c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
71c7a4214aSPierre Jolivet }
72c7a4214aSPierre Jolivet 
73c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */
749371c9d4SSatish Balay static PetscErrorCode MatMultAdd_Htool(Mat A, Vec v1, Vec v2, Vec v3) {
75c7a4214aSPierre Jolivet   Mat_Htool        *a = (Mat_Htool *)A->data;
76c7a4214aSPierre Jolivet   Vec               tmp;
77c7a4214aSPierre Jolivet   const PetscScalar scale = a->s;
78c7a4214aSPierre Jolivet 
79c7a4214aSPierre Jolivet   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(v2, &tmp));
819566063dSJacob Faibussowitsch   PetscCall(VecCopy(v2, v3)); /* no-op in MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]) since VecCopy() checks for x == y */
82c7a4214aSPierre Jolivet   a->s = 1.0;                 /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */
839566063dSJacob Faibussowitsch   PetscCall(MatMult_Htool(A, v1, tmp));
849566063dSJacob Faibussowitsch   PetscCall(VecAXPY(v3, scale, tmp));
859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp));
86c7a4214aSPierre Jolivet   a->s = scale; /* set s back to its original value */
87c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
88c7a4214aSPierre Jolivet }
89c7a4214aSPierre Jolivet 
909371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y) {
918b8ddd11SPierre Jolivet   Mat_Htool         *a = (Mat_Htool *)A->data;
928b8ddd11SPierre Jolivet   const PetscScalar *in;
938b8ddd11SPierre Jolivet   PetscScalar       *out;
948b8ddd11SPierre Jolivet 
958b8ddd11SPierre Jolivet   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &in));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(y, &out));
988b8ddd11SPierre Jolivet   a->hmatrix->mvprod_transp_local_to_local(in, out);
999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &in));
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(y, &out));
1019566063dSJacob Faibussowitsch   PetscCall(VecScale(y, a->s));
1028b8ddd11SPierre Jolivet   PetscFunctionReturn(0);
1038b8ddd11SPierre Jolivet }
1048b8ddd11SPierre Jolivet 
1059371c9d4SSatish Balay static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov) {
106c7a4214aSPierre Jolivet   std::set<PetscInt> set;
107c7a4214aSPierre Jolivet   const PetscInt    *idx;
1088f308287SPierre Jolivet   PetscInt          *oidx, size, bs[2];
109c7a4214aSPierre Jolivet   PetscMPIInt        csize;
110c7a4214aSPierre Jolivet 
111c7a4214aSPierre Jolivet   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, bs, bs + 1));
1138f308287SPierre Jolivet   if (bs[0] != bs[1]) bs[0] = 1;
114c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < is_max; ++i) {
115c7a4214aSPierre Jolivet     /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */
116c7a4214aSPierre Jolivet     /* needed to avoid subdomain matrices to replicate A since it is dense                           */
1179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize));
1185f80ce2aSJacob Faibussowitsch     PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported parallel IS");
1199566063dSJacob Faibussowitsch     PetscCall(ISGetSize(is[i], &size));
1209566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[i], &idx));
121c7a4214aSPierre Jolivet     for (PetscInt j = 0; j < size; ++j) {
122c7a4214aSPierre Jolivet       set.insert(idx[j]);
123c7a4214aSPierre Jolivet       for (PetscInt k = 1; k <= ov; ++k) {                                              /* for each layer of overlap      */
124c7a4214aSPierre Jolivet         if (idx[j] - k >= 0) set.insert(idx[j] - k);                                    /* do not insert negative indices */
125c7a4214aSPierre 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 */
126c7a4214aSPierre Jolivet       }
127c7a4214aSPierre Jolivet     }
1289566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[i], &idx));
1299566063dSJacob Faibussowitsch     PetscCall(ISDestroy(is + i));
1308f308287SPierre Jolivet     if (bs[0] > 1) {
1318f308287SPierre Jolivet       for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) {
1328f308287SPierre Jolivet         std::vector<PetscInt> block(bs[0]);
1338f308287SPierre Jolivet         std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]);
1348f308287SPierre Jolivet         set.insert(block.cbegin(), block.cend());
1358f308287SPierre Jolivet       }
1368f308287SPierre Jolivet     }
137c7a4214aSPierre Jolivet     size = set.size(); /* size with overlap */
1389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &oidx));
139c7a4214aSPierre Jolivet     for (const PetscInt j : set) *oidx++ = j;
140c7a4214aSPierre Jolivet     oidx -= size;
1419566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i));
142c7a4214aSPierre Jolivet   }
143c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
144c7a4214aSPierre Jolivet }
145c7a4214aSPierre Jolivet 
1469371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) {
147c7a4214aSPierre Jolivet   Mat_Htool         *a = (Mat_Htool *)A->data;
148c7a4214aSPierre Jolivet   Mat                D, B, BT;
149c7a4214aSPierre Jolivet   const PetscScalar *copy;
150c7a4214aSPierre Jolivet   PetscScalar       *ptr;
151c7a4214aSPierre Jolivet   const PetscInt    *idxr, *idxc, *it;
152c7a4214aSPierre Jolivet   PetscInt           nrow, m, i;
153c7a4214aSPierre Jolivet   PetscBool          flg;
154c7a4214aSPierre Jolivet 
155c7a4214aSPierre Jolivet   PetscFunctionBegin;
156*48a46eb9SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat));
157c7a4214aSPierre Jolivet   for (i = 0; i < n; ++i) {
1589566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &nrow));
1599566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &m));
1609566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &idxr));
1619566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &idxc));
162*48a46eb9SPierre Jolivet     if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, NULL, (*submat) + i));
1639566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr));
164c7a4214aSPierre Jolivet     if (irow[i] == icol[i]) { /* same row and column IS? */
1659566063dSJacob Faibussowitsch       PetscCall(MatHasCongruentLayouts(A, &flg));
166c7a4214aSPierre Jolivet       if (flg) {
1679566063dSJacob Faibussowitsch         PetscCall(ISSorted(irow[i], &flg));
168c7a4214aSPierre Jolivet         if (flg) { /* sorted IS? */
169c7a4214aSPierre Jolivet           it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart);
170c7a4214aSPierre Jolivet           if (it != idxr + nrow && *it == A->rmap->rstart) {    /* rmap->rstart in IS? */
171c7a4214aSPierre Jolivet             if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */
1729371c9d4SSatish Balay               for (PetscInt j = 0; j < A->rmap->n && flg; ++j)
1739371c9d4SSatish Balay                 if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE;
174c7a4214aSPierre Jolivet               if (flg) { /* complete local diagonal block in IS? */
175c7a4214aSPierre Jolivet                 /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM
176c7a4214aSPierre Jolivet                  *      [   B   C   E   ]
177c7a4214aSPierre Jolivet                  *  A = [   B   D   E   ]
178c7a4214aSPierre Jolivet                  *      [   B   F   E   ]
179c7a4214aSPierre Jolivet                  */
180c7a4214aSPierre Jolivet                 m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */
1819566063dSJacob Faibussowitsch                 PetscCall(MatGetDiagonalBlock_Htool(A, &D));
1829566063dSJacob Faibussowitsch                 PetscCall(MatDenseGetArrayRead(D, &copy));
1839371c9d4SSatish 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 */ }
1849566063dSJacob Faibussowitsch                 PetscCall(MatDenseRestoreArrayRead(D, &copy));
185c7a4214aSPierre Jolivet                 if (m) {
186c7a4214aSPierre Jolivet                   a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */
187c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
188b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
1899566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B));
1909566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
1919566063dSJacob Faibussowitsch                     PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT));
1929566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
193b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
1949566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
195c7a4214aSPierre Jolivet                     } else {
1967fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
1979566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
198c7a4214aSPierre Jolivet                     }
1999566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2009566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
201c7a4214aSPierre Jolivet                   } else {
202c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */
203c7a4214aSPierre Jolivet                       a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow);
204c7a4214aSPierre Jolivet                     }
205c7a4214aSPierre Jolivet                   }
206c7a4214aSPierre Jolivet                 }
207c7a4214aSPierre Jolivet                 if (m + A->rmap->n != nrow) {
208c7a4214aSPierre 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 */
209c7a4214aSPierre Jolivet                   /* entry-wise assembly may be costly, so transpose already-computed entries when possible */
210b94d7dedSBarry Smith                   if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) {
2119566063dSJacob 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));
2129566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(B, nrow));
2139566063dSJacob 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));
2149566063dSJacob Faibussowitsch                     PetscCall(MatDenseSetLDA(BT, nrow));
215b94d7dedSBarry Smith                     if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) {
2169566063dSJacob Faibussowitsch                       PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT));
217c7a4214aSPierre Jolivet                     } else {
2187fb60732SBarry Smith                       PetscCall(MatTransposeSetPrecursor(B, BT));
2199566063dSJacob Faibussowitsch                       PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT));
220c7a4214aSPierre Jolivet                     }
2219566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&B));
2229566063dSJacob Faibussowitsch                     PetscCall(MatDestroy(&BT));
223c7a4214aSPierre Jolivet                   } else {
224c7a4214aSPierre Jolivet                     for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */
225c7a4214aSPierre 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);
226c7a4214aSPierre Jolivet                     }
227c7a4214aSPierre Jolivet                   }
228c7a4214aSPierre Jolivet                 }
229c7a4214aSPierre Jolivet               }                       /* complete local diagonal block not in IS */
230c7a4214aSPierre Jolivet             } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */
231c7a4214aSPierre Jolivet           } else flg = PETSC_FALSE;   /* rmap->rstart not in IS */
232c7a4214aSPierre Jolivet         }                             /* unsorted IS */
233c7a4214aSPierre Jolivet       }
234c7a4214aSPierre Jolivet     } else flg = PETSC_FALSE;                                       /* different row and column IS */
235c7a4214aSPierre Jolivet     if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */
2369566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &idxr));
2379566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &idxc));
2389566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr));
2399566063dSJacob Faibussowitsch     PetscCall(MatScale((*submat)[i], a->s));
240c7a4214aSPierre Jolivet   }
241c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
242c7a4214aSPierre Jolivet }
243c7a4214aSPierre Jolivet 
2449371c9d4SSatish Balay static PetscErrorCode MatDestroy_Htool(Mat A) {
245c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool *)A->data;
246c7a4214aSPierre Jolivet   PetscContainer           container;
247c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
248c7a4214aSPierre Jolivet 
249c7a4214aSPierre Jolivet   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
2519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", NULL));
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", NULL));
2539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", NULL));
2549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", NULL));
2559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", NULL));
2569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", NULL));
2579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", NULL));
2589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", NULL));
2599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", NULL));
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container));
261c7a4214aSPierre Jolivet   if (container) { /* created in MatTranspose_Htool() */
2629566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
2639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
2649566063dSJacob Faibussowitsch     PetscCall(PetscFree(kernelt));
2659566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&container));
2669566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", NULL));
267c7a4214aSPierre Jolivet   }
268*48a46eb9SPierre Jolivet   if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->gcoords_target));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->work_source, a->work_target));
271c7a4214aSPierre Jolivet   delete a->wrapper;
272c7a4214aSPierre Jolivet   delete a->hmatrix;
2739566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
274c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
275c7a4214aSPierre Jolivet }
276c7a4214aSPierre Jolivet 
2779371c9d4SSatish Balay static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv) {
278c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
279c7a4214aSPierre Jolivet   PetscBool  flg;
280c7a4214aSPierre Jolivet 
281c7a4214aSPierre Jolivet   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg));
283c7a4214aSPierre Jolivet   if (flg) {
2849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->hmatrix->get_symmetry_type()));
285c7a4214aSPierre Jolivet     if (PetscAbsScalar(a->s - 1.0) > PETSC_MACHINE_EPSILON) {
286c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX)
2879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(a->s), (double)PetscImaginaryPart(a->s)));
288c7a4214aSPierre Jolivet #else
2899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)a->s));
290c7a4214aSPierre Jolivet #endif
291c7a4214aSPierre Jolivet     }
2929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum cluster size: %" PetscInt_FMT "\n", a->bs[0]));
2939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "maximum block size: %" PetscInt_FMT "\n", a->bs[1]));
2949566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon));
2959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta));
2969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0]));
2979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1]));
2989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor]));
2999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering]));
3009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", a->hmatrix->get_infos("Compression_ratio").c_str()));
3019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", a->hmatrix->get_infos("Space_saving").c_str()));
3029566063dSJacob 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()));
3039371c9d4SSatish 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(),
3049371c9d4SSatish Balay                                      a->hmatrix->get_infos("Dense_block_size_max").c_str()));
3059371c9d4SSatish 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(),
3069371c9d4SSatish Balay                                      a->hmatrix->get_infos("Low_rank_block_size_max").c_str()));
3079566063dSJacob 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()));
308c7a4214aSPierre Jolivet   }
309c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
310c7a4214aSPierre Jolivet }
311c7a4214aSPierre Jolivet 
3129371c9d4SSatish Balay static PetscErrorCode MatScale_Htool(Mat A, PetscScalar s) {
313c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
314c7a4214aSPierre Jolivet 
315c7a4214aSPierre Jolivet   PetscFunctionBegin;
316c7a4214aSPierre Jolivet   a->s *= s;
317c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
318c7a4214aSPierre Jolivet }
319c7a4214aSPierre Jolivet 
320c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */
3219371c9d4SSatish Balay static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) {
322c7a4214aSPierre Jolivet   Mat_Htool   *a = (Mat_Htool *)A->data;
323c7a4214aSPierre Jolivet   PetscInt    *idxc;
324c7a4214aSPierre Jolivet   PetscBLASInt one = 1, bn;
325c7a4214aSPierre Jolivet 
326c7a4214aSPierre Jolivet   PetscFunctionBegin;
327c7a4214aSPierre Jolivet   if (nz) *nz = A->cmap->N;
328c7a4214aSPierre Jolivet   if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */
3299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, &idxc));
330c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i;
331c7a4214aSPierre Jolivet   }
332c7a4214aSPierre Jolivet   if (idx) *idx = idxc;
333c7a4214aSPierre Jolivet   if (v) {
3349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(A->cmap->N, v));
335c7a4214aSPierre Jolivet     if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
336cab00e0dSPierre Jolivet     else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v);
3379566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(A->cmap->N, &bn));
338792fecdfSBarry Smith     PetscCallBLAS("BLASscal", BLASscal_(&bn, &a->s, *v, &one));
339c7a4214aSPierre Jolivet   }
340*48a46eb9SPierre Jolivet   if (!idx) PetscCall(PetscFree(idxc));
341c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
342c7a4214aSPierre Jolivet }
343c7a4214aSPierre Jolivet 
3449371c9d4SSatish Balay static PetscErrorCode MatRestoreRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) {
345c7a4214aSPierre Jolivet   PetscFunctionBegin;
346c7a4214aSPierre Jolivet   if (nz) *nz = 0;
347*48a46eb9SPierre Jolivet   if (idx) PetscCall(PetscFree(*idx));
348*48a46eb9SPierre Jolivet   if (v) PetscCall(PetscFree(*v));
349c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
350c7a4214aSPierre Jolivet }
351c7a4214aSPierre Jolivet 
3529371c9d4SSatish Balay static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject) {
353c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
354c7a4214aSPierre Jolivet   PetscInt   n;
355c7a4214aSPierre Jolivet   PetscBool  flg;
356c7a4214aSPierre Jolivet 
357c7a4214aSPierre Jolivet   PetscFunctionBegin;
358d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Htool options");
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", NULL, a->bs[0], a->bs, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_max_block_size", "Maximal number of coefficients in a dense block", NULL, a->bs[1], a->bs + 1, NULL));
3619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", NULL, a->epsilon, &a->epsilon, NULL));
3629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", NULL, a->eta, &a->eta, NULL));
3639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", NULL, a->depth[0], a->depth, NULL));
3649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", NULL, a->depth[1], a->depth + 1, NULL));
365c7a4214aSPierre Jolivet   n = 0;
366dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg));
367c7a4214aSPierre Jolivet   if (flg) a->compressor = MatHtoolCompressorType(n);
36898e73e17SPierre Jolivet   n = 0;
369dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg));
37098e73e17SPierre Jolivet   if (flg) a->clustering = MatHtoolClusteringType(n);
371d0609cedSBarry Smith   PetscOptionsHeadEnd();
372c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
373c7a4214aSPierre Jolivet }
374c7a4214aSPierre Jolivet 
3759371c9d4SSatish Balay static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType type) {
376c7a4214aSPierre Jolivet   Mat_Htool                                                   *a = (Mat_Htool *)A->data;
377c7a4214aSPierre Jolivet   const PetscInt                                              *ranges;
378c7a4214aSPierre Jolivet   PetscInt                                                    *offset;
379c7a4214aSPierre Jolivet   PetscMPIInt                                                  size;
380b94d7dedSBarry Smith   char                                                         S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U';
381cab00e0dSPierre Jolivet   htool::VirtualGenerator<PetscScalar>                        *generator = nullptr;
38298e73e17SPierre Jolivet   std::shared_ptr<htool::VirtualCluster>                       t, s = nullptr;
3833b9338faSPierre Jolivet   std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr;
384c7a4214aSPierre Jolivet 
385c7a4214aSPierre Jolivet   PetscFunctionBegin;
3869566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite));
387c7a4214aSPierre Jolivet   delete a->wrapper;
388c7a4214aSPierre Jolivet   delete a->hmatrix;
3899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
3909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &offset));
3919566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
392c7a4214aSPierre Jolivet   for (PetscInt i = 0; i < size; ++i) {
393c7a4214aSPierre Jolivet     offset[2 * i]     = ranges[i];
394c7a4214aSPierre Jolivet     offset[2 * i + 1] = ranges[i + 1] - ranges[i];
395c7a4214aSPierre Jolivet   }
39698e73e17SPierre Jolivet   switch (a->clustering) {
3979371c9d4SSatish Balay   case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); break;
3989371c9d4SSatish Balay   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); break;
3999371c9d4SSatish Balay   case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); break;
4009371c9d4SSatish Balay   default: t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
40198e73e17SPierre Jolivet   }
402c7a4214aSPierre Jolivet   t->set_minclustersize(a->bs[0]);
40398e73e17SPierre Jolivet   t->build(A->rmap->N, a->gcoords_target, offset);
404c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx);
405c7a4214aSPierre Jolivet   else {
406c7a4214aSPierre Jolivet     a->wrapper = NULL;
407cab00e0dSPierre Jolivet     generator  = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx);
408c7a4214aSPierre Jolivet   }
409c7a4214aSPierre Jolivet   if (a->gcoords_target != a->gcoords_source) {
4109566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRangesColumn(A, &ranges));
411c7a4214aSPierre Jolivet     for (PetscInt i = 0; i < size; ++i) {
412c7a4214aSPierre Jolivet       offset[2 * i]     = ranges[i];
413c7a4214aSPierre Jolivet       offset[2 * i + 1] = ranges[i + 1] - ranges[i];
414c7a4214aSPierre Jolivet     }
41598e73e17SPierre Jolivet     switch (a->clustering) {
4169371c9d4SSatish Balay     case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); break;
4179371c9d4SSatish Balay     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); break;
4189371c9d4SSatish Balay     case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); break;
4199371c9d4SSatish Balay     default: s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim);
42098e73e17SPierre Jolivet     }
421c7a4214aSPierre Jolivet     s->set_minclustersize(a->bs[0]);
42298e73e17SPierre Jolivet     s->build(A->cmap->N, a->gcoords_source, offset);
423c7a4214aSPierre Jolivet     S = uplo = 'N';
424c7a4214aSPierre Jolivet   }
4259566063dSJacob Faibussowitsch   PetscCall(PetscFree(offset));
426c7a4214aSPierre Jolivet   switch (a->compressor) {
4279371c9d4SSatish Balay   case MAT_HTOOL_COMPRESSOR_FULL_ACA: compressor = std::make_shared<htool::fullACA<PetscScalar>>(); break;
4289371c9d4SSatish Balay   case MAT_HTOOL_COMPRESSOR_SVD: compressor = std::make_shared<htool::SVD<PetscScalar>>(); break;
4299371c9d4SSatish Balay   default: compressor = std::make_shared<htool::sympartialACA<PetscScalar>>();
430c7a4214aSPierre Jolivet   }
4313b9338faSPierre Jolivet   a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar> *>(new htool::HMatrix<PetscScalar>(t, s ? s : t, a->epsilon, a->eta, S, uplo));
4323b9338faSPierre Jolivet   a->hmatrix->set_compression(compressor);
433c7a4214aSPierre Jolivet   a->hmatrix->set_maxblocksize(a->bs[1]);
434c7a4214aSPierre Jolivet   a->hmatrix->set_mintargetdepth(a->depth[0]);
435c7a4214aSPierre Jolivet   a->hmatrix->set_minsourcedepth(a->depth[1]);
4363b9338faSPierre Jolivet   if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target, a->gcoords_source);
4373b9338faSPierre Jolivet   else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target);
438c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
439c7a4214aSPierre Jolivet }
440c7a4214aSPierre Jolivet 
4419371c9d4SSatish Balay static PetscErrorCode MatProductNumeric_Htool(Mat C) {
442c7a4214aSPierre Jolivet   Mat_Product       *product = C->product;
443c7a4214aSPierre Jolivet   Mat_Htool         *a       = (Mat_Htool *)product->A->data;
444c7a4214aSPierre Jolivet   const PetscScalar *in;
445c7a4214aSPierre Jolivet   PetscScalar       *out;
4468b8ddd11SPierre Jolivet   PetscInt           N, lda;
447c7a4214aSPierre Jolivet 
448c7a4214aSPierre Jolivet   PetscFunctionBegin;
449c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
4509566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, NULL, &N));
4519566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &lda));
4525f80ce2aSJacob Faibussowitsch   PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
4539566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(product->B, &in));
4549566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &out));
4558b8ddd11SPierre Jolivet   switch (product->type) {
4569371c9d4SSatish Balay   case MATPRODUCT_AB: a->hmatrix->mvprod_local_to_local(in, out, N); break;
4579371c9d4SSatish Balay   case MATPRODUCT_AtB: a->hmatrix->mvprod_transp_local_to_local(in, out, N); break;
4589371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]);
459c7a4214aSPierre Jolivet   }
4609566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &out));
4619566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(product->B, &in));
4629566063dSJacob Faibussowitsch   PetscCall(MatScale(C, a->s));
463c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
464c7a4214aSPierre Jolivet }
465c7a4214aSPierre Jolivet 
4669371c9d4SSatish Balay static PetscErrorCode MatProductSymbolic_Htool(Mat C) {
467c7a4214aSPierre Jolivet   Mat_Product *product = C->product;
468c7a4214aSPierre Jolivet   Mat          A, B;
469c7a4214aSPierre Jolivet   PetscBool    flg;
470c7a4214aSPierre Jolivet 
471c7a4214aSPierre Jolivet   PetscFunctionBegin;
472c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
473c7a4214aSPierre Jolivet   A = product->A;
474c7a4214aSPierre Jolivet   B = product->B;
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, ""));
4765f80ce2aSJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MatProduct_AB not supported for %s", ((PetscObject)product->B)->type_name);
477c7a4214aSPierre Jolivet   switch (product->type) {
478c7a4214aSPierre Jolivet   case MATPRODUCT_AB:
479*48a46eb9SPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
4808b8ddd11SPierre Jolivet     break;
4818b8ddd11SPierre Jolivet   case MATPRODUCT_AtB:
482*48a46eb9SPierre Jolivet     if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
4838b8ddd11SPierre Jolivet     break;
4849371c9d4SSatish Balay   default: SETERRQ(PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "ProductType %s is not supported", MatProductTypes[product->type]);
4858b8ddd11SPierre Jolivet   }
4869566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
4879566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
4889566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
4899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
4909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
491c7a4214aSPierre Jolivet   C->ops->productsymbolic = NULL;
492c7a4214aSPierre Jolivet   C->ops->productnumeric  = MatProductNumeric_Htool;
493c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
494c7a4214aSPierre Jolivet }
495c7a4214aSPierre Jolivet 
4969371c9d4SSatish Balay static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) {
497c7a4214aSPierre Jolivet   PetscFunctionBegin;
498c7a4214aSPierre Jolivet   MatCheckProduct(C, 1);
4998b8ddd11SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool;
500c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
501c7a4214aSPierre Jolivet }
502c7a4214aSPierre Jolivet 
5039371c9d4SSatish Balay static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) {
504c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
505c7a4214aSPierre Jolivet 
506c7a4214aSPierre Jolivet   PetscFunctionBegin;
507c7a4214aSPierre Jolivet   *hmatrix = a->hmatrix;
508c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
509c7a4214aSPierre Jolivet }
510c7a4214aSPierre Jolivet 
511c7a4214aSPierre Jolivet /*@C
512c7a4214aSPierre Jolivet      MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a MATHTOOL.
513c7a4214aSPierre Jolivet 
514c7a4214aSPierre Jolivet    Input Parameter:
515c7a4214aSPierre Jolivet .     A - hierarchical matrix
516c7a4214aSPierre Jolivet 
517c7a4214aSPierre Jolivet    Output Parameter:
518c7a4214aSPierre Jolivet .     hmatrix - opaque pointer to a Htool virtual matrix
519c7a4214aSPierre Jolivet 
520c7a4214aSPierre Jolivet    Level: advanced
521c7a4214aSPierre Jolivet 
522db781477SPatrick Sanan .seealso: `MATHTOOL`
523c7a4214aSPierre Jolivet @*/
5249371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) {
525c7a4214aSPierre Jolivet   PetscFunctionBegin;
526c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
527c7a4214aSPierre Jolivet   PetscValidPointer(hmatrix, 2);
528cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::VirtualHMatrix<PetscScalar> **), (A, hmatrix));
529c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
530c7a4214aSPierre Jolivet }
531c7a4214aSPierre Jolivet 
5329371c9d4SSatish Balay static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernel kernel, void *kernelctx) {
533c7a4214aSPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
534c7a4214aSPierre Jolivet 
535c7a4214aSPierre Jolivet   PetscFunctionBegin;
536c7a4214aSPierre Jolivet   a->kernel    = kernel;
537c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
538c7a4214aSPierre Jolivet   delete a->wrapper;
539c7a4214aSPierre Jolivet   if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx);
540c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
541c7a4214aSPierre Jolivet }
542c7a4214aSPierre Jolivet 
543c7a4214aSPierre Jolivet /*@C
54498e73e17SPierre Jolivet      MatHtoolSetKernel - Sets the kernel and context used for the assembly of a MATHTOOL.
545c7a4214aSPierre Jolivet 
546c7a4214aSPierre Jolivet    Input Parameters:
547c7a4214aSPierre Jolivet +     A - hierarchical matrix
548c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
549cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
550c7a4214aSPierre Jolivet 
551c7a4214aSPierre Jolivet    Level: advanced
552c7a4214aSPierre Jolivet 
553db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatCreateHtoolFromKernel()`
554c7a4214aSPierre Jolivet @*/
5559371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernel kernel, void *kernelctx) {
556c7a4214aSPierre Jolivet   PetscFunctionBegin;
557c7a4214aSPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
558c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 2);
559c7a4214aSPierre Jolivet   if (!kernel) PetscValidPointer(kernelctx, 3);
560cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernel, void *), (A, kernel, kernelctx));
561c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
562c7a4214aSPierre Jolivet }
563c7a4214aSPierre Jolivet 
5649371c9d4SSatish Balay static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) {
56598e73e17SPierre Jolivet   Mat_Htool            *a = (Mat_Htool *)A->data;
56698e73e17SPierre Jolivet   std::vector<PetscInt> source;
56798e73e17SPierre Jolivet 
56898e73e17SPierre Jolivet   PetscFunctionBegin;
569a9918087SPierre Jolivet   source = a->hmatrix->get_source_cluster()->get_local_perm();
5709566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), source.size(), source.data(), PETSC_COPY_VALUES, is));
5719566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
57298e73e17SPierre Jolivet   PetscFunctionReturn(0);
57398e73e17SPierre Jolivet }
57498e73e17SPierre Jolivet 
57598e73e17SPierre Jolivet /*@C
57698e73e17SPierre Jolivet      MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster.
57798e73e17SPierre Jolivet 
57898e73e17SPierre Jolivet    Input Parameter:
57998e73e17SPierre Jolivet .     A - hierarchical matrix
58098e73e17SPierre Jolivet 
58198e73e17SPierre Jolivet    Output Parameter:
58298e73e17SPierre Jolivet .     is - permutation
58398e73e17SPierre Jolivet 
58498e73e17SPierre Jolivet    Level: advanced
58598e73e17SPierre Jolivet 
586db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()`
58798e73e17SPierre Jolivet @*/
5889371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) {
58998e73e17SPierre Jolivet   PetscFunctionBegin;
59098e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
59198e73e17SPierre Jolivet   if (!is) PetscValidPointer(is, 2);
592cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is));
59398e73e17SPierre Jolivet   PetscFunctionReturn(0);
59498e73e17SPierre Jolivet }
59598e73e17SPierre Jolivet 
5969371c9d4SSatish Balay static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) {
59798e73e17SPierre Jolivet   Mat_Htool            *a = (Mat_Htool *)A->data;
59898e73e17SPierre Jolivet   std::vector<PetscInt> target;
59998e73e17SPierre Jolivet 
60098e73e17SPierre Jolivet   PetscFunctionBegin;
601a9918087SPierre Jolivet   target = a->hmatrix->get_target_cluster()->get_local_perm();
6029566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), target.size(), target.data(), PETSC_COPY_VALUES, is));
6039566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*is));
60498e73e17SPierre Jolivet   PetscFunctionReturn(0);
60598e73e17SPierre Jolivet }
60698e73e17SPierre Jolivet 
60798e73e17SPierre Jolivet /*@C
60898e73e17SPierre Jolivet      MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster.
60998e73e17SPierre Jolivet 
61098e73e17SPierre Jolivet    Input Parameter:
61198e73e17SPierre Jolivet .     A - hierarchical matrix
61298e73e17SPierre Jolivet 
61398e73e17SPierre Jolivet    Output Parameter:
61498e73e17SPierre Jolivet .     is - permutation
61598e73e17SPierre Jolivet 
61698e73e17SPierre Jolivet    Level: advanced
61798e73e17SPierre Jolivet 
618db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()`
61998e73e17SPierre Jolivet @*/
6209371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) {
62198e73e17SPierre Jolivet   PetscFunctionBegin;
62298e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
62398e73e17SPierre Jolivet   if (!is) PetscValidPointer(is, 2);
624cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is));
62598e73e17SPierre Jolivet   PetscFunctionReturn(0);
62698e73e17SPierre Jolivet }
62798e73e17SPierre Jolivet 
6289371c9d4SSatish Balay static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) {
62998e73e17SPierre Jolivet   Mat_Htool *a = (Mat_Htool *)A->data;
63098e73e17SPierre Jolivet 
63198e73e17SPierre Jolivet   PetscFunctionBegin;
63298e73e17SPierre Jolivet   a->hmatrix->set_use_permutation(use);
63398e73e17SPierre Jolivet   PetscFunctionReturn(0);
63498e73e17SPierre Jolivet }
63598e73e17SPierre Jolivet 
63698e73e17SPierre Jolivet /*@C
63798e73e17SPierre Jolivet      MatHtoolUsePermutation - Sets whether MATHTOOL should permute input (resp. output) vectors following its internal source (resp. target) permutation.
63898e73e17SPierre Jolivet 
63998e73e17SPierre Jolivet    Input Parameters:
64098e73e17SPierre Jolivet +     A - hierarchical matrix
64198e73e17SPierre Jolivet -     use - Boolean value
64298e73e17SPierre Jolivet 
64398e73e17SPierre Jolivet    Level: advanced
64498e73e17SPierre Jolivet 
645db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()`
64698e73e17SPierre Jolivet @*/
6479371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) {
64898e73e17SPierre Jolivet   PetscFunctionBegin;
64998e73e17SPierre Jolivet   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
65098e73e17SPierre Jolivet   PetscValidLogicalCollectiveBool(A, use, 2);
651cac4c232SBarry Smith   PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use));
65298e73e17SPierre Jolivet   PetscFunctionReturn(0);
65398e73e17SPierre Jolivet }
65498e73e17SPierre Jolivet 
6559371c9d4SSatish Balay static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType newtype, MatReuse reuse, Mat *B) {
656c7a4214aSPierre Jolivet   Mat          C;
657c7a4214aSPierre Jolivet   Mat_Htool   *a = (Mat_Htool *)A->data;
658c7a4214aSPierre Jolivet   PetscInt     lda;
659c7a4214aSPierre Jolivet   PetscScalar *array;
660c7a4214aSPierre Jolivet 
661c7a4214aSPierre Jolivet   PetscFunctionBegin;
662c7a4214aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
663c7a4214aSPierre Jolivet     C = *B;
6645f80ce2aSJacob Faibussowitsch     PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions");
6659566063dSJacob Faibussowitsch     PetscCall(MatDenseGetLDA(C, &lda));
6665f80ce2aSJacob Faibussowitsch     PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n);
667c7a4214aSPierre Jolivet   } else {
6689566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
6699566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
6709566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, MATDENSE));
6719566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
6729566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
673c7a4214aSPierre Jolivet   }
6749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
6759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
6769566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &array));
677c7a4214aSPierre Jolivet   a->hmatrix->copy_local_dense_perm(array);
6789566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(C, &array));
6799566063dSJacob Faibussowitsch   PetscCall(MatScale(C, a->s));
680c7a4214aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
6819566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &C));
682c7a4214aSPierre Jolivet   } else *B = C;
683c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
684c7a4214aSPierre Jolivet }
685c7a4214aSPierre Jolivet 
6869371c9d4SSatish Balay static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx) {
687c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx;
68898e73e17SPierre Jolivet   PetscScalar             *tmp;
68998e73e17SPierre Jolivet 
69098e73e17SPierre Jolivet   PetscFunctionBegin;
69198e73e17SPierre Jolivet   generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx);
6929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M * N, &tmp));
6939566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tmp, ptr, M * N));
69498e73e17SPierre Jolivet   for (PetscInt i = 0; i < M; ++i) {
69598e73e17SPierre Jolivet     for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N];
69698e73e17SPierre Jolivet   }
6979566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmp));
69898e73e17SPierre Jolivet   PetscFunctionReturn(0);
699c7a4214aSPierre Jolivet }
700c7a4214aSPierre Jolivet 
701c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */
7029371c9d4SSatish Balay static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) {
703c7a4214aSPierre Jolivet   Mat                      C;
704c7a4214aSPierre Jolivet   Mat_Htool               *a = (Mat_Htool *)A->data, *c;
705c7a4214aSPierre Jolivet   PetscInt                 M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n;
706c7a4214aSPierre Jolivet   PetscContainer           container;
707c7a4214aSPierre Jolivet   MatHtoolKernelTranspose *kernelt;
708c7a4214aSPierre Jolivet 
709c7a4214aSPierre Jolivet   PetscFunctionBegin;
7107fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
7115f80ce2aSJacob Faibussowitsch   PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported");
712c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
7139566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
7149566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, n, m, N, M));
7159566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, ((PetscObject)A)->type_name));
7169566063dSJacob Faibussowitsch     PetscCall(MatSetUp(C));
7179566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C), &container));
7189566063dSJacob Faibussowitsch     PetscCall(PetscNew(&kernelt));
7199566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(container, kernelt));
7209566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)C, "KernelTranspose", (PetscObject)container));
721c7a4214aSPierre Jolivet   } else {
722c7a4214aSPierre Jolivet     C = *B;
7239566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container));
7245f80ce2aSJacob Faibussowitsch     PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first");
7259566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(container, (void **)&kernelt));
726c7a4214aSPierre Jolivet   }
727c7a4214aSPierre Jolivet   c         = (Mat_Htool *)C->data;
728c7a4214aSPierre Jolivet   c->dim    = a->dim;
729c7a4214aSPierre Jolivet   c->s      = a->s;
73098e73e17SPierre Jolivet   c->kernel = GenEntriesTranspose;
731c7a4214aSPierre Jolivet   if (kernelt->A != A) {
7329566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&kernelt->A));
733c7a4214aSPierre Jolivet     kernelt->A = A;
7349566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
735c7a4214aSPierre Jolivet   }
736c7a4214aSPierre Jolivet   kernelt->kernel    = a->kernel;
737c7a4214aSPierre Jolivet   kernelt->kernelctx = a->kernelctx;
738c7a4214aSPierre Jolivet   c->kernelctx       = kernelt;
739c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) {
7409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target));
7419566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim));
742c7a4214aSPierre Jolivet     if (a->gcoords_target != a->gcoords_source) {
7439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source));
7449566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim));
745c7a4214aSPierre Jolivet     } else c->gcoords_source = c->gcoords_target;
7469566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target));
747c7a4214aSPierre Jolivet   }
7489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
7499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
750c7a4214aSPierre Jolivet   if (reuse == MAT_INITIAL_MATRIX) *B = C;
751c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
752c7a4214aSPierre Jolivet }
753c7a4214aSPierre Jolivet 
754c7a4214aSPierre Jolivet /*@C
755c7a4214aSPierre Jolivet      MatCreateHtoolFromKernel - Creates a MATHTOOL from a user-supplied kernel.
756c7a4214aSPierre Jolivet 
757c7a4214aSPierre Jolivet    Input Parameters:
758c7a4214aSPierre Jolivet +     comm - MPI communicator
759c7a4214aSPierre Jolivet .     m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
760c7a4214aSPierre Jolivet .     n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
761c7a4214aSPierre Jolivet .     M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
762c7a4214aSPierre Jolivet .     N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
763c7a4214aSPierre Jolivet .     spacedim - dimension of the space coordinates
764c7a4214aSPierre Jolivet .     coords_target - coordinates of the target
765c7a4214aSPierre Jolivet .     coords_source - coordinates of the source
766c7a4214aSPierre Jolivet .     kernel - computational kernel (or NULL)
767cab00e0dSPierre Jolivet -     kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*)
768c7a4214aSPierre Jolivet 
769c7a4214aSPierre Jolivet    Output Parameter:
770c7a4214aSPierre Jolivet .     B - matrix
771c7a4214aSPierre Jolivet 
772c7a4214aSPierre Jolivet    Options Database Keys:
773c7a4214aSPierre Jolivet +     -mat_htool_min_cluster_size <PetscInt> - minimal leaf size in cluster tree
774c7a4214aSPierre Jolivet .     -mat_htool_max_block_size <PetscInt> - maximal number of coefficients in a dense block
775c7a4214aSPierre Jolivet .     -mat_htool_epsilon <PetscReal> - relative error in Frobenius norm when approximating a block
776c7a4214aSPierre Jolivet .     -mat_htool_eta <PetscReal> - admissibility condition tolerance
777c7a4214aSPierre Jolivet .     -mat_htool_min_target_depth <PetscInt> - minimal cluster tree depth associated with the rows
778c7a4214aSPierre Jolivet .     -mat_htool_min_source_depth <PetscInt> - minimal cluster tree depth associated with the columns
77998e73e17SPierre Jolivet .     -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression
78098e73e17SPierre Jolivet -     -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering
781c7a4214aSPierre Jolivet 
782c7a4214aSPierre Jolivet    Level: intermediate
783c7a4214aSPierre Jolivet 
784db781477SPatrick Sanan .seealso: `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
785c7a4214aSPierre Jolivet @*/
7869371c9d4SSatish Balay PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords_target[], const PetscReal coords_source[], MatHtoolKernel kernel, void *kernelctx, Mat *B) {
787c7a4214aSPierre Jolivet   Mat        A;
788c7a4214aSPierre Jolivet   Mat_Htool *a;
789c7a4214aSPierre Jolivet 
790c7a4214aSPierre Jolivet   PetscFunctionBegin;
7919566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
792c7a4214aSPierre Jolivet   PetscValidLogicalCollectiveInt(A, spacedim, 6);
793c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_target, 7);
794c7a4214aSPierre Jolivet   PetscValidRealPointer(coords_source, 8);
795c7a4214aSPierre Jolivet   if (!kernelctx) PetscValidFunction(kernel, 9);
796c7a4214aSPierre Jolivet   if (!kernel) PetscValidPointer(kernelctx, 10);
7979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m, n, M, N));
7989566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATHTOOL));
7999566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
800c7a4214aSPierre Jolivet   a            = (Mat_Htool *)A->data;
801c7a4214aSPierre Jolivet   a->dim       = spacedim;
802c7a4214aSPierre Jolivet   a->s         = 1.0;
803c7a4214aSPierre Jolivet   a->kernel    = kernel;
804c7a4214aSPierre Jolivet   a->kernelctx = kernelctx;
8059566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target));
8069566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim));
8071c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */
808c7a4214aSPierre Jolivet   if (coords_target != coords_source) {
8099566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source));
8109566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim));
8111c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */
812c7a4214aSPierre Jolivet   } else a->gcoords_source = a->gcoords_target;
8139566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target));
814c7a4214aSPierre Jolivet   *B = A;
815c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
816c7a4214aSPierre Jolivet }
817c7a4214aSPierre Jolivet 
818c7a4214aSPierre Jolivet /*MC
819c7a4214aSPierre Jolivet      MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package.
820c7a4214aSPierre Jolivet 
821c7a4214aSPierre Jolivet   Use ./configure --download-htool to install PETSc to use Htool.
822c7a4214aSPierre Jolivet 
823c7a4214aSPierre Jolivet    Options Database Keys:
824e85a1b63SPatrick Sanan .     -mat_type htool - matrix type to "htool" during a call to MatSetFromOptions()
825c7a4214aSPierre Jolivet 
826c7a4214aSPierre Jolivet    Level: beginner
827c7a4214aSPierre Jolivet 
828db781477SPatrick Sanan .seealso: `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()`
829c7a4214aSPierre Jolivet M*/
8309371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) {
831c7a4214aSPierre Jolivet   Mat_Htool *a;
832c7a4214aSPierre Jolivet 
833c7a4214aSPierre Jolivet   PetscFunctionBegin;
8349566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(A, &a));
835c7a4214aSPierre Jolivet   A->data = (void *)a;
8369566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL));
8379566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
838c7a4214aSPierre Jolivet   A->ops->getdiagonal      = MatGetDiagonal_Htool;
839c7a4214aSPierre Jolivet   A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool;
840c7a4214aSPierre Jolivet   A->ops->mult             = MatMult_Htool;
841c7a4214aSPierre Jolivet   A->ops->multadd          = MatMultAdd_Htool;
8428b8ddd11SPierre Jolivet   A->ops->multtranspose    = MatMultTranspose_Htool;
8438b8ddd11SPierre Jolivet   if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool;
844c7a4214aSPierre Jolivet   A->ops->increaseoverlap   = MatIncreaseOverlap_Htool;
845c7a4214aSPierre Jolivet   A->ops->createsubmatrices = MatCreateSubMatrices_Htool;
846c7a4214aSPierre Jolivet   A->ops->transpose         = MatTranspose_Htool;
847c7a4214aSPierre Jolivet   A->ops->destroy           = MatDestroy_Htool;
848c7a4214aSPierre Jolivet   A->ops->view              = MatView_Htool;
849c7a4214aSPierre Jolivet   A->ops->setfromoptions    = MatSetFromOptions_Htool;
850c7a4214aSPierre Jolivet   A->ops->scale             = MatScale_Htool;
851c7a4214aSPierre Jolivet   A->ops->getrow            = MatGetRow_Htool;
852c7a4214aSPierre Jolivet   A->ops->restorerow        = MatRestoreRow_Htool;
853c7a4214aSPierre Jolivet   A->ops->assemblyend       = MatAssemblyEnd_Htool;
854c7a4214aSPierre Jolivet   a->dim                    = 0;
855c7a4214aSPierre Jolivet   a->gcoords_target         = NULL;
856c7a4214aSPierre Jolivet   a->gcoords_source         = NULL;
857c7a4214aSPierre Jolivet   a->s                      = 1.0;
858c7a4214aSPierre Jolivet   a->bs[0]                  = 10;
859c7a4214aSPierre Jolivet   a->bs[1]                  = 1000000;
860c7a4214aSPierre Jolivet   a->epsilon                = PetscSqrtReal(PETSC_SMALL);
861c7a4214aSPierre Jolivet   a->eta                    = 10.0;
862c7a4214aSPierre Jolivet   a->depth[0]               = 0;
863c7a4214aSPierre Jolivet   a->depth[1]               = 0;
864c7a4214aSPierre Jolivet   a->compressor             = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA;
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool));
8669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool));
8679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense));
8689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense));
8699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool));
8709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool));
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool));
8729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool));
8739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool));
874c7a4214aSPierre Jolivet   PetscFunctionReturn(0);
875c7a4214aSPierre Jolivet }
876