1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/ 2*c9a4fd69SAudic XU #include <../src/mat/impls/shell/shell.h> 3c7a4214aSPierre Jolivet #include <petscblaslapack.h> 4c7a4214aSPierre Jolivet #include <set> 5c7a4214aSPierre Jolivet 6c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = {"sympartialACA", "fullACA", "SVD"}; 798e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = {"PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric"}; 8c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n" 9c7a4214aSPierre Jolivet " Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n" 10c7a4214aSPierre Jolivet " Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n" 11c7a4214aSPierre Jolivet " Year = {2020},\n" 12c7a4214aSPierre Jolivet " Publisher = {Elsevier},\n" 13c7a4214aSPierre Jolivet " Journal = {Numerische Mathematik},\n" 14c7a4214aSPierre Jolivet " Volume = {146},\n" 15c7a4214aSPierre Jolivet " Pages = {597--628},\n" 16c7a4214aSPierre Jolivet " Url = {https://github.com/htool-ddm/htool}\n" 17c7a4214aSPierre Jolivet "}\n"; 18c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE; 19c7a4214aSPierre Jolivet 20d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Htool(Mat A, Vec v) 21d71ae5a4SJacob Faibussowitsch { 22*c9a4fd69SAudic XU Mat_Htool *a; 23c7a4214aSPierre Jolivet PetscScalar *x; 24c7a4214aSPierre Jolivet PetscBool flg; 25c7a4214aSPierre Jolivet 26c7a4214aSPierre Jolivet PetscFunctionBegin; 279566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 285f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported"); 29*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 309566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(v, &x)); 31c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal(x); 329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(v, &x)); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34c7a4214aSPierre Jolivet } 35c7a4214aSPierre Jolivet 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b) 37d71ae5a4SJacob Faibussowitsch { 38*c9a4fd69SAudic XU Mat_Htool *a; 39c7a4214aSPierre Jolivet Mat B; 40c7a4214aSPierre Jolivet PetscScalar *ptr; 41c7a4214aSPierre Jolivet PetscBool flg; 42c7a4214aSPierre Jolivet 43c7a4214aSPierre Jolivet PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 455f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported"); 46*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 479566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */ 48c7a4214aSPierre Jolivet if (!B) { 49f22e26b7SPierre Jolivet PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B)); 509566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(B, &ptr)); 51c7a4214aSPierre Jolivet a->hmatrix->copy_local_diagonal_block(ptr); 529566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(B, &ptr)); 539566063dSJacob Faibussowitsch PetscCall(MatPropagateSymmetryOptions(A, B)); 549566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B)); 55c7a4214aSPierre Jolivet *b = B; 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 57*c9a4fd69SAudic XU } else { 58*c9a4fd69SAudic XU PetscCheck(PetscAbsScalar(((Mat_Shell *)A->data)->vscale - 1.0) <= PETSC_MACHINE_EPSILON && PetscAbsScalar(((Mat_Shell *)A->data)->vshift) <= PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot reuse DiagonalBlock with non-trivial shift and scaling"); // TODO FIXME 59*c9a4fd69SAudic XU *b = B; 60*c9a4fd69SAudic XU } 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c7a4214aSPierre Jolivet } 63c7a4214aSPierre Jolivet 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y) 65d71ae5a4SJacob Faibussowitsch { 66*c9a4fd69SAudic XU Mat_Htool *a; 67c7a4214aSPierre Jolivet const PetscScalar *in; 68c7a4214aSPierre Jolivet PetscScalar *out; 69c7a4214aSPierre Jolivet 70c7a4214aSPierre Jolivet PetscFunctionBegin; 71*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 74c7a4214aSPierre Jolivet a->hmatrix->mvprod_local_to_local(in, out); 759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78c7a4214aSPierre Jolivet } 79c7a4214aSPierre Jolivet 80d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y) 81d71ae5a4SJacob Faibussowitsch { 82*c9a4fd69SAudic XU Mat_Htool *a; 838b8ddd11SPierre Jolivet const PetscScalar *in; 848b8ddd11SPierre Jolivet PetscScalar *out; 858b8ddd11SPierre Jolivet 868b8ddd11SPierre Jolivet PetscFunctionBegin; 87*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 899566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 908b8ddd11SPierre Jolivet a->hmatrix->mvprod_transp_local_to_local(in, out); 919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 948b8ddd11SPierre Jolivet } 958b8ddd11SPierre Jolivet 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov) 97d71ae5a4SJacob Faibussowitsch { 98c7a4214aSPierre Jolivet std::set<PetscInt> set; 99c7a4214aSPierre Jolivet const PetscInt *idx; 1008f308287SPierre Jolivet PetscInt *oidx, size, bs[2]; 101c7a4214aSPierre Jolivet PetscMPIInt csize; 102c7a4214aSPierre Jolivet 103c7a4214aSPierre Jolivet PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, bs, bs + 1)); 1058f308287SPierre Jolivet if (bs[0] != bs[1]) bs[0] = 1; 106c7a4214aSPierre Jolivet for (PetscInt i = 0; i < is_max; ++i) { 107c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 108c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 1099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize)); 1105f80ce2aSJacob Faibussowitsch PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported parallel IS"); 1119566063dSJacob Faibussowitsch PetscCall(ISGetSize(is[i], &size)); 1129566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx)); 113c7a4214aSPierre Jolivet for (PetscInt j = 0; j < size; ++j) { 114c7a4214aSPierre Jolivet set.insert(idx[j]); 115c7a4214aSPierre Jolivet for (PetscInt k = 1; k <= ov; ++k) { /* for each layer of overlap */ 116c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 117c7a4214aSPierre 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 */ 118c7a4214aSPierre Jolivet } 119c7a4214aSPierre Jolivet } 1209566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[i], &idx)); 1219566063dSJacob Faibussowitsch PetscCall(ISDestroy(is + i)); 1228f308287SPierre Jolivet if (bs[0] > 1) { 1238f308287SPierre Jolivet for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) { 1248f308287SPierre Jolivet std::vector<PetscInt> block(bs[0]); 1258f308287SPierre Jolivet std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]); 1268f308287SPierre Jolivet set.insert(block.cbegin(), block.cend()); 1278f308287SPierre Jolivet } 1288f308287SPierre Jolivet } 129c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 1309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &oidx)); 131c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 132c7a4214aSPierre Jolivet oidx -= size; 1339566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i)); 134c7a4214aSPierre Jolivet } 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136c7a4214aSPierre Jolivet } 137c7a4214aSPierre Jolivet 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 139d71ae5a4SJacob Faibussowitsch { 140*c9a4fd69SAudic XU Mat_Htool *a; 141c7a4214aSPierre Jolivet Mat D, B, BT; 142c7a4214aSPierre Jolivet const PetscScalar *copy; 143c7a4214aSPierre Jolivet PetscScalar *ptr; 144c7a4214aSPierre Jolivet const PetscInt *idxr, *idxc, *it; 145c7a4214aSPierre Jolivet PetscInt nrow, m, i; 146c7a4214aSPierre Jolivet PetscBool flg; 147c7a4214aSPierre Jolivet 148c7a4214aSPierre Jolivet PetscFunctionBegin; 149*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME 150*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatAXPY() has been called on the input Mat"); // TODO FIXME 151*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->dshift, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatCreateSubMatrices() if MatDiagonalSet() has been called on the input Mat"); // TODO FIXME 152*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 15348a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 154c7a4214aSPierre Jolivet for (i = 0; i < n; ++i) { 1559566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &nrow)); 1569566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &m)); 1579566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i], &idxr)); 1589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i], &idxc)); 159f22e26b7SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i)); 1609566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr)); 161c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 1629566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 163c7a4214aSPierre Jolivet if (flg) { 1649566063dSJacob Faibussowitsch PetscCall(ISSorted(irow[i], &flg)); 165c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 166c7a4214aSPierre Jolivet it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart); 167c7a4214aSPierre Jolivet if (it != idxr + nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 168c7a4214aSPierre Jolivet if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 1699371c9d4SSatish Balay for (PetscInt j = 0; j < A->rmap->n && flg; ++j) 1709371c9d4SSatish Balay if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE; 171c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 172c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 173c7a4214aSPierre Jolivet * [ B C E ] 174c7a4214aSPierre Jolivet * A = [ B D E ] 175c7a4214aSPierre Jolivet * [ B F E ] 176c7a4214aSPierre Jolivet */ 177c7a4214aSPierre Jolivet m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */ 178*c9a4fd69SAudic XU PetscCall(MatGetDiagonalBlock(A, &D)); 1799566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(D, ©)); 1809371c9d4SSatish Balay for (PetscInt k = 0; k < A->rmap->n; ++k) { PetscCall(PetscArraycpy(ptr + (m + k) * nrow + m, copy + k * A->rmap->n, A->rmap->n)); /* block D from above */ } 1819566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(D, ©)); 182c7a4214aSPierre Jolivet if (m) { 183c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */ 184c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 185b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 1869566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B)); 1879566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 1889566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT)); 1899566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 190b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 1919566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 192c7a4214aSPierre Jolivet } else { 1937fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 1949566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 195c7a4214aSPierre Jolivet } 1969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 198c7a4214aSPierre Jolivet } else { 199c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */ 200c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow); 201c7a4214aSPierre Jolivet } 202c7a4214aSPierre Jolivet } 203c7a4214aSPierre Jolivet } 204c7a4214aSPierre Jolivet if (m + A->rmap->n != nrow) { 205c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, std::distance(it + A->rmap->n, idxr + nrow), idxr, idxc + m + A->rmap->n, ptr + (m + A->rmap->n) * nrow); /* vertical block E from above */ 206c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 207b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2089566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), ptr + (m + A->rmap->n) * nrow + m, &B)); 2099566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 2109566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, ptr + m * nrow + m + A->rmap->n, &BT)); 2119566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 212b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2139566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 214c7a4214aSPierre Jolivet } else { 2157fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 2169566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 217c7a4214aSPierre Jolivet } 2189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 220c7a4214aSPierre Jolivet } else { 221c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */ 222c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(std::distance(it + A->rmap->n, idxr + nrow), 1, it + A->rmap->n, idxc + m + k, ptr + (m + k) * nrow + m + A->rmap->n); 223c7a4214aSPierre Jolivet } 224c7a4214aSPierre Jolivet } 225c7a4214aSPierre Jolivet } 226c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 227c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 228c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 229c7a4214aSPierre Jolivet } /* unsorted IS */ 230c7a4214aSPierre Jolivet } 231c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 232c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */ 2339566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i], &idxr)); 2349566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i], &idxc)); 2359566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr)); 236*c9a4fd69SAudic XU PetscCall(MatScale((*submat)[i], ((Mat_Shell *)A->data)->vscale)); 237*c9a4fd69SAudic XU PetscCall(MatShift((*submat)[i], ((Mat_Shell *)A->data)->vshift)); 238c7a4214aSPierre Jolivet } 2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 240c7a4214aSPierre Jolivet } 241c7a4214aSPierre Jolivet 242d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A) 243d71ae5a4SJacob Faibussowitsch { 244*c9a4fd69SAudic XU Mat_Htool *a; 245c7a4214aSPierre Jolivet PetscContainer container; 246c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 247c7a4214aSPierre Jolivet 248c7a4214aSPierre Jolivet PetscFunctionBegin; 249*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 250f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr)); 251f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr)); 252f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr)); 253f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr)); 254f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr)); 255f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr)); 256f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr)); 257f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr)); 258f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container)); 260c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 2619566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 2639566063dSJacob Faibussowitsch PetscCall(PetscFree(kernelt)); 2649566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&container)); 265f22e26b7SPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr)); 266c7a4214aSPierre Jolivet } 26748a46eb9SPierre Jolivet if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source)); 2689566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFree2(a->work_source, a->work_target)); 270c7a4214aSPierre Jolivet delete a->wrapper; 271c7a4214aSPierre Jolivet delete a->hmatrix; 272*c9a4fd69SAudic XU PetscCall(PetscFree(a)); 273*c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable() 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275c7a4214aSPierre Jolivet } 276c7a4214aSPierre Jolivet 277d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv) 278d71ae5a4SJacob Faibussowitsch { 279*c9a4fd69SAudic XU Mat_Htool *a; 280c7a4214aSPierre Jolivet PetscBool flg; 281c7a4214aSPierre Jolivet 282c7a4214aSPierre Jolivet PetscFunctionBegin; 283*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 2849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg)); 285c7a4214aSPierre Jolivet if (flg) { 2869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->hmatrix->get_symmetry_type())); 287*c9a4fd69SAudic XU if (PetscAbsScalar(((Mat_Shell *)A->data)->vscale - 1.0) > PETSC_MACHINE_EPSILON) { 288c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 289*c9a4fd69SAudic XU PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(((Mat_Shell *)A->data)->vscale), (double)PetscImaginaryPart(((Mat_Shell *)A->data)->vscale))); 290c7a4214aSPierre Jolivet #else 291*c9a4fd69SAudic XU PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)((Mat_Shell *)A->data)->vscale)); 292*c9a4fd69SAudic XU #endif 293*c9a4fd69SAudic XU } 294*c9a4fd69SAudic XU if (PetscAbsScalar(((Mat_Shell *)A->data)->vshift) > PETSC_MACHINE_EPSILON) { 295*c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX) 296*c9a4fd69SAudic XU PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(((Mat_Shell *)A->data)->vshift), (double)PetscImaginaryPart(((Mat_Shell *)A->data)->vshift))); 297*c9a4fd69SAudic XU #else 298*c9a4fd69SAudic XU PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)((Mat_Shell *)A->data)->vshift)); 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 321c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 322d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 323d71ae5a4SJacob Faibussowitsch { 324*c9a4fd69SAudic XU Mat_Htool *a; 325c7a4214aSPierre Jolivet PetscInt *idxc; 326c7a4214aSPierre Jolivet PetscBLASInt one = 1, bn; 327c7a4214aSPierre Jolivet 328c7a4214aSPierre Jolivet PetscFunctionBegin; 329*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatGetRow() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME 330*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatGetRow() if MatAXPY() has been called on the input Mat"); // TODO FIXME 331*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 332c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 333c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, &idxc)); 335c7a4214aSPierre Jolivet for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i; 336c7a4214aSPierre Jolivet } 337c7a4214aSPierre Jolivet if (idx) *idx = idxc; 338c7a4214aSPierre Jolivet if (v) { 3399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, v)); 340c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 341cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 3429566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N, &bn)); 343*c9a4fd69SAudic XU PetscCallBLAS("BLASscal", BLASscal_(&bn, &(((Mat_Shell *)A->data)->vscale), *v, &one)); 344*c9a4fd69SAudic XU if (row < A->cmap->N) (*v)[row] += ((Mat_Shell *)A->data)->vshift; 345c7a4214aSPierre Jolivet } 34648a46eb9SPierre Jolivet if (!idx) PetscCall(PetscFree(idxc)); 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348c7a4214aSPierre Jolivet } 349c7a4214aSPierre Jolivet 350f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v) 351d71ae5a4SJacob Faibussowitsch { 352c7a4214aSPierre Jolivet PetscFunctionBegin; 35348a46eb9SPierre Jolivet if (idx) PetscCall(PetscFree(*idx)); 35448a46eb9SPierre Jolivet if (v) PetscCall(PetscFree(*v)); 3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 356c7a4214aSPierre Jolivet } 357c7a4214aSPierre Jolivet 358d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems *PetscOptionsObject) 359d71ae5a4SJacob Faibussowitsch { 360*c9a4fd69SAudic XU Mat_Htool *a; 361c7a4214aSPierre Jolivet PetscInt n; 362c7a4214aSPierre Jolivet PetscBool flg; 363c7a4214aSPierre Jolivet 364c7a4214aSPierre Jolivet PetscFunctionBegin; 365*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 366d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Htool options"); 367f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", nullptr, a->bs[0], a->bs, nullptr)); 368f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_max_block_size", "Maximal number of coefficients in a dense block", nullptr, a->bs[1], a->bs + 1, nullptr)); 369f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr)); 370f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr)); 371f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr)); 372f22e26b7SPierre Jolivet PetscCall(PetscOptionsInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", nullptr, a->depth[1], a->depth + 1, nullptr)); 373c7a4214aSPierre Jolivet n = 0; 374dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg)); 375c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 37698e73e17SPierre Jolivet n = 0; 377dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg)); 37898e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 379d0609cedSBarry Smith PetscOptionsHeadEnd(); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381c7a4214aSPierre Jolivet } 382c7a4214aSPierre Jolivet 383cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType) 384d71ae5a4SJacob Faibussowitsch { 385*c9a4fd69SAudic XU Mat_Htool *a; 386c7a4214aSPierre Jolivet const PetscInt *ranges; 387c7a4214aSPierre Jolivet PetscInt *offset; 388c7a4214aSPierre Jolivet PetscMPIInt size; 389b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U'; 390cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 39198e73e17SPierre Jolivet std::shared_ptr<htool::VirtualCluster> t, s = nullptr; 3923b9338faSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor = nullptr; 393c7a4214aSPierre Jolivet 394c7a4214aSPierre Jolivet PetscFunctionBegin; 3959566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite)); 396*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 397c7a4214aSPierre Jolivet delete a->wrapper; 398c7a4214aSPierre Jolivet delete a->hmatrix; 3999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 4009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &offset)); 4019566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 402c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 403c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 404c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 405c7a4214aSPierre Jolivet } 40698e73e17SPierre Jolivet switch (a->clustering) { 407d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 408d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 409d71ae5a4SJacob Faibussowitsch break; 410d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 411d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 412d71ae5a4SJacob Faibussowitsch break; 413d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 414d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 415d71ae5a4SJacob Faibussowitsch break; 416d71ae5a4SJacob Faibussowitsch default: 417d71ae5a4SJacob Faibussowitsch t = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 41898e73e17SPierre Jolivet } 419c7a4214aSPierre Jolivet t->set_minclustersize(a->bs[0]); 4200429413eSPierre Jolivet t->build(A->rmap->N, a->gcoords_target, offset, -1, PetscObjectComm((PetscObject)A)); 421c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 422c7a4214aSPierre Jolivet else { 423f22e26b7SPierre Jolivet a->wrapper = nullptr; 424cab00e0dSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx); 425c7a4214aSPierre Jolivet } 426c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4279566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A, &ranges)); 428c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 429c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 430c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 431c7a4214aSPierre Jolivet } 43298e73e17SPierre Jolivet switch (a->clustering) { 433d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 434d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 435d71ae5a4SJacob Faibussowitsch break; 436d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 437d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::GeometricSplitting>>>(a->dim); 438d71ae5a4SJacob Faibussowitsch break; 439d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 440d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::BoundingBox1<htool::SplittingTypes::RegularSplitting>>>(a->dim); 441d71ae5a4SJacob Faibussowitsch break; 442d71ae5a4SJacob Faibussowitsch default: 443d71ae5a4SJacob Faibussowitsch s = std::make_shared<htool::Cluster<htool::PCA<htool::SplittingTypes::RegularSplitting>>>(a->dim); 44498e73e17SPierre Jolivet } 445c7a4214aSPierre Jolivet s->set_minclustersize(a->bs[0]); 4460429413eSPierre Jolivet s->build(A->cmap->N, a->gcoords_source, offset, -1, PetscObjectComm((PetscObject)A)); 447c7a4214aSPierre Jolivet S = uplo = 'N'; 448c7a4214aSPierre Jolivet } 4499566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 450c7a4214aSPierre Jolivet switch (a->compressor) { 451d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_FULL_ACA: 452d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 453d71ae5a4SJacob Faibussowitsch break; 454d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_SVD: 455d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::SVD<PetscScalar>>(); 456d71ae5a4SJacob Faibussowitsch break; 457d71ae5a4SJacob Faibussowitsch default: 458d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 459c7a4214aSPierre Jolivet } 4600429413eSPierre 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))); 4613b9338faSPierre Jolivet a->hmatrix->set_compression(compressor); 462c7a4214aSPierre Jolivet a->hmatrix->set_maxblocksize(a->bs[1]); 463c7a4214aSPierre Jolivet a->hmatrix->set_mintargetdepth(a->depth[0]); 464c7a4214aSPierre Jolivet a->hmatrix->set_minsourcedepth(a->depth[1]); 4653b9338faSPierre Jolivet if (s) a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target, a->gcoords_source); 4663b9338faSPierre Jolivet else a->hmatrix->build(a->wrapper ? *a->wrapper : *generator, a->gcoords_target); 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 468c7a4214aSPierre Jolivet } 469c7a4214aSPierre Jolivet 470d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C) 471d71ae5a4SJacob Faibussowitsch { 472c7a4214aSPierre Jolivet Mat_Product *product = C->product; 473*c9a4fd69SAudic XU Mat_Htool *a; 474c7a4214aSPierre Jolivet const PetscScalar *in; 475c7a4214aSPierre Jolivet PetscScalar *out; 4768b8ddd11SPierre Jolivet PetscInt N, lda; 477c7a4214aSPierre Jolivet 478c7a4214aSPierre Jolivet PetscFunctionBegin; 479c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 480f22e26b7SPierre Jolivet PetscCall(MatGetSize(C, nullptr, &N)); 4819566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 4825f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 4839566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B, &in)); 4849566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &out)); 485*c9a4fd69SAudic XU PetscCall(MatShellGetContext(product->A, &a)); 4868b8ddd11SPierre Jolivet switch (product->type) { 487d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AB: 488d71ae5a4SJacob Faibussowitsch a->hmatrix->mvprod_local_to_local(in, out, N); 489d71ae5a4SJacob Faibussowitsch break; 490d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AtB: 491d71ae5a4SJacob Faibussowitsch a->hmatrix->mvprod_transp_local_to_local(in, out, N); 492d71ae5a4SJacob Faibussowitsch break; 493d71ae5a4SJacob Faibussowitsch default: 494d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]); 495c7a4214aSPierre Jolivet } 4969566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &out)); 4979566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B, &in)); 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 499c7a4214aSPierre Jolivet } 500c7a4214aSPierre Jolivet 501d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C) 502d71ae5a4SJacob Faibussowitsch { 503c7a4214aSPierre Jolivet Mat_Product *product = C->product; 504c7a4214aSPierre Jolivet Mat A, B; 505c7a4214aSPierre Jolivet PetscBool flg; 506c7a4214aSPierre Jolivet 507c7a4214aSPierre Jolivet PetscFunctionBegin; 508c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 509c7a4214aSPierre Jolivet A = product->A; 510c7a4214aSPierre Jolivet B = product->B; 5119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, "")); 51297713f8cSPierre 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); 51397713f8cSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 51497713f8cSPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N)); 51597713f8cSPierre Jolivet else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N)); 5168b8ddd11SPierre Jolivet } 5179566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 5189566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5199566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 5209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 522f22e26b7SPierre Jolivet C->ops->productsymbolic = nullptr; 523c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 5243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 525c7a4214aSPierre Jolivet } 526c7a4214aSPierre Jolivet 527d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 528d71ae5a4SJacob Faibussowitsch { 529c7a4214aSPierre Jolivet PetscFunctionBegin; 530c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 5318b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 533c7a4214aSPierre Jolivet } 534c7a4214aSPierre Jolivet 535d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::VirtualHMatrix<PetscScalar> **hmatrix) 536d71ae5a4SJacob Faibussowitsch { 537*c9a4fd69SAudic XU Mat_Htool *a; 538c7a4214aSPierre Jolivet 539c7a4214aSPierre Jolivet PetscFunctionBegin; 540*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 541c7a4214aSPierre Jolivet *hmatrix = a->hmatrix; 5423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 543c7a4214aSPierre Jolivet } 544c7a4214aSPierre Jolivet 545c7a4214aSPierre Jolivet /*@C 54611a5261eSBarry Smith MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`. 547c7a4214aSPierre Jolivet 548cc4c1da9SBarry Smith No Fortran Support, No C Support 549cc4c1da9SBarry Smith 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 5698434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 570d71ae5a4SJacob Faibussowitsch { 571*c9a4fd69SAudic XU Mat_Htool *a; 572c7a4214aSPierre Jolivet 573c7a4214aSPierre Jolivet PetscFunctionBegin; 574*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 575c7a4214aSPierre Jolivet a->kernel = kernel; 576c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 577c7a4214aSPierre Jolivet delete a->wrapper; 578c7a4214aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(A->rmap->N, A->cmap->N, a->dim, a->kernel, a->kernelctx); 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 580c7a4214aSPierre Jolivet } 581c7a4214aSPierre Jolivet 582c7a4214aSPierre Jolivet /*@C 58311a5261eSBarry Smith MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`. 584c7a4214aSPierre Jolivet 585cc4c1da9SBarry Smith Collective, No Fortran Support 586cc4c1da9SBarry Smith 587c7a4214aSPierre Jolivet Input Parameters: 588c7a4214aSPierre Jolivet + A - hierarchical matrix 5892ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 5902ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 591c7a4214aSPierre Jolivet 592c7a4214aSPierre Jolivet Level: advanced 593c7a4214aSPierre Jolivet 5941cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()` 595c7a4214aSPierre Jolivet @*/ 596cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 597d71ae5a4SJacob Faibussowitsch { 598c7a4214aSPierre Jolivet PetscFunctionBegin; 599c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 600c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 2); 6014f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 3); 6028434afd1SBarry Smith PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx)); 6033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 604c7a4214aSPierre Jolivet } 605c7a4214aSPierre Jolivet 606d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) 607d71ae5a4SJacob Faibussowitsch { 608*c9a4fd69SAudic XU Mat_Htool *a; 60998e73e17SPierre Jolivet std::vector<PetscInt> source; 61098e73e17SPierre Jolivet 61198e73e17SPierre Jolivet PetscFunctionBegin; 612*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 613a9918087SPierre Jolivet source = a->hmatrix->get_source_cluster()->get_local_perm(); 6149566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), source.size(), source.data(), PETSC_COPY_VALUES, is)); 6159566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61798e73e17SPierre Jolivet } 61898e73e17SPierre Jolivet 619cc4c1da9SBarry Smith /*@ 62011a5261eSBarry Smith MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix. 62198e73e17SPierre Jolivet 62298e73e17SPierre Jolivet Input Parameter: 62398e73e17SPierre Jolivet . A - hierarchical matrix 62498e73e17SPierre Jolivet 62598e73e17SPierre Jolivet Output Parameter: 62698e73e17SPierre Jolivet . is - permutation 62798e73e17SPierre Jolivet 62898e73e17SPierre Jolivet Level: advanced 62998e73e17SPierre Jolivet 6301cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 63198e73e17SPierre Jolivet @*/ 632cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) 633d71ae5a4SJacob Faibussowitsch { 63498e73e17SPierre Jolivet PetscFunctionBegin; 63598e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6364f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 637cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is)); 6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63998e73e17SPierre Jolivet } 64098e73e17SPierre Jolivet 641d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) 642d71ae5a4SJacob Faibussowitsch { 643*c9a4fd69SAudic XU Mat_Htool *a; 64498e73e17SPierre Jolivet std::vector<PetscInt> target; 64598e73e17SPierre Jolivet 64698e73e17SPierre Jolivet PetscFunctionBegin; 647*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 648a9918087SPierre Jolivet target = a->hmatrix->get_target_cluster()->get_local_perm(); 6499566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), target.size(), target.data(), PETSC_COPY_VALUES, is)); 6509566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65298e73e17SPierre Jolivet } 65398e73e17SPierre Jolivet 654cc4c1da9SBarry Smith /*@ 65511a5261eSBarry Smith MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix. 65698e73e17SPierre Jolivet 65798e73e17SPierre Jolivet Input Parameter: 65898e73e17SPierre Jolivet . A - hierarchical matrix 65998e73e17SPierre Jolivet 66098e73e17SPierre Jolivet Output Parameter: 66198e73e17SPierre Jolivet . is - permutation 66298e73e17SPierre Jolivet 66398e73e17SPierre Jolivet Level: advanced 66498e73e17SPierre Jolivet 6651cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 66698e73e17SPierre Jolivet @*/ 667cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) 668d71ae5a4SJacob Faibussowitsch { 66998e73e17SPierre Jolivet PetscFunctionBegin; 67098e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6714f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 672cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is)); 6733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67498e73e17SPierre Jolivet } 67598e73e17SPierre Jolivet 676d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) 677d71ae5a4SJacob Faibussowitsch { 678*c9a4fd69SAudic XU Mat_Htool *a; 67998e73e17SPierre Jolivet 68098e73e17SPierre Jolivet PetscFunctionBegin; 681*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 68298e73e17SPierre Jolivet a->hmatrix->set_use_permutation(use); 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68498e73e17SPierre Jolivet } 68598e73e17SPierre Jolivet 686cc4c1da9SBarry Smith /*@ 68711a5261eSBarry Smith MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation. 68898e73e17SPierre Jolivet 68998e73e17SPierre Jolivet Input Parameters: 69098e73e17SPierre Jolivet + A - hierarchical matrix 69198e73e17SPierre Jolivet - use - Boolean value 69298e73e17SPierre Jolivet 69398e73e17SPierre Jolivet Level: advanced 69498e73e17SPierre Jolivet 6951cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 69698e73e17SPierre Jolivet @*/ 697cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) 698d71ae5a4SJacob Faibussowitsch { 69998e73e17SPierre Jolivet PetscFunctionBegin; 70098e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 70198e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 702cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use)); 7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70498e73e17SPierre Jolivet } 70598e73e17SPierre Jolivet 706cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 707d71ae5a4SJacob Faibussowitsch { 708c7a4214aSPierre Jolivet Mat C; 709*c9a4fd69SAudic XU Mat_Htool *a; 710c7a4214aSPierre Jolivet PetscInt lda; 711c7a4214aSPierre Jolivet PetscScalar *array; 712c7a4214aSPierre Jolivet 713c7a4214aSPierre Jolivet PetscFunctionBegin; 714*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 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)); 732*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatConvert() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME 733*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatConvert() if MatAXPY() has been called on the input Mat"); // TODO FIXME 734*c9a4fd69SAudic XU PetscCall(MatScale(C, ((Mat_Shell *)A->data)->vscale)); 735*c9a4fd69SAudic XU PetscCall(MatShift(C, ((Mat_Shell *)A->data)->vshift)); 736c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7379566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &C)); 738c7a4214aSPierre Jolivet } else *B = C; 7393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 740c7a4214aSPierre Jolivet } 741c7a4214aSPierre Jolivet 742d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx) 743d71ae5a4SJacob Faibussowitsch { 744c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx; 74598e73e17SPierre Jolivet PetscScalar *tmp; 74698e73e17SPierre Jolivet 74798e73e17SPierre Jolivet PetscFunctionBegin; 7483ba16761SJacob Faibussowitsch PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx)); 7499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * N, &tmp)); 7509566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, ptr, M * N)); 75198e73e17SPierre Jolivet for (PetscInt i = 0; i < M; ++i) { 75298e73e17SPierre Jolivet for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N]; 75398e73e17SPierre Jolivet } 7549566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 756c7a4214aSPierre Jolivet } 757c7a4214aSPierre Jolivet 758c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 759d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) 760d71ae5a4SJacob Faibussowitsch { 761c7a4214aSPierre Jolivet Mat C; 762*c9a4fd69SAudic XU Mat_Htool *a, *c; 763c7a4214aSPierre Jolivet PetscInt M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n; 764c7a4214aSPierre Jolivet PetscContainer container; 765c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 766c7a4214aSPierre Jolivet 767c7a4214aSPierre Jolivet PetscFunctionBegin; 768*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 7697fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 7705f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported"); 771c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 7729566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, n, m, N, M)); 7749566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 7759566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7769566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)C), &container)); 7779566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 7789566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(container, kernelt)); 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)C, "KernelTranspose", (PetscObject)container)); 780c7a4214aSPierre Jolivet } else { 781c7a4214aSPierre Jolivet C = *B; 7829566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container)); 7835f80ce2aSJacob Faibussowitsch PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 7849566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 785c7a4214aSPierre Jolivet } 786*c9a4fd69SAudic XU PetscCall(MatShellGetContext(C, &c)); 787c7a4214aSPierre Jolivet c->dim = a->dim; 788*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->zrows && !((Mat_Shell *)A->data)->zcols, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatTranspose() if MatZeroRows() or MatZeroRowsColumns() has been called on the input Mat"); // TODO FIXME 789*c9a4fd69SAudic XU PetscCheck(!((Mat_Shell *)A->data)->axpy, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot call MatTranspose() if MatAXPY() has been called on the input Mat"); // TODO FIXME 790*c9a4fd69SAudic XU ((Mat_Shell *)C->data)->vscale = ((Mat_Shell *)A->data)->vscale; 791*c9a4fd69SAudic XU ((Mat_Shell *)C->data)->vshift = ((Mat_Shell *)A->data)->vshift; 79298e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 793c7a4214aSPierre Jolivet if (kernelt->A != A) { 7949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 795c7a4214aSPierre Jolivet kernelt->A = A; 7969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 797c7a4214aSPierre Jolivet } 798c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 799c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 800c7a4214aSPierre Jolivet c->kernelctx = kernelt; 801c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target)); 8039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim)); 804c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 8059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source)); 8069566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim)); 807c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 8089566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target)); 809c7a4214aSPierre Jolivet } 8109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 8119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 812c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 8133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 814c7a4214aSPierre Jolivet } 815c7a4214aSPierre Jolivet 816c7a4214aSPierre Jolivet /*@C 81711a5261eSBarry Smith MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel. 818c7a4214aSPierre Jolivet 819cc4c1da9SBarry Smith Collective, No Fortran Support 820cc4c1da9SBarry Smith 821c7a4214aSPierre Jolivet Input Parameters: 822c7a4214aSPierre Jolivet + comm - MPI communicator 8232ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 8242ef1f0ffSBarry Smith . n - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given) 8252ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 8262ef1f0ffSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 827c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 828c7a4214aSPierre Jolivet . coords_target - coordinates of the target 829c7a4214aSPierre Jolivet . coords_source - coordinates of the source 8302ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 8312ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 832c7a4214aSPierre Jolivet 833c7a4214aSPierre Jolivet Output Parameter: 834c7a4214aSPierre Jolivet . B - matrix 835c7a4214aSPierre Jolivet 836c7a4214aSPierre Jolivet Options Database Keys: 83711a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`> - minimal leaf size in cluster tree 83811a5261eSBarry Smith . -mat_htool_max_block_size <`PetscInt`> - maximal number of coefficients in a dense block 83911a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`> - relative error in Frobenius norm when approximating a block 84011a5261eSBarry Smith . -mat_htool_eta <`PetscReal`> - admissibility condition tolerance 84111a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`> - minimal cluster tree depth associated with the rows 84211a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`> - minimal cluster tree depth associated with the columns 84398e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 84498e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 845c7a4214aSPierre Jolivet 846c7a4214aSPierre Jolivet Level: intermediate 847c7a4214aSPierre Jolivet 8481cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 849c7a4214aSPierre Jolivet @*/ 8508434afd1SBarry Smith PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords_target[], const PetscReal coords_source[], MatHtoolKernelFn *kernel, void *kernelctx, Mat *B) 851d71ae5a4SJacob Faibussowitsch { 852c7a4214aSPierre Jolivet Mat A; 853c7a4214aSPierre Jolivet Mat_Htool *a; 854c7a4214aSPierre Jolivet 855c7a4214aSPierre Jolivet PetscFunctionBegin; 8569566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 857c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A, spacedim, 6); 8584f572ea9SToby Isaac PetscAssertPointer(coords_target, 7); 8594f572ea9SToby Isaac PetscAssertPointer(coords_source, 8); 860c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 9); 8614f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 10); 8629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, M, N)); 8639566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATHTOOL)); 8649566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 865*c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 866c7a4214aSPierre Jolivet a->dim = spacedim; 867c7a4214aSPierre Jolivet a->kernel = kernel; 868c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 8699566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target)); 8709566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim)); 8711c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */ 872c7a4214aSPierre Jolivet if (coords_target != coords_source) { 8739566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source)); 8749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim)); 8751c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */ 876c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 8779566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target)); 878c7a4214aSPierre Jolivet *B = A; 8793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 880c7a4214aSPierre Jolivet } 881c7a4214aSPierre Jolivet 882c7a4214aSPierre Jolivet /*MC 883c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 884c7a4214aSPierre Jolivet 8852ef1f0ffSBarry Smith Use `./configure --download-htool` to install PETSc to use Htool. 886c7a4214aSPierre Jolivet 8872ef1f0ffSBarry Smith Options Database Key: 8882ef1f0ffSBarry Smith . -mat_type htool - matrix type to `MATHTOOL` 889c7a4214aSPierre Jolivet 890c7a4214aSPierre Jolivet Level: beginner 891c7a4214aSPierre Jolivet 8921cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 893c7a4214aSPierre Jolivet M*/ 894d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 895d71ae5a4SJacob Faibussowitsch { 896c7a4214aSPierre Jolivet Mat_Htool *a; 897c7a4214aSPierre Jolivet 898c7a4214aSPierre Jolivet PetscFunctionBegin; 899*c9a4fd69SAudic XU PetscCall(MatSetType(A, MATSHELL)); 9004dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 901*c9a4fd69SAudic XU PetscCall(MatShellSetContext(A, a)); 902*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Htool)); 903*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Htool)); 904*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Htool)); 905*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool)); 906*c9a4fd69SAudic XU if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool)); 907c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 908c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 909*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_Htool)); 910*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_Htool)); 911*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (void (*)(void))MatGetRow_Htool)); 912*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (void (*)(void))MatRestoreRow_Htool)); 913*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_Htool)); 914*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_Htool)); 915*c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_Htool)); 916c7a4214aSPierre Jolivet a->dim = 0; 917f22e26b7SPierre Jolivet a->gcoords_target = nullptr; 918f22e26b7SPierre Jolivet a->gcoords_source = nullptr; 919c7a4214aSPierre Jolivet a->bs[0] = 10; 920c7a4214aSPierre Jolivet a->bs[1] = 1000000; 921c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 922c7a4214aSPierre Jolivet a->eta = 10.0; 923c7a4214aSPierre Jolivet a->depth[0] = 0; 924c7a4214aSPierre Jolivet a->depth[1] = 0; 925c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 9269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool)); 9279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool)); 9289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense)); 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense)); 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool)); 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool)); 9329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool)); 9339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool)); 9349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool)); 935*c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable)); 936*c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 937*c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 938*c9a4fd69SAudic XU PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL)); 9393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 940c7a4214aSPierre Jolivet } 941