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 { 21c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 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"); 289566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(v, &x)); 29c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal(x); 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(v, &x)); 319566063dSJacob Faibussowitsch PetscCall(VecScale(v, a->s)); 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33c7a4214aSPierre Jolivet } 34c7a4214aSPierre Jolivet 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b) 36d71ae5a4SJacob Faibussowitsch { 37c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 38c7a4214aSPierre Jolivet Mat B; 39c7a4214aSPierre Jolivet PetscScalar *ptr; 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"); 459566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */ 46c7a4214aSPierre Jolivet if (!B) { 47f22e26b7SPierre Jolivet PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B)); 489566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(B, &ptr)); 49c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal_block(ptr); 509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(B, &ptr)); 519566063dSJacob Faibussowitsch PetscCall(MatPropagateSymmetryOptions(A, B)); 529566063dSJacob Faibussowitsch PetscCall(MatScale(B, a->s)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B)); 54c7a4214aSPierre Jolivet *b = B; 559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 56c7a4214aSPierre Jolivet } else *b = B; 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58c7a4214aSPierre Jolivet } 59c7a4214aSPierre Jolivet 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y) 61d71ae5a4SJacob Faibussowitsch { 62c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 63c7a4214aSPierre Jolivet const PetscScalar *in; 64c7a4214aSPierre Jolivet PetscScalar *out; 65c7a4214aSPierre Jolivet 66c7a4214aSPierre Jolivet PetscFunctionBegin; 679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 689566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 69c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in, out); 709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 729566063dSJacob Faibussowitsch PetscCall(VecScale(y, a->s)); 733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74c7a4214aSPierre Jolivet } 75c7a4214aSPierre Jolivet 76c7a4214aSPierre Jolivet /* naive implementation of MatMultAdd() needed for FEM-BEM coupling via MATNEST */ 77d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_Htool(Mat A, Vec v1, Vec v2, Vec v3) 78d71ae5a4SJacob Faibussowitsch { 79c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 80c7a4214aSPierre Jolivet Vec tmp; 81c7a4214aSPierre Jolivet const PetscScalar scale = a->s; 82c7a4214aSPierre Jolivet 83c7a4214aSPierre Jolivet PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(v2, &tmp)); 859566063dSJacob 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 */ 86c7a4214aSPierre Jolivet a->s = 1.0; /* set s to 1.0 since VecAXPY() may be used to scale the MatMult() output Vec */ 879566063dSJacob Faibussowitsch PetscCall(MatMult_Htool(A, v1, tmp)); 889566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3, scale, tmp)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp)); 90c7a4214aSPierre Jolivet a->s = scale; /* set s back to its original value */ 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92c7a4214aSPierre Jolivet } 93c7a4214aSPierre Jolivet 94d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y) 95d71ae5a4SJacob Faibussowitsch { 968b8ddd11SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 978b8ddd11SPierre Jolivet const PetscScalar *in; 988b8ddd11SPierre Jolivet PetscScalar *out; 998b8ddd11SPierre Jolivet 1008b8ddd11SPierre Jolivet PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 1038b8ddd11SPierre Jolivet a->hmatrix->mvprod_transp_local_to_local(in, out); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 1069566063dSJacob Faibussowitsch PetscCall(VecScale(y, a->s)); 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1088b8ddd11SPierre Jolivet } 1098b8ddd11SPierre Jolivet 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov) 111d71ae5a4SJacob Faibussowitsch { 112c7a4214aSPierre Jolivet std::set<PetscInt> set; 113c7a4214aSPierre Jolivet const PetscInt *idx; 1148f308287SPierre Jolivet PetscInt *oidx, size, bs[2]; 115c7a4214aSPierre Jolivet PetscMPIInt csize; 116c7a4214aSPierre Jolivet 117c7a4214aSPierre Jolivet PetscFunctionBegin; 1189566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, bs, bs + 1)); 1198f308287SPierre Jolivet if (bs[0] != bs[1]) bs[0] = 1; 120c7a4214aSPierre Jolivet for (PetscInt i = 0; i < is_max; ++i) { 121c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 122c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize)); 1245f80ce2aSJacob Faibussowitsch PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported parallel IS"); 1259566063dSJacob Faibussowitsch PetscCall(ISGetSize(is[i], &size)); 1269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx)); 127c7a4214aSPierre Jolivet for (PetscInt j = 0; j < size; ++j) { 128c7a4214aSPierre Jolivet set.insert(idx[j]); 129c7a4214aSPierre Jolivet for (PetscInt k = 1; k <= ov; ++k) { /* for each layer of overlap */ 130c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 131c7a4214aSPierre 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 */ 132c7a4214aSPierre Jolivet } 133c7a4214aSPierre Jolivet } 1349566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[i], &idx)); 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(is + i)); 1368f308287SPierre Jolivet if (bs[0] > 1) { 1378f308287SPierre Jolivet for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) { 1388f308287SPierre Jolivet std::vector<PetscInt> block(bs[0]); 1398f308287SPierre Jolivet std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]); 1408f308287SPierre Jolivet set.insert(block.cbegin(), block.cend()); 1418f308287SPierre Jolivet } 1428f308287SPierre Jolivet } 143c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &oidx)); 145c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 146c7a4214aSPierre Jolivet oidx -= size; 1479566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i)); 148c7a4214aSPierre Jolivet } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150c7a4214aSPierre Jolivet } 151c7a4214aSPierre Jolivet 152d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 153d71ae5a4SJacob Faibussowitsch { 154c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 155c7a4214aSPierre Jolivet Mat D, B, BT; 156c7a4214aSPierre Jolivet const PetscScalar *copy; 157c7a4214aSPierre Jolivet PetscScalar *ptr; 158c7a4214aSPierre Jolivet const PetscInt *idxr, *idxc, *it; 159c7a4214aSPierre Jolivet PetscInt nrow, m, i; 160c7a4214aSPierre Jolivet PetscBool flg; 161c7a4214aSPierre Jolivet 162c7a4214aSPierre Jolivet PetscFunctionBegin; 16348a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 164c7a4214aSPierre Jolivet for (i = 0; i < n; ++i) { 1659566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &nrow)); 1669566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &m)); 1679566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i], &idxr)); 1689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i], &idxc)); 169f22e26b7SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i)); 1709566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr)); 171c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 1729566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 173c7a4214aSPierre Jolivet if (flg) { 1749566063dSJacob Faibussowitsch PetscCall(ISSorted(irow[i], &flg)); 175c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 176c7a4214aSPierre Jolivet it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart); 177c7a4214aSPierre Jolivet if (it != idxr + nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 178c7a4214aSPierre Jolivet if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 1799371c9d4SSatish Balay for (PetscInt j = 0; j < A->rmap->n && flg; ++j) 1809371c9d4SSatish Balay if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE; 181c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 182c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 183c7a4214aSPierre Jolivet * [ B C E ] 184c7a4214aSPierre Jolivet * A = [ B D E ] 185c7a4214aSPierre Jolivet * [ B F E ] 186c7a4214aSPierre Jolivet */ 187c7a4214aSPierre Jolivet m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */ 1889566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock_Htool(A, &D)); 1899566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(D, ©)); 1909371c9d4SSatish 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 */ } 1919566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(D, ©)); 192c7a4214aSPierre Jolivet if (m) { 193c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */ 194c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 195b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 1969566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B)); 1979566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 1989566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT)); 1999566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 200b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2019566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 202c7a4214aSPierre Jolivet } else { 2037fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 2049566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 205c7a4214aSPierre Jolivet } 2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 208c7a4214aSPierre Jolivet } else { 209c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */ 210c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow); 211c7a4214aSPierre Jolivet } 212c7a4214aSPierre Jolivet } 213c7a4214aSPierre Jolivet } 214c7a4214aSPierre Jolivet if (m + A->rmap->n != nrow) { 215c7a4214aSPierre 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 */ 216c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 217b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2189566063dSJacob 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)); 2199566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 2209566063dSJacob 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)); 2219566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 222b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2239566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 224c7a4214aSPierre Jolivet } else { 2257fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 2269566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 227c7a4214aSPierre Jolivet } 2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 230c7a4214aSPierre Jolivet } else { 231c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */ 232c7a4214aSPierre 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); 233c7a4214aSPierre Jolivet } 234c7a4214aSPierre Jolivet } 235c7a4214aSPierre Jolivet } 236c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 237c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 238c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 239c7a4214aSPierre Jolivet } /* unsorted IS */ 240c7a4214aSPierre Jolivet } 241c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 242c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */ 2439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i], &idxr)); 2449566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i], &idxc)); 2459566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr)); 2469566063dSJacob Faibussowitsch PetscCall(MatScale((*submat)[i], a->s)); 247c7a4214aSPierre Jolivet } 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 249c7a4214aSPierre Jolivet } 250c7a4214aSPierre Jolivet 251d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A) 252d71ae5a4SJacob Faibussowitsch { 253c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 254c7a4214aSPierre Jolivet PetscContainer container; 255c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 256c7a4214aSPierre Jolivet 257c7a4214aSPierre Jolivet PetscFunctionBegin; 258f22e26b7SPierre Jolivet PetscCall(PetscObjectChangeTypeName((PetscObject)A, nullptr)); 259f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr)); 260f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr)); 261f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr)); 262f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr)); 263f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr)); 264f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr)); 265f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr)); 266f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr)); 267f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr)); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container)); 269c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 2709566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 2719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(kernelt)); 2739566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 274f22e26b7SPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr)); 275c7a4214aSPierre Jolivet } 27648a46eb9SPierre Jolivet if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree2(a->work_source, a->work_target)); 279c7a4214aSPierre Jolivet delete a->wrapper; 280c7a4214aSPierre Jolivet delete a->hmatrix; 2819566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 283c7a4214aSPierre Jolivet } 284c7a4214aSPierre Jolivet 285d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv) 286d71ae5a4SJacob Faibussowitsch { 287c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 288c7a4214aSPierre Jolivet PetscBool flg; 289c7a4214aSPierre Jolivet 290c7a4214aSPierre Jolivet PetscFunctionBegin; 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg)); 292c7a4214aSPierre Jolivet if (flg) { 2939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->hmatrix->get_symmetry_type())); 294c7a4214aSPierre Jolivet if (PetscAbsScalar(a->s - 1.0) > PETSC_MACHINE_EPSILON) { 295c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 2969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(a->s), (double)PetscImaginaryPart(a->s))); 297c7a4214aSPierre Jolivet #else 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)a->s)); 299c7a4214aSPierre Jolivet #endif 300c7a4214aSPierre Jolivet } 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum cluster size: %" PetscInt_FMT "\n", a->bs[0])); 3029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "maximum block size: %" PetscInt_FMT "\n", a->bs[1])); 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon)); 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta)); 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0])); 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1])); 3079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor])); 3089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering])); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", a->hmatrix->get_infos("Compression_ratio").c_str())); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", a->hmatrix->get_infos("Space_saving").c_str())); 3119566063dSJacob 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())); 3129371c9d4SSatish 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(), 3139371c9d4SSatish Balay a->hmatrix->get_infos("Dense_block_size_max").c_str())); 3149371c9d4SSatish 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(), 3159371c9d4SSatish Balay a->hmatrix->get_infos("Low_rank_block_size_max").c_str())); 3169566063dSJacob 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())); 317c7a4214aSPierre Jolivet } 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319c7a4214aSPierre Jolivet } 320c7a4214aSPierre Jolivet 321d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_Htool(Mat A, PetscScalar s) 322d71ae5a4SJacob Faibussowitsch { 323c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 324c7a4214aSPierre Jolivet 325c7a4214aSPierre Jolivet PetscFunctionBegin; 326c7a4214aSPierre Jolivet a->s *= s; 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328c7a4214aSPierre Jolivet } 329c7a4214aSPierre Jolivet 330c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 331d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 332d71ae5a4SJacob Faibussowitsch { 333c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 334c7a4214aSPierre Jolivet PetscInt *idxc; 335c7a4214aSPierre Jolivet PetscBLASInt one = 1, bn; 336c7a4214aSPierre Jolivet 337c7a4214aSPierre Jolivet PetscFunctionBegin; 338c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 339c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, &idxc)); 341c7a4214aSPierre Jolivet for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i; 342c7a4214aSPierre Jolivet } 343c7a4214aSPierre Jolivet if (idx) *idx = idxc; 344c7a4214aSPierre Jolivet if (v) { 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, v)); 346c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 347cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 3489566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N, &bn)); 349792fecdfSBarry Smith PetscCallBLAS("BLASscal", BLASscal_(&bn, &a->s, *v, &one)); 350c7a4214aSPierre Jolivet } 35148a46eb9SPierre Jolivet if (!idx) PetscCall(PetscFree(idxc)); 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 353c7a4214aSPierre Jolivet } 354c7a4214aSPierre Jolivet 355f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v) 356d71ae5a4SJacob Faibussowitsch { 357c7a4214aSPierre Jolivet PetscFunctionBegin; 35848a46eb9SPierre Jolivet if (idx) PetscCall(PetscFree(*idx)); 35948a46eb9SPierre Jolivet if (v) PetscCall(PetscFree(*v)); 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 361c7a4214aSPierre Jolivet } 362c7a4214aSPierre Jolivet 363d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject) 364d71ae5a4SJacob Faibussowitsch { 365c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 366c7a4214aSPierre Jolivet PetscInt n; 367c7a4214aSPierre Jolivet PetscBool flg; 368c7a4214aSPierre Jolivet 369c7a4214aSPierre Jolivet PetscFunctionBegin; 370d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Htool options"); 371f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", nullptr, a->bs[0], a->bs, nullptr)); 372f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_max_block_size", "Maximal number of coefficients in a dense block", nullptr, a->bs[1], a->bs + 1, nullptr)); 373f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr)); 374f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr)); 375f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr)); 376f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", nullptr, a->depth[1], a->depth + 1, nullptr)); 377c7a4214aSPierre Jolivet n = 0; 378dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg)); 379c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 38098e73e17SPierre Jolivet n = 0; 381dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg)); 38298e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 383d0609cedSBarry Smith PetscOptionsHeadEnd(); 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 385c7a4214aSPierre Jolivet } 386c7a4214aSPierre Jolivet 387cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType) 388d71ae5a4SJacob Faibussowitsch { 389c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 390c7a4214aSPierre Jolivet const PetscInt *ranges; 391c7a4214aSPierre Jolivet PetscInt *offset; 392c7a4214aSPierre Jolivet PetscMPIInt size; 393b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U'; 394cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 39598e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t, s = nullptr; 3963b9338faSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr; 397c7a4214aSPierre Jolivet 398c7a4214aSPierre Jolivet PetscFunctionBegin; 3999566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite)); 400c7a4214aSPierre Jolivet delete a->wrapper; 401c7a4214aSPierre Jolivet delete a->hmatrix; 4029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 4039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &offset)); 4049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 405c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 406c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 407c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 408c7a4214aSPierre Jolivet } 40998e73e17SPierre Jolivet switch (a->clustering) { 410d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 411d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 412d71ae5a4SJacob Faibussowitsch break; 413d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 414d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 415d71ae5a4SJacob Faibussowitsch break; 416d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 417d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 418d71ae5a4SJacob Faibussowitsch break; 419d71ae5a4SJacob Faibussowitsch default: 420d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 42198e73e17SPierre Jolivet } 422c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 4230429413eSPierre Jolivet t->build(A->rmap->N, a->gcoords_target, offset, -1, PetscObjectComm((PetscObject)A)); 424c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 425c7a4214aSPierre Jolivet else { 426f22e26b7SPierre Jolivet a->wrapper = nullptr; 427cab00e0dSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx); 428c7a4214aSPierre Jolivet } 429c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4309566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A, &ranges)); 431c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 432c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 433c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 434c7a4214aSPierre Jolivet } 43598e73e17SPierre Jolivet switch (a->clustering) { 436d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 437d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 438d71ae5a4SJacob Faibussowitsch break; 439d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 440d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 441d71ae5a4SJacob Faibussowitsch break; 442d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 443d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 444d71ae5a4SJacob Faibussowitsch break; 445d71ae5a4SJacob Faibussowitsch default: 446d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 44798e73e17SPierre Jolivet } 448c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 4490429413eSPierre Jolivet s->build(A->cmap->N, a->gcoords_source, offset, -1, PetscObjectComm((PetscObject)A)); 450c7a4214aSPierre Jolivet S = uplo = 'N'; 451c7a4214aSPierre Jolivet } 4529566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 453c7a4214aSPierre Jolivet switch (a->compressor) { 454d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_FULL_ACA: 455d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 456d71ae5a4SJacob Faibussowitsch break; 457d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_SVD: 458d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::SVD<PetscScalar>>(); 459d71ae5a4SJacob Faibussowitsch break; 460d71ae5a4SJacob Faibussowitsch default: 461d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 462c7a4214aSPierre Jolivet } 4630429413eSPierre 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))); 4643b9338faSPierre Jolivet a->hmatrix->set_compression(compressor); 465c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 466c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 467c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 4683b9338faSPierre Jolivet if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target, a->gcoords_source); 4693b9338faSPierre Jolivet else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target); 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471c7a4214aSPierre Jolivet } 472c7a4214aSPierre Jolivet 473d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C) 474d71ae5a4SJacob Faibussowitsch { 475c7a4214aSPierre Jolivet Mat_Product *product = C->product; 476c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)product->A->data; 477c7a4214aSPierre Jolivet const PetscScalar *in; 478c7a4214aSPierre Jolivet PetscScalar *out; 4798b8ddd11SPierre Jolivet PetscInt N, lda; 480c7a4214aSPierre Jolivet 481c7a4214aSPierre Jolivet PetscFunctionBegin; 482c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 483f22e26b7SPierre Jolivet PetscCall(MatGetSize(C, nullptr, &N)); 4849566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 4855f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 4869566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B, &in)); 4879566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &out)); 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)); 5009566063dSJacob Faibussowitsch PetscCall(MatScale(C, a->s)); 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 502c7a4214aSPierre Jolivet } 503c7a4214aSPierre Jolivet 504d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C) 505d71ae5a4SJacob Faibussowitsch { 506c7a4214aSPierre Jolivet Mat_Product *product = C->product; 507c7a4214aSPierre Jolivet Mat A, B; 508c7a4214aSPierre Jolivet PetscBool flg; 509c7a4214aSPierre Jolivet 510c7a4214aSPierre Jolivet PetscFunctionBegin; 511c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 512c7a4214aSPierre Jolivet A = product->A; 513c7a4214aSPierre Jolivet B = product->B; 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, "")); 515*97713f8cSPierre 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); 516*97713f8cSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 517*97713f8cSPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N)); 518*97713f8cSPierre Jolivet else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N)); 5198b8ddd11SPierre Jolivet } 5209566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 5219566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5229566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 5239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 525f22e26b7SPierre Jolivet C->ops->productsymbolic = nullptr; 526c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 5273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 528c7a4214aSPierre Jolivet } 529c7a4214aSPierre Jolivet 530d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 531d71ae5a4SJacob Faibussowitsch { 532c7a4214aSPierre Jolivet PetscFunctionBegin; 533c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 5348b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 5353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 536c7a4214aSPierre Jolivet } 537c7a4214aSPierre Jolivet 538d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) 539d71ae5a4SJacob Faibussowitsch { 540c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 541c7a4214aSPierre Jolivet 542c7a4214aSPierre Jolivet PetscFunctionBegin; 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 550c7a4214aSPierre Jolivet Input Parameter: 551c7a4214aSPierre Jolivet . A - hierarchical matrix 552c7a4214aSPierre Jolivet 553c7a4214aSPierre Jolivet Output Parameter: 554c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 555c7a4214aSPierre Jolivet 556c7a4214aSPierre Jolivet Level: advanced 557c7a4214aSPierre Jolivet 5581cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL` 559c7a4214aSPierre Jolivet @*/ 560d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) 561d71ae5a4SJacob Faibussowitsch { 562c7a4214aSPierre Jolivet PetscFunctionBegin; 563c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5644f572ea9SToby Isaac PetscAssertPointer(hmatrix, 2); 565cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::VirtualHMatrix<PetscScalar> **), (A, hmatrix)); 5663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 567c7a4214aSPierre Jolivet } 568c7a4214aSPierre Jolivet 569d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernel kernel, void *kernelctx) 570d71ae5a4SJacob Faibussowitsch { 571c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 572c7a4214aSPierre Jolivet 573c7a4214aSPierre Jolivet PetscFunctionBegin; 574c7a4214aSPierre Jolivet a->kernel = kernel; 575c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 576c7a4214aSPierre Jolivet delete a->wrapper; 577c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 5783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 579c7a4214aSPierre Jolivet } 580c7a4214aSPierre Jolivet 581c7a4214aSPierre Jolivet /*@C 58211a5261eSBarry Smith MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`. 583c7a4214aSPierre Jolivet 584c7a4214aSPierre Jolivet Input Parameters: 585c7a4214aSPierre Jolivet + A - hierarchical matrix 5862ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 5872ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 588c7a4214aSPierre Jolivet 589c7a4214aSPierre Jolivet Level: advanced 590c7a4214aSPierre Jolivet 5911cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()` 592c7a4214aSPierre Jolivet @*/ 593d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernel kernel, void *kernelctx) 594d71ae5a4SJacob Faibussowitsch { 595c7a4214aSPierre Jolivet PetscFunctionBegin; 596c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 597c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 2); 5984f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 3); 599cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernel, void *), (A, kernel, kernelctx)); 6003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 601c7a4214aSPierre Jolivet } 602c7a4214aSPierre Jolivet 603d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) 604d71ae5a4SJacob Faibussowitsch { 60598e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 60698e73e17SPierre Jolivet std::vector<PetscInt> source; 60798e73e17SPierre Jolivet 60898e73e17SPierre Jolivet PetscFunctionBegin; 609a9918087SPierre Jolivet source = a->hmatrix->get_source_cluster()->get_local_perm(); 6109566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), source.size(), source.data(), PETSC_COPY_VALUES, is)); 6119566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61398e73e17SPierre Jolivet } 61498e73e17SPierre Jolivet 61598e73e17SPierre Jolivet /*@C 61611a5261eSBarry Smith MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix. 61798e73e17SPierre Jolivet 61898e73e17SPierre Jolivet Input Parameter: 61998e73e17SPierre Jolivet . A - hierarchical matrix 62098e73e17SPierre Jolivet 62198e73e17SPierre Jolivet Output Parameter: 62298e73e17SPierre Jolivet . is - permutation 62398e73e17SPierre Jolivet 62498e73e17SPierre Jolivet Level: advanced 62598e73e17SPierre Jolivet 6261cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 62798e73e17SPierre Jolivet @*/ 628d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) 629d71ae5a4SJacob Faibussowitsch { 63098e73e17SPierre Jolivet PetscFunctionBegin; 63198e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6324f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 633cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is)); 6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63598e73e17SPierre Jolivet } 63698e73e17SPierre Jolivet 637d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) 638d71ae5a4SJacob Faibussowitsch { 63998e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 64098e73e17SPierre Jolivet std::vector<PetscInt> target; 64198e73e17SPierre Jolivet 64298e73e17SPierre Jolivet PetscFunctionBegin; 643a9918087SPierre Jolivet target = a->hmatrix->get_target_cluster()->get_local_perm(); 6449566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), target.size(), target.data(), PETSC_COPY_VALUES, is)); 6459566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64798e73e17SPierre Jolivet } 64898e73e17SPierre Jolivet 64998e73e17SPierre Jolivet /*@C 65011a5261eSBarry Smith MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix. 65198e73e17SPierre Jolivet 65298e73e17SPierre Jolivet Input Parameter: 65398e73e17SPierre Jolivet . A - hierarchical matrix 65498e73e17SPierre Jolivet 65598e73e17SPierre Jolivet Output Parameter: 65698e73e17SPierre Jolivet . is - permutation 65798e73e17SPierre Jolivet 65898e73e17SPierre Jolivet Level: advanced 65998e73e17SPierre Jolivet 6601cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 66198e73e17SPierre Jolivet @*/ 662d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) 663d71ae5a4SJacob Faibussowitsch { 66498e73e17SPierre Jolivet PetscFunctionBegin; 66598e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6664f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 667cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is)); 6683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66998e73e17SPierre Jolivet } 67098e73e17SPierre Jolivet 671d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) 672d71ae5a4SJacob Faibussowitsch { 67398e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 67498e73e17SPierre Jolivet 67598e73e17SPierre Jolivet PetscFunctionBegin; 67698e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67898e73e17SPierre Jolivet } 67998e73e17SPierre Jolivet 68098e73e17SPierre Jolivet /*@C 68111a5261eSBarry Smith MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation. 68298e73e17SPierre Jolivet 68398e73e17SPierre Jolivet Input Parameters: 68498e73e17SPierre Jolivet + A - hierarchical matrix 68598e73e17SPierre Jolivet - use - Boolean value 68698e73e17SPierre Jolivet 68798e73e17SPierre Jolivet Level: advanced 68898e73e17SPierre Jolivet 6891cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 69098e73e17SPierre Jolivet @*/ 691d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) 692d71ae5a4SJacob Faibussowitsch { 69398e73e17SPierre Jolivet PetscFunctionBegin; 69498e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 69598e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 696cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use)); 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69898e73e17SPierre Jolivet } 69998e73e17SPierre Jolivet 700cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 701d71ae5a4SJacob Faibussowitsch { 702c7a4214aSPierre Jolivet Mat C; 703c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 704c7a4214aSPierre Jolivet PetscInt lda; 705c7a4214aSPierre Jolivet PetscScalar *array; 706c7a4214aSPierre Jolivet 707c7a4214aSPierre Jolivet PetscFunctionBegin; 708c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 709c7a4214aSPierre Jolivet C = *B; 7105f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions"); 7119566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 7125f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 713c7a4214aSPierre Jolivet } else { 7149566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7159566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 7169566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 7179566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7189566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 719c7a4214aSPierre Jolivet } 7209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 7219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 7229566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 723c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 7249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 7259566063dSJacob Faibussowitsch PetscCall(MatScale(C, a->s)); 726c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7279566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &C)); 728c7a4214aSPierre Jolivet } else *B = C; 7293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 730c7a4214aSPierre Jolivet } 731c7a4214aSPierre Jolivet 732d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx) 733d71ae5a4SJacob Faibussowitsch { 734c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx; 73598e73e17SPierre Jolivet PetscScalar *tmp; 73698e73e17SPierre Jolivet 73798e73e17SPierre Jolivet PetscFunctionBegin; 7383ba16761SJacob Faibussowitsch PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx)); 7399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * N, &tmp)); 7409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, ptr, M * N)); 74198e73e17SPierre Jolivet for (PetscInt i = 0; i < M; ++i) { 74298e73e17SPierre Jolivet for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N]; 74398e73e17SPierre Jolivet } 7449566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 7453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 746c7a4214aSPierre Jolivet } 747c7a4214aSPierre Jolivet 748c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 749d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) 750d71ae5a4SJacob Faibussowitsch { 751c7a4214aSPierre Jolivet Mat C; 752c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data, *c; 753c7a4214aSPierre Jolivet PetscInt M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n; 754c7a4214aSPierre Jolivet PetscContainer container; 755c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 756c7a4214aSPierre Jolivet 757c7a4214aSPierre Jolivet PetscFunctionBegin; 7587fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 7595f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported"); 760c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7619566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, n, m, N, M)); 7639566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 7649566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7659566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C), &container)); 7669566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 7679566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container, kernelt)); 7689566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)C, "KernelTranspose", (PetscObject)container)); 769c7a4214aSPierre Jolivet } else { 770c7a4214aSPierre Jolivet C = *B; 7719566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container)); 7725f80ce2aSJacob Faibussowitsch PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 7739566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 774c7a4214aSPierre Jolivet } 775c7a4214aSPierre Jolivet c = (Mat_Htool *)C->data; 776c7a4214aSPierre Jolivet c->dim = a->dim; 777c7a4214aSPierre Jolivet c->s = a->s; 77898e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 779c7a4214aSPierre Jolivet if (kernelt->A != A) { 7809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 781c7a4214aSPierre Jolivet kernelt->A = A; 7829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 783c7a4214aSPierre Jolivet } 784c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 785c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 786c7a4214aSPierre Jolivet c->kernelctx = kernelt; 787c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target)); 7899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim)); 790c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 7919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source)); 7929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim)); 793c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 7949566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target)); 795c7a4214aSPierre Jolivet } 7969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 7979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 798c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 800c7a4214aSPierre Jolivet } 801c7a4214aSPierre Jolivet 802c7a4214aSPierre Jolivet /*@C 80311a5261eSBarry Smith MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel. 804c7a4214aSPierre Jolivet 805c7a4214aSPierre Jolivet Input Parameters: 806c7a4214aSPierre Jolivet + comm - MPI communicator 8072ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 8082ef1f0ffSBarry Smith . n - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given) 8092ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 8102ef1f0ffSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 811c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 812c7a4214aSPierre Jolivet . coords_target - coordinates of the target 813c7a4214aSPierre Jolivet . coords_source - coordinates of the source 8142ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 8152ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 816c7a4214aSPierre Jolivet 817c7a4214aSPierre Jolivet Output Parameter: 818c7a4214aSPierre Jolivet . B - matrix 819c7a4214aSPierre Jolivet 820c7a4214aSPierre Jolivet Options Database Keys: 82111a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`> - minimal leaf size in cluster tree 82211a5261eSBarry Smith . -mat_htool_max_block_size <`PetscInt`> - maximal number of coefficients in a dense block 82311a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`> - relative error in Frobenius norm when approximating a block 82411a5261eSBarry Smith . -mat_htool_eta <`PetscReal`> - admissibility condition tolerance 82511a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`> - minimal cluster tree depth associated with the rows 82611a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`> - minimal cluster tree depth associated with the columns 82798e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 82898e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 829c7a4214aSPierre Jolivet 830c7a4214aSPierre Jolivet Level: intermediate 831c7a4214aSPierre Jolivet 8321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 833c7a4214aSPierre Jolivet @*/ 834d71ae5a4SJacob Faibussowitsch 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) 835d71ae5a4SJacob Faibussowitsch { 836c7a4214aSPierre Jolivet Mat A; 837c7a4214aSPierre Jolivet Mat_Htool *a; 838c7a4214aSPierre Jolivet 839c7a4214aSPierre Jolivet PetscFunctionBegin; 8409566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 841c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A, spacedim, 6); 8424f572ea9SToby Isaac PetscAssertPointer(coords_target, 7); 8434f572ea9SToby Isaac PetscAssertPointer(coords_source, 8); 844c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 9); 8454f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 10); 8469566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, M, N)); 8479566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATHTOOL)); 8489566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 849c7a4214aSPierre Jolivet a = (Mat_Htool *)A->data; 850c7a4214aSPierre Jolivet a->dim = spacedim; 851c7a4214aSPierre Jolivet a->s = 1.0; 852c7a4214aSPierre Jolivet a->kernel = kernel; 853c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 8549566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target)); 8559566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim)); 8561c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */ 857c7a4214aSPierre Jolivet if (coords_target != coords_source) { 8589566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source)); 8599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim)); 8601c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */ 861c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 8629566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target)); 863c7a4214aSPierre Jolivet *B = A; 8643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 865c7a4214aSPierre Jolivet } 866c7a4214aSPierre Jolivet 867c7a4214aSPierre Jolivet /*MC 868c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 869c7a4214aSPierre Jolivet 8702ef1f0ffSBarry Smith Use `./configure --download-htool` to install PETSc to use Htool. 871c7a4214aSPierre Jolivet 8722ef1f0ffSBarry Smith Options Database Key: 8732ef1f0ffSBarry Smith . -mat_type htool - matrix type to `MATHTOOL` 874c7a4214aSPierre Jolivet 875c7a4214aSPierre Jolivet Level: beginner 876c7a4214aSPierre Jolivet 8771cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 878c7a4214aSPierre Jolivet M*/ 879d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 880d71ae5a4SJacob Faibussowitsch { 881c7a4214aSPierre Jolivet Mat_Htool *a; 882c7a4214aSPierre Jolivet 883c7a4214aSPierre Jolivet PetscFunctionBegin; 8844dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 885c7a4214aSPierre Jolivet A->data = (void *)a; 8869566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL)); 8879566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps))); 888c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 889c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 890c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 891c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 8928b8ddd11SPierre Jolivet A->ops->multtranspose = MatMultTranspose_Htool; 8938b8ddd11SPierre Jolivet if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool; 894c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 895c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 896c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 897c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 898c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 899c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 900c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 901c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 902c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 903c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 904c7a4214aSPierre Jolivet a->dim = 0; 905f22e26b7SPierre Jolivet a->gcoords_target = nullptr; 906f22e26b7SPierre Jolivet a->gcoords_source = nullptr; 907c7a4214aSPierre Jolivet a->s = 1.0; 908c7a4214aSPierre Jolivet a->bs[0] = 10; 909c7a4214aSPierre Jolivet a->bs[1] = 1000000; 910c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 911c7a4214aSPierre Jolivet a->eta = 10.0; 912c7a4214aSPierre Jolivet a->depth[0] = 0; 913c7a4214aSPierre Jolivet a->depth[1] = 0; 914c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 9159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool)); 9169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool)); 9179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense)); 9189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense)); 9199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool)); 9209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool)); 9219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool)); 9229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool)); 9239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool)); 9243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 925c7a4214aSPierre Jolivet } 926