11e07b27eSBarry Smith 21e07b27eSBarry Smith #include <petsc/private/pcimpl.h> 31e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/ 41e07b27eSBarry Smith #include <petscdm.h> 51e07b27eSBarry Smith #include <petscdmda.h> 61e07b27eSBarry Smith 71e07b27eSBarry Smith #include "telescope.h" 81e07b27eSBarry Smith 91e07b27eSBarry Smith #undef __FUNCT__ 101e07b27eSBarry Smith #define __FUNCT__ "_DMDADetermineRankFromGlobalIJK" 111e07b27eSBarry Smith PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim,PetscInt i,PetscInt j,PetscInt k, 121e07b27eSBarry Smith PetscInt Mp,PetscInt Np,PetscInt Pp, 131e07b27eSBarry Smith PetscInt start_i[],PetscInt start_j[],PetscInt start_k[], 141e07b27eSBarry Smith PetscInt span_i[],PetscInt span_j[],PetscInt span_k[], 151e07b27eSBarry Smith PetscMPIInt *_pi,PetscMPIInt *_pj,PetscMPIInt *_pk,PetscMPIInt *rank_re) 161e07b27eSBarry Smith { 171e07b27eSBarry Smith PetscInt pi,pj,pk,n; 181e07b27eSBarry Smith 191e07b27eSBarry Smith PetscFunctionBegin; 201e07b27eSBarry Smith pi = pj = pk = -1; 211e07b27eSBarry Smith if (_pi) { 221e07b27eSBarry Smith for (n=0; n<Mp; n++) { 231e07b27eSBarry Smith if ( (i >= start_i[n]) && (i < start_i[n]+span_i[n]) ) { 241e07b27eSBarry Smith pi = n; 251e07b27eSBarry Smith break; 261e07b27eSBarry Smith } 271e07b27eSBarry Smith } 281e07b27eSBarry Smith if (pi == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pi cannot be determined : range %D, val %D",Mp,i); 291e07b27eSBarry Smith *_pi = pi; 301e07b27eSBarry Smith } 311e07b27eSBarry Smith 321e07b27eSBarry Smith if (_pj) { 331e07b27eSBarry Smith for (n=0; n<Np; n++) { 341e07b27eSBarry Smith if ( (j >= start_j[n]) && (j < start_j[n]+span_j[n]) ) { 351e07b27eSBarry Smith pj = n; 361e07b27eSBarry Smith break; 371e07b27eSBarry Smith } 381e07b27eSBarry Smith } 391e07b27eSBarry Smith if (pj == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pj cannot be determined : range %D, val %D",Np,j); 401e07b27eSBarry Smith *_pj = pj; 411e07b27eSBarry Smith } 421e07b27eSBarry Smith 431e07b27eSBarry Smith if (_pk) { 441e07b27eSBarry Smith for (n=0; n<Pp; n++) { 451e07b27eSBarry Smith if ( (k >= start_k[n]) && (k < start_k[n]+span_k[n]) ) { 461e07b27eSBarry Smith pk = n; 471e07b27eSBarry Smith break; 481e07b27eSBarry Smith } 491e07b27eSBarry Smith } 501e07b27eSBarry Smith if (pk == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pk cannot be determined : range %D, val %D",Pp,k); 511e07b27eSBarry Smith *_pk = pk; 521e07b27eSBarry Smith } 531e07b27eSBarry Smith 541e07b27eSBarry Smith switch (dim) { 551e07b27eSBarry Smith case 1: 561e07b27eSBarry Smith *rank_re = pi; 571e07b27eSBarry Smith break; 581e07b27eSBarry Smith case 2: 591e07b27eSBarry Smith *rank_re = pi + pj * Mp; 601e07b27eSBarry Smith break; 611e07b27eSBarry Smith case 3: 621e07b27eSBarry Smith *rank_re = pi + pj * Mp + pk * (Mp*Np); 631e07b27eSBarry Smith break; 641e07b27eSBarry Smith } 651e07b27eSBarry Smith PetscFunctionReturn(0); 661e07b27eSBarry Smith } 671e07b27eSBarry Smith 681e07b27eSBarry Smith #undef __FUNCT__ 691e07b27eSBarry Smith #define __FUNCT__ "_DMDADetermineGlobalS0" 701e07b27eSBarry Smith PetscErrorCode _DMDADetermineGlobalS0(PetscInt dim,PetscMPIInt rank_re,PetscInt Mp_re,PetscInt Np_re,PetscInt Pp_re, 711e07b27eSBarry Smith PetscInt range_i_re[],PetscInt range_j_re[],PetscInt range_k_re[],PetscInt *s0) 721e07b27eSBarry Smith { 731e07b27eSBarry Smith PetscInt i,j,k,start_IJK; 741e07b27eSBarry Smith PetscInt rank_ijk; 751e07b27eSBarry Smith 761e07b27eSBarry Smith PetscFunctionBegin; 771e07b27eSBarry Smith switch (dim) { 781e07b27eSBarry Smith case 1: 791e07b27eSBarry Smith start_IJK = 0; 801e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 811e07b27eSBarry Smith rank_ijk = i; 821e07b27eSBarry Smith if (rank_ijk < rank_re) { 831e07b27eSBarry Smith start_IJK += range_i_re[i]; 841e07b27eSBarry Smith } 851e07b27eSBarry Smith } 861e07b27eSBarry Smith break; 871e07b27eSBarry Smith case 2: 881e07b27eSBarry Smith start_IJK = 0; 891e07b27eSBarry Smith for (j=0; j<Np_re; j++) { 901e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 911e07b27eSBarry Smith rank_ijk = i + j*Mp_re; 921e07b27eSBarry Smith if (rank_ijk < rank_re) { 931e07b27eSBarry Smith start_IJK += range_i_re[i]*range_j_re[j]; 941e07b27eSBarry Smith } 951e07b27eSBarry Smith } 961e07b27eSBarry Smith } 971e07b27eSBarry Smith break; 981e07b27eSBarry Smith case 3: 991e07b27eSBarry Smith start_IJK = 0; 1001e07b27eSBarry Smith for (k=0; k<Pp_re; k++) { 1011e07b27eSBarry Smith for (j=0; j<Np_re; j++) { 1021e07b27eSBarry Smith for (i=0; i<Mp_re; i++) { 1031e07b27eSBarry Smith rank_ijk = i + j*Mp_re + k*Mp_re*Np_re; 1041e07b27eSBarry Smith if (rank_ijk < rank_re) { 1051e07b27eSBarry Smith start_IJK += range_i_re[i]*range_j_re[j]*range_k_re[k]; 1061e07b27eSBarry Smith } 1071e07b27eSBarry Smith } 1081e07b27eSBarry Smith } 1091e07b27eSBarry Smith } 1101e07b27eSBarry Smith break; 1111e07b27eSBarry Smith } 1121e07b27eSBarry Smith *s0 = start_IJK; 1131e07b27eSBarry Smith PetscFunctionReturn(0); 1141e07b27eSBarry Smith } 1151e07b27eSBarry Smith 1161e07b27eSBarry Smith #undef __FUNCT__ 1171e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_repart_coors2d" 1181e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PetscSubcomm psubcomm,DM dm,DM subdm) 1191e07b27eSBarry Smith { 1201e07b27eSBarry Smith PetscErrorCode ierr; 1211e07b27eSBarry Smith DM cdm; 1221e07b27eSBarry Smith Vec coor,coor_natural,perm_coors; 1231e07b27eSBarry Smith PetscInt i,j,si,sj,ni,nj,M,N,Ml,Nl,c,nidx; 1241e07b27eSBarry Smith PetscInt *fine_indices; 1251e07b27eSBarry Smith IS is_fine,is_local; 1261e07b27eSBarry Smith VecScatter sctx; 1271e07b27eSBarry Smith 1281e07b27eSBarry Smith PetscFunctionBegin; 1291e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 1301e07b27eSBarry Smith if (!coor) return(0); 1311e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 1321e07b27eSBarry Smith ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1331e07b27eSBarry Smith } 1341e07b27eSBarry Smith /* Get the coordinate vector from the distributed array */ 1351e07b27eSBarry Smith ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 1361e07b27eSBarry Smith ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 1371e07b27eSBarry Smith 1381e07b27eSBarry Smith ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 1391e07b27eSBarry Smith ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 1401e07b27eSBarry Smith 1411e07b27eSBarry Smith /* get indices of the guys I want to grab */ 1421e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1431e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 1441e07b27eSBarry Smith ierr = DMDAGetCorners(subdm,&si,&sj,NULL,&ni,&nj,NULL);CHKERRQ(ierr); 1451e07b27eSBarry Smith 1461e07b27eSBarry Smith Ml = ni; 1471e07b27eSBarry Smith Nl = nj; 1481e07b27eSBarry Smith } else { 1491e07b27eSBarry Smith Ml = Nl = 1; 1501e07b27eSBarry Smith } 1511e07b27eSBarry Smith 1521e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*Ml*Nl*2,&fine_indices);CHKERRQ(ierr); 1531e07b27eSBarry Smith c = 0; 1541e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 1551e07b27eSBarry Smith for (j=sj; j<sj+nj; j++) { 1561e07b27eSBarry Smith for (i=si; i<si+ni; i++) { 1571e07b27eSBarry Smith nidx = (i) + (j)*M; 1581e07b27eSBarry Smith fine_indices[c ] = 2 * nidx ; 1591e07b27eSBarry Smith fine_indices[c+1] = 2 * nidx + 1 ; 1601e07b27eSBarry Smith c = c + 2; 1611e07b27eSBarry Smith } 1621e07b27eSBarry Smith } 1631e07b27eSBarry Smith } else { 1641e07b27eSBarry Smith i = si; 1651e07b27eSBarry Smith j = sj; 1661e07b27eSBarry Smith nidx = (i) + (j)*M; 1671e07b27eSBarry Smith fine_indices[0] = 2 * nidx ; 1681e07b27eSBarry Smith fine_indices[1] = 2 * nidx + 1 ; 1691e07b27eSBarry Smith } 1701e07b27eSBarry Smith 1711e07b27eSBarry Smith /* generate scatter */ 1721e07b27eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*2,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 1731e07b27eSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*2,0,1,&is_local);CHKERRQ(ierr); 1741e07b27eSBarry Smith 1751e07b27eSBarry Smith /* scatter */ 1761e07b27eSBarry Smith ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 1771e07b27eSBarry Smith ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*2);CHKERRQ(ierr); 1781e07b27eSBarry Smith ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 1791e07b27eSBarry Smith 1801e07b27eSBarry Smith ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 1811e07b27eSBarry Smith ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1821e07b27eSBarry Smith ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1831e07b27eSBarry Smith /* access */ 1841e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 1851e07b27eSBarry Smith Vec _coors; 1861e07b27eSBarry Smith const PetscScalar *LA_perm; 1871e07b27eSBarry Smith PetscScalar *LA_coors; 1881e07b27eSBarry Smith 1891e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 1901e07b27eSBarry Smith ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 1911e07b27eSBarry Smith ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 1921e07b27eSBarry Smith for (i=0; i<Ml*Nl*2; i++) { 1931e07b27eSBarry Smith LA_coors[i] = LA_perm[i]; 1941e07b27eSBarry Smith } 1951e07b27eSBarry Smith ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 1961e07b27eSBarry Smith ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 1971e07b27eSBarry Smith } 1981e07b27eSBarry Smith 1991e07b27eSBarry Smith /* update local coords */ 2001e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 2011e07b27eSBarry Smith DM _dmc; 2021e07b27eSBarry Smith Vec _coors,_coors_local; 2031e07b27eSBarry Smith ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 2041e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 2051e07b27eSBarry Smith ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 2061e07b27eSBarry Smith ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 2071e07b27eSBarry Smith ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 2081e07b27eSBarry Smith } 2091e07b27eSBarry Smith ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 2101e07b27eSBarry Smith ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 2111e07b27eSBarry Smith ierr = PetscFree(fine_indices);CHKERRQ(ierr); 2121e07b27eSBarry Smith ierr = ISDestroy(&is_local);CHKERRQ(ierr); 2131e07b27eSBarry Smith ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 2141e07b27eSBarry Smith ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 2151e07b27eSBarry Smith PetscFunctionReturn(0); 2161e07b27eSBarry Smith } 2171e07b27eSBarry Smith 2181e07b27eSBarry Smith #undef __FUNCT__ 2191e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_repart_coors3d" 2201e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PetscSubcomm psubcomm,DM dm,DM subdm) 2211e07b27eSBarry Smith { 2221e07b27eSBarry Smith PetscErrorCode ierr; 2231e07b27eSBarry Smith DM cdm; 2241e07b27eSBarry Smith Vec coor,coor_natural,perm_coors; 2251e07b27eSBarry Smith PetscInt i,j,k,si,sj,sk,ni,nj,nk,M,N,P,Ml,Nl,Pl,c,nidx; 2261e07b27eSBarry Smith PetscInt *fine_indices; 2271e07b27eSBarry Smith IS is_fine,is_local; 2281e07b27eSBarry Smith VecScatter sctx; 2291e07b27eSBarry Smith 2301e07b27eSBarry Smith PetscFunctionBegin; 2311e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 2321e07b27eSBarry Smith if (!coor) PetscFunctionReturn(0); 2331e07b27eSBarry Smith 2341e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 2351e07b27eSBarry Smith ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 2361e07b27eSBarry Smith } 2371e07b27eSBarry Smith 2381e07b27eSBarry Smith /* Get the coordinate vector from the distributed array */ 2391e07b27eSBarry Smith ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 2401e07b27eSBarry Smith ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 2411e07b27eSBarry Smith ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 2421e07b27eSBarry Smith ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 2431e07b27eSBarry Smith 2441e07b27eSBarry Smith /* get indices of the guys I want to grab */ 2451e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2461e07b27eSBarry Smith 2471e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 2481e07b27eSBarry Smith ierr = DMDAGetCorners(subdm,&si,&sj,&sk,&ni,&nj,&nk);CHKERRQ(ierr); 2491e07b27eSBarry Smith 2501e07b27eSBarry Smith Ml = ni; 2511e07b27eSBarry Smith Nl = nj; 2521e07b27eSBarry Smith Pl = nk; 2531e07b27eSBarry Smith } else { 2541e07b27eSBarry Smith Ml = Nl = Pl = 1; 2551e07b27eSBarry Smith } 2561e07b27eSBarry Smith 2571e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*Ml*Nl*Pl*3,&fine_indices);CHKERRQ(ierr); 2581e07b27eSBarry Smith 2591e07b27eSBarry Smith c = 0; 2601e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 2611e07b27eSBarry Smith for (k=sk; k<sk+nk; k++) { 2621e07b27eSBarry Smith for (j=sj; j<sj+nj; j++) { 2631e07b27eSBarry Smith for (i=si; i<si+ni; i++) { 2641e07b27eSBarry Smith nidx = (i) + (j)*M + (k)*M*N; 2651e07b27eSBarry Smith fine_indices[c ] = 3 * nidx ; 2661e07b27eSBarry Smith fine_indices[c+1] = 3 * nidx + 1 ; 2671e07b27eSBarry Smith fine_indices[c+2] = 3 * nidx + 2 ; 2681e07b27eSBarry Smith c = c + 3; 2691e07b27eSBarry Smith } 2701e07b27eSBarry Smith } 2711e07b27eSBarry Smith } 2721e07b27eSBarry Smith } else { 2731e07b27eSBarry Smith i = si; 2741e07b27eSBarry Smith j = sj; 2751e07b27eSBarry Smith k = sk; 2761e07b27eSBarry Smith nidx = (i) + (j)*M + (k)*M*N; 2771e07b27eSBarry Smith fine_indices[0] = 3 * nidx ; 2781e07b27eSBarry Smith fine_indices[1] = 3 * nidx + 1 ; 2791e07b27eSBarry Smith fine_indices[2] = 3 * nidx + 2 ; 2801e07b27eSBarry Smith } 2811e07b27eSBarry Smith 2821e07b27eSBarry Smith /* generate scatter */ 2831e07b27eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*Pl*3,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 2841e07b27eSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*Pl*3,0,1,&is_local);CHKERRQ(ierr); 2851e07b27eSBarry Smith 2861e07b27eSBarry Smith /* scatter */ 2871e07b27eSBarry Smith ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 2881e07b27eSBarry Smith ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*Pl*3);CHKERRQ(ierr); 2891e07b27eSBarry Smith ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 2901e07b27eSBarry Smith ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 2911e07b27eSBarry Smith ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2921e07b27eSBarry Smith ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2931e07b27eSBarry Smith 2941e07b27eSBarry Smith /* access */ 2951e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 2961e07b27eSBarry Smith Vec _coors; 2971e07b27eSBarry Smith const PetscScalar *LA_perm; 2981e07b27eSBarry Smith PetscScalar *LA_coors; 2991e07b27eSBarry Smith 3001e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 3011e07b27eSBarry Smith ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 3021e07b27eSBarry Smith ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 3031e07b27eSBarry Smith for (i=0; i<Ml*Nl*Pl*3; i++) { 3041e07b27eSBarry Smith LA_coors[i] = LA_perm[i]; 3051e07b27eSBarry Smith } 3061e07b27eSBarry Smith ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 3071e07b27eSBarry Smith ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 3081e07b27eSBarry Smith } 3091e07b27eSBarry Smith 3101e07b27eSBarry Smith /* update local coords */ 3111e07b27eSBarry Smith if (isActiveRank(psubcomm)) { 3121e07b27eSBarry Smith DM _dmc; 3131e07b27eSBarry Smith Vec _coors,_coors_local; 3141e07b27eSBarry Smith 3151e07b27eSBarry Smith ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 3161e07b27eSBarry Smith ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 3171e07b27eSBarry Smith ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 3181e07b27eSBarry Smith ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 3191e07b27eSBarry Smith ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 3201e07b27eSBarry Smith } 3211e07b27eSBarry Smith 3221e07b27eSBarry Smith ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 3231e07b27eSBarry Smith ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 3241e07b27eSBarry Smith ierr = PetscFree(fine_indices);CHKERRQ(ierr); 3251e07b27eSBarry Smith ierr = ISDestroy(&is_local);CHKERRQ(ierr); 3261e07b27eSBarry Smith ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 3271e07b27eSBarry Smith ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 3281e07b27eSBarry Smith PetscFunctionReturn(0); 3291e07b27eSBarry Smith } 3301e07b27eSBarry Smith 3311e07b27eSBarry Smith #undef __FUNCT__ 3321e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_repart_coors" 3331e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_repart_coors(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 3341e07b27eSBarry Smith { 3351e07b27eSBarry Smith PetscInt dim; 3361e07b27eSBarry Smith DM dm,subdm; 3371e07b27eSBarry Smith PetscSubcomm psubcomm; 3381e07b27eSBarry Smith MPI_Comm comm; 3391e07b27eSBarry Smith Vec coor; 3401e07b27eSBarry Smith PetscErrorCode ierr; 3411e07b27eSBarry Smith 3421e07b27eSBarry Smith PetscFunctionBegin; 3431e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3441e07b27eSBarry Smith ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 3451e07b27eSBarry Smith if (!coor) PetscFunctionReturn(0); 3461e07b27eSBarry Smith psubcomm = sred->psubcomm; 3471e07b27eSBarry Smith comm = PetscSubcommParent(psubcomm); 3481e07b27eSBarry Smith subdm = ctx->dmrepart; 3491e07b27eSBarry Smith 3501e07b27eSBarry Smith 3511e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the coordinates (DMDA)\n");CHKERRQ(ierr); 3521e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 3531e07b27eSBarry Smith switch (dim) { 3541e07b27eSBarry Smith case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 3551e07b27eSBarry Smith break; 3561e07b27eSBarry Smith case 2: PCTelescopeSetUp_dmda_repart_coors2d(psubcomm,dm,subdm); 3571e07b27eSBarry Smith break; 3581e07b27eSBarry Smith case 3: PCTelescopeSetUp_dmda_repart_coors3d(psubcomm,dm,subdm); 3591e07b27eSBarry Smith break; 3601e07b27eSBarry Smith } 3611e07b27eSBarry Smith PetscFunctionReturn(0); 3621e07b27eSBarry Smith } 3631e07b27eSBarry Smith 3641e07b27eSBarry Smith /* setup repartitioned dm */ 3651e07b27eSBarry Smith #undef __FUNCT__ 3661e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_repart" 3671e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 3681e07b27eSBarry Smith { 3691e07b27eSBarry Smith PetscErrorCode ierr; 3701e07b27eSBarry Smith DM dm; 3711e07b27eSBarry Smith PetscInt dim,nx,ny,nz,ndof,nsw,sum,k; 3721e07b27eSBarry Smith DMBoundaryType bx,by,bz; 3731e07b27eSBarry Smith DMDAStencilType stencil; 3741e07b27eSBarry Smith const PetscInt *_range_i_re; 3751e07b27eSBarry Smith const PetscInt *_range_j_re; 3761e07b27eSBarry Smith const PetscInt *_range_k_re; 3771e07b27eSBarry Smith DMDAInterpolationType itype; 3781e07b27eSBarry Smith PetscInt refine_x,refine_y,refine_z; 3791e07b27eSBarry Smith MPI_Comm comm,subcomm; 3801e07b27eSBarry Smith const char *prefix; 3811e07b27eSBarry Smith 3821e07b27eSBarry Smith PetscFunctionBegin; 3831e07b27eSBarry Smith comm = PetscSubcommParent(sred->psubcomm); 3841e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 3851e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 3861e07b27eSBarry Smith 3871e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,&nsw,&bx,&by,&bz,&stencil);CHKERRQ(ierr); 3881e07b27eSBarry Smith ierr = DMDAGetInterpolationType(dm,&itype);CHKERRQ(ierr); 3891e07b27eSBarry Smith ierr = DMDAGetRefinementFactor(dm,&refine_x,&refine_y,&refine_z);CHKERRQ(ierr); 3901e07b27eSBarry Smith 3911e07b27eSBarry Smith ctx->dmrepart = NULL; 3921e07b27eSBarry Smith _range_i_re = _range_j_re = _range_k_re = NULL; 3931e07b27eSBarry Smith /* Create DMDA on the child communicator */ 3941e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 3951e07b27eSBarry Smith switch (dim) { 3961e07b27eSBarry Smith case 1: 3971e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (1D)\n");CHKERRQ(ierr); 3981e07b27eSBarry Smith /*ierr = DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 3991e07b27eSBarry Smith ny = nz = 1; 4001e07b27eSBarry Smith by = bz = DM_BOUNDARY_NONE; 4011e07b27eSBarry Smith break; 4021e07b27eSBarry Smith case 2: 4031e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (2D)\n");CHKERRQ(ierr); 4041e07b27eSBarry Smith /*ierr = DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 4051e07b27eSBarry Smith nz = 1; 4061e07b27eSBarry Smith bz = DM_BOUNDARY_NONE; 4071e07b27eSBarry Smith break; 4081e07b27eSBarry Smith case 3: 4091e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (3D)\n");CHKERRQ(ierr); 4101e07b27eSBarry 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);*/ 4111e07b27eSBarry Smith break; 4121e07b27eSBarry Smith } 4131e07b27eSBarry Smith /* 4141e07b27eSBarry Smith The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append 4151e07b27eSBarry Smith a unique option prefix for the DM, thus I prefer to expose the contents of these API's here. 4161e07b27eSBarry Smith This allows users to control the partitioning of the subDM. 4171e07b27eSBarry Smith */ 4181e07b27eSBarry Smith ierr = DMDACreate(subcomm,&ctx->dmrepart);CHKERRQ(ierr); 4191e07b27eSBarry Smith /* Set unique option prefix name */ 4201e07b27eSBarry Smith ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 4211e07b27eSBarry Smith ierr = DMSetOptionsPrefix(ctx->dmrepart,prefix);CHKERRQ(ierr); 4221e07b27eSBarry Smith ierr = DMAppendOptionsPrefix(ctx->dmrepart,"repart_");CHKERRQ(ierr); 4231e07b27eSBarry Smith /* standard setup from DMDACreate{1,2,3}d() */ 4241e07b27eSBarry Smith ierr = DMSetDimension(ctx->dmrepart,dim);CHKERRQ(ierr); 4251e07b27eSBarry Smith ierr = DMDASetSizes(ctx->dmrepart,nx,ny,nz);CHKERRQ(ierr); 4261e07b27eSBarry Smith ierr = DMDASetNumProcs(ctx->dmrepart,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 4271e07b27eSBarry Smith ierr = DMDASetBoundaryType(ctx->dmrepart,bx,by,bz);CHKERRQ(ierr); 4281e07b27eSBarry Smith ierr = DMDASetDof(ctx->dmrepart,ndof);CHKERRQ(ierr); 4291e07b27eSBarry Smith ierr = DMDASetStencilType(ctx->dmrepart,stencil);CHKERRQ(ierr); 4301e07b27eSBarry Smith ierr = DMDASetStencilWidth(ctx->dmrepart,nsw);CHKERRQ(ierr); 4311e07b27eSBarry Smith ierr = DMDASetOwnershipRanges(ctx->dmrepart,NULL,NULL,NULL);CHKERRQ(ierr); 4321e07b27eSBarry Smith ierr = DMSetFromOptions(ctx->dmrepart);CHKERRQ(ierr); 4331e07b27eSBarry Smith ierr = DMSetUp(ctx->dmrepart);CHKERRQ(ierr); 4341e07b27eSBarry Smith /* Set refinement factors and interpolation type from the partent */ 4351e07b27eSBarry Smith ierr = DMDASetRefinementFactor(ctx->dmrepart,refine_x,refine_y,refine_z);CHKERRQ(ierr); 4361e07b27eSBarry Smith ierr = DMDASetInterpolationType(ctx->dmrepart,itype);CHKERRQ(ierr); 4371e07b27eSBarry Smith 4381e07b27eSBarry 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); 4391e07b27eSBarry Smith ierr = DMDAGetOwnershipRanges(ctx->dmrepart,&_range_i_re,&_range_j_re,&_range_k_re);CHKERRQ(ierr); 4401e07b27eSBarry Smith } 4411e07b27eSBarry Smith 4421e07b27eSBarry Smith /* generate ranges for repartitioned dm */ 4431e07b27eSBarry Smith /* note - assume rank 0 always participates */ 4441e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Mp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4451e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Np_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4461e07b27eSBarry Smith ierr = MPI_Bcast(&ctx->Pp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 4471e07b27eSBarry Smith 4481e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Mp_re,&ctx->range_i_re);CHKERRQ(ierr); 4491e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Np_re,&ctx->range_j_re);CHKERRQ(ierr); 4501e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Pp_re,&ctx->range_k_re);CHKERRQ(ierr); 4511e07b27eSBarry Smith 4521e07b27eSBarry Smith if (_range_i_re != NULL) {ierr = PetscMemcpy(ctx->range_i_re,_range_i_re,sizeof(PetscInt)*ctx->Mp_re);CHKERRQ(ierr);} 4531e07b27eSBarry Smith if (_range_j_re != NULL) {ierr = PetscMemcpy(ctx->range_j_re,_range_j_re,sizeof(PetscInt)*ctx->Np_re);CHKERRQ(ierr);} 4541e07b27eSBarry Smith if (_range_k_re != NULL) {ierr = PetscMemcpy(ctx->range_k_re,_range_k_re,sizeof(PetscInt)*ctx->Pp_re);CHKERRQ(ierr);} 4551e07b27eSBarry Smith 4561e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4571e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4581e07b27eSBarry Smith ierr = MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 4591e07b27eSBarry Smith 4601e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Mp_re,&ctx->start_i_re);CHKERRQ(ierr); 4611e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Np_re,&ctx->start_j_re);CHKERRQ(ierr); 4621e07b27eSBarry Smith ierr = PetscMalloc(sizeof(PetscInt)*ctx->Pp_re,&ctx->start_k_re);CHKERRQ(ierr); 4631e07b27eSBarry Smith 4641e07b27eSBarry Smith sum = 0; 4651e07b27eSBarry Smith for (k=0; k<ctx->Mp_re; k++) { 4661e07b27eSBarry Smith ctx->start_i_re[k] = sum; 4671e07b27eSBarry Smith sum += ctx->range_i_re[k]; 4681e07b27eSBarry Smith } 4691e07b27eSBarry Smith 4701e07b27eSBarry Smith sum = 0; 4711e07b27eSBarry Smith for (k=0; k<ctx->Np_re; k++) { 4721e07b27eSBarry Smith ctx->start_j_re[k] = sum; 4731e07b27eSBarry Smith sum += ctx->range_j_re[k]; 4741e07b27eSBarry Smith } 4751e07b27eSBarry Smith 4761e07b27eSBarry Smith sum = 0; 4771e07b27eSBarry Smith for (k=0; k<ctx->Pp_re; k++) { 4781e07b27eSBarry Smith ctx->start_k_re[k] = sum; 4791e07b27eSBarry Smith sum += ctx->range_k_re[k]; 4801e07b27eSBarry Smith } 4811e07b27eSBarry Smith 4821e07b27eSBarry Smith /* attach dm to ksp on sub communicator */ 4831e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 4841e07b27eSBarry Smith ierr = KSPSetDM(sred->ksp,ctx->dmrepart);CHKERRQ(ierr); 4851e07b27eSBarry Smith ierr = KSPSetDMActive(sred->ksp,PETSC_FALSE);CHKERRQ(ierr); 4861e07b27eSBarry Smith } 4871e07b27eSBarry Smith PetscFunctionReturn(0); 4881e07b27eSBarry Smith } 4891e07b27eSBarry Smith 4901e07b27eSBarry Smith #undef __FUNCT__ 4911e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_permutation_3d" 4921e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 4931e07b27eSBarry Smith { 4941e07b27eSBarry Smith PetscErrorCode ierr; 4951e07b27eSBarry Smith DM dm; 4961e07b27eSBarry Smith MPI_Comm comm; 4971e07b27eSBarry Smith Mat Pscalar,P; 4981e07b27eSBarry Smith PetscInt ndof; 4991e07b27eSBarry Smith PetscInt i,j,k,location,startI[3],endI[3],lenI[3],nx,ny,nz; 5001e07b27eSBarry Smith PetscInt sr,er,Mr; 5011e07b27eSBarry Smith Vec V; 5021e07b27eSBarry Smith 5031e07b27eSBarry Smith PetscFunctionBegin; 5041e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-3D)\n");CHKERRQ(ierr); 5051e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5061e07b27eSBarry Smith 5071e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 5081e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 5091e07b27eSBarry Smith 5101e07b27eSBarry Smith ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 5111e07b27eSBarry Smith ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 5121e07b27eSBarry Smith ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 5131e07b27eSBarry Smith ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 5141e07b27eSBarry Smith sr = sr / ndof; 5151e07b27eSBarry Smith er = er / ndof; 5161e07b27eSBarry Smith Mr = Mr / ndof; 5171e07b27eSBarry Smith 5181e07b27eSBarry Smith ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 5191e07b27eSBarry Smith ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 5201e07b27eSBarry Smith ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 5211e07b27eSBarry Smith ierr = MatSeqAIJSetPreallocation(Pscalar,2,NULL);CHKERRQ(ierr); 5221e07b27eSBarry Smith ierr = MatMPIAIJSetPreallocation(Pscalar,2,NULL,2,NULL);CHKERRQ(ierr); 5231e07b27eSBarry Smith 5241e07b27eSBarry Smith ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]);CHKERRQ(ierr); 5251e07b27eSBarry Smith ierr = DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]);CHKERRQ(ierr); 5261e07b27eSBarry Smith endI[0] += startI[0]; 5271e07b27eSBarry Smith endI[1] += startI[1]; 5281e07b27eSBarry Smith endI[2] += startI[2]; 5291e07b27eSBarry Smith 5301e07b27eSBarry Smith for (k=startI[2]; k<endI[2]; k++) { 5311e07b27eSBarry Smith for (j=startI[1]; j<endI[1]; j++) { 5321e07b27eSBarry Smith for (i=startI[0]; i<endI[0]; i++) { 5331e07b27eSBarry Smith PetscMPIInt rank_ijk_re,rank_reI[3]; 5341e07b27eSBarry Smith PetscInt s0_re; 5351e07b27eSBarry Smith PetscInt ii,jj,kk,local_ijk_re,mapped_ijk,natural_ijk; 5361e07b27eSBarry Smith PetscInt lenI_re[3]; 5371e07b27eSBarry Smith 5381e07b27eSBarry Smith location = (i - startI[0]) + (j - startI[1])*lenI[0] + (k - startI[2])*lenI[0]*lenI[1]; 5391e07b27eSBarry Smith ierr = _DMDADetermineRankFromGlobalIJK(3,i,j,k, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 5401e07b27eSBarry Smith ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 5411e07b27eSBarry Smith ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 5421e07b27eSBarry Smith &rank_reI[0],&rank_reI[1],&rank_reI[2],&rank_ijk_re);CHKERRQ(ierr); 5431e07b27eSBarry 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); 5441e07b27eSBarry Smith ii = i - ctx->start_i_re[ rank_reI[0] ]; 5451e07b27eSBarry Smith if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error ii"); 5461e07b27eSBarry Smith jj = j - ctx->start_j_re[ rank_reI[1] ]; 5471e07b27eSBarry Smith if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error jj"); 5481e07b27eSBarry Smith kk = k - ctx->start_k_re[ rank_reI[2] ]; 5491e07b27eSBarry Smith if (kk < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error kk"); 5501e07b27eSBarry Smith lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 5511e07b27eSBarry Smith lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 5521e07b27eSBarry Smith lenI_re[2] = ctx->range_k_re[ rank_reI[2] ]; 5531e07b27eSBarry Smith local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1]; 5541e07b27eSBarry Smith mapped_ijk = s0_re + local_ijk_re; 5551e07b27eSBarry Smith natural_ijk = i + j*nx + k*nx*ny; 5561e07b27eSBarry Smith ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 5571e07b27eSBarry Smith } 5581e07b27eSBarry Smith } 5591e07b27eSBarry Smith } 5601e07b27eSBarry Smith ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5611e07b27eSBarry Smith ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5621e07b27eSBarry Smith ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 5631e07b27eSBarry Smith ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 5641e07b27eSBarry Smith ctx->permutation = P; 5651e07b27eSBarry Smith PetscFunctionReturn(0); 5661e07b27eSBarry Smith } 5671e07b27eSBarry Smith 5681e07b27eSBarry Smith #undef __FUNCT__ 5691e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_permutation_2d" 5701e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 5711e07b27eSBarry Smith { 5721e07b27eSBarry Smith PetscErrorCode ierr; 5731e07b27eSBarry Smith DM dm; 5741e07b27eSBarry Smith MPI_Comm comm; 5751e07b27eSBarry Smith Mat Pscalar,P; 5761e07b27eSBarry Smith PetscInt ndof; 5771e07b27eSBarry Smith PetscInt i,j,location,startI[2],endI[2],lenI[2],nx,ny,nz; 5781e07b27eSBarry Smith PetscInt sr,er,Mr; 5791e07b27eSBarry Smith Vec V; 5801e07b27eSBarry Smith 5811e07b27eSBarry Smith PetscFunctionBegin; 5821e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-2D)\n");CHKERRQ(ierr); 5831e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 5841e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 5851e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 5861e07b27eSBarry Smith ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 5871e07b27eSBarry Smith ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 5881e07b27eSBarry Smith ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 5891e07b27eSBarry Smith ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 5901e07b27eSBarry Smith sr = sr / ndof; 5911e07b27eSBarry Smith er = er / ndof; 5921e07b27eSBarry Smith Mr = Mr / ndof; 5931e07b27eSBarry Smith 5941e07b27eSBarry Smith ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 5951e07b27eSBarry Smith ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 5961e07b27eSBarry Smith ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 5971e07b27eSBarry Smith ierr = MatSeqAIJSetPreallocation(Pscalar,2,NULL);CHKERRQ(ierr); 5981e07b27eSBarry Smith ierr = MatMPIAIJSetPreallocation(Pscalar,2,NULL,2,NULL);CHKERRQ(ierr); 5991e07b27eSBarry Smith 6001e07b27eSBarry Smith ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL);CHKERRQ(ierr); 6011e07b27eSBarry Smith ierr = DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL);CHKERRQ(ierr); 6021e07b27eSBarry Smith endI[0] += startI[0]; 6031e07b27eSBarry Smith endI[1] += startI[1]; 6041e07b27eSBarry Smith 6051e07b27eSBarry Smith for (j=startI[1]; j<endI[1]; j++) { 6061e07b27eSBarry Smith for (i=startI[0]; i<endI[0]; i++) { 6071e07b27eSBarry Smith PetscMPIInt rank_ijk_re,rank_reI[3]; 6081e07b27eSBarry Smith PetscInt s0_re; 6091e07b27eSBarry Smith PetscInt ii,jj,local_ijk_re,mapped_ijk,natural_ijk; 6101e07b27eSBarry Smith PetscInt lenI_re[3]; 6111e07b27eSBarry Smith 6121e07b27eSBarry Smith location = (i - startI[0]) + (j - startI[1])*lenI[0]; 6131e07b27eSBarry Smith ierr = _DMDADetermineRankFromGlobalIJK(2,i,j,0, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 6141e07b27eSBarry Smith ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 6151e07b27eSBarry Smith ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 6161e07b27eSBarry Smith &rank_reI[0],&rank_reI[1],NULL,&rank_ijk_re);CHKERRQ(ierr); 6171e07b27eSBarry Smith 6181e07b27eSBarry 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); 6191e07b27eSBarry Smith 6201e07b27eSBarry Smith ii = i - ctx->start_i_re[ rank_reI[0] ]; 6211e07b27eSBarry Smith if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error ii"); 6221e07b27eSBarry Smith jj = j - ctx->start_j_re[ rank_reI[1] ]; 6231e07b27eSBarry Smith if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error jj"); 6241e07b27eSBarry Smith 6251e07b27eSBarry Smith lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 6261e07b27eSBarry Smith lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 6271e07b27eSBarry Smith local_ijk_re = ii + jj * lenI_re[0]; 6281e07b27eSBarry Smith mapped_ijk = s0_re + local_ijk_re; 6291e07b27eSBarry Smith natural_ijk = i + j*nx; 6301e07b27eSBarry Smith ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 6311e07b27eSBarry Smith } 6321e07b27eSBarry Smith } 6331e07b27eSBarry Smith ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6341e07b27eSBarry Smith ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6351e07b27eSBarry Smith ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 6361e07b27eSBarry Smith ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 6371e07b27eSBarry Smith ctx->permutation = P; 6381e07b27eSBarry Smith PetscFunctionReturn(0); 6391e07b27eSBarry Smith } 6401e07b27eSBarry Smith 6411e07b27eSBarry Smith #undef __FUNCT__ 6421e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda_scatters" 6431e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 6441e07b27eSBarry Smith { 6451e07b27eSBarry Smith PetscErrorCode ierr; 6461e07b27eSBarry Smith Vec xred,yred,xtmp,x,xp; 6471e07b27eSBarry Smith VecScatter scatter; 6481e07b27eSBarry Smith IS isin; 6491e07b27eSBarry Smith Mat B; 6501e07b27eSBarry Smith PetscInt m,bs,st,ed; 6511e07b27eSBarry Smith MPI_Comm comm; 6521e07b27eSBarry Smith 6531e07b27eSBarry Smith PetscFunctionBegin; 6541e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 6551e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 6561e07b27eSBarry Smith ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 6571e07b27eSBarry Smith ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 6581e07b27eSBarry Smith ierr = VecDuplicate(x,&xp);CHKERRQ(ierr); 6591e07b27eSBarry Smith m = bs; 6601e07b27eSBarry Smith xred = NULL; 6611e07b27eSBarry Smith yred = NULL; 6621e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 6631e07b27eSBarry Smith ierr = DMCreateGlobalVector(ctx->dmrepart,&xred);CHKERRQ(ierr); 6641e07b27eSBarry Smith ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 6651e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 6661e07b27eSBarry Smith ierr = ISCreateStride(comm,ed-st,st,1,&isin);CHKERRQ(ierr); 6671e07b27eSBarry Smith ierr = VecGetLocalSize(xred,&m); 6681e07b27eSBarry Smith } else { 6691e07b27eSBarry Smith /* fetch some local owned data - just to deal with avoiding zero length ownership on range */ 6701e07b27eSBarry Smith ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 6711e07b27eSBarry Smith ierr = ISCreateStride(comm,bs,st,1,&isin);CHKERRQ(ierr); 6721e07b27eSBarry Smith } 6731e07b27eSBarry Smith ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 6741e07b27eSBarry Smith ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 6751e07b27eSBarry Smith ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 6761e07b27eSBarry Smith ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 6771e07b27eSBarry Smith ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 6781e07b27eSBarry Smith ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 6791e07b27eSBarry Smith sred->xred = xred; 6801e07b27eSBarry Smith sred->yred = yred; 6811e07b27eSBarry Smith sred->isin = isin; 6821e07b27eSBarry Smith sred->scatter = scatter; 6831e07b27eSBarry Smith sred->xtmp = xtmp; 6841e07b27eSBarry Smith 6851e07b27eSBarry Smith ctx->xp = xp; 6861e07b27eSBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 6871e07b27eSBarry Smith PetscFunctionReturn(0); 6881e07b27eSBarry Smith } 6891e07b27eSBarry Smith 6901e07b27eSBarry Smith #undef __FUNCT__ 6911e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_dmda" 6921e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_dmda(PC pc,PC_Telescope sred) 6931e07b27eSBarry Smith { 6941e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 6951e07b27eSBarry Smith PetscInt dim; 6961e07b27eSBarry Smith DM dm; 6971e07b27eSBarry Smith MPI_Comm comm; 6981e07b27eSBarry Smith PetscErrorCode ierr; 6991e07b27eSBarry Smith 7001e07b27eSBarry Smith PetscFunctionBegin; 7011e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: setup (DMDA)\n");CHKERRQ(ierr); 7021e07b27eSBarry Smith PetscMalloc(sizeof(PC_Telescope_DMDACtx),&ctx); 7031e07b27eSBarry Smith PetscMemzero(ctx,sizeof(PC_Telescope_DMDACtx)); 7041e07b27eSBarry Smith sred->dm_ctx = (void*)ctx; 7051e07b27eSBarry Smith 7061e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7071e07b27eSBarry Smith ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 7081e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 7091e07b27eSBarry Smith 7101e07b27eSBarry Smith PCTelescopeSetUp_dmda_repart(pc,sred,ctx); 7111e07b27eSBarry Smith PCTelescopeSetUp_dmda_repart_coors(pc,sred,ctx); 7121e07b27eSBarry Smith switch (dim) { 7131e07b27eSBarry Smith case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 7141e07b27eSBarry Smith break; 7151e07b27eSBarry Smith case 2: ierr = PCTelescopeSetUp_dmda_permutation_2d(pc,sred,ctx);CHKERRQ(ierr); 7161e07b27eSBarry Smith break; 7171e07b27eSBarry Smith case 3: ierr = PCTelescopeSetUp_dmda_permutation_3d(pc,sred,ctx);CHKERRQ(ierr); 7181e07b27eSBarry Smith break; 7191e07b27eSBarry Smith } 7201e07b27eSBarry Smith ierr = PCTelescopeSetUp_dmda_scatters(pc,sred,ctx);CHKERRQ(ierr); 7211e07b27eSBarry Smith PetscFunctionReturn(0); 7221e07b27eSBarry Smith } 7231e07b27eSBarry Smith 7241e07b27eSBarry Smith #undef __FUNCT__ 7251e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatCreate_dmda" 7261e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_dmda(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 7271e07b27eSBarry Smith { 7281e07b27eSBarry Smith PetscErrorCode ierr; 7291e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 7301e07b27eSBarry Smith MPI_Comm comm,subcomm; 7311e07b27eSBarry Smith Mat Bperm,Bred,B,P; 7321e07b27eSBarry Smith PetscInt nr,nc; 7331e07b27eSBarry Smith IS isrow,iscol; 7341e07b27eSBarry Smith Mat Blocal,*_Blocal; 7351e07b27eSBarry Smith 7361e07b27eSBarry Smith PetscFunctionBegin; 7371e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (DMDA)\n");CHKERRQ(ierr); 7381e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7391e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 7401e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 7411e07b27eSBarry Smith 7421e07b27eSBarry Smith ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 7431e07b27eSBarry Smith ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 7441e07b27eSBarry Smith 7451e07b27eSBarry Smith P = ctx->permutation; 7461e07b27eSBarry Smith ierr = MatPtAP(B,P,MAT_INITIAL_MATRIX,1.1,&Bperm);CHKERRQ(ierr); 7471e07b27eSBarry Smith 7481e07b27eSBarry Smith /* Get submatrices */ 7491e07b27eSBarry Smith isrow = sred->isin; 7501e07b27eSBarry Smith ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 7511e07b27eSBarry Smith 7521e07b27eSBarry Smith ierr = MatGetSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 7531e07b27eSBarry Smith Blocal = *_Blocal; 7541e07b27eSBarry Smith Bred = NULL; 7551e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 7561e07b27eSBarry Smith PetscInt mm; 7571e07b27eSBarry Smith 7581e07b27eSBarry Smith if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;} 7591e07b27eSBarry Smith ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 760*bfd6bcc6SSatish Balay /* ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred);CHKERRQ(ierr); */ 7611e07b27eSBarry Smith ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 7621e07b27eSBarry Smith } 7631e07b27eSBarry Smith *A = Bred; 7641e07b27eSBarry Smith 7651e07b27eSBarry Smith ierr = ISDestroy(&iscol);CHKERRQ(ierr); 7661e07b27eSBarry Smith ierr = MatDestroy(&Bperm);CHKERRQ(ierr); 7671e07b27eSBarry Smith ierr = MatDestroyMatrices(1,&_Blocal);CHKERRQ(ierr); 7681e07b27eSBarry Smith PetscFunctionReturn(0); 7691e07b27eSBarry Smith } 7701e07b27eSBarry Smith 7711e07b27eSBarry Smith #undef __FUNCT__ 7721e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_dmda" 7731e07b27eSBarry Smith PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat) 7741e07b27eSBarry Smith { 7751e07b27eSBarry Smith PetscErrorCode ierr; 7761e07b27eSBarry Smith MatNullSpace nullspace,sub_nullspace; 7771e07b27eSBarry Smith Mat A,B; 7781e07b27eSBarry Smith PetscBool has_const; 7791e07b27eSBarry Smith PetscInt i,k,n; 7801e07b27eSBarry Smith const Vec *vecs; 7811e07b27eSBarry Smith Vec *sub_vecs; 7821e07b27eSBarry Smith MPI_Comm subcomm; 7831e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 7841e07b27eSBarry Smith 7851e07b27eSBarry Smith PetscFunctionBegin; 7861e07b27eSBarry Smith ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr); 7871e07b27eSBarry Smith ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 7881e07b27eSBarry Smith if (!nullspace) return(0); 7891e07b27eSBarry Smith 7901e07b27eSBarry Smith ierr = PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n");CHKERRQ(ierr); 7911e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 7921e07b27eSBarry Smith subcomm = PetscSubcommChild(sred->psubcomm); 7931e07b27eSBarry Smith ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 7941e07b27eSBarry Smith 7951e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 7961e07b27eSBarry Smith sub_vecs = NULL; 7971e07b27eSBarry Smith /* create new vectors */ 7981e07b27eSBarry Smith if (n != 0) { 7991e07b27eSBarry Smith PetscMalloc(sizeof(Vec)*n,&sub_vecs); 8001e07b27eSBarry Smith for (k=0; k<n; k++) { 8011e07b27eSBarry Smith ierr = VecDuplicate(sred->xred,&sub_vecs[k]);CHKERRQ(ierr); 8021e07b27eSBarry Smith } 8031e07b27eSBarry Smith } 8041e07b27eSBarry Smith } 8051e07b27eSBarry Smith 8061e07b27eSBarry Smith /* copy entries */ 8071e07b27eSBarry Smith for (k=0; k<n; k++) { 8081e07b27eSBarry Smith const PetscScalar *x_array; 8091e07b27eSBarry Smith PetscScalar *LA_sub_vec; 8101e07b27eSBarry Smith PetscInt st,ed,bs; 8111e07b27eSBarry Smith 8121e07b27eSBarry Smith /* permute vector into ordering associated with re-partitioned dmda */ 8131e07b27eSBarry Smith ierr = MatMultTranspose(ctx->permutation,vecs[k],ctx->xp);CHKERRQ(ierr); 8141e07b27eSBarry Smith 8151e07b27eSBarry Smith /* pull in vector x->xtmp */ 8161e07b27eSBarry Smith ierr = VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8171e07b27eSBarry Smith ierr = VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8181e07b27eSBarry Smith 8191e07b27eSBarry Smith /* copy vector entires into xred */ 8201e07b27eSBarry Smith ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr); 8211e07b27eSBarry Smith ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 8221e07b27eSBarry Smith if (sub_vecs[k]) { 8231e07b27eSBarry Smith ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 8241e07b27eSBarry Smith ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 8251e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 8261e07b27eSBarry Smith LA_sub_vec[i] = x_array[i]; 8271e07b27eSBarry Smith } 8281e07b27eSBarry Smith ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 8291e07b27eSBarry Smith } 8301e07b27eSBarry Smith ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 8311e07b27eSBarry Smith } 8321e07b27eSBarry Smith 8331e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 8341e07b27eSBarry Smith /* create new nullspace for redundant object */ 8351e07b27eSBarry Smith ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr); 8361e07b27eSBarry Smith 8371e07b27eSBarry Smith /* attach redundant nullspace to Bred */ 8381e07b27eSBarry Smith ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 8391e07b27eSBarry Smith 8401e07b27eSBarry Smith for (k=0; k<n; k++) { 8411e07b27eSBarry Smith ierr = VecDestroy(&sub_vecs[k]);CHKERRQ(ierr); 8421e07b27eSBarry Smith } 8431e07b27eSBarry Smith ierr = PetscFree(sub_vecs);CHKERRQ(ierr); 8441e07b27eSBarry Smith } 8451e07b27eSBarry Smith PetscFunctionReturn(0); 8461e07b27eSBarry Smith } 8471e07b27eSBarry Smith 8481e07b27eSBarry Smith #undef __FUNCT__ 8491e07b27eSBarry Smith #define __FUNCT__ "PCApply_Telescope_dmda" 8501e07b27eSBarry Smith PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y) 8511e07b27eSBarry Smith { 8521e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 8531e07b27eSBarry Smith PetscErrorCode ierr; 8541e07b27eSBarry Smith Mat perm; 8551e07b27eSBarry Smith Vec xtmp,xp,xred,yred; 8561e07b27eSBarry Smith PetscInt i,st,ed,bs; 8571e07b27eSBarry Smith VecScatter scatter; 8581e07b27eSBarry Smith PetscScalar *array; 8591e07b27eSBarry Smith const PetscScalar *x_array; 8601e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 8611e07b27eSBarry Smith 8621e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 8631e07b27eSBarry Smith xtmp = sred->xtmp; 8641e07b27eSBarry Smith scatter = sred->scatter; 8651e07b27eSBarry Smith xred = sred->xred; 8661e07b27eSBarry Smith yred = sred->yred; 8671e07b27eSBarry Smith perm = ctx->permutation; 8681e07b27eSBarry Smith xp = ctx->xp; 8691e07b27eSBarry Smith 8701e07b27eSBarry Smith PetscFunctionBegin; 8711e07b27eSBarry Smith /* permute vector into ordering associated with re-partitioned dmda */ 8721e07b27eSBarry Smith ierr = MatMultTranspose(perm,x,xp);CHKERRQ(ierr); 8731e07b27eSBarry Smith 8741e07b27eSBarry Smith /* pull in vector x->xtmp */ 8751e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8761e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8771e07b27eSBarry Smith 8781e07b27eSBarry Smith /* copy vector entires into xred */ 8791e07b27eSBarry Smith ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 8801e07b27eSBarry Smith ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 8811e07b27eSBarry Smith if (xred) { 8821e07b27eSBarry Smith PetscScalar *LA_xred; 8831e07b27eSBarry Smith ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 8841e07b27eSBarry Smith 8851e07b27eSBarry Smith ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 8861e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 8871e07b27eSBarry Smith LA_xred[i] = x_array[i]; 8881e07b27eSBarry Smith } 8891e07b27eSBarry Smith ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 8901e07b27eSBarry Smith } 8911e07b27eSBarry Smith ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 8921e07b27eSBarry Smith 8931e07b27eSBarry Smith /* solve */ 8941e07b27eSBarry Smith if (isActiveRank(sred->psubcomm)) { 8951e07b27eSBarry Smith ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 8961e07b27eSBarry Smith } 8971e07b27eSBarry Smith 8981e07b27eSBarry Smith /* return vector */ 8991e07b27eSBarry Smith ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 9001e07b27eSBarry Smith ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 9011e07b27eSBarry Smith if (yred) { 9021e07b27eSBarry Smith const PetscScalar *LA_yred; 9031e07b27eSBarry Smith ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 9041e07b27eSBarry Smith ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 9051e07b27eSBarry Smith for (i=0; i<ed-st; i++) { 9061e07b27eSBarry Smith array[i] = LA_yred[i]; 9071e07b27eSBarry Smith } 9081e07b27eSBarry Smith ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 9091e07b27eSBarry Smith } else { 9101e07b27eSBarry Smith for (i=0; i<bs; i++) { 9111e07b27eSBarry Smith array[i] = 0.0; 9121e07b27eSBarry Smith } 9131e07b27eSBarry Smith } 9141e07b27eSBarry Smith ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 9151e07b27eSBarry Smith ierr = VecScatterBegin(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9161e07b27eSBarry Smith ierr = VecScatterEnd(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9171e07b27eSBarry Smith ierr = MatMult(perm,xp,y);CHKERRQ(ierr); 9181e07b27eSBarry Smith PetscFunctionReturn(0); 9191e07b27eSBarry Smith } 9201e07b27eSBarry Smith 9211e07b27eSBarry Smith #undef __FUNCT__ 9221e07b27eSBarry Smith #define __FUNCT__ "PCReset_Telescope_dmda" 9231e07b27eSBarry Smith PetscErrorCode PCReset_Telescope_dmda(PC pc) 9241e07b27eSBarry Smith { 9251e07b27eSBarry Smith PetscErrorCode ierr; 9261e07b27eSBarry Smith PC_Telescope sred = (PC_Telescope)pc->data; 9271e07b27eSBarry Smith PC_Telescope_DMDACtx *ctx; 9281e07b27eSBarry Smith 9291e07b27eSBarry Smith PetscFunctionBegin; 9301e07b27eSBarry Smith ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 9311e07b27eSBarry Smith ierr = VecDestroy(&ctx->xp);CHKERRQ(ierr); 9321e07b27eSBarry Smith ierr = MatDestroy(&ctx->permutation);CHKERRQ(ierr); 9331e07b27eSBarry Smith ierr = DMDestroy(&ctx->dmrepart);CHKERRQ(ierr); 9341e07b27eSBarry Smith ierr = PetscFree(ctx->range_i_re);CHKERRQ(ierr); 9351e07b27eSBarry Smith ierr = PetscFree(ctx->range_j_re);CHKERRQ(ierr); 9361e07b27eSBarry Smith ierr = PetscFree(ctx->range_k_re);CHKERRQ(ierr); 9371e07b27eSBarry Smith ierr = PetscFree(ctx->start_i_re);CHKERRQ(ierr); 9381e07b27eSBarry Smith ierr = PetscFree(ctx->start_j_re);CHKERRQ(ierr); 9391e07b27eSBarry Smith ierr = PetscFree(ctx->start_k_re);CHKERRQ(ierr); 9401e07b27eSBarry Smith PetscFunctionReturn(0); 9411e07b27eSBarry Smith } 9421e07b27eSBarry Smith 9431e07b27eSBarry Smith #undef __FUNCT__ 9441e07b27eSBarry Smith #define __FUNCT__ "DMView_DMDAShort_3d" 9451e07b27eSBarry Smith PetscErrorCode DMView_DMDAShort_3d(DM dm,PetscViewer v) 9461e07b27eSBarry Smith { 9471e07b27eSBarry Smith PetscInt M,N,P,m,n,p,ndof,nsw; 9481e07b27eSBarry Smith MPI_Comm comm; 9491e07b27eSBarry Smith PetscMPIInt size; 9501e07b27eSBarry Smith const char* prefix; 9511e07b27eSBarry Smith PetscErrorCode ierr; 9521e07b27eSBarry Smith 9531e07b27eSBarry Smith PetscFunctionBegin; 9541e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 9551e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 9561e07b27eSBarry Smith ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 9571e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 9581e07b27eSBarry Smith if (prefix) {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 9591e07b27eSBarry Smith else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 9601e07b27eSBarry 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); 9611e07b27eSBarry Smith PetscFunctionReturn(0); 9621e07b27eSBarry Smith } 9631e07b27eSBarry Smith 9641e07b27eSBarry Smith #undef __FUNCT__ 9651e07b27eSBarry Smith #define __FUNCT__ "DMView_DMDAShort_2d" 9661e07b27eSBarry Smith PetscErrorCode DMView_DMDAShort_2d(DM dm,PetscViewer v) 9671e07b27eSBarry Smith { 9681e07b27eSBarry Smith PetscInt M,N,m,n,ndof,nsw; 9691e07b27eSBarry Smith MPI_Comm comm; 9701e07b27eSBarry Smith PetscMPIInt size; 9711e07b27eSBarry Smith const char* prefix; 9721e07b27eSBarry Smith PetscErrorCode ierr; 9731e07b27eSBarry Smith 9741e07b27eSBarry Smith PetscFunctionBegin; 9751e07b27eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 9761e07b27eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 9771e07b27eSBarry Smith ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 9781e07b27eSBarry Smith ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 9791e07b27eSBarry Smith if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 9801e07b27eSBarry Smith else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 9811e07b27eSBarry 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); 9821e07b27eSBarry Smith PetscFunctionReturn(0); 9831e07b27eSBarry Smith } 9841e07b27eSBarry Smith 9851e07b27eSBarry Smith #undef __FUNCT__ 9861e07b27eSBarry Smith #define __FUNCT__ "DMView_DMDAShort" 9871e07b27eSBarry Smith PetscErrorCode DMView_DMDAShort(DM dm,PetscViewer v) 9881e07b27eSBarry Smith { 9891e07b27eSBarry Smith PetscErrorCode ierr; 9901e07b27eSBarry Smith PetscInt dim; 9911e07b27eSBarry Smith 9921e07b27eSBarry Smith PetscFunctionBegin; 9931e07b27eSBarry Smith ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 9941e07b27eSBarry Smith switch (dim) { 9951e07b27eSBarry Smith case 2: ierr = DMView_DMDAShort_2d(dm,v);CHKERRQ(ierr); 9961e07b27eSBarry Smith break; 9971e07b27eSBarry Smith case 3: ierr = DMView_DMDAShort_3d(dm,v);CHKERRQ(ierr); 9981e07b27eSBarry Smith break; 9991e07b27eSBarry Smith } 10001e07b27eSBarry Smith PetscFunctionReturn(0); 10011e07b27eSBarry Smith } 10021e07b27eSBarry Smith 1003