11e07b27eSBarry Smith 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" 20bf00f589SPatrick Sanan " url = {http://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 /* 301e07b27eSBarry Smith PCTelescopeSetUp_default() 311e07b27eSBarry Smith PCTelescopeMatCreate_default() 321e07b27eSBarry Smith 331e07b27eSBarry Smith default 341e07b27eSBarry Smith 351e07b27eSBarry Smith // scatter in 361e07b27eSBarry Smith x(comm) -> xtmp(comm) 371e07b27eSBarry Smith 381e07b27eSBarry Smith xred(subcomm) <- xtmp 391e07b27eSBarry Smith yred(subcomm) 401e07b27eSBarry Smith 411e07b27eSBarry Smith yred(subcomm) --> xtmp 421e07b27eSBarry Smith 431e07b27eSBarry Smith // scatter out 441e07b27eSBarry Smith xtmp(comm) -> y(comm) 451e07b27eSBarry Smith */ 461e07b27eSBarry Smith 471e07b27eSBarry Smith PetscBool isActiveRank(PetscSubcomm scomm) 481e07b27eSBarry Smith { 491e07b27eSBarry Smith if (scomm->color == 0) { return PETSC_TRUE; } 501e07b27eSBarry Smith else { return PETSC_FALSE; } 511e07b27eSBarry Smith } 521e07b27eSBarry Smith 531e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred) 541e07b27eSBarry Smith { 55c6a0d831SBarry Smith DM subdm = NULL; 561e07b27eSBarry Smith 571e07b27eSBarry Smith if (!isActiveRank(sred->psubcomm)) { subdm = NULL; } 581e07b27eSBarry Smith else { 591e07b27eSBarry Smith switch (sred->sr_type) { 601e07b27eSBarry Smith case TELESCOPE_DEFAULT: subdm = NULL; 611e07b27eSBarry Smith break; 621e07b27eSBarry Smith case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 631e07b27eSBarry Smith break; 641e07b27eSBarry Smith case TELESCOPE_DMPLEX: subdm = NULL; 651e07b27eSBarry Smith break; 661e07b27eSBarry Smith } 671e07b27eSBarry Smith } 681e07b27eSBarry Smith return(subdm); 691e07b27eSBarry Smith } 701e07b27eSBarry Smith 711e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 721e07b27eSBarry Smith { 731e07b27eSBarry Smith PetscErrorCode ierr; 741e07b27eSBarry Smith PetscInt m,M,bs,st,ed; 751e07b27eSBarry Smith Vec x,xred,yred,xtmp; 761e07b27eSBarry Smith Mat B; 771e07b27eSBarry Smith MPI_Comm comm,subcomm; 781e07b27eSBarry Smith VecScatter scatter; 791e07b27eSBarry Smith IS isin; 801e07b27eSBarry Smith 811e07b27eSBarry Smith PetscFunctionBegin; 821e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr); 831e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 841e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 851e07b27eSBarry Smith 861e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 871e07b27eSBarry Smith ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 881e07b27eSBarry Smith ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 891e07b27eSBarry Smith ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 901e07b27eSBarry Smith 911e07b27eSBarry Smith xred = NULL; 923ac26c5eSBarry Smith m = 0; 931e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 941e07b27eSBarry Smith ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr); 951e07b27eSBarry Smith ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr); 961e07b27eSBarry Smith ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr); 971e07b27eSBarry Smith ierr = VecSetFromOptions(xred);CHKERRQ(ierr); 98ca43db0aSBarry Smith ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr); 991e07b27eSBarry Smith } 1001e07b27eSBarry Smith 1011e07b27eSBarry Smith yred = NULL; 1021e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 1031e07b27eSBarry Smith ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 1041e07b27eSBarry Smith } 1051e07b27eSBarry Smith 1061e07b27eSBarry Smith ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 1071e07b27eSBarry Smith ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 1081e07b27eSBarry Smith ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 1091e07b27eSBarry Smith ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 1101e07b27eSBarry Smith 1111e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 1121e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 1131e07b27eSBarry Smith ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr); 1141e07b27eSBarry Smith } else { 1151e07b27eSBarry Smith ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 1163ac26c5eSBarry Smith ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 1171e07b27eSBarry Smith } 1181e07b27eSBarry Smith ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 1191e07b27eSBarry Smith 12035928de7SBarry Smith ierr = VecScatterCreateWithData(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 1211e07b27eSBarry Smith 1221e07b27eSBarry Smith sred->isin = isin; 1231e07b27eSBarry Smith sred->scatter = scatter; 1241e07b27eSBarry Smith sred->xred = xred; 1251e07b27eSBarry Smith sred->yred = yred; 1261e07b27eSBarry Smith sred->xtmp = xtmp; 1271e07b27eSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1281e07b27eSBarry Smith PetscFunctionReturn(0); 1291e07b27eSBarry Smith } 1301e07b27eSBarry Smith 1311e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 1321e07b27eSBarry Smith { 1331e07b27eSBarry Smith PetscErrorCode ierr; 1341e07b27eSBarry Smith MPI_Comm comm,subcomm; 1351e07b27eSBarry Smith Mat Bred,B; 1361e07b27eSBarry Smith PetscInt nr,nc; 1371e07b27eSBarry Smith IS isrow,iscol; 1381e07b27eSBarry Smith Mat Blocal,*_Blocal; 1391e07b27eSBarry Smith 1401e07b27eSBarry Smith PetscFunctionBegin; 1411e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr); 1421e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1431e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1441e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 1451e07b27eSBarry Smith ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 1461e07b27eSBarry Smith isrow = sred->isin; 1471e07b27eSBarry Smith ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 1487dae84e0SHong Zhang ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 1491e07b27eSBarry Smith Blocal = *_Blocal; 1501e07b27eSBarry Smith ierr = PetscFree(_Blocal);CHKERRQ(ierr); 1511e07b27eSBarry Smith Bred = NULL; 1521e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 1531e07b27eSBarry Smith PetscInt mm; 1541e07b27eSBarry Smith 1551e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 1561e07b27eSBarry Smith 1571e07b27eSBarry Smith ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 1581e07b27eSBarry Smith ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 1591e07b27eSBarry Smith } 1601e07b27eSBarry Smith *A = Bred; 1611e07b27eSBarry Smith ierr = ISDestroy(&iscol);CHKERRQ(ierr); 1621e07b27eSBarry Smith ierr = MatDestroy(&Blocal);CHKERRQ(ierr); 1631e07b27eSBarry Smith PetscFunctionReturn(0); 1641e07b27eSBarry Smith } 1651e07b27eSBarry Smith 166392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 1671e07b27eSBarry Smith { 1681e07b27eSBarry Smith PetscErrorCode ierr; 1691e07b27eSBarry Smith PetscBool has_const; 1701e07b27eSBarry Smith const Vec *vecs; 171c41e779fSDave May Vec *sub_vecs = NULL; 172392968a1SPatrick Sanan PetscInt i,k,n = 0; 1731e07b27eSBarry Smith MPI_Comm subcomm; 1741e07b27eSBarry Smith 1751e07b27eSBarry Smith PetscFunctionBegin; 1761e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1771e07b27eSBarry Smith ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 1781e07b27eSBarry Smith 1791e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 180e3acf2f7SBarry Smith if (n) { 181e3acf2f7SBarry Smith ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 1821e07b27eSBarry Smith } 1831e07b27eSBarry Smith } 1841e07b27eSBarry Smith 1851e07b27eSBarry Smith /* copy entries */ 1861e07b27eSBarry Smith for (k=0; k<n; k++) { 1871e07b27eSBarry Smith const PetscScalar *x_array; 1881e07b27eSBarry Smith PetscScalar *LA_sub_vec; 18913c30530SDave May PetscInt st,ed; 1901e07b27eSBarry Smith 1911e07b27eSBarry Smith /* pull in vector x->xtmp */ 1921e07b27eSBarry Smith ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1931e07b27eSBarry Smith ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 19447856c66SBarry Smith if (sub_vecs) { 195a04a6428SPatrick Sanan /* copy vector entries into xred */ 1961e07b27eSBarry Smith ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 197ea2b237eSDave May if (sub_vecs[k]) { 1981e07b27eSBarry Smith ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 1991e07b27eSBarry Smith ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 2001e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 2011e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 2021e07b27eSBarry Smith } 2031e07b27eSBarry Smith ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 2041e07b27eSBarry Smith } 2051e07b27eSBarry Smith ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 2061e07b27eSBarry Smith } 20747856c66SBarry Smith } 2081e07b27eSBarry Smith 2091e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 210d8b9d5b7SPatrick Sanan /* create new (near) nullspace for redundant object */ 211392968a1SPatrick Sanan ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr); 212392968a1SPatrick Sanan ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 213d8b9d5b7SPatrick Sanan if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 214d8b9d5b7SPatrick Sanan if (nullspace->rmctx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope"); 215d8b9d5b7SPatrick Sanan } 216392968a1SPatrick Sanan PetscFunctionReturn(0); 217392968a1SPatrick Sanan } 218392968a1SPatrick Sanan 219392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 220392968a1SPatrick Sanan { 221392968a1SPatrick Sanan PetscErrorCode ierr; 222392968a1SPatrick Sanan Mat B; 223392968a1SPatrick Sanan 224392968a1SPatrick Sanan PetscFunctionBegin; 225392968a1SPatrick Sanan ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 226392968a1SPatrick Sanan 227392968a1SPatrick Sanan /* Propagate the nullspace if it exists */ 228392968a1SPatrick Sanan { 229392968a1SPatrick Sanan MatNullSpace nullspace,sub_nullspace; 230392968a1SPatrick Sanan ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 231392968a1SPatrick Sanan if (nullspace) { 232392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr); 233392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr); 234392968a1SPatrick Sanan if (isActiveRank(sred->psubcomm)) { 235392968a1SPatrick Sanan ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 23641ff1ee9SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr); 2371e07b27eSBarry Smith } 238392968a1SPatrick Sanan } 239392968a1SPatrick Sanan } 240392968a1SPatrick Sanan 241392968a1SPatrick Sanan /* Propagate the near nullspace if it exists */ 242392968a1SPatrick Sanan { 243392968a1SPatrick Sanan MatNullSpace nearnullspace,sub_nearnullspace; 244392968a1SPatrick Sanan ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr); 245392968a1SPatrick Sanan if (nearnullspace) { 246392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr); 247392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr); 248392968a1SPatrick Sanan if (isActiveRank(sred->psubcomm)) { 249392968a1SPatrick Sanan ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr); 250392968a1SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr); 251392968a1SPatrick Sanan } 252392968a1SPatrick Sanan } 253392968a1SPatrick Sanan } 2541e07b27eSBarry Smith PetscFunctionReturn(0); 2551e07b27eSBarry Smith } 2561e07b27eSBarry Smith 2571e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 2581e07b27eSBarry Smith { 2591e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 2601e07b27eSBarry Smith PetscErrorCode ierr; 2611e07b27eSBarry Smith PetscBool iascii,isstring; 2621e07b27eSBarry Smith PetscViewer subviewer; 2631e07b27eSBarry Smith 2641e07b27eSBarry Smith PetscFunctionBegin; 2651e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2661e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 2671e07b27eSBarry Smith if (iascii) { 2681e07b27eSBarry Smith if (!sred->psubcomm) { 269efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," preconditioner not yet setup\n");CHKERRQ(ierr); 2701e07b27eSBarry Smith } else { 2711e07b27eSBarry Smith MPI_Comm comm,subcomm; 2721e07b27eSBarry Smith PetscMPIInt comm_size,subcomm_size; 2731e07b27eSBarry Smith DM dm,subdm; 2741e07b27eSBarry Smith 2751e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 2761e07b27eSBarry Smith subdm = private_PCTelescopeGetSubDM(sred); 2771e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 2781e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2791e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 2801e07b27eSBarry Smith ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 2811e07b27eSBarry Smith 282efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 283efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 28448a10b22SPatrick Sanan switch (sred->subcommtype) { 28548a10b22SPatrick Sanan case PETSC_SUBCOMM_INTERLACED : 286efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," subcomm type: interlaced\n",sred->subcommtype);CHKERRQ(ierr); 28748a10b22SPatrick Sanan break; 28848a10b22SPatrick Sanan case PETSC_SUBCOMM_CONTIGUOUS : 289efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," subcomm type: contiguous\n",sred->subcommtype);CHKERRQ(ierr); 29048a10b22SPatrick Sanan break; 29148a10b22SPatrick Sanan default : 29248a10b22SPatrick Sanan SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope"); 29348a10b22SPatrick Sanan } 2941e07b27eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 2951e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 2961e07b27eSBarry Smith ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 2971e07b27eSBarry Smith 2981e07b27eSBarry Smith if (dm && sred->ignore_dm) { 299efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," ignoring DM\n");CHKERRQ(ierr); 3001e07b27eSBarry Smith } 3017c5279cbSDave May if (sred->ignore_kspcomputeoperators) { 302efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," ignoring KSPComputeOperators\n");CHKERRQ(ierr); 3037c5279cbSDave May } 3041e07b27eSBarry Smith switch (sred->sr_type) { 3051e07b27eSBarry Smith case TELESCOPE_DEFAULT: 306efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," using default setup\n");CHKERRQ(ierr); 3071e07b27eSBarry Smith break; 3081e07b27eSBarry Smith case TELESCOPE_DMDA: 309efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," DMDA detected\n");CHKERRQ(ierr); 3108ef9ca65SPatrick Sanan ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr); 3111e07b27eSBarry Smith break; 3121e07b27eSBarry Smith case TELESCOPE_DMPLEX: 313efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," DMPLEX detected\n");CHKERRQ(ierr); 3141e07b27eSBarry Smith break; 3151e07b27eSBarry Smith } 3161e07b27eSBarry Smith 3171e07b27eSBarry Smith ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 3181e07b27eSBarry Smith ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 3191e07b27eSBarry Smith } 3201e07b27eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 3211e07b27eSBarry Smith } 3221e07b27eSBarry Smith } 3231e07b27eSBarry Smith PetscFunctionReturn(0); 3241e07b27eSBarry Smith } 3251e07b27eSBarry Smith 3261e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc) 3271e07b27eSBarry Smith { 3281e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 3291e07b27eSBarry Smith PetscErrorCode ierr; 330bd49479cSSatish Balay MPI_Comm comm,subcomm=0; 3311e07b27eSBarry Smith PCTelescopeType sr_type; 3321e07b27eSBarry Smith 3331e07b27eSBarry Smith PetscFunctionBegin; 3341e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 3351e07b27eSBarry Smith 3361e07b27eSBarry Smith /* subcomm definition */ 3371e07b27eSBarry Smith if (!pc->setupcalled) { 3381e07b27eSBarry Smith if (!sred->psubcomm) { 3391e07b27eSBarry Smith ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 3401e07b27eSBarry Smith ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 34148a10b22SPatrick Sanan ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr); 3421e07b27eSBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3431e07b27eSBarry Smith } 3441e07b27eSBarry Smith } 3450f6f40a7SSatish Balay subcomm = PetscSubcommChild(sred->psubcomm); 3461e07b27eSBarry Smith 3471e07b27eSBarry Smith /* internal KSP */ 3481e07b27eSBarry Smith if (!pc->setupcalled) { 3491e07b27eSBarry Smith const char *prefix; 3501e07b27eSBarry Smith 3511e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 3521e07b27eSBarry Smith ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 3531e07b27eSBarry Smith ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 3541e07b27eSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3551e07b27eSBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 3561e07b27eSBarry Smith ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 3571e07b27eSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3581e07b27eSBarry Smith ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 3591e07b27eSBarry Smith ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 3601e07b27eSBarry Smith } 3611e07b27eSBarry Smith } 3621e07b27eSBarry Smith /* Determine type of setup/update */ 3631e07b27eSBarry Smith if (!pc->setupcalled) { 3641e07b27eSBarry Smith PetscBool has_dm,same; 3651e07b27eSBarry Smith DM dm; 3661e07b27eSBarry Smith 3671e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 3681e07b27eSBarry Smith has_dm = PETSC_FALSE; 3691e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3701e07b27eSBarry Smith if (dm) { has_dm = PETSC_TRUE; } 3711e07b27eSBarry Smith if (has_dm) { 3721e07b27eSBarry Smith /* check for dmda */ 3731e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 3741e07b27eSBarry Smith if (same) { 3751e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 3761e07b27eSBarry Smith sr_type = TELESCOPE_DMDA; 3771e07b27eSBarry Smith } 3781e07b27eSBarry Smith /* check for dmplex */ 3791e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 3801e07b27eSBarry Smith if (same) { 381994fe344SLisandro Dalcin ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr); 3821e07b27eSBarry Smith sr_type = TELESCOPE_DMPLEX; 3831e07b27eSBarry Smith } 3841e07b27eSBarry Smith } 3851e07b27eSBarry Smith 3861e07b27eSBarry Smith if (sred->ignore_dm) { 387994fe344SLisandro Dalcin ierr = PetscInfo(pc,"PCTelescope: ignore DM\n");CHKERRQ(ierr); 3881e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 3891e07b27eSBarry Smith } 3901e07b27eSBarry Smith sred->sr_type = sr_type; 3911e07b27eSBarry Smith } else { 3921e07b27eSBarry Smith sr_type = sred->sr_type; 3931e07b27eSBarry Smith } 3941e07b27eSBarry Smith 395d8b9d5b7SPatrick Sanan /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */ 3961e07b27eSBarry Smith switch (sr_type) { 3971e07b27eSBarry Smith case TELESCOPE_DEFAULT: 3981e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 3991e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 4001e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 4011e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 4021e07b27eSBarry Smith break; 4031e07b27eSBarry Smith case TELESCOPE_DMDA: 4041e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope_dmda; 405f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_dmda; 4061e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 4071e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 4081e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 4091e07b27eSBarry Smith sred->pctelescope_reset_type = PCReset_Telescope_dmda; 4101e07b27eSBarry Smith break; 411d8b9d5b7SPatrick Sanan case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available"); 4121e07b27eSBarry Smith break; 413a04a6428SPatrick Sanan default: SETERRQ(comm,PETSC_ERR_SUP,"Only support for repartitioning DMDA is provided"); 4141e07b27eSBarry Smith break; 4151e07b27eSBarry Smith } 4161e07b27eSBarry Smith 4171e07b27eSBarry Smith /* setup */ 4181e07b27eSBarry Smith if (sred->pctelescope_setup_type) { 4191e07b27eSBarry Smith ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 4201e07b27eSBarry Smith } 4211e07b27eSBarry Smith /* update */ 4221e07b27eSBarry Smith if (!pc->setupcalled) { 4231e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 4241e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 4251e07b27eSBarry Smith } 4261e07b27eSBarry Smith if (sred->pctelescope_matnullspacecreate_type) { 427392968a1SPatrick Sanan ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 4281e07b27eSBarry Smith } 4291e07b27eSBarry Smith } else { 4301e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 4311e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 4321e07b27eSBarry Smith } 4331e07b27eSBarry Smith } 4341e07b27eSBarry Smith 4351e07b27eSBarry Smith /* common - no construction */ 4361e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 4371e07b27eSBarry Smith ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 4381e07b27eSBarry Smith if (pc->setfromoptionscalled && !pc->setupcalled){ 4391e07b27eSBarry Smith ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 4401e07b27eSBarry Smith } 4411e07b27eSBarry Smith } 4421e07b27eSBarry Smith PetscFunctionReturn(0); 4431e07b27eSBarry Smith } 4441e07b27eSBarry Smith 4451e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 4461e07b27eSBarry Smith { 4471e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 4481e07b27eSBarry Smith PetscErrorCode ierr; 4491e07b27eSBarry Smith Vec xtmp,xred,yred; 45013c30530SDave May PetscInt i,st,ed; 4511e07b27eSBarry Smith VecScatter scatter; 4521e07b27eSBarry Smith PetscScalar *array; 4531e07b27eSBarry Smith const PetscScalar *x_array; 4541e07b27eSBarry Smith 4551e07b27eSBarry Smith PetscFunctionBegin; 456bf00f589SPatrick Sanan ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 457bf00f589SPatrick Sanan 4581e07b27eSBarry Smith xtmp = sred->xtmp; 4591e07b27eSBarry Smith scatter = sred->scatter; 4601e07b27eSBarry Smith xred = sred->xred; 4611e07b27eSBarry Smith yred = sred->yred; 4621e07b27eSBarry Smith 4631e07b27eSBarry Smith /* pull in vector x->xtmp */ 4641e07b27eSBarry Smith ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4651e07b27eSBarry Smith ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4661e07b27eSBarry Smith 467bf00f589SPatrick Sanan /* copy vector entries into xred */ 4681e07b27eSBarry Smith ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 4691e07b27eSBarry Smith if (xred) { 4701e07b27eSBarry Smith PetscScalar *LA_xred; 4711e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 4721e07b27eSBarry Smith ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 4731e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 4741e07b27eSBarry Smith LA_xred[i] = x_array[i]; 4751e07b27eSBarry Smith } 4761e07b27eSBarry Smith ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 4771e07b27eSBarry Smith } 4781e07b27eSBarry Smith ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 4791e07b27eSBarry Smith /* solve */ 4801e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 4811e07b27eSBarry Smith ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 482*c0decd05SBarry Smith ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr); 4831e07b27eSBarry Smith } 4841e07b27eSBarry Smith /* return vector */ 4851e07b27eSBarry Smith ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 4861e07b27eSBarry Smith if (yred) { 4871e07b27eSBarry Smith const PetscScalar *LA_yred; 4881e07b27eSBarry Smith ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 4891e07b27eSBarry Smith ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 4901e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 4911e07b27eSBarry Smith array[i] = LA_yred[i]; 4921e07b27eSBarry Smith } 4931e07b27eSBarry Smith ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 4941e07b27eSBarry Smith } 4951e07b27eSBarry Smith ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 4961e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4971e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4981e07b27eSBarry Smith PetscFunctionReturn(0); 4991e07b27eSBarry Smith } 5001e07b27eSBarry Smith 501f650675bSDave 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) 502f650675bSDave May { 503f650675bSDave May PC_Telescope sred = (PC_Telescope)pc->data; 504f650675bSDave May PetscErrorCode ierr; 505a1d91a28SDave May Vec xtmp,yred; 506f650675bSDave May PetscInt i,st,ed; 507f650675bSDave May VecScatter scatter; 508f650675bSDave May const PetscScalar *x_array; 509f650675bSDave May PetscBool default_init_guess_value; 510f650675bSDave May 511f650675bSDave May PetscFunctionBegin; 512f650675bSDave May xtmp = sred->xtmp; 513f650675bSDave May scatter = sred->scatter; 514f650675bSDave May yred = sred->yred; 515f650675bSDave May 516f650675bSDave May if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1"); 517f650675bSDave May *reason = (PCRichardsonConvergedReason)0; 518f650675bSDave May 519f650675bSDave May if (!zeroguess) { 520f650675bSDave May ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr); 521f650675bSDave May /* pull in vector y->xtmp */ 522f650675bSDave May ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 523f650675bSDave May ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 524f650675bSDave May 525bf00f589SPatrick Sanan /* copy vector entries into xred */ 526f650675bSDave May ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 527f650675bSDave May if (yred) { 528f650675bSDave May PetscScalar *LA_yred; 529f650675bSDave May ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 530f650675bSDave May ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 531f650675bSDave May for (i=0; i<ed-st; i++) { 532f650675bSDave May LA_yred[i] = x_array[i]; 533f650675bSDave May } 534f650675bSDave May ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 535f650675bSDave May } 536f650675bSDave May ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 537f650675bSDave May } 538f650675bSDave May 539f650675bSDave May if (isActiveRank(sred->psubcomm)) { 540f650675bSDave May ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr); 541f650675bSDave May if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 542f650675bSDave May } 543f650675bSDave May 544f650675bSDave May ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr); 545f650675bSDave May 546f650675bSDave May if (isActiveRank(sred->psubcomm)) { 547f650675bSDave May ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr); 548f650675bSDave May } 549f650675bSDave May 550f650675bSDave May if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 551f650675bSDave May *outits = 1; 552f650675bSDave May PetscFunctionReturn(0); 553f650675bSDave May } 554f650675bSDave May 5551e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc) 5561e07b27eSBarry Smith { 5571e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5581e07b27eSBarry Smith PetscErrorCode ierr; 5591e07b27eSBarry Smith 5601e07b27eSBarry Smith ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 5611e07b27eSBarry Smith ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 562e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); 563e3acf2f7SBarry Smith ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); 564e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); 565e3acf2f7SBarry Smith ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); 566e3acf2f7SBarry Smith ierr = KSPReset(sred->ksp);CHKERRQ(ierr); 5671e07b27eSBarry Smith if (sred->pctelescope_reset_type) { 5681e07b27eSBarry Smith ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 5691e07b27eSBarry Smith } 5701e07b27eSBarry Smith PetscFunctionReturn(0); 5711e07b27eSBarry Smith } 5721e07b27eSBarry Smith 5731e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc) 5741e07b27eSBarry Smith { 5751e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5761e07b27eSBarry Smith PetscErrorCode ierr; 5771e07b27eSBarry Smith 5781e07b27eSBarry Smith PetscFunctionBegin; 5791e07b27eSBarry Smith ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 580e3acf2f7SBarry Smith ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); 5811e07b27eSBarry Smith ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 582e3acf2f7SBarry Smith ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr); 583e3acf2f7SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 5841e07b27eSBarry Smith PetscFunctionReturn(0); 5851e07b27eSBarry Smith } 5861e07b27eSBarry Smith 5874416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc) 5881e07b27eSBarry Smith { 5891e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5901e07b27eSBarry Smith PetscErrorCode ierr; 5911e07b27eSBarry Smith MPI_Comm comm; 5921e07b27eSBarry Smith PetscMPIInt size; 59348a10b22SPatrick Sanan PetscBool flg; 59448a10b22SPatrick Sanan PetscSubcommType subcommtype; 5951e07b27eSBarry Smith 5961e07b27eSBarry Smith PetscFunctionBegin; 5971e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5981e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 5991e07b27eSBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 60048a10b22SPatrick Sanan ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr); 60148a10b22SPatrick Sanan if (flg) { 60248a10b22SPatrick Sanan ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr); 60348a10b22SPatrick Sanan } 6041e07b27eSBarry Smith ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 6051e07b27eSBarry Smith if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 6061e07b27eSBarry Smith ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 6077c5279cbSDave May ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr); 6081e07b27eSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6091e07b27eSBarry Smith PetscFunctionReturn(0); 6101e07b27eSBarry Smith } 6111e07b27eSBarry Smith 6121e07b27eSBarry Smith /* PC simplementation specific API's */ 6131e07b27eSBarry Smith 6141e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 6151e07b27eSBarry Smith { 6161e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 617bd49479cSSatish Balay PetscFunctionBegin; 6181e07b27eSBarry Smith if (ksp) *ksp = red->ksp; 619bd49479cSSatish Balay PetscFunctionReturn(0); 6201e07b27eSBarry Smith } 6211e07b27eSBarry Smith 62248a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype) 62348a10b22SPatrick Sanan { 62448a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 62548a10b22SPatrick Sanan PetscFunctionBegin; 62648a10b22SPatrick Sanan if (subcommtype) *subcommtype = red->subcommtype; 62748a10b22SPatrick Sanan PetscFunctionReturn(0); 62848a10b22SPatrick Sanan } 62948a10b22SPatrick Sanan 63048a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype) 63148a10b22SPatrick Sanan { 63248a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 63348a10b22SPatrick Sanan 63448a10b22SPatrick Sanan PetscFunctionBegin; 63548a10b22SPatrick Sanan if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"You cannot change the subcommunicator type for PCTelescope after it has been set up."); 63648a10b22SPatrick Sanan red->subcommtype = subcommtype; 63748a10b22SPatrick Sanan PetscFunctionReturn(0); 63848a10b22SPatrick Sanan } 63948a10b22SPatrick Sanan 6401e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 6411e07b27eSBarry Smith { 6421e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 643bd49479cSSatish Balay PetscFunctionBegin; 6441e07b27eSBarry Smith if (fact) *fact = red->redfactor; 645bd49479cSSatish Balay PetscFunctionReturn(0); 6461e07b27eSBarry Smith } 6471e07b27eSBarry Smith 6481e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 6491e07b27eSBarry Smith { 6501e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 6511e07b27eSBarry Smith PetscMPIInt size; 6521e07b27eSBarry Smith PetscErrorCode ierr; 6531e07b27eSBarry Smith 654bd49479cSSatish Balay PetscFunctionBegin; 6551e07b27eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 6561e07b27eSBarry Smith if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 6571e07b27eSBarry Smith if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 6581e07b27eSBarry Smith red->redfactor = fact; 659bd49479cSSatish Balay PetscFunctionReturn(0); 6601e07b27eSBarry Smith } 6611e07b27eSBarry Smith 6621e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 6631e07b27eSBarry Smith { 6641e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 665bd49479cSSatish Balay PetscFunctionBegin; 6661e07b27eSBarry Smith if (v) *v = red->ignore_dm; 667bd49479cSSatish Balay PetscFunctionReturn(0); 6681e07b27eSBarry Smith } 66948a10b22SPatrick Sanan 6701e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 6711e07b27eSBarry Smith { 6721e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 673bd49479cSSatish Balay PetscFunctionBegin; 6741e07b27eSBarry Smith red->ignore_dm = v; 675bd49479cSSatish Balay PetscFunctionReturn(0); 6761e07b27eSBarry Smith } 6771e07b27eSBarry Smith 6780ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v) 6790ae7c45bSDave May { 6800ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 6810ae7c45bSDave May PetscFunctionBegin; 6820ae7c45bSDave May if (v) *v = red->ignore_kspcomputeoperators; 6830ae7c45bSDave May PetscFunctionReturn(0); 6840ae7c45bSDave May } 68548a10b22SPatrick Sanan 6860ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v) 6870ae7c45bSDave May { 6880ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 6890ae7c45bSDave May PetscFunctionBegin; 6900ae7c45bSDave May red->ignore_kspcomputeoperators = v; 6910ae7c45bSDave May PetscFunctionReturn(0); 6920ae7c45bSDave May } 6930ae7c45bSDave May 6941e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 6951e07b27eSBarry Smith { 6961e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 697bd49479cSSatish Balay PetscFunctionBegin; 6981e07b27eSBarry Smith *dm = private_PCTelescopeGetSubDM(red); 699bd49479cSSatish Balay PetscFunctionReturn(0); 7001e07b27eSBarry Smith } 7011e07b27eSBarry Smith 7021e07b27eSBarry Smith /*@ 7031e07b27eSBarry Smith PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 7041e07b27eSBarry Smith 7051e07b27eSBarry Smith Not Collective 7061e07b27eSBarry Smith 7071e07b27eSBarry Smith Input Parameter: 7081e07b27eSBarry Smith . pc - the preconditioner context 7091e07b27eSBarry Smith 7101e07b27eSBarry Smith Output Parameter: 7111e07b27eSBarry Smith . subksp - the KSP defined the smaller set of processes 7121e07b27eSBarry Smith 7131e07b27eSBarry Smith Level: advanced 7141e07b27eSBarry Smith 7151e07b27eSBarry Smith .keywords: PC, telescopting solve 7161e07b27eSBarry Smith @*/ 7171e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 7181e07b27eSBarry Smith { 719bd49479cSSatish Balay PetscErrorCode ierr; 720bd49479cSSatish Balay PetscFunctionBegin; 721163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 722bd49479cSSatish Balay PetscFunctionReturn(0); 7231e07b27eSBarry Smith } 7241e07b27eSBarry Smith 7251e07b27eSBarry Smith /*@ 7261e07b27eSBarry Smith PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 7271e07b27eSBarry Smith 7281e07b27eSBarry Smith Not Collective 7291e07b27eSBarry Smith 7301e07b27eSBarry Smith Input Parameter: 7311e07b27eSBarry Smith . pc - the preconditioner context 7321e07b27eSBarry Smith 7331e07b27eSBarry Smith Output Parameter: 7341e07b27eSBarry Smith . fact - the reduction factor 7351e07b27eSBarry Smith 7361e07b27eSBarry Smith Level: advanced 7371e07b27eSBarry Smith 7381e07b27eSBarry Smith .keywords: PC, telescoping solve 7391e07b27eSBarry Smith @*/ 7401e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 7411e07b27eSBarry Smith { 742bd49479cSSatish Balay PetscErrorCode ierr; 743bd49479cSSatish Balay PetscFunctionBegin; 744163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 745bd49479cSSatish Balay PetscFunctionReturn(0); 7461e07b27eSBarry Smith } 7471e07b27eSBarry Smith 7481e07b27eSBarry Smith /*@ 7491e07b27eSBarry Smith PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 7501e07b27eSBarry Smith 7511e07b27eSBarry Smith Not Collective 7521e07b27eSBarry Smith 7531e07b27eSBarry Smith Input Parameter: 7541e07b27eSBarry Smith . pc - the preconditioner context 7551e07b27eSBarry Smith 7561e07b27eSBarry Smith Output Parameter: 7571e07b27eSBarry Smith . fact - the reduction factor 7581e07b27eSBarry Smith 7591e07b27eSBarry Smith Level: advanced 7601e07b27eSBarry Smith 7611e07b27eSBarry Smith .keywords: PC, telescoping solve 7621e07b27eSBarry Smith @*/ 7631e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 7641e07b27eSBarry Smith { 765bd49479cSSatish Balay PetscErrorCode ierr; 766bd49479cSSatish Balay PetscFunctionBegin; 767bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 768bd49479cSSatish Balay PetscFunctionReturn(0); 7691e07b27eSBarry Smith } 7701e07b27eSBarry Smith 7711e07b27eSBarry Smith /*@ 7721e07b27eSBarry Smith PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 7731e07b27eSBarry Smith 7741e07b27eSBarry Smith Not Collective 7751e07b27eSBarry Smith 7761e07b27eSBarry Smith Input Parameter: 7771e07b27eSBarry Smith . pc - the preconditioner context 7781e07b27eSBarry Smith 7791e07b27eSBarry Smith Output Parameter: 7801e07b27eSBarry Smith . v - the flag 7811e07b27eSBarry Smith 7821e07b27eSBarry Smith Level: advanced 7831e07b27eSBarry Smith 7841e07b27eSBarry Smith .keywords: PC, telescoping solve 7851e07b27eSBarry Smith @*/ 7861e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 7871e07b27eSBarry Smith { 788bd49479cSSatish Balay PetscErrorCode ierr; 789bd49479cSSatish Balay PetscFunctionBegin; 790163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 791bd49479cSSatish Balay PetscFunctionReturn(0); 7921e07b27eSBarry Smith } 7931e07b27eSBarry Smith 7941e07b27eSBarry Smith /*@ 7951e07b27eSBarry Smith PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 7961e07b27eSBarry Smith 7971e07b27eSBarry Smith Not Collective 7981e07b27eSBarry Smith 7991e07b27eSBarry Smith Input Parameter: 8001e07b27eSBarry Smith . pc - the preconditioner context 8011e07b27eSBarry Smith 8021e07b27eSBarry Smith Output Parameter: 8031e07b27eSBarry Smith . v - Use PETSC_TRUE to ignore any DM 8041e07b27eSBarry Smith 8051e07b27eSBarry Smith Level: advanced 8061e07b27eSBarry Smith 8071e07b27eSBarry Smith .keywords: PC, telescoping solve 8081e07b27eSBarry Smith @*/ 809bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 8101e07b27eSBarry Smith { 811bd49479cSSatish Balay PetscErrorCode ierr; 812bd49479cSSatish Balay PetscFunctionBegin; 813bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 814bd49479cSSatish Balay PetscFunctionReturn(0); 8151e07b27eSBarry Smith } 8161e07b27eSBarry Smith 8171e07b27eSBarry Smith /*@ 8180ae7c45bSDave May PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used. 8190ae7c45bSDave May 8200ae7c45bSDave May Not Collective 8210ae7c45bSDave May 8220ae7c45bSDave May Input Parameter: 8230ae7c45bSDave May . pc - the preconditioner context 8240ae7c45bSDave May 8250ae7c45bSDave May Output Parameter: 8260ae7c45bSDave May . v - the flag 8270ae7c45bSDave May 8280ae7c45bSDave May Level: advanced 8290ae7c45bSDave May 8300ae7c45bSDave May .keywords: PC, telescoping solve 8310ae7c45bSDave May @*/ 8320ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v) 8330ae7c45bSDave May { 8340ae7c45bSDave May PetscErrorCode ierr; 8350ae7c45bSDave May PetscFunctionBegin; 836163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 8370ae7c45bSDave May PetscFunctionReturn(0); 8380ae7c45bSDave May } 8390ae7c45bSDave May 8400ae7c45bSDave May /*@ 8410ae7c45bSDave May PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators. 8420ae7c45bSDave May 8430ae7c45bSDave May Not Collective 8440ae7c45bSDave May 8450ae7c45bSDave May Input Parameter: 8460ae7c45bSDave May . pc - the preconditioner context 8470ae7c45bSDave May 8480ae7c45bSDave May Output Parameter: 849a954d8f4SDave May . v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc 8500ae7c45bSDave May 8510ae7c45bSDave May Level: advanced 8520ae7c45bSDave May 8530ae7c45bSDave May .keywords: PC, telescoping solve 8540ae7c45bSDave May @*/ 8550ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v) 8560ae7c45bSDave May { 8570ae7c45bSDave May PetscErrorCode ierr; 8580ae7c45bSDave May PetscFunctionBegin; 8590ae7c45bSDave May ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 8600ae7c45bSDave May PetscFunctionReturn(0); 8610ae7c45bSDave May } 8620ae7c45bSDave May 8630ae7c45bSDave May /*@ 8641e07b27eSBarry Smith PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 8651e07b27eSBarry Smith 8661e07b27eSBarry Smith Not Collective 8671e07b27eSBarry Smith 8681e07b27eSBarry Smith Input Parameter: 8691e07b27eSBarry Smith . pc - the preconditioner context 8701e07b27eSBarry Smith 8711e07b27eSBarry Smith Output Parameter: 8721e07b27eSBarry Smith . subdm - The re-partitioned DM 8731e07b27eSBarry Smith 8741e07b27eSBarry Smith Level: advanced 8751e07b27eSBarry Smith 8761e07b27eSBarry Smith .keywords: PC, telescoping solve 8771e07b27eSBarry Smith @*/ 8781e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 8791e07b27eSBarry Smith { 880bd49479cSSatish Balay PetscErrorCode ierr; 881bd49479cSSatish Balay PetscFunctionBegin; 882163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 883bd49479cSSatish Balay PetscFunctionReturn(0); 8841e07b27eSBarry Smith } 8851e07b27eSBarry Smith 88648a10b22SPatrick Sanan /*@ 88748a10b22SPatrick Sanan PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous) 88848a10b22SPatrick Sanan 88948a10b22SPatrick Sanan Logically Collective 89048a10b22SPatrick Sanan 89148a10b22SPatrick Sanan Input Parameter: 8921dae98e4SBarry Smith + pc - the preconditioner context 8931dae98e4SBarry Smith - subcommtype - the subcommunicator type (see PetscSubcommType) 89448a10b22SPatrick Sanan 89548a10b22SPatrick Sanan Level: advanced 89648a10b22SPatrick Sanan 89748a10b22SPatrick Sanan .keywords: PC, telescoping solve 89848a10b22SPatrick Sanan 89948a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE 90048a10b22SPatrick Sanan @*/ 90148a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) 90248a10b22SPatrick Sanan { 90348a10b22SPatrick Sanan PetscErrorCode ierr; 90448a10b22SPatrick Sanan PetscFunctionBegin; 90548a10b22SPatrick Sanan ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr); 90648a10b22SPatrick Sanan PetscFunctionReturn(0); 90748a10b22SPatrick Sanan } 90848a10b22SPatrick Sanan 90948a10b22SPatrick Sanan /*@ 91048a10b22SPatrick Sanan PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous) 91148a10b22SPatrick Sanan 91248a10b22SPatrick Sanan Not Collective 91348a10b22SPatrick Sanan 91448a10b22SPatrick Sanan Input Parameter: 91548a10b22SPatrick Sanan . pc - the preconditioner context 91648a10b22SPatrick Sanan 91748a10b22SPatrick Sanan Output Parameter: 91848a10b22SPatrick Sanan . subcommtype - the subcommunicator type (see PetscSubcommType) 91948a10b22SPatrick Sanan 92048a10b22SPatrick Sanan Level: advanced 92148a10b22SPatrick Sanan 92248a10b22SPatrick Sanan .keywords: PC, telescoping solve 92348a10b22SPatrick Sanan 9241dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE 92548a10b22SPatrick Sanan @*/ 9261dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) 92748a10b22SPatrick Sanan { 92848a10b22SPatrick Sanan PetscErrorCode ierr; 92948a10b22SPatrick Sanan PetscFunctionBegin; 93048a10b22SPatrick Sanan ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr); 93148a10b22SPatrick Sanan PetscFunctionReturn(0); 93248a10b22SPatrick Sanan } 93348a10b22SPatrick Sanan 9341e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/ 9351e07b27eSBarry Smith /*MC 9361e07b27eSBarry Smith PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 9371e07b27eSBarry Smith 9381e07b27eSBarry Smith Options Database: 939a04a6428SPatrick Sanan + -pc_telescope_reduction_factor <r> - factor to use communicator size by. e.g. with 64 MPI processes and r=4, the new sub-communicator will have 64/4 = 16 ranks. 940a04a6428SPatrick Sanan - -pc_telescope_ignore_dm - flag to indicate whether an attached DM should be ignored 941a04a6428SPatrick Sanan - -pc_telescope_subcomm_type <interlaced,contiguous> - how to define the reduced communicator. see PetscSubcomm for more. 9421e07b27eSBarry Smith 9431e07b27eSBarry Smith Level: advanced 9441e07b27eSBarry Smith 9451e07b27eSBarry Smith Notes: 9466fc41876SBarry Smith The preconditioner is deemed telescopic as it only calls KSPSolve() on a single 9478439623fSPatrick Sanan sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 9488439623fSPatrick Sanan This means there will be MPI processes which will be idle during the application of this preconditioner. 9496fc41876SBarry Smith 9501e07b27eSBarry Smith The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 9511e07b27eSBarry Smith Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 9521e07b27eSBarry Smith Currently only support for re-partitioning a DMDA is provided. 9538439623fSPatrick Sanan Any nullspace attached to the original Bmat operator is extracted, re-partitioned and set on the repartitioned Bmat operator. 9541e07b27eSBarry Smith KSPSetComputeOperators() is not propagated to the sub KSP. 9551e07b27eSBarry Smith Currently there is no support for the flag -pc_use_amat 9561e07b27eSBarry Smith 9576fc41876SBarry Smith Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 9588439623fSPatrick Sanan creates a child sub-communicator (c') containing fewer MPI processes than the original parent preconditioner (PC). 9596fc41876SBarry Smith 9606fc41876SBarry Smith Developer Notes: 9616fc41876SBarry Smith During PCSetup, the B operator is scattered onto c'. 9626fc41876SBarry Smith Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 9638439623fSPatrick Sanan Then, KSPSolve() is executed on the c' communicator. 9646fc41876SBarry Smith 9656fc41876SBarry Smith The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 966a04a6428SPatrick 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. 9676fc41876SBarry Smith 968005d9f20SPatrick Sanan The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator. 9698439623fSPatrick Sanan In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into 970005d9f20SPatrick Sanan a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c') 9716fc41876SBarry Smith 9726fc41876SBarry Smith The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D - 9736fc41876SBarry Smith 1D support for 1D DMDAs is not provided), a new DMDA is created on c' (e.g. it is re-partitioned), and this new DM 9746fc41876SBarry Smith is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support 9758439623fSPatrick Sanan for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 9766fc41876SBarry Smith 9776fc41876SBarry Smith By default, B' is defined by simply fusing rows from different MPI processes 9786fc41876SBarry Smith 9798439623fSPatrick Sanan When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B 9807dae84e0SHong Zhang into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the 9816fc41876SBarry Smith locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 9826fc41876SBarry Smith 9838439623fSPatrick Sanan Limitations/improvements include the following. 9848439623fSPatrick Sanan VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage. 9856fc41876SBarry Smith 9866fc41876SBarry Smith The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P, 9878439623fSPatrick 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 9888439623fSPatrick Sanan VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a 9896fc41876SBarry Smith consistent manner. 9906fc41876SBarry Smith 9918439623fSPatrick Sanan Mapping of vectors is performed in the following way. 9928439623fSPatrick Sanan Suppose the parent comm size was 4, and we set a reduction factor of 2; this would give a comm size on c' of 2. 9938439623fSPatrick Sanan Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2. 9948439623fSPatrick Sanan We perform the scatter to the sub-comm in the following way. 9956fc41876SBarry Smith [1] Given a vector x defined on comm c 9966fc41876SBarry Smith 9976fc41876SBarry Smith rank(c) : _________ 0 ______ ________ 1 _______ ________ 2 _____________ ___________ 3 __________ 9986fc41876SBarry Smith x : [0, 1, 2, 3, 4, 5] [6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17] [18, 19, 20, 21, 22, 23] 9996fc41876SBarry Smith 10006fc41876SBarry Smith scatter to xtmp defined also on comm c so that we have the following values 10016fc41876SBarry Smith 10026fc41876SBarry Smith rank(c) : ___________________ 0 ________________ _1_ ______________________ 2 _______________________ __3_ 10036fc41876SBarry Smith xtmp : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [ ] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] [ ] 10046fc41876SBarry Smith 10056fc41876SBarry Smith The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 10066fc41876SBarry Smith 10076fc41876SBarry Smith 10086fc41876SBarry Smith [2] Copy the value from rank 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'. 10096fc41876SBarry Smith Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 10106fc41876SBarry Smith 10116fc41876SBarry Smith rank(c') : ___________________ 0 _______________ ______________________ 1 _____________________ 10126fc41876SBarry Smith xred : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] 10136fc41876SBarry Smith 10141e07b27eSBarry Smith 10151e07b27eSBarry Smith Contributed by Dave May 10161e07b27eSBarry Smith 1017bf00f589SPatrick Sanan Reference: 1018bf00f589SPatrick 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 1019bf00f589SPatrick Sanan 10206fc41876SBarry Smith .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 10211e07b27eSBarry Smith M*/ 10221e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 10231e07b27eSBarry Smith { 10241e07b27eSBarry Smith PetscErrorCode ierr; 10251e07b27eSBarry Smith struct _PC_Telescope *sred; 10261e07b27eSBarry Smith 10271e07b27eSBarry Smith PetscFunctionBegin; 10281e07b27eSBarry Smith ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 102948a10b22SPatrick Sanan sred->subcommtype = PETSC_SUBCOMM_INTERLACED; 10301e07b27eSBarry Smith sred->redfactor = 1; 10311e07b27eSBarry Smith sred->ignore_dm = PETSC_FALSE; 10327c5279cbSDave May sred->ignore_kspcomputeoperators = PETSC_FALSE; 10331e07b27eSBarry Smith pc->data = (void*)sred; 10341e07b27eSBarry Smith 10351e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope; 10361e07b27eSBarry Smith pc->ops->applytranspose = NULL; 1037f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope; 10381e07b27eSBarry Smith pc->ops->setup = PCSetUp_Telescope; 10391e07b27eSBarry Smith pc->ops->destroy = PCDestroy_Telescope; 10401e07b27eSBarry Smith pc->ops->reset = PCReset_Telescope; 10411e07b27eSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Telescope; 10421e07b27eSBarry Smith pc->ops->view = PCView_Telescope; 10431e07b27eSBarry Smith 10441e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 10451e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 10461e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 10471e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 10481e07b27eSBarry Smith 10491e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 105048a10b22SPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr); 105148a10b22SPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr); 10521e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 10531e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 10541e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 10551e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 10560ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 10570ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 10581e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 10591e07b27eSBarry Smith PetscFunctionReturn(0); 10601e07b27eSBarry Smith } 1061