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; 9bf00f589SPatrick Sanan static const char citation[] = 10bf00f589SPatrick Sanan "@inproceedings{MaySananRuppKnepleySmith2016,\n" 11bf00f589SPatrick Sanan " title = {Extreme-Scale Multigrid Components within PETSc},\n" 12bf00f589SPatrick Sanan " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 13bf00f589SPatrick Sanan " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 14bf00f589SPatrick Sanan " series = {PASC '16},\n" 15bf00f589SPatrick Sanan " isbn = {978-1-4503-4126-4},\n" 16bf00f589SPatrick Sanan " location = {Lausanne, Switzerland},\n" 17bf00f589SPatrick Sanan " pages = {5:1--5:12},\n" 18bf00f589SPatrick Sanan " articleno = {5},\n" 19bf00f589SPatrick Sanan " numpages = {12},\n" 20a8d69d7bSBarry Smith " url = {https://doi.acm.org/10.1145/2929908.2929913},\n" 21bf00f589SPatrick Sanan " doi = {10.1145/2929908.2929913},\n" 22bf00f589SPatrick Sanan " acmid = {2929913},\n" 23bf00f589SPatrick Sanan " publisher = {ACM},\n" 24bf00f589SPatrick Sanan " address = {New York, NY, USA},\n" 25bf00f589SPatrick Sanan " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 26bf00f589SPatrick Sanan " year = {2016}\n" 27bf00f589SPatrick Sanan "}\n"; 28bf00f589SPatrick Sanan 291e07b27eSBarry Smith /* 30c5083d92SDave May default setup mode 311e07b27eSBarry Smith 32c5083d92SDave May [1a] scatter to (FORWARD) 331e07b27eSBarry Smith x(comm) -> xtmp(comm) 34c5083d92SDave May [1b] local copy (to) ranks with color = 0 351e07b27eSBarry Smith xred(subcomm) <- xtmp 361e07b27eSBarry Smith 37c5083d92SDave May [2] solve on sub KSP to obtain yred(subcomm) 38c5083d92SDave May 39c5083d92SDave May [3a] local copy (from) ranks with color = 0 401e07b27eSBarry Smith yred(subcomm) --> xtmp 41c5083d92SDave May [2b] scatter from (REVERSE) 421e07b27eSBarry Smith xtmp(comm) -> y(comm) 431e07b27eSBarry Smith */ 441e07b27eSBarry Smith 458d9f7141SDave May /* 46d083f849SBarry Smith Collective[comm_f] 478d9f7141SDave May Notes 488d9f7141SDave May * Using comm_f = MPI_COMM_NULL will result in an error 498d9f7141SDave May * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid. 508d9f7141SDave May * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid. 518d9f7141SDave May */ 528d9f7141SDave May PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool *isvalid) 538d9f7141SDave May { 5457f12427SDave May PetscInt valid = 1; 558d9f7141SDave May MPI_Group group_f,group_c; 568d9f7141SDave May PetscMPIInt count,k,size_f = 0,size_c = 0,size_c_sum = 0; 575bd1e576SStefano Zampini PetscMPIInt *ranks_f,*ranks_c; 588d9f7141SDave May 5957f12427SDave May PetscFunctionBegin; 60*08401ef6SPierre Jolivet PetscCheck(comm_f != MPI_COMM_NULL,PETSC_COMM_SELF,PETSC_ERR_SUP,"comm_f cannot be MPI_COMM_NULL"); 618d9f7141SDave May 629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(comm_f,&group_f)); 638d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(comm_c,&group_c)); 658d9f7141SDave May } 668d9f7141SDave May 679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_f,&size_f)); 688d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_c,&size_c)); 708d9f7141SDave May } 718d9f7141SDave May 728d9f7141SDave May /* check not all comm_c's are NULL */ 738d9f7141SDave May size_c_sum = size_c; 749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE,&size_c_sum,1,MPI_INT,MPI_SUM,comm_f)); 755bd1e576SStefano Zampini if (size_c_sum == 0) valid = 0; 768d9f7141SDave May 778d9f7141SDave May /* check we can map at least 1 rank in comm_c to comm_f */ 789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size_f,&ranks_f)); 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size_c,&ranks_c)); 805bd1e576SStefano Zampini for (k=0; k<size_f; k++) ranks_f[k] = MPI_UNDEFINED; 815bd1e576SStefano Zampini for (k=0; k<size_c; k++) ranks_c[k] = k; 828d9f7141SDave May 8357f12427SDave May /* 8457f12427SDave May MPI_Group_translate_ranks() returns a non-zero exit code if any rank cannot be translated. 8557f12427SDave May I do not want the code to terminate immediately if this occurs, rather I want to throw 8657f12427SDave May the error later (during PCSetUp_Telescope()) via SETERRQ() with a message indicating 8757f12427SDave May that comm_c is not a valid sub-communicator. 889566063dSJacob Faibussowitsch Hence I purposefully do not call PetscCall() after MPI_Group_translate_ranks(). 8957f12427SDave May */ 908d9f7141SDave May count = 0; 918d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 9266b79024SDave May (void)MPI_Group_translate_ranks(group_c,size_c,ranks_c,group_f,ranks_f); 938d9f7141SDave May for (k=0; k<size_f; k++) { 948d9f7141SDave May if (ranks_f[k] == MPI_UNDEFINED) { 958d9f7141SDave May count++; 968d9f7141SDave May } 978d9f7141SDave May } 988d9f7141SDave May } 995bd1e576SStefano Zampini if (count == size_f) valid = 0; 1008d9f7141SDave May 1019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE,&valid,1,MPIU_INT,MPI_MIN,comm_f)); 1025bd1e576SStefano Zampini if (valid == 1) *isvalid = PETSC_TRUE; 1035bd1e576SStefano Zampini else *isvalid = PETSC_FALSE; 1048d9f7141SDave May 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks_f)); 1069566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks_c)); 1079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group_f)); 1088d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 1099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group_c)); 1108d9f7141SDave May } 1118d9f7141SDave May PetscFunctionReturn(0); 1128d9f7141SDave May } 1138d9f7141SDave May 1141e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred) 1151e07b27eSBarry Smith { 116c6a0d831SBarry Smith DM subdm = NULL; 1171e07b27eSBarry Smith 11857f12427SDave May if (!PCTelescope_isActiveRank(sred)) { subdm = NULL; } 1191e07b27eSBarry Smith else { 1201e07b27eSBarry Smith switch (sred->sr_type) { 1211e07b27eSBarry Smith case TELESCOPE_DEFAULT: subdm = NULL; 1221e07b27eSBarry Smith break; 1231e07b27eSBarry Smith case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 1241e07b27eSBarry Smith break; 1251e07b27eSBarry Smith case TELESCOPE_DMPLEX: subdm = NULL; 1261e07b27eSBarry Smith break; 1278d9f7141SDave May case TELESCOPE_COARSEDM: if (sred->ksp) { KSPGetDM(sred->ksp,&subdm); } 1288d9f7141SDave May break; 1291e07b27eSBarry Smith } 1301e07b27eSBarry Smith } 1311e07b27eSBarry Smith return(subdm); 1321e07b27eSBarry Smith } 1331e07b27eSBarry Smith 1341e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 1351e07b27eSBarry Smith { 1361e07b27eSBarry Smith PetscInt m,M,bs,st,ed; 1371e07b27eSBarry Smith Vec x,xred,yred,xtmp; 1381e07b27eSBarry Smith Mat B; 1391e07b27eSBarry Smith MPI_Comm comm,subcomm; 1401e07b27eSBarry Smith VecScatter scatter; 1411e07b27eSBarry Smith IS isin; 14254cd43dcSJunchao Zhang VecType vectype; 1431e07b27eSBarry Smith 1441e07b27eSBarry Smith PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: setup (default)\n")); 1461e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 1471e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1481e07b27eSBarry Smith 1499566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc,NULL,&B)); 1509566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,&M,NULL)); 1519566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(B,&bs)); 1529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B,&x,NULL)); 1539566063dSJacob Faibussowitsch PetscCall(MatGetVecType(B,&vectype)); 1541e07b27eSBarry Smith 1551e07b27eSBarry Smith xred = NULL; 1563ac26c5eSBarry Smith m = 0; 15757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1589566063dSJacob Faibussowitsch PetscCall(VecCreate(subcomm,&xred)); 1599566063dSJacob Faibussowitsch PetscCall(VecSetSizes(xred,PETSC_DECIDE,M)); 1609566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(xred,bs)); 1619566063dSJacob Faibussowitsch PetscCall(VecSetType(xred,vectype)); /* Use the preconditioner matrix's vectype by default */ 1629566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(xred)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xred,&m)); 1641e07b27eSBarry Smith } 1651e07b27eSBarry Smith 1661e07b27eSBarry Smith yred = NULL; 16757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xred,&yred)); 1691e07b27eSBarry Smith } 1701e07b27eSBarry Smith 1719566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&xtmp)); 1729566063dSJacob Faibussowitsch PetscCall(VecSetSizes(xtmp,m,PETSC_DECIDE)); 1739566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(xtmp,bs)); 1749566063dSJacob Faibussowitsch PetscCall(VecSetType(xtmp,vectype)); 1751e07b27eSBarry Smith 17657f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1779566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xred,&st,&ed)); 1789566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm,(ed-st),st,1,&isin)); 1791e07b27eSBarry Smith } else { 1809566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&st,&ed)); 1819566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm,0,st,1,&isin)); 1821e07b27eSBarry Smith } 1839566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(isin,bs)); 1841e07b27eSBarry Smith 1859566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,isin,xtmp,NULL,&scatter)); 1861e07b27eSBarry Smith 1871e07b27eSBarry Smith sred->isin = isin; 1881e07b27eSBarry Smith sred->scatter = scatter; 1891e07b27eSBarry Smith sred->xred = xred; 1901e07b27eSBarry Smith sred->yred = yred; 1911e07b27eSBarry Smith sred->xtmp = xtmp; 1929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1931e07b27eSBarry Smith PetscFunctionReturn(0); 1941e07b27eSBarry Smith } 1951e07b27eSBarry Smith 1961e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 1971e07b27eSBarry Smith { 1981e07b27eSBarry Smith MPI_Comm comm,subcomm; 1991e07b27eSBarry Smith Mat Bred,B; 2005624f943SPierre Jolivet PetscInt nr,nc,bs; 2011e07b27eSBarry Smith IS isrow,iscol; 2021e07b27eSBarry Smith Mat Blocal,*_Blocal; 2031e07b27eSBarry Smith 2041e07b27eSBarry Smith PetscFunctionBegin; 2059566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n")); 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 2071e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2089566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc,NULL,&B)); 2099566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,&nr,&nc)); 2101e07b27eSBarry Smith isrow = sred->isin; 2119566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,nc,0,1,&iscol)); 2129566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(iscol)); 2139566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(B,NULL,&bs)); 2149566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(iscol,bs)); 2159566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_SUBMAT_SINGLEIS,PETSC_TRUE)); 2169566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal)); 2171e07b27eSBarry Smith Blocal = *_Blocal; 2189566063dSJacob Faibussowitsch PetscCall(PetscFree(_Blocal)); 2191e07b27eSBarry Smith Bred = NULL; 22057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2211e07b27eSBarry Smith PetscInt mm; 2221e07b27eSBarry Smith 2231e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 2241e07b27eSBarry Smith 2259566063dSJacob Faibussowitsch PetscCall(MatGetSize(Blocal,&mm,NULL)); 2269566063dSJacob Faibussowitsch PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred)); 2271e07b27eSBarry Smith } 2281e07b27eSBarry Smith *A = Bred; 2299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Blocal)); 2311e07b27eSBarry Smith PetscFunctionReturn(0); 2321e07b27eSBarry Smith } 2331e07b27eSBarry Smith 234392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 2351e07b27eSBarry Smith { 2361e07b27eSBarry Smith PetscBool has_const; 2371e07b27eSBarry Smith const Vec *vecs; 238c41e779fSDave May Vec *sub_vecs = NULL; 239392968a1SPatrick Sanan PetscInt i,k,n = 0; 2401e07b27eSBarry Smith MPI_Comm subcomm; 2411e07b27eSBarry Smith 2421e07b27eSBarry Smith PetscFunctionBegin; 2431e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2449566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs)); 2451e07b27eSBarry Smith 24657f12427SDave May if (PCTelescope_isActiveRank(sred)) { 247e3acf2f7SBarry Smith if (n) { 2489566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(sred->xred,n,&sub_vecs)); 2491e07b27eSBarry Smith } 2501e07b27eSBarry Smith } 2511e07b27eSBarry Smith 2521e07b27eSBarry Smith /* copy entries */ 2531e07b27eSBarry Smith for (k=0; k<n; k++) { 2541e07b27eSBarry Smith const PetscScalar *x_array; 2551e07b27eSBarry Smith PetscScalar *LA_sub_vec; 25613c30530SDave May PetscInt st,ed; 2571e07b27eSBarry Smith 2581e07b27eSBarry Smith /* pull in vector x->xtmp */ 2599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD)); 2609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD)); 26147856c66SBarry Smith if (sub_vecs) { 262a04a6428SPatrick Sanan /* copy vector entries into xred */ 2639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(sred->xtmp,&x_array)); 264ea2b237eSDave May if (sub_vecs[k]) { 2659566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(sub_vecs[k],&st,&ed)); 2669566063dSJacob Faibussowitsch PetscCall(VecGetArray(sub_vecs[k],&LA_sub_vec)); 2671e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 2681e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 2691e07b27eSBarry Smith } 2709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(sub_vecs[k],&LA_sub_vec)); 2711e07b27eSBarry Smith } 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(sred->xtmp,&x_array)); 2731e07b27eSBarry Smith } 27447856c66SBarry Smith } 2751e07b27eSBarry Smith 27657f12427SDave May if (PCTelescope_isActiveRank(sred)) { 277d8b9d5b7SPatrick Sanan /* create new (near) nullspace for redundant object */ 2789566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(n,&sub_vecs)); 28028b400f6SJacob Faibussowitsch PetscCheck(!nullspace->remove,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 28128b400f6SJacob Faibussowitsch PetscCheck(!nullspace->rmctx,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope"); 282d8b9d5b7SPatrick Sanan } 283392968a1SPatrick Sanan PetscFunctionReturn(0); 284392968a1SPatrick Sanan } 285392968a1SPatrick Sanan 286392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 287392968a1SPatrick Sanan { 288392968a1SPatrick Sanan Mat B; 289392968a1SPatrick Sanan 290392968a1SPatrick Sanan PetscFunctionBegin; 2919566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc,NULL,&B)); 292392968a1SPatrick Sanan /* Propagate the nullspace if it exists */ 293392968a1SPatrick Sanan { 294392968a1SPatrick Sanan MatNullSpace nullspace,sub_nullspace; 2959566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(B,&nullspace)); 296392968a1SPatrick Sanan if (nullspace) { 2979566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: generating nullspace (default)\n")); 2989566063dSJacob Faibussowitsch PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace)); 29957f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3009566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(sub_mat,sub_nullspace)); 3019566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sub_nullspace)); 3021e07b27eSBarry Smith } 303392968a1SPatrick Sanan } 304392968a1SPatrick Sanan } 305392968a1SPatrick Sanan /* Propagate the near nullspace if it exists */ 306392968a1SPatrick Sanan { 307392968a1SPatrick Sanan MatNullSpace nearnullspace,sub_nearnullspace; 3089566063dSJacob Faibussowitsch PetscCall(MatGetNearNullSpace(B,&nearnullspace)); 309392968a1SPatrick Sanan if (nearnullspace) { 3109566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n")); 3119566063dSJacob Faibussowitsch PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace)); 31257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3139566063dSJacob Faibussowitsch PetscCall(MatSetNearNullSpace(sub_mat,sub_nearnullspace)); 3149566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sub_nearnullspace)); 315392968a1SPatrick Sanan } 316392968a1SPatrick Sanan } 317392968a1SPatrick Sanan } 3181e07b27eSBarry Smith PetscFunctionReturn(0); 3191e07b27eSBarry Smith } 3201e07b27eSBarry Smith 3211e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 3221e07b27eSBarry Smith { 3231e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 3241e07b27eSBarry Smith PetscBool iascii,isstring; 3251e07b27eSBarry Smith PetscViewer subviewer; 3261e07b27eSBarry Smith 3271e07b27eSBarry Smith PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 3301e07b27eSBarry Smith if (iascii) { 3318d9f7141SDave May { 3321e07b27eSBarry Smith MPI_Comm comm,subcomm; 3331e07b27eSBarry Smith PetscMPIInt comm_size,subcomm_size; 3348d9f7141SDave May DM dm = NULL,subdm = NULL; 3351e07b27eSBarry Smith 3369566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc,&dm)); 3371e07b27eSBarry Smith subdm = private_PCTelescopeGetSubDM(sred); 3388d9f7141SDave May 3398d9f7141SDave May if (sred->psubcomm) { 3401e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 3411e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 3429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&comm_size)); 3439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm,&subcomm_size)); 3441e07b27eSBarry Smith 3459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent comm size reduction factor = %D\n",sred->redfactor)); 3479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size)); 34848a10b22SPatrick Sanan switch (sred->subcommtype) { 34948a10b22SPatrick Sanan case PETSC_SUBCOMM_INTERLACED : 3509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"petsc subcomm: type = interlaced\n",sred->subcommtype)); 35148a10b22SPatrick Sanan break; 35248a10b22SPatrick Sanan case PETSC_SUBCOMM_CONTIGUOUS : 3539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"petsc subcomm type = contiguous\n",sred->subcommtype)); 35448a10b22SPatrick Sanan break; 35548a10b22SPatrick Sanan default : 35648a10b22SPatrick Sanan SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope"); 35748a10b22SPatrick Sanan } 3589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3598d9f7141SDave May } else { 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 3618d9f7141SDave May subcomm = sred->subcomm; 36257f12427SDave May if (!PCTelescope_isActiveRank(sred)) { 3638d9f7141SDave May subcomm = PETSC_COMM_SELF; 3648d9f7141SDave May } 3658d9f7141SDave May 3669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"subcomm: using user provided sub-communicator\n")); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3698d9f7141SDave May } 3708d9f7141SDave May 3719566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,subcomm,&subviewer)); 37257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(subviewer)); 3741e07b27eSBarry Smith 3751e07b27eSBarry Smith if (dm && sred->ignore_dm) { 3769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"ignoring DM\n")); 3771e07b27eSBarry Smith } 3787c5279cbSDave May if (sred->ignore_kspcomputeoperators) { 3799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"ignoring KSPComputeOperators\n")); 3807c5279cbSDave May } 3811e07b27eSBarry Smith switch (sred->sr_type) { 3821e07b27eSBarry Smith case TELESCOPE_DEFAULT: 3839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"setup type: default\n")); 3841e07b27eSBarry Smith break; 3851e07b27eSBarry Smith case TELESCOPE_DMDA: 3869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"setup type: DMDA auto-repartitioning\n")); 3879566063dSJacob Faibussowitsch PetscCall(DMView_DA_Short(subdm,subviewer)); 3881e07b27eSBarry Smith break; 3891e07b27eSBarry Smith case TELESCOPE_DMPLEX: 3909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"setup type: DMPLEX auto-repartitioning\n")); 3911e07b27eSBarry Smith break; 3928d9f7141SDave May case TELESCOPE_COARSEDM: 3939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"setup type: coarse DM\n")); 3948d9f7141SDave May break; 3958d9f7141SDave May } 3968d9f7141SDave May 3978d9f7141SDave May if (dm) { 3988d9f7141SDave May PetscObject obj = (PetscObject)dm; 3999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Parent DM object:")); 4008d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 4018d9f7141SDave May if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 4028d9f7141SDave May if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 4038d9f7141SDave May if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 4049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"\n")); 4058d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 4068d9f7141SDave May } else { 4079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Parent DM object: NULL\n")); 4088d9f7141SDave May } 4098d9f7141SDave May if (subdm) { 4108d9f7141SDave May PetscObject obj = (PetscObject)subdm; 4119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Sub DM object:")); 4128d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 4138d9f7141SDave May if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 4148d9f7141SDave May if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 4158d9f7141SDave May if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 4169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"\n")); 4178d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 4188d9f7141SDave May } else { 4199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(subviewer,"Sub DM object: NULL\n")); 4201e07b27eSBarry Smith } 4211e07b27eSBarry Smith 4229566063dSJacob Faibussowitsch PetscCall(KSPView(sred->ksp,subviewer)); 4239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(subviewer)); 4241e07b27eSBarry Smith } 4259566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer)); 4261e07b27eSBarry Smith } 4271e07b27eSBarry Smith } 4281e07b27eSBarry Smith PetscFunctionReturn(0); 4291e07b27eSBarry Smith } 4301e07b27eSBarry Smith 4311e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc) 4321e07b27eSBarry Smith { 4331e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 434bd49479cSSatish Balay MPI_Comm comm,subcomm=0; 4351e07b27eSBarry Smith PCTelescopeType sr_type; 4361e07b27eSBarry Smith 4371e07b27eSBarry Smith PetscFunctionBegin; 4389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 4391e07b27eSBarry Smith 4401e07b27eSBarry Smith /* Determine type of setup/update */ 4411e07b27eSBarry Smith if (!pc->setupcalled) { 4421e07b27eSBarry Smith PetscBool has_dm,same; 4431e07b27eSBarry Smith DM dm; 4441e07b27eSBarry Smith 4451e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4461e07b27eSBarry Smith has_dm = PETSC_FALSE; 4479566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc,&dm)); 4481e07b27eSBarry Smith if (dm) { has_dm = PETSC_TRUE; } 4491e07b27eSBarry Smith if (has_dm) { 4501e07b27eSBarry Smith /* check for dmda */ 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMDA,&same)); 4521e07b27eSBarry Smith if (same) { 4539566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: found DMDA\n")); 4541e07b27eSBarry Smith sr_type = TELESCOPE_DMDA; 4551e07b27eSBarry Smith } 4561e07b27eSBarry Smith /* check for dmplex */ 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same)); 4581e07b27eSBarry Smith if (same) { 4599566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: found DMPLEX\n")); 4601e07b27eSBarry Smith sr_type = TELESCOPE_DMPLEX; 4611e07b27eSBarry Smith } 4628d9f7141SDave May 4638d9f7141SDave May if (sred->use_coarse_dm) { 4649566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: using coarse DM\n")); 4658d9f7141SDave May sr_type = TELESCOPE_COARSEDM; 4661e07b27eSBarry Smith } 4671e07b27eSBarry Smith 4681e07b27eSBarry Smith if (sred->ignore_dm) { 4699566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: ignoring DM\n")); 4701e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4711e07b27eSBarry Smith } 4728d9f7141SDave May } 4731e07b27eSBarry Smith sred->sr_type = sr_type; 4741e07b27eSBarry Smith } else { 4751e07b27eSBarry Smith sr_type = sred->sr_type; 4761e07b27eSBarry Smith } 4771e07b27eSBarry Smith 478d8b9d5b7SPatrick Sanan /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */ 4791e07b27eSBarry Smith switch (sr_type) { 4801e07b27eSBarry Smith case TELESCOPE_DEFAULT: 4811e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 4821e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 4831e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 4841e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 4851e07b27eSBarry Smith break; 4861e07b27eSBarry Smith case TELESCOPE_DMDA: 4871e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope_dmda; 488f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_dmda; 4891e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 4901e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 4911e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 4921e07b27eSBarry Smith sred->pctelescope_reset_type = PCReset_Telescope_dmda; 4931e07b27eSBarry Smith break; 494b458e8f1SJose E. Roman case TELESCOPE_DMPLEX: 495b458e8f1SJose E. Roman SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available"); 4968d9f7141SDave May case TELESCOPE_COARSEDM: 4978d9f7141SDave May pc->ops->apply = PCApply_Telescope_CoarseDM; 4988d9f7141SDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_CoarseDM; 4998d9f7141SDave May sred->pctelescope_setup_type = PCTelescopeSetUp_CoarseDM; 5008d9f7141SDave May sred->pctelescope_matcreate_type = NULL; 5018d9f7141SDave May sred->pctelescope_matnullspacecreate_type = NULL; /* PCTelescopeMatNullSpaceCreate_CoarseDM; */ 5028d9f7141SDave May sred->pctelescope_reset_type = PCReset_Telescope_CoarseDM; 5031e07b27eSBarry Smith break; 504b458e8f1SJose E. Roman default: 505b458e8f1SJose E. Roman SETERRQ(comm,PETSC_ERR_SUP,"Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM"); 5068d9f7141SDave May } 5078d9f7141SDave May 5088d9f7141SDave May /* subcomm definition */ 5098d9f7141SDave May if (!pc->setupcalled) { 5108d9f7141SDave May if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) { 5118d9f7141SDave May if (!sred->psubcomm) { 5129566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(comm,&sred->psubcomm)); 5139566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(sred->psubcomm,sred->redfactor)); 5149566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(sred->psubcomm,sred->subcommtype)); 5159566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm))); 5168d9f7141SDave May sred->subcomm = PetscSubcommChild(sred->psubcomm); 5178d9f7141SDave May } 5188d9f7141SDave May } else { /* query PC for DM, check communicators */ 5198d9f7141SDave May DM dm,dm_coarse_partition = NULL; 5208d9f7141SDave May MPI_Comm comm_fine,comm_coarse_partition = MPI_COMM_NULL; 5218d9f7141SDave May PetscMPIInt csize_fine=0,csize_coarse_partition=0,cs[2],csg[2],cnt=0; 5228d9f7141SDave May PetscBool isvalidsubcomm; 5238d9f7141SDave May 5249566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc,&dm)); 5258d9f7141SDave May comm_fine = PetscObjectComm((PetscObject)dm); 5269566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm,&dm_coarse_partition)); 5278d9f7141SDave May if (dm_coarse_partition) { cnt = 1; } 5289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE,&cnt,1,MPI_INT,MPI_SUM,comm_fine)); 529*08401ef6SPierre Jolivet PetscCheck(cnt != 0,comm_fine,PETSC_ERR_SUP,"Zero instances of a coarse DM were found"); 5308d9f7141SDave May 5319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_fine,&csize_fine)); 5328d9f7141SDave May if (dm_coarse_partition) { 5338d9f7141SDave May comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition); 5349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm_coarse_partition,&csize_coarse_partition)); 5358d9f7141SDave May } 5368d9f7141SDave May 5378d9f7141SDave May cs[0] = csize_fine; 5388d9f7141SDave May cs[1] = csize_coarse_partition; 5399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(cs,csg,2,MPI_INT,MPI_MAX,comm_fine)); 540*08401ef6SPierre 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"); 5418d9f7141SDave May 5429566063dSJacob Faibussowitsch PetscCall(PCTelescopeTestValidSubcomm(comm_fine,comm_coarse_partition,&isvalidsubcomm)); 54328b400f6SJacob Faibussowitsch PetscCheck(isvalidsubcomm,comm_fine,PETSC_ERR_SUP,"Coarse DM communicator is not a sub-communicator of parentDM->comm"); 5448d9f7141SDave May sred->subcomm = comm_coarse_partition; 5458d9f7141SDave May } 5468d9f7141SDave May } 5478d9f7141SDave May subcomm = sred->subcomm; 5488d9f7141SDave May 5498d9f7141SDave May /* internal KSP */ 5508d9f7141SDave May if (!pc->setupcalled) { 5518d9f7141SDave May const char *prefix; 5528d9f7141SDave May 55357f12427SDave May if (PCTelescope_isActiveRank(sred)) { 5549566063dSJacob Faibussowitsch PetscCall(KSPCreate(subcomm,&sred->ksp)); 5559566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure)); 5569566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1)); 5579566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp)); 5589566063dSJacob Faibussowitsch PetscCall(KSPSetType(sred->ksp,KSPPREONLY)); 5599566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 5609566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(sred->ksp,prefix)); 5619566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(sred->ksp,"telescope_")); 5628d9f7141SDave May } 5631e07b27eSBarry Smith } 5641e07b27eSBarry Smith 5651e07b27eSBarry Smith /* setup */ 566aaa7a805SDave May if (!pc->setupcalled && sred->pctelescope_setup_type) { 5679566063dSJacob Faibussowitsch PetscCall(sred->pctelescope_setup_type(pc,sred)); 5681e07b27eSBarry Smith } 5691e07b27eSBarry Smith /* update */ 5701e07b27eSBarry Smith if (!pc->setupcalled) { 5711e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 5729566063dSJacob Faibussowitsch PetscCall(sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred)); 5731e07b27eSBarry Smith } 5741e07b27eSBarry Smith if (sred->pctelescope_matnullspacecreate_type) { 5759566063dSJacob Faibussowitsch PetscCall(sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred)); 5761e07b27eSBarry Smith } 5771e07b27eSBarry Smith } else { 5781e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 5799566063dSJacob Faibussowitsch PetscCall(sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred)); 5801e07b27eSBarry Smith } 5811e07b27eSBarry Smith } 5821e07b27eSBarry Smith 5831e07b27eSBarry Smith /* common - no construction */ 58457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 5859566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(sred->ksp,sred->Bred,sred->Bred)); 5861e07b27eSBarry Smith if (pc->setfromoptionscalled && !pc->setupcalled) { 5879566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(sred->ksp)); 5881e07b27eSBarry Smith } 5891e07b27eSBarry Smith } 5901e07b27eSBarry Smith PetscFunctionReturn(0); 5911e07b27eSBarry Smith } 5921e07b27eSBarry Smith 5931e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 5941e07b27eSBarry Smith { 5951e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5961e07b27eSBarry Smith Vec xtmp,xred,yred; 59713c30530SDave May PetscInt i,st,ed; 5981e07b27eSBarry Smith VecScatter scatter; 5991e07b27eSBarry Smith PetscScalar *array; 6001e07b27eSBarry Smith const PetscScalar *x_array; 6011e07b27eSBarry Smith 6021e07b27eSBarry Smith PetscFunctionBegin; 6039566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation,&cited)); 604bf00f589SPatrick Sanan 6051e07b27eSBarry Smith xtmp = sred->xtmp; 6061e07b27eSBarry Smith scatter = sred->scatter; 6071e07b27eSBarry Smith xred = sred->xred; 6081e07b27eSBarry Smith yred = sred->yred; 6091e07b27eSBarry Smith 6101e07b27eSBarry Smith /* pull in vector x->xtmp */ 6119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD)); 6131e07b27eSBarry Smith 614bf00f589SPatrick Sanan /* copy vector entries into xred */ 6159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xtmp,&x_array)); 6161e07b27eSBarry Smith if (xred) { 6171e07b27eSBarry Smith PetscScalar *LA_xred; 6189566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(xred,&st,&ed)); 6199566063dSJacob Faibussowitsch PetscCall(VecGetArray(xred,&LA_xred)); 6201e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 6211e07b27eSBarry Smith LA_xred[i] = x_array[i]; 6221e07b27eSBarry Smith } 6239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xred,&LA_xred)); 6241e07b27eSBarry Smith } 6259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xtmp,&x_array)); 6261e07b27eSBarry Smith /* solve */ 62757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 6289566063dSJacob Faibussowitsch PetscCall(KSPSolve(sred->ksp,xred,yred)); 6299566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(sred->ksp,pc,yred)); 6301e07b27eSBarry Smith } 6311e07b27eSBarry Smith /* return vector */ 6329566063dSJacob Faibussowitsch PetscCall(VecGetArray(xtmp,&array)); 6331e07b27eSBarry Smith if (yred) { 6341e07b27eSBarry Smith const PetscScalar *LA_yred; 6359566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(yred,&st,&ed)); 6369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(yred,&LA_yred)); 6371e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 6381e07b27eSBarry Smith array[i] = LA_yred[i]; 6391e07b27eSBarry Smith } 6409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(yred,&LA_yred)); 6411e07b27eSBarry Smith } 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xtmp,&array)); 6439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE)); 6449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE)); 6451e07b27eSBarry Smith PetscFunctionReturn(0); 6461e07b27eSBarry Smith } 6471e07b27eSBarry Smith 648f650675bSDave May 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) 649f650675bSDave May { 650f650675bSDave May PC_Telescope sred = (PC_Telescope)pc->data; 651a1d91a28SDave May Vec xtmp,yred; 652f650675bSDave May PetscInt i,st,ed; 653f650675bSDave May VecScatter scatter; 654f650675bSDave May const PetscScalar *x_array; 655f650675bSDave May PetscBool default_init_guess_value; 656f650675bSDave May 657f650675bSDave May PetscFunctionBegin; 658f650675bSDave May xtmp = sred->xtmp; 659f650675bSDave May scatter = sred->scatter; 660f650675bSDave May yred = sred->yred; 661f650675bSDave May 662*08401ef6SPierre Jolivet PetscCheck(its <= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1"); 663f650675bSDave May *reason = (PCRichardsonConvergedReason)0; 664f650675bSDave May 665f650675bSDave May if (!zeroguess) { 6669566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n")); 667f650675bSDave May /* pull in vector y->xtmp */ 6689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD)); 6699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD)); 670f650675bSDave May 671bf00f589SPatrick Sanan /* copy vector entries into xred */ 6729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xtmp,&x_array)); 673f650675bSDave May if (yred) { 674f650675bSDave May PetscScalar *LA_yred; 6759566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(yred,&st,&ed)); 6769566063dSJacob Faibussowitsch PetscCall(VecGetArray(yred,&LA_yred)); 677f650675bSDave May for (i=0; i<ed-st; i++) { 678f650675bSDave May LA_yred[i] = x_array[i]; 679f650675bSDave May } 6809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(yred,&LA_yred)); 681f650675bSDave May } 6829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xtmp,&x_array)); 683f650675bSDave May } 684f650675bSDave May 68557f12427SDave May if (PCTelescope_isActiveRank(sred)) { 6869566063dSJacob Faibussowitsch PetscCall(KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value)); 6879566063dSJacob Faibussowitsch if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE)); 688f650675bSDave May } 689f650675bSDave May 6909566063dSJacob Faibussowitsch PetscCall(PCApply_Telescope(pc,x,y)); 691f650675bSDave May 69257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 6939566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value)); 694f650675bSDave May } 695f650675bSDave May 696f650675bSDave May if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 697f650675bSDave May *outits = 1; 698f650675bSDave May PetscFunctionReturn(0); 699f650675bSDave May } 700f650675bSDave May 7011e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc) 7021e07b27eSBarry Smith { 7031e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7041e07b27eSBarry Smith 705362febeeSStefano Zampini PetscFunctionBegin; 7069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&sred->isin)); 7079566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&sred->scatter)); 7089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->xred)); 7099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->yred)); 7109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sred->xtmp)); 7119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sred->Bred)); 7129566063dSJacob Faibussowitsch PetscCall(KSPReset(sred->ksp)); 7131e07b27eSBarry Smith if (sred->pctelescope_reset_type) { 7149566063dSJacob Faibussowitsch PetscCall(sred->pctelescope_reset_type(pc)); 7151e07b27eSBarry Smith } 7161e07b27eSBarry Smith PetscFunctionReturn(0); 7171e07b27eSBarry Smith } 7181e07b27eSBarry Smith 7191e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc) 7201e07b27eSBarry Smith { 7211e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7221e07b27eSBarry Smith 7231e07b27eSBarry Smith PetscFunctionBegin; 7249566063dSJacob Faibussowitsch PetscCall(PCReset_Telescope(pc)); 7259566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&sred->ksp)); 7269566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&sred->psubcomm)); 7279566063dSJacob Faibussowitsch PetscCall(PetscFree(sred->dm_ctx)); 7289566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 7291e07b27eSBarry Smith PetscFunctionReturn(0); 7301e07b27eSBarry Smith } 7311e07b27eSBarry Smith 7324416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc) 7331e07b27eSBarry Smith { 7341e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7351e07b27eSBarry Smith MPI_Comm comm; 7361e07b27eSBarry Smith PetscMPIInt size; 73748a10b22SPatrick Sanan PetscBool flg; 73848a10b22SPatrick Sanan PetscSubcommType subcommtype; 7391e07b27eSBarry Smith 7401e07b27eSBarry Smith PetscFunctionBegin; 7419566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 7429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 7439566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"Telescope options")); 7449566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg)); 74548a10b22SPatrick Sanan if (flg) { 7469566063dSJacob Faibussowitsch PetscCall(PCTelescopeSetSubcommType(pc,subcommtype)); 74748a10b22SPatrick Sanan } 7489566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,NULL)); 749*08401ef6SPierre Jolivet PetscCheck(sred->redfactor <= size,comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 7509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,NULL)); 7519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,NULL)); 7529566063dSJacob 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)); 7539566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 7541e07b27eSBarry Smith PetscFunctionReturn(0); 7551e07b27eSBarry Smith } 7561e07b27eSBarry Smith 7571e07b27eSBarry Smith /* PC simplementation specific API's */ 7581e07b27eSBarry Smith 7591e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 7601e07b27eSBarry Smith { 7611e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 762bd49479cSSatish Balay PetscFunctionBegin; 7631e07b27eSBarry Smith if (ksp) *ksp = red->ksp; 764bd49479cSSatish Balay PetscFunctionReturn(0); 7651e07b27eSBarry Smith } 7661e07b27eSBarry Smith 76748a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype) 76848a10b22SPatrick Sanan { 76948a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 77048a10b22SPatrick Sanan PetscFunctionBegin; 77148a10b22SPatrick Sanan if (subcommtype) *subcommtype = red->subcommtype; 77248a10b22SPatrick Sanan PetscFunctionReturn(0); 77348a10b22SPatrick Sanan } 77448a10b22SPatrick Sanan 77548a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype) 77648a10b22SPatrick Sanan { 77748a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 77848a10b22SPatrick Sanan 77948a10b22SPatrick Sanan PetscFunctionBegin; 78028b400f6SJacob 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."); 78148a10b22SPatrick Sanan red->subcommtype = subcommtype; 78248a10b22SPatrick Sanan PetscFunctionReturn(0); 78348a10b22SPatrick Sanan } 78448a10b22SPatrick Sanan 7851e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 7861e07b27eSBarry Smith { 7871e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 788bd49479cSSatish Balay PetscFunctionBegin; 7891e07b27eSBarry Smith if (fact) *fact = red->redfactor; 790bd49479cSSatish Balay PetscFunctionReturn(0); 7911e07b27eSBarry Smith } 7921e07b27eSBarry Smith 7931e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 7941e07b27eSBarry Smith { 7951e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 7961e07b27eSBarry Smith PetscMPIInt size; 7971e07b27eSBarry Smith 798bd49479cSSatish Balay PetscFunctionBegin; 7999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 800*08401ef6SPierre Jolivet PetscCheck(fact > 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 801*08401ef6SPierre Jolivet PetscCheck(fact <= size,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 8021e07b27eSBarry Smith red->redfactor = fact; 803bd49479cSSatish Balay PetscFunctionReturn(0); 8041e07b27eSBarry Smith } 8051e07b27eSBarry Smith 8061e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 8071e07b27eSBarry Smith { 8081e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 809bd49479cSSatish Balay PetscFunctionBegin; 8101e07b27eSBarry Smith if (v) *v = red->ignore_dm; 811bd49479cSSatish Balay PetscFunctionReturn(0); 8121e07b27eSBarry Smith } 81348a10b22SPatrick Sanan 8141e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 8151e07b27eSBarry Smith { 8161e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 817bd49479cSSatish Balay PetscFunctionBegin; 8181e07b27eSBarry Smith red->ignore_dm = v; 819bd49479cSSatish Balay PetscFunctionReturn(0); 8201e07b27eSBarry Smith } 8211e07b27eSBarry Smith 8228d9f7141SDave May static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool *v) 8238d9f7141SDave May { 8248d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 8258d9f7141SDave May PetscFunctionBegin; 8268d9f7141SDave May if (v) *v = red->use_coarse_dm; 8278d9f7141SDave May PetscFunctionReturn(0); 8288d9f7141SDave May } 8298d9f7141SDave May 8308d9f7141SDave May static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v) 8318d9f7141SDave May { 8328d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 8338d9f7141SDave May PetscFunctionBegin; 8348d9f7141SDave May red->use_coarse_dm = v; 8358d9f7141SDave May PetscFunctionReturn(0); 8368d9f7141SDave May } 8378d9f7141SDave May 8380ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v) 8390ae7c45bSDave May { 8400ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8410ae7c45bSDave May PetscFunctionBegin; 8420ae7c45bSDave May if (v) *v = red->ignore_kspcomputeoperators; 8430ae7c45bSDave May PetscFunctionReturn(0); 8440ae7c45bSDave May } 84548a10b22SPatrick Sanan 8460ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v) 8470ae7c45bSDave May { 8480ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8490ae7c45bSDave May PetscFunctionBegin; 8500ae7c45bSDave May red->ignore_kspcomputeoperators = v; 8510ae7c45bSDave May PetscFunctionReturn(0); 8520ae7c45bSDave May } 8530ae7c45bSDave May 8541e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 8551e07b27eSBarry Smith { 8561e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 857bd49479cSSatish Balay PetscFunctionBegin; 8581e07b27eSBarry Smith *dm = private_PCTelescopeGetSubDM(red); 859bd49479cSSatish Balay PetscFunctionReturn(0); 8601e07b27eSBarry Smith } 8611e07b27eSBarry Smith 8621e07b27eSBarry Smith /*@ 8631e07b27eSBarry Smith PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 8641e07b27eSBarry Smith 8651e07b27eSBarry Smith Not Collective 8661e07b27eSBarry Smith 8671e07b27eSBarry Smith Input Parameter: 8681e07b27eSBarry Smith . pc - the preconditioner context 8691e07b27eSBarry Smith 8701e07b27eSBarry Smith Output Parameter: 8711e07b27eSBarry Smith . subksp - the KSP defined the smaller set of processes 8721e07b27eSBarry Smith 8731e07b27eSBarry Smith Level: advanced 8741e07b27eSBarry Smith 8751e07b27eSBarry Smith @*/ 8761e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 8771e07b27eSBarry Smith { 878bd49479cSSatish Balay PetscFunctionBegin; 879cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp)); 880bd49479cSSatish Balay PetscFunctionReturn(0); 8811e07b27eSBarry Smith } 8821e07b27eSBarry Smith 8831e07b27eSBarry Smith /*@ 8841e07b27eSBarry Smith PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 8851e07b27eSBarry Smith 8861e07b27eSBarry Smith Not Collective 8871e07b27eSBarry Smith 8881e07b27eSBarry Smith Input Parameter: 8891e07b27eSBarry Smith . pc - the preconditioner context 8901e07b27eSBarry Smith 8911e07b27eSBarry Smith Output Parameter: 8921e07b27eSBarry Smith . fact - the reduction factor 8931e07b27eSBarry Smith 8941e07b27eSBarry Smith Level: advanced 8951e07b27eSBarry Smith 8961e07b27eSBarry Smith @*/ 8971e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 8981e07b27eSBarry Smith { 899bd49479cSSatish Balay PetscFunctionBegin; 900cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact)); 901bd49479cSSatish Balay PetscFunctionReturn(0); 9021e07b27eSBarry Smith } 9031e07b27eSBarry Smith 9041e07b27eSBarry Smith /*@ 9051e07b27eSBarry Smith PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 9061e07b27eSBarry Smith 9071e07b27eSBarry Smith Not Collective 9081e07b27eSBarry Smith 9091e07b27eSBarry Smith Input Parameter: 9101e07b27eSBarry Smith . pc - the preconditioner context 9111e07b27eSBarry Smith 9121e07b27eSBarry Smith Output Parameter: 9131e07b27eSBarry Smith . fact - the reduction factor 9141e07b27eSBarry Smith 9151e07b27eSBarry Smith Level: advanced 9161e07b27eSBarry Smith 9171e07b27eSBarry Smith @*/ 9181e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 9191e07b27eSBarry Smith { 920bd49479cSSatish Balay PetscFunctionBegin; 921cac4c232SBarry Smith PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact)); 922bd49479cSSatish Balay PetscFunctionReturn(0); 9231e07b27eSBarry Smith } 9241e07b27eSBarry Smith 9251e07b27eSBarry Smith /*@ 9261e07b27eSBarry Smith PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 9271e07b27eSBarry Smith 9281e07b27eSBarry Smith Not Collective 9291e07b27eSBarry Smith 9301e07b27eSBarry Smith Input Parameter: 9311e07b27eSBarry Smith . pc - the preconditioner context 9321e07b27eSBarry Smith 9331e07b27eSBarry Smith Output Parameter: 9341e07b27eSBarry Smith . v - the flag 9351e07b27eSBarry Smith 9361e07b27eSBarry Smith Level: advanced 9371e07b27eSBarry Smith 9381e07b27eSBarry Smith @*/ 9391e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 9401e07b27eSBarry Smith { 941bd49479cSSatish Balay PetscFunctionBegin; 942cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v)); 943bd49479cSSatish Balay PetscFunctionReturn(0); 9441e07b27eSBarry Smith } 9451e07b27eSBarry Smith 9461e07b27eSBarry Smith /*@ 9471e07b27eSBarry Smith PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 9481e07b27eSBarry Smith 9491e07b27eSBarry Smith Not Collective 9501e07b27eSBarry Smith 9511e07b27eSBarry Smith Input Parameter: 9521e07b27eSBarry Smith . pc - the preconditioner context 9531e07b27eSBarry Smith 9541e07b27eSBarry Smith Output Parameter: 9551e07b27eSBarry Smith . v - Use PETSC_TRUE to ignore any DM 9561e07b27eSBarry Smith 9571e07b27eSBarry Smith Level: advanced 9581e07b27eSBarry Smith 9591e07b27eSBarry Smith @*/ 960bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 9611e07b27eSBarry Smith { 962bd49479cSSatish Balay PetscFunctionBegin; 963cac4c232SBarry Smith PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v)); 964bd49479cSSatish Balay PetscFunctionReturn(0); 9651e07b27eSBarry Smith } 9661e07b27eSBarry Smith 9671e07b27eSBarry Smith /*@ 9688d9f7141SDave May PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used. 9698d9f7141SDave May 9708d9f7141SDave May Not Collective 9718d9f7141SDave May 9728d9f7141SDave May Input Parameter: 9738d9f7141SDave May . pc - the preconditioner context 9748d9f7141SDave May 9758d9f7141SDave May Output Parameter: 9768d9f7141SDave May . v - the flag 9778d9f7141SDave May 9788d9f7141SDave May Level: advanced 9798d9f7141SDave May 9808d9f7141SDave May @*/ 9818d9f7141SDave May PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc,PetscBool *v) 9828d9f7141SDave May { 9838d9f7141SDave May PetscFunctionBegin; 984cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetUseCoarseDM_C",(PC,PetscBool*),(pc,v)); 9858d9f7141SDave May PetscFunctionReturn(0); 9868d9f7141SDave May } 9878d9f7141SDave May 9888d9f7141SDave May /*@ 9898d9f7141SDave May PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM 9908d9f7141SDave May 9918d9f7141SDave May Not Collective 9928d9f7141SDave May 9938d9f7141SDave May Input Parameter: 9948d9f7141SDave May . pc - the preconditioner context 9958d9f7141SDave May 9968d9f7141SDave May Output Parameter: 997aaa7a805SDave May . v - Use PETSC_FALSE to ignore any coarse DM 9988d9f7141SDave May 9998d9f7141SDave May Notes: 10008d9f7141SDave May When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope 10018d9f7141SDave May will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and 10028d9f7141SDave May -pc_telescope_subcomm_type will no longer have any meaning. 10038d9f7141SDave May It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes. 10048d9f7141SDave May An error will occur of the size of the communicator associated with the coarse DM 10058d9f7141SDave May is the same as that of the parent DM. 10068d9f7141SDave May Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent. 10078d9f7141SDave May This will be checked at the time the preconditioner is setup and an error will occur if 10088d9f7141SDave May the coarse DM does not define a sub-communicator of that used by the parent DM. 10098d9f7141SDave May 10108d9f7141SDave May The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of 10118d9f7141SDave May the DM used (e.g. it supports DMSHELL, DMPLEX, etc). 10128d9f7141SDave May 10138d9f7141SDave May Support is currently only provided for the case when you are using KSPSetComputeOperators() 10148d9f7141SDave May 10158d9f7141SDave May 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 DMs. 10168d9f7141SDave May In the user code, this is achieved via 10178d9f7141SDave May .vb 10188d9f7141SDave May { 10198d9f7141SDave May DM dm_fine; 10208d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method); 10218d9f7141SDave May } 10228d9f7141SDave May .ve 10238d9f7141SDave May The signature of the user provided field scatter method is 10248d9f7141SDave May .vb 10258d9f7141SDave May PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse); 10268d9f7141SDave May .ve 10278d9f7141SDave May The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE. 10288d9f7141SDave May SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM. 10298d9f7141SDave May 10308d9f7141SDave May Optionally, the user may also compose a function with the parent DM to facilitate the transfer 10318d9f7141SDave May of state variables between the fine and coarse DMs. 10328d9f7141SDave May In the context of a finite element discretization, an example state variable might be 10338d9f7141SDave May values associated with quadrature points within each element. 10348d9f7141SDave May A user provided state scatter method is composed via 10358d9f7141SDave May .vb 10368d9f7141SDave May { 10378d9f7141SDave May DM dm_fine; 10388d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method); 10398d9f7141SDave May } 10408d9f7141SDave May .ve 10418d9f7141SDave May The signature of the user provided state scatter method is 10428d9f7141SDave May .vb 10438d9f7141SDave May PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse); 10448d9f7141SDave May .ve 10458d9f7141SDave May SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM. 10468d9f7141SDave May The user is only required to support mode = SCATTER_FORWARD. 10478d9f7141SDave May No assumption is made about the data type of the state variables. 10488d9f7141SDave May These must be managed by the user and must be accessible from the DM. 10498d9f7141SDave May 10508d9f7141SDave May Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be 10518d9f7141SDave May associated with the sub-KSP residing within PCTelescope. 10528d9f7141SDave May In general, PCTelescope assumes that the context on the fine and coarse DM used with 10538d9f7141SDave May KSPSetComputeOperators() should be "similar" in type or origin. 10548d9f7141SDave May Specifically the following rules are used to infer what context on the sub-KSP should be. 10558d9f7141SDave May 10568d9f7141SDave May First the contexts from the KSP and the fine and coarse DMs are retrieved. 10578d9f7141SDave May Note that the special case of a DMSHELL context is queried. 10588d9f7141SDave May 10598d9f7141SDave May .vb 10608d9f7141SDave May DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx); 10618d9f7141SDave May DMGetApplicationContext(dm_fine,&dmfine_appctx); 10628d9f7141SDave May DMShellGetContext(dm_fine,&dmfine_shellctx); 10638d9f7141SDave May 10648d9f7141SDave May DMGetApplicationContext(dm_coarse,&dmcoarse_appctx); 10658d9f7141SDave May DMShellGetContext(dm_coarse,&dmcoarse_shellctx); 10668d9f7141SDave May .ve 10678d9f7141SDave May 10688d9f7141SDave May The following rules are then enforced: 10698d9f7141SDave May 10708d9f7141SDave May 1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP: 10718d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL); 10728d9f7141SDave May 10738d9f7141SDave May 2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx, 10748d9f7141SDave May check that dmcoarse_appctx is also non-NULL. If this is true, then: 10758d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx); 10768d9f7141SDave May 10778d9f7141SDave May 3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx, 10788d9f7141SDave May check that dmcoarse_shellctx is also non-NULL. If this is true, then: 10798d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx); 10808d9f7141SDave May 10818d9f7141SDave May If neither of the above three tests passed, then PCTelescope cannot safely determine what 10828d9f7141SDave May context should be provided to KSPSetComputeOperators() for use with the sub-KSP. 10838d9f7141SDave May In this case, an additional mechanism is provided via a composed function which will return 10848d9f7141SDave May the actual context to be used. To use this feature you must compose the "getter" function 10858d9f7141SDave May with the coarse DM, e.g. 10868d9f7141SDave May .vb 10878d9f7141SDave May { 10888d9f7141SDave May DM dm_coarse; 10898d9f7141SDave May PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter); 10908d9f7141SDave May } 10918d9f7141SDave May .ve 10928d9f7141SDave May The signature of the user provided method is 10938d9f7141SDave May .vb 10948d9f7141SDave May PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext); 10958d9f7141SDave May .ve 10968d9f7141SDave May 10978d9f7141SDave May Level: advanced 10988d9f7141SDave May 10998d9f7141SDave May @*/ 11008d9f7141SDave May PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc,PetscBool v) 11018d9f7141SDave May { 11028d9f7141SDave May PetscFunctionBegin; 1103cac4c232SBarry Smith PetscTryMethod(pc,"PCTelescopeSetUseCoarseDM_C",(PC,PetscBool),(pc,v)); 11048d9f7141SDave May PetscFunctionReturn(0); 11058d9f7141SDave May } 11068d9f7141SDave May 11078d9f7141SDave May /*@ 11080ae7c45bSDave May PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used. 11090ae7c45bSDave May 11100ae7c45bSDave May Not Collective 11110ae7c45bSDave May 11120ae7c45bSDave May Input Parameter: 11130ae7c45bSDave May . pc - the preconditioner context 11140ae7c45bSDave May 11150ae7c45bSDave May Output Parameter: 11160ae7c45bSDave May . v - the flag 11170ae7c45bSDave May 11180ae7c45bSDave May Level: advanced 11190ae7c45bSDave May 11200ae7c45bSDave May @*/ 11210ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v) 11220ae7c45bSDave May { 11230ae7c45bSDave May PetscFunctionBegin; 1124cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v)); 11250ae7c45bSDave May PetscFunctionReturn(0); 11260ae7c45bSDave May } 11270ae7c45bSDave May 11280ae7c45bSDave May /*@ 11290ae7c45bSDave May PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators. 11300ae7c45bSDave May 11310ae7c45bSDave May Not Collective 11320ae7c45bSDave May 11330ae7c45bSDave May Input Parameter: 11340ae7c45bSDave May . pc - the preconditioner context 11350ae7c45bSDave May 11360ae7c45bSDave May Output Parameter: 1137a954d8f4SDave May . v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc 11380ae7c45bSDave May 11390ae7c45bSDave May Level: advanced 11400ae7c45bSDave May 11410ae7c45bSDave May @*/ 11420ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v) 11430ae7c45bSDave May { 11440ae7c45bSDave May PetscFunctionBegin; 1145cac4c232SBarry Smith PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v)); 11460ae7c45bSDave May PetscFunctionReturn(0); 11470ae7c45bSDave May } 11480ae7c45bSDave May 11490ae7c45bSDave May /*@ 11501e07b27eSBarry Smith PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 11511e07b27eSBarry Smith 11521e07b27eSBarry Smith Not Collective 11531e07b27eSBarry Smith 11541e07b27eSBarry Smith Input Parameter: 11551e07b27eSBarry Smith . pc - the preconditioner context 11561e07b27eSBarry Smith 11571e07b27eSBarry Smith Output Parameter: 11581e07b27eSBarry Smith . subdm - The re-partitioned DM 11591e07b27eSBarry Smith 11601e07b27eSBarry Smith Level: advanced 11611e07b27eSBarry Smith 11621e07b27eSBarry Smith @*/ 11631e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 11641e07b27eSBarry Smith { 1165bd49479cSSatish Balay PetscFunctionBegin; 1166cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm)); 1167bd49479cSSatish Balay PetscFunctionReturn(0); 11681e07b27eSBarry Smith } 11691e07b27eSBarry Smith 117048a10b22SPatrick Sanan /*@ 117148a10b22SPatrick Sanan PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous) 117248a10b22SPatrick Sanan 117348a10b22SPatrick Sanan Logically Collective 117448a10b22SPatrick Sanan 1175d8d19677SJose E. Roman Input Parameters: 11761dae98e4SBarry Smith + pc - the preconditioner context 11771dae98e4SBarry Smith - subcommtype - the subcommunicator type (see PetscSubcommType) 117848a10b22SPatrick Sanan 117948a10b22SPatrick Sanan Level: advanced 118048a10b22SPatrick Sanan 118148a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE 118248a10b22SPatrick Sanan @*/ 118348a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) 118448a10b22SPatrick Sanan { 118548a10b22SPatrick Sanan PetscFunctionBegin; 1186cac4c232SBarry Smith PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype)); 118748a10b22SPatrick Sanan PetscFunctionReturn(0); 118848a10b22SPatrick Sanan } 118948a10b22SPatrick Sanan 119048a10b22SPatrick Sanan /*@ 119148a10b22SPatrick Sanan PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous) 119248a10b22SPatrick Sanan 119348a10b22SPatrick Sanan Not Collective 119448a10b22SPatrick Sanan 119548a10b22SPatrick Sanan Input Parameter: 119648a10b22SPatrick Sanan . pc - the preconditioner context 119748a10b22SPatrick Sanan 119848a10b22SPatrick Sanan Output Parameter: 119948a10b22SPatrick Sanan . subcommtype - the subcommunicator type (see PetscSubcommType) 120048a10b22SPatrick Sanan 120148a10b22SPatrick Sanan Level: advanced 120248a10b22SPatrick Sanan 12031dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE 120448a10b22SPatrick Sanan @*/ 12051dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) 120648a10b22SPatrick Sanan { 120748a10b22SPatrick Sanan PetscFunctionBegin; 1208cac4c232SBarry Smith PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype)); 120948a10b22SPatrick Sanan PetscFunctionReturn(0); 121048a10b22SPatrick Sanan } 121148a10b22SPatrick Sanan 12121e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/ 12131e07b27eSBarry Smith /*MC 121400fea0ebSDave May PCTELESCOPE - Runs a KSP solver on a sub-communicator. MPI ranks not in the sub-communicator are idle during the solve. 12151e07b27eSBarry Smith 12161e07b27eSBarry Smith Options Database: 121700fea0ebSDave 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. 121800fea0ebSDave May . -pc_telescope_ignore_dm - flag to indicate whether an attached DM should be ignored. 121900fea0ebSDave May . -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI ranks on the sub-communicator. see PetscSubcomm for more information. 122000fea0ebSDave May . -pc_telescope_ignore_kspcomputeoperators - flag to indicate whether KSPSetComputeOperators should be used on the sub-KSP. 122100fea0ebSDave May - -pc_telescope_use_coarse_dm - flag to indicate whether the coarse DM should be used to define the sub-communicator. 12221e07b27eSBarry Smith 12231e07b27eSBarry Smith Level: advanced 12241e07b27eSBarry Smith 12251e07b27eSBarry Smith Notes: 122600fea0ebSDave May Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 122700fea0ebSDave May creates a child sub-communicator (c') containing fewer MPI ranks than the original parent preconditioner (PC). 12286fc41876SBarry Smith The preconditioner is deemed telescopic as it only calls KSPSolve() on a single 12298439623fSPatrick Sanan sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 123000fea0ebSDave May This means there will be MPI ranks which will be idle during the application of this preconditioner. 123100fea0ebSDave May Additionally, in comparison with PCREDUNDANT, PCTELESCOPE can utilize an attached DM. 12326fc41876SBarry Smith 123300fea0ebSDave May The default type of the sub KSP (the KSP defined on c') is PREONLY. 123400fea0ebSDave May 123500fea0ebSDave May There are three setup mechanisms for PCTelescope. Features support by each type are described below. 123600fea0ebSDave May In the following, we will refer to the operators B and B', these are the Bmat provided to the KSP on the 123700fea0ebSDave May communicators c and c' respectively. 123800fea0ebSDave May 123900fea0ebSDave May [1] Default setup 124000fea0ebSDave May The sub-communicator c' is created via PetscSubcommCreate(). 1241a5b23f4aSJose E. Roman Explicitly defined nullspace and near nullspace vectors will be propagated from B to B'. 124200fea0ebSDave May Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()). 124300fea0ebSDave May No support is provided for KSPSetComputeOperators(). 124400fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 124500fea0ebSDave May 124600fea0ebSDave May [2] DM aware setup 124700fea0ebSDave May If a DM is attached to the PC, it is re-partitioned on the sub-communicator c'. 124800fea0ebSDave May c' is created via PetscSubcommCreate(). 12491e07b27eSBarry Smith Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 12501e07b27eSBarry Smith Currently only support for re-partitioning a DMDA is provided. 125100fea0ebSDave 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'). 125200fea0ebSDave May Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()). 125300fea0ebSDave May Support is provided for KSPSetComputeOperators(). The user provided function and context is propagated to the sub KSP. 125400fea0ebSDave May This is fragile since the user must ensure that their user context is valid for use on c'. 125500fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 12561e07b27eSBarry Smith 125700fea0ebSDave May [3] Coarse DM setup 125800fea0ebSDave May If a DM (dmfine) is attached to the PC, dmfine is queried for a "coarse" DM (call this dmcoarse) via DMGetCoarseDM(). 125900fea0ebSDave May PCTELESCOPE will interpret the coarse DM as being defined on a sub-communicator of c. 126000fea0ebSDave May The communicator associated with dmcoarse will define the c' to be used within PCTELESCOPE. 126100fea0ebSDave May PCTELESCOPE will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported. 126200fea0ebSDave May The intention of this setup type is that PCTELESCOPE will use an existing (e.g. user defined) communicator hierarchy, say as would be 126300fea0ebSDave May available with using multi-grid on unstructured meshes. 126400fea0ebSDave May This setup will not use the command line options -pc_telescope_reduction_factor or -pc_telescope_subcomm_type. 126500fea0ebSDave 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'. 126600fea0ebSDave May Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()). 126700fea0ebSDave May There is no general method to permute field orderings, hence only KSPSetComputeOperators() is supported. 126800fea0ebSDave May The user must use PetscObjectComposeFunction() with dmfine to define the method to scatter fields from dmfine to dmcoarse. 1269a5b23f4aSJose E. Roman 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(). 127000fea0ebSDave May Currently there is no support for the flag -pc_use_amat. 127100fea0ebSDave May This setup can be invoked by the option -pc_telescope_use_coarse_dm or by calling PCTelescopeSetUseCoarseDM(pc,PETSC_TRUE); 127200fea0ebSDave May Further information about the user-provided methods required by this setup type are described here PCTelescopeSetUseCoarseDM(). 12736fc41876SBarry Smith 12746fc41876SBarry Smith Developer Notes: 12756fc41876SBarry Smith During PCSetup, the B operator is scattered onto c'. 12766fc41876SBarry Smith Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 12778439623fSPatrick Sanan Then, KSPSolve() is executed on the c' communicator. 12786fc41876SBarry Smith 12796fc41876SBarry Smith The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 1280a04a6428SPatrick Sanan 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. 12816fc41876SBarry Smith 1282005d9f20SPatrick Sanan The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator. 12838439623fSPatrick Sanan In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into 1284005d9f20SPatrick Sanan a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c') 12856fc41876SBarry Smith 128600fea0ebSDave May The telescoping preconditioner can re-partition an attached DM if it is a DMDA (2D or 3D - 12874c500f23SPierre Jolivet support for 1D DMDAs is not provided). If a DMDA is found, a topologically equivalent DMDA is created on c' 128800fea0ebSDave May and this new DM is attached the sub KSP. The design of telescope is such that it should be possible to extend support 12898439623fSPatrick Sanan for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 129000fea0ebSDave May Alternatively, user-provided re-partitioned DMs can be used via -pc_telescope_use_coarse_dm. 12916fc41876SBarry Smith 129200fea0ebSDave May With the default setup mode, B' is defined by fusing rows (in order) associated with MPI ranks common to c and c'. 12936fc41876SBarry Smith 12948439623fSPatrick Sanan When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B 12957dae84e0SHong Zhang into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the 12966fc41876SBarry Smith locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 12976fc41876SBarry Smith 12988439623fSPatrick Sanan Limitations/improvements include the following. 12998439623fSPatrick Sanan VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage. 130000fea0ebSDave May A unified mechanism to query for user contexts as required by KSPSetComputeOperators() and MatNullSpaceSetFunction(). 13016fc41876SBarry Smith 13020705ae88SPierre Jolivet The symmetric permutation used when a DMDA is encountered is performed via explicitly assembling a permutation matrix P, 13038439623fSPatrick Sanan 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 13040705ae88SPierre Jolivet VecPermute() does not support the use case required here. By computing P, one can permute both the operator and RHS in a 13056fc41876SBarry Smith consistent manner. 13066fc41876SBarry Smith 130700fea0ebSDave May Mapping of vectors (default setup mode) is performed in the following way. 130800fea0ebSDave 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. 13098439623fSPatrick Sanan Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2. 131000fea0ebSDave May We perform the scatter to the sub-communicator in the following way. 1311514db83aSDave May [1] Given a vector x defined on communicator c 13126fc41876SBarry Smith 1313514db83aSDave May .vb 1314514db83aSDave May rank(c) local values of x 1315514db83aSDave May ------- ---------------------------------------- 1316514db83aSDave May 0 [ 0.0, 1.0, 2.0, 3.0, 4.0, 5.0 ] 1317514db83aSDave May 1 [ 6.0, 7.0, 8.0, 9.0, 10.0, 11.0 ] 1318514db83aSDave May 2 [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ] 1319514db83aSDave May 3 [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ] 1320514db83aSDave May .ve 13216fc41876SBarry Smith 1322514db83aSDave May scatter into xtmp defined also on comm c, so that we have the following values 13236fc41876SBarry Smith 1324514db83aSDave May .vb 1325514db83aSDave May rank(c) local values of xtmp 1326514db83aSDave May ------- ---------------------------------------------------------------------------- 1327514db83aSDave 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 ] 1328514db83aSDave May 1 [ ] 1329514db83aSDave 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 ] 1330514db83aSDave May 3 [ ] 1331514db83aSDave May .ve 13326fc41876SBarry Smith 13336fc41876SBarry Smith The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 13346fc41876SBarry Smith 1335514db83aSDave 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'. 13366fc41876SBarry Smith Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 13376fc41876SBarry Smith 1338514db83aSDave May .vb 1339514db83aSDave May rank(c') local values of xred 1340514db83aSDave May -------- ---------------------------------------------------------------------------- 1341514db83aSDave 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 ] 1342514db83aSDave 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 ] 1343514db83aSDave May .ve 13441e07b27eSBarry Smith 13451e07b27eSBarry Smith Contributed by Dave May 13461e07b27eSBarry Smith 1347bf00f589SPatrick Sanan Reference: 1348bf00f589SPatrick 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 1349bf00f589SPatrick Sanan 13506fc41876SBarry Smith .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 13511e07b27eSBarry Smith M*/ 13521e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 13531e07b27eSBarry Smith { 13541e07b27eSBarry Smith struct _PC_Telescope *sred; 13551e07b27eSBarry Smith 13561e07b27eSBarry Smith PetscFunctionBegin; 13579566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&sred)); 13582a22aa42SDave May sred->psubcomm = NULL; 135948a10b22SPatrick Sanan sred->subcommtype = PETSC_SUBCOMM_INTERLACED; 13602a22aa42SDave May sred->subcomm = MPI_COMM_NULL; 13611e07b27eSBarry Smith sred->redfactor = 1; 13621e07b27eSBarry Smith sred->ignore_dm = PETSC_FALSE; 13637c5279cbSDave May sred->ignore_kspcomputeoperators = PETSC_FALSE; 13648d9f7141SDave May sred->use_coarse_dm = PETSC_FALSE; 13651e07b27eSBarry Smith pc->data = (void*)sred; 13661e07b27eSBarry Smith 13671e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope; 13681e07b27eSBarry Smith pc->ops->applytranspose = NULL; 1369f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope; 13701e07b27eSBarry Smith pc->ops->setup = PCSetUp_Telescope; 13711e07b27eSBarry Smith pc->ops->destroy = PCDestroy_Telescope; 13721e07b27eSBarry Smith pc->ops->reset = PCReset_Telescope; 13731e07b27eSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Telescope; 13741e07b27eSBarry Smith pc->ops->view = PCView_Telescope; 13751e07b27eSBarry Smith 13761e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 13771e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 13781e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 13791e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 13801e07b27eSBarry Smith 13819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope)); 13829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope)); 13839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope)); 13849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope)); 13859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope)); 13869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope)); 13879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope)); 13889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope)); 13899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope)); 13909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope)); 13919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetUseCoarseDM_C",PCTelescopeGetUseCoarseDM_Telescope)); 13929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetUseCoarseDM_C",PCTelescopeSetUseCoarseDM_Telescope)); 13931e07b27eSBarry Smith PetscFunctionReturn(0); 13941e07b27eSBarry Smith } 1395