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)); 791dd4f53aSPierre Jolivet a->distributed_operator_holder->distributed_operator.vector_product_local_to_local(in, out, nullptr); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83c7a4214aSPierre Jolivet } 84c7a4214aSPierre Jolivet 85d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Htool(Mat A, Vec x, Vec y) 86d71ae5a4SJacob Faibussowitsch { 87c9a4fd69SAudic XU Mat_Htool *a; 888b8ddd11SPierre Jolivet const PetscScalar *in; 898b8ddd11SPierre Jolivet PetscScalar *out; 908b8ddd11SPierre Jolivet 918b8ddd11SPierre Jolivet PetscFunctionBegin; 92c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &in)); 949566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y, &out)); 951dd4f53aSPierre Jolivet a->distributed_operator_holder->distributed_operator.vector_product_transp_local_to_local(in, out, nullptr); 969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &in)); 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y, &out)); 983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 998b8ddd11SPierre Jolivet } 1008b8ddd11SPierre Jolivet 101d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_Htool(Mat A, PetscInt is_max, IS is[], PetscInt ov) 102d71ae5a4SJacob Faibussowitsch { 103c7a4214aSPierre Jolivet std::set<PetscInt> set; 104c7a4214aSPierre Jolivet const PetscInt *idx; 1058f308287SPierre Jolivet PetscInt *oidx, size, bs[2]; 106c7a4214aSPierre Jolivet PetscMPIInt csize; 107c7a4214aSPierre Jolivet 108c7a4214aSPierre Jolivet PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, bs, bs + 1)); 1108f308287SPierre Jolivet if (bs[0] != bs[1]) bs[0] = 1; 111c7a4214aSPierre Jolivet for (PetscInt i = 0; i < is_max; ++i) { 112c7a4214aSPierre Jolivet /* basic implementation that adds indices by shifting an IS by -ov, -ov+1..., -1, 1..., ov-1, ov */ 113c7a4214aSPierre Jolivet /* needed to avoid subdomain matrices to replicate A since it is dense */ 1149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is[i]), &csize)); 1151dd4f53aSPierre Jolivet PetscCheck(csize == 1, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel IS"); 1169566063dSJacob Faibussowitsch PetscCall(ISGetSize(is[i], &size)); 1179566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx)); 118c7a4214aSPierre Jolivet for (PetscInt j = 0; j < size; ++j) { 119c7a4214aSPierre Jolivet set.insert(idx[j]); 120c7a4214aSPierre Jolivet for (PetscInt k = 1; k <= ov; ++k) { /* for each layer of overlap */ 121c7a4214aSPierre Jolivet if (idx[j] - k >= 0) set.insert(idx[j] - k); /* do not insert negative indices */ 122c7a4214aSPierre 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 */ 123c7a4214aSPierre Jolivet } 124c7a4214aSPierre Jolivet } 1259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[i], &idx)); 1269566063dSJacob Faibussowitsch PetscCall(ISDestroy(is + i)); 1278f308287SPierre Jolivet if (bs[0] > 1) { 1288f308287SPierre Jolivet for (std::set<PetscInt>::iterator it = set.cbegin(); it != set.cend(); it++) { 1298f308287SPierre Jolivet std::vector<PetscInt> block(bs[0]); 1308f308287SPierre Jolivet std::iota(block.begin(), block.end(), (*it / bs[0]) * bs[0]); 1318f308287SPierre Jolivet set.insert(block.cbegin(), block.cend()); 1328f308287SPierre Jolivet } 1338f308287SPierre Jolivet } 134c7a4214aSPierre Jolivet size = set.size(); /* size with overlap */ 1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &oidx)); 136c7a4214aSPierre Jolivet for (const PetscInt j : set) *oidx++ = j; 137c7a4214aSPierre Jolivet oidx -= size; 1389566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, oidx, PETSC_OWN_POINTER, is + i)); 139c7a4214aSPierre Jolivet } 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141c7a4214aSPierre Jolivet } 142c7a4214aSPierre Jolivet 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_Htool(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 144d71ae5a4SJacob Faibussowitsch { 145c9a4fd69SAudic XU Mat_Htool *a; 146c7a4214aSPierre Jolivet Mat D, B, BT; 147c7a4214aSPierre Jolivet const PetscScalar *copy; 14866b4df27SPierre Jolivet PetscScalar *ptr, shift, scale; 149c7a4214aSPierre Jolivet const PetscInt *idxr, *idxc, *it; 150c7a4214aSPierre Jolivet PetscInt nrow, m, i; 151c7a4214aSPierre Jolivet PetscBool flg; 152c7a4214aSPierre Jolivet 153c7a4214aSPierre Jolivet PetscFunctionBegin; 15466b4df27SPierre 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)); 155c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 15648a46eb9SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(PetscCalloc1(n, submat)); 157c7a4214aSPierre Jolivet for (i = 0; i < n; ++i) { 1589566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &nrow)); 1599566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &m)); 1609566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i], &idxr)); 1619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i], &idxc)); 162f22e26b7SPierre Jolivet if (scall != MAT_REUSE_MATRIX) PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow, m, nrow, m, nullptr, (*submat) + i)); 1639566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite((*submat)[i], &ptr)); 164c7a4214aSPierre Jolivet if (irow[i] == icol[i]) { /* same row and column IS? */ 1659566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A, &flg)); 166c7a4214aSPierre Jolivet if (flg) { 1679566063dSJacob Faibussowitsch PetscCall(ISSorted(irow[i], &flg)); 168c7a4214aSPierre Jolivet if (flg) { /* sorted IS? */ 169c7a4214aSPierre Jolivet it = std::lower_bound(idxr, idxr + nrow, A->rmap->rstart); 170c7a4214aSPierre Jolivet if (it != idxr + nrow && *it == A->rmap->rstart) { /* rmap->rstart in IS? */ 171c7a4214aSPierre Jolivet if (std::distance(idxr, it) + A->rmap->n <= nrow) { /* long enough IS to store the local diagonal block? */ 1729371c9d4SSatish Balay for (PetscInt j = 0; j < A->rmap->n && flg; ++j) 1739371c9d4SSatish Balay if (PetscUnlikely(it[j] != A->rmap->rstart + j)) flg = PETSC_FALSE; 174c7a4214aSPierre Jolivet if (flg) { /* complete local diagonal block in IS? */ 175c7a4214aSPierre Jolivet /* fast extraction when the local diagonal block is part of the submatrix, e.g., for PCASM or PCHPDDM 176c7a4214aSPierre Jolivet * [ B C E ] 177c7a4214aSPierre Jolivet * A = [ B D E ] 178c7a4214aSPierre Jolivet * [ B F E ] 179c7a4214aSPierre Jolivet */ 180c7a4214aSPierre Jolivet m = std::distance(idxr, it); /* shift of the coefficient (0,0) of block D from above */ 181c9a4fd69SAudic XU PetscCall(MatGetDiagonalBlock(A, &D)); 1829566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(D, ©)); 1839371c9d4SSatish Balay for (PetscInt k = 0; k < A->rmap->n; ++k) { PetscCall(PetscArraycpy(ptr + (m + k) * nrow + m, copy + k * A->rmap->n, A->rmap->n)); /* block D from above */ } 1849566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(D, ©)); 185c7a4214aSPierre Jolivet if (m) { 186c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* vertical block B from above */ 187c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 188b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 1899566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, m, A->rmap->n, m, ptr + m, &B)); 1909566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 1919566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, m, A->rmap->n, m, A->rmap->n, ptr + m * nrow, &BT)); 1929566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 193b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 1949566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 195c7a4214aSPierre Jolivet } else { 1967fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 1979566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 198c7a4214aSPierre Jolivet } 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 201c7a4214aSPierre Jolivet } else { 202c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block C from above */ 203c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(m, 1, idxr, idxc + m + k, ptr + (m + k) * nrow); 204c7a4214aSPierre Jolivet } 205c7a4214aSPierre Jolivet } 206c7a4214aSPierre Jolivet } 207c7a4214aSPierre Jolivet if (m + A->rmap->n != nrow) { 208c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(nrow, std::distance(it + A->rmap->n, idxr + nrow), idxr, idxc + m + A->rmap->n, ptr + (m + A->rmap->n) * nrow); /* vertical block E from above */ 209c7a4214aSPierre Jolivet /* entry-wise assembly may be costly, so transpose already-computed entries when possible */ 210b94d7dedSBarry Smith if (A->symmetric == PETSC_BOOL3_TRUE || A->hermitian == PETSC_BOOL3_TRUE) { 2119566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), ptr + (m + A->rmap->n) * nrow + m, &B)); 2129566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, nrow)); 2139566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF, nrow - (m + A->rmap->n), A->rmap->n, nrow - (m + A->rmap->n), A->rmap->n, ptr + m * nrow + m + A->rmap->n, &BT)); 2149566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(BT, nrow)); 215b94d7dedSBarry Smith if (A->hermitian == PETSC_BOOL3_TRUE && PetscDefined(USE_COMPLEX)) { 2169566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(B, MAT_REUSE_MATRIX, &BT)); 217c7a4214aSPierre Jolivet } else { 2187fb60732SBarry Smith PetscCall(MatTransposeSetPrecursor(B, BT)); 2199566063dSJacob Faibussowitsch PetscCall(MatTranspose(B, MAT_REUSE_MATRIX, &BT)); 220c7a4214aSPierre Jolivet } 2219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BT)); 223c7a4214aSPierre Jolivet } else { 224c7a4214aSPierre Jolivet for (PetscInt k = 0; k < A->rmap->n; ++k) { /* block F from above */ 225c7a4214aSPierre Jolivet a->wrapper->copy_submatrix(std::distance(it + A->rmap->n, idxr + nrow), 1, it + A->rmap->n, idxc + m + k, ptr + (m + k) * nrow + m + A->rmap->n); 226c7a4214aSPierre Jolivet } 227c7a4214aSPierre Jolivet } 228c7a4214aSPierre Jolivet } 229c7a4214aSPierre Jolivet } /* complete local diagonal block not in IS */ 230c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* IS not long enough to store the local diagonal block */ 231c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* rmap->rstart not in IS */ 232c7a4214aSPierre Jolivet } /* unsorted IS */ 233c7a4214aSPierre Jolivet } 234c7a4214aSPierre Jolivet } else flg = PETSC_FALSE; /* different row and column IS */ 235c7a4214aSPierre Jolivet if (!flg) a->wrapper->copy_submatrix(nrow, m, idxr, idxc, ptr); /* reassemble everything */ 2369566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i], &idxr)); 2379566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i], &idxc)); 2389566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite((*submat)[i], &ptr)); 23966b4df27SPierre Jolivet PetscCall(MatShift((*submat)[i], shift)); 24066b4df27SPierre Jolivet PetscCall(MatScale((*submat)[i], scale)); 241c7a4214aSPierre Jolivet } 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 243c7a4214aSPierre Jolivet } 244c7a4214aSPierre Jolivet 245d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Htool(Mat A) 246d71ae5a4SJacob Faibussowitsch { 247c9a4fd69SAudic XU Mat_Htool *a; 248c7a4214aSPierre Jolivet PetscContainer container; 249c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 250c7a4214aSPierre Jolivet 251c7a4214aSPierre Jolivet PetscFunctionBegin; 252c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 253f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", nullptr)); 254f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", nullptr)); 255f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", nullptr)); 256f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", nullptr)); 257f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", nullptr)); 258f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", nullptr)); 259f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", nullptr)); 260f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", nullptr)); 261f22e26b7SPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", nullptr)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "KernelTranspose", (PetscObject *)&container)); 263c7a4214aSPierre Jolivet if (container) { /* created in MatTranspose_Htool() */ 2649566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 2659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 266f22e26b7SPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)A, "KernelTranspose", nullptr)); 267c7a4214aSPierre Jolivet } 26848a46eb9SPierre Jolivet if (a->gcoords_source != a->gcoords_target) PetscCall(PetscFree(a->gcoords_source)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFree(a->gcoords_target)); 2709566063dSJacob Faibussowitsch PetscCall(PetscFree2(a->work_source, a->work_target)); 271c7a4214aSPierre Jolivet delete a->wrapper; 2721dd4f53aSPierre Jolivet a->target_cluster.reset(); 2731dd4f53aSPierre Jolivet a->source_cluster.reset(); 2741dd4f53aSPierre Jolivet a->distributed_operator_holder.reset(); 275c9a4fd69SAudic XU PetscCall(PetscFree(a)); 276c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", nullptr)); // needed to avoid a call to MatShellSetContext_Immutable() 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278c7a4214aSPierre Jolivet } 279c7a4214aSPierre Jolivet 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Htool(Mat A, PetscViewer pv) 281d71ae5a4SJacob Faibussowitsch { 282c9a4fd69SAudic XU Mat_Htool *a; 28366b4df27SPierre Jolivet PetscScalar shift, scale; 284c7a4214aSPierre Jolivet PetscBool flg; 2851dd4f53aSPierre Jolivet std::map<std::string, std::string> hmatrix_information; 286c7a4214aSPierre Jolivet 287c7a4214aSPierre Jolivet PetscFunctionBegin; 288c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 2891dd4f53aSPierre Jolivet hmatrix_information = htool::get_distributed_hmatrix_information(a->distributed_operator_holder->hmatrix, PetscObjectComm((PetscObject)A)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pv, PETSCVIEWERASCII, &flg)); 291c7a4214aSPierre Jolivet if (flg) { 29266b4df27SPierre 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)); 2931dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "symmetry: %c\n", a->distributed_operator_holder->distributed_operator.get_symmetry_type())); 29466b4df27SPierre Jolivet if (PetscAbsScalar(scale - 1.0) > PETSC_MACHINE_EPSILON) { 295c7a4214aSPierre Jolivet #if defined(PETSC_USE_COMPLEX) 29666b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g+%gi\n", (double)PetscRealPart(scale), (double)PetscImaginaryPart(scale))); 297c7a4214aSPierre Jolivet #else 29866b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "scaling: %g\n", (double)scale)); 299c9a4fd69SAudic XU #endif 300c9a4fd69SAudic XU } 30166b4df27SPierre Jolivet if (PetscAbsScalar(shift) > PETSC_MACHINE_EPSILON) { 302c9a4fd69SAudic XU #if defined(PETSC_USE_COMPLEX) 30366b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g+%gi\n", (double)PetscRealPart(shift), (double)PetscImaginaryPart(shift))); 304c9a4fd69SAudic XU #else 30566b4df27SPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "shift: %g\n", (double)shift)); 306c7a4214aSPierre Jolivet #endif 307c7a4214aSPierre Jolivet } 3081dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "minimum cluster size: %" PetscInt_FMT "\n", a->min_cluster_size)); 3099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "epsilon: %g\n", (double)a->epsilon)); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "eta: %g\n", (double)a->eta)); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum target depth: %" PetscInt_FMT "\n", a->depth[0])); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "minimum source depth: %" PetscInt_FMT "\n", a->depth[1])); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "compressor: %s\n", MatHtoolCompressorTypes[a->compressor])); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(pv, "clustering: %s\n", MatHtoolClusteringTypes[a->clustering])); 3151dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "compression ratio: %s\n", hmatrix_information["Compression_ratio"].c_str())); 3161dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "space saving: %s\n", hmatrix_information["Space_saving"].c_str())); 3171dd4f53aSPierre Jolivet PetscCall(PetscViewerASCIIPrintf(pv, "block tree consistency: %s\n", PetscBools[a->distributed_operator_holder->hmatrix.is_block_tree_consistent()])); 3181dd4f53aSPierre 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())); 3191dd4f53aSPierre Jolivet PetscCall( 3201dd4f53aSPierre 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())); 3211dd4f53aSPierre 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(), 3221dd4f53aSPierre Jolivet hmatrix_information["Low_rank_block_size_max"].c_str())); 3231dd4f53aSPierre 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())); 324c7a4214aSPierre Jolivet } 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 326c7a4214aSPierre Jolivet } 327c7a4214aSPierre Jolivet 328c7a4214aSPierre Jolivet /* naive implementation of MatGetRow() needed for MatConvert_Nest_AIJ() */ 329d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_Htool(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 330d71ae5a4SJacob Faibussowitsch { 331c9a4fd69SAudic XU Mat_Htool *a; 33266b4df27SPierre Jolivet PetscScalar shift, scale; 333c7a4214aSPierre Jolivet PetscInt *idxc; 334c7a4214aSPierre Jolivet PetscBLASInt one = 1, bn; 335c7a4214aSPierre Jolivet 336c7a4214aSPierre Jolivet PetscFunctionBegin; 33766b4df27SPierre 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)); 338c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 339c7a4214aSPierre Jolivet if (nz) *nz = A->cmap->N; 340c7a4214aSPierre Jolivet if (idx || v) { /* even if !idx, need to set idxc for htool::copy_submatrix() */ 3419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, &idxc)); 342c7a4214aSPierre Jolivet for (PetscInt i = 0; i < A->cmap->N; ++i) idxc[i] = i; 343c7a4214aSPierre Jolivet } 344c7a4214aSPierre Jolivet if (idx) *idx = idxc; 345c7a4214aSPierre Jolivet if (v) { 3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->cmap->N, v)); 347c7a4214aSPierre Jolivet if (a->wrapper) a->wrapper->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 348cab00e0dSPierre Jolivet else reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx)->copy_submatrix(1, A->cmap->N, &row, idxc, *v); 3499566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(A->cmap->N, &bn)); 350*ebd42cecSPierre Jolivet PetscCallCXX(htool::Blas<PetscScalar>::scal(&bn, &scale, *v, &one)); 35166b4df27SPierre Jolivet if (row < A->cmap->N) (*v)[row] += shift; 352c7a4214aSPierre Jolivet } 35348a46eb9SPierre Jolivet if (!idx) PetscCall(PetscFree(idxc)); 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355c7a4214aSPierre Jolivet } 356c7a4214aSPierre Jolivet 357f9178db3SPierre Jolivet static PetscErrorCode MatRestoreRow_Htool(Mat, PetscInt, PetscInt *, PetscInt **idx, PetscScalar **v) 358d71ae5a4SJacob Faibussowitsch { 359c7a4214aSPierre Jolivet PetscFunctionBegin; 36048a46eb9SPierre Jolivet if (idx) PetscCall(PetscFree(*idx)); 36148a46eb9SPierre Jolivet if (v) PetscCall(PetscFree(*v)); 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 363c7a4214aSPierre Jolivet } 364c7a4214aSPierre Jolivet 365ce78bad3SBarry Smith static PetscErrorCode MatSetFromOptions_Htool(Mat A, PetscOptionItems PetscOptionsObject) 366d71ae5a4SJacob Faibussowitsch { 367c9a4fd69SAudic XU Mat_Htool *a; 368c7a4214aSPierre Jolivet PetscInt n; 369c7a4214aSPierre Jolivet PetscBool flg; 370c7a4214aSPierre Jolivet 371c7a4214aSPierre Jolivet PetscFunctionBegin; 372c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 373d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Htool options"); 3741dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedInt("-mat_htool_min_cluster_size", "Minimal leaf size in cluster tree", nullptr, a->min_cluster_size, &a->min_cluster_size, nullptr, 0)); 3751dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedReal("-mat_htool_epsilon", "Relative error in Frobenius norm when approximating a block", nullptr, a->epsilon, &a->epsilon, nullptr, 0.0)); 376f22e26b7SPierre Jolivet PetscCall(PetscOptionsReal("-mat_htool_eta", "Admissibility condition tolerance", nullptr, a->eta, &a->eta, nullptr)); 3771dd4f53aSPierre Jolivet PetscCall(PetscOptionsBoundedInt("-mat_htool_min_target_depth", "Minimal cluster tree depth associated with the rows", nullptr, a->depth[0], a->depth, nullptr, 0)); 3781dd4f53aSPierre 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)); 3791dd4f53aSPierre Jolivet PetscCall(PetscOptionsBool("-mat_htool_block_tree_consistency", "Block tree consistency", nullptr, a->block_tree_consistency, &a->block_tree_consistency, nullptr)); 3801dd4f53aSPierre Jolivet 381c7a4214aSPierre Jolivet n = 0; 382dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_compressor", "Type of compression", "MatHtoolCompressorType", MatHtoolCompressorTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolCompressorTypes), MatHtoolCompressorTypes[MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA], &n, &flg)); 383c7a4214aSPierre Jolivet if (flg) a->compressor = MatHtoolCompressorType(n); 38498e73e17SPierre Jolivet n = 0; 385dd39110bSPierre Jolivet PetscCall(PetscOptionsEList("-mat_htool_clustering", "Type of clustering", "MatHtoolClusteringType", MatHtoolClusteringTypes, PETSC_STATIC_ARRAY_LENGTH(MatHtoolClusteringTypes), MatHtoolClusteringTypes[MAT_HTOOL_CLUSTERING_PCA_REGULAR], &n, &flg)); 38698e73e17SPierre Jolivet if (flg) a->clustering = MatHtoolClusteringType(n); 387d0609cedSBarry Smith PetscOptionsHeadEnd(); 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 389c7a4214aSPierre Jolivet } 390c7a4214aSPierre Jolivet 391cedd07caSPierre Jolivet static PetscErrorCode MatAssemblyEnd_Htool(Mat A, MatAssemblyType) 392d71ae5a4SJacob Faibussowitsch { 393c9a4fd69SAudic XU Mat_Htool *a; 394c7a4214aSPierre Jolivet const PetscInt *ranges; 395c7a4214aSPierre Jolivet PetscInt *offset; 3961dd4f53aSPierre Jolivet PetscMPIInt size, rank; 397b94d7dedSBarry Smith char S = PetscDefined(USE_COMPLEX) && A->hermitian == PETSC_BOOL3_TRUE ? 'H' : (A->symmetric == PETSC_BOOL3_TRUE ? 'S' : 'N'), uplo = S == 'N' ? 'N' : 'U'; 398cab00e0dSPierre Jolivet htool::VirtualGenerator<PetscScalar> *generator = nullptr; 3991dd4f53aSPierre Jolivet htool::ClusterTreeBuilder<PetscReal> recursive_build_strategy; 4001dd4f53aSPierre Jolivet htool::Cluster<PetscReal> *source_cluster; 4011dd4f53aSPierre Jolivet std::shared_ptr<htool::VirtualLowRankGenerator<PetscScalar>> compressor; 402c7a4214aSPierre Jolivet 403c7a4214aSPierre Jolivet PetscFunctionBegin; 4049566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(HtoolCitation, &HtoolCite)); 405c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 406c7a4214aSPierre Jolivet delete a->wrapper; 4071dd4f53aSPierre Jolivet a->target_cluster.reset(); 4081dd4f53aSPierre Jolivet a->source_cluster.reset(); 4091dd4f53aSPierre Jolivet a->distributed_operator_holder.reset(); 4101dd4f53aSPierre Jolivet // clustering 4119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 4129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &offset)); 4139566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 414c7a4214aSPierre Jolivet for (PetscInt i = 0; i < size; ++i) { 415c7a4214aSPierre Jolivet offset[2 * i] = ranges[i]; 416c7a4214aSPierre Jolivet offset[2 * i + 1] = ranges[i + 1] - ranges[i]; 417c7a4214aSPierre Jolivet } 41898e73e17SPierre Jolivet switch (a->clustering) { 419d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC: 4201dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>()); 4211dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>()); 422d71ae5a4SJacob Faibussowitsch break; 423d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 4241dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>()); 4251dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>()); 426d71ae5a4SJacob Faibussowitsch break; 427d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 4281dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>()); 4291dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>()); 430d71ae5a4SJacob Faibussowitsch break; 431d71ae5a4SJacob Faibussowitsch default: 4321dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>()); 4331dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>()); 43498e73e17SPierre Jolivet } 4351dd4f53aSPierre Jolivet recursive_build_strategy.set_minclustersize(a->min_cluster_size); 4361dd4f53aSPierre Jolivet a->target_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree(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: 4451dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>()); 4461dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>()); 447d71ae5a4SJacob Faibussowitsch break; 448d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC: 4491dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>()); 4501dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::GeometricSplitting<PetscReal>>()); 451d71ae5a4SJacob Faibussowitsch break; 452d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR: 4531dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeBoundingBox<PetscReal>>()); 4541dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>()); 455d71ae5a4SJacob Faibussowitsch break; 456d71ae5a4SJacob Faibussowitsch default: 4571dd4f53aSPierre Jolivet recursive_build_strategy.set_direction_computation_strategy(std::make_shared<htool::ComputeLargestExtent<PetscReal>>()); 4581dd4f53aSPierre Jolivet recursive_build_strategy.set_splitting_strategy(std::make_shared<htool::RegularSplitting<PetscReal>>()); 45998e73e17SPierre Jolivet } 4601dd4f53aSPierre Jolivet recursive_build_strategy.set_minclustersize(a->min_cluster_size); 4611dd4f53aSPierre Jolivet a->source_cluster = std::make_unique<htool::Cluster<PetscReal>>(recursive_build_strategy.create_cluster_tree(A->cmap->N, a->dim, a->gcoords_source, 2, size, offset)); 462c7a4214aSPierre Jolivet S = uplo = 'N'; 4631dd4f53aSPierre Jolivet source_cluster = a->source_cluster.get(); 4641dd4f53aSPierre Jolivet } else source_cluster = a->target_cluster.get(); 4659566063dSJacob Faibussowitsch PetscCall(PetscFree(offset)); 4661dd4f53aSPierre Jolivet // generator 4671dd4f53aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx); 4681dd4f53aSPierre Jolivet else { 4691dd4f53aSPierre Jolivet a->wrapper = nullptr; 4701dd4f53aSPierre Jolivet generator = reinterpret_cast<htool::VirtualGenerator<PetscScalar> *>(a->kernelctx); 4711dd4f53aSPierre Jolivet } 4721dd4f53aSPierre Jolivet // compressor 473c7a4214aSPierre Jolivet switch (a->compressor) { 474d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_FULL_ACA: 475d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::fullACA<PetscScalar>>(); 476d71ae5a4SJacob Faibussowitsch break; 477d71ae5a4SJacob Faibussowitsch case MAT_HTOOL_COMPRESSOR_SVD: 478d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::SVD<PetscScalar>>(); 479d71ae5a4SJacob Faibussowitsch break; 480d71ae5a4SJacob Faibussowitsch default: 481d71ae5a4SJacob Faibussowitsch compressor = std::make_shared<htool::sympartialACA<PetscScalar>>(); 482c7a4214aSPierre Jolivet } 4831dd4f53aSPierre Jolivet // local hierarchical matrix 4841dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 4851dd4f53aSPierre Jolivet auto hmatrix_builder = htool::HMatrixTreeBuilder<PetscScalar>(*a->target_cluster, *source_cluster, a->epsilon, a->eta, S, uplo, -1, rank, rank); 4861dd4f53aSPierre Jolivet hmatrix_builder.set_low_rank_generator(compressor); 4871dd4f53aSPierre Jolivet hmatrix_builder.set_minimal_target_depth(a->depth[0]); 4881dd4f53aSPierre Jolivet hmatrix_builder.set_minimal_source_depth(a->depth[1]); 4891dd4f53aSPierre 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"); 4901dd4f53aSPierre Jolivet hmatrix_builder.set_block_tree_consistency(a->block_tree_consistency); 4911dd4f53aSPierre Jolivet a->distributed_operator_holder = std::make_unique<htool::DistributedOperatorFromHMatrix<PetscScalar>>(a->wrapper ? *a->wrapper : *generator, *a->target_cluster, *source_cluster, hmatrix_builder, PetscObjectComm((PetscObject)A)); 4923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 493c7a4214aSPierre Jolivet } 494c7a4214aSPierre Jolivet 495d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_Htool(Mat C) 496d71ae5a4SJacob Faibussowitsch { 497c7a4214aSPierre Jolivet Mat_Product *product = C->product; 498c9a4fd69SAudic XU Mat_Htool *a; 499c7a4214aSPierre Jolivet const PetscScalar *in; 500c7a4214aSPierre Jolivet PetscScalar *out; 5018b8ddd11SPierre Jolivet PetscInt N, lda; 502c7a4214aSPierre Jolivet 503c7a4214aSPierre Jolivet PetscFunctionBegin; 504c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 505f22e26b7SPierre Jolivet PetscCall(MatGetSize(C, nullptr, &N)); 5069566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 5075f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 5089566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(product->B, &in)); 5099566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &out)); 510c9a4fd69SAudic XU PetscCall(MatShellGetContext(product->A, &a)); 5118b8ddd11SPierre Jolivet switch (product->type) { 512d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AB: 5131dd4f53aSPierre Jolivet a->distributed_operator_holder->distributed_operator.matrix_product_local_to_local(in, out, N, nullptr); 514d71ae5a4SJacob Faibussowitsch break; 515d71ae5a4SJacob Faibussowitsch case MATPRODUCT_AtB: 5161dd4f53aSPierre Jolivet a->distributed_operator_holder->distributed_operator.matrix_product_transp_local_to_local(in, out, N, nullptr); 517d71ae5a4SJacob Faibussowitsch break; 518d71ae5a4SJacob Faibussowitsch default: 519d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType %s is not supported", MatProductTypes[product->type]); 520c7a4214aSPierre Jolivet } 5219566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &out)); 5229566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(product->B, &in)); 5233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 524c7a4214aSPierre Jolivet } 525c7a4214aSPierre Jolivet 526d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_Htool(Mat C) 527d71ae5a4SJacob Faibussowitsch { 528c7a4214aSPierre Jolivet Mat_Product *product = C->product; 529c7a4214aSPierre Jolivet Mat A, B; 530c7a4214aSPierre Jolivet PetscBool flg; 531c7a4214aSPierre Jolivet 532c7a4214aSPierre Jolivet PetscFunctionBegin; 533c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 534c7a4214aSPierre Jolivet A = product->A; 535c7a4214aSPierre Jolivet B = product->B; 5369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, "")); 53797713f8cSPierre 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); 53897713f8cSPierre Jolivet if (C->rmap->n == PETSC_DECIDE || C->cmap->n == PETSC_DECIDE || C->rmap->N == PETSC_DECIDE || C->cmap->N == PETSC_DECIDE) { 53997713f8cSPierre Jolivet if (product->type == MATPRODUCT_AB) PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N)); 54097713f8cSPierre Jolivet else PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N)); 5418b8ddd11SPierre Jolivet } 5429566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 5439566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 5449566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 5459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 5469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 547f22e26b7SPierre Jolivet C->ops->productsymbolic = nullptr; 548c7a4214aSPierre Jolivet C->ops->productnumeric = MatProductNumeric_Htool; 5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 550c7a4214aSPierre Jolivet } 551c7a4214aSPierre Jolivet 552d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Htool(Mat C) 553d71ae5a4SJacob Faibussowitsch { 554c7a4214aSPierre Jolivet PetscFunctionBegin; 555c7a4214aSPierre Jolivet MatCheckProduct(C, 1); 5568b8ddd11SPierre Jolivet if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_Htool; 5573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 558c7a4214aSPierre Jolivet } 559c7a4214aSPierre Jolivet 5601dd4f53aSPierre Jolivet static PetscErrorCode MatHtoolGetHierarchicalMat_Htool(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator) 561d71ae5a4SJacob Faibussowitsch { 562c9a4fd69SAudic XU Mat_Htool *a; 563c7a4214aSPierre Jolivet 564c7a4214aSPierre Jolivet PetscFunctionBegin; 565c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 5661dd4f53aSPierre Jolivet *distributed_operator = &a->distributed_operator_holder->distributed_operator; 5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 568c7a4214aSPierre Jolivet } 569c7a4214aSPierre Jolivet 570c7a4214aSPierre Jolivet /*@C 57111a5261eSBarry Smith MatHtoolGetHierarchicalMat - Retrieves the opaque pointer to a Htool virtual matrix stored in a `MATHTOOL`. 572c7a4214aSPierre Jolivet 573cc4c1da9SBarry Smith No Fortran Support, No C Support 574cc4c1da9SBarry Smith 575c7a4214aSPierre Jolivet Input Parameter: 576c7a4214aSPierre Jolivet . A - hierarchical matrix 577c7a4214aSPierre Jolivet 578c7a4214aSPierre Jolivet Output Parameter: 5791dd4f53aSPierre Jolivet . distributed_operator - opaque pointer to a Htool virtual matrix 580c7a4214aSPierre Jolivet 581c7a4214aSPierre Jolivet Level: advanced 582c7a4214aSPierre Jolivet 5831cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL` 584c7a4214aSPierre Jolivet @*/ 5851dd4f53aSPierre Jolivet PETSC_EXTERN PetscErrorCode MatHtoolGetHierarchicalMat(Mat A, const htool::DistributedOperator<PetscScalar> **distributed_operator) 586d71ae5a4SJacob Faibussowitsch { 587c7a4214aSPierre Jolivet PetscFunctionBegin; 588c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5891dd4f53aSPierre Jolivet PetscAssertPointer(distributed_operator, 2); 5901dd4f53aSPierre Jolivet PetscTryMethod(A, "MatHtoolGetHierarchicalMat_C", (Mat, const htool::DistributedOperator<PetscScalar> **), (A, distributed_operator)); 5913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 592c7a4214aSPierre Jolivet } 593c7a4214aSPierre Jolivet 5948434afd1SBarry Smith static PetscErrorCode MatHtoolSetKernel_Htool(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 595d71ae5a4SJacob Faibussowitsch { 596c9a4fd69SAudic XU Mat_Htool *a; 597c7a4214aSPierre Jolivet 598c7a4214aSPierre Jolivet PetscFunctionBegin; 599c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 600c7a4214aSPierre Jolivet a->kernel = kernel; 601c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 602c7a4214aSPierre Jolivet delete a->wrapper; 6031dd4f53aSPierre Jolivet if (a->kernel) a->wrapper = new WrapperHtool(a->dim, a->kernel, a->kernelctx); 6043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 605c7a4214aSPierre Jolivet } 606c7a4214aSPierre Jolivet 607c7a4214aSPierre Jolivet /*@C 60811a5261eSBarry Smith MatHtoolSetKernel - Sets the kernel and context used for the assembly of a `MATHTOOL`. 609c7a4214aSPierre Jolivet 610cc4c1da9SBarry Smith Collective, No Fortran Support 611cc4c1da9SBarry Smith 612c7a4214aSPierre Jolivet Input Parameters: 613c7a4214aSPierre Jolivet + A - hierarchical matrix 6142ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 6152ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 616c7a4214aSPierre Jolivet 617c7a4214aSPierre Jolivet Level: advanced 618c7a4214aSPierre Jolivet 6191cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatCreateHtoolFromKernel()` 620c7a4214aSPierre Jolivet @*/ 621cc4c1da9SBarry Smith PetscErrorCode MatHtoolSetKernel(Mat A, MatHtoolKernelFn *kernel, void *kernelctx) 622d71ae5a4SJacob Faibussowitsch { 623c7a4214aSPierre Jolivet PetscFunctionBegin; 624c7a4214aSPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 625c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 2); 6264f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 3); 6278434afd1SBarry Smith PetscTryMethod(A, "MatHtoolSetKernel_C", (Mat, MatHtoolKernelFn *, void *), (A, kernel, kernelctx)); 6283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 629c7a4214aSPierre Jolivet } 630c7a4214aSPierre Jolivet 631d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationSource_Htool(Mat A, IS *is) 632d71ae5a4SJacob Faibussowitsch { 633c9a4fd69SAudic XU Mat_Htool *a; 6341dd4f53aSPierre Jolivet PetscMPIInt rank; 6351dd4f53aSPierre Jolivet const std::vector<PetscInt> *source; 6361dd4f53aSPierre Jolivet const htool::Cluster<PetscReal> *local_source_cluster; 63798e73e17SPierre Jolivet 63898e73e17SPierre Jolivet PetscFunctionBegin; 639c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 6401dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 6411dd4f53aSPierre Jolivet local_source_cluster = a->source_cluster ? &a->source_cluster->get_cluster_on_partition(rank) : &a->target_cluster->get_cluster_on_partition(rank); 6421dd4f53aSPierre Jolivet source = &local_source_cluster->get_permutation(); 6431dd4f53aSPierre Jolivet PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)A), local_source_cluster->get_size(), source->data() + local_source_cluster->get_offset(), PETSC_COPY_VALUES, is)); 6449566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64698e73e17SPierre Jolivet } 64798e73e17SPierre Jolivet 648cc4c1da9SBarry Smith /*@ 64911a5261eSBarry Smith MatHtoolGetPermutationSource - Gets the permutation associated to the source cluster for a `MATHTOOL` matrix. 65098e73e17SPierre Jolivet 65198e73e17SPierre Jolivet Input Parameter: 65298e73e17SPierre Jolivet . A - hierarchical matrix 65398e73e17SPierre Jolivet 65498e73e17SPierre Jolivet Output Parameter: 65598e73e17SPierre Jolivet . is - permutation 65698e73e17SPierre Jolivet 65798e73e17SPierre Jolivet Level: advanced 65898e73e17SPierre Jolivet 6591cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationTarget()`, `MatHtoolUsePermutation()` 66098e73e17SPierre Jolivet @*/ 661cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationSource(Mat A, IS *is) 662d71ae5a4SJacob Faibussowitsch { 66398e73e17SPierre Jolivet PetscFunctionBegin; 66498e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6654f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 666cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationSource_C", (Mat, IS *), (A, is)); 6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66898e73e17SPierre Jolivet } 66998e73e17SPierre Jolivet 670d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolGetPermutationTarget_Htool(Mat A, IS *is) 671d71ae5a4SJacob Faibussowitsch { 672c9a4fd69SAudic XU Mat_Htool *a; 6731dd4f53aSPierre Jolivet const std::vector<PetscInt> *target; 6741dd4f53aSPierre Jolivet PetscMPIInt rank; 67598e73e17SPierre Jolivet 67698e73e17SPierre Jolivet PetscFunctionBegin; 677c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 6781dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 6791dd4f53aSPierre Jolivet target = &a->target_cluster->get_permutation(); 6801dd4f53aSPierre 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)); 6819566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*is)); 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68398e73e17SPierre Jolivet } 68498e73e17SPierre Jolivet 685cc4c1da9SBarry Smith /*@ 68611a5261eSBarry Smith MatHtoolGetPermutationTarget - Gets the permutation associated to the target cluster for a `MATHTOOL` matrix. 68798e73e17SPierre Jolivet 68898e73e17SPierre Jolivet Input Parameter: 68998e73e17SPierre Jolivet . A - hierarchical matrix 69098e73e17SPierre Jolivet 69198e73e17SPierre Jolivet Output Parameter: 69298e73e17SPierre Jolivet . is - permutation 69398e73e17SPierre Jolivet 69498e73e17SPierre Jolivet Level: advanced 69598e73e17SPierre Jolivet 6961cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolUsePermutation()` 69798e73e17SPierre Jolivet @*/ 698cc4c1da9SBarry Smith PetscErrorCode MatHtoolGetPermutationTarget(Mat A, IS *is) 699d71ae5a4SJacob Faibussowitsch { 70098e73e17SPierre Jolivet PetscFunctionBegin; 70198e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 7024f572ea9SToby Isaac if (!is) PetscAssertPointer(is, 2); 703cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolGetPermutationTarget_C", (Mat, IS *), (A, is)); 7043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70598e73e17SPierre Jolivet } 70698e73e17SPierre Jolivet 707d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatHtoolUsePermutation_Htool(Mat A, PetscBool use) 708d71ae5a4SJacob Faibussowitsch { 709c9a4fd69SAudic XU Mat_Htool *a; 71098e73e17SPierre Jolivet 71198e73e17SPierre Jolivet PetscFunctionBegin; 712c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 7131dd4f53aSPierre Jolivet a->distributed_operator_holder->distributed_operator.use_permutation() = use; 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71598e73e17SPierre Jolivet } 71698e73e17SPierre Jolivet 717cc4c1da9SBarry Smith /*@ 71811a5261eSBarry Smith MatHtoolUsePermutation - Sets whether a `MATHTOOL` matrix should permute input (resp. output) vectors following its internal source (resp. target) permutation. 71998e73e17SPierre Jolivet 72098e73e17SPierre Jolivet Input Parameters: 72198e73e17SPierre Jolivet + A - hierarchical matrix 72298e73e17SPierre Jolivet - use - Boolean value 72398e73e17SPierre Jolivet 72498e73e17SPierre Jolivet Level: advanced 72598e73e17SPierre Jolivet 7261cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATHTOOL`, `MatHtoolGetPermutationSource()`, `MatHtoolGetPermutationTarget()` 72798e73e17SPierre Jolivet @*/ 728cc4c1da9SBarry Smith PetscErrorCode MatHtoolUsePermutation(Mat A, PetscBool use) 729d71ae5a4SJacob Faibussowitsch { 73098e73e17SPierre Jolivet PetscFunctionBegin; 73198e73e17SPierre Jolivet PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 73298e73e17SPierre Jolivet PetscValidLogicalCollectiveBool(A, use, 2); 733cac4c232SBarry Smith PetscTryMethod(A, "MatHtoolUsePermutation_C", (Mat, PetscBool), (A, use)); 7343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73598e73e17SPierre Jolivet } 73698e73e17SPierre Jolivet 737cedd07caSPierre Jolivet static PetscErrorCode MatConvert_Htool_Dense(Mat A, MatType, MatReuse reuse, Mat *B) 738d71ae5a4SJacob Faibussowitsch { 739c7a4214aSPierre Jolivet Mat C; 740c9a4fd69SAudic XU Mat_Htool *a; 74166b4df27SPierre Jolivet PetscScalar *array, shift, scale; 742c7a4214aSPierre Jolivet PetscInt lda; 743c7a4214aSPierre Jolivet 744c7a4214aSPierre Jolivet PetscFunctionBegin; 74566b4df27SPierre 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)); 746c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 747c7a4214aSPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 748c7a4214aSPierre Jolivet C = *B; 7495f80ce2aSJacob Faibussowitsch PetscCheck(C->rmap->n == A->rmap->n && C->cmap->N == A->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible dimensions"); 7509566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &lda)); 7515f80ce2aSJacob Faibussowitsch PetscCheck(lda == C->rmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported leading dimension (%" PetscInt_FMT " != %" PetscInt_FMT ")", lda, C->rmap->n); 752c7a4214aSPierre Jolivet } else { 7539566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 7549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 7559566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 7569566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 7579566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 758c7a4214aSPierre Jolivet } 7599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 7609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 7619566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &array)); 7621dd4f53aSPierre Jolivet htool::copy_to_dense_in_user_numbering(a->distributed_operator_holder->hmatrix, array); 7639566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C, &array)); 76466b4df27SPierre Jolivet PetscCall(MatShift(C, shift)); 76566b4df27SPierre Jolivet PetscCall(MatScale(C, scale)); 766c7a4214aSPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 7679566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &C)); 768c7a4214aSPierre Jolivet } else *B = C; 7693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 770c7a4214aSPierre Jolivet } 771c7a4214aSPierre Jolivet 772d71ae5a4SJacob Faibussowitsch static PetscErrorCode GenEntriesTranspose(PetscInt sdim, PetscInt M, PetscInt N, const PetscInt *rows, const PetscInt *cols, PetscScalar *ptr, void *ctx) 773d71ae5a4SJacob Faibussowitsch { 774c7a4214aSPierre Jolivet MatHtoolKernelTranspose *generator = (MatHtoolKernelTranspose *)ctx; 77598e73e17SPierre Jolivet PetscScalar *tmp; 77698e73e17SPierre Jolivet 77798e73e17SPierre Jolivet PetscFunctionBegin; 7783ba16761SJacob Faibussowitsch PetscCall(generator->kernel(sdim, N, M, cols, rows, ptr, generator->kernelctx)); 7799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * N, &tmp)); 7809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, ptr, M * N)); 78198e73e17SPierre Jolivet for (PetscInt i = 0; i < M; ++i) { 78298e73e17SPierre Jolivet for (PetscInt j = 0; j < N; ++j) ptr[i + j * M] = tmp[j + i * N]; 78398e73e17SPierre Jolivet } 7849566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 7853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 786c7a4214aSPierre Jolivet } 787c7a4214aSPierre Jolivet 788c7a4214aSPierre Jolivet /* naive implementation which keeps a reference to the original Mat */ 789d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Htool(Mat A, MatReuse reuse, Mat *B) 790d71ae5a4SJacob Faibussowitsch { 791c7a4214aSPierre Jolivet Mat C; 792c9a4fd69SAudic XU Mat_Htool *a, *c; 79366b4df27SPierre Jolivet PetscScalar shift, scale; 794c7a4214aSPierre Jolivet PetscInt M = A->rmap->N, N = A->cmap->N, m = A->rmap->n, n = A->cmap->n; 795c7a4214aSPierre Jolivet PetscContainer container; 796c7a4214aSPierre Jolivet MatHtoolKernelTranspose *kernelt; 797c7a4214aSPierre Jolivet 798c7a4214aSPierre Jolivet PetscFunctionBegin; 79966b4df27SPierre 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)); 800c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 8017fb60732SBarry Smith if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B)); 8025f80ce2aSJacob Faibussowitsch PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatTranspose() with MAT_INPLACE_MATRIX not supported"); 803c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8049566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C)); 8059566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, n, m, N, M)); 8069566063dSJacob Faibussowitsch PetscCall(MatSetType(C, ((PetscObject)A)->type_name)); 8079566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 8089566063dSJacob Faibussowitsch PetscCall(PetscNew(&kernelt)); 80949abdd8aSBarry Smith PetscCall(PetscObjectContainerCompose((PetscObject)C, "KernelTranspose", kernelt, PetscCtxDestroyDefault)); 810c7a4214aSPierre Jolivet } else { 811c7a4214aSPierre Jolivet C = *B; 8129566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)C, "KernelTranspose", (PetscObject *)&container)); 8135f80ce2aSJacob Faibussowitsch PetscCheck(container, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatTranspose() with MAT_INITIAL_MATRIX first"); 8149566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(container, (void **)&kernelt)); 815c7a4214aSPierre Jolivet } 816c9a4fd69SAudic XU PetscCall(MatShellGetContext(C, &c)); 817c7a4214aSPierre Jolivet c->dim = a->dim; 81866b4df27SPierre Jolivet PetscCall(MatShift(C, shift)); 81966b4df27SPierre Jolivet PetscCall(MatScale(C, scale)); 82098e73e17SPierre Jolivet c->kernel = GenEntriesTranspose; 821c7a4214aSPierre Jolivet if (kernelt->A != A) { 8229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&kernelt->A)); 823c7a4214aSPierre Jolivet kernelt->A = A; 8249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 825c7a4214aSPierre Jolivet } 826c7a4214aSPierre Jolivet kernelt->kernel = a->kernel; 827c7a4214aSPierre Jolivet kernelt->kernelctx = a->kernelctx; 828c7a4214aSPierre Jolivet c->kernelctx = kernelt; 829c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) { 8309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * c->dim, &c->gcoords_target)); 8319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_target, a->gcoords_source, N * c->dim)); 832c7a4214aSPierre Jolivet if (a->gcoords_target != a->gcoords_source) { 8339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M * c->dim, &c->gcoords_source)); 8349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(c->gcoords_source, a->gcoords_target, M * c->dim)); 835c7a4214aSPierre Jolivet } else c->gcoords_source = c->gcoords_target; 8369566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(M, &c->work_source, N, &c->work_target)); 837c7a4214aSPierre Jolivet } 8389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 8399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 840c7a4214aSPierre Jolivet if (reuse == MAT_INITIAL_MATRIX) *B = C; 8413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 842c7a4214aSPierre Jolivet } 843c7a4214aSPierre Jolivet 8441dd4f53aSPierre Jolivet static PetscErrorCode MatDestroy_Factor(Mat F) 8451dd4f53aSPierre Jolivet { 8461dd4f53aSPierre Jolivet PetscContainer container; 8471dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *A; 8481dd4f53aSPierre Jolivet 8491dd4f53aSPierre Jolivet PetscFunctionBegin; 8501dd4f53aSPierre Jolivet PetscCall(PetscObjectQuery((PetscObject)F, "HMatrix", (PetscObject *)&container)); 8511dd4f53aSPierre Jolivet if (container) { 8521dd4f53aSPierre Jolivet PetscCall(PetscContainerGetPointer(container, (void **)&A)); 8531dd4f53aSPierre Jolivet delete A; 8541dd4f53aSPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)F, "HMatrix", nullptr)); 8551dd4f53aSPierre Jolivet } 8561dd4f53aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", nullptr)); 8571dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 8581dd4f53aSPierre Jolivet } 8591dd4f53aSPierre Jolivet 8601dd4f53aSPierre Jolivet static PetscErrorCode MatFactorGetSolverType_Htool(Mat, MatSolverType *type) 8611dd4f53aSPierre Jolivet { 8621dd4f53aSPierre Jolivet PetscFunctionBegin; 8631dd4f53aSPierre Jolivet *type = MATSOLVERHTOOL; 8641dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 8651dd4f53aSPierre Jolivet } 8661dd4f53aSPierre Jolivet 8671dd4f53aSPierre Jolivet template <char trans> 8681dd4f53aSPierre Jolivet static inline PetscErrorCode MatSolve_Private(Mat A, htool::Matrix<PetscScalar> &X) 8691dd4f53aSPierre Jolivet { 8701dd4f53aSPierre Jolivet PetscContainer container; 8711dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *B; 8721dd4f53aSPierre Jolivet 8731dd4f53aSPierre Jolivet PetscFunctionBegin; 8741dd4f53aSPierre 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"); 8751dd4f53aSPierre Jolivet PetscCall(PetscObjectQuery((PetscObject)A, "HMatrix", (PetscObject *)&container)); 8761dd4f53aSPierre 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"); 8771dd4f53aSPierre Jolivet PetscCall(PetscContainerGetPointer(container, (void **)&B)); 8781dd4f53aSPierre Jolivet if (A->factortype == MAT_FACTOR_LU) htool::lu_solve(trans, *B, X); 8791dd4f53aSPierre Jolivet else htool::cholesky_solve('L', *B, X); 8801dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 8811dd4f53aSPierre Jolivet } 8821dd4f53aSPierre Jolivet 8831dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Vec>::value>::type * = nullptr> 8841dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type b, Type x) 8851dd4f53aSPierre Jolivet { 8861dd4f53aSPierre Jolivet PetscInt n; 8871dd4f53aSPierre Jolivet htool::Matrix<PetscScalar> v; 8881dd4f53aSPierre Jolivet PetscScalar *array; 8891dd4f53aSPierre Jolivet 8901dd4f53aSPierre Jolivet PetscFunctionBegin; 8911dd4f53aSPierre Jolivet PetscCall(VecGetLocalSize(b, &n)); 8921dd4f53aSPierre Jolivet PetscCall(VecCopy(b, x)); 8931dd4f53aSPierre Jolivet PetscCall(VecGetArrayWrite(x, &array)); 8941dd4f53aSPierre Jolivet v.assign(n, 1, array, false); 8951dd4f53aSPierre Jolivet PetscCall(VecRestoreArrayWrite(x, &array)); 8961dd4f53aSPierre Jolivet PetscCall(MatSolve_Private<trans>(A, v)); 8971dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 8981dd4f53aSPierre Jolivet } 8991dd4f53aSPierre Jolivet 9001dd4f53aSPierre Jolivet template <char trans, class Type, typename std::enable_if<std::is_same<Type, Mat>::value>::type * = nullptr> 9011dd4f53aSPierre Jolivet static PetscErrorCode MatSolve_Htool(Mat A, Type B, Type X) 9021dd4f53aSPierre Jolivet { 9031dd4f53aSPierre Jolivet PetscInt m, N; 9041dd4f53aSPierre Jolivet htool::Matrix<PetscScalar> v; 9051dd4f53aSPierre Jolivet PetscScalar *array; 9061dd4f53aSPierre Jolivet 9071dd4f53aSPierre Jolivet PetscFunctionBegin; 9081dd4f53aSPierre Jolivet PetscCall(MatGetLocalSize(B, &m, nullptr)); 9091dd4f53aSPierre Jolivet PetscCall(MatGetLocalSize(B, nullptr, &N)); 9101dd4f53aSPierre Jolivet PetscCall(MatCopy(B, X, SAME_NONZERO_PATTERN)); 9111dd4f53aSPierre Jolivet PetscCall(MatDenseGetArrayWrite(X, &array)); 9121dd4f53aSPierre Jolivet v.assign(m, N, array, false); 9131dd4f53aSPierre Jolivet PetscCall(MatDenseRestoreArrayWrite(X, &array)); 9141dd4f53aSPierre Jolivet PetscCall(MatSolve_Private<trans>(A, v)); 9151dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9161dd4f53aSPierre Jolivet } 9171dd4f53aSPierre Jolivet 9181dd4f53aSPierre Jolivet template <MatFactorType ftype> 9191dd4f53aSPierre Jolivet static PetscErrorCode MatFactorNumeric_Htool(Mat F, Mat A, const MatFactorInfo *) 9201dd4f53aSPierre Jolivet { 9211dd4f53aSPierre Jolivet Mat_Htool *a; 9221dd4f53aSPierre Jolivet htool::HMatrix<PetscScalar> *B; 9231dd4f53aSPierre Jolivet 9241dd4f53aSPierre Jolivet PetscFunctionBegin; 9251dd4f53aSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 9261dd4f53aSPierre Jolivet B = new htool::HMatrix<PetscScalar>(a->distributed_operator_holder->hmatrix); 9271dd4f53aSPierre Jolivet if (ftype == MAT_FACTOR_LU) htool::lu_factorization(*B); 9281dd4f53aSPierre Jolivet else htool::cholesky_factorization('L', *B); 9291dd4f53aSPierre Jolivet PetscCall(PetscObjectContainerCompose((PetscObject)F, "HMatrix", B, nullptr)); 9301dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9311dd4f53aSPierre Jolivet } 9321dd4f53aSPierre Jolivet 9331dd4f53aSPierre Jolivet template <MatFactorType ftype> 9341dd4f53aSPierre Jolivet PetscErrorCode MatFactorSymbolic_Htool(Mat F, Mat) 9351dd4f53aSPierre Jolivet { 9361dd4f53aSPierre Jolivet PetscFunctionBegin; 9371dd4f53aSPierre Jolivet F->preallocated = PETSC_TRUE; 9381dd4f53aSPierre Jolivet F->assembled = PETSC_TRUE; 9391dd4f53aSPierre Jolivet F->ops->solve = MatSolve_Htool<'N', Vec>; 9401dd4f53aSPierre Jolivet F->ops->matsolve = MatSolve_Htool<'N', Mat>; 9411dd4f53aSPierre Jolivet if (!PetscDefined(USE_COMPLEX) || ftype == MAT_FACTOR_LU) { 9421dd4f53aSPierre Jolivet F->ops->solvetranspose = MatSolve_Htool<'T', Vec>; 9431dd4f53aSPierre Jolivet F->ops->matsolvetranspose = MatSolve_Htool<'T', Mat>; 9441dd4f53aSPierre Jolivet } 9451dd4f53aSPierre Jolivet F->ops->destroy = MatDestroy_Factor; 9461dd4f53aSPierre Jolivet if (ftype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_LU>; 9471dd4f53aSPierre Jolivet else F->ops->choleskyfactornumeric = MatFactorNumeric_Htool<MAT_FACTOR_CHOLESKY>; 9481dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9491dd4f53aSPierre Jolivet } 9501dd4f53aSPierre Jolivet 9511dd4f53aSPierre Jolivet static PetscErrorCode MatLUFactorSymbolic_Htool(Mat F, Mat A, IS, IS, const MatFactorInfo *) 9521dd4f53aSPierre Jolivet { 9531dd4f53aSPierre Jolivet PetscFunctionBegin; 9541dd4f53aSPierre Jolivet PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_LU>(F, A)); 9551dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9561dd4f53aSPierre Jolivet } 9571dd4f53aSPierre Jolivet 9581dd4f53aSPierre Jolivet static PetscErrorCode MatCholeskyFactorSymbolic_Htool(Mat F, Mat A, IS, const MatFactorInfo *) 9591dd4f53aSPierre Jolivet { 9601dd4f53aSPierre Jolivet PetscFunctionBegin; 9611dd4f53aSPierre Jolivet PetscCall(MatFactorSymbolic_Htool<MAT_FACTOR_CHOLESKY>(F, A)); 9621dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9631dd4f53aSPierre Jolivet } 9641dd4f53aSPierre Jolivet 9651dd4f53aSPierre Jolivet static PetscErrorCode MatGetFactor_htool_htool(Mat A, MatFactorType ftype, Mat *F) 9661dd4f53aSPierre Jolivet { 9671dd4f53aSPierre Jolivet Mat B; 9681dd4f53aSPierre Jolivet Mat_Htool *a; 9691dd4f53aSPierre Jolivet PetscMPIInt size; 9701dd4f53aSPierre Jolivet 9711dd4f53aSPierre Jolivet PetscFunctionBegin; 9721dd4f53aSPierre Jolivet PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 9731dd4f53aSPierre Jolivet PetscCall(MatShellGetContext(A, &a)); 9741dd4f53aSPierre Jolivet PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_WRONG_MPI_SIZE, "Unsupported parallel MatGetFactor()"); 9751dd4f53aSPierre Jolivet PetscCheck(a->block_tree_consistency, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot factor a MatHtool with inconsistent block tree"); 9761dd4f53aSPierre Jolivet PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 9771dd4f53aSPierre Jolivet PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 9781dd4f53aSPierre Jolivet PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &((PetscObject)B)->type_name)); 9791dd4f53aSPierre Jolivet PetscCall(MatSetUp(B)); 9801dd4f53aSPierre Jolivet 9811dd4f53aSPierre Jolivet B->ops->getinfo = MatGetInfo_External; 9821dd4f53aSPierre Jolivet B->factortype = ftype; 9831dd4f53aSPierre Jolivet B->trivialsymbolic = PETSC_TRUE; 9841dd4f53aSPierre Jolivet 9851dd4f53aSPierre Jolivet if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_Htool; 9861dd4f53aSPierre Jolivet else if (ftype == MAT_FACTOR_CHOLESKY) B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_Htool; 9871dd4f53aSPierre Jolivet 9881dd4f53aSPierre Jolivet PetscCall(PetscFree(B->solvertype)); 9891dd4f53aSPierre Jolivet PetscCall(PetscStrallocpy(MATSOLVERHTOOL, &B->solvertype)); 9901dd4f53aSPierre Jolivet 9911dd4f53aSPierre Jolivet PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_Htool)); 9921dd4f53aSPierre Jolivet *F = B; 9931dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 9941dd4f53aSPierre Jolivet } 9951dd4f53aSPierre Jolivet 9961dd4f53aSPierre Jolivet PETSC_INTERN PetscErrorCode MatSolverTypeRegister_Htool(void) 9971dd4f53aSPierre Jolivet { 9981dd4f53aSPierre Jolivet PetscFunctionBegin; 9991dd4f53aSPierre Jolivet PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_LU, MatGetFactor_htool_htool)); 10001dd4f53aSPierre Jolivet PetscCall(MatSolverTypeRegister(MATSOLVERHTOOL, MATHTOOL, MAT_FACTOR_CHOLESKY, MatGetFactor_htool_htool)); 10011dd4f53aSPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 10021dd4f53aSPierre Jolivet } 10031dd4f53aSPierre Jolivet 1004c7a4214aSPierre Jolivet /*@C 100511a5261eSBarry Smith MatCreateHtoolFromKernel - Creates a `MATHTOOL` from a user-supplied kernel. 1006c7a4214aSPierre Jolivet 1007cc4c1da9SBarry Smith Collective, No Fortran Support 1008cc4c1da9SBarry Smith 1009c7a4214aSPierre Jolivet Input Parameters: 1010c7a4214aSPierre Jolivet + comm - MPI communicator 10112ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 10122ef1f0ffSBarry Smith . n - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given) 10132ef1f0ffSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given) 10142ef1f0ffSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given) 1015c7a4214aSPierre Jolivet . spacedim - dimension of the space coordinates 1016c7a4214aSPierre Jolivet . coords_target - coordinates of the target 1017c7a4214aSPierre Jolivet . coords_source - coordinates of the source 10182ef1f0ffSBarry Smith . kernel - computational kernel (or `NULL`) 10192ef1f0ffSBarry Smith - kernelctx - kernel context (if kernel is `NULL`, the pointer must be of type htool::VirtualGenerator<PetscScalar>*) 1020c7a4214aSPierre Jolivet 1021c7a4214aSPierre Jolivet Output Parameter: 1022c7a4214aSPierre Jolivet . B - matrix 1023c7a4214aSPierre Jolivet 1024c7a4214aSPierre Jolivet Options Database Keys: 102511a5261eSBarry Smith + -mat_htool_min_cluster_size <`PetscInt`> - minimal leaf size in cluster tree 102611a5261eSBarry Smith . -mat_htool_epsilon <`PetscReal`> - relative error in Frobenius norm when approximating a block 102711a5261eSBarry Smith . -mat_htool_eta <`PetscReal`> - admissibility condition tolerance 102811a5261eSBarry Smith . -mat_htool_min_target_depth <`PetscInt`> - minimal cluster tree depth associated with the rows 102911a5261eSBarry Smith . -mat_htool_min_source_depth <`PetscInt`> - minimal cluster tree depth associated with the columns 10301dd4f53aSPierre Jolivet . -mat_htool_block_tree_consistency <`PetscBool`> - block tree consistency 103198e73e17SPierre Jolivet . -mat_htool_compressor <sympartialACA, fullACA, SVD> - type of compression 103298e73e17SPierre Jolivet - -mat_htool_clustering <PCARegular, PCAGeometric, BounbingBox1Regular, BoundingBox1Geometric> - type of clustering 1033c7a4214aSPierre Jolivet 1034c7a4214aSPierre Jolivet Level: intermediate 1035c7a4214aSPierre Jolivet 10361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATHTOOL`, `PCSetCoordinates()`, `MatHtoolSetKernel()`, `MatHtoolCompressorType`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()` 1037c7a4214aSPierre Jolivet @*/ 10388434afd1SBarry 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) 1039d71ae5a4SJacob Faibussowitsch { 1040c7a4214aSPierre Jolivet Mat A; 1041c7a4214aSPierre Jolivet Mat_Htool *a; 1042c7a4214aSPierre Jolivet 1043c7a4214aSPierre Jolivet PetscFunctionBegin; 10449566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 1045c7a4214aSPierre Jolivet PetscValidLogicalCollectiveInt(A, spacedim, 6); 10464f572ea9SToby Isaac PetscAssertPointer(coords_target, 7); 10474f572ea9SToby Isaac PetscAssertPointer(coords_source, 8); 1048c7a4214aSPierre Jolivet if (!kernelctx) PetscValidFunction(kernel, 9); 10494f572ea9SToby Isaac if (!kernel) PetscAssertPointer(kernelctx, 10); 10509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, n, M, N)); 10519566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATHTOOL)); 10529566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1053c9a4fd69SAudic XU PetscCall(MatShellGetContext(A, &a)); 1054c7a4214aSPierre Jolivet a->dim = spacedim; 1055c7a4214aSPierre Jolivet a->kernel = kernel; 1056c7a4214aSPierre Jolivet a->kernelctx = kernelctx; 10579566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->rmap->N * spacedim, &a->gcoords_target)); 10589566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_target + A->rmap->rstart * spacedim, coords_target, A->rmap->n * spacedim)); 1059462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_target, A->rmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global target coordinates */ 1060c7a4214aSPierre Jolivet if (coords_target != coords_source) { 10619566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(A->cmap->N * spacedim, &a->gcoords_source)); 10629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a->gcoords_source + A->cmap->rstart * spacedim, coords_source, A->cmap->n * spacedim)); 1063462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, a->gcoords_source, A->cmap->N * spacedim, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)A))); /* global source coordinates */ 1064c7a4214aSPierre Jolivet } else a->gcoords_source = a->gcoords_target; 10659566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(A->cmap->N, &a->work_source, A->rmap->N, &a->work_target)); 1066c7a4214aSPierre Jolivet *B = A; 10673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1068c7a4214aSPierre Jolivet } 1069c7a4214aSPierre Jolivet 1070c7a4214aSPierre Jolivet /*MC 1071c7a4214aSPierre Jolivet MATHTOOL = "htool" - A matrix type for hierarchical matrices using the Htool package. 1072c7a4214aSPierre Jolivet 10732ef1f0ffSBarry Smith Use `./configure --download-htool` to install PETSc to use Htool. 1074c7a4214aSPierre Jolivet 10752ef1f0ffSBarry Smith Options Database Key: 10762ef1f0ffSBarry Smith . -mat_type htool - matrix type to `MATHTOOL` 1077c7a4214aSPierre Jolivet 1078c7a4214aSPierre Jolivet Level: beginner 1079c7a4214aSPierre Jolivet 10801cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATDENSE`, `MatCreateHtoolFromKernel()`, `MatHtoolSetKernel()` 1081c7a4214aSPierre Jolivet M*/ 1082d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Htool(Mat A) 1083d71ae5a4SJacob Faibussowitsch { 1084c7a4214aSPierre Jolivet Mat_Htool *a; 1085c7a4214aSPierre Jolivet 1086c7a4214aSPierre Jolivet PetscFunctionBegin; 1087c9a4fd69SAudic XU PetscCall(MatSetType(A, MATSHELL)); 10884dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&a)); 1089c9a4fd69SAudic XU PetscCall(MatShellSetContext(A, a)); 1090c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Htool)); 1091c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Htool)); 1092c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))MatMult_Htool)); 1093c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool)); 1094c9a4fd69SAudic XU if (!PetscDefined(USE_COMPLEX)) PetscCall(MatShellSetOperation(A, MATOP_MULT_HERMITIAN_TRANSPOSE, (void (*)(void))MatMultTranspose_Htool)); 1095c7a4214aSPierre Jolivet A->ops->increaseoverlap = MatIncreaseOverlap_Htool; 1096c7a4214aSPierre Jolivet A->ops->createsubmatrices = MatCreateSubMatrices_Htool; 1097c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_Htool)); 1098c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_SET_FROM_OPTIONS, (void (*)(void))MatSetFromOptions_Htool)); 1099c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_GET_ROW, (void (*)(void))MatGetRow_Htool)); 1100c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_RESTORE_ROW, (void (*)(void))MatRestoreRow_Htool)); 1101c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_ASSEMBLY_END, (void (*)(void))MatAssemblyEnd_Htool)); 1102c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_TRANSPOSE, (void (*)(void))MatTranspose_Htool)); 1103c9a4fd69SAudic XU PetscCall(MatShellSetOperation(A, MATOP_DESTROY, (void (*)(void))MatDestroy_Htool)); 1104c7a4214aSPierre Jolivet a->dim = 0; 1105f22e26b7SPierre Jolivet a->gcoords_target = nullptr; 1106f22e26b7SPierre Jolivet a->gcoords_source = nullptr; 11071dd4f53aSPierre Jolivet a->min_cluster_size = 10; 1108c7a4214aSPierre Jolivet a->epsilon = PetscSqrtReal(PETSC_SMALL); 1109c7a4214aSPierre Jolivet a->eta = 10.0; 1110c7a4214aSPierre Jolivet a->depth[0] = 0; 1111c7a4214aSPierre Jolivet a->depth[1] = 0; 11121dd4f53aSPierre Jolivet a->block_tree_consistency = PETSC_TRUE; 1113c7a4214aSPierre Jolivet a->compressor = MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA; 11149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_seqdense_C", MatProductSetFromOptions_Htool)); 11159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_htool_mpidense_C", MatProductSetFromOptions_Htool)); 11169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_seqdense_C", MatConvert_Htool_Dense)); 11179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_htool_mpidense_C", MatConvert_Htool_Dense)); 11189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetHierarchicalMat_C", MatHtoolGetHierarchicalMat_Htool)); 11199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolSetKernel_C", MatHtoolSetKernel_Htool)); 11209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationSource_C", MatHtoolGetPermutationSource_Htool)); 11219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolGetPermutationTarget_C", MatHtoolGetPermutationTarget_Htool)); 11229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatHtoolUsePermutation_C", MatHtoolUsePermutation_Htool)); 1123c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContext_C", MatShellSetContext_Immutable)); 1124c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetContextDestroy_C", MatShellSetContextDestroy_Immutable)); 1125c9a4fd69SAudic XU PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatShellSetManageScalingShifts_C", MatShellSetManageScalingShifts_Immutable)); 1126c9a4fd69SAudic XU PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATHTOOL)); 11273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1128c7a4214aSPierre Jolivet } 1129