1c7a4214aSPierre Jolivet #include <../src/mat/impls/htool/htool.hpp> /*I "petscmat.h" I*/ 2c7a4214aSPierre Jolivet #include <set> 3c7a4214aSPierre Jolivet 4c7a4214aSPierre Jolivet const char *const MatHtoolCompressorTypes[] = {"sympartialACA", "fullACA", "SVD"}; 598e73e17SPierre Jolivet const char *const MatHtoolClusteringTypes[] = {"PCARegular", "PCAGeometric", "BoundingBox1Regular", "BoundingBox1Geometric"}; 6c7a4214aSPierre Jolivet const char HtoolCitation[] = "@article{marchand2020two,\n" 7c7a4214aSPierre Jolivet " Author = {Marchand, Pierre and Claeys, Xavier and Jolivet, Pierre and Nataf, Fr\\'ed\\'eric and Tournier, Pierre-Henri},\n" 8c7a4214aSPierre Jolivet " Title = {Two-level preconditioning for $h$-version boundary element approximation of hypersingular operator with {GenEO}},\n" 9c7a4214aSPierre Jolivet " Year = {2020},\n" 10c7a4214aSPierre Jolivet " Publisher = {Elsevier},\n" 11c7a4214aSPierre Jolivet " Journal = {Numerische Mathematik},\n" 12c7a4214aSPierre Jolivet " Volume = {146},\n" 13c7a4214aSPierre Jolivet " Pages = {597--628},\n" 14c7a4214aSPierre Jolivet " Url = {https://github.com/htool-ddm/htool}\n" 15c7a4214aSPierre Jolivet "}\n"; 16c7a4214aSPierre Jolivet static PetscBool HtoolCite = PETSC_FALSE; 17c7a4214aSPierre Jolivet 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Htool(Mat A, Vec v) 19d71ae5a4SJacob Faibussowitsch { 20c9a4fd69SAudic XU Mat_Htool *a; 21c7a4214aSPierre Jolivet PetscScalar *x; 22c7a4214aSPierre Jolivet PetscBool flg; 23c7a4214aSPierre Jolivet 24c7a4214aSPierre Jolivet PetscFunctionBegin; 259566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 265f80ce2aSJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only congruent layouts supported"); 27c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 289566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(v, &x)); 291dd4f53aSPierre Jolivet PetscStackCallExternalVoid("copy_diagonal_in_user_numbering", htool::copy_diagonal_in_user_numbering(a->distributed_operator_holder->hmatrix, x)); 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(v, &x)); 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32c7a4214aSPierre Jolivet } 33c7a4214aSPierre Jolivet 34d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_Htool(Mat A, Mat *b) 35d71ae5a4SJacob Faibussowitsch { 36c9a4fd69SAudic XU Mat_Htool *a; 37c7a4214aSPierre Jolivet Mat B; 381dd4f53aSPierre Jolivet PetscScalar *ptr, shift, scale; 39c7a4214aSPierre Jolivet PetscBool flg; 401dd4f53aSPierre Jolivet PetscMPIInt rank; 411dd4f53aSPierre Jolivet htool::Cluster<PetscReal> *source_cluster = nullptr; 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"); 46c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 479566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "DiagonalBlock", (PetscObject *)&B)); /* same logic as in MatGetDiagonalBlock_MPIDense() */ 48c7a4214aSPierre Jolivet if (!B) { 4966b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 50f22e26b7SPierre Jolivet PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, A->rmap->n, A->rmap->n, A->rmap->n, nullptr, &B)); 519566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(B, &ptr)); 521dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 531dd4f53aSPierre Jolivet source_cluster = a->source_cluster ? a->source_cluster.get() : a->target_cluster.get(); 541dd4f53aSPierre Jolivet PetscStackCallExternalVoid("copy_to_dense_in_user_numbering", htool::copy_to_dense_in_user_numbering(*a->distributed_operator_holder->hmatrix.get_sub_hmatrix(a->target_cluster->get_cluster_on_partition(rank), source_cluster->get_cluster_on_partition(rank)), ptr)); 559566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(B, &ptr)); 569566063dSJacob Faibussowitsch PetscCall(MatPropagateSymmetryOptions(A, B)); 579566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "DiagonalBlock", (PetscObject)B)); 58c7a4214aSPierre Jolivet *b = B; 599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 6066b4df27SPierre Jolivet PetscCall(MatShift(*b, shift)); 6166b4df27SPierre Jolivet PetscCall(MatScale(*b, scale)); 62c9a4fd69SAudic XU } else { 6366b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, (PetscScalar *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 64c9a4fd69SAudic XU *b = B; 65c9a4fd69SAudic XU } 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67c7a4214aSPierre Jolivet } 68c7a4214aSPierre Jolivet 69d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Htool(Mat A, Vec x, Vec y) 70d71ae5a4SJacob Faibussowitsch { 71c9a4fd69SAudic XU Mat_Htool *a; 72c7a4214aSPierre Jolivet const PetscScalar *in; 73c7a4214aSPierre Jolivet PetscScalar *out; 74c7a4214aSPierre Jolivet 75c7a4214aSPierre Jolivet PetscFunctionBegin; 76c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 789566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 79c095613aSPierre Jolivet if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_vector_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr); 80c095613aSPierre Jolivet else htool::internal_add_distributed_operator_vector_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84c7a4214aSPierre Jolivet } 85c7a4214aSPierre Jolivet 86d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y) 87d71ae5a4SJacob Faibussowitsch { 88c9a4fd69SAudic XU Mat_Htool *a; 898b8ddd11SPierre Jolivet const PetscScalar *in; 908b8ddd11SPierre Jolivet PetscScalar *out; 918b8ddd11SPierre Jolivet 928b8ddd11SPierre Jolivet PetscFunctionBegin; 93c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 959566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 96c095613aSPierre Jolivet if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_vector_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr); 97c095613aSPierre Jolivet else htool::internal_add_distributed_operator_vector_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, nullptr); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1018b8ddd11SPierre Jolivet } 1028b8ddd11SPierre Jolivet 103d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov) 104d71ae5a4SJacob Faibussowitsch { 105c7a4214aSPierre Jolivet std::set<PetscInt> set; 106c7a4214aSPierre Jolivet const PetscInt *idx; 1078f308287SPierre Jolivet PetscInt *oidx, size, bs[2]; 108c7a4214aSPierre Jolivet PetscMPIInt csize; 109c7a4214aSPierre Jolivet 110c7a4214aSPierre Jolivet PetscFunctionBegin; 1119566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, bs, bs + 1)); 1128f308287SPierre Jolivet if (bs[0] != bs[1]) bs[0] = 1; 113c7a4214aSPierre Jolivet for (PetscInt i = 0; i < is_max; ++i) { 114c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 115c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 1169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize)); 1171dd4f53aSPierre Jolivet PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel IS"); 1189566063dSJacob Faibussowitsch PetscCall(ISGetSize(is[i], &size)); 1199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx)); 120c7a4214aSPierre Jolivet for (PetscInt j = 0; j < size; ++j) { 121c7a4214aSPierre Jolivet set.insert(idx[j]); 122c7a4214aSPierre Jolivet for (PetscInt k = 1; k <= ov; ++k) { /* for each layer of overlap */ 123c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 124c7a4214aSPierre 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 */ 125c7a4214aSPierre Jolivet } 126c7a4214aSPierre Jolivet } 1279566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[i], &idx)); 1289566063dSJacob Faibussowitsch PetscCall(ISDestroy(is + i)); 1298f308287SPierre Jolivet if (bs[0] > 1) { 1308f308287SPierre Jolivet for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) { 1318f308287SPierre Jolivet std::vector<PetscInt> block(bs[0]); 1328f308287SPierre Jolivet std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]); 1338f308287SPierre Jolivet set.insert(block.cbegin(), block.cend()); 1348f308287SPierre Jolivet } 1358f308287SPierre Jolivet } 136c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &oidx)); 138c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 139c7a4214aSPierre Jolivet oidx -= size; 1409566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i)); 141c7a4214aSPierre Jolivet } 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143c7a4214aSPierre Jolivet } 144c7a4214aSPierre Jolivet 145d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 146d71ae5a4SJacob Faibussowitsch { 147c9a4fd69SAudic XU Mat_Htool *a; 148c7a4214aSPierre Jolivet Mat D, B, BT; 149c7a4214aSPierre Jolivet const PetscScalar *copy; 15066b4df27SPierre Jolivet PetscScalar *ptr, shift, scale; 151c7a4214aSPierre Jolivet const PetscInt *idxr, *idxc, *it; 152c7a4214aSPierre Jolivet PetscInt nrow, m, i; 153c7a4214aSPierre Jolivet PetscBool flg; 154c7a4214aSPierre Jolivet 155c7a4214aSPierre Jolivet PetscFunctionBegin; 15666b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 157c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 15848a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 159c7a4214aSPierre Jolivet for (i = 0; i < n; ++i) { 1609566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &nrow)); 1619566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &m)); 1629566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i], &idxr)); 1639566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i], &idxc)); 164f22e26b7SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i)); 1659566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr)); 166c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 1679566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 168c7a4214aSPierre Jolivet if (flg) { 1699566063dSJacob Faibussowitsch PetscCall(ISSorted(irow[i], &flg)); 170c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 171c7a4214aSPierre Jolivet it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart); 172c7a4214aSPierre Jolivet if (it != idxr + nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 173c7a4214aSPierre Jolivet if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 1749371c9d4SSatish Balay for (PetscInt j = 0; j < A->rmap->n && flg; ++j) 1759371c9d4SSatish Balay if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE; 176c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 177c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 178c7a4214aSPierre Jolivet * [ B C E ] 179c7a4214aSPierre Jolivet * A = [ B D E ] 180c7a4214aSPierre Jolivet * [ B F E ] 181c7a4214aSPierre Jolivet */ 182c7a4214aSPierre Jolivet m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */ 183c9a4fd69SAudic XU PetscCall(MatGetDiagonalBlock(A, &D)); 1849566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(D, ©)); 185ac530a7eSPierre Jolivet 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 */ 1869566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(D, ©)); 187c7a4214aSPierre Jolivet if (m) { 188c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */ 189c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 190b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 1919566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B)); 1929566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 1939566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT)); 1949566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 195b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 1969566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 197c7a4214aSPierre Jolivet } else { 1987fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 1999566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 200c7a4214aSPierre Jolivet } 2019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 203c7a4214aSPierre Jolivet } else { 204c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */ 205c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow); 206c7a4214aSPierre Jolivet } 207c7a4214aSPierre Jolivet } 208c7a4214aSPierre Jolivet } 209c7a4214aSPierre Jolivet if (m + A->rmap->n != nrow) { 210c7a4214aSPierre 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 */ 211c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 212b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2139566063dSJacob 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)); 2149566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 2159566063dSJacob 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)); 2169566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 217b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2189566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 219c7a4214aSPierre Jolivet } else { 2207fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 2219566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 222c7a4214aSPierre Jolivet } 2239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 225c7a4214aSPierre Jolivet } else { 226c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */ 227c7a4214aSPierre 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); 228c7a4214aSPierre Jolivet } 229c7a4214aSPierre Jolivet } 230c7a4214aSPierre Jolivet } 231c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 232c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 233c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 234c7a4214aSPierre Jolivet } /* unsorted IS */ 235c7a4214aSPierre Jolivet } 236c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 237c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */ 2389566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i], &idxr)); 2399566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i], &idxc)); 2409566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr)); 24166b4df27SPierre Jolivet PetscCall(MatShift((*submat)[i], shift)); 24266b4df27SPierre Jolivet PetscCall(MatScale((*submat)[i], scale)); 243c7a4214aSPierre Jolivet } 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245c7a4214aSPierre Jolivet } 246c7a4214aSPierre Jolivet 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A) 248d71ae5a4SJacob Faibussowitsch { 249c9a4fd69SAudic XU Mat_Htool *a; 250c7a4214aSPierre Jolivet PetscContainer container; 251c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 252c7a4214aSPierre Jolivet 253c7a4214aSPierre Jolivet PetscFunctionBegin; 254c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 255f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr)); 256f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr)); 257f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr)); 258f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr)); 259f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr)); 260f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr)); 261f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr)); 262f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr)); 263f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr)); 264c095613aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUseRecompression_C", nullptr)); 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container)); 266c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 267*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, &kernelt)); 2689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 269f22e26b7SPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr)); 270c7a4214aSPierre Jolivet } 27148a46eb9SPierre Jolivet if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 273c7a4214aSPierre Jolivet delete a->wrapper; 2741dd4f53aSPierre Jolivet a->target_cluster.reset(); 2751dd4f53aSPierre Jolivet a->source_cluster.reset(); 2761dd4f53aSPierre Jolivet a->distributed_operator_holder.reset(); 277c9a4fd69SAudic XU PetscCall(PetscFree(a)); 278c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable() 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 280c7a4214aSPierre Jolivet } 281c7a4214aSPierre Jolivet 282d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv) 283d71ae5a4SJacob Faibussowitsch { 284c9a4fd69SAudic XU Mat_Htool *a; 28566b4df27SPierre Jolivet PetscScalar shift, scale; 286c7a4214aSPierre Jolivet PetscBool flg; 2871dd4f53aSPierre Jolivet std::map<std::string, std::string> hmatrix_information; 288c7a4214aSPierre Jolivet 289c7a4214aSPierre Jolivet PetscFunctionBegin; 290c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 2911dd4f53aSPierre Jolivet hmatrix_information = htool::get_distributed_hmatrix_information(a->distributed_operator_holder->hmatrix, PetscObjectComm((PetscObject)A)); 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg)); 293c7a4214aSPierre Jolivet if (flg) { 29466b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 295c095613aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->distributed_operator_holder->block_diagonal_hmatrix->get_symmetry())); 29666b4df27SPierre Jolivet if (PetscAbsScalar(scale - 1.0) > PETSC_MACHINE_EPSILON) { 297c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 29866b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(scale), (double)PetscImaginaryPart(scale))); 299c7a4214aSPierre Jolivet #else 30066b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)scale)); 301c9a4fd69SAudic XU #endif 302c9a4fd69SAudic XU } 30366b4df27SPierre Jolivet if (PetscAbsScalar(shift) > PETSC_MACHINE_EPSILON) { 304c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX) 30566b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(shift), (double)PetscImaginaryPart(shift))); 306c9a4fd69SAudic XU #else 30766b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)shift)); 308c7a4214aSPierre Jolivet #endif 309c7a4214aSPierre Jolivet } 310c095613aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "maximal cluster leaf size: %" PetscInt_FMT "\n", a->max_cluster_leaf_size)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon)); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta)); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0])); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1])); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor])); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering])); 3171dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", hmatrix_information["Compression_ratio"].c_str())); 3181dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", hmatrix_information["Space_saving"].c_str())); 3191dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "block tree consistency: %s\n", PetscBools[a->distributed_operator_holder->hmatrix.is_block_tree_consistent()])); 320c095613aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "recompression: %s\n", PetscBools[a->recompression])); 3211dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "number of dense (resp. low rank) matrices: %s (resp. %s)\n", hmatrix_information["Number_of_dense_blocks"].c_str(), hmatrix_information["Number_of_low_rank_blocks"].c_str())); 3221dd4f53aSPierre Jolivet PetscCall( 3231dd4f53aSPierre Jolivet PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) dense block sizes: (%s, %s, %s)\n", hmatrix_information["Dense_block_size_min"].c_str(), hmatrix_information["Dense_block_size_mean"].c_str(), hmatrix_information["Dense_block_size_max"].c_str())); 3241dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) low rank block sizes: (%s, %s, %s)\n", hmatrix_information["Low_rank_block_size_min"].c_str(), hmatrix_information["Low_rank_block_size_mean"].c_str(), 3251dd4f53aSPierre Jolivet hmatrix_information["Low_rank_block_size_max"].c_str())); 3261dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "(minimum, mean, maximum) ranks: (%s, %s, %s)\n", hmatrix_information["Rank_min"].c_str(), hmatrix_information["Rank_mean"].c_str(), hmatrix_information["Rank_max"].c_str())); 327c7a4214aSPierre Jolivet } 3283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 329c7a4214aSPierre Jolivet } 330c7a4214aSPierre Jolivet 331c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 332d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 333d71ae5a4SJacob Faibussowitsch { 334c9a4fd69SAudic XU Mat_Htool *a; 33566b4df27SPierre Jolivet PetscScalar shift, scale; 336c7a4214aSPierre Jolivet PetscInt *idxc; 337c7a4214aSPierre Jolivet PetscBLASInt one = 1, bn; 338c7a4214aSPierre Jolivet 339c7a4214aSPierre Jolivet PetscFunctionBegin; 34066b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 341c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 342c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 343c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, &idxc)); 345c7a4214aSPierre Jolivet for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i; 346c7a4214aSPierre Jolivet } 347c7a4214aSPierre Jolivet if (idx) *idx = idxc; 348c7a4214aSPierre Jolivet if (v) { 3499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, v)); 350c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 351cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 3529566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N, &bn)); 353ebd42cecSPierre Jolivet PetscCallCXX(htool::Blas<PetscScalar>::scal(&bn, &scale, *v, &one)); 35466b4df27SPierre Jolivet if (row < A->cmap->N) (*v)[row] += shift; 355c7a4214aSPierre Jolivet } 35648a46eb9SPierre Jolivet if (!idx) PetscCall(PetscFree(idxc)); 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 358c7a4214aSPierre Jolivet } 359c7a4214aSPierre Jolivet 360f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v) 361d71ae5a4SJacob Faibussowitsch { 362c7a4214aSPierre Jolivet PetscFunctionBegin; 36348a46eb9SPierre Jolivet if (idx) PetscCall(PetscFree(*idx)); 36448a46eb9SPierre Jolivet if (v) PetscCall(PetscFree(*v)); 3653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 366c7a4214aSPierre Jolivet } 367c7a4214aSPierre Jolivet 368ce78bad3SBarry Smith static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems PetscOptionsObject) 369d71ae5a4SJacob Faibussowitsch { 370c9a4fd69SAudic XU Mat_Htool *a; 371c7a4214aSPierre Jolivet PetscInt n; 372c7a4214aSPierre Jolivet PetscBool flg; 373c7a4214aSPierre Jolivet 374c7a4214aSPierre Jolivet PetscFunctionBegin; 375c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 376d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Htool options"); 377c095613aSPierre Jolivet PetscCall(PetscOptionsBoundedInt("-mat_htool_max_cluster_leaf_size", "Maximal leaf size in cluster tree", nullptr, a->max_cluster_leaf_size, &a->max_cluster_leaf_size, nullptr, 0)); 3781dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr, 0.0)); 379f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr)); 3801dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr, 0)); 3811dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedInt("-mat_htool_min_source_depth", "Minimal cluster tree depth associated with the columns", nullptr, a->depth[1], a->depth + 1, nullptr, 0)); 3821dd4f53aSPierre Jolivet PetscCall(PetscOptionsBool("-mat_htool_block_tree_consistency", "Block tree consistency", nullptr, a->block_tree_consistency, &a->block_tree_consistency, nullptr)); 383c095613aSPierre Jolivet PetscCall(PetscOptionsBool("-mat_htool_recompression", "Use recompression", nullptr, a->recompression, &a->recompression, nullptr)); 3841dd4f53aSPierre Jolivet 385c7a4214aSPierre Jolivet n = 0; 386dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg)); 387c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 38898e73e17SPierre Jolivet n = 0; 389dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg)); 39098e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 391d0609cedSBarry Smith PetscOptionsHeadEnd(); 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 393c7a4214aSPierre Jolivet } 394c7a4214aSPierre Jolivet 395cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType) 396d71ae5a4SJacob Faibussowitsch { 397c9a4fd69SAudic XU Mat_Htool *a; 398c7a4214aSPierre Jolivet const PetscInt *ranges; 399c7a4214aSPierre Jolivet PetscInt *offset; 4001dd4f53aSPierre Jolivet PetscMPIInt size, rank; 401b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U'; 402cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 4031dd4f53aSPierre Jolivet htool::ClusterTreeBuilder<PetscReal> recursive_build_strategy; 4041dd4f53aSPierre Jolivet htool::Cluster<PetscReal> *source_cluster; 405c095613aSPierre Jolivet std::shared_ptr<htool::VirtualInternalLowRankGenerator<PetscScalar>> compressor; 406c7a4214aSPierre Jolivet 407c7a4214aSPierre Jolivet PetscFunctionBegin; 4089566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite)); 409c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 410c7a4214aSPierre Jolivet delete a->wrapper; 4111dd4f53aSPierre Jolivet a->target_cluster.reset(); 4121dd4f53aSPierre Jolivet a->source_cluster.reset(); 4131dd4f53aSPierre Jolivet a->distributed_operator_holder.reset(); 4141dd4f53aSPierre Jolivet // clustering 4159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 4169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &offset)); 4179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 418c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 419c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 420c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 421c7a4214aSPierre Jolivet } 42298e73e17SPierre Jolivet switch (a->clustering) { 423d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 424c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::GeometricSplitting<PetscReal>>>()); 425d71ae5a4SJacob Faibussowitsch break; 426d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 427c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::GeometricSplitting<PetscReal>>>()); 428d71ae5a4SJacob Faibussowitsch break; 429d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 430c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::RegularSplitting<PetscReal>>>()); 431d71ae5a4SJacob Faibussowitsch break; 432d71ae5a4SJacob Faibussowitsch default: 433c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::RegularSplitting<PetscReal>>>()); 43498e73e17SPierre Jolivet } 435c095613aSPierre Jolivet recursive_build_strategy.set_maximal_leaf_size(a->max_cluster_leaf_size); 436c095613aSPierre Jolivet a->target_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree_from_local_partition(A->rmap->N, a->dim, a->gcoords_target, 2, size, offset)); 437c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 4389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangesColumn(A, &ranges)); 439c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 440c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 441c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 442c7a4214aSPierre Jolivet } 44398e73e17SPierre Jolivet switch (a->clustering) { 444d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 445c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::GeometricSplitting<PetscReal>>>()); 446d71ae5a4SJacob Faibussowitsch break; 447d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 448c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::GeometricSplitting<PetscReal>>>()); 449d71ae5a4SJacob Faibussowitsch break; 450d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 451c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeBoundingBox<PetscReal>, htool::RegularSplitting<PetscReal>>>()); 452d71ae5a4SJacob Faibussowitsch break; 453d71ae5a4SJacob Faibussowitsch default: 454c095613aSPierre Jolivet recursive_build_strategy.set_partitioning_strategy(std::make_shared<htool::Partitioning<PetscReal, htool::ComputeLargestExtent<PetscReal>, htool::RegularSplitting<PetscReal>>>()); 45598e73e17SPierre Jolivet } 456c095613aSPierre Jolivet recursive_build_strategy.set_maximal_leaf_size(a->max_cluster_leaf_size); 457c095613aSPierre Jolivet a->source_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree_from_local_partition(A->cmap->N, a->dim, a->gcoords_source, 2, size, offset)); 458c7a4214aSPierre Jolivet S = uplo = 'N'; 4591dd4f53aSPierre Jolivet source_cluster = a->source_cluster.get(); 4601dd4f53aSPierre Jolivet } else source_cluster = a->target_cluster.get(); 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 4621dd4f53aSPierre Jolivet // generator 4631dd4f53aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx); 4641dd4f53aSPierre Jolivet else { 4651dd4f53aSPierre Jolivet a->wrapper = nullptr; 4661dd4f53aSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx); 4671dd4f53aSPierre Jolivet } 4681dd4f53aSPierre Jolivet // compressor 469c7a4214aSPierre Jolivet switch (a->compressor) { 470d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_FULL_ACA: 471c095613aSPierre Jolivet compressor = std::make_shared<htool::fullACA<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data()); 472d71ae5a4SJacob Faibussowitsch break; 473d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_SVD: 474c095613aSPierre Jolivet compressor = std::make_shared<htool::SVD<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data()); 475d71ae5a4SJacob Faibussowitsch break; 476d71ae5a4SJacob Faibussowitsch default: 477c095613aSPierre Jolivet compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, a->target_cluster->get_permutation().data(), source_cluster->get_permutation().data()); 478c7a4214aSPierre Jolivet } 4791dd4f53aSPierre Jolivet // local hierarchical matrix 4801dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 481c095613aSPierre Jolivet auto hmatrix_builder = htool::HMatrixTreeBuilder<PetscScalar>(a->epsilon, a->eta, S, uplo); 482c095613aSPierre Jolivet if (a->recompression) { 483c095613aSPierre Jolivet std::shared_ptr<htool::VirtualInternalLowRankGenerator<PetscScalar>> RecompressedLowRankGenerator = std::make_shared<htool::RecompressedLowRankGenerator<PetscScalar>>(*compressor, std::function<void(htool::LowRankMatrix<PetscScalar> &)>(htool::SVD_recompression<PetscScalar>)); 484c095613aSPierre Jolivet hmatrix_builder.set_low_rank_generator(RecompressedLowRankGenerator); 485c095613aSPierre Jolivet } else { 4861dd4f53aSPierre Jolivet hmatrix_builder.set_low_rank_generator(compressor); 487c095613aSPierre Jolivet } 4881dd4f53aSPierre Jolivet hmatrix_builder.set_minimal_target_depth(a->depth[0]); 4891dd4f53aSPierre Jolivet hmatrix_builder.set_minimal_source_depth(a->depth[1]); 4901dd4f53aSPierre Jolivet PetscCheck(a->block_tree_consistency || (!a->block_tree_consistency && !(A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE)), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot have a MatHtool with inconsistent block tree which is either symmetric or Hermitian"); 4911dd4f53aSPierre Jolivet hmatrix_builder.set_block_tree_consistency(a->block_tree_consistency); 492c095613aSPierre Jolivet a->distributed_operator_holder = std::make_unique<htool::DefaultApproximationBuilder<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, *a->target_cluster, *source_cluster, hmatrix_builder, PetscObjectComm((PetscObject)A)); 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 494c7a4214aSPierre Jolivet } 495c7a4214aSPierre Jolivet 496d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C) 497d71ae5a4SJacob Faibussowitsch { 498c7a4214aSPierre Jolivet Mat_Product *product = C->product; 499c9a4fd69SAudic XU Mat_Htool *a; 500c7a4214aSPierre Jolivet const PetscScalar *in; 501c7a4214aSPierre Jolivet PetscScalar *out; 5028b8ddd11SPierre Jolivet PetscInt N, lda; 503c7a4214aSPierre Jolivet 504c7a4214aSPierre Jolivet PetscFunctionBegin; 505c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 506f22e26b7SPierre Jolivet PetscCall(MatGetSize(C, nullptr, &N)); 5079566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 5085f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 5099566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B, &in)); 5109566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &out)); 511c9a4fd69SAudic XU PetscCall(MatShellGetContext(product->A, &a)); 5128b8ddd11SPierre Jolivet switch (product->type) { 513d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AB: 514c095613aSPierre Jolivet if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_matrix_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr); 515c095613aSPierre Jolivet else htool::internal_add_distributed_operator_matrix_product_local_to_local<PetscScalar>('N', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr); 516d71ae5a4SJacob Faibussowitsch break; 517d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AtB: 518c095613aSPierre Jolivet if (a->permutation == PETSC_TRUE) htool::add_distributed_operator_matrix_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr); 519c095613aSPierre Jolivet else htool::internal_add_distributed_operator_matrix_product_local_to_local<PetscScalar>('T', 1.0, a->distributed_operator_holder->distributed_operator, in, 0.0, out, N, nullptr); 520d71ae5a4SJacob Faibussowitsch break; 521d71ae5a4SJacob Faibussowitsch default: 522d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]); 523c7a4214aSPierre Jolivet } 5249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &out)); 5259566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B, &in)); 5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 527c7a4214aSPierre Jolivet } 528c7a4214aSPierre Jolivet 529d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C) 530d71ae5a4SJacob Faibussowitsch { 531c7a4214aSPierre Jolivet Mat_Product *product = C->product; 532c7a4214aSPierre Jolivet Mat A, B; 533c7a4214aSPierre Jolivet PetscBool flg; 534c7a4214aSPierre Jolivet 535c7a4214aSPierre Jolivet PetscFunctionBegin; 536c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 537c7a4214aSPierre Jolivet A = product->A; 538c7a4214aSPierre Jolivet B = product->B; 5399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, "")); 54097713f8cSPierre 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); 54197713f8cSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 54297713f8cSPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N)); 54397713f8cSPierre Jolivet else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N)); 5448b8ddd11SPierre Jolivet } 5459566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 5469566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5479566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 5489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 550f22e26b7SPierre Jolivet C->ops->productsymbolic = nullptr; 551c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 553c7a4214aSPierre Jolivet } 554c7a4214aSPierre Jolivet 555d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 556d71ae5a4SJacob Faibussowitsch { 557c7a4214aSPierre Jolivet PetscFunctionBegin; 558c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 5598b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 561c7a4214aSPierre Jolivet } 562c7a4214aSPierre Jolivet 5631dd4f53aSPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator) 564d71ae5a4SJacob Faibussowitsch { 565c9a4fd69SAudic XU Mat_Htool *a; 566c7a4214aSPierre Jolivet 567c7a4214aSPierre Jolivet PetscFunctionBegin; 568c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 5691dd4f53aSPierre Jolivet *distributed_operator = &a->distributed_operator_holder->distributed_operator; 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 571c7a4214aSPierre Jolivet } 572c7a4214aSPierre Jolivet 573c7a4214aSPierre Jolivet /*@C 57411a5261eSBarry Smith MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`. 575c7a4214aSPierre Jolivet 576cc4c1da9SBarry Smith No Fortran Support, No C Support 577cc4c1da9SBarry Smith 578c7a4214aSPierre Jolivet Input Parameter: 579c7a4214aSPierre Jolivet . A - hierarchical matrix 580c7a4214aSPierre Jolivet 581c7a4214aSPierre Jolivet Output Parameter: 5821dd4f53aSPierre Jolivet . distributed_operator - opaque pointer to a Htool virtual matrix 583c7a4214aSPierre Jolivet 584c7a4214aSPierre Jolivet Level: advanced 585c7a4214aSPierre Jolivet 5861cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL` 587c7a4214aSPierre Jolivet @*/ 5881dd4f53aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator) 589d71ae5a4SJacob Faibussowitsch { 590c7a4214aSPierre Jolivet PetscFunctionBegin; 591c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5921dd4f53aSPierre Jolivet PetscAssertPointer(distributed_operator, 2); 5931dd4f53aSPierre Jolivet PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::DistributedOperator<PetscScalar> **), (A, distributed_operator)); 5943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 595c7a4214aSPierre Jolivet } 596c7a4214aSPierre Jolivet 5978434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 598d71ae5a4SJacob Faibussowitsch { 599c9a4fd69SAudic XU Mat_Htool *a; 600c7a4214aSPierre Jolivet 601c7a4214aSPierre Jolivet PetscFunctionBegin; 602c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 603c7a4214aSPierre Jolivet a->kernel = kernel; 604c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 605c7a4214aSPierre Jolivet delete a->wrapper; 6061dd4f53aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx); 6073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 608c7a4214aSPierre Jolivet } 609c7a4214aSPierre Jolivet 610c7a4214aSPierre Jolivet /*@C 61111a5261eSBarry Smith MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`. 612c7a4214aSPierre Jolivet 613cc4c1da9SBarry Smith Collective, No Fortran Support 614cc4c1da9SBarry Smith 615c7a4214aSPierre Jolivet Input Parameters: 616c7a4214aSPierre Jolivet + A - hierarchical matrix 6172ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 6182ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 619c7a4214aSPierre Jolivet 620c7a4214aSPierre Jolivet Level: advanced 621c7a4214aSPierre Jolivet 6221cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()` 623c7a4214aSPierre Jolivet @*/ 624cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 625d71ae5a4SJacob Faibussowitsch { 626c7a4214aSPierre Jolivet PetscFunctionBegin; 627c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 628c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 2); 6294f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 3); 6308434afd1SBarry Smith PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx)); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 632c7a4214aSPierre Jolivet } 633c7a4214aSPierre Jolivet 634d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) 635d71ae5a4SJacob Faibussowitsch { 636c9a4fd69SAudic XU Mat_Htool *a; 6371dd4f53aSPierre Jolivet PetscMPIInt rank; 6381dd4f53aSPierre Jolivet const std::vector<PetscInt> *source; 6391dd4f53aSPierre Jolivet const htool::Cluster<PetscReal> *local_source_cluster; 64098e73e17SPierre Jolivet 64198e73e17SPierre Jolivet PetscFunctionBegin; 642c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 6431dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 6441dd4f53aSPierre Jolivet local_source_cluster = a->source_cluster ? &a->source_cluster->get_cluster_on_partition(rank) : &a->target_cluster->get_cluster_on_partition(rank); 6451dd4f53aSPierre Jolivet source = &local_source_cluster->get_permutation(); 6461dd4f53aSPierre Jolivet PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), local_source_cluster->get_size(), source->data() + local_source_cluster->get_offset(), PETSC_COPY_VALUES, is)); 6479566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64998e73e17SPierre Jolivet } 65098e73e17SPierre Jolivet 651cc4c1da9SBarry Smith /*@ 65211a5261eSBarry Smith MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix. 65398e73e17SPierre Jolivet 65498e73e17SPierre Jolivet Input Parameter: 65598e73e17SPierre Jolivet . A - hierarchical matrix 65698e73e17SPierre Jolivet 65798e73e17SPierre Jolivet Output Parameter: 65898e73e17SPierre Jolivet . is - permutation 65998e73e17SPierre Jolivet 66098e73e17SPierre Jolivet Level: advanced 66198e73e17SPierre Jolivet 6621cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 66398e73e17SPierre Jolivet @*/ 664cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) 665d71ae5a4SJacob Faibussowitsch { 66698e73e17SPierre Jolivet PetscFunctionBegin; 66798e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6684f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 669cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is)); 6703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67198e73e17SPierre Jolivet } 67298e73e17SPierre Jolivet 673d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) 674d71ae5a4SJacob Faibussowitsch { 675c9a4fd69SAudic XU Mat_Htool *a; 6761dd4f53aSPierre Jolivet const std::vector<PetscInt> *target; 6771dd4f53aSPierre Jolivet PetscMPIInt rank; 67898e73e17SPierre Jolivet 67998e73e17SPierre Jolivet PetscFunctionBegin; 680c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 6811dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 6821dd4f53aSPierre Jolivet target = &a->target_cluster->get_permutation(); 6831dd4f53aSPierre Jolivet PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), a->target_cluster->get_cluster_on_partition(rank).get_size(), target->data() + a->target_cluster->get_cluster_on_partition(rank).get_offset(), PETSC_COPY_VALUES, is)); 6849566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68698e73e17SPierre Jolivet } 68798e73e17SPierre Jolivet 688cc4c1da9SBarry Smith /*@ 68911a5261eSBarry Smith MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix. 69098e73e17SPierre Jolivet 69198e73e17SPierre Jolivet Input Parameter: 69298e73e17SPierre Jolivet . A - hierarchical matrix 69398e73e17SPierre Jolivet 69498e73e17SPierre Jolivet Output Parameter: 69598e73e17SPierre Jolivet . is - permutation 69698e73e17SPierre Jolivet 69798e73e17SPierre Jolivet Level: advanced 69898e73e17SPierre Jolivet 6991cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 70098e73e17SPierre Jolivet @*/ 701cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) 702d71ae5a4SJacob Faibussowitsch { 70398e73e17SPierre Jolivet PetscFunctionBegin; 70498e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 7054f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 706cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is)); 7073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70898e73e17SPierre Jolivet } 70998e73e17SPierre Jolivet 710d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) 711d71ae5a4SJacob Faibussowitsch { 712c9a4fd69SAudic XU Mat_Htool *a; 71398e73e17SPierre Jolivet 71498e73e17SPierre Jolivet PetscFunctionBegin; 715c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 716c095613aSPierre Jolivet a->permutation = use; 7173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71898e73e17SPierre Jolivet } 71998e73e17SPierre Jolivet 720cc4c1da9SBarry Smith /*@ 72111a5261eSBarry Smith MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation. 72298e73e17SPierre Jolivet 72398e73e17SPierre Jolivet Input Parameters: 72498e73e17SPierre Jolivet + A - hierarchical matrix 72598e73e17SPierre Jolivet - use - Boolean value 72698e73e17SPierre Jolivet 72798e73e17SPierre Jolivet Level: advanced 72898e73e17SPierre Jolivet 7291cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 73098e73e17SPierre Jolivet @*/ 731cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) 732d71ae5a4SJacob Faibussowitsch { 73398e73e17SPierre Jolivet PetscFunctionBegin; 73498e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 73598e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 736cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use)); 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73898e73e17SPierre Jolivet } 73998e73e17SPierre Jolivet 740c095613aSPierre Jolivet static PetscErrorCode MatHtoolUseRecompression_Htool(Mat A, PetscBool use) 741c095613aSPierre Jolivet { 742c095613aSPierre Jolivet Mat_Htool *a; 743c095613aSPierre Jolivet 744c095613aSPierre Jolivet PetscFunctionBegin; 745c095613aSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 746c095613aSPierre Jolivet a->recompression = use; 747c095613aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 748c095613aSPierre Jolivet } 749c095613aSPierre Jolivet 750c095613aSPierre Jolivet /*@ 751c095613aSPierre Jolivet MatHtoolUseRecompression - Sets whether a `MATHTOOL` matrix should use recompression. 752c095613aSPierre Jolivet 753c095613aSPierre Jolivet Input Parameters: 754c095613aSPierre Jolivet + A - hierarchical matrix 755c095613aSPierre Jolivet - use - Boolean value 756c095613aSPierre Jolivet 757c095613aSPierre Jolivet Level: advanced 758c095613aSPierre Jolivet 759c095613aSPierre Jolivet .seealso: [](ch_matrices), `Mat`, `MATHTOOL` 760c095613aSPierre Jolivet @*/ 761c095613aSPierre Jolivet PetscErrorCode MatHtoolUseRecompression(Mat A, PetscBool use) 762c095613aSPierre Jolivet { 763c095613aSPierre Jolivet PetscFunctionBegin; 764c095613aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 765c095613aSPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 766c095613aSPierre Jolivet PetscTryMethod(A, "MatHtoolUseRecompression_C", (Mat, PetscBool), (A, use)); 767c095613aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 768c095613aSPierre Jolivet } 769c095613aSPierre Jolivet 770cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 771d71ae5a4SJacob Faibussowitsch { 772c7a4214aSPierre Jolivet Mat C; 773c9a4fd69SAudic XU Mat_Htool *a; 77466b4df27SPierre Jolivet PetscScalar *array, shift, scale; 775c7a4214aSPierre Jolivet PetscInt lda; 776c7a4214aSPierre Jolivet 777c7a4214aSPierre Jolivet PetscFunctionBegin; 77866b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 779c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 780c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 781c7a4214aSPierre Jolivet C = *B; 7825f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions"); 7839566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 7845f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 785c7a4214aSPierre Jolivet } else { 7869566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7879566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 7889566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 7899566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7909566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 791c7a4214aSPierre Jolivet } 7929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 7939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 7949566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 7951dd4f53aSPierre Jolivet htool::copy_to_dense_in_user_numbering(a->distributed_operator_holder->hmatrix, array); 7969566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 79766b4df27SPierre Jolivet PetscCall(MatShift(C, shift)); 79866b4df27SPierre Jolivet PetscCall(MatScale(C, scale)); 799ac530a7eSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) PetscCall(MatHeaderReplace(A, &C)); 800ac530a7eSPierre Jolivet else *B = C; 8013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 802c7a4214aSPierre Jolivet } 803c7a4214aSPierre Jolivet 804*2a8381b2SBarry Smith static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, PetscCtx ctx) 805d71ae5a4SJacob Faibussowitsch { 806c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx; 80798e73e17SPierre Jolivet PetscScalar *tmp; 80898e73e17SPierre Jolivet 80998e73e17SPierre Jolivet PetscFunctionBegin; 8103ba16761SJacob Faibussowitsch PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx)); 8119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * N, &tmp)); 8129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, ptr, M * N)); 81398e73e17SPierre Jolivet for (PetscInt i = 0; i < M; ++i) { 81498e73e17SPierre Jolivet for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N]; 81598e73e17SPierre Jolivet } 8169566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 8173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 818c7a4214aSPierre Jolivet } 819c7a4214aSPierre Jolivet 820c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 821d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) 822d71ae5a4SJacob Faibussowitsch { 823c7a4214aSPierre Jolivet Mat C; 824c9a4fd69SAudic XU Mat_Htool *a, *c; 82566b4df27SPierre Jolivet PetscScalar shift, scale; 826c7a4214aSPierre Jolivet PetscInt M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n; 827c7a4214aSPierre Jolivet PetscContainer container; 828c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 829c7a4214aSPierre Jolivet 830c7a4214aSPierre Jolivet PetscFunctionBegin; 83166b4df27SPierre Jolivet PetscCall(MatShellGetScalingShifts(A, &shift, &scale, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Vec *)MAT_SHELL_NOT_ALLOWED, (Mat *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED, (IS *)MAT_SHELL_NOT_ALLOWED)); 832c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 8337fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 8345f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported"); 835c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8369566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 8379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, n, m, N, M)); 8389566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 8399566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 8409566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 84149abdd8aSBarry Smith PetscCall(PetscObjectContainerCompose((PetscObject)C, "KernelTranspose", kernelt, PetscCtxDestroyDefault)); 842c7a4214aSPierre Jolivet } else { 843c7a4214aSPierre Jolivet C = *B; 8449566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container)); 8455f80ce2aSJacob Faibussowitsch PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 846*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, &kernelt)); 847c7a4214aSPierre Jolivet } 848c9a4fd69SAudic XU PetscCall(MatShellGetContext(C, &c)); 849c7a4214aSPierre Jolivet c->dim = a->dim; 85066b4df27SPierre Jolivet PetscCall(MatShift(C, shift)); 85166b4df27SPierre Jolivet PetscCall(MatScale(C, scale)); 85298e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 853c7a4214aSPierre Jolivet if (kernelt->A != A) { 8549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 855c7a4214aSPierre Jolivet kernelt->A = A; 8569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 857c7a4214aSPierre Jolivet } 858c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 859c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 860c7a4214aSPierre Jolivet c->kernelctx = kernelt; 861c095613aSPierre Jolivet c->max_cluster_leaf_size = a->max_cluster_leaf_size; 862de647125SPierre Marchand c->epsilon = a->epsilon; 863de647125SPierre Marchand c->eta = a->eta; 864de647125SPierre Marchand c->block_tree_consistency = a->block_tree_consistency; 865c095613aSPierre Jolivet c->permutation = a->permutation; 866c095613aSPierre Jolivet c->recompression = a->recompression; 867de647125SPierre Marchand c->compressor = a->compressor; 868de647125SPierre Marchand c->clustering = a->clustering; 869c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target)); 8719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim)); 872c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 8739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source)); 8749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim)); 875c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 876c7a4214aSPierre Jolivet } 8779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 8789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 879c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 8803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 881c7a4214aSPierre Jolivet } 882c7a4214aSPierre Jolivet 8831dd4f53aSPierre Jolivet static PetscErrorCode MatDestroy_Factor(Mat F) 8841dd4f53aSPierre Jolivet { 8851dd4f53aSPierre Jolivet PetscContainer container; 8861dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *A; 8871dd4f53aSPierre Jolivet 8881dd4f53aSPierre Jolivet PetscFunctionBegin; 8891dd4f53aSPierre Jolivet PetscCall(PetscObjectQuery((PetscObject)F, "HMatrix", (PetscObject *)&container)); 8901dd4f53aSPierre Jolivet if (container) { 891*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, &A)); 8921dd4f53aSPierre Jolivet delete A; 8931dd4f53aSPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)F, "HMatrix", nullptr)); 8941dd4f53aSPierre Jolivet } 8951dd4f53aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", nullptr)); 8961dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 8971dd4f53aSPierre Jolivet } 8981dd4f53aSPierre Jolivet 8991dd4f53aSPierre Jolivet static PetscErrorCode MatFactorGetSolverType_Htool(Mat, MatSolverType *type) 9001dd4f53aSPierre Jolivet { 9011dd4f53aSPierre Jolivet PetscFunctionBegin; 9021dd4f53aSPierre Jolivet *type = MATSOLVERHTOOL; 9031dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9041dd4f53aSPierre Jolivet } 9051dd4f53aSPierre Jolivet 9061dd4f53aSPierre Jolivet template <char trans> 9071dd4f53aSPierre Jolivet static inline PetscErrorCode MatSolve_Private(Mat A, htool::Matrix<PetscScalar> &X) 9081dd4f53aSPierre Jolivet { 9091dd4f53aSPierre Jolivet PetscContainer container; 9101dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *B; 9111dd4f53aSPierre Jolivet 9121dd4f53aSPierre Jolivet PetscFunctionBegin; 9131dd4f53aSPierre Jolivet PetscCheck(A->factortype == MAT_FACTOR_LU || A->factortype == MAT_FACTOR_CHOLESKY, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_UNKNOWN_TYPE, "Only MAT_LU_FACTOR and MAT_CHOLESKY_FACTOR are supported"); 9141dd4f53aSPierre Jolivet PetscCall(PetscObjectQuery((PetscObject)A, "HMatrix", (PetscObject *)&container)); 9151dd4f53aSPierre Jolivet PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call Mat%sFactorNumeric() before Mat%sSolve%s()", A->factortype == MAT_FACTOR_LU ? "LU" : "Cholesky", X.nb_cols() == 1 ? "" : "Mat", trans == 'N' ? "" : "Transpose"); 916*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(container, &B)); 9171dd4f53aSPierre Jolivet if (A->factortype == MAT_FACTOR_LU) htool::lu_solve(trans, *B, X); 9181dd4f53aSPierre Jolivet else htool::cholesky_solve('L', *B, X); 9191dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9201dd4f53aSPierre Jolivet } 9211dd4f53aSPierre Jolivet 9221dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type * = nullptr> 9231dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type b, Type x) 9241dd4f53aSPierre Jolivet { 9251dd4f53aSPierre Jolivet PetscInt n; 9261dd4f53aSPierre Jolivet htool::Matrix<PetscScalar> v; 9271dd4f53aSPierre Jolivet PetscScalar *array; 9281dd4f53aSPierre Jolivet 9291dd4f53aSPierre Jolivet PetscFunctionBegin; 9301dd4f53aSPierre Jolivet PetscCall(VecGetLocalSize(b, &n)); 9311dd4f53aSPierre Jolivet PetscCall(VecCopy(b, x)); 9321dd4f53aSPierre Jolivet PetscCall(VecGetArrayWrite(x, &array)); 9331dd4f53aSPierre Jolivet v.assign(n, 1, array, false); 9341dd4f53aSPierre Jolivet PetscCall(VecRestoreArrayWrite(x, &array)); 9351dd4f53aSPierre Jolivet PetscCall(MatSolve_Private<trans>(A, v)); 9361dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9371dd4f53aSPierre Jolivet } 9381dd4f53aSPierre Jolivet 9391dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type * = nullptr> 9401dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type B, Type X) 9411dd4f53aSPierre Jolivet { 9421dd4f53aSPierre Jolivet PetscInt m, N; 9431dd4f53aSPierre Jolivet htool::Matrix<PetscScalar> v; 9441dd4f53aSPierre Jolivet PetscScalar *array; 9451dd4f53aSPierre Jolivet 9461dd4f53aSPierre Jolivet PetscFunctionBegin; 9471dd4f53aSPierre Jolivet PetscCall(MatGetLocalSize(B, &m, nullptr)); 9481dd4f53aSPierre Jolivet PetscCall(MatGetLocalSize(B, nullptr, &N)); 9491dd4f53aSPierre Jolivet PetscCall(MatCopy(B, X, SAME_NONZERO_PATTERN)); 9501dd4f53aSPierre Jolivet PetscCall(MatDenseGetArrayWrite(X, &array)); 9511dd4f53aSPierre Jolivet v.assign(m, N, array, false); 9521dd4f53aSPierre Jolivet PetscCall(MatDenseRestoreArrayWrite(X, &array)); 9531dd4f53aSPierre Jolivet PetscCall(MatSolve_Private<trans>(A, v)); 9541dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9551dd4f53aSPierre Jolivet } 9561dd4f53aSPierre Jolivet 9571dd4f53aSPierre Jolivet template <MatFactorType ftype> 9581dd4f53aSPierre Jolivet static PetscErrorCode MatFactorNumeric_Htool(Mat F, Mat A, const MatFactorInfo *) 9591dd4f53aSPierre Jolivet { 9601dd4f53aSPierre Jolivet Mat_Htool *a; 9611dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *B; 9621dd4f53aSPierre Jolivet 9631dd4f53aSPierre Jolivet PetscFunctionBegin; 9641dd4f53aSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 9651dd4f53aSPierre Jolivet B = new htool::HMatrix<PetscScalar>(a->distributed_operator_holder->hmatrix); 966c095613aSPierre Jolivet if (ftype == MAT_FACTOR_LU) htool::sequential_lu_factorization(*B); 967c095613aSPierre Jolivet else htool::sequential_cholesky_factorization('L', *B); 9681dd4f53aSPierre Jolivet PetscCall(PetscObjectContainerCompose((PetscObject)F, "HMatrix", B, nullptr)); 9691dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9701dd4f53aSPierre Jolivet } 9711dd4f53aSPierre Jolivet 9721dd4f53aSPierre Jolivet template <MatFactorType ftype> 9731dd4f53aSPierre Jolivet PetscErrorCode MatFactorSymbolic_Htool(Mat F, Mat) 9741dd4f53aSPierre Jolivet { 9751dd4f53aSPierre Jolivet PetscFunctionBegin; 9761dd4f53aSPierre Jolivet F->preallocated = PETSC_TRUE; 9771dd4f53aSPierre Jolivet F->assembled = PETSC_TRUE; 9781dd4f53aSPierre Jolivet F->ops->solve = MatSolve_Htool<'N', Vec>; 9791dd4f53aSPierre Jolivet F->ops->matsolve = MatSolve_Htool<'N', Mat>; 9801dd4f53aSPierre Jolivet if (!PetscDefined(USE_COMPLEX) || ftype == MAT_FACTOR_LU) { 9811dd4f53aSPierre Jolivet F->ops->solvetranspose = MatSolve_Htool<'T', Vec>; 9821dd4f53aSPierre Jolivet F->ops->matsolvetranspose = MatSolve_Htool<'T', Mat>; 9831dd4f53aSPierre Jolivet } 9841dd4f53aSPierre Jolivet F->ops->destroy = MatDestroy_Factor; 9851dd4f53aSPierre Jolivet if (ftype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_LU>; 9861dd4f53aSPierre Jolivet else F->ops->choleskyfactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_CHOLESKY>; 9871dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9881dd4f53aSPierre Jolivet } 9891dd4f53aSPierre Jolivet 9901dd4f53aSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Htool(Mat F, Mat A, IS, IS, const MatFactorInfo *) 9911dd4f53aSPierre Jolivet { 9921dd4f53aSPierre Jolivet PetscFunctionBegin; 9931dd4f53aSPierre Jolivet PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_LU>(F, A)); 9941dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9951dd4f53aSPierre Jolivet } 9961dd4f53aSPierre Jolivet 9971dd4f53aSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Htool(Mat F, Mat A, IS, const MatFactorInfo *) 9981dd4f53aSPierre Jolivet { 9991dd4f53aSPierre Jolivet PetscFunctionBegin; 10001dd4f53aSPierre Jolivet PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_CHOLESKY>(F, A)); 10011dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 10021dd4f53aSPierre Jolivet } 10031dd4f53aSPierre Jolivet 10041dd4f53aSPierre Jolivet static PetscErrorCode MatGetFactor_htool_htool(Mat A, MatFactorType ftype, Mat *F) 10051dd4f53aSPierre Jolivet { 10061dd4f53aSPierre Jolivet Mat B; 10071dd4f53aSPierre Jolivet Mat_Htool *a; 10081dd4f53aSPierre Jolivet PetscMPIInt size; 10091dd4f53aSPierre Jolivet 10101dd4f53aSPierre Jolivet PetscFunctionBegin; 10111dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 10121dd4f53aSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 10131dd4f53aSPierre Jolivet PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel MatGetFactor()"); 10141dd4f53aSPierre Jolivet PetscCheck(a->block_tree_consistency, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot factor a MatHtool with inconsistent block tree"); 10151dd4f53aSPierre Jolivet PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 10161dd4f53aSPierre Jolivet PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 10171dd4f53aSPierre Jolivet PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &((PetscObject)B)->type_name)); 10181dd4f53aSPierre Jolivet PetscCall(MatSetUp(B)); 10191dd4f53aSPierre Jolivet 10201dd4f53aSPierre Jolivet B->ops->getinfo = MatGetInfo_External; 10211dd4f53aSPierre Jolivet B->factortype = ftype; 10221dd4f53aSPierre Jolivet B->trivialsymbolic = PETSC_TRUE; 10231dd4f53aSPierre Jolivet 10241dd4f53aSPierre Jolivet if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_Htool; 10251dd4f53aSPierre Jolivet else if (ftype == MAT_FACTOR_CHOLESKY) B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_Htool; 10261dd4f53aSPierre Jolivet 10271dd4f53aSPierre Jolivet PetscCall(PetscFree(B->solvertype)); 10281dd4f53aSPierre Jolivet PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &B->solvertype)); 10291dd4f53aSPierre Jolivet 10301dd4f53aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_Htool)); 10311dd4f53aSPierre Jolivet *F = B; 10321dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 10331dd4f53aSPierre Jolivet } 10341dd4f53aSPierre Jolivet 10351dd4f53aSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Htool(void) 10361dd4f53aSPierre Jolivet { 10371dd4f53aSPierre Jolivet PetscFunctionBegin; 10381dd4f53aSPierre Jolivet PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_LU, MatGetFactor_htool_htool)); 10391dd4f53aSPierre Jolivet PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_CHOLESKY, MatGetFactor_htool_htool)); 10401dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 10411dd4f53aSPierre Jolivet } 10421dd4f53aSPierre Jolivet 1043c7a4214aSPierre Jolivet /*@C 104411a5261eSBarry Smith MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel. 1045c7a4214aSPierre Jolivet 1046cc4c1da9SBarry Smith Collective, No Fortran Support 1047cc4c1da9SBarry Smith 1048c7a4214aSPierre Jolivet Input Parameters: 1049c7a4214aSPierre Jolivet + comm - MPI communicator 10502ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 10512ef1f0ffSBarry Smith . n - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given) 10522ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 10532ef1f0ffSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 1054c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 1055c7a4214aSPierre Jolivet . coords_target - coordinates of the target 1056c7a4214aSPierre Jolivet . coords_source - coordinates of the source 10572ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 10582ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 1059c7a4214aSPierre Jolivet 1060c7a4214aSPierre Jolivet Output Parameter: 1061c7a4214aSPierre Jolivet . B - matrix 1062c7a4214aSPierre Jolivet 1063c7a4214aSPierre Jolivet Options Database Keys: 1064c095613aSPierre Jolivet + -mat_htool_max_cluster_leaf_size <`PetscInt`> - maximal leaf size in cluster tree 106511a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`> - relative error in Frobenius norm when approximating a block 106611a5261eSBarry Smith . -mat_htool_eta <`PetscReal`> - admissibility condition tolerance 106711a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`> - minimal cluster tree depth associated with the rows 106811a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`> - minimal cluster tree depth associated with the columns 10691dd4f53aSPierre Jolivet . -mat_htool_block_tree_consistency <`PetscBool`> - block tree consistency 1070c095613aSPierre Jolivet . -mat_htool_recompression <`PetscBool`> - use recompression 107198e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 107298e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 1073c7a4214aSPierre Jolivet 1074c7a4214aSPierre Jolivet Level: intermediate 1075c7a4214aSPierre Jolivet 10761cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 1077c7a4214aSPierre Jolivet @*/ 10788434afd1SBarry 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) 1079d71ae5a4SJacob Faibussowitsch { 1080c7a4214aSPierre Jolivet Mat A; 1081c7a4214aSPierre Jolivet Mat_Htool *a; 1082c7a4214aSPierre Jolivet 1083c7a4214aSPierre Jolivet PetscFunctionBegin; 10849566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 1085c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A, spacedim, 6); 10864f572ea9SToby Isaac PetscAssertPointer(coords_target, 7); 10874f572ea9SToby Isaac PetscAssertPointer(coords_source, 8); 1088c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 9); 10894f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 10); 10909566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, M, N)); 10919566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATHTOOL)); 10929566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1093c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 1094c7a4214aSPierre Jolivet a->dim = spacedim; 1095c7a4214aSPierre Jolivet a->kernel = kernel; 1096c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 10979566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target)); 10989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim)); 1099462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */ 1100c7a4214aSPierre Jolivet if (coords_target != coords_source) { 11019566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source)); 11029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim)); 1103462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */ 1104c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 1105c7a4214aSPierre Jolivet *B = A; 11063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1107c7a4214aSPierre Jolivet } 1108c7a4214aSPierre Jolivet 1109c7a4214aSPierre Jolivet /*MC 1110c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 1111c7a4214aSPierre Jolivet 11122ef1f0ffSBarry Smith Use `./configure --download-htool` to install PETSc to use Htool. 1113c7a4214aSPierre Jolivet 11142ef1f0ffSBarry Smith Options Database Key: 11152ef1f0ffSBarry Smith . -mat_type htool - matrix type to `MATHTOOL` 1116c7a4214aSPierre Jolivet 1117c7a4214aSPierre Jolivet Level: beginner 1118c7a4214aSPierre Jolivet 11191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 1120c7a4214aSPierre Jolivet M*/ 1121d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 1122d71ae5a4SJacob Faibussowitsch { 1123c7a4214aSPierre Jolivet Mat_Htool *a; 1124c7a4214aSPierre Jolivet 1125c7a4214aSPierre Jolivet PetscFunctionBegin; 1126c9a4fd69SAudic XU PetscCall(MatSetType(A, MATSHELL)); 11274dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 1128c9a4fd69SAudic XU PetscCall(MatShellSetContext(A, a)); 112957d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (PetscErrorCodeFn *)MatGetDiagonal_Htool)); 113057d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (PetscErrorCodeFn *)MatGetDiagonalBlock_Htool)); 113157d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_MULT, (PetscErrorCodeFn *)MatMult_Htool)); 113257d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Htool)); 113357d50842SBarry Smith if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_Htool)); 1134c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 1135c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 113657d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_VIEW, (PetscErrorCodeFn *)MatView_Htool)); 113757d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (PetscErrorCodeFn *)MatSetFromOptions_Htool)); 113857d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (PetscErrorCodeFn *)MatGetRow_Htool)); 113957d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (PetscErrorCodeFn *)MatRestoreRow_Htool)); 114057d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (PetscErrorCodeFn *)MatAssemblyEnd_Htool)); 114157d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (PetscErrorCodeFn *)MatTranspose_Htool)); 114257d50842SBarry Smith PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_Htool)); 1143c7a4214aSPierre Jolivet a->dim = 0; 1144f22e26b7SPierre Jolivet a->gcoords_target = nullptr; 1145f22e26b7SPierre Jolivet a->gcoords_source = nullptr; 1146c095613aSPierre Jolivet a->max_cluster_leaf_size = 10; 1147c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 1148c7a4214aSPierre Jolivet a->eta = 10.0; 1149c7a4214aSPierre Jolivet a->depth[0] = 0; 1150c7a4214aSPierre Jolivet a->depth[1] = 0; 11511dd4f53aSPierre Jolivet a->block_tree_consistency = PETSC_TRUE; 1152c095613aSPierre Jolivet a->permutation = PETSC_TRUE; 1153c095613aSPierre Jolivet a->recompression = PETSC_FALSE; 1154c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 11559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool)); 11569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool)); 11579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense)); 11589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense)); 11599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool)); 11609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool)); 11619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool)); 11629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool)); 11639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool)); 1164c095613aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUseRecompression_C", MatHtoolUseRecompression_Htool)); 1165c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable)); 1166c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 1167c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 1168c9a4fd69SAudic XU PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL)); 11693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1170c7a4214aSPierre Jolivet } 1171