11e07b27eSBarry Smith 2*8d9f7141SDave May #include <petsc/private/petscimpl.h> 3120bdd93SDave May #include <petsc/private/matimpl.h> 46fc41876SBarry Smith #include <petsc/private/pcimpl.h> 51e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/ 61e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/ 7575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h" 81e07b27eSBarry Smith 9bf00f589SPatrick Sanan static PetscBool cited = PETSC_FALSE; 10bf00f589SPatrick Sanan static const char citation[] = 11bf00f589SPatrick Sanan "@inproceedings{MaySananRuppKnepleySmith2016,\n" 12bf00f589SPatrick Sanan " title = {Extreme-Scale Multigrid Components within PETSc},\n" 13bf00f589SPatrick Sanan " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 14bf00f589SPatrick Sanan " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 15bf00f589SPatrick Sanan " series = {PASC '16},\n" 16bf00f589SPatrick Sanan " isbn = {978-1-4503-4126-4},\n" 17bf00f589SPatrick Sanan " location = {Lausanne, Switzerland},\n" 18bf00f589SPatrick Sanan " pages = {5:1--5:12},\n" 19bf00f589SPatrick Sanan " articleno = {5},\n" 20bf00f589SPatrick Sanan " numpages = {12},\n" 21bf00f589SPatrick Sanan " url = {http://doi.acm.org/10.1145/2929908.2929913},\n" 22bf00f589SPatrick Sanan " doi = {10.1145/2929908.2929913},\n" 23bf00f589SPatrick Sanan " acmid = {2929913},\n" 24bf00f589SPatrick Sanan " publisher = {ACM},\n" 25bf00f589SPatrick Sanan " address = {New York, NY, USA},\n" 26bf00f589SPatrick Sanan " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 27bf00f589SPatrick Sanan " year = {2016}\n" 28bf00f589SPatrick Sanan "}\n"; 29bf00f589SPatrick Sanan 301e07b27eSBarry Smith /* 311e07b27eSBarry Smith PCTelescopeSetUp_default() 321e07b27eSBarry Smith PCTelescopeMatCreate_default() 331e07b27eSBarry Smith 341e07b27eSBarry Smith default 351e07b27eSBarry Smith 361e07b27eSBarry Smith // scatter in 371e07b27eSBarry Smith x(comm) -> xtmp(comm) 381e07b27eSBarry Smith 391e07b27eSBarry Smith xred(subcomm) <- xtmp 401e07b27eSBarry Smith yred(subcomm) 411e07b27eSBarry Smith 421e07b27eSBarry Smith yred(subcomm) --> xtmp 431e07b27eSBarry Smith 441e07b27eSBarry Smith // scatter out 451e07b27eSBarry Smith xtmp(comm) -> y(comm) 461e07b27eSBarry Smith */ 471e07b27eSBarry Smith 482a22aa42SDave May PetscBool PetscSubComm_isActiveRank(PetscSubcomm scomm) 491e07b27eSBarry Smith { 501e07b27eSBarry Smith if (scomm->color == 0) { return PETSC_TRUE; } 511e07b27eSBarry Smith else { return PETSC_FALSE; } 521e07b27eSBarry Smith } 531e07b27eSBarry Smith 542a22aa42SDave May PetscBool isActiveRank(PC_Telescope sred) 552a22aa42SDave May { 562a22aa42SDave May if (sred->psubcomm) { 572a22aa42SDave May return(PetscSubComm_isActiveRank(sred->psubcomm)); 582a22aa42SDave May } else { 592a22aa42SDave May if (sred->subcomm != MPI_COMM_NULL) { return PETSC_TRUE; } 602a22aa42SDave May else { return PETSC_FALSE; } 612a22aa42SDave May } 622a22aa42SDave May } 632a22aa42SDave May 64*8d9f7141SDave May /* 65*8d9f7141SDave May Collective on MPI_Comm[comm_f] 66*8d9f7141SDave May Notes 67*8d9f7141SDave May * Using comm_f = MPI_COMM_NULL will result in an error 68*8d9f7141SDave May * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid. 69*8d9f7141SDave May * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid. 70*8d9f7141SDave May */ 71*8d9f7141SDave May PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool *isvalid) 72*8d9f7141SDave May { 73*8d9f7141SDave May int valid = 1; 74*8d9f7141SDave May MPI_Group group_f,group_c; 75*8d9f7141SDave May PetscErrorCode ierr; 76*8d9f7141SDave May int errorcode; 77*8d9f7141SDave May PetscMPIInt count,k,size_f = 0,size_c = 0,size_c_sum = 0; 78*8d9f7141SDave May int *ranks_f = NULL,*ranks_c = NULL; 79*8d9f7141SDave May 80*8d9f7141SDave May if (comm_f == MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"comm_f cannot be MPI_COMM_NULL"); 81*8d9f7141SDave May 82*8d9f7141SDave May ierr = MPI_Comm_group(comm_f,&group_f);CHKERRQ(ierr); 83*8d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 84*8d9f7141SDave May ierr = MPI_Comm_group(comm_c,&group_c);CHKERRQ(ierr); 85*8d9f7141SDave May } 86*8d9f7141SDave May 87*8d9f7141SDave May ierr = MPI_Comm_size(comm_f,&size_f);CHKERRQ(ierr); 88*8d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 89*8d9f7141SDave May ierr = MPI_Comm_size(comm_c,&size_c);CHKERRQ(ierr); 90*8d9f7141SDave May } 91*8d9f7141SDave May 92*8d9f7141SDave May /* check not all comm_c's are NULL */ 93*8d9f7141SDave May size_c_sum = size_c; 94*8d9f7141SDave May ierr = MPI_Allreduce(MPI_IN_PLACE,&size_c_sum,1,MPI_INT,MPI_SUM,comm_f);CHKERRQ(ierr); 95*8d9f7141SDave May if (size_c_sum == 0) { 96*8d9f7141SDave May valid = 0; 97*8d9f7141SDave May } 98*8d9f7141SDave May 99*8d9f7141SDave May /* check we can map at least 1 rank in comm_c to comm_f */ 100*8d9f7141SDave May ierr = PetscMalloc1(size_f,&ranks_f);CHKERRQ(ierr); 101*8d9f7141SDave May ierr = PetscMalloc1(size_c,&ranks_c);CHKERRQ(ierr); 102*8d9f7141SDave May for (k=0; k<size_f; k++) { 103*8d9f7141SDave May ranks_f[k] = MPI_UNDEFINED; 104*8d9f7141SDave May } 105*8d9f7141SDave May for (k=0; k<size_c; k++) { 106*8d9f7141SDave May ranks_c[k] = (int)k; 107*8d9f7141SDave May } 108*8d9f7141SDave May 109*8d9f7141SDave May count = 0; 110*8d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 111*8d9f7141SDave May errorcode = MPI_Group_translate_ranks(group_c,size_c,ranks_c,group_f,ranks_f); 112*8d9f7141SDave May for (k=0; k<size_f; k++) { 113*8d9f7141SDave May if (ranks_f[k] == MPI_UNDEFINED) { 114*8d9f7141SDave May count++; 115*8d9f7141SDave May } 116*8d9f7141SDave May } 117*8d9f7141SDave May } 118*8d9f7141SDave May if (count == size_f) { 119*8d9f7141SDave May valid = 0; 120*8d9f7141SDave May } 121*8d9f7141SDave May 122*8d9f7141SDave May ierr = MPI_Allreduce(MPI_IN_PLACE,&valid,1,MPI_INT,MPI_MIN,comm_f);CHKERRQ(ierr); 123*8d9f7141SDave May if (valid == 1) { *isvalid = PETSC_TRUE; } 124*8d9f7141SDave May else { *isvalid = PETSC_FALSE; } 125*8d9f7141SDave May 126*8d9f7141SDave May ierr = PetscFree(ranks_f);CHKERRQ(ierr); 127*8d9f7141SDave May ierr = PetscFree(ranks_c);CHKERRQ(ierr); 128*8d9f7141SDave May ierr = MPI_Group_free(&group_f);CHKERRQ(ierr); 129*8d9f7141SDave May if (comm_c != MPI_COMM_NULL) { 130*8d9f7141SDave May ierr = MPI_Group_free(&group_c);CHKERRQ(ierr); 131*8d9f7141SDave May } 132*8d9f7141SDave May PetscFunctionReturn(0); 133*8d9f7141SDave May } 134*8d9f7141SDave May 1351e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred) 1361e07b27eSBarry Smith { 137c6a0d831SBarry Smith DM subdm = NULL; 1381e07b27eSBarry Smith 1392a22aa42SDave May if (!isActiveRank(sred)) { subdm = NULL; } 1401e07b27eSBarry Smith else { 1411e07b27eSBarry Smith switch (sred->sr_type) { 1421e07b27eSBarry Smith case TELESCOPE_DEFAULT: subdm = NULL; 1431e07b27eSBarry Smith break; 1441e07b27eSBarry Smith case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 1451e07b27eSBarry Smith break; 1461e07b27eSBarry Smith case TELESCOPE_DMPLEX: subdm = NULL; 1471e07b27eSBarry Smith break; 148*8d9f7141SDave May case TELESCOPE_COARSEDM: if (sred->ksp) { KSPGetDM(sred->ksp,&subdm); } 149*8d9f7141SDave May break; 1501e07b27eSBarry Smith } 1511e07b27eSBarry Smith } 1521e07b27eSBarry Smith return(subdm); 1531e07b27eSBarry Smith } 1541e07b27eSBarry Smith 1551e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 1561e07b27eSBarry Smith { 1571e07b27eSBarry Smith PetscErrorCode ierr; 1581e07b27eSBarry Smith PetscInt m,M,bs,st,ed; 1591e07b27eSBarry Smith Vec x,xred,yred,xtmp; 1601e07b27eSBarry Smith Mat B; 1611e07b27eSBarry Smith MPI_Comm comm,subcomm; 1621e07b27eSBarry Smith VecScatter scatter; 1631e07b27eSBarry Smith IS isin; 1641e07b27eSBarry Smith 1651e07b27eSBarry Smith PetscFunctionBegin; 1661e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr); 1671e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 1681e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 1691e07b27eSBarry Smith 1701e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 1711e07b27eSBarry Smith ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 1721e07b27eSBarry Smith ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 1731e07b27eSBarry Smith ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 1741e07b27eSBarry Smith 1751e07b27eSBarry Smith xred = NULL; 1763ac26c5eSBarry Smith m = 0; 1772a22aa42SDave May if (isActiveRank(sred)) { 1781e07b27eSBarry Smith ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr); 1791e07b27eSBarry Smith ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr); 1801e07b27eSBarry Smith ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr); 1811e07b27eSBarry Smith ierr = VecSetFromOptions(xred);CHKERRQ(ierr); 182ca43db0aSBarry Smith ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr); 1831e07b27eSBarry Smith } 1841e07b27eSBarry Smith 1851e07b27eSBarry Smith yred = NULL; 1862a22aa42SDave May if (isActiveRank(sred)) { 1871e07b27eSBarry Smith ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 1881e07b27eSBarry Smith } 1891e07b27eSBarry Smith 1901e07b27eSBarry Smith ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 1911e07b27eSBarry Smith ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 1921e07b27eSBarry Smith ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 1931e07b27eSBarry Smith ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 1941e07b27eSBarry Smith 1952a22aa42SDave May if (isActiveRank(sred)) { 1961e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 1971e07b27eSBarry Smith ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr); 1981e07b27eSBarry Smith } else { 1991e07b27eSBarry Smith ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 2003ac26c5eSBarry Smith ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 2011e07b27eSBarry Smith } 2021e07b27eSBarry Smith ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 2031e07b27eSBarry Smith 20435928de7SBarry Smith ierr = VecScatterCreateWithData(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 2051e07b27eSBarry Smith 2061e07b27eSBarry Smith sred->isin = isin; 2071e07b27eSBarry Smith sred->scatter = scatter; 2081e07b27eSBarry Smith sred->xred = xred; 2091e07b27eSBarry Smith sred->yred = yred; 2101e07b27eSBarry Smith sred->xtmp = xtmp; 2111e07b27eSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 2121e07b27eSBarry Smith PetscFunctionReturn(0); 2131e07b27eSBarry Smith } 2141e07b27eSBarry Smith 2151e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 2161e07b27eSBarry Smith { 2171e07b27eSBarry Smith PetscErrorCode ierr; 2181e07b27eSBarry Smith MPI_Comm comm,subcomm; 2191e07b27eSBarry Smith Mat Bred,B; 2201e07b27eSBarry Smith PetscInt nr,nc; 2211e07b27eSBarry Smith IS isrow,iscol; 2221e07b27eSBarry Smith Mat Blocal,*_Blocal; 2231e07b27eSBarry Smith 2241e07b27eSBarry Smith PetscFunctionBegin; 2251e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr); 2261e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 2271e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2281e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 2291e07b27eSBarry Smith ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 2301e07b27eSBarry Smith isrow = sred->isin; 2311e07b27eSBarry Smith ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 2327dae84e0SHong Zhang ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 2331e07b27eSBarry Smith Blocal = *_Blocal; 2341e07b27eSBarry Smith ierr = PetscFree(_Blocal);CHKERRQ(ierr); 2351e07b27eSBarry Smith Bred = NULL; 2362a22aa42SDave May if (isActiveRank(sred)) { 2371e07b27eSBarry Smith PetscInt mm; 2381e07b27eSBarry Smith 2391e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 2401e07b27eSBarry Smith 2411e07b27eSBarry Smith ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 2421e07b27eSBarry Smith ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 2431e07b27eSBarry Smith } 2441e07b27eSBarry Smith *A = Bred; 2451e07b27eSBarry Smith ierr = ISDestroy(&iscol);CHKERRQ(ierr); 2461e07b27eSBarry Smith ierr = MatDestroy(&Blocal);CHKERRQ(ierr); 2471e07b27eSBarry Smith PetscFunctionReturn(0); 2481e07b27eSBarry Smith } 2491e07b27eSBarry Smith 250392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 2511e07b27eSBarry Smith { 2521e07b27eSBarry Smith PetscErrorCode ierr; 2531e07b27eSBarry Smith PetscBool has_const; 2541e07b27eSBarry Smith const Vec *vecs; 255c41e779fSDave May Vec *sub_vecs = NULL; 256392968a1SPatrick Sanan PetscInt i,k,n = 0; 2571e07b27eSBarry Smith MPI_Comm subcomm; 2581e07b27eSBarry Smith 2591e07b27eSBarry Smith PetscFunctionBegin; 2601e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 2611e07b27eSBarry Smith ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 2621e07b27eSBarry Smith 2632a22aa42SDave May if (isActiveRank(sred)) { 264e3acf2f7SBarry Smith if (n) { 265e3acf2f7SBarry Smith ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 2661e07b27eSBarry Smith } 2671e07b27eSBarry Smith } 2681e07b27eSBarry Smith 2691e07b27eSBarry Smith /* copy entries */ 2701e07b27eSBarry Smith for (k=0; k<n; k++) { 2711e07b27eSBarry Smith const PetscScalar *x_array; 2721e07b27eSBarry Smith PetscScalar *LA_sub_vec; 27313c30530SDave May PetscInt st,ed; 2741e07b27eSBarry Smith 2751e07b27eSBarry Smith /* pull in vector x->xtmp */ 2761e07b27eSBarry Smith ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2771e07b27eSBarry Smith ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 27847856c66SBarry Smith if (sub_vecs) { 279a04a6428SPatrick Sanan /* copy vector entries into xred */ 2801e07b27eSBarry Smith ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 281ea2b237eSDave May if (sub_vecs[k]) { 2821e07b27eSBarry Smith ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 2831e07b27eSBarry Smith ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 2841e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 2851e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 2861e07b27eSBarry Smith } 2871e07b27eSBarry Smith ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 2881e07b27eSBarry Smith } 2891e07b27eSBarry Smith ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 2901e07b27eSBarry Smith } 29147856c66SBarry Smith } 2921e07b27eSBarry Smith 2932a22aa42SDave May if (isActiveRank(sred)) { 294d8b9d5b7SPatrick Sanan /* create new (near) nullspace for redundant object */ 295392968a1SPatrick Sanan ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr); 296392968a1SPatrick Sanan ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 297d8b9d5b7SPatrick Sanan if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 298d8b9d5b7SPatrick 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"); 299d8b9d5b7SPatrick Sanan } 300392968a1SPatrick Sanan PetscFunctionReturn(0); 301392968a1SPatrick Sanan } 302392968a1SPatrick Sanan 303392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 304392968a1SPatrick Sanan { 305392968a1SPatrick Sanan PetscErrorCode ierr; 306392968a1SPatrick Sanan Mat B; 307392968a1SPatrick Sanan 308392968a1SPatrick Sanan PetscFunctionBegin; 309392968a1SPatrick Sanan ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 310392968a1SPatrick Sanan 311392968a1SPatrick Sanan /* Propagate the nullspace if it exists */ 312392968a1SPatrick Sanan { 313392968a1SPatrick Sanan MatNullSpace nullspace,sub_nullspace; 314392968a1SPatrick Sanan ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 315392968a1SPatrick Sanan if (nullspace) { 316392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr); 317392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr); 3182a22aa42SDave May if (isActiveRank(sred)) { 319392968a1SPatrick Sanan ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 32041ff1ee9SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr); 3211e07b27eSBarry Smith } 322392968a1SPatrick Sanan } 323392968a1SPatrick Sanan } 324392968a1SPatrick Sanan 325392968a1SPatrick Sanan /* Propagate the near nullspace if it exists */ 326392968a1SPatrick Sanan { 327392968a1SPatrick Sanan MatNullSpace nearnullspace,sub_nearnullspace; 328392968a1SPatrick Sanan ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr); 329392968a1SPatrick Sanan if (nearnullspace) { 330392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr); 331392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr); 3322a22aa42SDave May if (isActiveRank(sred)) { 333392968a1SPatrick Sanan ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr); 334392968a1SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr); 335392968a1SPatrick Sanan } 336392968a1SPatrick Sanan } 337392968a1SPatrick Sanan } 3381e07b27eSBarry Smith PetscFunctionReturn(0); 3391e07b27eSBarry Smith } 3401e07b27eSBarry Smith 3411e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 3421e07b27eSBarry Smith { 3431e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 3441e07b27eSBarry Smith PetscErrorCode ierr; 3451e07b27eSBarry Smith PetscBool iascii,isstring; 3461e07b27eSBarry Smith PetscViewer subviewer; 3471e07b27eSBarry Smith 3481e07b27eSBarry Smith PetscFunctionBegin; 3491e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3501e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 3511e07b27eSBarry Smith if (iascii) { 352*8d9f7141SDave May { 3531e07b27eSBarry Smith MPI_Comm comm,subcomm; 3541e07b27eSBarry Smith PetscMPIInt comm_size,subcomm_size; 355*8d9f7141SDave May DM dm = NULL,subdm = NULL; 3561e07b27eSBarry Smith 3571e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3581e07b27eSBarry Smith subdm = private_PCTelescopeGetSubDM(sred); 359*8d9f7141SDave May 360*8d9f7141SDave May if (sred->psubcomm) { 3611e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 3621e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 3631e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 3641e07b27eSBarry Smith ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 3651e07b27eSBarry Smith 366*8d9f7141SDave May ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 367*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 368*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 36948a10b22SPatrick Sanan switch (sred->subcommtype) { 37048a10b22SPatrick Sanan case PETSC_SUBCOMM_INTERLACED : 371*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: type = interlaced\n",sred->subcommtype);CHKERRQ(ierr); 37248a10b22SPatrick Sanan break; 37348a10b22SPatrick Sanan case PETSC_SUBCOMM_CONTIGUOUS : 374*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm type = contiguous\n",sred->subcommtype);CHKERRQ(ierr); 37548a10b22SPatrick Sanan break; 37648a10b22SPatrick Sanan default : 37748a10b22SPatrick Sanan SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope"); 37848a10b22SPatrick Sanan } 379*8d9f7141SDave May ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 380*8d9f7141SDave May } else { 381*8d9f7141SDave May ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 382*8d9f7141SDave May subcomm = sred->subcomm; 383*8d9f7141SDave May if (!isActiveRank(sred)) { 384*8d9f7141SDave May subcomm = PETSC_COMM_SELF; 385*8d9f7141SDave May } 386*8d9f7141SDave May 387*8d9f7141SDave May ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 388*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(viewer,"subcomm: using user provided sub-communicator\n");CHKERRQ(ierr); 389*8d9f7141SDave May ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 390*8d9f7141SDave May } 391*8d9f7141SDave May 3921e07b27eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 3932a22aa42SDave May if (isActiveRank(sred)) { 3941e07b27eSBarry Smith ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 3951e07b27eSBarry Smith 3961e07b27eSBarry Smith if (dm && sred->ignore_dm) { 397efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer,"ignoring DM\n");CHKERRQ(ierr); 3981e07b27eSBarry Smith } 3997c5279cbSDave May if (sred->ignore_kspcomputeoperators) { 400efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(subviewer,"ignoring KSPComputeOperators\n");CHKERRQ(ierr); 4017c5279cbSDave May } 4021e07b27eSBarry Smith switch (sred->sr_type) { 4031e07b27eSBarry Smith case TELESCOPE_DEFAULT: 404*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"setup type: default\n");CHKERRQ(ierr); 4051e07b27eSBarry Smith break; 4061e07b27eSBarry Smith case TELESCOPE_DMDA: 407*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMDA auto-repartitioning\n");CHKERRQ(ierr); 4088ef9ca65SPatrick Sanan ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr); 4091e07b27eSBarry Smith break; 4101e07b27eSBarry Smith case TELESCOPE_DMPLEX: 411*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMPLEX auto-repartitioning\n");CHKERRQ(ierr); 4121e07b27eSBarry Smith break; 413*8d9f7141SDave May case TELESCOPE_COARSEDM: 414*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"setup type: coarse DM\n");CHKERRQ(ierr); 415*8d9f7141SDave May break; 416*8d9f7141SDave May } 417*8d9f7141SDave May 418*8d9f7141SDave May if (dm) { 419*8d9f7141SDave May PetscObject obj = (PetscObject)dm; 420*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object:");CHKERRQ(ierr); 421*8d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 422*8d9f7141SDave May if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 423*8d9f7141SDave May if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 424*8d9f7141SDave May if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 425*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr); 426*8d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 427*8d9f7141SDave May } else { 428*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object: NULL\n");CHKERRQ(ierr); 429*8d9f7141SDave May } 430*8d9f7141SDave May if (subdm) { 431*8d9f7141SDave May PetscObject obj = (PetscObject)subdm; 432*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object:");CHKERRQ(ierr); 433*8d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 434*8d9f7141SDave May if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 435*8d9f7141SDave May if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 436*8d9f7141SDave May if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 437*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr); 438*8d9f7141SDave May PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 439*8d9f7141SDave May } else { 440*8d9f7141SDave May ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object: NULL\n");CHKERRQ(ierr); 4411e07b27eSBarry Smith } 4421e07b27eSBarry Smith 4431e07b27eSBarry Smith ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 4441e07b27eSBarry Smith ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 4451e07b27eSBarry Smith } 4461e07b27eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 4471e07b27eSBarry Smith } 4481e07b27eSBarry Smith } 4491e07b27eSBarry Smith PetscFunctionReturn(0); 4501e07b27eSBarry Smith } 4511e07b27eSBarry Smith 4521e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc) 4531e07b27eSBarry Smith { 4541e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 4551e07b27eSBarry Smith PetscErrorCode ierr; 456bd49479cSSatish Balay MPI_Comm comm,subcomm=0; 4571e07b27eSBarry Smith PCTelescopeType sr_type; 4581e07b27eSBarry Smith 4591e07b27eSBarry Smith PetscFunctionBegin; 4601e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 4611e07b27eSBarry Smith 4621e07b27eSBarry Smith /* Determine type of setup/update */ 4631e07b27eSBarry Smith if (!pc->setupcalled) { 4641e07b27eSBarry Smith PetscBool has_dm,same; 4651e07b27eSBarry Smith DM dm; 4661e07b27eSBarry Smith 4671e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4681e07b27eSBarry Smith has_dm = PETSC_FALSE; 4691e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 4701e07b27eSBarry Smith if (dm) { has_dm = PETSC_TRUE; } 4711e07b27eSBarry Smith if (has_dm) { 4721e07b27eSBarry Smith /* check for dmda */ 4731e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 4741e07b27eSBarry Smith if (same) { 4751e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 4761e07b27eSBarry Smith sr_type = TELESCOPE_DMDA; 4771e07b27eSBarry Smith } 4781e07b27eSBarry Smith /* check for dmplex */ 4791e07b27eSBarry Smith ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 4801e07b27eSBarry Smith if (same) { 481994fe344SLisandro Dalcin ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr); 4821e07b27eSBarry Smith sr_type = TELESCOPE_DMPLEX; 4831e07b27eSBarry Smith } 484*8d9f7141SDave May 485*8d9f7141SDave May if (sred->use_coarse_dm) { 486*8d9f7141SDave May ierr = PetscInfo(pc,"PCTelescope: using coarse DM\n");CHKERRQ(ierr); 487*8d9f7141SDave May sr_type = TELESCOPE_COARSEDM; 4881e07b27eSBarry Smith } 4891e07b27eSBarry Smith 4901e07b27eSBarry Smith if (sred->ignore_dm) { 491*8d9f7141SDave May ierr = PetscInfo(pc,"PCTelescope: ignoring DM\n");CHKERRQ(ierr); 4921e07b27eSBarry Smith sr_type = TELESCOPE_DEFAULT; 4931e07b27eSBarry Smith } 494*8d9f7141SDave May } 4951e07b27eSBarry Smith sred->sr_type = sr_type; 4961e07b27eSBarry Smith } else { 4971e07b27eSBarry Smith sr_type = sred->sr_type; 4981e07b27eSBarry Smith } 4991e07b27eSBarry Smith 500d8b9d5b7SPatrick Sanan /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */ 5011e07b27eSBarry Smith switch (sr_type) { 5021e07b27eSBarry Smith case TELESCOPE_DEFAULT: 5031e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 5041e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 5051e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 5061e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 5071e07b27eSBarry Smith break; 5081e07b27eSBarry Smith case TELESCOPE_DMDA: 5091e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope_dmda; 510f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_dmda; 5111e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 5121e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 5131e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 5141e07b27eSBarry Smith sred->pctelescope_reset_type = PCReset_Telescope_dmda; 5151e07b27eSBarry Smith break; 516d8b9d5b7SPatrick Sanan case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available"); 5171e07b27eSBarry Smith break; 518*8d9f7141SDave May case TELESCOPE_COARSEDM: 519*8d9f7141SDave May pc->ops->apply = PCApply_Telescope_CoarseDM; 520*8d9f7141SDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope_CoarseDM; 521*8d9f7141SDave May sred->pctelescope_setup_type = PCTelescopeSetUp_CoarseDM; 522*8d9f7141SDave May sred->pctelescope_matcreate_type = NULL; 523*8d9f7141SDave May sred->pctelescope_matnullspacecreate_type = NULL;/*PCTelescopeMatNullSpaceCreate_CoarseDM;*/ 524*8d9f7141SDave May sred->pctelescope_reset_type = PCReset_Telescope_CoarseDM; 5251e07b27eSBarry Smith break; 526*8d9f7141SDave May default: SETERRQ(comm,PETSC_ERR_SUP,"Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM"); 527*8d9f7141SDave May break; 528*8d9f7141SDave May } 529*8d9f7141SDave May 530*8d9f7141SDave May /* subcomm definition */ 531*8d9f7141SDave May if (!pc->setupcalled) { 532*8d9f7141SDave May if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) { 533*8d9f7141SDave May if (!sred->psubcomm) { 534*8d9f7141SDave May ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 535*8d9f7141SDave May ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 536*8d9f7141SDave May ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr); 537*8d9f7141SDave May ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 538*8d9f7141SDave May sred->subcomm = PetscSubcommChild(sred->psubcomm); 539*8d9f7141SDave May } 540*8d9f7141SDave May } else { /* query PC for DM, check communicators */ 541*8d9f7141SDave May DM dm,dm_coarse_partition = NULL; 542*8d9f7141SDave May MPI_Comm comm_fine,comm_coarse_partition = MPI_COMM_NULL; 543*8d9f7141SDave May PetscMPIInt csize_fine=0,csize_coarse_partition=0,cs[2],csg[2],cnt=0; 544*8d9f7141SDave May PetscBool isvalidsubcomm; 545*8d9f7141SDave May 546*8d9f7141SDave May ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 547*8d9f7141SDave May comm_fine = PetscObjectComm((PetscObject)dm); 548*8d9f7141SDave May ierr = DMGetCoarseDM(dm,&dm_coarse_partition);CHKERRQ(ierr); 549*8d9f7141SDave May if (dm_coarse_partition) { cnt = 1; } 550*8d9f7141SDave May ierr = MPI_Allreduce(MPI_IN_PLACE,&cnt,1,MPI_INT,MPI_SUM,comm_fine);CHKERRQ(ierr); 551*8d9f7141SDave May if (cnt == 0) SETERRQ(comm_fine,PETSC_ERR_SUP,"Zero instances of a coarse DM were found"); 552*8d9f7141SDave May 553*8d9f7141SDave May ierr = MPI_Comm_size(comm_fine,&csize_fine);CHKERRQ(ierr); 554*8d9f7141SDave May if (dm_coarse_partition) { 555*8d9f7141SDave May comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition); 556*8d9f7141SDave May ierr = MPI_Comm_size(comm_coarse_partition,&csize_coarse_partition);CHKERRQ(ierr); 557*8d9f7141SDave May } 558*8d9f7141SDave May 559*8d9f7141SDave May cs[0] = csize_fine; 560*8d9f7141SDave May cs[1] = csize_coarse_partition; 561*8d9f7141SDave May ierr = MPI_Allreduce(cs,csg,2,MPI_INT,MPI_MAX,comm_fine);CHKERRQ(ierr); 562*8d9f7141SDave May if (csg[0] == csg[1]) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM uses the same size communicator as the parent DM attached to the PC"); 563*8d9f7141SDave May 564*8d9f7141SDave May ierr = PCTelescopeTestValidSubcomm(comm_fine,comm_coarse_partition,&isvalidsubcomm);CHKERRQ(ierr); 565*8d9f7141SDave May if (!isvalidsubcomm) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM communicator is not a sub-communicator of parentDM->comm"); 566*8d9f7141SDave May sred->subcomm = comm_coarse_partition; 567*8d9f7141SDave May } 568*8d9f7141SDave May } 569*8d9f7141SDave May subcomm = sred->subcomm; 570*8d9f7141SDave May 571*8d9f7141SDave May /* internal KSP */ 572*8d9f7141SDave May if (!pc->setupcalled) { 573*8d9f7141SDave May const char *prefix; 574*8d9f7141SDave May 575*8d9f7141SDave May if (isActiveRank(sred)) { 576*8d9f7141SDave May ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 577*8d9f7141SDave May ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 578*8d9f7141SDave May ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 579*8d9f7141SDave May ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 580*8d9f7141SDave May ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 581*8d9f7141SDave May ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 582*8d9f7141SDave May ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 583*8d9f7141SDave May ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 584*8d9f7141SDave May } 5851e07b27eSBarry Smith } 5861e07b27eSBarry Smith 5871e07b27eSBarry Smith /* setup */ 5881e07b27eSBarry Smith if (sred->pctelescope_setup_type) { 5891e07b27eSBarry Smith ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 5901e07b27eSBarry Smith } 5911e07b27eSBarry Smith /* update */ 5921e07b27eSBarry Smith if (!pc->setupcalled) { 5931e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 5941e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 5951e07b27eSBarry Smith } 5961e07b27eSBarry Smith if (sred->pctelescope_matnullspacecreate_type) { 597392968a1SPatrick Sanan ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 5981e07b27eSBarry Smith } 5991e07b27eSBarry Smith } else { 6001e07b27eSBarry Smith if (sred->pctelescope_matcreate_type) { 6011e07b27eSBarry Smith ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 6021e07b27eSBarry Smith } 6031e07b27eSBarry Smith } 6041e07b27eSBarry Smith 6051e07b27eSBarry Smith /* common - no construction */ 6062a22aa42SDave May if (isActiveRank(sred)) { 6071e07b27eSBarry Smith ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 6081e07b27eSBarry Smith if (pc->setfromoptionscalled && !pc->setupcalled){ 6091e07b27eSBarry Smith ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 6101e07b27eSBarry Smith } 6111e07b27eSBarry Smith } 612*8d9f7141SDave May 613*8d9f7141SDave May #if 0 614*8d9f7141SDave May /* we perform this last as Bred is not available with KSPSetComputeOperators() until KSPSetUp has been called */ 615*8d9f7141SDave May if (!pc->setupcalled) { 616*8d9f7141SDave May if (isActiveRank(sred)) { 617*8d9f7141SDave May ierr = KSPSetUp(sred->ksp);CHKERRQ(ierr); 618*8d9f7141SDave May } 619*8d9f7141SDave May if (sred->pctelescope_matnullspacecreate_type) { 620*8d9f7141SDave May ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 621*8d9f7141SDave May } 622*8d9f7141SDave May } 623*8d9f7141SDave May #endif 624*8d9f7141SDave May 6251e07b27eSBarry Smith PetscFunctionReturn(0); 6261e07b27eSBarry Smith } 6271e07b27eSBarry Smith 6281e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 6291e07b27eSBarry Smith { 6301e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 6311e07b27eSBarry Smith PetscErrorCode ierr; 6321e07b27eSBarry Smith Vec xtmp,xred,yred; 63313c30530SDave May PetscInt i,st,ed; 6341e07b27eSBarry Smith VecScatter scatter; 6351e07b27eSBarry Smith PetscScalar *array; 6361e07b27eSBarry Smith const PetscScalar *x_array; 6371e07b27eSBarry Smith 6381e07b27eSBarry Smith PetscFunctionBegin; 639bf00f589SPatrick Sanan ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 640bf00f589SPatrick Sanan 6411e07b27eSBarry Smith xtmp = sred->xtmp; 6421e07b27eSBarry Smith scatter = sred->scatter; 6431e07b27eSBarry Smith xred = sred->xred; 6441e07b27eSBarry Smith yred = sred->yred; 6451e07b27eSBarry Smith 6461e07b27eSBarry Smith /* pull in vector x->xtmp */ 6471e07b27eSBarry Smith ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6481e07b27eSBarry Smith ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6491e07b27eSBarry Smith 650bf00f589SPatrick Sanan /* copy vector entries into xred */ 6511e07b27eSBarry Smith ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 6521e07b27eSBarry Smith if (xred) { 6531e07b27eSBarry Smith PetscScalar *LA_xred; 6541e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 6551e07b27eSBarry Smith ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 6561e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 6571e07b27eSBarry Smith LA_xred[i] = x_array[i]; 6581e07b27eSBarry Smith } 6591e07b27eSBarry Smith ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 6601e07b27eSBarry Smith } 6611e07b27eSBarry Smith ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 6621e07b27eSBarry Smith /* solve */ 6632a22aa42SDave May if (isActiveRank(sred)) { 6641e07b27eSBarry Smith ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 665c0decd05SBarry Smith ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr); 6661e07b27eSBarry Smith } 6671e07b27eSBarry Smith /* return vector */ 6681e07b27eSBarry Smith ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 6691e07b27eSBarry Smith if (yred) { 6701e07b27eSBarry Smith const PetscScalar *LA_yred; 6711e07b27eSBarry Smith ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 6721e07b27eSBarry Smith ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 6731e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 6741e07b27eSBarry Smith array[i] = LA_yred[i]; 6751e07b27eSBarry Smith } 6761e07b27eSBarry Smith ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 6771e07b27eSBarry Smith } 6781e07b27eSBarry Smith ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 6791e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6801e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6811e07b27eSBarry Smith PetscFunctionReturn(0); 6821e07b27eSBarry Smith } 6831e07b27eSBarry Smith 684f650675bSDave 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) 685f650675bSDave May { 686f650675bSDave May PC_Telescope sred = (PC_Telescope)pc->data; 687f650675bSDave May PetscErrorCode ierr; 688a1d91a28SDave May Vec xtmp,yred; 689f650675bSDave May PetscInt i,st,ed; 690f650675bSDave May VecScatter scatter; 691f650675bSDave May const PetscScalar *x_array; 692f650675bSDave May PetscBool default_init_guess_value; 693f650675bSDave May 694f650675bSDave May PetscFunctionBegin; 695f650675bSDave May xtmp = sred->xtmp; 696f650675bSDave May scatter = sred->scatter; 697f650675bSDave May yred = sred->yred; 698f650675bSDave May 699f650675bSDave May if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1"); 700f650675bSDave May *reason = (PCRichardsonConvergedReason)0; 701f650675bSDave May 702f650675bSDave May if (!zeroguess) { 703f650675bSDave May ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr); 704f650675bSDave May /* pull in vector y->xtmp */ 705f650675bSDave May ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 706f650675bSDave May ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 707f650675bSDave May 708bf00f589SPatrick Sanan /* copy vector entries into xred */ 709f650675bSDave May ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 710f650675bSDave May if (yred) { 711f650675bSDave May PetscScalar *LA_yred; 712f650675bSDave May ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 713f650675bSDave May ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 714f650675bSDave May for (i=0; i<ed-st; i++) { 715f650675bSDave May LA_yred[i] = x_array[i]; 716f650675bSDave May } 717f650675bSDave May ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 718f650675bSDave May } 719f650675bSDave May ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 720f650675bSDave May } 721f650675bSDave May 7222a22aa42SDave May if (isActiveRank(sred)) { 723f650675bSDave May ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr); 724f650675bSDave May if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 725f650675bSDave May } 726f650675bSDave May 727f650675bSDave May ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr); 728f650675bSDave May 7292a22aa42SDave May if (isActiveRank(sred)) { 730f650675bSDave May ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr); 731f650675bSDave May } 732f650675bSDave May 733f650675bSDave May if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 734f650675bSDave May *outits = 1; 735f650675bSDave May PetscFunctionReturn(0); 736f650675bSDave May } 737f650675bSDave May 7381e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc) 7391e07b27eSBarry Smith { 7401e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7411e07b27eSBarry Smith PetscErrorCode ierr; 7421e07b27eSBarry Smith 7431e07b27eSBarry Smith ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 7441e07b27eSBarry Smith ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 745e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); 746e3acf2f7SBarry Smith ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); 747e3acf2f7SBarry Smith ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); 748e3acf2f7SBarry Smith ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); 749e3acf2f7SBarry Smith ierr = KSPReset(sred->ksp);CHKERRQ(ierr); 7501e07b27eSBarry Smith if (sred->pctelescope_reset_type) { 7511e07b27eSBarry Smith ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 7521e07b27eSBarry Smith } 7531e07b27eSBarry Smith PetscFunctionReturn(0); 7541e07b27eSBarry Smith } 7551e07b27eSBarry Smith 7561e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc) 7571e07b27eSBarry Smith { 7581e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7591e07b27eSBarry Smith PetscErrorCode ierr; 7601e07b27eSBarry Smith 7611e07b27eSBarry Smith PetscFunctionBegin; 7621e07b27eSBarry Smith ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 763e3acf2f7SBarry Smith ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); 7641e07b27eSBarry Smith ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 765e3acf2f7SBarry Smith ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr); 766e3acf2f7SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 7671e07b27eSBarry Smith PetscFunctionReturn(0); 7681e07b27eSBarry Smith } 7691e07b27eSBarry Smith 7704416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc) 7711e07b27eSBarry Smith { 7721e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 7731e07b27eSBarry Smith PetscErrorCode ierr; 7741e07b27eSBarry Smith MPI_Comm comm; 7751e07b27eSBarry Smith PetscMPIInt size; 77648a10b22SPatrick Sanan PetscBool flg; 77748a10b22SPatrick Sanan PetscSubcommType subcommtype; 7781e07b27eSBarry Smith 7791e07b27eSBarry Smith PetscFunctionBegin; 7801e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7811e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 7821e07b27eSBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 78348a10b22SPatrick Sanan ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr); 78448a10b22SPatrick Sanan if (flg) { 78548a10b22SPatrick Sanan ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr); 78648a10b22SPatrick Sanan } 7871e07b27eSBarry Smith ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 7881e07b27eSBarry Smith if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 7891e07b27eSBarry Smith ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 7907c5279cbSDave May ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr); 791*8d9f7141SDave May ierr = PetscOptionsBool("-pc_telescope_use_coarse_dm","Define sub-communicator from the coarse DM","PCTelescopeSetUseCoarseDM",sred->use_coarse_dm,&sred->use_coarse_dm,0);CHKERRQ(ierr); 7921e07b27eSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7931e07b27eSBarry Smith PetscFunctionReturn(0); 7941e07b27eSBarry Smith } 7951e07b27eSBarry Smith 7961e07b27eSBarry Smith /* PC simplementation specific API's */ 7971e07b27eSBarry Smith 7981e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 7991e07b27eSBarry Smith { 8001e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 801bd49479cSSatish Balay PetscFunctionBegin; 8021e07b27eSBarry Smith if (ksp) *ksp = red->ksp; 803bd49479cSSatish Balay PetscFunctionReturn(0); 8041e07b27eSBarry Smith } 8051e07b27eSBarry Smith 80648a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype) 80748a10b22SPatrick Sanan { 80848a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 80948a10b22SPatrick Sanan PetscFunctionBegin; 81048a10b22SPatrick Sanan if (subcommtype) *subcommtype = red->subcommtype; 81148a10b22SPatrick Sanan PetscFunctionReturn(0); 81248a10b22SPatrick Sanan } 81348a10b22SPatrick Sanan 81448a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype) 81548a10b22SPatrick Sanan { 81648a10b22SPatrick Sanan PC_Telescope red = (PC_Telescope)pc->data; 81748a10b22SPatrick Sanan 81848a10b22SPatrick Sanan PetscFunctionBegin; 81948a10b22SPatrick 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."); 82048a10b22SPatrick Sanan red->subcommtype = subcommtype; 82148a10b22SPatrick Sanan PetscFunctionReturn(0); 82248a10b22SPatrick Sanan } 82348a10b22SPatrick Sanan 8241e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 8251e07b27eSBarry Smith { 8261e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 827bd49479cSSatish Balay PetscFunctionBegin; 8281e07b27eSBarry Smith if (fact) *fact = red->redfactor; 829bd49479cSSatish Balay PetscFunctionReturn(0); 8301e07b27eSBarry Smith } 8311e07b27eSBarry Smith 8321e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 8331e07b27eSBarry Smith { 8341e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 8351e07b27eSBarry Smith PetscMPIInt size; 8361e07b27eSBarry Smith PetscErrorCode ierr; 8371e07b27eSBarry Smith 838bd49479cSSatish Balay PetscFunctionBegin; 8391e07b27eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 8401e07b27eSBarry Smith if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 8411e07b27eSBarry Smith if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 8421e07b27eSBarry Smith red->redfactor = fact; 843bd49479cSSatish Balay PetscFunctionReturn(0); 8441e07b27eSBarry Smith } 8451e07b27eSBarry Smith 8461e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 8471e07b27eSBarry Smith { 8481e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 849bd49479cSSatish Balay PetscFunctionBegin; 8501e07b27eSBarry Smith if (v) *v = red->ignore_dm; 851bd49479cSSatish Balay PetscFunctionReturn(0); 8521e07b27eSBarry Smith } 85348a10b22SPatrick Sanan 8541e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 8551e07b27eSBarry Smith { 8561e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 857bd49479cSSatish Balay PetscFunctionBegin; 8581e07b27eSBarry Smith red->ignore_dm = v; 859bd49479cSSatish Balay PetscFunctionReturn(0); 8601e07b27eSBarry Smith } 8611e07b27eSBarry Smith 862*8d9f7141SDave May static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool *v) 863*8d9f7141SDave May { 864*8d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 865*8d9f7141SDave May PetscFunctionBegin; 866*8d9f7141SDave May if (v) *v = red->use_coarse_dm; 867*8d9f7141SDave May PetscFunctionReturn(0); 868*8d9f7141SDave May } 869*8d9f7141SDave May 870*8d9f7141SDave May static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v) 871*8d9f7141SDave May { 872*8d9f7141SDave May PC_Telescope red = (PC_Telescope)pc->data; 873*8d9f7141SDave May PetscFunctionBegin; 874*8d9f7141SDave May red->use_coarse_dm = v; 875*8d9f7141SDave May PetscFunctionReturn(0); 876*8d9f7141SDave May } 877*8d9f7141SDave May 8780ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v) 8790ae7c45bSDave May { 8800ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8810ae7c45bSDave May PetscFunctionBegin; 8820ae7c45bSDave May if (v) *v = red->ignore_kspcomputeoperators; 8830ae7c45bSDave May PetscFunctionReturn(0); 8840ae7c45bSDave May } 88548a10b22SPatrick Sanan 8860ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v) 8870ae7c45bSDave May { 8880ae7c45bSDave May PC_Telescope red = (PC_Telescope)pc->data; 8890ae7c45bSDave May PetscFunctionBegin; 8900ae7c45bSDave May red->ignore_kspcomputeoperators = v; 8910ae7c45bSDave May PetscFunctionReturn(0); 8920ae7c45bSDave May } 8930ae7c45bSDave May 8941e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 8951e07b27eSBarry Smith { 8961e07b27eSBarry Smith PC_Telescope red = (PC_Telescope)pc->data; 897bd49479cSSatish Balay PetscFunctionBegin; 8981e07b27eSBarry Smith *dm = private_PCTelescopeGetSubDM(red); 899bd49479cSSatish Balay PetscFunctionReturn(0); 9001e07b27eSBarry Smith } 9011e07b27eSBarry Smith 9021e07b27eSBarry Smith /*@ 9031e07b27eSBarry Smith PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 9041e07b27eSBarry Smith 9051e07b27eSBarry Smith Not Collective 9061e07b27eSBarry Smith 9071e07b27eSBarry Smith Input Parameter: 9081e07b27eSBarry Smith . pc - the preconditioner context 9091e07b27eSBarry Smith 9101e07b27eSBarry Smith Output Parameter: 9111e07b27eSBarry Smith . subksp - the KSP defined the smaller set of processes 9121e07b27eSBarry Smith 9131e07b27eSBarry Smith Level: advanced 9141e07b27eSBarry Smith 9151e07b27eSBarry Smith .keywords: PC, telescopting solve 9161e07b27eSBarry Smith @*/ 9171e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 9181e07b27eSBarry Smith { 919bd49479cSSatish Balay PetscErrorCode ierr; 920bd49479cSSatish Balay PetscFunctionBegin; 921163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 922bd49479cSSatish Balay PetscFunctionReturn(0); 9231e07b27eSBarry Smith } 9241e07b27eSBarry Smith 9251e07b27eSBarry Smith /*@ 9261e07b27eSBarry Smith PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 9271e07b27eSBarry Smith 9281e07b27eSBarry Smith Not Collective 9291e07b27eSBarry Smith 9301e07b27eSBarry Smith Input Parameter: 9311e07b27eSBarry Smith . pc - the preconditioner context 9321e07b27eSBarry Smith 9331e07b27eSBarry Smith Output Parameter: 9341e07b27eSBarry Smith . fact - the reduction factor 9351e07b27eSBarry Smith 9361e07b27eSBarry Smith Level: advanced 9371e07b27eSBarry Smith 9381e07b27eSBarry Smith .keywords: PC, telescoping solve 9391e07b27eSBarry Smith @*/ 9401e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 9411e07b27eSBarry Smith { 942bd49479cSSatish Balay PetscErrorCode ierr; 943bd49479cSSatish Balay PetscFunctionBegin; 944163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 945bd49479cSSatish Balay PetscFunctionReturn(0); 9461e07b27eSBarry Smith } 9471e07b27eSBarry Smith 9481e07b27eSBarry Smith /*@ 9491e07b27eSBarry Smith PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 9501e07b27eSBarry Smith 9511e07b27eSBarry Smith Not Collective 9521e07b27eSBarry Smith 9531e07b27eSBarry Smith Input Parameter: 9541e07b27eSBarry Smith . pc - the preconditioner context 9551e07b27eSBarry Smith 9561e07b27eSBarry Smith Output Parameter: 9571e07b27eSBarry Smith . fact - the reduction factor 9581e07b27eSBarry Smith 9591e07b27eSBarry Smith Level: advanced 9601e07b27eSBarry Smith 9611e07b27eSBarry Smith .keywords: PC, telescoping solve 9621e07b27eSBarry Smith @*/ 9631e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 9641e07b27eSBarry Smith { 965bd49479cSSatish Balay PetscErrorCode ierr; 966bd49479cSSatish Balay PetscFunctionBegin; 967bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 968bd49479cSSatish Balay PetscFunctionReturn(0); 9691e07b27eSBarry Smith } 9701e07b27eSBarry Smith 9711e07b27eSBarry Smith /*@ 9721e07b27eSBarry Smith PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 9731e07b27eSBarry Smith 9741e07b27eSBarry Smith Not Collective 9751e07b27eSBarry Smith 9761e07b27eSBarry Smith Input Parameter: 9771e07b27eSBarry Smith . pc - the preconditioner context 9781e07b27eSBarry Smith 9791e07b27eSBarry Smith Output Parameter: 9801e07b27eSBarry Smith . v - the flag 9811e07b27eSBarry Smith 9821e07b27eSBarry Smith Level: advanced 9831e07b27eSBarry Smith 9841e07b27eSBarry Smith .keywords: PC, telescoping solve 9851e07b27eSBarry Smith @*/ 9861e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 9871e07b27eSBarry Smith { 988bd49479cSSatish Balay PetscErrorCode ierr; 989bd49479cSSatish Balay PetscFunctionBegin; 990163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 991bd49479cSSatish Balay PetscFunctionReturn(0); 9921e07b27eSBarry Smith } 9931e07b27eSBarry Smith 9941e07b27eSBarry Smith /*@ 9951e07b27eSBarry Smith PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 9961e07b27eSBarry Smith 9971e07b27eSBarry Smith Not Collective 9981e07b27eSBarry Smith 9991e07b27eSBarry Smith Input Parameter: 10001e07b27eSBarry Smith . pc - the preconditioner context 10011e07b27eSBarry Smith 10021e07b27eSBarry Smith Output Parameter: 10031e07b27eSBarry Smith . v - Use PETSC_TRUE to ignore any DM 10041e07b27eSBarry Smith 10051e07b27eSBarry Smith Level: advanced 10061e07b27eSBarry Smith 10071e07b27eSBarry Smith .keywords: PC, telescoping solve 10081e07b27eSBarry Smith @*/ 1009bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 10101e07b27eSBarry Smith { 1011bd49479cSSatish Balay PetscErrorCode ierr; 1012bd49479cSSatish Balay PetscFunctionBegin; 1013bd49479cSSatish Balay ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 1014bd49479cSSatish Balay PetscFunctionReturn(0); 10151e07b27eSBarry Smith } 10161e07b27eSBarry Smith 10171e07b27eSBarry Smith /*@ 1018*8d9f7141SDave May PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used. 1019*8d9f7141SDave May 1020*8d9f7141SDave May Not Collective 1021*8d9f7141SDave May 1022*8d9f7141SDave May Input Parameter: 1023*8d9f7141SDave May . pc - the preconditioner context 1024*8d9f7141SDave May 1025*8d9f7141SDave May Output Parameter: 1026*8d9f7141SDave May . v - the flag 1027*8d9f7141SDave May 1028*8d9f7141SDave May Level: advanced 1029*8d9f7141SDave May 1030*8d9f7141SDave May .keywords: PC, telescoping solve 1031*8d9f7141SDave May @*/ 1032*8d9f7141SDave May PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc,PetscBool *v) 1033*8d9f7141SDave May { 1034*8d9f7141SDave May PetscErrorCode ierr; 1035*8d9f7141SDave May PetscFunctionBegin; 1036*8d9f7141SDave May ierr = PetscUseMethod(pc,"PCTelescopeGetUseCoarseDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 1037*8d9f7141SDave May PetscFunctionReturn(0); 1038*8d9f7141SDave May } 1039*8d9f7141SDave May 1040*8d9f7141SDave May /*@ 1041*8d9f7141SDave May PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM 1042*8d9f7141SDave May 1043*8d9f7141SDave May Not Collective 1044*8d9f7141SDave May 1045*8d9f7141SDave May Input Parameter: 1046*8d9f7141SDave May . pc - the preconditioner context 1047*8d9f7141SDave May 1048*8d9f7141SDave May Output Parameter: 1049*8d9f7141SDave May . v - Use PETSC_TRUE to ignore any DM 1050*8d9f7141SDave May 1051*8d9f7141SDave May Notes: 1052*8d9f7141SDave May When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope 1053*8d9f7141SDave May will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and 1054*8d9f7141SDave May -pc_telescope_subcomm_type will no longer have any meaning. 1055*8d9f7141SDave May It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes. 1056*8d9f7141SDave May An error will occur of the size of the communicator associated with the coarse DM 1057*8d9f7141SDave May is the same as that of the parent DM. 1058*8d9f7141SDave May Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent. 1059*8d9f7141SDave May This will be checked at the time the preconditioner is setup and an error will occur if 1060*8d9f7141SDave May the coarse DM does not define a sub-communicator of that used by the parent DM. 1061*8d9f7141SDave May 1062*8d9f7141SDave May The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of 1063*8d9f7141SDave May the DM used (e.g. it supports DMSHELL, DMPLEX, etc). 1064*8d9f7141SDave May 1065*8d9f7141SDave May Support is currently only provided for the case when you are using KSPSetComputeOperators() 1066*8d9f7141SDave May 1067*8d9f7141SDave May The user is required to compose a function with the parent DM to facilitate the transfer of fields (Vec) between the different decompositions defined by the fine and coarse DMs. 1068*8d9f7141SDave May In the user code, this is achieved via 1069*8d9f7141SDave May .vb 1070*8d9f7141SDave May { 1071*8d9f7141SDave May DM dm_fine; 1072*8d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method); 1073*8d9f7141SDave May } 1074*8d9f7141SDave May .ve 1075*8d9f7141SDave May The signature of the user provided field scatter method is 1076*8d9f7141SDave May .vb 1077*8d9f7141SDave May PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse); 1078*8d9f7141SDave May .ve 1079*8d9f7141SDave May The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE. 1080*8d9f7141SDave May SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM. 1081*8d9f7141SDave May 1082*8d9f7141SDave May Optionally, the user may also compose a function with the parent DM to facilitate the transfer 1083*8d9f7141SDave May of state variables between the fine and coarse DMs. 1084*8d9f7141SDave May In the context of a finite element discretization, an example state variable might be 1085*8d9f7141SDave May values associated with quadrature points within each element. 1086*8d9f7141SDave May A user provided state scatter method is composed via 1087*8d9f7141SDave May .vb 1088*8d9f7141SDave May { 1089*8d9f7141SDave May DM dm_fine; 1090*8d9f7141SDave May PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method); 1091*8d9f7141SDave May } 1092*8d9f7141SDave May .ve 1093*8d9f7141SDave May The signature of the user provided state scatter method is 1094*8d9f7141SDave May .vb 1095*8d9f7141SDave May PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse); 1096*8d9f7141SDave May .ve 1097*8d9f7141SDave May SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM. 1098*8d9f7141SDave May The user is only required to support mode = SCATTER_FORWARD. 1099*8d9f7141SDave May No assumption is made about the data type of the state variables. 1100*8d9f7141SDave May These must be managed by the user and must be accessible from the DM. 1101*8d9f7141SDave May 1102*8d9f7141SDave May Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be 1103*8d9f7141SDave May associated with the sub-KSP residing within PCTelescope. 1104*8d9f7141SDave May In general, PCTelescope assumes that the context on the fine and coarse DM used with 1105*8d9f7141SDave May KSPSetComputeOperators() should be "similar" in type or origin. 1106*8d9f7141SDave May Specifically the following rules are used to infer what context on the sub-KSP should be. 1107*8d9f7141SDave May 1108*8d9f7141SDave May First the contexts from the KSP and the fine and coarse DMs are retrieved. 1109*8d9f7141SDave May Note that the special case of a DMSHELL context is queried. 1110*8d9f7141SDave May 1111*8d9f7141SDave May .vb 1112*8d9f7141SDave May DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx); 1113*8d9f7141SDave May DMGetApplicationContext(dm_fine,&dmfine_appctx); 1114*8d9f7141SDave May DMShellGetContext(dm_fine,&dmfine_shellctx); 1115*8d9f7141SDave May 1116*8d9f7141SDave May DMGetApplicationContext(dm_coarse,&dmcoarse_appctx); 1117*8d9f7141SDave May DMShellGetContext(dm_coarse,&dmcoarse_shellctx); 1118*8d9f7141SDave May .ve 1119*8d9f7141SDave May 1120*8d9f7141SDave May The following rules are then enforced: 1121*8d9f7141SDave May 1122*8d9f7141SDave May 1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP: 1123*8d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL); 1124*8d9f7141SDave May 1125*8d9f7141SDave May 2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx, 1126*8d9f7141SDave May check that dmcoarse_appctx is also non-NULL. If this is true, then: 1127*8d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx); 1128*8d9f7141SDave May 1129*8d9f7141SDave May 3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx, 1130*8d9f7141SDave May check that dmcoarse_shellctx is also non-NULL. If this is true, then: 1131*8d9f7141SDave May KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx); 1132*8d9f7141SDave May 1133*8d9f7141SDave May If neither of the above three tests passed, then PCTelescope cannot safely determine what 1134*8d9f7141SDave May context should be provided to KSPSetComputeOperators() for use with the sub-KSP. 1135*8d9f7141SDave May In this case, an additional mechanism is provided via a composed function which will return 1136*8d9f7141SDave May the actual context to be used. To use this feature you must compose the "getter" function 1137*8d9f7141SDave May with the coarse DM, e.g. 1138*8d9f7141SDave May .vb 1139*8d9f7141SDave May { 1140*8d9f7141SDave May DM dm_coarse; 1141*8d9f7141SDave May PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter); 1142*8d9f7141SDave May } 1143*8d9f7141SDave May .ve 1144*8d9f7141SDave May The signature of the user provided method is 1145*8d9f7141SDave May .vb 1146*8d9f7141SDave May PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext); 1147*8d9f7141SDave May .ve 1148*8d9f7141SDave May 1149*8d9f7141SDave May Level: advanced 1150*8d9f7141SDave May 1151*8d9f7141SDave May .keywords: PC, telescoping solve 1152*8d9f7141SDave May @*/ 1153*8d9f7141SDave May PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc,PetscBool v) 1154*8d9f7141SDave May { 1155*8d9f7141SDave May PetscErrorCode ierr; 1156*8d9f7141SDave May PetscFunctionBegin; 1157*8d9f7141SDave May ierr = PetscTryMethod(pc,"PCTelescopeSetUseCoarseDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 1158*8d9f7141SDave May PetscFunctionReturn(0); 1159*8d9f7141SDave May } 1160*8d9f7141SDave May 1161*8d9f7141SDave May /*@ 11620ae7c45bSDave May PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used. 11630ae7c45bSDave May 11640ae7c45bSDave May Not Collective 11650ae7c45bSDave May 11660ae7c45bSDave May Input Parameter: 11670ae7c45bSDave May . pc - the preconditioner context 11680ae7c45bSDave May 11690ae7c45bSDave May Output Parameter: 11700ae7c45bSDave May . v - the flag 11710ae7c45bSDave May 11720ae7c45bSDave May Level: advanced 11730ae7c45bSDave May 11740ae7c45bSDave May .keywords: PC, telescoping solve 11750ae7c45bSDave May @*/ 11760ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v) 11770ae7c45bSDave May { 11780ae7c45bSDave May PetscErrorCode ierr; 11790ae7c45bSDave May PetscFunctionBegin; 1180163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 11810ae7c45bSDave May PetscFunctionReturn(0); 11820ae7c45bSDave May } 11830ae7c45bSDave May 11840ae7c45bSDave May /*@ 11850ae7c45bSDave May PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators. 11860ae7c45bSDave May 11870ae7c45bSDave May Not Collective 11880ae7c45bSDave May 11890ae7c45bSDave May Input Parameter: 11900ae7c45bSDave May . pc - the preconditioner context 11910ae7c45bSDave May 11920ae7c45bSDave May Output Parameter: 1193a954d8f4SDave May . v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc 11940ae7c45bSDave May 11950ae7c45bSDave May Level: advanced 11960ae7c45bSDave May 11970ae7c45bSDave May .keywords: PC, telescoping solve 11980ae7c45bSDave May @*/ 11990ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v) 12000ae7c45bSDave May { 12010ae7c45bSDave May PetscErrorCode ierr; 12020ae7c45bSDave May PetscFunctionBegin; 12030ae7c45bSDave May ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 12040ae7c45bSDave May PetscFunctionReturn(0); 12050ae7c45bSDave May } 12060ae7c45bSDave May 12070ae7c45bSDave May /*@ 12081e07b27eSBarry Smith PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 12091e07b27eSBarry Smith 12101e07b27eSBarry Smith Not Collective 12111e07b27eSBarry Smith 12121e07b27eSBarry Smith Input Parameter: 12131e07b27eSBarry Smith . pc - the preconditioner context 12141e07b27eSBarry Smith 12151e07b27eSBarry Smith Output Parameter: 12161e07b27eSBarry Smith . subdm - The re-partitioned DM 12171e07b27eSBarry Smith 12181e07b27eSBarry Smith Level: advanced 12191e07b27eSBarry Smith 12201e07b27eSBarry Smith .keywords: PC, telescoping solve 12211e07b27eSBarry Smith @*/ 12221e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 12231e07b27eSBarry Smith { 1224bd49479cSSatish Balay PetscErrorCode ierr; 1225bd49479cSSatish Balay PetscFunctionBegin; 1226163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 1227bd49479cSSatish Balay PetscFunctionReturn(0); 12281e07b27eSBarry Smith } 12291e07b27eSBarry Smith 123048a10b22SPatrick Sanan /*@ 123148a10b22SPatrick Sanan PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous) 123248a10b22SPatrick Sanan 123348a10b22SPatrick Sanan Logically Collective 123448a10b22SPatrick Sanan 123548a10b22SPatrick Sanan Input Parameter: 12361dae98e4SBarry Smith + pc - the preconditioner context 12371dae98e4SBarry Smith - subcommtype - the subcommunicator type (see PetscSubcommType) 123848a10b22SPatrick Sanan 123948a10b22SPatrick Sanan Level: advanced 124048a10b22SPatrick Sanan 124148a10b22SPatrick Sanan .keywords: PC, telescoping solve 124248a10b22SPatrick Sanan 124348a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE 124448a10b22SPatrick Sanan @*/ 124548a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) 124648a10b22SPatrick Sanan { 124748a10b22SPatrick Sanan PetscErrorCode ierr; 124848a10b22SPatrick Sanan PetscFunctionBegin; 124948a10b22SPatrick Sanan ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr); 125048a10b22SPatrick Sanan PetscFunctionReturn(0); 125148a10b22SPatrick Sanan } 125248a10b22SPatrick Sanan 125348a10b22SPatrick Sanan /*@ 125448a10b22SPatrick Sanan PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous) 125548a10b22SPatrick Sanan 125648a10b22SPatrick Sanan Not Collective 125748a10b22SPatrick Sanan 125848a10b22SPatrick Sanan Input Parameter: 125948a10b22SPatrick Sanan . pc - the preconditioner context 126048a10b22SPatrick Sanan 126148a10b22SPatrick Sanan Output Parameter: 126248a10b22SPatrick Sanan . subcommtype - the subcommunicator type (see PetscSubcommType) 126348a10b22SPatrick Sanan 126448a10b22SPatrick Sanan Level: advanced 126548a10b22SPatrick Sanan 126648a10b22SPatrick Sanan .keywords: PC, telescoping solve 126748a10b22SPatrick Sanan 12681dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE 126948a10b22SPatrick Sanan @*/ 12701dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) 127148a10b22SPatrick Sanan { 127248a10b22SPatrick Sanan PetscErrorCode ierr; 127348a10b22SPatrick Sanan PetscFunctionBegin; 127448a10b22SPatrick Sanan ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr); 127548a10b22SPatrick Sanan PetscFunctionReturn(0); 127648a10b22SPatrick Sanan } 127748a10b22SPatrick Sanan 12781e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/ 12791e07b27eSBarry Smith /*MC 12801e07b27eSBarry Smith PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 12811e07b27eSBarry Smith 12821e07b27eSBarry Smith Options Database: 1283a04a6428SPatrick 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. 1284a04a6428SPatrick Sanan - -pc_telescope_ignore_dm - flag to indicate whether an attached DM should be ignored 1285a04a6428SPatrick Sanan - -pc_telescope_subcomm_type <interlaced,contiguous> - how to define the reduced communicator. see PetscSubcomm for more. 12861e07b27eSBarry Smith 12871e07b27eSBarry Smith Level: advanced 12881e07b27eSBarry Smith 12891e07b27eSBarry Smith Notes: 12906fc41876SBarry Smith The preconditioner is deemed telescopic as it only calls KSPSolve() on a single 12918439623fSPatrick Sanan sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 12928439623fSPatrick Sanan This means there will be MPI processes which will be idle during the application of this preconditioner. 12936fc41876SBarry Smith 12941e07b27eSBarry Smith The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 12951e07b27eSBarry Smith Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 12961e07b27eSBarry Smith Currently only support for re-partitioning a DMDA is provided. 12978439623fSPatrick Sanan Any nullspace attached to the original Bmat operator is extracted, re-partitioned and set on the repartitioned Bmat operator. 12981e07b27eSBarry Smith KSPSetComputeOperators() is not propagated to the sub KSP. 12991e07b27eSBarry Smith Currently there is no support for the flag -pc_use_amat 13001e07b27eSBarry Smith 13016fc41876SBarry Smith Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 13028439623fSPatrick Sanan creates a child sub-communicator (c') containing fewer MPI processes than the original parent preconditioner (PC). 13036fc41876SBarry Smith 13046fc41876SBarry Smith Developer Notes: 13056fc41876SBarry Smith During PCSetup, the B operator is scattered onto c'. 13066fc41876SBarry Smith Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 13078439623fSPatrick Sanan Then, KSPSolve() is executed on the c' communicator. 13086fc41876SBarry Smith 13096fc41876SBarry Smith The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 1310a04a6428SPatrick 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. 13116fc41876SBarry Smith 1312005d9f20SPatrick Sanan The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator. 13138439623fSPatrick Sanan In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into 1314005d9f20SPatrick Sanan a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c') 13156fc41876SBarry Smith 13166fc41876SBarry Smith The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D - 13176fc41876SBarry 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 13186fc41876SBarry Smith is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support 13198439623fSPatrick Sanan for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 13206fc41876SBarry Smith 13216fc41876SBarry Smith By default, B' is defined by simply fusing rows from different MPI processes 13226fc41876SBarry Smith 13238439623fSPatrick Sanan When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B 13247dae84e0SHong Zhang into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the 13256fc41876SBarry Smith locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 13266fc41876SBarry Smith 13278439623fSPatrick Sanan Limitations/improvements include the following. 13288439623fSPatrick Sanan VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage. 13296fc41876SBarry Smith 13306fc41876SBarry Smith The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P, 13318439623fSPatrick 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 13328439623fSPatrick Sanan VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a 13336fc41876SBarry Smith consistent manner. 13346fc41876SBarry Smith 13358439623fSPatrick Sanan Mapping of vectors is performed in the following way. 13368439623fSPatrick 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. 13378439623fSPatrick Sanan Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2. 13388439623fSPatrick Sanan We perform the scatter to the sub-comm in the following way. 13396fc41876SBarry Smith [1] Given a vector x defined on comm c 13406fc41876SBarry Smith 13416fc41876SBarry Smith rank(c) : _________ 0 ______ ________ 1 _______ ________ 2 _____________ ___________ 3 __________ 13426fc41876SBarry 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] 13436fc41876SBarry Smith 13446fc41876SBarry Smith scatter to xtmp defined also on comm c so that we have the following values 13456fc41876SBarry Smith 13466fc41876SBarry Smith rank(c) : ___________________ 0 ________________ _1_ ______________________ 2 _______________________ __3_ 13476fc41876SBarry 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] [ ] 13486fc41876SBarry Smith 13496fc41876SBarry Smith The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 13506fc41876SBarry Smith 13516fc41876SBarry Smith 13526fc41876SBarry 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'. 13536fc41876SBarry Smith Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 13546fc41876SBarry Smith 13556fc41876SBarry Smith rank(c') : ___________________ 0 _______________ ______________________ 1 _____________________ 13566fc41876SBarry 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] 13576fc41876SBarry Smith 13581e07b27eSBarry Smith 13591e07b27eSBarry Smith Contributed by Dave May 13601e07b27eSBarry Smith 1361bf00f589SPatrick Sanan Reference: 1362bf00f589SPatrick 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 1363bf00f589SPatrick Sanan 13646fc41876SBarry Smith .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 13651e07b27eSBarry Smith M*/ 13661e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 13671e07b27eSBarry Smith { 13681e07b27eSBarry Smith PetscErrorCode ierr; 13691e07b27eSBarry Smith struct _PC_Telescope *sred; 13701e07b27eSBarry Smith 13711e07b27eSBarry Smith PetscFunctionBegin; 13721e07b27eSBarry Smith ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 13732a22aa42SDave May sred->psubcomm = NULL; 137448a10b22SPatrick Sanan sred->subcommtype = PETSC_SUBCOMM_INTERLACED; 13752a22aa42SDave May sred->subcomm = MPI_COMM_NULL; 13761e07b27eSBarry Smith sred->redfactor = 1; 13771e07b27eSBarry Smith sred->ignore_dm = PETSC_FALSE; 13787c5279cbSDave May sred->ignore_kspcomputeoperators = PETSC_FALSE; 1379*8d9f7141SDave May sred->use_coarse_dm = PETSC_FALSE; 13801e07b27eSBarry Smith pc->data = (void*)sred; 13811e07b27eSBarry Smith 13821e07b27eSBarry Smith pc->ops->apply = PCApply_Telescope; 13831e07b27eSBarry Smith pc->ops->applytranspose = NULL; 1384f650675bSDave May pc->ops->applyrichardson = PCApplyRichardson_Telescope; 13851e07b27eSBarry Smith pc->ops->setup = PCSetUp_Telescope; 13861e07b27eSBarry Smith pc->ops->destroy = PCDestroy_Telescope; 13871e07b27eSBarry Smith pc->ops->reset = PCReset_Telescope; 13881e07b27eSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Telescope; 13891e07b27eSBarry Smith pc->ops->view = PCView_Telescope; 13901e07b27eSBarry Smith 13911e07b27eSBarry Smith sred->pctelescope_setup_type = PCTelescopeSetUp_default; 13921e07b27eSBarry Smith sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 13931e07b27eSBarry Smith sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 13941e07b27eSBarry Smith sred->pctelescope_reset_type = NULL; 13951e07b27eSBarry Smith 13961e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 139748a10b22SPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr); 139848a10b22SPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr); 13991e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 14001e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 14011e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 14021e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 14030ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 14040ae7c45bSDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 14051e07b27eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 1406*8d9f7141SDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetUseCoarseDM_C",PCTelescopeGetUseCoarseDM_Telescope);CHKERRQ(ierr); 1407*8d9f7141SDave May ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetUseCoarseDM_C",PCTelescopeSetUseCoarseDM_Telescope);CHKERRQ(ierr); 14081e07b27eSBarry Smith PetscFunctionReturn(0); 14091e07b27eSBarry Smith } 1410