11e07b27eSBarry Smith 2120bdd93SDave May 3120bdd93SDave May #include <petsc/private/matimpl.h> 41e07b27eSBarry Smith #include <petsc/private/pcimpl.h> 55e897e82SDave May #include <petsc/private/dmimpl.h> 61e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/ 71e07b27eSBarry Smith #include <petscdm.h> 81e07b27eSBarry Smith #include <petscdmda.h> 91e07b27eSBarry Smith 10575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h" 111e07b27eSBarry Smith 12bf00f589SPatrick Sanan static PetscBool cited = PETSC_FALSE; 13bf00f589SPatrick Sanan static const char citation[] = 14bf00f589SPatrick Sanan "@inproceedings{MaySananRuppKnepleySmith2016,\n" 15bf00f589SPatrick Sanan " title = {Extreme-Scale Multigrid Components within PETSc},\n" 16bf00f589SPatrick Sanan " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 17bf00f589SPatrick Sanan " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 18bf00f589SPatrick Sanan " series = {PASC '16},\n" 19bf00f589SPatrick Sanan " isbn = {978-1-4503-4126-4},\n" 20bf00f589SPatrick Sanan " location = {Lausanne, Switzerland},\n" 21bf00f589SPatrick Sanan " pages = {5:1--5:12},\n" 22bf00f589SPatrick Sanan " articleno = {5},\n" 23bf00f589SPatrick Sanan " numpages = {12},\n" 24a8d69d7bSBarry Smith " url = {https://doi.acm.org/10.1145/2929908.2929913},\n" 25bf00f589SPatrick Sanan " doi = {10.1145/2929908.2929913},\n" 26bf00f589SPatrick Sanan " acmid = {2929913},\n" 27bf00f589SPatrick Sanan " publisher = {ACM},\n" 28bf00f589SPatrick Sanan " address = {New York, NY, USA},\n" 29bf00f589SPatrick Sanan " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 30bf00f589SPatrick Sanan " year = {2016}\n" 31bf00f589SPatrick Sanan "}\n"; 32bf00f589SPatrick Sanan 3357f12427SDave May static PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim,PetscInt i,PetscInt j,PetscInt k, 341e07b27eSBarry Smith PetscInt Mp,PetscInt Np,PetscInt Pp, 351e07b27eSBarry Smith PetscInt start_i[],PetscInt start_j[],PetscInt start_k[], 361e07b27eSBarry Smith PetscInt span_i[],PetscInt span_j[],PetscInt span_k[], 371e07b27eSBarry Smith PetscMPIInt *_pi,PetscMPIInt *_pj,PetscMPIInt *_pk,PetscMPIInt *rank_re) 381e07b27eSBarry Smith { 391e07b27eSBarry Smith PetscInt pi,pj,pk,n; 401e07b27eSBarry Smith 411e07b27eSBarry Smith PetscFunctionBegin; 42137d0469SJed Brown *rank_re = -1; 43137d0469SJed Brown if (_pi) *_pi = -1; 44137d0469SJed Brown if (_pj) *_pj = -1; 45137d0469SJed Brown if (_pk) *_pk = -1; 461e07b27eSBarry Smith pi = pj = pk = -1; 471e07b27eSBarry Smith if (_pi) { 481e07b27eSBarry Smith for (n=0; n<Mp; n++) { 491e07b27eSBarry Smith if ( (i >= start_i[n]) && (i < start_i[n]+span_i[n]) ) { 501e07b27eSBarry Smith pi = n; 511e07b27eSBarry Smith break; 521e07b27eSBarry Smith } 531e07b27eSBarry Smith } 541e07b27eSBarry Smith if (pi == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pi cannot be determined : range %D, val %D",Mp,i); 551e07b27eSBarry Smith *_pi = pi; 561e07b27eSBarry Smith } 571e07b27eSBarry Smith 581e07b27eSBarry Smith if (_pj) { 591e07b27eSBarry Smith for (n=0; n<Np; n++) { 601e07b27eSBarry Smith if ( (j >= start_j[n]) && (j < start_j[n]+span_j[n]) ) { 611e07b27eSBarry Smith pj = n; 621e07b27eSBarry Smith break; 631e07b27eSBarry Smith } 641e07b27eSBarry Smith } 651e07b27eSBarry Smith if (pj == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pj cannot be determined : range %D, val %D",Np,j); 661e07b27eSBarry Smith *_pj = pj; 671e07b27eSBarry Smith } 681e07b27eSBarry Smith 691e07b27eSBarry Smith if (_pk) { 701e07b27eSBarry Smith for (n=0; n<Pp; n++) { 711e07b27eSBarry Smith if ( (k >= start_k[n]) && (k < start_k[n]+span_k[n]) ) { 721e07b27eSBarry Smith pk = n; 731e07b27eSBarry Smith break; 741e07b27eSBarry Smith } 751e07b27eSBarry Smith } 761e07b27eSBarry Smith if (pk == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pk cannot be determined : range %D, val %D",Pp,k); 771e07b27eSBarry Smith *_pk = pk; 781e07b27eSBarry Smith } 791e07b27eSBarry Smith 801e07b27eSBarry Smith switch (dim) { 811e07b27eSBarry Smith case 1: 821e07b27eSBarry Smith *rank_re = pi; 831e07b27eSBarry Smith break; 841e07b27eSBarry Smith case 2: 851e07b27eSBarry Smith *rank_re = pi + pj * Mp; 861e07b27eSBarry Smith break; 871e07b27eSBarry Smith case 3: 881e07b27eSBarry Smith *rank_re = pi + pj * Mp + pk * (Mp*Np); 891e07b27eSBarry Smith break; 901e07b27eSBarry Smith } 911e07b27eSBarry Smith PetscFunctionReturn(0); 921e07b27eSBarry Smith } 931e07b27eSBarry Smith 9457f12427SDave May static PetscErrorCode _DMDADetermineGlobalS0(PetscInt dim,PetscMPIInt rank_re,PetscInt Mp_re,PetscInt Np_re,PetscInt Pp_re, 951e07b27eSBarry Smith PetscInt range_i_re[],PetscInt range_j_re[],PetscInt range_k_re[],PetscInt *s0) 961e07b27eSBarry Smith { 97c6a0d831SBarry Smith PetscInt i,j,k,start_IJK = 0; 981e07b27eSBarry Smith PetscInt rank_ijk; 991e07b27eSBarry Smith 1001e07b27eSBarry Smith PetscFunctionBegin; 1011e07b27eSBarry Smith switch (dim) { 1021e07b27eSBarry Smith case 1: 1031e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 1041e07b27eSBarry Smith rank_ijk = i; 1051e07b27eSBarry Smith if (rank_ijk < rank_re) { 1061e07b27eSBarry Smith start_IJK += range_i_re[i]; 1071e07b27eSBarry Smith } 1081e07b27eSBarry Smith } 1091e07b27eSBarry Smith break; 1101e07b27eSBarry Smith case 2: 1111e07b27eSBarry Smith for (j=0; j<Np_re; j++) { 1121e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 1131e07b27eSBarry Smith rank_ijk = i + j*Mp_re; 1141e07b27eSBarry Smith if (rank_ijk < rank_re) { 1151e07b27eSBarry Smith start_IJK += range_i_re[i]*range_j_re[j]; 1161e07b27eSBarry Smith } 1171e07b27eSBarry Smith } 1181e07b27eSBarry Smith } 1191e07b27eSBarry Smith break; 1201e07b27eSBarry Smith case 3: 1211e07b27eSBarry Smith for (k=0; k<Pp_re; k++) { 1221e07b27eSBarry Smith for (j=0; j<Np_re; j++) { 1231e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 1241e07b27eSBarry Smith rank_ijk = i + j*Mp_re + k*Mp_re*Np_re; 1251e07b27eSBarry Smith if (rank_ijk < rank_re) { 1261e07b27eSBarry Smith start_IJK += range_i_re[i]*range_j_re[j]*range_k_re[k]; 1271e07b27eSBarry Smith } 1281e07b27eSBarry Smith } 1291e07b27eSBarry Smith } 1301e07b27eSBarry Smith } 1311e07b27eSBarry Smith break; 1321e07b27eSBarry Smith } 1331e07b27eSBarry Smith *s0 = start_IJK; 1341e07b27eSBarry Smith PetscFunctionReturn(0); 1351e07b27eSBarry Smith } 1361e07b27eSBarry Smith 13757f12427SDave May static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PC_Telescope sred,DM dm,DM subdm) 1381e07b27eSBarry Smith { 1391e07b27eSBarry Smith PetscErrorCode ierr; 1401e07b27eSBarry Smith DM cdm; 1411e07b27eSBarry Smith Vec coor,coor_natural,perm_coors; 1421e07b27eSBarry Smith PetscInt i,j,si,sj,ni,nj,M,N,Ml,Nl,c,nidx; 1431e07b27eSBarry Smith PetscInt *fine_indices; 1441e07b27eSBarry Smith IS is_fine,is_local; 1451e07b27eSBarry Smith VecScatter sctx; 1461e07b27eSBarry Smith 1471e07b27eSBarry Smith PetscFunctionBegin; 1481e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 1491e07b27eSBarry Smith if (!coor) return(0); 15057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1511e07b27eSBarry Smith ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1521e07b27eSBarry Smith } 1531e07b27eSBarry Smith /* Get the coordinate vector from the distributed array */ 1541e07b27eSBarry Smith ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 1551e07b27eSBarry Smith ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 1561e07b27eSBarry Smith 1571e07b27eSBarry Smith ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 1581e07b27eSBarry Smith ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 1591e07b27eSBarry Smith 1601e07b27eSBarry Smith /* get indices of the guys I want to grab */ 1611e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 16257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1631e07b27eSBarry Smith ierr = DMDAGetCorners(subdm,&si,&sj,NULL,&ni,&nj,NULL);CHKERRQ(ierr); 16415dd08bcSBarry Smith Ml = ni; 16515dd08bcSBarry Smith Nl = nj; 1661e07b27eSBarry Smith } else { 167c41e779fSDave May si = sj = 0; 168c41e779fSDave May ni = nj = 0; 1693ac26c5eSBarry Smith Ml = Nl = 0; 1701e07b27eSBarry Smith } 1711e07b27eSBarry Smith 172e3acf2f7SBarry Smith ierr = PetscMalloc1(Ml*Nl*2,&fine_indices);CHKERRQ(ierr); 1731e07b27eSBarry Smith c = 0; 17457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1751e07b27eSBarry Smith for (j=sj; j<sj+nj; j++) { 1761e07b27eSBarry Smith for (i=si; i<si+ni; i++) { 1771e07b27eSBarry Smith nidx = (i) + (j)*M; 1781e07b27eSBarry Smith fine_indices[c ] = 2 * nidx ; 1791e07b27eSBarry Smith fine_indices[c+1] = 2 * nidx + 1 ; 1801e07b27eSBarry Smith c = c + 2; 1811e07b27eSBarry Smith } 1821e07b27eSBarry Smith } 183c2be2c73SBarry Smith if (c != Ml*Nl*2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c %D should equal 2 * Ml %D * Nl %D",c,Ml,Nl); 1841e07b27eSBarry Smith } 1851e07b27eSBarry Smith 1861e07b27eSBarry Smith /* generate scatter */ 1871e07b27eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*2,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 1881e07b27eSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*2,0,1,&is_local);CHKERRQ(ierr); 1891e07b27eSBarry Smith 1901e07b27eSBarry Smith /* scatter */ 1911e07b27eSBarry Smith ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 1921e07b27eSBarry Smith ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*2);CHKERRQ(ierr); 1931e07b27eSBarry Smith ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 1941e07b27eSBarry Smith 1959448b7f1SJunchao Zhang ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 1961e07b27eSBarry Smith ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1971e07b27eSBarry Smith ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1981e07b27eSBarry Smith /* access */ 19957f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2001e07b27eSBarry Smith Vec _coors; 2011e07b27eSBarry Smith const PetscScalar *LA_perm; 2021e07b27eSBarry Smith PetscScalar *LA_coors; 2031e07b27eSBarry Smith 2041e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 2051e07b27eSBarry Smith ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 2061e07b27eSBarry Smith ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 2071e07b27eSBarry Smith for (i=0; i<Ml*Nl*2; i++) { 2081e07b27eSBarry Smith LA_coors[i] = LA_perm[i]; 2091e07b27eSBarry Smith } 2101e07b27eSBarry Smith ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 2111e07b27eSBarry Smith ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 2121e07b27eSBarry Smith } 2131e07b27eSBarry Smith 2141e07b27eSBarry Smith /* update local coords */ 21557f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2161e07b27eSBarry Smith DM _dmc; 2171e07b27eSBarry Smith Vec _coors,_coors_local; 2181e07b27eSBarry Smith ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 2191e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 2201e07b27eSBarry Smith ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 2211e07b27eSBarry Smith ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 2221e07b27eSBarry Smith ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 2231e07b27eSBarry Smith } 2241e07b27eSBarry Smith ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 2251e07b27eSBarry Smith ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 2261e07b27eSBarry Smith ierr = PetscFree(fine_indices);CHKERRQ(ierr); 2271e07b27eSBarry Smith ierr = ISDestroy(&is_local);CHKERRQ(ierr); 2281e07b27eSBarry Smith ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 2291e07b27eSBarry Smith ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 2301e07b27eSBarry Smith PetscFunctionReturn(0); 2311e07b27eSBarry Smith } 2321e07b27eSBarry Smith 23357f12427SDave May static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PC_Telescope sred,DM dm,DM subdm) 2341e07b27eSBarry Smith { 2351e07b27eSBarry Smith PetscErrorCode ierr; 2361e07b27eSBarry Smith DM cdm; 2371e07b27eSBarry Smith Vec coor,coor_natural,perm_coors; 2381e07b27eSBarry Smith PetscInt i,j,k,si,sj,sk,ni,nj,nk,M,N,P,Ml,Nl,Pl,c,nidx; 2391e07b27eSBarry Smith PetscInt *fine_indices; 2401e07b27eSBarry Smith IS is_fine,is_local; 2411e07b27eSBarry Smith VecScatter sctx; 2421e07b27eSBarry Smith 2431e07b27eSBarry Smith PetscFunctionBegin; 2441e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 2451e07b27eSBarry Smith if (!coor) PetscFunctionReturn(0); 2461e07b27eSBarry Smith 24757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2481e07b27eSBarry Smith ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 2491e07b27eSBarry Smith } 2501e07b27eSBarry Smith 2511e07b27eSBarry Smith /* Get the coordinate vector from the distributed array */ 2521e07b27eSBarry Smith ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 2531e07b27eSBarry Smith ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 2541e07b27eSBarry Smith ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 2551e07b27eSBarry Smith ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 2561e07b27eSBarry Smith 2571e07b27eSBarry Smith /* get indices of the guys I want to grab */ 2581e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2591e07b27eSBarry Smith 26057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2611e07b27eSBarry Smith ierr = DMDAGetCorners(subdm,&si,&sj,&sk,&ni,&nj,&nk);CHKERRQ(ierr); 262553d0ae9SBarry Smith Ml = ni; 263553d0ae9SBarry Smith Nl = nj; 264553d0ae9SBarry Smith Pl = nk; 2651e07b27eSBarry Smith } else { 266c41e779fSDave May si = sj = sk = 0; 267c41e779fSDave May ni = nj = nk = 0; 2683ac26c5eSBarry Smith Ml = Nl = Pl = 0; 2691e07b27eSBarry Smith } 2701e07b27eSBarry Smith 271e3acf2f7SBarry Smith ierr = PetscMalloc1(Ml*Nl*Pl*3,&fine_indices);CHKERRQ(ierr); 2721e07b27eSBarry Smith 2731e07b27eSBarry Smith c = 0; 27457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 2751e07b27eSBarry Smith for (k=sk; k<sk+nk; k++) { 2761e07b27eSBarry Smith for (j=sj; j<sj+nj; j++) { 2771e07b27eSBarry Smith for (i=si; i<si+ni; i++) { 2781e07b27eSBarry Smith nidx = (i) + (j)*M + (k)*M*N; 2791e07b27eSBarry Smith fine_indices[c ] = 3 * nidx ; 2801e07b27eSBarry Smith fine_indices[c+1] = 3 * nidx + 1 ; 2811e07b27eSBarry Smith fine_indices[c+2] = 3 * nidx + 2 ; 2821e07b27eSBarry Smith c = c + 3; 2831e07b27eSBarry Smith } 2841e07b27eSBarry Smith } 2851e07b27eSBarry Smith } 2861e07b27eSBarry Smith } 2871e07b27eSBarry Smith 2881e07b27eSBarry Smith /* generate scatter */ 2891e07b27eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*Pl*3,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 2901e07b27eSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*Pl*3,0,1,&is_local);CHKERRQ(ierr); 2911e07b27eSBarry Smith 2921e07b27eSBarry Smith /* scatter */ 2931e07b27eSBarry Smith ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 2941e07b27eSBarry Smith ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*Pl*3);CHKERRQ(ierr); 2951e07b27eSBarry Smith ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 2969448b7f1SJunchao Zhang ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 2971e07b27eSBarry Smith ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2981e07b27eSBarry Smith ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2991e07b27eSBarry Smith 3001e07b27eSBarry Smith /* access */ 30157f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3021e07b27eSBarry Smith Vec _coors; 3031e07b27eSBarry Smith const PetscScalar *LA_perm; 3041e07b27eSBarry Smith PetscScalar *LA_coors; 3051e07b27eSBarry Smith 3061e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 3071e07b27eSBarry Smith ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 3081e07b27eSBarry Smith ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 3091e07b27eSBarry Smith for (i=0; i<Ml*Nl*Pl*3; i++) { 3101e07b27eSBarry Smith LA_coors[i] = LA_perm[i]; 3111e07b27eSBarry Smith } 3121e07b27eSBarry Smith ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 3131e07b27eSBarry Smith ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 3141e07b27eSBarry Smith } 3151e07b27eSBarry Smith 3161e07b27eSBarry Smith /* update local coords */ 31757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3181e07b27eSBarry Smith DM _dmc; 3191e07b27eSBarry Smith Vec _coors,_coors_local; 3201e07b27eSBarry Smith 3211e07b27eSBarry Smith ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 3221e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 3231e07b27eSBarry Smith ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 3241e07b27eSBarry Smith ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 3251e07b27eSBarry Smith ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 3261e07b27eSBarry Smith } 3271e07b27eSBarry Smith 3281e07b27eSBarry Smith ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 3291e07b27eSBarry Smith ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 3301e07b27eSBarry Smith ierr = PetscFree(fine_indices);CHKERRQ(ierr); 3311e07b27eSBarry Smith ierr = ISDestroy(&is_local);CHKERRQ(ierr); 3321e07b27eSBarry Smith ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 3331e07b27eSBarry Smith ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 3341e07b27eSBarry Smith PetscFunctionReturn(0); 3351e07b27eSBarry Smith } 3361e07b27eSBarry Smith 33757f12427SDave May static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 3381e07b27eSBarry Smith { 3391e07b27eSBarry Smith PetscInt dim; 3401e07b27eSBarry Smith DM dm,subdm; 3411e07b27eSBarry Smith PetscSubcomm psubcomm; 3421e07b27eSBarry Smith MPI_Comm comm; 3431e07b27eSBarry Smith Vec coor; 3441e07b27eSBarry Smith PetscErrorCode ierr; 3451e07b27eSBarry Smith 3461e07b27eSBarry Smith PetscFunctionBegin; 3471e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3481e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 3491e07b27eSBarry Smith if (!coor) PetscFunctionReturn(0); 3501e07b27eSBarry Smith psubcomm = sred->psubcomm; 3511e07b27eSBarry Smith comm = PetscSubcommParent(psubcomm); 3521e07b27eSBarry Smith subdm = ctx->dmrepart; 3531e07b27eSBarry Smith 3541e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the coordinates (DMDA)\n");CHKERRQ(ierr); 3551e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 3561e07b27eSBarry Smith switch (dim) { 3571e07b27eSBarry Smith case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 3581e07b27eSBarry Smith break; 35957f12427SDave May case 2: ierr = PCTelescopeSetUp_dmda_repart_coors2d(sred,dm,subdm);CHKERRQ(ierr); 3601e07b27eSBarry Smith break; 36157f12427SDave May case 3: ierr = PCTelescopeSetUp_dmda_repart_coors3d(sred,dm,subdm);CHKERRQ(ierr); 3621e07b27eSBarry Smith break; 3631e07b27eSBarry Smith } 3641e07b27eSBarry Smith PetscFunctionReturn(0); 3651e07b27eSBarry Smith } 3661e07b27eSBarry Smith 3671e07b27eSBarry Smith /* setup repartitioned dm */ 3681e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 3691e07b27eSBarry Smith { 3701e07b27eSBarry Smith PetscErrorCode ierr; 3711e07b27eSBarry Smith DM dm; 3721e07b27eSBarry Smith PetscInt dim,nx,ny,nz,ndof,nsw,sum,k; 3731e07b27eSBarry Smith DMBoundaryType bx,by,bz; 3741e07b27eSBarry Smith DMDAStencilType stencil; 3751e07b27eSBarry Smith const PetscInt *_range_i_re; 3761e07b27eSBarry Smith const PetscInt *_range_j_re; 3771e07b27eSBarry Smith const PetscInt *_range_k_re; 3781e07b27eSBarry Smith DMDAInterpolationType itype; 3791e07b27eSBarry Smith PetscInt refine_x,refine_y,refine_z; 3801e07b27eSBarry Smith MPI_Comm comm,subcomm; 3811e07b27eSBarry Smith const char *prefix; 3821e07b27eSBarry Smith 3831e07b27eSBarry Smith PetscFunctionBegin; 3841e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 3851e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 3861e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3871e07b27eSBarry Smith 3881e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,&nsw,&bx,&by,&bz,&stencil);CHKERRQ(ierr); 3891e07b27eSBarry Smith ierr = DMDAGetInterpolationType(dm,&itype);CHKERRQ(ierr); 3901e07b27eSBarry Smith ierr = DMDAGetRefinementFactor(dm,&refine_x,&refine_y,&refine_z);CHKERRQ(ierr); 3911e07b27eSBarry Smith 3921e07b27eSBarry Smith ctx->dmrepart = NULL; 3931e07b27eSBarry Smith _range_i_re = _range_j_re = _range_k_re = NULL; 3941e07b27eSBarry Smith /* Create DMDA on the child communicator */ 39557f12427SDave May if (PCTelescope_isActiveRank(sred)) { 3961e07b27eSBarry Smith switch (dim) { 3971e07b27eSBarry Smith case 1: 3981e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (1D)\n");CHKERRQ(ierr); 3991e07b27eSBarry Smith /*ierr = DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 4001e07b27eSBarry Smith ny = nz = 1; 4011e07b27eSBarry Smith by = bz = DM_BOUNDARY_NONE; 4021e07b27eSBarry Smith break; 4031e07b27eSBarry Smith case 2: 4041e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (2D)\n");CHKERRQ(ierr); 4051e07b27eSBarry Smith /*ierr = DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 4061e07b27eSBarry Smith nz = 1; 4071e07b27eSBarry Smith bz = DM_BOUNDARY_NONE; 4081e07b27eSBarry Smith break; 4091e07b27eSBarry Smith case 3: 4101e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (3D)\n");CHKERRQ(ierr); 4111e07b27eSBarry Smith /*ierr = DMDACreate3d(subcomm,bx,by,bz,stencil,nx,ny,nz, PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 4121e07b27eSBarry Smith break; 4131e07b27eSBarry Smith } 4141e07b27eSBarry Smith /* 4151e07b27eSBarry Smith The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append 4161e07b27eSBarry Smith a unique option prefix for the DM, thus I prefer to expose the contents of these API's here. 4171e07b27eSBarry Smith This allows users to control the partitioning of the subDM. 4181e07b27eSBarry Smith */ 4191e07b27eSBarry Smith ierr = DMDACreate(subcomm,&ctx->dmrepart);CHKERRQ(ierr); 4201e07b27eSBarry Smith /* Set unique option prefix name */ 4217c5279cbSDave May ierr = KSPGetOptionsPrefix(sred->ksp,&prefix);CHKERRQ(ierr); 4221e07b27eSBarry Smith ierr = DMSetOptionsPrefix(ctx->dmrepart,prefix);CHKERRQ(ierr); 4231e07b27eSBarry Smith ierr = DMAppendOptionsPrefix(ctx->dmrepart,"repart_");CHKERRQ(ierr); 4241e07b27eSBarry Smith /* standard setup from DMDACreate{1,2,3}d() */ 4251e07b27eSBarry Smith ierr = DMSetDimension(ctx->dmrepart,dim);CHKERRQ(ierr); 4261e07b27eSBarry Smith ierr = DMDASetSizes(ctx->dmrepart,nx,ny,nz);CHKERRQ(ierr); 4271e07b27eSBarry Smith ierr = DMDASetNumProcs(ctx->dmrepart,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 4281e07b27eSBarry Smith ierr = DMDASetBoundaryType(ctx->dmrepart,bx,by,bz);CHKERRQ(ierr); 4291e07b27eSBarry Smith ierr = DMDASetDof(ctx->dmrepart,ndof);CHKERRQ(ierr); 4301e07b27eSBarry Smith ierr = DMDASetStencilType(ctx->dmrepart,stencil);CHKERRQ(ierr); 4311e07b27eSBarry Smith ierr = DMDASetStencilWidth(ctx->dmrepart,nsw);CHKERRQ(ierr); 4321e07b27eSBarry Smith ierr = DMDASetOwnershipRanges(ctx->dmrepart,NULL,NULL,NULL);CHKERRQ(ierr); 4331e07b27eSBarry Smith ierr = DMSetFromOptions(ctx->dmrepart);CHKERRQ(ierr); 4341e07b27eSBarry Smith ierr = DMSetUp(ctx->dmrepart);CHKERRQ(ierr); 4351e07b27eSBarry Smith /* Set refinement factors and interpolation type from the partent */ 4361e07b27eSBarry Smith ierr = DMDASetRefinementFactor(ctx->dmrepart,refine_x,refine_y,refine_z);CHKERRQ(ierr); 4371e07b27eSBarry Smith ierr = DMDASetInterpolationType(ctx->dmrepart,itype);CHKERRQ(ierr); 4381e07b27eSBarry Smith 4391e07b27eSBarry Smith ierr = DMDAGetInfo(ctx->dmrepart,NULL,NULL,NULL,NULL,&ctx->Mp_re,&ctx->Np_re,&ctx->Pp_re,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 4401e07b27eSBarry Smith ierr = DMDAGetOwnershipRanges(ctx->dmrepart,&_range_i_re,&_range_j_re,&_range_k_re);CHKERRQ(ierr); 4415e897e82SDave May 4425e897e82SDave May ctx->dmrepart->ops->creatematrix = dm->ops->creatematrix; 4435e897e82SDave May ctx->dmrepart->ops->createdomaindecomposition = dm->ops->createdomaindecomposition; 4441e07b27eSBarry Smith } 4451e07b27eSBarry Smith 4461e07b27eSBarry Smith /* generate ranges for repartitioned dm */ 4471e07b27eSBarry Smith /* note - assume rank 0 always participates */ 448*071fcb05SBarry Smith /* TODO: use a single MPI call */ 4491e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Mp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4501e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Np_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4511e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Pp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4521e07b27eSBarry Smith 453*071fcb05SBarry Smith ierr = PetscCalloc3(ctx->Mp_re,&ctx->range_i_re,ctx->Np_re,&ctx->range_j_re,ctx->Pp_re,&ctx->range_k_re);CHKERRQ(ierr); 4541e07b27eSBarry Smith 455580bdb30SBarry Smith if (_range_i_re) {ierr = PetscArraycpy(ctx->range_i_re,_range_i_re,ctx->Mp_re);CHKERRQ(ierr);} 456580bdb30SBarry Smith if (_range_j_re) {ierr = PetscArraycpy(ctx->range_j_re,_range_j_re,ctx->Np_re);CHKERRQ(ierr);} 457580bdb30SBarry Smith if (_range_k_re) {ierr = PetscArraycpy(ctx->range_k_re,_range_k_re,ctx->Pp_re);CHKERRQ(ierr);} 4581e07b27eSBarry Smith 459*071fcb05SBarry Smith /* TODO: use a single MPI call */ 4601e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4611e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4621e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4631e07b27eSBarry Smith 464*071fcb05SBarry Smith ierr = PetscMalloc3(ctx->Mp_re,&ctx->start_i_re,ctx->Np_re,&ctx->start_j_re,ctx->Pp_re,&ctx->start_k_re);CHKERRQ(ierr); 4651e07b27eSBarry Smith 4661e07b27eSBarry Smith sum = 0; 4671e07b27eSBarry Smith for (k=0; k<ctx->Mp_re; k++) { 4681e07b27eSBarry Smith ctx->start_i_re[k] = sum; 4691e07b27eSBarry Smith sum += ctx->range_i_re[k]; 4701e07b27eSBarry Smith } 4711e07b27eSBarry Smith 4721e07b27eSBarry Smith sum = 0; 4731e07b27eSBarry Smith for (k=0; k<ctx->Np_re; k++) { 4741e07b27eSBarry Smith ctx->start_j_re[k] = sum; 4751e07b27eSBarry Smith sum += ctx->range_j_re[k]; 4761e07b27eSBarry Smith } 4771e07b27eSBarry Smith 4781e07b27eSBarry Smith sum = 0; 4791e07b27eSBarry Smith for (k=0; k<ctx->Pp_re; k++) { 4801e07b27eSBarry Smith ctx->start_k_re[k] = sum; 4811e07b27eSBarry Smith sum += ctx->range_k_re[k]; 4821e07b27eSBarry Smith } 4831e07b27eSBarry Smith 484ba1c3560SDave May /* attach repartitioned dm to child ksp */ 485ba1c3560SDave May { 486ba1c3560SDave May PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*); 487ba1c3560SDave May void *dmksp_ctx; 488ba1c3560SDave May 489ba1c3560SDave May ierr = DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);CHKERRQ(ierr); 490ba1c3560SDave May 4911e07b27eSBarry Smith /* attach dm to ksp on sub communicator */ 49257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 4931e07b27eSBarry Smith ierr = KSPSetDM(sred->ksp,ctx->dmrepart);CHKERRQ(ierr); 494ba1c3560SDave May 495c5db1f53SDave May if (!dmksp_func || sred->ignore_kspcomputeoperators) { 4961e07b27eSBarry Smith ierr = KSPSetDMActive(sred->ksp,PETSC_FALSE);CHKERRQ(ierr); 497ba1c3560SDave May } else { 498ba1c3560SDave May /* sub ksp inherits dmksp_func and context provided by user */ 499ba1c3560SDave May ierr = KSPSetComputeOperators(sred->ksp,dmksp_func,dmksp_ctx);CHKERRQ(ierr); 500ba1c3560SDave May ierr = KSPSetDMActive(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 501ba1c3560SDave May } 502ba1c3560SDave May } 5031e07b27eSBarry Smith } 5041e07b27eSBarry Smith PetscFunctionReturn(0); 5051e07b27eSBarry Smith } 5061e07b27eSBarry Smith 5071e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 5081e07b27eSBarry Smith { 5091e07b27eSBarry Smith PetscErrorCode ierr; 5101e07b27eSBarry Smith DM dm; 5111e07b27eSBarry Smith MPI_Comm comm; 5121e07b27eSBarry Smith Mat Pscalar,P; 5131e07b27eSBarry Smith PetscInt ndof; 5141e07b27eSBarry Smith PetscInt i,j,k,location,startI[3],endI[3],lenI[3],nx,ny,nz; 5151e07b27eSBarry Smith PetscInt sr,er,Mr; 5161e07b27eSBarry Smith Vec V; 5171e07b27eSBarry Smith 5181e07b27eSBarry Smith PetscFunctionBegin; 5191e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-3D)\n");CHKERRQ(ierr); 5201e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5211e07b27eSBarry Smith 5221e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 5231e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 5241e07b27eSBarry Smith 5251e07b27eSBarry Smith ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 5261e07b27eSBarry Smith ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 5271e07b27eSBarry Smith ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 5281e07b27eSBarry Smith ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 5291e07b27eSBarry Smith sr = sr / ndof; 5301e07b27eSBarry Smith er = er / ndof; 5311e07b27eSBarry Smith Mr = Mr / ndof; 5321e07b27eSBarry Smith 5331e07b27eSBarry Smith ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 5341e07b27eSBarry Smith ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 5351e07b27eSBarry Smith ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 536308c2a70SDave May ierr = MatSeqAIJSetPreallocation(Pscalar,1,NULL);CHKERRQ(ierr); 537308c2a70SDave May ierr = MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);CHKERRQ(ierr); 5381e07b27eSBarry Smith 5391e07b27eSBarry Smith ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]);CHKERRQ(ierr); 5401e07b27eSBarry Smith ierr = DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]);CHKERRQ(ierr); 5411e07b27eSBarry Smith endI[0] += startI[0]; 5421e07b27eSBarry Smith endI[1] += startI[1]; 5431e07b27eSBarry Smith endI[2] += startI[2]; 5441e07b27eSBarry Smith 5451e07b27eSBarry Smith for (k=startI[2]; k<endI[2]; k++) { 5461e07b27eSBarry Smith for (j=startI[1]; j<endI[1]; j++) { 5471e07b27eSBarry Smith for (i=startI[0]; i<endI[0]; i++) { 5481e07b27eSBarry Smith PetscMPIInt rank_ijk_re,rank_reI[3]; 5491e07b27eSBarry Smith PetscInt s0_re; 550c6a0d831SBarry Smith PetscInt ii,jj,kk,local_ijk_re,mapped_ijk; 5511e07b27eSBarry Smith PetscInt lenI_re[3]; 5521e07b27eSBarry Smith 5531e07b27eSBarry Smith location = (i - startI[0]) + (j - startI[1])*lenI[0] + (k - startI[2])*lenI[0]*lenI[1]; 5541e07b27eSBarry Smith ierr = _DMDADetermineRankFromGlobalIJK(3,i,j,k, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 5551e07b27eSBarry Smith ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 5561e07b27eSBarry Smith ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 5571e07b27eSBarry Smith &rank_reI[0],&rank_reI[1],&rank_reI[2],&rank_ijk_re);CHKERRQ(ierr); 5581e07b27eSBarry Smith ierr = _DMDADetermineGlobalS0(3,rank_ijk_re, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, &s0_re);CHKERRQ(ierr); 5591e07b27eSBarry Smith ii = i - ctx->start_i_re[ rank_reI[0] ]; 5601e07b27eSBarry Smith if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error ii"); 5611e07b27eSBarry Smith jj = j - ctx->start_j_re[ rank_reI[1] ]; 5621e07b27eSBarry Smith if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error jj"); 5631e07b27eSBarry Smith kk = k - ctx->start_k_re[ rank_reI[2] ]; 5641e07b27eSBarry Smith if (kk < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error kk"); 5651e07b27eSBarry Smith lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 5661e07b27eSBarry Smith lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 5671e07b27eSBarry Smith lenI_re[2] = ctx->range_k_re[ rank_reI[2] ]; 5681e07b27eSBarry Smith local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1]; 5691e07b27eSBarry Smith mapped_ijk = s0_re + local_ijk_re; 5701e07b27eSBarry Smith ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 5711e07b27eSBarry Smith } 5721e07b27eSBarry Smith } 5731e07b27eSBarry Smith } 5741e07b27eSBarry Smith ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5751e07b27eSBarry Smith ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5761e07b27eSBarry Smith ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 5771e07b27eSBarry Smith ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 5781e07b27eSBarry Smith ctx->permutation = P; 5791e07b27eSBarry Smith PetscFunctionReturn(0); 5801e07b27eSBarry Smith } 5811e07b27eSBarry Smith 5821e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 5831e07b27eSBarry Smith { 5841e07b27eSBarry Smith PetscErrorCode ierr; 5851e07b27eSBarry Smith DM dm; 5861e07b27eSBarry Smith MPI_Comm comm; 5871e07b27eSBarry Smith Mat Pscalar,P; 5881e07b27eSBarry Smith PetscInt ndof; 5891e07b27eSBarry Smith PetscInt i,j,location,startI[2],endI[2],lenI[2],nx,ny,nz; 5901e07b27eSBarry Smith PetscInt sr,er,Mr; 5911e07b27eSBarry Smith Vec V; 5921e07b27eSBarry Smith 5931e07b27eSBarry Smith PetscFunctionBegin; 5941e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-2D)\n");CHKERRQ(ierr); 5951e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5961e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 5971e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 5981e07b27eSBarry Smith ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 5991e07b27eSBarry Smith ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 6001e07b27eSBarry Smith ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 6011e07b27eSBarry Smith ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 6021e07b27eSBarry Smith sr = sr / ndof; 6031e07b27eSBarry Smith er = er / ndof; 6041e07b27eSBarry Smith Mr = Mr / ndof; 6051e07b27eSBarry Smith 6061e07b27eSBarry Smith ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 6071e07b27eSBarry Smith ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 6081e07b27eSBarry Smith ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 609308c2a70SDave May ierr = MatSeqAIJSetPreallocation(Pscalar,1,NULL);CHKERRQ(ierr); 610308c2a70SDave May ierr = MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);CHKERRQ(ierr); 6111e07b27eSBarry Smith 6121e07b27eSBarry Smith ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL);CHKERRQ(ierr); 6131e07b27eSBarry Smith ierr = DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL);CHKERRQ(ierr); 6141e07b27eSBarry Smith endI[0] += startI[0]; 6151e07b27eSBarry Smith endI[1] += startI[1]; 6161e07b27eSBarry Smith 6171e07b27eSBarry Smith for (j=startI[1]; j<endI[1]; j++) { 6181e07b27eSBarry Smith for (i=startI[0]; i<endI[0]; i++) { 6191e07b27eSBarry Smith PetscMPIInt rank_ijk_re,rank_reI[3]; 6201e07b27eSBarry Smith PetscInt s0_re; 621c6a0d831SBarry Smith PetscInt ii,jj,local_ijk_re,mapped_ijk; 6221e07b27eSBarry Smith PetscInt lenI_re[3]; 6231e07b27eSBarry Smith 6241e07b27eSBarry Smith location = (i - startI[0]) + (j - startI[1])*lenI[0]; 6251e07b27eSBarry Smith ierr = _DMDADetermineRankFromGlobalIJK(2,i,j,0, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 6261e07b27eSBarry Smith ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 6271e07b27eSBarry Smith ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 6281e07b27eSBarry Smith &rank_reI[0],&rank_reI[1],NULL,&rank_ijk_re);CHKERRQ(ierr); 6291e07b27eSBarry Smith 6301e07b27eSBarry Smith ierr = _DMDADetermineGlobalS0(2,rank_ijk_re, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, &s0_re);CHKERRQ(ierr); 6311e07b27eSBarry Smith 6321e07b27eSBarry Smith ii = i - ctx->start_i_re[ rank_reI[0] ]; 6331e07b27eSBarry Smith if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error ii"); 6341e07b27eSBarry Smith jj = j - ctx->start_j_re[ rank_reI[1] ]; 6351e07b27eSBarry Smith if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error jj"); 6361e07b27eSBarry Smith 6371e07b27eSBarry Smith lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 6381e07b27eSBarry Smith lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 6391e07b27eSBarry Smith local_ijk_re = ii + jj * lenI_re[0]; 6401e07b27eSBarry Smith mapped_ijk = s0_re + local_ijk_re; 6411e07b27eSBarry Smith ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 6421e07b27eSBarry Smith } 6431e07b27eSBarry Smith } 6441e07b27eSBarry Smith ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6451e07b27eSBarry Smith ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6461e07b27eSBarry Smith ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 6471e07b27eSBarry Smith ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 6481e07b27eSBarry Smith ctx->permutation = P; 6491e07b27eSBarry Smith PetscFunctionReturn(0); 6501e07b27eSBarry Smith } 6511e07b27eSBarry Smith 6521e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 6531e07b27eSBarry Smith { 6541e07b27eSBarry Smith PetscErrorCode ierr; 6551e07b27eSBarry Smith Vec xred,yred,xtmp,x,xp; 6561e07b27eSBarry Smith VecScatter scatter; 6571e07b27eSBarry Smith IS isin; 6581e07b27eSBarry Smith Mat B; 6591e07b27eSBarry Smith PetscInt m,bs,st,ed; 6601e07b27eSBarry Smith MPI_Comm comm; 6611e07b27eSBarry Smith 6621e07b27eSBarry Smith PetscFunctionBegin; 6631e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 6641e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 6651e07b27eSBarry Smith ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 6661e07b27eSBarry Smith ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 6671e07b27eSBarry Smith ierr = VecDuplicate(x,&xp);CHKERRQ(ierr); 6683ac26c5eSBarry Smith m = 0; 6691e07b27eSBarry Smith xred = NULL; 6701e07b27eSBarry Smith yred = NULL; 67157f12427SDave May if (PCTelescope_isActiveRank(sred)) { 6721e07b27eSBarry Smith ierr = DMCreateGlobalVector(ctx->dmrepart,&xred);CHKERRQ(ierr); 6731e07b27eSBarry Smith ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 6741e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 6751e07b27eSBarry Smith ierr = ISCreateStride(comm,ed-st,st,1,&isin);CHKERRQ(ierr); 676ca43db0aSBarry Smith ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr); 6771e07b27eSBarry Smith } else { 6781e07b27eSBarry Smith ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 6793ac26c5eSBarry Smith ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 6801e07b27eSBarry Smith } 6811e07b27eSBarry Smith ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 6821e07b27eSBarry Smith ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 6831e07b27eSBarry Smith ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 6841e07b27eSBarry Smith ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 6851e07b27eSBarry Smith ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 6869448b7f1SJunchao Zhang ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 6871e07b27eSBarry Smith sred->xred = xred; 6881e07b27eSBarry Smith sred->yred = yred; 6891e07b27eSBarry Smith sred->isin = isin; 6901e07b27eSBarry Smith sred->scatter = scatter; 6911e07b27eSBarry Smith sred->xtmp = xtmp; 6921e07b27eSBarry Smith 6931e07b27eSBarry Smith ctx->xp = xp; 6941e07b27eSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 6951e07b27eSBarry Smith PetscFunctionReturn(0); 6961e07b27eSBarry Smith } 6971e07b27eSBarry Smith 6981e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda(PC pc,PC_Telescope sred) 6991e07b27eSBarry Smith { 7001e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 7011e07b27eSBarry Smith PetscInt dim; 7021e07b27eSBarry Smith DM dm; 7031e07b27eSBarry Smith MPI_Comm comm; 7041e07b27eSBarry Smith PetscErrorCode ierr; 7051e07b27eSBarry Smith 7061e07b27eSBarry Smith PetscFunctionBegin; 7071e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setup (DMDA)\n");CHKERRQ(ierr); 708580bdb30SBarry Smith ierr = PetscNew(&ctx);CHKERRQ(ierr); 7091e07b27eSBarry Smith sred->dm_ctx = (void*)ctx; 7101e07b27eSBarry Smith 7111e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7121e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 7131e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 7141e07b27eSBarry Smith 7151e07b27eSBarry Smith PCTelescopeSetUp_dmda_repart(pc,sred,ctx); 7161e07b27eSBarry Smith PCTelescopeSetUp_dmda_repart_coors(pc,sred,ctx); 7171e07b27eSBarry Smith switch (dim) { 7181e07b27eSBarry Smith case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 7191e07b27eSBarry Smith break; 7201e07b27eSBarry Smith case 2: ierr = PCTelescopeSetUp_dmda_permutation_2d(pc,sred,ctx);CHKERRQ(ierr); 7211e07b27eSBarry Smith break; 7221e07b27eSBarry Smith case 3: ierr = PCTelescopeSetUp_dmda_permutation_3d(pc,sred,ctx);CHKERRQ(ierr); 7231e07b27eSBarry Smith break; 7241e07b27eSBarry Smith } 7251e07b27eSBarry Smith ierr = PCTelescopeSetUp_dmda_scatters(pc,sred,ctx);CHKERRQ(ierr); 7261e07b27eSBarry Smith PetscFunctionReturn(0); 7271e07b27eSBarry Smith } 7281e07b27eSBarry Smith 729ba1c3560SDave May PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 7301e07b27eSBarry Smith { 7311e07b27eSBarry Smith PetscErrorCode ierr; 7321e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 7331e07b27eSBarry Smith MPI_Comm comm,subcomm; 7341e07b27eSBarry Smith Mat Bperm,Bred,B,P; 7351e07b27eSBarry Smith PetscInt nr,nc; 7361e07b27eSBarry Smith IS isrow,iscol; 7371e07b27eSBarry Smith Mat Blocal,*_Blocal; 7381e07b27eSBarry Smith 7391e07b27eSBarry Smith PetscFunctionBegin; 7401e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (DMDA)\n");CHKERRQ(ierr); 7411e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7421e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 7431e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 7441e07b27eSBarry Smith 7451e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 7461e07b27eSBarry Smith ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 7471e07b27eSBarry Smith 7481e07b27eSBarry Smith P = ctx->permutation; 7491e07b27eSBarry Smith ierr = MatPtAP(B,P,MAT_INITIAL_MATRIX,1.1,&Bperm);CHKERRQ(ierr); 7501e07b27eSBarry Smith 7511e07b27eSBarry Smith /* Get submatrices */ 7521e07b27eSBarry Smith isrow = sred->isin; 7531e07b27eSBarry Smith ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 7541e07b27eSBarry Smith 7557dae84e0SHong Zhang ierr = MatCreateSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 7561e07b27eSBarry Smith Blocal = *_Blocal; 7571e07b27eSBarry Smith Bred = NULL; 75857f12427SDave May if (PCTelescope_isActiveRank(sred)) { 7591e07b27eSBarry Smith PetscInt mm; 7601e07b27eSBarry Smith 7611e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;} 7621e07b27eSBarry Smith ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 763bfd6bcc6SSatish Balay /* ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred);CHKERRQ(ierr); */ 7641e07b27eSBarry Smith ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 7651e07b27eSBarry Smith } 7661e07b27eSBarry Smith *A = Bred; 7671e07b27eSBarry Smith 7681e07b27eSBarry Smith ierr = ISDestroy(&iscol);CHKERRQ(ierr); 7691e07b27eSBarry Smith ierr = MatDestroy(&Bperm);CHKERRQ(ierr); 7701e07b27eSBarry Smith ierr = MatDestroyMatrices(1,&_Blocal);CHKERRQ(ierr); 7711e07b27eSBarry Smith PetscFunctionReturn(0); 7721e07b27eSBarry Smith } 7731e07b27eSBarry Smith 774ba1c3560SDave May PetscErrorCode PCTelescopeMatCreate_dmda(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 775ba1c3560SDave May { 776ba1c3560SDave May PetscErrorCode ierr; 777ba1c3560SDave May DM dm; 778ba1c3560SDave May PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*); 779ba1c3560SDave May void *dmksp_ctx; 780ba1c3560SDave May 781ba1c3560SDave May PetscFunctionBegin; 782ba1c3560SDave May ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 783ba1c3560SDave May ierr = DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);CHKERRQ(ierr); 784dc9ee9fdSDave May /* We assume that dmksp_func = NULL, is equivalent to dmActive = PETSC_FALSE */ 7857c5279cbSDave May if (dmksp_func && !sred->ignore_kspcomputeoperators) { 786ba1c3560SDave May DM dmrepart; 78728323a89SDave May Mat Ak; 788ba1c3560SDave May 789ba1c3560SDave May *A = NULL; 79057f12427SDave May if (PCTelescope_isActiveRank(sred)) { 791ba1c3560SDave May ierr = KSPGetDM(sred->ksp,&dmrepart);CHKERRQ(ierr); 792ba1c3560SDave May if (reuse == MAT_INITIAL_MATRIX) { 793ba1c3560SDave May ierr = DMCreateMatrix(dmrepart,&Ak);CHKERRQ(ierr); 794ba1c3560SDave May *A = Ak; 795ba1c3560SDave May } else if (reuse == MAT_REUSE_MATRIX) { 796ba1c3560SDave May Ak = *A; 797ba1c3560SDave May } 7985c5dbb1cSDave May /* 7995c5dbb1cSDave May There is no need to explicitly assemble the operator now, 8005c5dbb1cSDave May the sub-KSP will call the method provided to KSPSetComputeOperators() during KSPSetUp() 8015c5dbb1cSDave May */ 802ba1c3560SDave May } 803ba1c3560SDave May } else { 804ba1c3560SDave May ierr = PCTelescopeMatCreate_dmda_dmactivefalse(pc,sred,reuse,A);CHKERRQ(ierr); 805ba1c3560SDave May } 806ba1c3560SDave May PetscFunctionReturn(0); 807ba1c3560SDave May } 808ba1c3560SDave May 809392968a1SPatrick Sanan PetscErrorCode PCTelescopeSubNullSpaceCreate_dmda_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 8101e07b27eSBarry Smith { 8111e07b27eSBarry Smith PetscErrorCode ierr; 8121e07b27eSBarry Smith PetscBool has_const; 813a947c41eSDave May PetscInt i,k,n = 0; 8141e07b27eSBarry Smith const Vec *vecs; 815c41e779fSDave May Vec *sub_vecs = NULL; 8161e07b27eSBarry Smith MPI_Comm subcomm; 8171e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 8181e07b27eSBarry Smith 8191e07b27eSBarry Smith PetscFunctionBegin; 8201e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 8211e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 8221e07b27eSBarry Smith ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 8231e07b27eSBarry Smith 82457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 8251e07b27eSBarry Smith /* create new vectors */ 826e3acf2f7SBarry Smith if (n) { 827e3acf2f7SBarry Smith ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 8281e07b27eSBarry Smith } 8291e07b27eSBarry Smith } 8301e07b27eSBarry Smith 8311e07b27eSBarry Smith /* copy entries */ 8321e07b27eSBarry Smith for (k=0; k<n; k++) { 8331e07b27eSBarry Smith const PetscScalar *x_array; 8341e07b27eSBarry Smith PetscScalar *LA_sub_vec; 83513c30530SDave May PetscInt st,ed; 8361e07b27eSBarry Smith 8371e07b27eSBarry Smith /* permute vector into ordering associated with re-partitioned dmda */ 8381e07b27eSBarry Smith ierr = MatMultTranspose(ctx->permutation,vecs[k],ctx->xp);CHKERRQ(ierr); 8391e07b27eSBarry Smith 8401e07b27eSBarry Smith /* pull in vector x->xtmp */ 8411e07b27eSBarry Smith ierr = VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8421e07b27eSBarry Smith ierr = VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8431e07b27eSBarry Smith 844392968a1SPatrick Sanan /* copy vector entries into xred */ 8451e07b27eSBarry Smith ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 846ea2b237eSDave May if (sub_vecs) { 847ea2b237eSDave May if (sub_vecs[k]) { 8481e07b27eSBarry Smith ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 8491e07b27eSBarry Smith ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 8501e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 8511e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 8521e07b27eSBarry Smith } 8531e07b27eSBarry Smith ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 8541e07b27eSBarry Smith } 855ea2b237eSDave May } 8561e07b27eSBarry Smith ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 8571e07b27eSBarry Smith } 8581e07b27eSBarry Smith 85957f12427SDave May if (PCTelescope_isActiveRank(sred)) { 860d8b9d5b7SPatrick Sanan /* create new (near) nullspace for redundant object */ 861392968a1SPatrick Sanan ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr); 862392968a1SPatrick Sanan ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 863d8b9d5b7SPatrick Sanan if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 864d8b9d5b7SPatrick 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"); 865d8b9d5b7SPatrick Sanan } 866392968a1SPatrick Sanan PetscFunctionReturn(0); 867392968a1SPatrick Sanan } 868392968a1SPatrick Sanan 869392968a1SPatrick Sanan PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat) 870392968a1SPatrick Sanan { 871392968a1SPatrick Sanan PetscErrorCode ierr; 872392968a1SPatrick Sanan Mat B; 873392968a1SPatrick Sanan 874392968a1SPatrick Sanan PetscFunctionBegin; 875392968a1SPatrick Sanan ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 876392968a1SPatrick Sanan { 877392968a1SPatrick Sanan MatNullSpace nullspace,sub_nullspace; 878392968a1SPatrick Sanan ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 879392968a1SPatrick Sanan if (nullspace) { 880392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n");CHKERRQ(ierr); 881392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr); 88257f12427SDave May if (PCTelescope_isActiveRank(sred)) { 883392968a1SPatrick Sanan ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 88441ff1ee9SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr); 8851e07b27eSBarry Smith } 886392968a1SPatrick Sanan } 887392968a1SPatrick Sanan } 888392968a1SPatrick Sanan { 889392968a1SPatrick Sanan MatNullSpace nearnullspace,sub_nearnullspace; 89064f1816dSDave May ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr); 891392968a1SPatrick Sanan if (nearnullspace) { 892392968a1SPatrick Sanan ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (DMDA)\n");CHKERRQ(ierr); 893392968a1SPatrick Sanan ierr = PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr); 89457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 89564f1816dSDave May ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr); 896392968a1SPatrick Sanan ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr); 897392968a1SPatrick Sanan } 898392968a1SPatrick Sanan } 899392968a1SPatrick Sanan } 9001e07b27eSBarry Smith PetscFunctionReturn(0); 9011e07b27eSBarry Smith } 9021e07b27eSBarry Smith 9031e07b27eSBarry Smith PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y) 9041e07b27eSBarry Smith { 9051e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 9061e07b27eSBarry Smith PetscErrorCode ierr; 9071e07b27eSBarry Smith Mat perm; 9081e07b27eSBarry Smith Vec xtmp,xp,xred,yred; 90913c30530SDave May PetscInt i,st,ed; 9101e07b27eSBarry Smith VecScatter scatter; 9111e07b27eSBarry Smith PetscScalar *array; 9121e07b27eSBarry Smith const PetscScalar *x_array; 9131e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 9141e07b27eSBarry Smith 9151e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 9161e07b27eSBarry Smith xtmp = sred->xtmp; 9171e07b27eSBarry Smith scatter = sred->scatter; 9181e07b27eSBarry Smith xred = sred->xred; 9191e07b27eSBarry Smith yred = sred->yred; 9201e07b27eSBarry Smith perm = ctx->permutation; 9211e07b27eSBarry Smith xp = ctx->xp; 9221e07b27eSBarry Smith 9231e07b27eSBarry Smith PetscFunctionBegin; 924bf00f589SPatrick Sanan ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 92514c9fce5SDave May 9261e07b27eSBarry Smith /* permute vector into ordering associated with re-partitioned dmda */ 9271e07b27eSBarry Smith ierr = MatMultTranspose(perm,x,xp);CHKERRQ(ierr); 9281e07b27eSBarry Smith 9291e07b27eSBarry Smith /* pull in vector x->xtmp */ 9301e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9311e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9321e07b27eSBarry Smith 9331e07b27eSBarry Smith /* copy vector entires into xred */ 9341e07b27eSBarry Smith ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 9351e07b27eSBarry Smith if (xred) { 9361e07b27eSBarry Smith PetscScalar *LA_xred; 9371e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 9381e07b27eSBarry Smith 9391e07b27eSBarry Smith ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 9401e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 9411e07b27eSBarry Smith LA_xred[i] = x_array[i]; 9421e07b27eSBarry Smith } 9431e07b27eSBarry Smith ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 9441e07b27eSBarry Smith } 9451e07b27eSBarry Smith ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 9461e07b27eSBarry Smith 9471e07b27eSBarry Smith /* solve */ 94857f12427SDave May if (PCTelescope_isActiveRank(sred)) { 9491e07b27eSBarry Smith ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 950c0decd05SBarry Smith ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr); 9511e07b27eSBarry Smith } 9521e07b27eSBarry Smith 9531e07b27eSBarry Smith /* return vector */ 9541e07b27eSBarry Smith ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 9551e07b27eSBarry Smith if (yred) { 9561e07b27eSBarry Smith const PetscScalar *LA_yred; 9571e07b27eSBarry Smith ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 9581e07b27eSBarry Smith ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 9591e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 9601e07b27eSBarry Smith array[i] = LA_yred[i]; 9611e07b27eSBarry Smith } 9621e07b27eSBarry Smith ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 9631e07b27eSBarry Smith } 9641e07b27eSBarry Smith ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 9651e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9661e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9671e07b27eSBarry Smith ierr = MatMult(perm,xp,y);CHKERRQ(ierr); 9681e07b27eSBarry Smith PetscFunctionReturn(0); 9691e07b27eSBarry Smith } 9701e07b27eSBarry Smith 971f650675bSDave May PetscErrorCode PCApplyRichardson_Telescope_dmda(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 972f650675bSDave May { 973f650675bSDave May PC_Telescope sred = (PC_Telescope)pc->data; 974f650675bSDave May PetscErrorCode ierr; 975f650675bSDave May Mat perm; 976a1d91a28SDave May Vec xtmp,xp,yred; 977f650675bSDave May PetscInt i,st,ed; 978f650675bSDave May VecScatter scatter; 979f650675bSDave May const PetscScalar *x_array; 980c41e779fSDave May PetscBool default_init_guess_value = PETSC_FALSE; 981f650675bSDave May PC_Telescope_DMDACtx *ctx; 982f650675bSDave May 98357f12427SDave May PetscFunctionBegin; 984f650675bSDave May ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 985f650675bSDave May xtmp = sred->xtmp; 986f650675bSDave May scatter = sred->scatter; 987f650675bSDave May yred = sred->yred; 988f650675bSDave May perm = ctx->permutation; 989f650675bSDave May xp = ctx->xp; 990f650675bSDave May 991f650675bSDave May if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope_dmda only supports max_it = 1"); 992f650675bSDave May *reason = (PCRichardsonConvergedReason)0; 993f650675bSDave May 994f650675bSDave May if (!zeroguess) { 995f650675bSDave May ierr = PetscInfo(pc,"PCTelescopeDMDA: Scattering y for non-zero-initial guess\n");CHKERRQ(ierr); 996f650675bSDave May /* permute vector into ordering associated with re-partitioned dmda */ 997f650675bSDave May ierr = MatMultTranspose(perm,y,xp);CHKERRQ(ierr); 998f650675bSDave May 999f650675bSDave May /* pull in vector x->xtmp */ 1000f650675bSDave May ierr = VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1001f650675bSDave May ierr = VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1002f650675bSDave May 1003f650675bSDave May /* copy vector entires into xred */ 1004f650675bSDave May ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 1005f650675bSDave May if (yred) { 1006f650675bSDave May PetscScalar *LA_yred; 1007f650675bSDave May ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 1008f650675bSDave May ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 1009f650675bSDave May for (i=0; i<ed-st; i++) { 1010f650675bSDave May LA_yred[i] = x_array[i]; 1011f650675bSDave May } 1012f650675bSDave May ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 1013f650675bSDave May } 1014f650675bSDave May ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 1015f650675bSDave May } 1016f650675bSDave May 101757f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1018f650675bSDave May ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr); 1019f650675bSDave May if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 1020f650675bSDave May } 1021f650675bSDave May 1022f650675bSDave May ierr = PCApply_Telescope_dmda(pc,x,y);CHKERRQ(ierr); 1023f650675bSDave May 102457f12427SDave May if (PCTelescope_isActiveRank(sred)) { 1025f650675bSDave May ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr); 1026f650675bSDave May } 1027f650675bSDave May 1028f650675bSDave May if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1029f650675bSDave May *outits = 1; 1030f650675bSDave May PetscFunctionReturn(0); 1031f650675bSDave May } 1032f650675bSDave May 10331e07b27eSBarry Smith PetscErrorCode PCReset_Telescope_dmda(PC pc) 10341e07b27eSBarry Smith { 10351e07b27eSBarry Smith PetscErrorCode ierr; 10361e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 10371e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 10381e07b27eSBarry Smith 10391e07b27eSBarry Smith PetscFunctionBegin; 10401e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 10411e07b27eSBarry Smith ierr = VecDestroy(&ctx->xp);CHKERRQ(ierr); 10421e07b27eSBarry Smith ierr = MatDestroy(&ctx->permutation);CHKERRQ(ierr); 10431e07b27eSBarry Smith ierr = DMDestroy(&ctx->dmrepart);CHKERRQ(ierr); 1044*071fcb05SBarry Smith ierr = PetscFree3(ctx->range_i_re,ctx->range_j_re,ctx->range_k_re);CHKERRQ(ierr); 1045*071fcb05SBarry Smith ierr = PetscFree3(ctx->start_i_re,ctx->start_j_re,ctx->start_k_re);CHKERRQ(ierr); 10461e07b27eSBarry Smith PetscFunctionReturn(0); 10471e07b27eSBarry Smith } 10481e07b27eSBarry Smith 10498ef9ca65SPatrick Sanan PetscErrorCode DMView_DA_Short_3d(DM dm,PetscViewer v) 10501e07b27eSBarry Smith { 10511e07b27eSBarry Smith PetscInt M,N,P,m,n,p,ndof,nsw; 10521e07b27eSBarry Smith MPI_Comm comm; 10531e07b27eSBarry Smith PetscMPIInt size; 10541e07b27eSBarry Smith const char* prefix; 10551e07b27eSBarry Smith PetscErrorCode ierr; 10561e07b27eSBarry Smith 10571e07b27eSBarry Smith PetscFunctionBegin; 10581e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 10591e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 10601e07b27eSBarry Smith ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 10611e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 10621e07b27eSBarry Smith if (prefix) {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 10631e07b27eSBarry Smith else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 10641e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(v," M %D N %D P %D m %D n %D p %D dof %D overlap %D\n",M,N,P,m,n,p,ndof,nsw);CHKERRQ(ierr); 10651e07b27eSBarry Smith PetscFunctionReturn(0); 10661e07b27eSBarry Smith } 10671e07b27eSBarry Smith 10688ef9ca65SPatrick Sanan PetscErrorCode DMView_DA_Short_2d(DM dm,PetscViewer v) 10691e07b27eSBarry Smith { 10701e07b27eSBarry Smith PetscInt M,N,m,n,ndof,nsw; 10711e07b27eSBarry Smith MPI_Comm comm; 10721e07b27eSBarry Smith PetscMPIInt size; 10731e07b27eSBarry Smith const char* prefix; 10741e07b27eSBarry Smith PetscErrorCode ierr; 10751e07b27eSBarry Smith 10761e07b27eSBarry Smith PetscFunctionBegin; 10771e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 10781e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 10791e07b27eSBarry Smith ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 10801e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 10811e07b27eSBarry Smith if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 10821e07b27eSBarry Smith else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 10831e07b27eSBarry Smith ierr = PetscViewerASCIIPrintf(v," M %D N %D m %D n %D dof %D overlap %D\n",M,N,m,n,ndof,nsw);CHKERRQ(ierr); 10841e07b27eSBarry Smith PetscFunctionReturn(0); 10851e07b27eSBarry Smith } 10861e07b27eSBarry Smith 10878ef9ca65SPatrick Sanan PetscErrorCode DMView_DA_Short(DM dm,PetscViewer v) 10881e07b27eSBarry Smith { 10891e07b27eSBarry Smith PetscErrorCode ierr; 10901e07b27eSBarry Smith PetscInt dim; 10911e07b27eSBarry Smith 10921e07b27eSBarry Smith PetscFunctionBegin; 10931e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 10941e07b27eSBarry Smith switch (dim) { 10958ef9ca65SPatrick Sanan case 2: ierr = DMView_DA_Short_2d(dm,v);CHKERRQ(ierr); 10961e07b27eSBarry Smith break; 10978ef9ca65SPatrick Sanan case 3: ierr = DMView_DA_Short_3d(dm,v);CHKERRQ(ierr); 10981e07b27eSBarry Smith break; 10991e07b27eSBarry Smith } 11001e07b27eSBarry Smith PetscFunctionReturn(0); 11011e07b27eSBarry Smith } 11021e07b27eSBarry Smith 1103