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)); 32c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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) { 479566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, NULL, &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; 57c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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)); 73c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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 */ 91c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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)); 1078b8ddd11SPierre Jolivet PetscFunctionReturn(0); 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 } 149c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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)); 16948a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, NULL, (*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 } 248c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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; 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", NULL)); 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", NULL)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", NULL)); 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", NULL)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", NULL)); 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)); 2749566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", NULL)); 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)); 282c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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 } 318c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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; 327c7a4214aSPierre Jolivet PetscFunctionReturn(0); 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)); 352c7a4214aSPierre Jolivet PetscFunctionReturn(0); 353c7a4214aSPierre Jolivet } 354c7a4214aSPierre Jolivet 355*cedd07caSPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *nz, PetscInt **idx, PetscScalar **v) 356d71ae5a4SJacob Faibussowitsch { 357c7a4214aSPierre Jolivet PetscFunctionBegin; 358c7a4214aSPierre Jolivet if (nz) *nz = 0; 35948a46eb9SPierre Jolivet if (idx) PetscCall(PetscFree(*idx)); 36048a46eb9SPierre Jolivet if (v) PetscCall(PetscFree(*v)); 361c7a4214aSPierre Jolivet PetscFunctionReturn(0); 362c7a4214aSPierre Jolivet } 363c7a4214aSPierre Jolivet 364d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject) 365d71ae5a4SJacob Faibussowitsch { 366c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 367c7a4214aSPierre Jolivet PetscInt n; 368c7a4214aSPierre Jolivet PetscBool flg; 369c7a4214aSPierre Jolivet 370c7a4214aSPierre Jolivet PetscFunctionBegin; 371d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Htool options"); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", NULL, a->bs[0], a->bs, NULL)); 3739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_max_block_size", "Maximal number of coefficients in a dense block", NULL, a->bs[1], a->bs + 1, NULL)); 3749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", NULL, a->epsilon, &a->epsilon, NULL)); 3759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", NULL, a->eta, &a->eta, NULL)); 3769566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", NULL, a->depth[0], a->depth, NULL)); 3779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", NULL, a->depth[1], a->depth + 1, NULL)); 378c7a4214aSPierre Jolivet n = 0; 379dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg)); 380c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 38198e73e17SPierre Jolivet n = 0; 382dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg)); 38398e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 384d0609cedSBarry Smith PetscOptionsHeadEnd(); 385c7a4214aSPierre Jolivet PetscFunctionReturn(0); 386c7a4214aSPierre Jolivet } 387c7a4214aSPierre Jolivet 388*cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType) 389d71ae5a4SJacob Faibussowitsch { 390c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 391c7a4214aSPierre Jolivet const PetscInt *ranges; 392c7a4214aSPierre Jolivet PetscInt *offset; 393c7a4214aSPierre Jolivet PetscMPIInt size; 394b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U'; 395cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 39698e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t, s = nullptr; 3973b9338faSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr; 398c7a4214aSPierre Jolivet 399c7a4214aSPierre Jolivet PetscFunctionBegin; 4009566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite)); 401c7a4214aSPierre Jolivet delete a->wrapper; 402c7a4214aSPierre Jolivet delete a->hmatrix; 4039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 4049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &offset)); 4059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 406c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 407c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 408c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 409c7a4214aSPierre Jolivet } 41098e73e17SPierre Jolivet switch (a->clustering) { 411d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 412d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 413d71ae5a4SJacob Faibussowitsch break; 414d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 415d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 416d71ae5a4SJacob Faibussowitsch break; 417d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 418d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 419d71ae5a4SJacob Faibussowitsch break; 420d71ae5a4SJacob Faibussowitsch default: 421d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 42298e73e17SPierre Jolivet } 423c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 42498e73e17SPierre Jolivet t->build(A->rmap->N, a->gcoords_target, offset); 425c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 426c7a4214aSPierre Jolivet else { 427c7a4214aSPierre Jolivet a->wrapper = NULL; 428cab00e0dSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx); 429c7a4214aSPierre Jolivet } 430c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A, &ranges)); 432c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 433c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 434c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 435c7a4214aSPierre Jolivet } 43698e73e17SPierre Jolivet switch (a->clustering) { 437d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 438d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 439d71ae5a4SJacob Faibussowitsch break; 440d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 441d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 442d71ae5a4SJacob Faibussowitsch break; 443d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 444d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 445d71ae5a4SJacob Faibussowitsch break; 446d71ae5a4SJacob Faibussowitsch default: 447d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 44898e73e17SPierre Jolivet } 449c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 45098e73e17SPierre Jolivet s->build(A->cmap->N, a->gcoords_source, offset); 451c7a4214aSPierre Jolivet S = uplo = 'N'; 452c7a4214aSPierre Jolivet } 4539566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 454c7a4214aSPierre Jolivet switch (a->compressor) { 455d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_FULL_ACA: 456d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 457d71ae5a4SJacob Faibussowitsch break; 458d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_SVD: 459d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::SVD<PetscScalar>>(); 460d71ae5a4SJacob Faibussowitsch break; 461d71ae5a4SJacob Faibussowitsch default: 462d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 463c7a4214aSPierre Jolivet } 4643b9338faSPierre Jolivet a->hmatrix = dynamic_cast<htool::VirtualHMatrix<PetscScalar> *>(new htool::HMatrix<PetscScalar>(t, s ? s : t, a->epsilon, a->eta, S, uplo)); 4653b9338faSPierre Jolivet a->hmatrix->set_compression(compressor); 466c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 467c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 468c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 4693b9338faSPierre Jolivet if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target, a->gcoords_source); 4703b9338faSPierre Jolivet else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target); 471c7a4214aSPierre Jolivet PetscFunctionReturn(0); 472c7a4214aSPierre Jolivet } 473c7a4214aSPierre Jolivet 474d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C) 475d71ae5a4SJacob Faibussowitsch { 476c7a4214aSPierre Jolivet Mat_Product *product = C->product; 477c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)product->A->data; 478c7a4214aSPierre Jolivet const PetscScalar *in; 479c7a4214aSPierre Jolivet PetscScalar *out; 4808b8ddd11SPierre Jolivet PetscInt N, lda; 481c7a4214aSPierre Jolivet 482c7a4214aSPierre Jolivet PetscFunctionBegin; 483c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 4849566063dSJacob Faibussowitsch PetscCall(MatGetSize(C, NULL, &N)); 4859566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 4865f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 4879566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B, &in)); 4889566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &out)); 4898b8ddd11SPierre Jolivet switch (product->type) { 490d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AB: 491d71ae5a4SJacob Faibussowitsch a->hmatrix->mvprod_local_to_local(in, out, N); 492d71ae5a4SJacob Faibussowitsch break; 493d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AtB: 494d71ae5a4SJacob Faibussowitsch a->hmatrix->mvprod_transp_local_to_local(in, out, N); 495d71ae5a4SJacob Faibussowitsch break; 496d71ae5a4SJacob Faibussowitsch default: 497d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]); 498c7a4214aSPierre Jolivet } 4999566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &out)); 5009566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B, &in)); 5019566063dSJacob Faibussowitsch PetscCall(MatScale(C, a->s)); 502c7a4214aSPierre Jolivet PetscFunctionReturn(0); 503c7a4214aSPierre Jolivet } 504c7a4214aSPierre Jolivet 505d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C) 506d71ae5a4SJacob Faibussowitsch { 507c7a4214aSPierre Jolivet Mat_Product *product = C->product; 508c7a4214aSPierre Jolivet Mat A, B; 509c7a4214aSPierre Jolivet PetscBool flg; 510c7a4214aSPierre Jolivet 511c7a4214aSPierre Jolivet PetscFunctionBegin; 512c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 513c7a4214aSPierre Jolivet A = product->A; 514c7a4214aSPierre Jolivet B = product->B; 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, "")); 5165f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MatProduct_AB not supported for %s", ((PetscObject)product->B)->type_name); 517c7a4214aSPierre Jolivet switch (product->type) { 518c7a4214aSPierre Jolivet case MATPRODUCT_AB: 51948a46eb9SPierre 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)); 5208b8ddd11SPierre Jolivet break; 5218b8ddd11SPierre Jolivet case MATPRODUCT_AtB: 52248a46eb9SPierre 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)); 5238b8ddd11SPierre Jolivet break; 524d71ae5a4SJacob Faibussowitsch default: 525d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "ProductType %s is not supported", MatProductTypes[product->type]); 5268b8ddd11SPierre Jolivet } 5279566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 5289566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5299566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 5309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 532c7a4214aSPierre Jolivet C->ops->productsymbolic = NULL; 533c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 534c7a4214aSPierre Jolivet PetscFunctionReturn(0); 535c7a4214aSPierre Jolivet } 536c7a4214aSPierre Jolivet 537d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 538d71ae5a4SJacob Faibussowitsch { 539c7a4214aSPierre Jolivet PetscFunctionBegin; 540c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 5418b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 542c7a4214aSPierre Jolivet PetscFunctionReturn(0); 543c7a4214aSPierre Jolivet } 544c7a4214aSPierre Jolivet 545d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) 546d71ae5a4SJacob Faibussowitsch { 547c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 548c7a4214aSPierre Jolivet 549c7a4214aSPierre Jolivet PetscFunctionBegin; 550c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 551c7a4214aSPierre Jolivet PetscFunctionReturn(0); 552c7a4214aSPierre Jolivet } 553c7a4214aSPierre Jolivet 554c7a4214aSPierre Jolivet /*@C 55511a5261eSBarry Smith MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`. 556c7a4214aSPierre Jolivet 557c7a4214aSPierre Jolivet Input Parameter: 558c7a4214aSPierre Jolivet . A - hierarchical matrix 559c7a4214aSPierre Jolivet 560c7a4214aSPierre Jolivet Output Parameter: 561c7a4214aSPierre Jolivet . hmatrix - opaque pointer to a Htool virtual matrix 562c7a4214aSPierre Jolivet 563c7a4214aSPierre Jolivet Level: advanced 564c7a4214aSPierre Jolivet 565db781477SPatrick Sanan .seealso: `MATHTOOL` 566c7a4214aSPierre Jolivet @*/ 567d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) 568d71ae5a4SJacob Faibussowitsch { 569c7a4214aSPierre Jolivet PetscFunctionBegin; 570c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 571c7a4214aSPierre Jolivet PetscValidPointer(hmatrix, 2); 572cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::VirtualHMatrix<PetscScalar> **), (A, hmatrix)); 573c7a4214aSPierre Jolivet PetscFunctionReturn(0); 574c7a4214aSPierre Jolivet } 575c7a4214aSPierre Jolivet 576d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernel kernel, void *kernelctx) 577d71ae5a4SJacob Faibussowitsch { 578c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 579c7a4214aSPierre Jolivet 580c7a4214aSPierre Jolivet PetscFunctionBegin; 581c7a4214aSPierre Jolivet a->kernel = kernel; 582c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 583c7a4214aSPierre Jolivet delete a->wrapper; 584c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 585c7a4214aSPierre Jolivet PetscFunctionReturn(0); 586c7a4214aSPierre Jolivet } 587c7a4214aSPierre Jolivet 588c7a4214aSPierre Jolivet /*@C 58911a5261eSBarry Smith MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`. 590c7a4214aSPierre Jolivet 591c7a4214aSPierre Jolivet Input Parameters: 592c7a4214aSPierre Jolivet + A - hierarchical matrix 593c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 594cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 595c7a4214aSPierre Jolivet 596c7a4214aSPierre Jolivet Level: advanced 597c7a4214aSPierre Jolivet 598db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatCreateHtoolFromKernel()` 599c7a4214aSPierre Jolivet @*/ 600d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernel kernel, void *kernelctx) 601d71ae5a4SJacob Faibussowitsch { 602c7a4214aSPierre Jolivet PetscFunctionBegin; 603c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 604c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 2); 605c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx, 3); 606cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernel, void *), (A, kernel, kernelctx)); 607c7a4214aSPierre Jolivet PetscFunctionReturn(0); 608c7a4214aSPierre Jolivet } 609c7a4214aSPierre Jolivet 610d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) 611d71ae5a4SJacob Faibussowitsch { 61298e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 61398e73e17SPierre Jolivet std::vector<PetscInt> source; 61498e73e17SPierre Jolivet 61598e73e17SPierre Jolivet PetscFunctionBegin; 616a9918087SPierre Jolivet source = a->hmatrix->get_source_cluster()->get_local_perm(); 6179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), source.size(), source.data(), PETSC_COPY_VALUES, is)); 6189566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 61998e73e17SPierre Jolivet PetscFunctionReturn(0); 62098e73e17SPierre Jolivet } 62198e73e17SPierre Jolivet 62298e73e17SPierre Jolivet /*@C 62311a5261eSBarry Smith MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix. 62498e73e17SPierre Jolivet 62598e73e17SPierre Jolivet Input Parameter: 62698e73e17SPierre Jolivet . A - hierarchical matrix 62798e73e17SPierre Jolivet 62898e73e17SPierre Jolivet Output Parameter: 62998e73e17SPierre Jolivet . is - permutation 63098e73e17SPierre Jolivet 63198e73e17SPierre Jolivet Level: advanced 63298e73e17SPierre Jolivet 633db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 63498e73e17SPierre Jolivet @*/ 635d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) 636d71ae5a4SJacob Faibussowitsch { 63798e73e17SPierre Jolivet PetscFunctionBegin; 63898e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 63998e73e17SPierre Jolivet if (!is) PetscValidPointer(is, 2); 640cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is)); 64198e73e17SPierre Jolivet PetscFunctionReturn(0); 64298e73e17SPierre Jolivet } 64398e73e17SPierre Jolivet 644d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) 645d71ae5a4SJacob Faibussowitsch { 64698e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 64798e73e17SPierre Jolivet std::vector<PetscInt> target; 64898e73e17SPierre Jolivet 64998e73e17SPierre Jolivet PetscFunctionBegin; 650a9918087SPierre Jolivet target = a->hmatrix->get_target_cluster()->get_local_perm(); 6519566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), target.size(), target.data(), PETSC_COPY_VALUES, is)); 6529566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 65398e73e17SPierre Jolivet PetscFunctionReturn(0); 65498e73e17SPierre Jolivet } 65598e73e17SPierre Jolivet 65698e73e17SPierre Jolivet /*@C 65711a5261eSBarry Smith MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix. 65898e73e17SPierre Jolivet 65998e73e17SPierre Jolivet Input Parameter: 66098e73e17SPierre Jolivet . A - hierarchical matrix 66198e73e17SPierre Jolivet 66298e73e17SPierre Jolivet Output Parameter: 66398e73e17SPierre Jolivet . is - permutation 66498e73e17SPierre Jolivet 66598e73e17SPierre Jolivet Level: advanced 66698e73e17SPierre Jolivet 667db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 66898e73e17SPierre Jolivet @*/ 669d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) 670d71ae5a4SJacob Faibussowitsch { 67198e73e17SPierre Jolivet PetscFunctionBegin; 67298e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 67398e73e17SPierre Jolivet if (!is) PetscValidPointer(is, 2); 674cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is)); 67598e73e17SPierre Jolivet PetscFunctionReturn(0); 67698e73e17SPierre Jolivet } 67798e73e17SPierre Jolivet 678d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) 679d71ae5a4SJacob Faibussowitsch { 68098e73e17SPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 68198e73e17SPierre Jolivet 68298e73e17SPierre Jolivet PetscFunctionBegin; 68398e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 68498e73e17SPierre Jolivet PetscFunctionReturn(0); 68598e73e17SPierre Jolivet } 68698e73e17SPierre Jolivet 68798e73e17SPierre Jolivet /*@C 68811a5261eSBarry Smith MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation. 68998e73e17SPierre Jolivet 69098e73e17SPierre Jolivet Input Parameters: 69198e73e17SPierre Jolivet + A - hierarchical matrix 69298e73e17SPierre Jolivet - use - Boolean value 69398e73e17SPierre Jolivet 69498e73e17SPierre Jolivet Level: advanced 69598e73e17SPierre Jolivet 696db781477SPatrick Sanan .seealso: `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 69798e73e17SPierre Jolivet @*/ 698d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) 699d71ae5a4SJacob Faibussowitsch { 70098e73e17SPierre Jolivet PetscFunctionBegin; 70198e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 70298e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 703cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use)); 70498e73e17SPierre Jolivet PetscFunctionReturn(0); 70598e73e17SPierre Jolivet } 70698e73e17SPierre Jolivet 707*cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 708d71ae5a4SJacob Faibussowitsch { 709c7a4214aSPierre Jolivet Mat C; 710c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data; 711c7a4214aSPierre Jolivet PetscInt lda; 712c7a4214aSPierre Jolivet PetscScalar *array; 713c7a4214aSPierre Jolivet 714c7a4214aSPierre Jolivet PetscFunctionBegin; 715c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 716c7a4214aSPierre Jolivet C = *B; 7175f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions"); 7189566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 7195f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 720c7a4214aSPierre Jolivet } else { 7219566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 7239566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 7249566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7259566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 726c7a4214aSPierre Jolivet } 7279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 7289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 7299566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 730c7a4214aSPierre Jolivet a->hmatrix->copy_local_dense_perm(array); 7319566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 7329566063dSJacob Faibussowitsch PetscCall(MatScale(C, a->s)); 733c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7349566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &C)); 735c7a4214aSPierre Jolivet } else *B = C; 736c7a4214aSPierre Jolivet PetscFunctionReturn(0); 737c7a4214aSPierre Jolivet } 738c7a4214aSPierre Jolivet 739d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx) 740d71ae5a4SJacob Faibussowitsch { 741c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx; 74298e73e17SPierre Jolivet PetscScalar *tmp; 74398e73e17SPierre Jolivet 74498e73e17SPierre Jolivet PetscFunctionBegin; 74598e73e17SPierre Jolivet generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx); 7469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * N, &tmp)); 7479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, ptr, M * N)); 74898e73e17SPierre Jolivet for (PetscInt i = 0; i < M; ++i) { 74998e73e17SPierre Jolivet for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N]; 75098e73e17SPierre Jolivet } 7519566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 75298e73e17SPierre Jolivet PetscFunctionReturn(0); 753c7a4214aSPierre Jolivet } 754c7a4214aSPierre Jolivet 755c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 756d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) 757d71ae5a4SJacob Faibussowitsch { 758c7a4214aSPierre Jolivet Mat C; 759c7a4214aSPierre Jolivet Mat_Htool *a = (Mat_Htool *)A->data, *c; 760c7a4214aSPierre Jolivet PetscInt M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n; 761c7a4214aSPierre Jolivet PetscContainer container; 762c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 763c7a4214aSPierre Jolivet 764c7a4214aSPierre Jolivet PetscFunctionBegin; 7657fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 7665f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported"); 767c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7689566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7699566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, n, m, N, M)); 7709566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 7719566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7729566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C), &container)); 7739566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 7749566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container, kernelt)); 7759566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)C, "KernelTranspose", (PetscObject)container)); 776c7a4214aSPierre Jolivet } else { 777c7a4214aSPierre Jolivet C = *B; 7789566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container)); 7795f80ce2aSJacob Faibussowitsch PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 7809566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 781c7a4214aSPierre Jolivet } 782c7a4214aSPierre Jolivet c = (Mat_Htool *)C->data; 783c7a4214aSPierre Jolivet c->dim = a->dim; 784c7a4214aSPierre Jolivet c->s = a->s; 78598e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 786c7a4214aSPierre Jolivet if (kernelt->A != A) { 7879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 788c7a4214aSPierre Jolivet kernelt->A = A; 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 790c7a4214aSPierre Jolivet } 791c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 792c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 793c7a4214aSPierre Jolivet c->kernelctx = kernelt; 794c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target)); 7969566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim)); 797c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 7989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source)); 7999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim)); 800c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 8019566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target)); 802c7a4214aSPierre Jolivet } 8039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 8049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 805c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 806c7a4214aSPierre Jolivet PetscFunctionReturn(0); 807c7a4214aSPierre Jolivet } 808c7a4214aSPierre Jolivet 809c7a4214aSPierre Jolivet /*@C 81011a5261eSBarry Smith MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel. 811c7a4214aSPierre Jolivet 812c7a4214aSPierre Jolivet Input Parameters: 813c7a4214aSPierre Jolivet + comm - MPI communicator 81411a5261eSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given) 81511a5261eSBarry Smith . n - number of local columns (or `PETSC_DECIDE` to have calculated if N is given) 81611a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given) 81711a5261eSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given) 818c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 819c7a4214aSPierre Jolivet . coords_target - coordinates of the target 820c7a4214aSPierre Jolivet . coords_source - coordinates of the source 821c7a4214aSPierre Jolivet . kernel - computational kernel (or NULL) 822cab00e0dSPierre Jolivet - kernelctx - kernel context (if kernel is NULL, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 823c7a4214aSPierre Jolivet 824c7a4214aSPierre Jolivet Output Parameter: 825c7a4214aSPierre Jolivet . B - matrix 826c7a4214aSPierre Jolivet 827c7a4214aSPierre Jolivet Options Database Keys: 82811a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`> - minimal leaf size in cluster tree 82911a5261eSBarry Smith . -mat_htool_max_block_size <`PetscInt`> - maximal number of coefficients in a dense block 83011a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`> - relative error in Frobenius norm when approximating a block 83111a5261eSBarry Smith . -mat_htool_eta <`PetscReal`> - admissibility condition tolerance 83211a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`> - minimal cluster tree depth associated with the rows 83311a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`> - minimal cluster tree depth associated with the columns 83498e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 83598e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 836c7a4214aSPierre Jolivet 837c7a4214aSPierre Jolivet Level: intermediate 838c7a4214aSPierre Jolivet 839db781477SPatrick Sanan .seealso: `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 840c7a4214aSPierre Jolivet @*/ 841d71ae5a4SJacob 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) 842d71ae5a4SJacob Faibussowitsch { 843c7a4214aSPierre Jolivet Mat A; 844c7a4214aSPierre Jolivet Mat_Htool *a; 845c7a4214aSPierre Jolivet 846c7a4214aSPierre Jolivet PetscFunctionBegin; 8479566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 848c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A, spacedim, 6); 849c7a4214aSPierre Jolivet PetscValidRealPointer(coords_target, 7); 850c7a4214aSPierre Jolivet PetscValidRealPointer(coords_source, 8); 851c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 9); 852c7a4214aSPierre Jolivet if (!kernel) PetscValidPointer(kernelctx, 10); 8539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, M, N)); 8549566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATHTOOL)); 8559566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 856c7a4214aSPierre Jolivet a = (Mat_Htool *)A->data; 857c7a4214aSPierre Jolivet a->dim = spacedim; 858c7a4214aSPierre Jolivet a->s = 1.0; 859c7a4214aSPierre Jolivet a->kernel = kernel; 860c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 8619566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target)); 8629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim)); 8631c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */ 864c7a4214aSPierre Jolivet if (coords_target != coords_source) { 8659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source)); 8669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim)); 8671c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */ 868c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 8699566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target)); 870c7a4214aSPierre Jolivet *B = A; 871c7a4214aSPierre Jolivet PetscFunctionReturn(0); 872c7a4214aSPierre Jolivet } 873c7a4214aSPierre Jolivet 874c7a4214aSPierre Jolivet /*MC 875c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 876c7a4214aSPierre Jolivet 877c7a4214aSPierre Jolivet Use ./configure --download-htool to install PETSc to use Htool. 878c7a4214aSPierre Jolivet 879c7a4214aSPierre Jolivet Options Database Keys: 88011a5261eSBarry Smith . -mat_type htool - matrix type to `MATHTOOL` during a call to `MatSetFromOptions()` 881c7a4214aSPierre Jolivet 882c7a4214aSPierre Jolivet Level: beginner 883c7a4214aSPierre Jolivet 884db781477SPatrick Sanan .seealso: `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 885c7a4214aSPierre Jolivet M*/ 886d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 887d71ae5a4SJacob Faibussowitsch { 888c7a4214aSPierre Jolivet Mat_Htool *a; 889c7a4214aSPierre Jolivet 890c7a4214aSPierre Jolivet PetscFunctionBegin; 8914dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 892c7a4214aSPierre Jolivet A->data = (void *)a; 8939566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL)); 8949566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps))); 895c7a4214aSPierre Jolivet A->ops->getdiagonal = MatGetDiagonal_Htool; 896c7a4214aSPierre Jolivet A->ops->getdiagonalblock = MatGetDiagonalBlock_Htool; 897c7a4214aSPierre Jolivet A->ops->mult = MatMult_Htool; 898c7a4214aSPierre Jolivet A->ops->multadd = MatMultAdd_Htool; 8998b8ddd11SPierre Jolivet A->ops->multtranspose = MatMultTranspose_Htool; 9008b8ddd11SPierre Jolivet if (!PetscDefined(USE_COMPLEX)) A->ops->multhermitiantranspose = MatMultTranspose_Htool; 901c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 902c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 903c7a4214aSPierre Jolivet A->ops->transpose = MatTranspose_Htool; 904c7a4214aSPierre Jolivet A->ops->destroy = MatDestroy_Htool; 905c7a4214aSPierre Jolivet A->ops->view = MatView_Htool; 906c7a4214aSPierre Jolivet A->ops->setfromoptions = MatSetFromOptions_Htool; 907c7a4214aSPierre Jolivet A->ops->scale = MatScale_Htool; 908c7a4214aSPierre Jolivet A->ops->getrow = MatGetRow_Htool; 909c7a4214aSPierre Jolivet A->ops->restorerow = MatRestoreRow_Htool; 910c7a4214aSPierre Jolivet A->ops->assemblyend = MatAssemblyEnd_Htool; 911c7a4214aSPierre Jolivet a->dim = 0; 912c7a4214aSPierre Jolivet a->gcoords_target = NULL; 913c7a4214aSPierre Jolivet a->gcoords_source = NULL; 914c7a4214aSPierre Jolivet a->s = 1.0; 915c7a4214aSPierre Jolivet a->bs[0] = 10; 916c7a4214aSPierre Jolivet a->bs[1] = 1000000; 917c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 918c7a4214aSPierre Jolivet a->eta = 10.0; 919c7a4214aSPierre Jolivet a->depth[0] = 0; 920c7a4214aSPierre Jolivet a->depth[1] = 0; 921c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 9229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool)); 9239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool)); 9249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense)); 9259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense)); 9269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool)); 9279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool)); 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool)); 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool)); 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool)); 931c7a4214aSPierre Jolivet PetscFunctionReturn(0); 932c7a4214aSPierre Jolivet } 933