18d9f7141SDave May #include <petsc/private/petscimpl.h> 2120bdd93SDave May #include <petsc/private/matimpl.h> 36fc41876SBarry Smith #include <petsc/private/pcimpl.h> 41e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/ 51e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/ 6575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h" 71e07b27eSBarry Smith 8bf00f589SPatrick Sanan static PetscBool cited = PETSC_FALSE; 99371c9d4SSatish Balay static const char citation[] = "@inproceedings{MaySananRuppKnepleySmith2016,\n" 10bf00f589SPatrick Sanan " title = {Extreme-Scale Multigrid Components within PETSc},\n" 11bf00f589SPatrick Sanan " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 12bf00f589SPatrick Sanan " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 13bf00f589SPatrick Sanan " series = {PASC '16},\n" 14bf00f589SPatrick Sanan " isbn = {978-1-4503-4126-4},\n" 15bf00f589SPatrick Sanan " location = {Lausanne, Switzerland},\n" 16bf00f589SPatrick Sanan " pages = {5:1--5:12},\n" 17bf00f589SPatrick Sanan " articleno = {5},\n" 18bf00f589SPatrick Sanan " numpages = {12},\n" 19a8d69d7bSBarry Smith " url = {https://doi.acm.org/10.1145/2929908.2929913},\n" 20bf00f589SPatrick Sanan " doi = {10.1145/2929908.2929913},\n" 21bf00f589SPatrick Sanan " acmid = {2929913},\n" 22bf00f589SPatrick Sanan " publisher = {ACM},\n" 23bf00f589SPatrick Sanan " address = {New York, NY, USA},\n" 24bf00f589SPatrick Sanan " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 25bf00f589SPatrick Sanan " year = {2016}\n" 26bf00f589SPatrick Sanan "}\n"; 27bf00f589SPatrick Sanan 281e07b27eSBarry Smith /* 29c5083d92SDave May default setup mode 301e07b27eSBarry Smith 31c5083d92SDave May [1a] scatter to (FORWARD) 321e07b27eSBarry Smith x(comm) -> xtmp(comm) 33c5083d92SDave May [1b] local copy (to) ranks with color = 0 341e07b27eSBarry Smith xred(subcomm) <- xtmp 351e07b27eSBarry Smith 36c5083d92SDave May [2] solve on sub KSP to obtain yred(subcomm) 37c5083d92SDave May 38c5083d92SDave May [3a] local copy (from) ranks with color = 0 391e07b27eSBarry Smith yred(subcomm) --> xtmp 40c5083d92SDave May [2b] scatter from (REVERSE) 411e07b27eSBarry Smith xtmp(comm) -> y(comm) 421e07b27eSBarry Smith */ 431e07b27eSBarry Smith 448d9f7141SDave May /* 45d083f849SBarry Smith Collective[comm_f] 468d9f7141SDave May Notes 478d9f7141SDave May * Using comm_f = MPI_COMM_NULL will result in an error 488d9f7141SDave May * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid. 498d9f7141SDave May * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid. 508d9f7141SDave May */ 51d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f, MPI_Comm comm_c, PetscBool *isvalid) 52d71ae5a4SJacob Faibussowitsch { 5357f12427SDave May PetscInt valid = 1; 548d9f7141SDave May MPI_Group group_f, group_c; 558d9f7141SDave May PetscMPIInt count, k, size_f = 0, size_c = 0, size_c_sum = 0; 565bd1e576SStefano Zampini PetscMPIInt *ranks_f, *ranks_c; 578d9f7141SDave May 5857f12427SDave May PetscFunctionBegin; 5908401ef6SPierre Jolivet PetscCheck(comm_f != MPI_COMM_NULL, PETSC_COMM_SELF, PETSC_ERR_SUP, "comm_f cannot be MPI_COMM_NULL"); 608d9f7141SDave May 619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(comm_f, &group_f)); 6248a46eb9SPierre Jolivet if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_group(comm_c, &group_c)); 638d9f7141SDave May 649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_f, &size_f)); 6548a46eb9SPierre Jolivet if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_size(comm_c, &size_c)); 668d9f7141SDave May 678d9f7141SDave May /* check not all comm_c's are NULL */ 688d9f7141SDave May size_c_sum = size_c; 699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &size_c_sum, 1, MPI_INT, MPI_SUM, comm_f)); 705bd1e576SStefano Zampini if (size_c_sum == 0) valid = 0; 718d9f7141SDave May 728d9f7141SDave May /* check we can map at least 1 rank in comm_c to comm_f */ 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size_f, &ranks_f)); 749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size_c, &ranks_c)); 755bd1e576SStefano Zampini for (k = 0; k < size_f; k++) ranks_f[k] = MPI_UNDEFINED; 765bd1e576SStefano Zampini for (k = 0; k < size_c; k++) ranks_c[k] = k; 778d9f7141SDave May 7857f12427SDave May /* 7957f12427SDave May MPI_Group_translate_ranks() returns a non-zero exit code if any rank cannot be translated. 8057f12427SDave May I do not want the code to terminate immediately if this occurs, rather I want to throw 8157f12427SDave May the error later (during PCSetUp_Telescope()) via SETERRQ() with a message indicating 8257f12427SDave May that comm_c is not a valid sub-communicator. 839566063dSJacob Faibussowitsch Hence I purposefully do not call PetscCall() after MPI_Group_translate_ranks(). 8457f12427SDave May */ 858d9f7141SDave May count = 0; 868d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 8766b79024SDave May (void)MPI_Group_translate_ranks(group_c, size_c, ranks_c, group_f, ranks_f); 888d9f7141SDave May for (k = 0; k < size_f; k++) { 89ad540459SPierre Jolivet if (ranks_f[k] == MPI_UNDEFINED) count++; 908d9f7141SDave May } 918d9f7141SDave May } 925bd1e576SStefano Zampini if (count == size_f) valid = 0; 938d9f7141SDave May 94712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &valid, 1, MPIU_INT, MPI_MIN, comm_f)); 955bd1e576SStefano Zampini if (valid == 1) *isvalid = PETSC_TRUE; 965bd1e576SStefano Zampini else *isvalid = PETSC_FALSE; 978d9f7141SDave May 989566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks_f)); 999566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks_c)); 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group_f)); 10148a46eb9SPierre Jolivet if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Group_free(&group_c)); 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1038d9f7141SDave May } 1048d9f7141SDave May 105d71ae5a4SJacob Faibussowitsch DM private_PCTelescopeGetSubDM(PC_Telescope sred) 106d71ae5a4SJacob Faibussowitsch { 107c6a0d831SBarry Smith DM subdm = NULL; 1081e07b27eSBarry Smith 1099371c9d4SSatish Balay if (!PCTelescope_isActiveRank(sred)) { 1109371c9d4SSatish Balay subdm = NULL; 1119371c9d4SSatish Balay } else { 1121e07b27eSBarry Smith switch (sred->sr_type) { 113d71ae5a4SJacob Faibussowitsch case TELESCOPE_DEFAULT: 114d71ae5a4SJacob Faibussowitsch subdm = NULL; 115d71ae5a4SJacob Faibussowitsch break; 116d71ae5a4SJacob Faibussowitsch case TELESCOPE_DMDA: 117d71ae5a4SJacob Faibussowitsch subdm = ((PC_Telescope_DMDACtx *)sred->dm_ctx)->dmrepart; 118d71ae5a4SJacob Faibussowitsch break; 119d71ae5a4SJacob Faibussowitsch case TELESCOPE_DMPLEX: 120d71ae5a4SJacob Faibussowitsch subdm = NULL; 121d71ae5a4SJacob Faibussowitsch break; 1229371c9d4SSatish Balay case TELESCOPE_COARSEDM: 1233ba16761SJacob Faibussowitsch if (sred->ksp) PetscCallAbort(PETSC_COMM_SELF, KSPGetDM(sred->ksp, &subdm)); 1248d9f7141SDave May break; 1251e07b27eSBarry Smith } 1261e07b27eSBarry Smith } 1273ba16761SJacob Faibussowitsch return subdm; 1281e07b27eSBarry Smith } 1291e07b27eSBarry Smith 130d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetUp_default(PC pc, PC_Telescope sred) 131d71ae5a4SJacob Faibussowitsch { 1321e07b27eSBarry Smith PetscInt m, M, bs, st, ed; 1331e07b27eSBarry Smith Vec x, xred, yred, xtmp; 1341e07b27eSBarry Smith Mat B; 1351e07b27eSBarry Smith MPI_Comm comm, subcomm; 1361e07b27eSBarry Smith VecScatter scatter; 1371e07b27eSBarry Smith IS isin; 13854cd43dcSJunchao Zhang VecType vectype; 1391e07b27eSBarry Smith 1401e07b27eSBarry Smith PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: setup (default)\n")); 1421e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 1431e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1441e07b27eSBarry Smith 1459566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &B)); 1469566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &M, NULL)); 1479566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(B, &bs)); 1489566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B, &x, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(MatGetVecType(B, &vectype)); 1501e07b27eSBarry Smith 1511e07b27eSBarry Smith xred = NULL; 1523ac26c5eSBarry Smith m = 0; 15357f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1549566063dSJacob Faibussowitsch PetscCall(VecCreate(subcomm, &xred)); 1559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(xred, PETSC_DECIDE, M)); 1569566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(xred, bs)); 1579566063dSJacob Faibussowitsch PetscCall(VecSetType(xred, vectype)); /* Use the preconditioner matrix's vectype by default */ 1589566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(xred)); 1599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xred, &m)); 1601e07b27eSBarry Smith } 1611e07b27eSBarry Smith 1621e07b27eSBarry Smith yred = NULL; 16348a46eb9SPierre Jolivet if (PCTelescope_isActiveRank(sred)) PetscCall(VecDuplicate(xred, &yred)); 1641e07b27eSBarry Smith 1659566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &xtmp)); 1669566063dSJacob Faibussowitsch PetscCall(VecSetSizes(xtmp, m, PETSC_DECIDE)); 1679566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(xtmp, bs)); 1689566063dSJacob Faibussowitsch PetscCall(VecSetType(xtmp, vectype)); 1691e07b27eSBarry Smith 17057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1719566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xred, &st, &ed)); 1729566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (ed - st), st, 1, &isin)); 1731e07b27eSBarry Smith } else { 1749566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &st, &ed)); 1759566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, 0, st, 1, &isin)); 1761e07b27eSBarry Smith } 1779566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(isin, bs)); 1781e07b27eSBarry Smith 1799566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x, isin, xtmp, NULL, &scatter)); 1801e07b27eSBarry Smith 1811e07b27eSBarry Smith sred->isin = isin; 1821e07b27eSBarry Smith sred->scatter = scatter; 1831e07b27eSBarry Smith sred->xred = xred; 1841e07b27eSBarry Smith sred->yred = yred; 1851e07b27eSBarry Smith sred->xtmp = xtmp; 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1881e07b27eSBarry Smith } 1891e07b27eSBarry Smith 190d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeMatCreate_default(PC pc, PC_Telescope sred, MatReuse reuse, Mat *A) 191d71ae5a4SJacob Faibussowitsch { 1921e07b27eSBarry Smith MPI_Comm comm, subcomm; 1931e07b27eSBarry Smith Mat Bred, B; 1945624f943SPierre Jolivet PetscInt nr, nc, bs; 1951e07b27eSBarry Smith IS isrow, iscol; 1961e07b27eSBarry Smith Mat Blocal, *_Blocal; 1971e07b27eSBarry Smith 1981e07b27eSBarry Smith PetscFunctionBegin; 1999566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: updating the redundant preconditioned operator (default)\n")); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 2011e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2029566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &B)); 2039566063dSJacob Faibussowitsch PetscCall(MatGetSize(B, &nr, &nc)); 2041e07b27eSBarry Smith isrow = sred->isin; 2059566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, nc, 0, 1, &iscol)); 2069566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(iscol)); 2079566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(B, NULL, &bs)); 2089566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(iscol, bs)); 2099566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_SUBMAT_SINGLEIS, PETSC_TRUE)); 2109566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &_Blocal)); 2111e07b27eSBarry Smith Blocal = *_Blocal; 2129566063dSJacob Faibussowitsch PetscCall(PetscFree(_Blocal)); 2131e07b27eSBarry Smith Bred = NULL; 21457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2151e07b27eSBarry Smith PetscInt mm; 2161e07b27eSBarry Smith 217ad540459SPierre Jolivet if (reuse != MAT_INITIAL_MATRIX) Bred = *A; 2181e07b27eSBarry Smith 2199566063dSJacob Faibussowitsch PetscCall(MatGetSize(Blocal, &mm, NULL)); 2209566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, Blocal, mm, reuse, &Bred)); 2211e07b27eSBarry Smith } 2221e07b27eSBarry Smith *A = Bred; 2239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Blocal)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2261e07b27eSBarry Smith } 2271e07b27eSBarry Smith 228d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc, PC_Telescope sred, MatNullSpace nullspace, MatNullSpace *sub_nullspace) 229d71ae5a4SJacob Faibussowitsch { 2301e07b27eSBarry Smith PetscBool has_const; 2311e07b27eSBarry Smith const Vec *vecs; 232c41e779fSDave May Vec *sub_vecs = NULL; 233392968a1SPatrick Sanan PetscInt i, k, n = 0; 2341e07b27eSBarry Smith MPI_Comm subcomm; 2351e07b27eSBarry Smith 2361e07b27eSBarry Smith PetscFunctionBegin; 2371e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2389566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(nullspace, &has_const, &n, &vecs)); 2391e07b27eSBarry Smith 24057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 24148a46eb9SPierre Jolivet if (n) PetscCall(VecDuplicateVecs(sred->xred, n, &sub_vecs)); 2421e07b27eSBarry Smith } 2431e07b27eSBarry Smith 2441e07b27eSBarry Smith /* copy entries */ 2451e07b27eSBarry Smith for (k = 0; k < n; k++) { 2461e07b27eSBarry Smith const PetscScalar *x_array; 2471e07b27eSBarry Smith PetscScalar *LA_sub_vec; 24813c30530SDave May PetscInt st, ed; 2491e07b27eSBarry Smith 2501e07b27eSBarry Smith /* pull in vector x->xtmp */ 2519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD)); 2529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD)); 25347856c66SBarry Smith if (sub_vecs) { 254a04a6428SPatrick Sanan /* copy vector entries into xred */ 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(sred->xtmp, &x_array)); 256ea2b237eSDave May if (sub_vecs[k]) { 2579566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(sub_vecs[k], &st, &ed)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArray(sub_vecs[k], &LA_sub_vec)); 259ad540459SPierre Jolivet for (i = 0; i < ed - st; i++) LA_sub_vec[i] = x_array[i]; 2609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sub_vecs[k], &LA_sub_vec)); 2611e07b27eSBarry Smith } 2629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(sred->xtmp, &x_array)); 2631e07b27eSBarry Smith } 26447856c66SBarry Smith } 2651e07b27eSBarry Smith 26657f12427SDave May if (PCTelescope_isActiveRank(sred)) { 267d8b9d5b7SPatrick Sanan /* create new (near) nullspace for redundant object */ 2689566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(subcomm, has_const, n, sub_vecs, sub_nullspace)); 2699566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(n, &sub_vecs)); 27028b400f6SJacob Faibussowitsch PetscCheck(!nullspace->remove, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 27128b400f6SJacob Faibussowitsch PetscCheck(!nullspace->rmctx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope"); 272d8b9d5b7SPatrick Sanan } 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274392968a1SPatrick Sanan } 275392968a1SPatrick Sanan 276d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc, PC_Telescope sred, Mat sub_mat) 277d71ae5a4SJacob Faibussowitsch { 278392968a1SPatrick Sanan Mat B; 279392968a1SPatrick Sanan 280392968a1SPatrick Sanan PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &B)); 282392968a1SPatrick Sanan /* Propagate the nullspace if it exists */ 283392968a1SPatrick Sanan { 284392968a1SPatrick Sanan MatNullSpace nullspace, sub_nullspace; 2859566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(B, &nullspace)); 286392968a1SPatrick Sanan if (nullspace) { 2879566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: generating nullspace (default)\n")); 2889566063dSJacob Faibussowitsch PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nullspace, &sub_nullspace)); 28957f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2909566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(sub_mat, sub_nullspace)); 2919566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sub_nullspace)); 2921e07b27eSBarry Smith } 293392968a1SPatrick Sanan } 294392968a1SPatrick Sanan } 295392968a1SPatrick Sanan /* Propagate the near nullspace if it exists */ 296392968a1SPatrick Sanan { 297392968a1SPatrick Sanan MatNullSpace nearnullspace, sub_nearnullspace; 2989566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(B, &nearnullspace)); 299392968a1SPatrick Sanan if (nearnullspace) { 3009566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: generating near nullspace (default)\n")); 3019566063dSJacob Faibussowitsch PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nearnullspace, &sub_nearnullspace)); 30257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3039566063dSJacob Faibussowitsch PetscCall(MatSetNearNullSpace(sub_mat, sub_nearnullspace)); 3049566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sub_nearnullspace)); 305392968a1SPatrick Sanan } 306392968a1SPatrick Sanan } 307392968a1SPatrick Sanan } 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3091e07b27eSBarry Smith } 3101e07b27eSBarry Smith 311d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Telescope(PC pc, PetscViewer viewer) 312d71ae5a4SJacob Faibussowitsch { 3131e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 3141e07b27eSBarry Smith PetscBool iascii, isstring; 3151e07b27eSBarry Smith PetscViewer subviewer; 3161e07b27eSBarry Smith 3171e07b27eSBarry Smith PetscFunctionBegin; 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 3199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 3201e07b27eSBarry Smith if (iascii) { 3218d9f7141SDave May { 3221e07b27eSBarry Smith MPI_Comm comm, subcomm; 3231e07b27eSBarry Smith PetscMPIInt comm_size, subcomm_size; 3248d9f7141SDave May DM dm = NULL, subdm = NULL; 3251e07b27eSBarry Smith 3269566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 3271e07b27eSBarry Smith subdm = private_PCTelescopeGetSubDM(sred); 3288d9f7141SDave May 3298d9f7141SDave May if (sred->psubcomm) { 3301e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 3311e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 3329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 3339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm, &subcomm_size)); 3341e07b27eSBarry Smith 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 33663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: parent comm size reduction factor = %" PetscInt_FMT "\n", sred->redfactor)); 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: parent_size = %d , subcomm_size = %d\n", (int)comm_size, (int)subcomm_size)); 33848a10b22SPatrick Sanan switch (sred->subcommtype) { 339d71ae5a4SJacob Faibussowitsch case PETSC_SUBCOMM_INTERLACED: 340d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: type = %s\n", PetscSubcommTypes[sred->subcommtype])); 341d71ae5a4SJacob Faibussowitsch break; 342d71ae5a4SJacob Faibussowitsch case PETSC_SUBCOMM_CONTIGUOUS: 343d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm type = %s\n", PetscSubcommTypes[sred->subcommtype])); 344d71ae5a4SJacob Faibussowitsch break; 345d71ae5a4SJacob Faibussowitsch default: 346d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "General subcomm type not supported by PCTelescope"); 34748a10b22SPatrick Sanan } 3489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3498d9f7141SDave May } else { 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 3518d9f7141SDave May subcomm = sred->subcomm; 352ad540459SPierre Jolivet if (!PCTelescope_isActiveRank(sred)) subcomm = PETSC_COMM_SELF; 3538d9f7141SDave May 3549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "subcomm: using user provided sub-communicator\n")); 3569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3578d9f7141SDave May } 3588d9f7141SDave May 3599566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, subcomm, &subviewer)); 36057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(subviewer)); 3621e07b27eSBarry Smith 36348a46eb9SPierre Jolivet if (dm && sred->ignore_dm) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring DM\n")); 36448a46eb9SPierre Jolivet if (sred->ignore_kspcomputeoperators) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring KSPComputeOperators\n")); 3651e07b27eSBarry Smith switch (sred->sr_type) { 366d71ae5a4SJacob Faibussowitsch case TELESCOPE_DEFAULT: 367d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: default\n")); 368d71ae5a4SJacob Faibussowitsch break; 3691e07b27eSBarry Smith case TELESCOPE_DMDA: 3709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMDA auto-repartitioning\n")); 3719566063dSJacob Faibussowitsch PetscCall(DMView_DA_Short(subdm, subviewer)); 3721e07b27eSBarry Smith break; 373d71ae5a4SJacob Faibussowitsch case TELESCOPE_DMPLEX: 374d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMPLEX auto-repartitioning\n")); 375d71ae5a4SJacob Faibussowitsch break; 376d71ae5a4SJacob Faibussowitsch case TELESCOPE_COARSEDM: 377d71ae5a4SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: coarse DM\n")); 378d71ae5a4SJacob Faibussowitsch break; 3798d9f7141SDave May } 3808d9f7141SDave May 3818d9f7141SDave May if (dm) { 3828d9f7141SDave May PetscObject obj = (PetscObject)dm; 3839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object:")); 3843ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE)); 3853ba16761SJacob Faibussowitsch if (obj->type_name) PetscCall(PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name)); 3863ba16761SJacob Faibussowitsch if (obj->name) PetscCall(PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name)); 3873ba16761SJacob Faibussowitsch if (obj->prefix) PetscCall(PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix)); 3889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "\n")); 3893ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE)); 3908d9f7141SDave May } else { 3919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object: NULL\n")); 3928d9f7141SDave May } 3938d9f7141SDave May if (subdm) { 3948d9f7141SDave May PetscObject obj = (PetscObject)subdm; 3959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object:")); 3963ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE)); 3973ba16761SJacob Faibussowitsch if (obj->type_name) PetscCall(PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name)); 3983ba16761SJacob Faibussowitsch if (obj->name) PetscCall(PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name)); 3993ba16761SJacob Faibussowitsch if (obj->prefix) PetscCall(PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix)); 4009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "\n")); 4013ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE)); 4028d9f7141SDave May } else { 4039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object: NULL\n")); 4041e07b27eSBarry Smith } 4051e07b27eSBarry Smith 4069566063dSJacob Faibussowitsch PetscCall(KSPView(sred->ksp, subviewer)); 4079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(subviewer)); 4081e07b27eSBarry Smith } 4099566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, subcomm, &subviewer)); 4101e07b27eSBarry Smith } 4111e07b27eSBarry Smith } 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4131e07b27eSBarry Smith } 4141e07b27eSBarry Smith 415d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Telescope(PC pc) 416d71ae5a4SJacob Faibussowitsch { 4171e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 418bd49479cSSatish Balay MPI_Comm comm, subcomm = 0; 4191e07b27eSBarry Smith PCTelescopeType sr_type; 4201e07b27eSBarry Smith 4211e07b27eSBarry Smith PetscFunctionBegin; 4229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 4231e07b27eSBarry Smith 4241e07b27eSBarry Smith /* Determine type of setup/update */ 4251e07b27eSBarry Smith if (!pc->setupcalled) { 4261e07b27eSBarry Smith PetscBool has_dm, same; 4271e07b27eSBarry Smith DM dm; 4281e07b27eSBarry Smith 4291e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4301e07b27eSBarry Smith has_dm = PETSC_FALSE; 4319566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 432ad540459SPierre Jolivet if (dm) has_dm = PETSC_TRUE; 4331e07b27eSBarry Smith if (has_dm) { 4341e07b27eSBarry Smith /* check for dmda */ 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &same)); 4361e07b27eSBarry Smith if (same) { 4379566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: found DMDA\n")); 4381e07b27eSBarry Smith sr_type = TELESCOPE_DMDA; 4391e07b27eSBarry Smith } 4401e07b27eSBarry Smith /* check for dmplex */ 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &same)); 4421e07b27eSBarry Smith if (same) { 4439566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: found DMPLEX\n")); 4441e07b27eSBarry Smith sr_type = TELESCOPE_DMPLEX; 4451e07b27eSBarry Smith } 4468d9f7141SDave May 4478d9f7141SDave May if (sred->use_coarse_dm) { 4489566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: using coarse DM\n")); 4498d9f7141SDave May sr_type = TELESCOPE_COARSEDM; 4501e07b27eSBarry Smith } 4511e07b27eSBarry Smith 4521e07b27eSBarry Smith if (sred->ignore_dm) { 4539566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: ignoring DM\n")); 4541e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4551e07b27eSBarry Smith } 4568d9f7141SDave May } 4571e07b27eSBarry Smith sred->sr_type = sr_type; 4581e07b27eSBarry Smith } else { 4591e07b27eSBarry Smith sr_type = sred->sr_type; 4601e07b27eSBarry Smith } 4611e07b27eSBarry Smith 462d8b9d5b7SPatrick Sanan /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */ 4631e07b27eSBarry Smith switch (sr_type) { 4641e07b27eSBarry Smith case TELESCOPE_DEFAULT: 4651e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 4661e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 4671e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 4681e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 4691e07b27eSBarry Smith break; 4701e07b27eSBarry Smith case TELESCOPE_DMDA: 4711e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope_dmda; 472f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_dmda; 4731e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 4741e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 4751e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 4761e07b27eSBarry Smith sred->pctelescope_reset_type = PCReset_Telescope_dmda; 4771e07b27eSBarry Smith break; 478d71ae5a4SJacob Faibussowitsch case TELESCOPE_DMPLEX: 479d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_SUP, "Support for DMPLEX is currently not available"); 4808d9f7141SDave May case TELESCOPE_COARSEDM: 4818d9f7141SDave May pc->ops->apply = PCApply_Telescope_CoarseDM; 4828d9f7141SDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_CoarseDM; 4838d9f7141SDave May sred->pctelescope_setup_type = PCTelescopeSetUp_CoarseDM; 4848d9f7141SDave May sred->pctelescope_matcreate_type = NULL; 4858d9f7141SDave May sred->pctelescope_matnullspacecreate_type = NULL; /* PCTelescopeMatNullSpaceCreate_CoarseDM; */ 4868d9f7141SDave May sred->pctelescope_reset_type = PCReset_Telescope_CoarseDM; 4871e07b27eSBarry Smith break; 488d71ae5a4SJacob Faibussowitsch default: 489d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_SUP, "Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM"); 4908d9f7141SDave May } 4918d9f7141SDave May 4928d9f7141SDave May /* subcomm definition */ 4938d9f7141SDave May if (!pc->setupcalled) { 4948d9f7141SDave May if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) { 4958d9f7141SDave May if (!sred->psubcomm) { 4969566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(comm, &sred->psubcomm)); 4979566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(sred->psubcomm, sred->redfactor)); 4989566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(sred->psubcomm, sred->subcommtype)); 4998d9f7141SDave May sred->subcomm = PetscSubcommChild(sred->psubcomm); 5008d9f7141SDave May } 5018d9f7141SDave May } else { /* query PC for DM, check communicators */ 5028d9f7141SDave May DM dm, dm_coarse_partition = NULL; 5038d9f7141SDave May MPI_Comm comm_fine, comm_coarse_partition = MPI_COMM_NULL; 5048d9f7141SDave May PetscMPIInt csize_fine = 0, csize_coarse_partition = 0, cs[2], csg[2], cnt = 0; 5058d9f7141SDave May PetscBool isvalidsubcomm; 5068d9f7141SDave May 5079566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc, &dm)); 5088d9f7141SDave May comm_fine = PetscObjectComm((PetscObject)dm); 5099566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm, &dm_coarse_partition)); 510ad540459SPierre Jolivet if (dm_coarse_partition) cnt = 1; 5119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &cnt, 1, MPI_INT, MPI_SUM, comm_fine)); 51208401ef6SPierre Jolivet PetscCheck(cnt != 0, comm_fine, PETSC_ERR_SUP, "Zero instances of a coarse DM were found"); 5138d9f7141SDave May 5149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_fine, &csize_fine)); 5158d9f7141SDave May if (dm_coarse_partition) { 5168d9f7141SDave May comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition); 5179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_coarse_partition, &csize_coarse_partition)); 5188d9f7141SDave May } 5198d9f7141SDave May 5208d9f7141SDave May cs[0] = csize_fine; 5218d9f7141SDave May cs[1] = csize_coarse_partition; 5229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(cs, csg, 2, MPI_INT, MPI_MAX, comm_fine)); 52308401ef6SPierre Jolivet PetscCheck(csg[0] != csg[1], comm_fine, PETSC_ERR_SUP, "Coarse DM uses the same size communicator as the parent DM attached to the PC"); 5248d9f7141SDave May 5259566063dSJacob Faibussowitsch PetscCall(PCTelescopeTestValidSubcomm(comm_fine, comm_coarse_partition, &isvalidsubcomm)); 52628b400f6SJacob Faibussowitsch PetscCheck(isvalidsubcomm, comm_fine, PETSC_ERR_SUP, "Coarse DM communicator is not a sub-communicator of parentDM->comm"); 5278d9f7141SDave May sred->subcomm = comm_coarse_partition; 5288d9f7141SDave May } 5298d9f7141SDave May } 5308d9f7141SDave May subcomm = sred->subcomm; 5318d9f7141SDave May 5328d9f7141SDave May /* internal KSP */ 5338d9f7141SDave May if (!pc->setupcalled) { 5348d9f7141SDave May const char *prefix; 5358d9f7141SDave May 53657f12427SDave May if (PCTelescope_isActiveRank(sred)) { 5379566063dSJacob Faibussowitsch PetscCall(KSPCreate(subcomm, &sred->ksp)); 538*3821be0aSBarry Smith PetscCall(KSPSetNestLevel(sred->ksp, pc->kspnestlevel)); 5399566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(sred->ksp, pc->erroriffailure)); 5409566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)sred->ksp, (PetscObject)pc, 1)); 5419566063dSJacob Faibussowitsch PetscCall(KSPSetType(sred->ksp, KSPPREONLY)); 5429566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 5439566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(sred->ksp, prefix)); 5449566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(sred->ksp, "telescope_")); 5458d9f7141SDave May } 5461e07b27eSBarry Smith } 5471e07b27eSBarry Smith 5481e07b27eSBarry Smith /* setup */ 54948a46eb9SPierre Jolivet if (!pc->setupcalled && sred->pctelescope_setup_type) PetscCall(sred->pctelescope_setup_type(pc, sred)); 5501e07b27eSBarry Smith /* update */ 5511e07b27eSBarry Smith if (!pc->setupcalled) { 55248a46eb9SPierre Jolivet if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_INITIAL_MATRIX, &sred->Bred)); 5531baa6e33SBarry Smith if (sred->pctelescope_matnullspacecreate_type) PetscCall(sred->pctelescope_matnullspacecreate_type(pc, sred, sred->Bred)); 5541e07b27eSBarry Smith } else { 55548a46eb9SPierre Jolivet if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_REUSE_MATRIX, &sred->Bred)); 5561e07b27eSBarry Smith } 5571e07b27eSBarry Smith 5581e07b27eSBarry Smith /* common - no construction */ 55957f12427SDave May if (PCTelescope_isActiveRank(sred)) { 5609566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(sred->ksp, sred->Bred, sred->Bred)); 56148a46eb9SPierre Jolivet if (pc->setfromoptionscalled && !pc->setupcalled) PetscCall(KSPSetFromOptions(sred->ksp)); 5621e07b27eSBarry Smith } 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5641e07b27eSBarry Smith } 5651e07b27eSBarry Smith 566d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Telescope(PC pc, Vec x, Vec y) 567d71ae5a4SJacob Faibussowitsch { 5681e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5691e07b27eSBarry Smith Vec xtmp, xred, yred; 57013c30530SDave May PetscInt i, st, ed; 5711e07b27eSBarry Smith VecScatter scatter; 5721e07b27eSBarry Smith PetscScalar *array; 5731e07b27eSBarry Smith const PetscScalar *x_array; 5741e07b27eSBarry Smith 5751e07b27eSBarry Smith PetscFunctionBegin; 5769566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 577bf00f589SPatrick Sanan 5781e07b27eSBarry Smith xtmp = sred->xtmp; 5791e07b27eSBarry Smith scatter = sred->scatter; 5801e07b27eSBarry Smith xred = sred->xred; 5811e07b27eSBarry Smith yred = sred->yred; 5821e07b27eSBarry Smith 5831e07b27eSBarry Smith /* pull in vector x->xtmp */ 5849566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD)); 5859566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD)); 5861e07b27eSBarry Smith 587bf00f589SPatrick Sanan /* copy vector entries into xred */ 5889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xtmp, &x_array)); 5891e07b27eSBarry Smith if (xred) { 5901e07b27eSBarry Smith PetscScalar *LA_xred; 5919566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xred, &st, &ed)); 5929566063dSJacob Faibussowitsch PetscCall(VecGetArray(xred, &LA_xred)); 593ad540459SPierre Jolivet for (i = 0; i < ed - st; i++) LA_xred[i] = x_array[i]; 5949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xred, &LA_xred)); 5951e07b27eSBarry Smith } 5969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xtmp, &x_array)); 5971e07b27eSBarry Smith /* solve */ 59857f12427SDave May if (PCTelescope_isActiveRank(sred)) { 5999566063dSJacob Faibussowitsch PetscCall(KSPSolve(sred->ksp, xred, yred)); 6009566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(sred->ksp, pc, yred)); 6011e07b27eSBarry Smith } 6021e07b27eSBarry Smith /* return vector */ 6039566063dSJacob Faibussowitsch PetscCall(VecGetArray(xtmp, &array)); 6041e07b27eSBarry Smith if (yred) { 6051e07b27eSBarry Smith const PetscScalar *LA_yred; 6069566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(yred, &st, &ed)); 6079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(yred, &LA_yred)); 608ad540459SPierre Jolivet for (i = 0; i < ed - st; i++) array[i] = LA_yred[i]; 6099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(yred, &LA_yred)); 6101e07b27eSBarry Smith } 6119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xtmp, &array)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE)); 6139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE)); 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6151e07b27eSBarry Smith } 6161e07b27eSBarry Smith 617d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyRichardson_Telescope(PC pc, Vec x, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool zeroguess, PetscInt *outits, PCRichardsonConvergedReason *reason) 618d71ae5a4SJacob Faibussowitsch { 619f650675bSDave May PC_Telescope sred = (PC_Telescope)pc->data; 620a1d91a28SDave May Vec xtmp, yred; 621f650675bSDave May PetscInt i, st, ed; 622f650675bSDave May VecScatter scatter; 623f650675bSDave May const PetscScalar *x_array; 624f650675bSDave May PetscBool default_init_guess_value; 625f650675bSDave May 626f650675bSDave May PetscFunctionBegin; 627f650675bSDave May xtmp = sred->xtmp; 628f650675bSDave May scatter = sred->scatter; 629f650675bSDave May yred = sred->yred; 630f650675bSDave May 63108401ef6SPierre Jolivet PetscCheck(its <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PCApplyRichardson_Telescope only supports max_it = 1"); 632f650675bSDave May *reason = (PCRichardsonConvergedReason)0; 633f650675bSDave May 634f650675bSDave May if (!zeroguess) { 6359566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PCTelescope: Scattering y for non-zero initial guess\n")); 636f650675bSDave May /* pull in vector y->xtmp */ 6379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD)); 6389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD)); 639f650675bSDave May 640bf00f589SPatrick Sanan /* copy vector entries into xred */ 6419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xtmp, &x_array)); 642f650675bSDave May if (yred) { 643f650675bSDave May PetscScalar *LA_yred; 6449566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(yred, &st, &ed)); 6459566063dSJacob Faibussowitsch PetscCall(VecGetArray(yred, &LA_yred)); 646ad540459SPierre Jolivet for (i = 0; i < ed - st; i++) LA_yred[i] = x_array[i]; 6479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(yred, &LA_yred)); 648f650675bSDave May } 6499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xtmp, &x_array)); 650f650675bSDave May } 651f650675bSDave May 65257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 6539566063dSJacob Faibussowitsch PetscCall(KSPGetInitialGuessNonzero(sred->ksp, &default_init_guess_value)); 6549566063dSJacob Faibussowitsch if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, PETSC_TRUE)); 655f650675bSDave May } 656f650675bSDave May 6579566063dSJacob Faibussowitsch PetscCall(PCApply_Telescope(pc, x, y)); 658f650675bSDave May 65948a46eb9SPierre Jolivet if (PCTelescope_isActiveRank(sred)) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, default_init_guess_value)); 660f650675bSDave May 661f650675bSDave May if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 662f650675bSDave May *outits = 1; 6633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 664f650675bSDave May } 665f650675bSDave May 666d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Telescope(PC pc) 667d71ae5a4SJacob Faibussowitsch { 6681e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 6691e07b27eSBarry Smith 670362febeeSStefano Zampini PetscFunctionBegin; 6719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sred->isin)); 6729566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sred->scatter)); 6739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->xred)); 6749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->yred)); 6759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->xtmp)); 6769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sred->Bred)); 6779566063dSJacob Faibussowitsch PetscCall(KSPReset(sred->ksp)); 6781baa6e33SBarry Smith if (sred->pctelescope_reset_type) PetscCall(sred->pctelescope_reset_type(pc)); 6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6801e07b27eSBarry Smith } 6811e07b27eSBarry Smith 682d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Telescope(PC pc) 683d71ae5a4SJacob Faibussowitsch { 6841e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 6851e07b27eSBarry Smith 6861e07b27eSBarry Smith PetscFunctionBegin; 6879566063dSJacob Faibussowitsch PetscCall(PCReset_Telescope(pc)); 6889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&sred->ksp)); 6899566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&sred->psubcomm)); 6909566063dSJacob Faibussowitsch PetscCall(PetscFree(sred->dm_ctx)); 6912e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", NULL)); 6922e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", NULL)); 6932e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", NULL)); 6942e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", NULL)); 6952e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", NULL)); 6962e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", NULL)); 6972e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", NULL)); 6982e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", NULL)); 6992e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", NULL)); 7002e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", NULL)); 7012e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", NULL)); 7022e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", NULL)); 7039566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 7043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7051e07b27eSBarry Smith } 7061e07b27eSBarry Smith 707d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Telescope(PC pc, PetscOptionItems *PetscOptionsObject) 708d71ae5a4SJacob Faibussowitsch { 7091e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7101e07b27eSBarry Smith MPI_Comm comm; 7111e07b27eSBarry Smith PetscMPIInt size; 71248a10b22SPatrick Sanan PetscBool flg; 71348a10b22SPatrick Sanan PetscSubcommType subcommtype; 7141e07b27eSBarry Smith 7151e07b27eSBarry Smith PetscFunctionBegin; 7169566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 7179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 718d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Telescope options"); 7199566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_telescope_subcomm_type", "Subcomm type (interlaced or contiguous)", "PCTelescopeSetSubcommType", PetscSubcommTypes, (PetscEnum)sred->subcommtype, (PetscEnum *)&subcommtype, &flg)); 7201baa6e33SBarry Smith if (flg) PetscCall(PCTelescopeSetSubcommType(pc, subcommtype)); 7219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_telescope_reduction_factor", "Factor to reduce comm size by", "PCTelescopeSetReductionFactor", sred->redfactor, &sred->redfactor, NULL)); 72208401ef6SPierre Jolivet PetscCheck(sred->redfactor <= size, comm, PETSC_ERR_ARG_WRONG, "-pc_telescope_reduction_factor <= comm size"); 7239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_telescope_ignore_dm", "Ignore any DM attached to the PC", "PCTelescopeSetIgnoreDM", sred->ignore_dm, &sred->ignore_dm, NULL)); 7249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators", "Ignore method used to compute A", "PCTelescopeSetIgnoreKSPComputeOperators", sred->ignore_kspcomputeoperators, &sred->ignore_kspcomputeoperators, NULL)); 7259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_telescope_use_coarse_dm", "Define sub-communicator from the coarse DM", "PCTelescopeSetUseCoarseDM", sred->use_coarse_dm, &sred->use_coarse_dm, NULL)); 726d0609cedSBarry Smith PetscOptionsHeadEnd(); 7273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7281e07b27eSBarry Smith } 7291e07b27eSBarry Smith 7301e07b27eSBarry Smith /* PC simplementation specific API's */ 7311e07b27eSBarry Smith 732d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc, KSP *ksp) 733d71ae5a4SJacob Faibussowitsch { 7341e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 735bd49479cSSatish Balay PetscFunctionBegin; 7361e07b27eSBarry Smith if (ksp) *ksp = red->ksp; 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7381e07b27eSBarry Smith } 7391e07b27eSBarry Smith 740d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc, PetscSubcommType *subcommtype) 741d71ae5a4SJacob Faibussowitsch { 74248a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 74348a10b22SPatrick Sanan PetscFunctionBegin; 74448a10b22SPatrick Sanan if (subcommtype) *subcommtype = red->subcommtype; 7453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74648a10b22SPatrick Sanan } 74748a10b22SPatrick Sanan 748d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc, PetscSubcommType subcommtype) 749d71ae5a4SJacob Faibussowitsch { 75048a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 75148a10b22SPatrick Sanan 75248a10b22SPatrick Sanan PetscFunctionBegin; 75328b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "You cannot change the subcommunicator type for PCTelescope after it has been set up."); 75448a10b22SPatrick Sanan red->subcommtype = subcommtype; 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75648a10b22SPatrick Sanan } 75748a10b22SPatrick Sanan 758d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc, PetscInt *fact) 759d71ae5a4SJacob Faibussowitsch { 7601e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 761bd49479cSSatish Balay PetscFunctionBegin; 7621e07b27eSBarry Smith if (fact) *fact = red->redfactor; 7633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7641e07b27eSBarry Smith } 7651e07b27eSBarry Smith 766d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc, PetscInt fact) 767d71ae5a4SJacob Faibussowitsch { 7681e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 7691e07b27eSBarry Smith PetscMPIInt size; 7701e07b27eSBarry Smith 771bd49479cSSatish Balay PetscFunctionBegin; 7729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 77363a3b9bcSJacob Faibussowitsch PetscCheck(fact > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be positive", fact); 77463a3b9bcSJacob Faibussowitsch PetscCheck(fact <= size, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be <= comm.size", fact); 7751e07b27eSBarry Smith red->redfactor = fact; 7763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7771e07b27eSBarry Smith } 7781e07b27eSBarry Smith 779d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc, PetscBool *v) 780d71ae5a4SJacob Faibussowitsch { 7811e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 782bd49479cSSatish Balay PetscFunctionBegin; 7831e07b27eSBarry Smith if (v) *v = red->ignore_dm; 7843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7851e07b27eSBarry Smith } 78648a10b22SPatrick Sanan 787d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc, PetscBool v) 788d71ae5a4SJacob Faibussowitsch { 7891e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 790bd49479cSSatish Balay PetscFunctionBegin; 7911e07b27eSBarry Smith red->ignore_dm = v; 7923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7931e07b27eSBarry Smith } 7941e07b27eSBarry Smith 795d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc, PetscBool *v) 796d71ae5a4SJacob Faibussowitsch { 7978d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 7988d9f7141SDave May PetscFunctionBegin; 7998d9f7141SDave May if (v) *v = red->use_coarse_dm; 8003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8018d9f7141SDave May } 8028d9f7141SDave May 803d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc, PetscBool v) 804d71ae5a4SJacob Faibussowitsch { 8058d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 8068d9f7141SDave May PetscFunctionBegin; 8078d9f7141SDave May red->use_coarse_dm = v; 8083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8098d9f7141SDave May } 8108d9f7141SDave May 811d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool *v) 812d71ae5a4SJacob Faibussowitsch { 8130ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8140ae7c45bSDave May PetscFunctionBegin; 8150ae7c45bSDave May if (v) *v = red->ignore_kspcomputeoperators; 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8170ae7c45bSDave May } 81848a10b22SPatrick Sanan 819d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool v) 820d71ae5a4SJacob Faibussowitsch { 8210ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8220ae7c45bSDave May PetscFunctionBegin; 8230ae7c45bSDave May red->ignore_kspcomputeoperators = v; 8243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8250ae7c45bSDave May } 8260ae7c45bSDave May 827d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc, DM *dm) 828d71ae5a4SJacob Faibussowitsch { 8291e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 830bd49479cSSatish Balay PetscFunctionBegin; 8311e07b27eSBarry Smith *dm = private_PCTelescopeGetSubDM(red); 8323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8331e07b27eSBarry Smith } 8341e07b27eSBarry Smith 8351e07b27eSBarry Smith /*@ 836f1580f4eSBarry Smith PCTelescopeGetKSP - Gets the `KSP` created by the telescoping `PC`. 8371e07b27eSBarry Smith 8381e07b27eSBarry Smith Not Collective 8391e07b27eSBarry Smith 8401e07b27eSBarry Smith Input Parameter: 8411e07b27eSBarry Smith . pc - the preconditioner context 8421e07b27eSBarry Smith 8431e07b27eSBarry Smith Output Parameter: 844f1580f4eSBarry Smith . subksp - the `KSP` defined the smaller set of processes 8451e07b27eSBarry Smith 8461e07b27eSBarry Smith Level: advanced 8471e07b27eSBarry Smith 848f1580f4eSBarry Smith .seealso: `PCTELESCOPE` 8491e07b27eSBarry Smith @*/ 850d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetKSP(PC pc, KSP *subksp) 851d71ae5a4SJacob Faibussowitsch { 852bd49479cSSatish Balay PetscFunctionBegin; 853cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetKSP_C", (PC, KSP *), (pc, subksp)); 8543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8551e07b27eSBarry Smith } 8561e07b27eSBarry Smith 8571e07b27eSBarry Smith /*@ 858f1580f4eSBarry Smith PCTelescopeGetReductionFactor - Gets the factor by which the original number of MPI ranks has been reduced by. 8591e07b27eSBarry Smith 8601e07b27eSBarry Smith Not Collective 8611e07b27eSBarry Smith 8621e07b27eSBarry Smith Input Parameter: 8631e07b27eSBarry Smith . pc - the preconditioner context 8641e07b27eSBarry Smith 8651e07b27eSBarry Smith Output Parameter: 8661e07b27eSBarry Smith . fact - the reduction factor 8671e07b27eSBarry Smith 8681e07b27eSBarry Smith Level: advanced 8691e07b27eSBarry Smith 870f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetReductionFactor()` 8711e07b27eSBarry Smith @*/ 872d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetReductionFactor(PC pc, PetscInt *fact) 873d71ae5a4SJacob Faibussowitsch { 874bd49479cSSatish Balay PetscFunctionBegin; 875cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetReductionFactor_C", (PC, PetscInt *), (pc, fact)); 8763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8771e07b27eSBarry Smith } 8781e07b27eSBarry Smith 8791e07b27eSBarry Smith /*@ 880f1580f4eSBarry Smith PCTelescopeSetReductionFactor - Sets the factor by which the original number of MPI ranks will been reduced by. 8811e07b27eSBarry Smith 8821e07b27eSBarry Smith Not Collective 8831e07b27eSBarry Smith 8841e07b27eSBarry Smith Input Parameter: 8851e07b27eSBarry Smith . pc - the preconditioner context 8861e07b27eSBarry Smith 8871e07b27eSBarry Smith Output Parameter: 8881e07b27eSBarry Smith . fact - the reduction factor 8891e07b27eSBarry Smith 8901e07b27eSBarry Smith Level: advanced 8911e07b27eSBarry Smith 892f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeGetReductionFactor()` 8931e07b27eSBarry Smith @*/ 894d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetReductionFactor(PC pc, PetscInt fact) 895d71ae5a4SJacob Faibussowitsch { 896bd49479cSSatish Balay PetscFunctionBegin; 897cac4c232SBarry Smith PetscTryMethod(pc, "PCTelescopeSetReductionFactor_C", (PC, PetscInt), (pc, fact)); 8983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8991e07b27eSBarry Smith } 9001e07b27eSBarry Smith 9011e07b27eSBarry Smith /*@ 902f1580f4eSBarry Smith PCTelescopeGetIgnoreDM - Get the flag indicating if any `DM` attached to the `PC` will be used. 9031e07b27eSBarry Smith 9041e07b27eSBarry Smith Not Collective 9051e07b27eSBarry Smith 9061e07b27eSBarry Smith Input Parameter: 9071e07b27eSBarry Smith . pc - the preconditioner context 9081e07b27eSBarry Smith 9091e07b27eSBarry Smith Output Parameter: 9101e07b27eSBarry Smith . v - the flag 9111e07b27eSBarry Smith 9121e07b27eSBarry Smith Level: advanced 9131e07b27eSBarry Smith 914f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()` 9151e07b27eSBarry Smith @*/ 916d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetIgnoreDM(PC pc, PetscBool *v) 917d71ae5a4SJacob Faibussowitsch { 918bd49479cSSatish Balay PetscFunctionBegin; 919cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetIgnoreDM_C", (PC, PetscBool *), (pc, v)); 9203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9211e07b27eSBarry Smith } 9221e07b27eSBarry Smith 9231e07b27eSBarry Smith /*@ 9241e07b27eSBarry Smith PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 9251e07b27eSBarry Smith 9261e07b27eSBarry Smith Not Collective 9271e07b27eSBarry Smith 9281e07b27eSBarry Smith Input Parameter: 9291e07b27eSBarry Smith . pc - the preconditioner context 9301e07b27eSBarry Smith 9311e07b27eSBarry Smith Output Parameter: 9321e07b27eSBarry Smith . v - Use PETSC_TRUE to ignore any DM 9331e07b27eSBarry Smith 9341e07b27eSBarry Smith Level: advanced 9351e07b27eSBarry Smith 936f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeGetIgnoreDM()` 9371e07b27eSBarry Smith @*/ 938d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetIgnoreDM(PC pc, PetscBool v) 939d71ae5a4SJacob Faibussowitsch { 940bd49479cSSatish Balay PetscFunctionBegin; 941cac4c232SBarry Smith PetscTryMethod(pc, "PCTelescopeSetIgnoreDM_C", (PC, PetscBool), (pc, v)); 9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9431e07b27eSBarry Smith } 9441e07b27eSBarry Smith 9451e07b27eSBarry Smith /*@ 946f1580f4eSBarry Smith PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse `DM` attached to `DM` associated with the `PC` will be used. 9478d9f7141SDave May 9488d9f7141SDave May Not Collective 9498d9f7141SDave May 9508d9f7141SDave May Input Parameter: 9518d9f7141SDave May . pc - the preconditioner context 9528d9f7141SDave May 9538d9f7141SDave May Output Parameter: 9548d9f7141SDave May . v - the flag 9558d9f7141SDave May 9568d9f7141SDave May Level: advanced 9578d9f7141SDave May 958f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()` 9598d9f7141SDave May @*/ 960d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc, PetscBool *v) 961d71ae5a4SJacob Faibussowitsch { 9628d9f7141SDave May PetscFunctionBegin; 963cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetUseCoarseDM_C", (PC, PetscBool *), (pc, v)); 9643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9658d9f7141SDave May } 9668d9f7141SDave May 9678d9f7141SDave May /*@ 968f1580f4eSBarry Smith PCTelescopeSetUseCoarseDM - Set a flag to query the `DM` attached to the `PC` if it also has a coarse `DM` 9698d9f7141SDave May 9708d9f7141SDave May Not Collective 9718d9f7141SDave May 9728d9f7141SDave May Input Parameter: 9738d9f7141SDave May . pc - the preconditioner context 9748d9f7141SDave May 9758d9f7141SDave May Output Parameter: 976f1580f4eSBarry Smith . v - Use `PETSC_FALSE` to ignore any coarse `DM` 9778d9f7141SDave May 9788d9f7141SDave May Notes: 979f1580f4eSBarry Smith When you have specified to use a coarse `DM`, the communicator used to create the sub-KSP within `PCTELESCOPE` 980f1580f4eSBarry Smith will be that of the coarse `DM`. Hence the flags -pc_telescope_reduction_factor and 9818d9f7141SDave May -pc_telescope_subcomm_type will no longer have any meaning. 982f1580f4eSBarry Smith It is required that the communicator associated with the parent (fine) and the coarse `DM` are of different sizes. 983f1580f4eSBarry Smith An error will occur of the size of the communicator associated with the coarse `DM` 984f1580f4eSBarry Smith is the same as that of the parent `DM`. 9858d9f7141SDave May Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent. 9868d9f7141SDave May This will be checked at the time the preconditioner is setup and an error will occur if 9878d9f7141SDave May the coarse DM does not define a sub-communicator of that used by the parent DM. 9888d9f7141SDave May 9898d9f7141SDave May The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of 990f1580f4eSBarry Smith the `DM` used (e.g. it supports `DMSHELL`, `DMPLEX`, etc). 9918d9f7141SDave May 992f1580f4eSBarry Smith Support is currently only provided for the case when you are using `KSPSetComputeOperators()` 9938d9f7141SDave May 994f1580f4eSBarry Smith The user is required to compose a function with the parent DM to facilitate the transfer of fields (`Vec`) between the different decompositions defined by the fine and coarse `DM`s. 9958d9f7141SDave May In the user code, this is achieved via 9968d9f7141SDave May .vb 9978d9f7141SDave May { 9988d9f7141SDave May DM dm_fine; 9998d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method); 10008d9f7141SDave May } 10018d9f7141SDave May .ve 10028d9f7141SDave May The signature of the user provided field scatter method is 10038d9f7141SDave May .vb 10048d9f7141SDave May PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse); 10058d9f7141SDave May .ve 1006f1580f4eSBarry Smith The user must provide support for both mode = `SCATTER_FORWARD` and mode = `SCATTER_REVERSE`. 1007f1580f4eSBarry Smith `SCATTER_FORWARD` implies the direction of transfer is from the parent (fine) `DM` to the coarse `DM`. 10088d9f7141SDave May 10098d9f7141SDave May Optionally, the user may also compose a function with the parent DM to facilitate the transfer 1010f1580f4eSBarry Smith of state variables between the fine and coarse `DM`s. 10118d9f7141SDave May In the context of a finite element discretization, an example state variable might be 10128d9f7141SDave May values associated with quadrature points within each element. 10138d9f7141SDave May A user provided state scatter method is composed via 10148d9f7141SDave May .vb 10158d9f7141SDave May { 10168d9f7141SDave May DM dm_fine; 10178d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method); 10188d9f7141SDave May } 10198d9f7141SDave May .ve 10208d9f7141SDave May The signature of the user provided state scatter method is 10218d9f7141SDave May .vb 10228d9f7141SDave May PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse); 10238d9f7141SDave May .ve 1024f1580f4eSBarry Smith `SCATTER_FORWARD` implies the direction of transfer is from the fine `DM` to the coarse `DM`. 1025f1580f4eSBarry Smith The user is only required to support mode = `SCATTER_FORWARD`. 10268d9f7141SDave May No assumption is made about the data type of the state variables. 1027f1580f4eSBarry Smith These must be managed by the user and must be accessible from the `DM`. 10288d9f7141SDave May 1029f1580f4eSBarry Smith Care must be taken in defining the user context passed to `KSPSetComputeOperators()` which is to be 1030f1580f4eSBarry Smith associated with the sub-`KSP` residing within `PCTELESCOPE`. 1031f1580f4eSBarry Smith In general, `PCTELESCOPE` assumes that the context on the fine and coarse `DM` used with 1032f1580f4eSBarry Smith `KSPSetComputeOperators()` should be "similar" in type or origin. 1033f1580f4eSBarry Smith Specifically the following rules are used to infer what context on the sub-`KSP` should be. 10348d9f7141SDave May 1035f1580f4eSBarry Smith First the contexts from the `KSP` and the fine and coarse `DM`s are retrieved. 1036f1580f4eSBarry Smith Note that the special case of a `DMSHELL` context is queried. 10378d9f7141SDave May 10388d9f7141SDave May .vb 10398d9f7141SDave May DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx); 10408d9f7141SDave May DMGetApplicationContext(dm_fine,&dmfine_appctx); 10418d9f7141SDave May DMShellGetContext(dm_fine,&dmfine_shellctx); 10428d9f7141SDave May 10438d9f7141SDave May DMGetApplicationContext(dm_coarse,&dmcoarse_appctx); 10448d9f7141SDave May DMShellGetContext(dm_coarse,&dmcoarse_shellctx); 10458d9f7141SDave May .ve 10468d9f7141SDave May 10478d9f7141SDave May The following rules are then enforced: 10488d9f7141SDave May 10498d9f7141SDave May 1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP: 1050f1580f4eSBarry Smith `KSPSetComputeOperators`(sub_ksp,dmfine_kspfunc,NULL); 10518d9f7141SDave May 10528d9f7141SDave May 2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx, 1053feefa0e1SJacob Faibussowitsch 10548d9f7141SDave May check that dmcoarse_appctx is also non-NULL. If this is true, then: 1055f1580f4eSBarry Smith `KSPSetComputeOperators`(sub_ksp,dmfine_kspfunc,dmcoarse_appctx); 10568d9f7141SDave May 10578d9f7141SDave May 3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx, 1058feefa0e1SJacob Faibussowitsch 10598d9f7141SDave May check that dmcoarse_shellctx is also non-NULL. If this is true, then: 1060f1580f4eSBarry Smith `KSPSetComputeOperators`(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx); 10618d9f7141SDave May 1062f1580f4eSBarry Smith If neither of the above three tests passed, then `PCTELESCOPE` cannot safely determine what 1063f1580f4eSBarry Smith context should be provided to `KSPSetComputeOperators()` for use with the sub-`KSP`. 10648d9f7141SDave May In this case, an additional mechanism is provided via a composed function which will return 10658d9f7141SDave May the actual context to be used. To use this feature you must compose the "getter" function 1066f1580f4eSBarry Smith with the coarse `DM`, e.g. 10678d9f7141SDave May .vb 10688d9f7141SDave May { 10698d9f7141SDave May DM dm_coarse; 10708d9f7141SDave May PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter); 10718d9f7141SDave May } 10728d9f7141SDave May .ve 10738d9f7141SDave May The signature of the user provided method is 10748d9f7141SDave May .vb 10758d9f7141SDave May PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext); 10768d9f7141SDave May .ve 10778d9f7141SDave May 10788d9f7141SDave May Level: advanced 10798d9f7141SDave May 108042747ad1SJacob Faibussowitsch .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()` 10818d9f7141SDave May @*/ 1082d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc, PetscBool v) 1083d71ae5a4SJacob Faibussowitsch { 10848d9f7141SDave May PetscFunctionBegin; 1085cac4c232SBarry Smith PetscTryMethod(pc, "PCTelescopeSetUseCoarseDM_C", (PC, PetscBool), (pc, v)); 10863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10878d9f7141SDave May } 10888d9f7141SDave May 10898d9f7141SDave May /*@ 1090f1580f4eSBarry Smith PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if `KSPComputeOperators()` will be used. 10910ae7c45bSDave May 10920ae7c45bSDave May Not Collective 10930ae7c45bSDave May 10940ae7c45bSDave May Input Parameter: 10950ae7c45bSDave May . pc - the preconditioner context 10960ae7c45bSDave May 10970ae7c45bSDave May Output Parameter: 10980ae7c45bSDave May . v - the flag 10990ae7c45bSDave May 11000ae7c45bSDave May Level: advanced 11010ae7c45bSDave May 1102f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeSetIgnoreKSPComputeOperators()` 11030ae7c45bSDave May @*/ 1104d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc, PetscBool *v) 1105d71ae5a4SJacob Faibussowitsch { 11060ae7c45bSDave May PetscFunctionBegin; 1107cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", (PC, PetscBool *), (pc, v)); 11083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11090ae7c45bSDave May } 11100ae7c45bSDave May 11110ae7c45bSDave May /*@ 1112f1580f4eSBarry Smith PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore `KSPComputeOperators()`. 11130ae7c45bSDave May 11140ae7c45bSDave May Not Collective 11150ae7c45bSDave May 11160ae7c45bSDave May Input Parameter: 11170ae7c45bSDave May . pc - the preconditioner context 11180ae7c45bSDave May 11190ae7c45bSDave May Output Parameter: 1120f1580f4eSBarry Smith . v - Use `PETSC_TRUE` to ignore the method (if defined) set via `KSPSetComputeOperators()` on pc 11210ae7c45bSDave May 11220ae7c45bSDave May Level: advanced 11230ae7c45bSDave May 1124f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeGetIgnoreKSPComputeOperators()` 11250ae7c45bSDave May @*/ 1126d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc, PetscBool v) 1127d71ae5a4SJacob Faibussowitsch { 11280ae7c45bSDave May PetscFunctionBegin; 1129cac4c232SBarry Smith PetscTryMethod(pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", (PC, PetscBool), (pc, v)); 11303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11310ae7c45bSDave May } 11320ae7c45bSDave May 11330ae7c45bSDave May /*@ 1134f1580f4eSBarry Smith PCTelescopeGetDM - Get the re-partitioned `DM` attached to the sub-`KSP`. 11351e07b27eSBarry Smith 11361e07b27eSBarry Smith Not Collective 11371e07b27eSBarry Smith 11381e07b27eSBarry Smith Input Parameter: 11391e07b27eSBarry Smith . pc - the preconditioner context 11401e07b27eSBarry Smith 11411e07b27eSBarry Smith Output Parameter: 11421e07b27eSBarry Smith . subdm - The re-partitioned DM 11431e07b27eSBarry Smith 11441e07b27eSBarry Smith Level: advanced 11451e07b27eSBarry Smith 1146f1580f4eSBarry Smith .seealso: `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeGetIgnoreKSPComputeOperators()` 11471e07b27eSBarry Smith @*/ 1148d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetDM(PC pc, DM *subdm) 1149d71ae5a4SJacob Faibussowitsch { 1150bd49479cSSatish Balay PetscFunctionBegin; 1151cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetDM_C", (PC, DM *), (pc, subdm)); 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11531e07b27eSBarry Smith } 11541e07b27eSBarry Smith 115548a10b22SPatrick Sanan /*@ 115648a10b22SPatrick Sanan PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous) 115748a10b22SPatrick Sanan 115848a10b22SPatrick Sanan Logically Collective 115948a10b22SPatrick Sanan 1160d8d19677SJose E. Roman Input Parameters: 11611dae98e4SBarry Smith + pc - the preconditioner context 1162f1580f4eSBarry Smith - subcommtype - the subcommunicator type (see `PetscSubcommType`) 116348a10b22SPatrick Sanan 116448a10b22SPatrick Sanan Level: advanced 116548a10b22SPatrick Sanan 1166db781477SPatrick Sanan .seealso: `PetscSubcommType`, `PetscSubcomm`, `PCTELESCOPE` 116748a10b22SPatrick Sanan @*/ 1168d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) 1169d71ae5a4SJacob Faibussowitsch { 117048a10b22SPatrick Sanan PetscFunctionBegin; 1171cac4c232SBarry Smith PetscTryMethod(pc, "PCTelescopeSetSubcommType_C", (PC, PetscSubcommType), (pc, subcommtype)); 11723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117348a10b22SPatrick Sanan } 117448a10b22SPatrick Sanan 117548a10b22SPatrick Sanan /*@ 117648a10b22SPatrick Sanan PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous) 117748a10b22SPatrick Sanan 117848a10b22SPatrick Sanan Not Collective 117948a10b22SPatrick Sanan 118048a10b22SPatrick Sanan Input Parameter: 118148a10b22SPatrick Sanan . pc - the preconditioner context 118248a10b22SPatrick Sanan 118348a10b22SPatrick Sanan Output Parameter: 1184f1580f4eSBarry Smith . subcommtype - the subcommunicator type (see `PetscSubcommType`) 118548a10b22SPatrick Sanan 118648a10b22SPatrick Sanan Level: advanced 118748a10b22SPatrick Sanan 1188db781477SPatrick Sanan .seealso: `PetscSubcomm`, `PetscSubcommType`, `PCTELESCOPE` 118948a10b22SPatrick Sanan @*/ 1190d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) 1191d71ae5a4SJacob Faibussowitsch { 119248a10b22SPatrick Sanan PetscFunctionBegin; 1193cac4c232SBarry Smith PetscUseMethod(pc, "PCTelescopeGetSubcommType_C", (PC, PetscSubcommType *), (pc, subcommtype)); 11943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119548a10b22SPatrick Sanan } 119648a10b22SPatrick Sanan 11971e07b27eSBarry Smith /*MC 1198f1580f4eSBarry Smith PCTELESCOPE - Runs a `KSP` solver on a sub-communicator. MPI ranks not in the sub-communicator are idle during the solve. 11991e07b27eSBarry Smith 1200f1580f4eSBarry Smith Options Database Keys: 120100fea0ebSDave May + -pc_telescope_reduction_factor <r> - factor to reduce the communicator size by. e.g. with 64 MPI ranks and r=4, the new sub-communicator will have 64/4 = 16 ranks. 120200fea0ebSDave May . -pc_telescope_ignore_dm - flag to indicate whether an attached DM should be ignored. 120300fea0ebSDave May . -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI ranks on the sub-communicator. see PetscSubcomm for more information. 1204f1580f4eSBarry Smith . -pc_telescope_ignore_kspcomputeoperators - flag to indicate whether `KSPSetComputeOperators()` should be used on the sub-KSP. 1205f1580f4eSBarry Smith - -pc_telescope_use_coarse_dm - flag to indicate whether the coarse `DM` should be used to define the sub-communicator. 12061e07b27eSBarry Smith 12071e07b27eSBarry Smith Level: advanced 12081e07b27eSBarry Smith 12091e07b27eSBarry Smith Notes: 1210f1580f4eSBarry Smith Assuming that the parent preconditioner `PC` is defined on a communicator c, this implementation 1211f1580f4eSBarry Smith creates a child sub-communicator (c') containing fewer MPI ranks than the original parent preconditioner `PC`. 1212f1580f4eSBarry Smith The preconditioner is deemed telescopic as it only calls `KSPSolve()` on a single 1213f1580f4eSBarry Smith sub-communicator, in contrast with `PCREDUNDANT` which calls `KSPSolve()` on N sub-communicators. 121400fea0ebSDave May This means there will be MPI ranks which will be idle during the application of this preconditioner. 1215f1580f4eSBarry Smith Additionally, in comparison with `PCREDUNDANT`, `PCTELESCOPE` can utilize an attached `DM`. 12166fc41876SBarry Smith 1217f1580f4eSBarry Smith The default type of the sub `KSP` (the `KSP` defined on c') is `KSPPREONLY`. 121800fea0ebSDave May 1219f1580f4eSBarry Smith There are three setup mechanisms for `PCTELESCOPE`. Features support by each type are described below. 1220f1580f4eSBarry Smith In the following, we will refer to the operators B and B', these are the Bmat provided to the `KSP` on the 122100fea0ebSDave May communicators c and c' respectively. 122200fea0ebSDave May 122300fea0ebSDave May [1] Default setup 1224f1580f4eSBarry Smith The sub-communicator c' is created via `PetscSubcommCreate()`. 1225a5b23f4aSJose E. Roman Explicitly defined nullspace and near nullspace vectors will be propagated from B to B'. 1226f1580f4eSBarry Smith Currently there is no support define nullspaces via a user supplied method (e.g. as passed to `MatNullSpaceSetFunction()`). 1227f1580f4eSBarry Smith No support is provided for `KSPSetComputeOperators()`. 122800fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 122900fea0ebSDave May 1230f1580f4eSBarry Smith [2] `DM` aware setup 1231f1580f4eSBarry Smith If a `DM` is attached to the `PC`, it is re-partitioned on the sub-communicator c'. 1232f1580f4eSBarry Smith c' is created via `PetscSubcommCreate()`. 1233f1580f4eSBarry Smith Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned `DM`. 1234f1580f4eSBarry Smith Currently only support for re-partitioning a `DMDA` is provided. 123500fea0ebSDave May Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, re-partitioned and set on the re-partitioned Bmat operator (B'). 1236f1580f4eSBarry Smith Currently there is no support define nullspaces via a user supplied method (e.g. as passed to `MatNullSpaceSetFunction()`). 1237f1580f4eSBarry Smith Support is provided for `KSPSetComputeOperators()`. The user provided function and context is propagated to the sub `KSP`. 123800fea0ebSDave May This is fragile since the user must ensure that their user context is valid for use on c'. 123900fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 12401e07b27eSBarry Smith 1241f1580f4eSBarry Smith [3] Coarse `DM` setup 1242f1580f4eSBarry Smith If a `DM` (dmfine) is attached to the `PC`, dmfine is queried for a "coarse" `DM` (call this dmcoarse) via `DMGetCoarseDM()`. 1243f1580f4eSBarry Smith `PCTELESCOPE` will interpret the coarse `DM` as being defined on a sub-communicator of c. 1244f1580f4eSBarry Smith The communicator associated with dmcoarse will define the c' to be used within `PCTELESCOPE`. 1245f1580f4eSBarry Smith `PCTELESCOPE` will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported. 1246f1580f4eSBarry Smith The intention of this setup type is that `PCTELESCOPE` will use an existing (e.g. user defined) communicator hierarchy, say as would be 124700fea0ebSDave May available with using multi-grid on unstructured meshes. 124800fea0ebSDave May This setup will not use the command line options -pc_telescope_reduction_factor or -pc_telescope_subcomm_type. 124900fea0ebSDave May Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, scattered into the correct ordering consistent with dmcoarse and set on B'. 1250f1580f4eSBarry Smith Currently there is no support define nullspaces via a user supplied method (e.g. as passed to `MatNullSpaceSetFunction()`). 1251f1580f4eSBarry Smith There is no general method to permute field orderings, hence only `KSPSetComputeOperators()` is supported. 1252f1580f4eSBarry Smith The user must use `PetscObjectComposeFunction()` with dmfine to define the method to scatter fields from dmfine to dmcoarse. 1253f1580f4eSBarry Smith Propagation of the user context for `KSPSetComputeOperators()` on the sub `KSP` is attempted by querying the `DM` contexts associated with dmfine and dmcoarse. Alternatively, the user may use `PetscObjectComposeFunction()` with dmcoarse to define a method which will return the appropriate user context for `KSPSetComputeOperators()`. 125400fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 1255f1580f4eSBarry Smith This setup can be invoked by the option -pc_telescope_use_coarse_dm or by calling `PCTelescopeSetUseCoarseDM`(pc,`PETSC_TRUE`); 1256f1580f4eSBarry Smith Further information about the user-provided methods required by this setup type are described here `PCTelescopeSetUseCoarseDM()`. 12576fc41876SBarry Smith 12586fc41876SBarry Smith Developer Notes: 1259f1580f4eSBarry Smith During `PCSetup()`, the B operator is scattered onto c'. 1260f1580f4eSBarry Smith Within `PCApply()`, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 1261f1580f4eSBarry Smith Then, `KSPSolve()` is executed on the c' communicator. 12626fc41876SBarry Smith 1263f1580f4eSBarry Smith The communicator used within the telescoping preconditioner is defined by a `PetscSubcomm` using the INTERLACED 1264f1580f4eSBarry Smith creation routine by default (this can be changed with -pc_telescope_subcomm_type). We run the sub `KSP` on only the ranks within the communicator which have a color equal to zero. 12656fc41876SBarry Smith 1266005d9f20SPatrick Sanan The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator. 12678439623fSPatrick Sanan In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into 1268005d9f20SPatrick Sanan a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c') 12696fc41876SBarry Smith 1270f1580f4eSBarry Smith The telescoping preconditioner can re-partition an attached DM if it is a `DMDA` (2D or 3D - 1271f1580f4eSBarry Smith support for 1D `DMDA`s is not provided). If a `DMDA` is found, a topologically equivalent `DMDA` is created on c' 1272f1580f4eSBarry Smith and this new `DM` is attached the sub `KSP`. The design of telescope is such that it should be possible to extend support 1273f1580f4eSBarry Smith for re-partitioning other to DM's (e.g. `DMPLEX`). The user can supply a flag to ignore attached DMs. 127400fea0ebSDave May Alternatively, user-provided re-partitioned DMs can be used via -pc_telescope_use_coarse_dm. 12756fc41876SBarry Smith 127600fea0ebSDave May With the default setup mode, B' is defined by fusing rows (in order) associated with MPI ranks common to c and c'. 12776fc41876SBarry Smith 1278f1580f4eSBarry Smith When a `DMDA` is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B 1279f1580f4eSBarry Smith into the ordering defined by the `DMDA` on c', (ii) extracting the local chunks via `MatCreateSubMatrices()`, (iii) fusing the 1280f1580f4eSBarry Smith locally (sequential) matrices defined on the ranks common to c and c' into B' using `MatCreateMPIMatConcatenateSeqMat()` 12816fc41876SBarry Smith 12828439623fSPatrick Sanan Limitations/improvements include the following. 1283f1580f4eSBarry Smith `VecPlaceArray()` could be used within `PCApply()` to improve efficiency and reduce memory usage. 1284f1580f4eSBarry Smith A unified mechanism to query for user contexts as required by `KSPSetComputeOperators()` and `MatNullSpaceSetFunction()`. 12856fc41876SBarry Smith 1286f1580f4eSBarry Smith The symmetric permutation used when a `DMDA` is encountered is performed via explicitly assembling a permutation matrix P, 1287f1580f4eSBarry Smith and performing P^T.A.P. Possibly it might be more efficient to use `MatPermute()`. We opted to use P^T.A.P as it appears 1288f1580f4eSBarry Smith `VecPermute()` does not support the use case required here. By computing P, one can permute both the operator and RHS in a 12896fc41876SBarry Smith consistent manner. 12906fc41876SBarry Smith 129100fea0ebSDave May Mapping of vectors (default setup mode) is performed in the following way. 129200fea0ebSDave May Suppose the parent communicator size was 4, and we set a reduction factor of 2; this would give a comm size on c' of 2. 12938439623fSPatrick Sanan Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2. 129400fea0ebSDave May We perform the scatter to the sub-communicator in the following way. 1295514db83aSDave May [1] Given a vector x defined on communicator c 12966fc41876SBarry Smith 1297514db83aSDave May .vb 1298514db83aSDave May rank(c) local values of x 1299514db83aSDave May ------- ---------------------------------------- 1300514db83aSDave May 0 [ 0.0, 1.0, 2.0, 3.0, 4.0, 5.0 ] 1301514db83aSDave May 1 [ 6.0, 7.0, 8.0, 9.0, 10.0, 11.0 ] 1302514db83aSDave May 2 [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ] 1303514db83aSDave May 3 [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ] 1304514db83aSDave May .ve 13056fc41876SBarry Smith 1306514db83aSDave May scatter into xtmp defined also on comm c, so that we have the following values 13076fc41876SBarry Smith 1308514db83aSDave May .vb 1309514db83aSDave May rank(c) local values of xtmp 1310514db83aSDave May ------- ---------------------------------------------------------------------------- 1311514db83aSDave May 0 [ 0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0, 10.0, 11.0 ] 1312514db83aSDave May 1 [ ] 1313514db83aSDave May 2 [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0, 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ] 1314514db83aSDave May 3 [ ] 1315514db83aSDave May .ve 13166fc41876SBarry Smith 13176fc41876SBarry Smith The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 13186fc41876SBarry Smith 1319514db83aSDave May [2] Copy the values from ranks 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'. 13206fc41876SBarry Smith Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 13216fc41876SBarry Smith 1322514db83aSDave May .vb 1323514db83aSDave May rank(c') local values of xred 1324514db83aSDave May -------- ---------------------------------------------------------------------------- 1325514db83aSDave May 0 [ 0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0, 10.0, 11.0 ] 1326514db83aSDave May 1 [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0, 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ] 1327514db83aSDave May .ve 13281e07b27eSBarry Smith 13291e07b27eSBarry Smith Contributed by Dave May 13301e07b27eSBarry Smith 1331bf00f589SPatrick Sanan Reference: 1332bf00f589SPatrick Sanan Dave A. May, Patrick Sanan, Karl Rupp, Matthew G. Knepley, and Barry F. Smith, "Extreme-Scale Multigrid Components within PETSc". 2016. In Proceedings of the Platform for Advanced Scientific Computing Conference (PASC '16). DOI: 10.1145/2929908.2929913 1333bf00f589SPatrick Sanan 1334db781477SPatrick Sanan .seealso: `PCTelescopeGetKSP()`, `PCTelescopeGetDM()`, `PCTelescopeGetReductionFactor()`, `PCTelescopeSetReductionFactor()`, `PCTelescopeGetIgnoreDM()`, `PCTelescopeSetIgnoreDM()`, `PCREDUNDANT` 13351e07b27eSBarry Smith M*/ 1336d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 1337d71ae5a4SJacob Faibussowitsch { 13381e07b27eSBarry Smith struct _PC_Telescope *sred; 13391e07b27eSBarry Smith 13401e07b27eSBarry Smith PetscFunctionBegin; 13414dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&sred)); 13422a22aa42SDave May sred->psubcomm = NULL; 134348a10b22SPatrick Sanan sred->subcommtype = PETSC_SUBCOMM_INTERLACED; 13442a22aa42SDave May sred->subcomm = MPI_COMM_NULL; 13451e07b27eSBarry Smith sred->redfactor = 1; 13461e07b27eSBarry Smith sred->ignore_dm = PETSC_FALSE; 13477c5279cbSDave May sred->ignore_kspcomputeoperators = PETSC_FALSE; 13488d9f7141SDave May sred->use_coarse_dm = PETSC_FALSE; 13491e07b27eSBarry Smith pc->data = (void *)sred; 13501e07b27eSBarry Smith 13511e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope; 13521e07b27eSBarry Smith pc->ops->applytranspose = NULL; 1353f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope; 13541e07b27eSBarry Smith pc->ops->setup = PCSetUp_Telescope; 13551e07b27eSBarry Smith pc->ops->destroy = PCDestroy_Telescope; 13561e07b27eSBarry Smith pc->ops->reset = PCReset_Telescope; 13571e07b27eSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Telescope; 13581e07b27eSBarry Smith pc->ops->view = PCView_Telescope; 13591e07b27eSBarry Smith 13601e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 13611e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 13621e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 13631e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 13641e07b27eSBarry Smith 13659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", PCTelescopeGetKSP_Telescope)); 13669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", PCTelescopeGetSubcommType_Telescope)); 13679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", PCTelescopeSetSubcommType_Telescope)); 13689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", PCTelescopeGetReductionFactor_Telescope)); 13699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", PCTelescopeSetReductionFactor_Telescope)); 13709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", PCTelescopeGetIgnoreDM_Telescope)); 13719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", PCTelescopeSetIgnoreDM_Telescope)); 13729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", PCTelescopeGetIgnoreKSPComputeOperators_Telescope)); 13739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", PCTelescopeSetIgnoreKSPComputeOperators_Telescope)); 13749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", PCTelescopeGetDM_Telescope)); 13759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", PCTelescopeGetUseCoarseDM_Telescope)); 13769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", PCTelescopeSetUseCoarseDM_Telescope)); 13773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13781e07b27eSBarry Smith } 1379