11e07b27eSBarry Smith 21e07b27eSBarry Smith 31e07b27eSBarry Smith 4120bdd93SDave May #include <petsc/private/matimpl.h> 56fc41876SBarry Smith #include <petsc/private/pcimpl.h> 61e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/ 71e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/ 81e07b27eSBarry Smith 91e07b27eSBarry Smith #include "telescope.h" 101e07b27eSBarry Smith 111e07b27eSBarry Smith /* 121e07b27eSBarry Smith PCTelescopeSetUp_default() 131e07b27eSBarry Smith PCTelescopeMatCreate_default() 141e07b27eSBarry Smith 151e07b27eSBarry Smith default 161e07b27eSBarry Smith 171e07b27eSBarry Smith // scatter in 181e07b27eSBarry Smith x(comm) -> xtmp(comm) 191e07b27eSBarry Smith 201e07b27eSBarry Smith xred(subcomm) <- xtmp 211e07b27eSBarry Smith yred(subcomm) 221e07b27eSBarry Smith 231e07b27eSBarry Smith yred(subcomm) --> xtmp 241e07b27eSBarry Smith 251e07b27eSBarry Smith // scatter out 261e07b27eSBarry Smith xtmp(comm) -> y(comm) 271e07b27eSBarry Smith */ 281e07b27eSBarry Smith 291e07b27eSBarry Smith PetscBool isActiveRank(PetscSubcomm scomm) 301e07b27eSBarry Smith { 311e07b27eSBarry Smith if (scomm->color == 0) { return PETSC_TRUE; } 321e07b27eSBarry Smith else { return PETSC_FALSE; } 331e07b27eSBarry Smith } 341e07b27eSBarry Smith 351e07b27eSBarry Smith #undef __FUNCT__ 361e07b27eSBarry Smith #define __FUNCT__ "private_PCTelescopeGetSubDM" 371e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred) 381e07b27eSBarry Smith { 39c6a0d831SBarry Smith DM subdm = NULL; 401e07b27eSBarry Smith 411e07b27eSBarry Smith if (!isActiveRank(sred->psubcomm)) { subdm = NULL; } 421e07b27eSBarry Smith else { 431e07b27eSBarry Smith switch (sred->sr_type) { 441e07b27eSBarry Smith case TELESCOPE_DEFAULT: subdm = NULL; 451e07b27eSBarry Smith break; 461e07b27eSBarry Smith case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 471e07b27eSBarry Smith break; 481e07b27eSBarry Smith case TELESCOPE_DMPLEX: subdm = NULL; 491e07b27eSBarry Smith break; 501e07b27eSBarry Smith } 511e07b27eSBarry Smith } 521e07b27eSBarry Smith return(subdm); 531e07b27eSBarry Smith } 541e07b27eSBarry Smith 551e07b27eSBarry Smith #undef __FUNCT__ 561e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_default" 571e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 581e07b27eSBarry Smith { 591e07b27eSBarry Smith PetscErrorCode ierr; 601e07b27eSBarry Smith PetscInt m,M,bs,st,ed; 611e07b27eSBarry Smith Vec x,xred,yred,xtmp; 621e07b27eSBarry Smith Mat B; 631e07b27eSBarry Smith MPI_Comm comm,subcomm; 641e07b27eSBarry Smith VecScatter scatter; 651e07b27eSBarry Smith IS isin; 661e07b27eSBarry Smith 671e07b27eSBarry Smith PetscFunctionBegin; 681e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr); 691e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 701e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 711e07b27eSBarry Smith 721e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 731e07b27eSBarry Smith ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 741e07b27eSBarry Smith ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 751e07b27eSBarry Smith ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 761e07b27eSBarry Smith 771e07b27eSBarry Smith xred = NULL; 783ac26c5eSBarry Smith m = 0; 791e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 801e07b27eSBarry Smith ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr); 811e07b27eSBarry Smith ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr); 821e07b27eSBarry Smith ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr); 831e07b27eSBarry Smith ierr = VecSetFromOptions(xred);CHKERRQ(ierr); 841e07b27eSBarry Smith ierr = VecGetLocalSize(xred,&m); 851e07b27eSBarry Smith } 861e07b27eSBarry Smith 871e07b27eSBarry Smith yred = NULL; 881e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 891e07b27eSBarry Smith ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 901e07b27eSBarry Smith } 911e07b27eSBarry Smith 921e07b27eSBarry Smith ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 931e07b27eSBarry Smith ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 941e07b27eSBarry Smith ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 951e07b27eSBarry Smith ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 961e07b27eSBarry Smith 971e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 981e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 991e07b27eSBarry Smith ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr); 1001e07b27eSBarry Smith } else { 1011e07b27eSBarry Smith ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 1023ac26c5eSBarry Smith ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 1031e07b27eSBarry Smith } 1041e07b27eSBarry Smith ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 1051e07b27eSBarry Smith 1061e07b27eSBarry Smith ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 1071e07b27eSBarry Smith 1081e07b27eSBarry Smith sred->isin = isin; 1091e07b27eSBarry Smith sred->scatter = scatter; 1101e07b27eSBarry Smith sred->xred = xred; 1111e07b27eSBarry Smith sred->yred = yred; 1121e07b27eSBarry Smith sred->xtmp = xtmp; 1131e07b27eSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1141e07b27eSBarry Smith PetscFunctionReturn(0); 1151e07b27eSBarry Smith } 1161e07b27eSBarry Smith 1171e07b27eSBarry Smith #undef __FUNCT__ 1181e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatCreate_default" 1191e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 1201e07b27eSBarry Smith { 1211e07b27eSBarry Smith PetscErrorCode ierr; 1221e07b27eSBarry Smith MPI_Comm comm,subcomm; 1231e07b27eSBarry Smith Mat Bred,B; 1241e07b27eSBarry Smith PetscInt nr,nc; 1251e07b27eSBarry Smith IS isrow,iscol; 1261e07b27eSBarry Smith Mat Blocal,*_Blocal; 1271e07b27eSBarry Smith 1281e07b27eSBarry Smith PetscFunctionBegin; 1291e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr); 1301e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1311e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1321e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 1331e07b27eSBarry Smith ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 1341e07b27eSBarry Smith isrow = sred->isin; 1351e07b27eSBarry Smith ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 1361e07b27eSBarry Smith ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 1371e07b27eSBarry Smith Blocal = *_Blocal; 1381e07b27eSBarry Smith ierr = PetscFree(_Blocal);CHKERRQ(ierr); 1391e07b27eSBarry Smith Bred = NULL; 1401e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 1411e07b27eSBarry Smith PetscInt mm; 1421e07b27eSBarry Smith 1431e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 1441e07b27eSBarry Smith 1451e07b27eSBarry Smith ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 1461e07b27eSBarry Smith ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 1471e07b27eSBarry Smith } 1481e07b27eSBarry Smith *A = Bred; 1491e07b27eSBarry Smith ierr = ISDestroy(&iscol);CHKERRQ(ierr); 1501e07b27eSBarry Smith ierr = MatDestroy(&Blocal);CHKERRQ(ierr); 1511e07b27eSBarry Smith PetscFunctionReturn(0); 1521e07b27eSBarry Smith } 1531e07b27eSBarry Smith 1541e07b27eSBarry Smith #undef __FUNCT__ 1551e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default" 1561e07b27eSBarry Smith PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 1571e07b27eSBarry Smith { 1581e07b27eSBarry Smith PetscErrorCode ierr; 1591e07b27eSBarry Smith MatNullSpace nullspace,sub_nullspace; 1601e07b27eSBarry Smith Mat A,B; 1611e07b27eSBarry Smith PetscBool has_const; 1621e07b27eSBarry Smith PetscInt i,k,n; 1631e07b27eSBarry Smith const Vec *vecs; 1641e07b27eSBarry Smith Vec *sub_vecs; 1651e07b27eSBarry Smith MPI_Comm subcomm; 1661e07b27eSBarry Smith 1671e07b27eSBarry Smith PetscFunctionBegin; 1681e07b27eSBarry Smith ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr); 1691e07b27eSBarry Smith ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 1701e07b27eSBarry Smith if (!nullspace) return(0); 1711e07b27eSBarry Smith 1721e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr); 1731e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1741e07b27eSBarry Smith ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 1751e07b27eSBarry Smith 1761e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 177e3acf2f7SBarry Smith if (n) { 178e3acf2f7SBarry Smith ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 1791e07b27eSBarry Smith } 1801e07b27eSBarry Smith } 1811e07b27eSBarry Smith 1821e07b27eSBarry Smith /* copy entries */ 1831e07b27eSBarry Smith for (k=0; k<n; k++) { 1841e07b27eSBarry Smith const PetscScalar *x_array; 1851e07b27eSBarry Smith PetscScalar *LA_sub_vec; 1861e07b27eSBarry Smith PetscInt st,ed,bs; 1871e07b27eSBarry Smith 1881e07b27eSBarry Smith /* pull in vector x->xtmp */ 1891e07b27eSBarry Smith ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1901e07b27eSBarry Smith ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1911e07b27eSBarry Smith /* copy vector entires into xred */ 1921e07b27eSBarry Smith ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr); 1931e07b27eSBarry Smith ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 1941e07b27eSBarry Smith if (sub_vecs[k]) { 1951e07b27eSBarry Smith ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 1961e07b27eSBarry Smith ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 1971e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 1981e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 1991e07b27eSBarry Smith } 2001e07b27eSBarry Smith ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 2011e07b27eSBarry Smith } 2021e07b27eSBarry Smith ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 2031e07b27eSBarry Smith } 2041e07b27eSBarry Smith 2051e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 2061e07b27eSBarry Smith /* create new nullspace for redundant object */ 2071e07b27eSBarry Smith ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr); 208120bdd93SDave May sub_nullspace->remove = nullspace->remove; 209120bdd93SDave May sub_nullspace->rmctx = nullspace->rmctx; 210120bdd93SDave May 2111e07b27eSBarry Smith /* attach redundant nullspace to Bred */ 2121e07b27eSBarry Smith ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 213e3acf2f7SBarry Smith ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 2141e07b27eSBarry Smith } 2151e07b27eSBarry Smith PetscFunctionReturn(0); 2161e07b27eSBarry Smith } 2171e07b27eSBarry Smith 2181e07b27eSBarry Smith #undef __FUNCT__ 2191e07b27eSBarry Smith #define __FUNCT__ "PCView_Telescope" 2201e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 2211e07b27eSBarry Smith { 2221e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 2231e07b27eSBarry Smith PetscErrorCode ierr; 2241e07b27eSBarry Smith PetscBool iascii,isstring; 2251e07b27eSBarry Smith PetscViewer subviewer; 2261e07b27eSBarry Smith 2271e07b27eSBarry Smith PetscFunctionBegin; 2281e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2291e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 2301e07b27eSBarry Smith if (iascii) { 2311e07b27eSBarry Smith if (!sred->psubcomm) { 2321e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Telescope: preconditioner not yet setup\n");CHKERRQ(ierr); 2331e07b27eSBarry Smith } else { 2341e07b27eSBarry Smith MPI_Comm comm,subcomm; 2351e07b27eSBarry Smith PetscMPIInt comm_size,subcomm_size; 2361e07b27eSBarry Smith DM dm,subdm; 2371e07b27eSBarry Smith 2381e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 2391e07b27eSBarry Smith subdm = private_PCTelescopeGetSubDM(sred); 2401e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 2411e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2421e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 2431e07b27eSBarry Smith ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 2441e07b27eSBarry Smith 2451e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 2461e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 2471e07b27eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 2481e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 2491e07b27eSBarry Smith ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 2501e07b27eSBarry Smith 2511e07b27eSBarry Smith if (dm && sred->ignore_dm) { 2521e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," Telescope: ignoring DM\n");CHKERRQ(ierr); 2531e07b27eSBarry Smith } 2547c5279cbSDave May if (sred->ignore_kspcomputeoperators) { 2557c5279cbSDave May ierr = PetscViewerASCIIPrintf(subviewer," Telescope: ignoring KSPComputeOperators\n");CHKERRQ(ierr); 2567c5279cbSDave May } 2571e07b27eSBarry Smith switch (sred->sr_type) { 2581e07b27eSBarry Smith case TELESCOPE_DEFAULT: 2591e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," Telescope: using default setup\n");CHKERRQ(ierr); 2601e07b27eSBarry Smith break; 2611e07b27eSBarry Smith case TELESCOPE_DMDA: 2621e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMDA detected\n");CHKERRQ(ierr); 2631e07b27eSBarry Smith ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr); 2641e07b27eSBarry Smith break; 2651e07b27eSBarry Smith case TELESCOPE_DMPLEX: 2661e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMPLEX detected\n");CHKERRQ(ierr); 2671e07b27eSBarry Smith break; 2681e07b27eSBarry Smith } 2691e07b27eSBarry Smith 2701e07b27eSBarry Smith ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 2711e07b27eSBarry Smith ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 2721e07b27eSBarry Smith } 2731e07b27eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 2741e07b27eSBarry Smith } 2751e07b27eSBarry Smith } 2761e07b27eSBarry Smith PetscFunctionReturn(0); 2771e07b27eSBarry Smith } 2781e07b27eSBarry Smith 2791e07b27eSBarry Smith #undef __FUNCT__ 2801e07b27eSBarry Smith #define __FUNCT__ "PCSetUp_Telescope" 2811e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc) 2821e07b27eSBarry Smith { 2831e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 2841e07b27eSBarry Smith PetscErrorCode ierr; 285bd49479cSSatish Balay MPI_Comm comm,subcomm=0; 2861e07b27eSBarry Smith PCTelescopeType sr_type; 2871e07b27eSBarry Smith 2881e07b27eSBarry Smith PetscFunctionBegin; 2891e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 2901e07b27eSBarry Smith 2911e07b27eSBarry Smith /* subcomm definition */ 2921e07b27eSBarry Smith if (!pc->setupcalled) { 2931e07b27eSBarry Smith if (!sred->psubcomm) { 2941e07b27eSBarry Smith ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 2951e07b27eSBarry Smith ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 2961e07b27eSBarry Smith ierr = PetscSubcommSetType(sred->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 2971e07b27eSBarry Smith /* disable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 2981e07b27eSBarry Smith /* ierr = PetscSubcommSetFromOptions(sred->psubcomm);CHKERRQ(ierr); */ 2991e07b27eSBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3001e07b27eSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3011e07b27eSBarry Smith } 3021e07b27eSBarry Smith } 3030f6f40a7SSatish Balay subcomm = PetscSubcommChild(sred->psubcomm); 3041e07b27eSBarry Smith 3051e07b27eSBarry Smith /* internal KSP */ 3061e07b27eSBarry Smith if (!pc->setupcalled) { 3071e07b27eSBarry Smith const char *prefix; 3081e07b27eSBarry Smith 3091e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 3101e07b27eSBarry Smith ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 3111e07b27eSBarry Smith ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 3121e07b27eSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3131e07b27eSBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 3141e07b27eSBarry Smith ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 31514c9fce5SDave May ierr = KSPSetInitialGuessNonzero(sred->ksp,pc->nonzero_guess);CHKERRQ(ierr); 3161e07b27eSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3171e07b27eSBarry Smith ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 3181e07b27eSBarry Smith ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 3191e07b27eSBarry Smith } 3201e07b27eSBarry Smith } 3211e07b27eSBarry Smith /* Determine type of setup/update */ 3221e07b27eSBarry Smith if (!pc->setupcalled) { 3231e07b27eSBarry Smith PetscBool has_dm,same; 3241e07b27eSBarry Smith DM dm; 3251e07b27eSBarry Smith 3261e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 3271e07b27eSBarry Smith has_dm = PETSC_FALSE; 3281e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3291e07b27eSBarry Smith if (dm) { has_dm = PETSC_TRUE; } 3301e07b27eSBarry Smith if (has_dm) { 3311e07b27eSBarry Smith /* check for dmda */ 3321e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 3331e07b27eSBarry Smith if (same) { 3341e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 3351e07b27eSBarry Smith sr_type = TELESCOPE_DMDA; 3361e07b27eSBarry Smith } 3371e07b27eSBarry Smith /* check for dmplex */ 3381e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 3391e07b27eSBarry Smith if (same) { 3401e07b27eSBarry Smith PetscInfo(pc,"PCTelescope: found DMPLEX\n"); 3411e07b27eSBarry Smith sr_type = TELESCOPE_DMPLEX; 3421e07b27eSBarry Smith } 3431e07b27eSBarry Smith } 3441e07b27eSBarry Smith 3451e07b27eSBarry Smith if (sred->ignore_dm) { 3461e07b27eSBarry Smith PetscInfo(pc,"PCTelescope: ignore DM\n"); 3471e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 3481e07b27eSBarry Smith } 3491e07b27eSBarry Smith sred->sr_type = sr_type; 3501e07b27eSBarry Smith } else { 3511e07b27eSBarry Smith sr_type = sred->sr_type; 3521e07b27eSBarry Smith } 3531e07b27eSBarry Smith 3541e07b27eSBarry Smith /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */ 3551e07b27eSBarry Smith switch (sr_type) { 3561e07b27eSBarry Smith case TELESCOPE_DEFAULT: 3571e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 3581e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 3591e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 3601e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 3611e07b27eSBarry Smith break; 3621e07b27eSBarry Smith case TELESCOPE_DMDA: 3631e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope_dmda; 3641e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 3651e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 3661e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 3671e07b27eSBarry Smith sred->pctelescope_reset_type = PCReset_Telescope_dmda; 3681e07b27eSBarry Smith break; 3691e07b27eSBarry Smith case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available"); 3701e07b27eSBarry Smith break; 3711e07b27eSBarry Smith default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided"); 3721e07b27eSBarry Smith break; 3731e07b27eSBarry Smith } 3741e07b27eSBarry Smith 3751e07b27eSBarry Smith /* setup */ 3761e07b27eSBarry Smith if (sred->pctelescope_setup_type) { 3771e07b27eSBarry Smith ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 3781e07b27eSBarry Smith } 3791e07b27eSBarry Smith /* update */ 3801e07b27eSBarry Smith if (!pc->setupcalled) { 3811e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 3821e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 3831e07b27eSBarry Smith } 3841e07b27eSBarry Smith if (sred->pctelescope_matnullspacecreate_type) { 3851e07b27eSBarry Smith ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 3861e07b27eSBarry Smith } 3871e07b27eSBarry Smith } else { 3881e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 3891e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 3901e07b27eSBarry Smith } 3911e07b27eSBarry Smith } 3921e07b27eSBarry Smith 3931e07b27eSBarry Smith /* common - no construction */ 3941e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 3951e07b27eSBarry Smith ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 3961e07b27eSBarry Smith if (pc->setfromoptionscalled && !pc->setupcalled){ 3971e07b27eSBarry Smith ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 3981e07b27eSBarry Smith } 3991e07b27eSBarry Smith } 4001e07b27eSBarry Smith PetscFunctionReturn(0); 4011e07b27eSBarry Smith } 4021e07b27eSBarry Smith 4031e07b27eSBarry Smith #undef __FUNCT__ 4041e07b27eSBarry Smith #define __FUNCT__ "PCApply_Telescope" 4051e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 4061e07b27eSBarry Smith { 4071e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 4081e07b27eSBarry Smith PetscErrorCode ierr; 4091e07b27eSBarry Smith Vec xtmp,xred,yred; 4101e07b27eSBarry Smith PetscInt i,st,ed,bs; 4111e07b27eSBarry Smith VecScatter scatter; 4121e07b27eSBarry Smith PetscScalar *array; 4131e07b27eSBarry Smith const PetscScalar *x_array; 4141e07b27eSBarry Smith 4151e07b27eSBarry Smith PetscFunctionBegin; 4161e07b27eSBarry Smith xtmp = sred->xtmp; 4171e07b27eSBarry Smith scatter = sred->scatter; 4181e07b27eSBarry Smith xred = sred->xred; 4191e07b27eSBarry Smith yred = sred->yred; 4201e07b27eSBarry Smith 42114c9fce5SDave May if (pc->nonzero_guess) { 42214c9fce5SDave May ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero-initial guess\n");CHKERRQ(ierr); 42314c9fce5SDave May /* pull in vector y->xtmp */ 42414c9fce5SDave May ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 42514c9fce5SDave May ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 42614c9fce5SDave May 42714c9fce5SDave May /* copy vector entires into xred */ 42814c9fce5SDave May ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 42914c9fce5SDave May ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 43014c9fce5SDave May if (yred) { 43114c9fce5SDave May PetscScalar *LA_yred; 43214c9fce5SDave May ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 43314c9fce5SDave May ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 43414c9fce5SDave May for (i=0; i<ed-st; i++) { 43514c9fce5SDave May LA_yred[i] = x_array[i]; 43614c9fce5SDave May } 43714c9fce5SDave May ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 43814c9fce5SDave May } 43914c9fce5SDave May ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 44014c9fce5SDave May } 44114c9fce5SDave May 4421e07b27eSBarry Smith /* pull in vector x->xtmp */ 4431e07b27eSBarry Smith ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4441e07b27eSBarry Smith ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4451e07b27eSBarry Smith 4461e07b27eSBarry Smith /* copy vector entires into xred */ 4471e07b27eSBarry Smith ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 4481e07b27eSBarry Smith ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 4491e07b27eSBarry Smith if (xred) { 4501e07b27eSBarry Smith PetscScalar *LA_xred; 4511e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 4521e07b27eSBarry Smith ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 4531e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 4541e07b27eSBarry Smith LA_xred[i] = x_array[i]; 4551e07b27eSBarry Smith } 4561e07b27eSBarry Smith ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 4571e07b27eSBarry Smith } 4581e07b27eSBarry Smith ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 4591e07b27eSBarry Smith /* solve */ 4601e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 4611e07b27eSBarry Smith ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 4621e07b27eSBarry Smith } 4631e07b27eSBarry Smith /* return vector */ 4641e07b27eSBarry Smith ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 4651e07b27eSBarry Smith ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 4661e07b27eSBarry Smith if (yred) { 4671e07b27eSBarry Smith const PetscScalar *LA_yred; 4681e07b27eSBarry Smith ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 4691e07b27eSBarry Smith ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 4701e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 4711e07b27eSBarry Smith array[i] = LA_yred[i]; 4721e07b27eSBarry Smith } 4731e07b27eSBarry Smith ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 4741e07b27eSBarry Smith } 4751e07b27eSBarry Smith ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 4761e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4771e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4781e07b27eSBarry Smith PetscFunctionReturn(0); 4791e07b27eSBarry Smith } 4801e07b27eSBarry Smith 4811e07b27eSBarry Smith #undef __FUNCT__ 4821e07b27eSBarry Smith #define __FUNCT__ "PCReset_Telescope" 4831e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc) 4841e07b27eSBarry Smith { 4851e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 4861e07b27eSBarry Smith PetscErrorCode ierr; 4871e07b27eSBarry Smith 4881e07b27eSBarry Smith ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 4891e07b27eSBarry Smith ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 490e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); 491e3acf2f7SBarry Smith ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); 492e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); 493e3acf2f7SBarry Smith ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); 494e3acf2f7SBarry Smith ierr = KSPReset(sred->ksp);CHKERRQ(ierr); 4951e07b27eSBarry Smith if (sred->pctelescope_reset_type) { 4961e07b27eSBarry Smith ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 4971e07b27eSBarry Smith } 4981e07b27eSBarry Smith PetscFunctionReturn(0); 4991e07b27eSBarry Smith } 5001e07b27eSBarry Smith 5011e07b27eSBarry Smith #undef __FUNCT__ 5021e07b27eSBarry Smith #define __FUNCT__ "PCDestroy_Telescope" 5031e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc) 5041e07b27eSBarry Smith { 5051e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5061e07b27eSBarry Smith PetscErrorCode ierr; 5071e07b27eSBarry Smith 5081e07b27eSBarry Smith PetscFunctionBegin; 5091e07b27eSBarry Smith ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 510e3acf2f7SBarry Smith ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); 5111e07b27eSBarry Smith ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 512e3acf2f7SBarry Smith ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr); 513e3acf2f7SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 5141e07b27eSBarry Smith PetscFunctionReturn(0); 5151e07b27eSBarry Smith } 5161e07b27eSBarry Smith 5171e07b27eSBarry Smith #undef __FUNCT__ 5181e07b27eSBarry Smith #define __FUNCT__ "PCSetFromOptions_Telescope" 5194416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc) 5201e07b27eSBarry Smith { 5211e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 5221e07b27eSBarry Smith PetscErrorCode ierr; 5231e07b27eSBarry Smith MPI_Comm comm; 5241e07b27eSBarry Smith PetscMPIInt size; 5251e07b27eSBarry Smith 5261e07b27eSBarry Smith PetscFunctionBegin; 5271e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5281e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 5291e07b27eSBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 5301e07b27eSBarry Smith ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 5311e07b27eSBarry Smith if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 5321e07b27eSBarry Smith ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 5337c5279cbSDave May ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr); 5341e07b27eSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5351e07b27eSBarry Smith PetscFunctionReturn(0); 5361e07b27eSBarry Smith } 5371e07b27eSBarry Smith 5381e07b27eSBarry Smith /* PC simplementation specific API's */ 5391e07b27eSBarry Smith 5401e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 5411e07b27eSBarry Smith { 5421e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 543bd49479cSSatish Balay PetscFunctionBegin; 5441e07b27eSBarry Smith if (ksp) *ksp = red->ksp; 545bd49479cSSatish Balay PetscFunctionReturn(0); 5461e07b27eSBarry Smith } 5471e07b27eSBarry Smith 5481e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 5491e07b27eSBarry Smith { 5501e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 551bd49479cSSatish Balay PetscFunctionBegin; 5521e07b27eSBarry Smith if (fact) *fact = red->redfactor; 553bd49479cSSatish Balay PetscFunctionReturn(0); 5541e07b27eSBarry Smith } 5551e07b27eSBarry Smith 5561e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 5571e07b27eSBarry Smith { 5581e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 5591e07b27eSBarry Smith PetscMPIInt size; 5601e07b27eSBarry Smith PetscErrorCode ierr; 5611e07b27eSBarry Smith 562bd49479cSSatish Balay PetscFunctionBegin; 5631e07b27eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 5641e07b27eSBarry Smith if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 5651e07b27eSBarry Smith if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 5661e07b27eSBarry Smith red->redfactor = fact; 567bd49479cSSatish Balay PetscFunctionReturn(0); 5681e07b27eSBarry Smith } 5691e07b27eSBarry Smith 5701e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 5711e07b27eSBarry Smith { 5721e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 573bd49479cSSatish Balay PetscFunctionBegin; 5741e07b27eSBarry Smith if (v) *v = red->ignore_dm; 575bd49479cSSatish Balay PetscFunctionReturn(0); 5761e07b27eSBarry Smith } 5771e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 5781e07b27eSBarry Smith { 5791e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 580bd49479cSSatish Balay PetscFunctionBegin; 5811e07b27eSBarry Smith red->ignore_dm = v; 582bd49479cSSatish Balay PetscFunctionReturn(0); 5831e07b27eSBarry Smith } 5841e07b27eSBarry Smith 585*0ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v) 586*0ae7c45bSDave May { 587*0ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 588*0ae7c45bSDave May PetscFunctionBegin; 589*0ae7c45bSDave May if (v) *v = red->ignore_kspcomputeoperators; 590*0ae7c45bSDave May PetscFunctionReturn(0); 591*0ae7c45bSDave May } 592*0ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v) 593*0ae7c45bSDave May { 594*0ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 595*0ae7c45bSDave May PetscFunctionBegin; 596*0ae7c45bSDave May red->ignore_kspcomputeoperators = v; 597*0ae7c45bSDave May PetscFunctionReturn(0); 598*0ae7c45bSDave May } 599*0ae7c45bSDave May 6001e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 6011e07b27eSBarry Smith { 6021e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 603bd49479cSSatish Balay PetscFunctionBegin; 6041e07b27eSBarry Smith *dm = private_PCTelescopeGetSubDM(red); 605bd49479cSSatish Balay PetscFunctionReturn(0); 6061e07b27eSBarry Smith } 6071e07b27eSBarry Smith 6081e07b27eSBarry Smith /*@ 6091e07b27eSBarry Smith PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 6101e07b27eSBarry Smith 6111e07b27eSBarry Smith Not Collective 6121e07b27eSBarry Smith 6131e07b27eSBarry Smith Input Parameter: 6141e07b27eSBarry Smith . pc - the preconditioner context 6151e07b27eSBarry Smith 6161e07b27eSBarry Smith Output Parameter: 6171e07b27eSBarry Smith . subksp - the KSP defined the smaller set of processes 6181e07b27eSBarry Smith 6191e07b27eSBarry Smith Level: advanced 6201e07b27eSBarry Smith 6211e07b27eSBarry Smith .keywords: PC, telescopting solve 6221e07b27eSBarry Smith @*/ 6231e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 6241e07b27eSBarry Smith { 625bd49479cSSatish Balay PetscErrorCode ierr; 626bd49479cSSatish Balay PetscFunctionBegin; 627bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 628bd49479cSSatish Balay PetscFunctionReturn(0); 6291e07b27eSBarry Smith } 6301e07b27eSBarry Smith 6311e07b27eSBarry Smith /*@ 6321e07b27eSBarry Smith PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 6331e07b27eSBarry Smith 6341e07b27eSBarry Smith Not Collective 6351e07b27eSBarry Smith 6361e07b27eSBarry Smith Input Parameter: 6371e07b27eSBarry Smith . pc - the preconditioner context 6381e07b27eSBarry Smith 6391e07b27eSBarry Smith Output Parameter: 6401e07b27eSBarry Smith . fact - the reduction factor 6411e07b27eSBarry Smith 6421e07b27eSBarry Smith Level: advanced 6431e07b27eSBarry Smith 6441e07b27eSBarry Smith .keywords: PC, telescoping solve 6451e07b27eSBarry Smith @*/ 6461e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 6471e07b27eSBarry Smith { 648bd49479cSSatish Balay PetscErrorCode ierr; 649bd49479cSSatish Balay PetscFunctionBegin; 650bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 651bd49479cSSatish Balay PetscFunctionReturn(0); 6521e07b27eSBarry Smith } 6531e07b27eSBarry Smith 6541e07b27eSBarry Smith /*@ 6551e07b27eSBarry Smith PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 6561e07b27eSBarry Smith 6571e07b27eSBarry Smith Not Collective 6581e07b27eSBarry Smith 6591e07b27eSBarry Smith Input Parameter: 6601e07b27eSBarry Smith . pc - the preconditioner context 6611e07b27eSBarry Smith 6621e07b27eSBarry Smith Output Parameter: 6631e07b27eSBarry Smith . fact - the reduction factor 6641e07b27eSBarry Smith 6651e07b27eSBarry Smith Level: advanced 6661e07b27eSBarry Smith 6671e07b27eSBarry Smith .keywords: PC, telescoping solve 6681e07b27eSBarry Smith @*/ 6691e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 6701e07b27eSBarry Smith { 671bd49479cSSatish Balay PetscErrorCode ierr; 672bd49479cSSatish Balay PetscFunctionBegin; 673bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 674bd49479cSSatish Balay PetscFunctionReturn(0); 6751e07b27eSBarry Smith } 6761e07b27eSBarry Smith 6771e07b27eSBarry Smith /*@ 6781e07b27eSBarry Smith PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 6791e07b27eSBarry Smith 6801e07b27eSBarry Smith Not Collective 6811e07b27eSBarry Smith 6821e07b27eSBarry Smith Input Parameter: 6831e07b27eSBarry Smith . pc - the preconditioner context 6841e07b27eSBarry Smith 6851e07b27eSBarry Smith Output Parameter: 6861e07b27eSBarry Smith . v - the flag 6871e07b27eSBarry Smith 6881e07b27eSBarry Smith Level: advanced 6891e07b27eSBarry Smith 6901e07b27eSBarry Smith .keywords: PC, telescoping solve 6911e07b27eSBarry Smith @*/ 6921e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 6931e07b27eSBarry Smith { 694bd49479cSSatish Balay PetscErrorCode ierr; 695bd49479cSSatish Balay PetscFunctionBegin; 696bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 697bd49479cSSatish Balay PetscFunctionReturn(0); 6981e07b27eSBarry Smith } 6991e07b27eSBarry Smith 7001e07b27eSBarry Smith /*@ 7011e07b27eSBarry Smith PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 7021e07b27eSBarry Smith 7031e07b27eSBarry Smith Not Collective 7041e07b27eSBarry Smith 7051e07b27eSBarry Smith Input Parameter: 7061e07b27eSBarry Smith . pc - the preconditioner context 7071e07b27eSBarry Smith 7081e07b27eSBarry Smith Output Parameter: 7091e07b27eSBarry Smith . v - Use PETSC_TRUE to ignore any DM 7101e07b27eSBarry Smith 7111e07b27eSBarry Smith Level: advanced 7121e07b27eSBarry Smith 7131e07b27eSBarry Smith .keywords: PC, telescoping solve 7141e07b27eSBarry Smith @*/ 715bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 7161e07b27eSBarry Smith { 717bd49479cSSatish Balay PetscErrorCode ierr; 718bd49479cSSatish Balay PetscFunctionBegin; 719bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 720bd49479cSSatish Balay PetscFunctionReturn(0); 7211e07b27eSBarry Smith } 7221e07b27eSBarry Smith 7231e07b27eSBarry Smith /*@ 724*0ae7c45bSDave May PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used. 725*0ae7c45bSDave May 726*0ae7c45bSDave May Not Collective 727*0ae7c45bSDave May 728*0ae7c45bSDave May Input Parameter: 729*0ae7c45bSDave May . pc - the preconditioner context 730*0ae7c45bSDave May 731*0ae7c45bSDave May Output Parameter: 732*0ae7c45bSDave May . v - the flag 733*0ae7c45bSDave May 734*0ae7c45bSDave May Level: advanced 735*0ae7c45bSDave May 736*0ae7c45bSDave May .keywords: PC, telescoping solve 737*0ae7c45bSDave May @*/ 738*0ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v) 739*0ae7c45bSDave May { 740*0ae7c45bSDave May PetscErrorCode ierr; 741*0ae7c45bSDave May PetscFunctionBegin; 742*0ae7c45bSDave May ierr = PetscTryMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 743*0ae7c45bSDave May PetscFunctionReturn(0); 744*0ae7c45bSDave May } 745*0ae7c45bSDave May 746*0ae7c45bSDave May /*@ 747*0ae7c45bSDave May PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators. 748*0ae7c45bSDave May 749*0ae7c45bSDave May Not Collective 750*0ae7c45bSDave May 751*0ae7c45bSDave May Input Parameter: 752*0ae7c45bSDave May . pc - the preconditioner context 753*0ae7c45bSDave May 754*0ae7c45bSDave May Output Parameter: 755*0ae7c45bSDave May . v - Use PETSC_TRUE to ignore any DM 756*0ae7c45bSDave May 757*0ae7c45bSDave May Level: advanced 758*0ae7c45bSDave May 759*0ae7c45bSDave May .keywords: PC, telescoping solve 760*0ae7c45bSDave May @*/ 761*0ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v) 762*0ae7c45bSDave May { 763*0ae7c45bSDave May PetscErrorCode ierr; 764*0ae7c45bSDave May PetscFunctionBegin; 765*0ae7c45bSDave May ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 766*0ae7c45bSDave May PetscFunctionReturn(0); 767*0ae7c45bSDave May } 768*0ae7c45bSDave May 769*0ae7c45bSDave May /*@ 7701e07b27eSBarry Smith PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 7711e07b27eSBarry Smith 7721e07b27eSBarry Smith Not Collective 7731e07b27eSBarry Smith 7741e07b27eSBarry Smith Input Parameter: 7751e07b27eSBarry Smith . pc - the preconditioner context 7761e07b27eSBarry Smith 7771e07b27eSBarry Smith Output Parameter: 7781e07b27eSBarry Smith . subdm - The re-partitioned DM 7791e07b27eSBarry Smith 7801e07b27eSBarry Smith Level: advanced 7811e07b27eSBarry Smith 7821e07b27eSBarry Smith .keywords: PC, telescoping solve 7831e07b27eSBarry Smith @*/ 7841e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 7851e07b27eSBarry Smith { 786bd49479cSSatish Balay PetscErrorCode ierr; 787bd49479cSSatish Balay PetscFunctionBegin; 788bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 789bd49479cSSatish Balay PetscFunctionReturn(0); 7901e07b27eSBarry Smith } 7911e07b27eSBarry Smith 7921e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/ 7931e07b27eSBarry Smith /*MC 7941e07b27eSBarry Smith PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 7951e07b27eSBarry Smith 7961e07b27eSBarry Smith Options Database: 7971e07b27eSBarry Smith + -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and 7981e07b27eSBarry Smith use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes 7991e07b27eSBarry Smith - -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored 8001e07b27eSBarry Smith 8011e07b27eSBarry Smith Level: advanced 8021e07b27eSBarry Smith 8031e07b27eSBarry Smith Notes: 8046fc41876SBarry Smith The preconditioner is deemed telescopic as it only calls KSPSolve() on a single 8056fc41876SBarry Smith sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 8066fc41876SBarry Smith This means there will be MPI processes within c, which will be idle during the application of this preconditioner. 8076fc41876SBarry Smith 8081e07b27eSBarry Smith The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 8091e07b27eSBarry Smith Both the B mat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 8101e07b27eSBarry Smith Currently only support for re-partitioning a DMDA is provided. 8111e07b27eSBarry Smith Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator. 8121e07b27eSBarry Smith KSPSetComputeOperators() is not propagated to the sub KSP. 8131e07b27eSBarry Smith Currently there is no support for the flag -pc_use_amat 8141e07b27eSBarry Smith 8156fc41876SBarry Smith Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 8166fc41876SBarry Smith creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC). 8176fc41876SBarry Smith 8186fc41876SBarry Smith Developer Notes: 8196fc41876SBarry Smith During PCSetup, the B operator is scattered onto c'. 8206fc41876SBarry Smith Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 8216fc41876SBarry Smith Then KSPSolve() is executed on the c' communicator. 8226fc41876SBarry Smith 8236fc41876SBarry Smith The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 8246fc41876SBarry Smith creation routine. We run the sub KSP on only the ranks within the communicator which have a color equal to zero. 8256fc41876SBarry Smith 8266fc41876SBarry Smith The telescoping preconditioner is aware of nullspaces which are attached to the only B operator. 8276fc41876SBarry Smith In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into 8286fc41876SBarry Smith a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c') 8296fc41876SBarry Smith 8306fc41876SBarry Smith The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D - 8316fc41876SBarry 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 8326fc41876SBarry Smith is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support 8336fc41876SBarry Smith for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 8346fc41876SBarry Smith 8356fc41876SBarry Smith By default, B' is defined by simply fusing rows from different MPI processes 8366fc41876SBarry Smith 8376fc41876SBarry Smith When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B 8386fc41876SBarry Smith into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the 8396fc41876SBarry Smith locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 8406fc41876SBarry Smith 8416fc41876SBarry Smith Limitations/improvements 8426fc41876SBarry Smith VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage. 8436fc41876SBarry Smith 8446fc41876SBarry Smith The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P, 8456fc41876SBarry Smith and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). I opted to use P^T.A.P as it appears 8466fc41876SBarry Smith VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a 8476fc41876SBarry Smith consistent manner. 8486fc41876SBarry Smith 8496fc41876SBarry Smith Mapping of vectors is performed this way 8506fc41876SBarry Smith Suppose the parent comm size was 4, and we set a reduction factor of 2, thus would give a comm size on c' of 2. 8516fc41876SBarry Smith Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2. 8526fc41876SBarry Smith We perform the scatter to the sub-comm in the following way, 8536fc41876SBarry Smith [1] Given a vector x defined on comm c 8546fc41876SBarry Smith 8556fc41876SBarry Smith rank(c) : _________ 0 ______ ________ 1 _______ ________ 2 _____________ ___________ 3 __________ 8566fc41876SBarry 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] 8576fc41876SBarry Smith 8586fc41876SBarry Smith scatter to xtmp defined also on comm c so that we have the following values 8596fc41876SBarry Smith 8606fc41876SBarry Smith rank(c) : ___________________ 0 ________________ _1_ ______________________ 2 _______________________ __3_ 8616fc41876SBarry 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] [ ] 8626fc41876SBarry Smith 8636fc41876SBarry Smith The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 8646fc41876SBarry Smith 8656fc41876SBarry Smith 8666fc41876SBarry 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'. 8676fc41876SBarry Smith Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 8686fc41876SBarry Smith 8696fc41876SBarry Smith rank(c') : ___________________ 0 _______________ ______________________ 1 _____________________ 8706fc41876SBarry 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] 8716fc41876SBarry Smith 8721e07b27eSBarry Smith 8731e07b27eSBarry Smith Contributed by Dave May 8741e07b27eSBarry Smith 8756fc41876SBarry Smith .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 8761e07b27eSBarry Smith M*/ 8771e07b27eSBarry Smith #undef __FUNCT__ 8781e07b27eSBarry Smith #define __FUNCT__ "PCCreate_Telescope" 8791e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 8801e07b27eSBarry Smith { 8811e07b27eSBarry Smith PetscErrorCode ierr; 8821e07b27eSBarry Smith struct _PC_Telescope *sred; 8831e07b27eSBarry Smith 8841e07b27eSBarry Smith PetscFunctionBegin; 8851e07b27eSBarry Smith ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 8861e07b27eSBarry Smith sred->redfactor = 1; 8871e07b27eSBarry Smith sred->ignore_dm = PETSC_FALSE; 8887c5279cbSDave May sred->ignore_kspcomputeoperators = PETSC_FALSE; 8891e07b27eSBarry Smith pc->data = (void*)sred; 8901e07b27eSBarry Smith 8911e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope; 8921e07b27eSBarry Smith pc->ops->applytranspose = NULL; 8931e07b27eSBarry Smith pc->ops->setup = PCSetUp_Telescope; 8941e07b27eSBarry Smith pc->ops->destroy = PCDestroy_Telescope; 8951e07b27eSBarry Smith pc->ops->reset = PCReset_Telescope; 8961e07b27eSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Telescope; 8971e07b27eSBarry Smith pc->ops->view = PCView_Telescope; 8981e07b27eSBarry Smith 8991e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 9001e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 9011e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 9021e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 9031e07b27eSBarry Smith 9041e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 9051e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 9061e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 9071e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 9081e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 909*0ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 910*0ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 9111e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 9121e07b27eSBarry Smith PetscFunctionReturn(0); 9131e07b27eSBarry Smith } 914