1 2 3 #include <petsc/private/matimpl.h> 4 #include <petsc/private/pcimpl.h> 5 #include <petsc/private/dmimpl.h> 6 #include <petscksp.h> /*I "petscksp.h" I*/ 7 #include <petscdm.h> 8 #include <petscdmda.h> 9 10 #include "../src/ksp/pc/impls/telescope/telescope.h" 11 12 static PetscBool cited = PETSC_FALSE; 13 static const char citation[] = 14 "@inproceedings{MaySananRuppKnepleySmith2016,\n" 15 " title = {Extreme-Scale Multigrid Components within PETSc},\n" 16 " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 17 " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 18 " series = {PASC '16},\n" 19 " isbn = {978-1-4503-4126-4},\n" 20 " location = {Lausanne, Switzerland},\n" 21 " pages = {5:1--5:12},\n" 22 " articleno = {5},\n" 23 " numpages = {12},\n" 24 " url = {https://doi.acm.org/10.1145/2929908.2929913},\n" 25 " doi = {10.1145/2929908.2929913},\n" 26 " acmid = {2929913},\n" 27 " publisher = {ACM},\n" 28 " address = {New York, NY, USA},\n" 29 " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 30 " year = {2016}\n" 31 "}\n"; 32 33 static PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim,PetscInt i,PetscInt j,PetscInt k, 34 PetscInt Mp,PetscInt Np,PetscInt Pp, 35 PetscInt start_i[],PetscInt start_j[],PetscInt start_k[], 36 PetscInt span_i[],PetscInt span_j[],PetscInt span_k[], 37 PetscMPIInt *_pi,PetscMPIInt *_pj,PetscMPIInt *_pk,PetscMPIInt *rank_re) 38 { 39 PetscInt pi,pj,pk,n; 40 41 PetscFunctionBegin; 42 *rank_re = -1; 43 if (_pi) *_pi = -1; 44 if (_pj) *_pj = -1; 45 if (_pk) *_pk = -1; 46 pi = pj = pk = -1; 47 if (_pi) { 48 for (n=0; n<Mp; n++) { 49 if ( (i >= start_i[n]) && (i < start_i[n]+span_i[n]) ) { 50 pi = n; 51 break; 52 } 53 } 54 if (pi == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pi cannot be determined : range %D, val %D",Mp,i); 55 *_pi = pi; 56 } 57 58 if (_pj) { 59 for (n=0; n<Np; n++) { 60 if ( (j >= start_j[n]) && (j < start_j[n]+span_j[n]) ) { 61 pj = n; 62 break; 63 } 64 } 65 if (pj == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pj cannot be determined : range %D, val %D",Np,j); 66 *_pj = pj; 67 } 68 69 if (_pk) { 70 for (n=0; n<Pp; n++) { 71 if ( (k >= start_k[n]) && (k < start_k[n]+span_k[n]) ) { 72 pk = n; 73 break; 74 } 75 } 76 if (pk == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pk cannot be determined : range %D, val %D",Pp,k); 77 *_pk = pk; 78 } 79 80 switch (dim) { 81 case 1: 82 *rank_re = pi; 83 break; 84 case 2: 85 *rank_re = pi + pj * Mp; 86 break; 87 case 3: 88 *rank_re = pi + pj * Mp + pk * (Mp*Np); 89 break; 90 } 91 PetscFunctionReturn(0); 92 } 93 94 static PetscErrorCode _DMDADetermineGlobalS0(PetscInt dim,PetscMPIInt rank_re,PetscInt Mp_re,PetscInt Np_re,PetscInt Pp_re, 95 PetscInt range_i_re[],PetscInt range_j_re[],PetscInt range_k_re[],PetscInt *s0) 96 { 97 PetscInt i,j,k,start_IJK = 0; 98 PetscInt rank_ijk; 99 100 PetscFunctionBegin; 101 switch (dim) { 102 case 1: 103 for (i=0; i<Mp_re; i++) { 104 rank_ijk = i; 105 if (rank_ijk < rank_re) { 106 start_IJK += range_i_re[i]; 107 } 108 } 109 break; 110 case 2: 111 for (j=0; j<Np_re; j++) { 112 for (i=0; i<Mp_re; i++) { 113 rank_ijk = i + j*Mp_re; 114 if (rank_ijk < rank_re) { 115 start_IJK += range_i_re[i]*range_j_re[j]; 116 } 117 } 118 } 119 break; 120 case 3: 121 for (k=0; k<Pp_re; k++) { 122 for (j=0; j<Np_re; j++) { 123 for (i=0; i<Mp_re; i++) { 124 rank_ijk = i + j*Mp_re + k*Mp_re*Np_re; 125 if (rank_ijk < rank_re) { 126 start_IJK += range_i_re[i]*range_j_re[j]*range_k_re[k]; 127 } 128 } 129 } 130 } 131 break; 132 } 133 *s0 = start_IJK; 134 PetscFunctionReturn(0); 135 } 136 137 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PC_Telescope sred,DM dm,DM subdm) 138 { 139 PetscErrorCode ierr; 140 DM cdm; 141 Vec coor,coor_natural,perm_coors; 142 PetscInt i,j,si,sj,ni,nj,M,N,Ml,Nl,c,nidx; 143 PetscInt *fine_indices; 144 IS is_fine,is_local; 145 VecScatter sctx; 146 147 PetscFunctionBegin; 148 ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 149 if (!coor) return(0); 150 if (PCTelescope_isActiveRank(sred)) { 151 ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 152 } 153 /* Get the coordinate vector from the distributed array */ 154 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 155 ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 156 157 ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 158 ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 159 160 /* get indices of the guys I want to grab */ 161 ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 162 if (PCTelescope_isActiveRank(sred)) { 163 ierr = DMDAGetCorners(subdm,&si,&sj,NULL,&ni,&nj,NULL);CHKERRQ(ierr); 164 Ml = ni; 165 Nl = nj; 166 } else { 167 si = sj = 0; 168 ni = nj = 0; 169 Ml = Nl = 0; 170 } 171 172 ierr = PetscMalloc1(Ml*Nl*2,&fine_indices);CHKERRQ(ierr); 173 c = 0; 174 if (PCTelescope_isActiveRank(sred)) { 175 for (j=sj; j<sj+nj; j++) { 176 for (i=si; i<si+ni; i++) { 177 nidx = (i) + (j)*M; 178 fine_indices[c ] = 2 * nidx ; 179 fine_indices[c+1] = 2 * nidx + 1 ; 180 c = c + 2; 181 } 182 } 183 if (c != Ml*Nl*2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c %D should equal 2 * Ml %D * Nl %D",c,Ml,Nl); 184 } 185 186 /* generate scatter */ 187 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*2,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 188 ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*2,0,1,&is_local);CHKERRQ(ierr); 189 190 /* scatter */ 191 ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 192 ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*2);CHKERRQ(ierr); 193 ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 194 195 ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 196 ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 197 ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 198 /* access */ 199 if (PCTelescope_isActiveRank(sred)) { 200 Vec _coors; 201 const PetscScalar *LA_perm; 202 PetscScalar *LA_coors; 203 204 ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 205 ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 206 ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 207 for (i=0; i<Ml*Nl*2; i++) { 208 LA_coors[i] = LA_perm[i]; 209 } 210 ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 211 ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 212 } 213 214 /* update local coords */ 215 if (PCTelescope_isActiveRank(sred)) { 216 DM _dmc; 217 Vec _coors,_coors_local; 218 ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 219 ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 220 ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 221 ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 222 ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 223 } 224 ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 225 ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 226 ierr = PetscFree(fine_indices);CHKERRQ(ierr); 227 ierr = ISDestroy(&is_local);CHKERRQ(ierr); 228 ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 229 ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 233 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PC_Telescope sred,DM dm,DM subdm) 234 { 235 PetscErrorCode ierr; 236 DM cdm; 237 Vec coor,coor_natural,perm_coors; 238 PetscInt i,j,k,si,sj,sk,ni,nj,nk,M,N,P,Ml,Nl,Pl,c,nidx; 239 PetscInt *fine_indices; 240 IS is_fine,is_local; 241 VecScatter sctx; 242 243 PetscFunctionBegin; 244 ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 245 if (!coor) PetscFunctionReturn(0); 246 247 if (PCTelescope_isActiveRank(sred)) { 248 ierr = DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 249 } 250 251 /* Get the coordinate vector from the distributed array */ 252 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 253 ierr = DMDACreateNaturalVector(cdm,&coor_natural);CHKERRQ(ierr); 254 ierr = DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 255 ierr = DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);CHKERRQ(ierr); 256 257 /* get indices of the guys I want to grab */ 258 ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 259 260 if (PCTelescope_isActiveRank(sred)) { 261 ierr = DMDAGetCorners(subdm,&si,&sj,&sk,&ni,&nj,&nk);CHKERRQ(ierr); 262 Ml = ni; 263 Nl = nj; 264 Pl = nk; 265 } else { 266 si = sj = sk = 0; 267 ni = nj = nk = 0; 268 Ml = Nl = Pl = 0; 269 } 270 271 ierr = PetscMalloc1(Ml*Nl*Pl*3,&fine_indices);CHKERRQ(ierr); 272 273 c = 0; 274 if (PCTelescope_isActiveRank(sred)) { 275 for (k=sk; k<sk+nk; k++) { 276 for (j=sj; j<sj+nj; j++) { 277 for (i=si; i<si+ni; i++) { 278 nidx = (i) + (j)*M + (k)*M*N; 279 fine_indices[c ] = 3 * nidx ; 280 fine_indices[c+1] = 3 * nidx + 1 ; 281 fine_indices[c+2] = 3 * nidx + 2 ; 282 c = c + 3; 283 } 284 } 285 } 286 } 287 288 /* generate scatter */ 289 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*Pl*3,fine_indices,PETSC_USE_POINTER,&is_fine);CHKERRQ(ierr); 290 ierr = ISCreateStride(PETSC_COMM_SELF,Ml*Nl*Pl*3,0,1,&is_local);CHKERRQ(ierr); 291 292 /* scatter */ 293 ierr = VecCreate(PETSC_COMM_SELF,&perm_coors);CHKERRQ(ierr); 294 ierr = VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*Pl*3);CHKERRQ(ierr); 295 ierr = VecSetType(perm_coors,VECSEQ);CHKERRQ(ierr); 296 ierr = VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);CHKERRQ(ierr); 297 ierr = VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 298 ierr = VecScatterEnd( sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 299 300 /* access */ 301 if (PCTelescope_isActiveRank(sred)) { 302 Vec _coors; 303 const PetscScalar *LA_perm; 304 PetscScalar *LA_coors; 305 306 ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 307 ierr = VecGetArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 308 ierr = VecGetArray(_coors,&LA_coors);CHKERRQ(ierr); 309 for (i=0; i<Ml*Nl*Pl*3; i++) { 310 LA_coors[i] = LA_perm[i]; 311 } 312 ierr = VecRestoreArray(_coors,&LA_coors);CHKERRQ(ierr); 313 ierr = VecRestoreArrayRead(perm_coors,&LA_perm);CHKERRQ(ierr); 314 } 315 316 /* update local coords */ 317 if (PCTelescope_isActiveRank(sred)) { 318 DM _dmc; 319 Vec _coors,_coors_local; 320 321 ierr = DMGetCoordinateDM(subdm,&_dmc);CHKERRQ(ierr); 322 ierr = DMGetCoordinates(subdm,&_coors);CHKERRQ(ierr); 323 ierr = DMGetCoordinatesLocal(subdm,&_coors_local);CHKERRQ(ierr); 324 ierr = DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 325 ierr = DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);CHKERRQ(ierr); 326 } 327 328 ierr = VecScatterDestroy(&sctx);CHKERRQ(ierr); 329 ierr = ISDestroy(&is_fine);CHKERRQ(ierr); 330 ierr = PetscFree(fine_indices);CHKERRQ(ierr); 331 ierr = ISDestroy(&is_local);CHKERRQ(ierr); 332 ierr = VecDestroy(&perm_coors);CHKERRQ(ierr); 333 ierr = VecDestroy(&coor_natural);CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 338 { 339 PetscInt dim; 340 DM dm,subdm; 341 PetscSubcomm psubcomm; 342 MPI_Comm comm; 343 Vec coor; 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 348 ierr = DMGetCoordinates(dm,&coor);CHKERRQ(ierr); 349 if (!coor) PetscFunctionReturn(0); 350 psubcomm = sred->psubcomm; 351 comm = PetscSubcommParent(psubcomm); 352 subdm = ctx->dmrepart; 353 354 ierr = PetscInfo(pc,"PCTelescope: setting up the coordinates (DMDA)\n");CHKERRQ(ierr); 355 ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 356 switch (dim) { 357 case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 358 break; 359 case 2: ierr = PCTelescopeSetUp_dmda_repart_coors2d(sred,dm,subdm);CHKERRQ(ierr); 360 break; 361 case 3: ierr = PCTelescopeSetUp_dmda_repart_coors3d(sred,dm,subdm);CHKERRQ(ierr); 362 break; 363 } 364 PetscFunctionReturn(0); 365 } 366 367 /* setup repartitioned dm */ 368 PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 369 { 370 PetscErrorCode ierr; 371 DM dm; 372 PetscInt dim,nx,ny,nz,ndof,nsw,sum,k; 373 DMBoundaryType bx,by,bz; 374 DMDAStencilType stencil; 375 const PetscInt *_range_i_re; 376 const PetscInt *_range_j_re; 377 const PetscInt *_range_k_re; 378 DMDAInterpolationType itype; 379 PetscInt refine_x,refine_y,refine_z; 380 MPI_Comm comm,subcomm; 381 const char *prefix; 382 383 PetscFunctionBegin; 384 comm = PetscSubcommParent(sred->psubcomm); 385 subcomm = PetscSubcommChild(sred->psubcomm); 386 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 387 388 ierr = DMDAGetInfo(dm,&dim,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,&nsw,&bx,&by,&bz,&stencil);CHKERRQ(ierr); 389 ierr = DMDAGetInterpolationType(dm,&itype);CHKERRQ(ierr); 390 ierr = DMDAGetRefinementFactor(dm,&refine_x,&refine_y,&refine_z);CHKERRQ(ierr); 391 392 ctx->dmrepart = NULL; 393 _range_i_re = _range_j_re = _range_k_re = NULL; 394 /* Create DMDA on the child communicator */ 395 if (PCTelescope_isActiveRank(sred)) { 396 switch (dim) { 397 case 1: 398 ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (1D)\n");CHKERRQ(ierr); 399 /*ierr = DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 400 ny = nz = 1; 401 by = bz = DM_BOUNDARY_NONE; 402 break; 403 case 2: 404 ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (2D)\n");CHKERRQ(ierr); 405 /*ierr = DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,&ctx->dmrepart);CHKERRQ(ierr);*/ 406 nz = 1; 407 bz = DM_BOUNDARY_NONE; 408 break; 409 case 3: 410 ierr = PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (3D)\n");CHKERRQ(ierr); 411 /*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);*/ 412 break; 413 } 414 /* 415 The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append 416 a unique option prefix for the DM, thus I prefer to expose the contents of these API's here. 417 This allows users to control the partitioning of the subDM. 418 */ 419 ierr = DMDACreate(subcomm,&ctx->dmrepart);CHKERRQ(ierr); 420 /* Set unique option prefix name */ 421 ierr = KSPGetOptionsPrefix(sred->ksp,&prefix);CHKERRQ(ierr); 422 ierr = DMSetOptionsPrefix(ctx->dmrepart,prefix);CHKERRQ(ierr); 423 ierr = DMAppendOptionsPrefix(ctx->dmrepart,"repart_");CHKERRQ(ierr); 424 /* standard setup from DMDACreate{1,2,3}d() */ 425 ierr = DMSetDimension(ctx->dmrepart,dim);CHKERRQ(ierr); 426 ierr = DMDASetSizes(ctx->dmrepart,nx,ny,nz);CHKERRQ(ierr); 427 ierr = DMDASetNumProcs(ctx->dmrepart,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 428 ierr = DMDASetBoundaryType(ctx->dmrepart,bx,by,bz);CHKERRQ(ierr); 429 ierr = DMDASetDof(ctx->dmrepart,ndof);CHKERRQ(ierr); 430 ierr = DMDASetStencilType(ctx->dmrepart,stencil);CHKERRQ(ierr); 431 ierr = DMDASetStencilWidth(ctx->dmrepart,nsw);CHKERRQ(ierr); 432 ierr = DMDASetOwnershipRanges(ctx->dmrepart,NULL,NULL,NULL);CHKERRQ(ierr); 433 ierr = DMSetFromOptions(ctx->dmrepart);CHKERRQ(ierr); 434 ierr = DMSetUp(ctx->dmrepart);CHKERRQ(ierr); 435 /* Set refinement factors and interpolation type from the partent */ 436 ierr = DMDASetRefinementFactor(ctx->dmrepart,refine_x,refine_y,refine_z);CHKERRQ(ierr); 437 ierr = DMDASetInterpolationType(ctx->dmrepart,itype);CHKERRQ(ierr); 438 439 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); 440 ierr = DMDAGetOwnershipRanges(ctx->dmrepart,&_range_i_re,&_range_j_re,&_range_k_re);CHKERRQ(ierr); 441 442 ctx->dmrepart->ops->creatematrix = dm->ops->creatematrix; 443 ctx->dmrepart->ops->createdomaindecomposition = dm->ops->createdomaindecomposition; 444 } 445 446 /* generate ranges for repartitioned dm */ 447 /* note - assume rank 0 always participates */ 448 /* TODO: use a single MPI call */ 449 ierr = MPI_Bcast(&ctx->Mp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 450 ierr = MPI_Bcast(&ctx->Np_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 451 ierr = MPI_Bcast(&ctx->Pp_re,1,MPIU_INT,0,comm);CHKERRQ(ierr); 452 453 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); 454 455 if (_range_i_re) {ierr = PetscArraycpy(ctx->range_i_re,_range_i_re,ctx->Mp_re);CHKERRQ(ierr);} 456 if (_range_j_re) {ierr = PetscArraycpy(ctx->range_j_re,_range_j_re,ctx->Np_re);CHKERRQ(ierr);} 457 if (_range_k_re) {ierr = PetscArraycpy(ctx->range_k_re,_range_k_re,ctx->Pp_re);CHKERRQ(ierr);} 458 459 /* TODO: use a single MPI call */ 460 ierr = MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 461 ierr = MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm);CHKERRQ(ierr); 462 ierr = MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm);CHKERRQ(ierr); 463 464 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); 465 466 sum = 0; 467 for (k=0; k<ctx->Mp_re; k++) { 468 ctx->start_i_re[k] = sum; 469 sum += ctx->range_i_re[k]; 470 } 471 472 sum = 0; 473 for (k=0; k<ctx->Np_re; k++) { 474 ctx->start_j_re[k] = sum; 475 sum += ctx->range_j_re[k]; 476 } 477 478 sum = 0; 479 for (k=0; k<ctx->Pp_re; k++) { 480 ctx->start_k_re[k] = sum; 481 sum += ctx->range_k_re[k]; 482 } 483 484 /* attach repartitioned dm to child ksp */ 485 { 486 PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*); 487 void *dmksp_ctx; 488 489 ierr = DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);CHKERRQ(ierr); 490 491 /* attach dm to ksp on sub communicator */ 492 if (PCTelescope_isActiveRank(sred)) { 493 ierr = KSPSetDM(sred->ksp,ctx->dmrepart);CHKERRQ(ierr); 494 495 if (!dmksp_func || sred->ignore_kspcomputeoperators) { 496 ierr = KSPSetDMActive(sred->ksp,PETSC_FALSE);CHKERRQ(ierr); 497 } else { 498 /* sub ksp inherits dmksp_func and context provided by user */ 499 ierr = KSPSetComputeOperators(sred->ksp,dmksp_func,dmksp_ctx);CHKERRQ(ierr); 500 ierr = KSPSetDMActive(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 501 } 502 } 503 } 504 PetscFunctionReturn(0); 505 } 506 507 PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 508 { 509 PetscErrorCode ierr; 510 DM dm; 511 MPI_Comm comm; 512 Mat Pscalar,P; 513 PetscInt ndof; 514 PetscInt i,j,k,location,startI[3],endI[3],lenI[3],nx,ny,nz; 515 PetscInt sr,er,Mr; 516 Vec V; 517 518 PetscFunctionBegin; 519 ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-3D)\n");CHKERRQ(ierr); 520 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 521 522 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 523 ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 524 525 ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 526 ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 527 ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 528 ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 529 sr = sr / ndof; 530 er = er / ndof; 531 Mr = Mr / ndof; 532 533 ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 534 ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 535 ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 536 ierr = MatSeqAIJSetPreallocation(Pscalar,1,NULL);CHKERRQ(ierr); 537 ierr = MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);CHKERRQ(ierr); 538 539 ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]);CHKERRQ(ierr); 540 ierr = DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]);CHKERRQ(ierr); 541 endI[0] += startI[0]; 542 endI[1] += startI[1]; 543 endI[2] += startI[2]; 544 545 for (k=startI[2]; k<endI[2]; k++) { 546 for (j=startI[1]; j<endI[1]; j++) { 547 for (i=startI[0]; i<endI[0]; i++) { 548 PetscMPIInt rank_ijk_re,rank_reI[3]; 549 PetscInt s0_re; 550 PetscInt ii,jj,kk,local_ijk_re,mapped_ijk; 551 PetscInt lenI_re[3]; 552 553 location = (i - startI[0]) + (j - startI[1])*lenI[0] + (k - startI[2])*lenI[0]*lenI[1]; 554 ierr = _DMDADetermineRankFromGlobalIJK(3,i,j,k, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 555 ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 556 ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 557 &rank_reI[0],&rank_reI[1],&rank_reI[2],&rank_ijk_re);CHKERRQ(ierr); 558 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); 559 ii = i - ctx->start_i_re[ rank_reI[0] ]; 560 if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error ii"); 561 jj = j - ctx->start_j_re[ rank_reI[1] ]; 562 if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error jj"); 563 kk = k - ctx->start_k_re[ rank_reI[2] ]; 564 if (kk < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error kk"); 565 lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 566 lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 567 lenI_re[2] = ctx->range_k_re[ rank_reI[2] ]; 568 local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1]; 569 mapped_ijk = s0_re + local_ijk_re; 570 ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 571 } 572 } 573 } 574 ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 575 ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 576 ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 577 ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 578 ctx->permutation = P; 579 PetscFunctionReturn(0); 580 } 581 582 PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 583 { 584 PetscErrorCode ierr; 585 DM dm; 586 MPI_Comm comm; 587 Mat Pscalar,P; 588 PetscInt ndof; 589 PetscInt i,j,location,startI[2],endI[2],lenI[2],nx,ny,nz; 590 PetscInt sr,er,Mr; 591 Vec V; 592 593 PetscFunctionBegin; 594 ierr = PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-2D)\n");CHKERRQ(ierr); 595 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 596 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 597 ierr = DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 598 ierr = DMGetGlobalVector(dm,&V);CHKERRQ(ierr); 599 ierr = VecGetSize(V,&Mr);CHKERRQ(ierr); 600 ierr = VecGetOwnershipRange(V,&sr,&er);CHKERRQ(ierr); 601 ierr = DMRestoreGlobalVector(dm,&V);CHKERRQ(ierr); 602 sr = sr / ndof; 603 er = er / ndof; 604 Mr = Mr / ndof; 605 606 ierr = MatCreate(comm,&Pscalar);CHKERRQ(ierr); 607 ierr = MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);CHKERRQ(ierr); 608 ierr = MatSetType(Pscalar,MATAIJ);CHKERRQ(ierr); 609 ierr = MatSeqAIJSetPreallocation(Pscalar,1,NULL);CHKERRQ(ierr); 610 ierr = MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);CHKERRQ(ierr); 611 612 ierr = DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL);CHKERRQ(ierr); 613 ierr = DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL);CHKERRQ(ierr); 614 endI[0] += startI[0]; 615 endI[1] += startI[1]; 616 617 for (j=startI[1]; j<endI[1]; j++) { 618 for (i=startI[0]; i<endI[0]; i++) { 619 PetscMPIInt rank_ijk_re,rank_reI[3]; 620 PetscInt s0_re; 621 PetscInt ii,jj,local_ijk_re,mapped_ijk; 622 PetscInt lenI_re[3]; 623 624 location = (i - startI[0]) + (j - startI[1])*lenI[0]; 625 ierr = _DMDADetermineRankFromGlobalIJK(2,i,j,0, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, 626 ctx->start_i_re,ctx->start_j_re,ctx->start_k_re, 627 ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, 628 &rank_reI[0],&rank_reI[1],NULL,&rank_ijk_re);CHKERRQ(ierr); 629 630 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); 631 632 ii = i - ctx->start_i_re[ rank_reI[0] ]; 633 if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error ii"); 634 jj = j - ctx->start_j_re[ rank_reI[1] ]; 635 if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error jj"); 636 637 lenI_re[0] = ctx->range_i_re[ rank_reI[0] ]; 638 lenI_re[1] = ctx->range_j_re[ rank_reI[1] ]; 639 local_ijk_re = ii + jj * lenI_re[0]; 640 mapped_ijk = s0_re + local_ijk_re; 641 ierr = MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);CHKERRQ(ierr); 642 } 643 } 644 ierr = MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 645 ierr = MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 646 ierr = MatCreateMAIJ(Pscalar,ndof,&P);CHKERRQ(ierr); 647 ierr = MatDestroy(&Pscalar);CHKERRQ(ierr); 648 ctx->permutation = P; 649 PetscFunctionReturn(0); 650 } 651 652 PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx) 653 { 654 PetscErrorCode ierr; 655 Vec xred,yred,xtmp,x,xp; 656 VecScatter scatter; 657 IS isin; 658 Mat B; 659 PetscInt m,bs,st,ed; 660 MPI_Comm comm; 661 662 PetscFunctionBegin; 663 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 664 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 665 ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 666 ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 667 ierr = VecDuplicate(x,&xp);CHKERRQ(ierr); 668 m = 0; 669 xred = NULL; 670 yred = NULL; 671 if (PCTelescope_isActiveRank(sred)) { 672 ierr = DMCreateGlobalVector(ctx->dmrepart,&xred);CHKERRQ(ierr); 673 ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 674 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 675 ierr = ISCreateStride(comm,ed-st,st,1,&isin);CHKERRQ(ierr); 676 ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr); 677 } else { 678 ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 679 ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 680 } 681 ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 682 ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 683 ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 684 ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 685 ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 686 ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 687 sred->xred = xred; 688 sred->yred = yred; 689 sred->isin = isin; 690 sred->scatter = scatter; 691 sred->xtmp = xtmp; 692 693 ctx->xp = xp; 694 ierr = VecDestroy(&x);CHKERRQ(ierr); 695 PetscFunctionReturn(0); 696 } 697 698 PetscErrorCode PCTelescopeSetUp_dmda(PC pc,PC_Telescope sred) 699 { 700 PC_Telescope_DMDACtx *ctx; 701 PetscInt dim; 702 DM dm; 703 MPI_Comm comm; 704 PetscErrorCode ierr; 705 706 PetscFunctionBegin; 707 ierr = PetscInfo(pc,"PCTelescope: setup (DMDA)\n");CHKERRQ(ierr); 708 ierr = PetscNew(&ctx);CHKERRQ(ierr); 709 sred->dm_ctx = (void*)ctx; 710 711 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 712 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 713 ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 714 715 PCTelescopeSetUp_dmda_repart(pc,sred,ctx); 716 PCTelescopeSetUp_dmda_repart_coors(pc,sred,ctx); 717 switch (dim) { 718 case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided"); 719 break; 720 case 2: ierr = PCTelescopeSetUp_dmda_permutation_2d(pc,sred,ctx);CHKERRQ(ierr); 721 break; 722 case 3: ierr = PCTelescopeSetUp_dmda_permutation_3d(pc,sred,ctx);CHKERRQ(ierr); 723 break; 724 } 725 ierr = PCTelescopeSetUp_dmda_scatters(pc,sred,ctx);CHKERRQ(ierr); 726 PetscFunctionReturn(0); 727 } 728 729 PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 730 { 731 PetscErrorCode ierr; 732 PC_Telescope_DMDACtx *ctx; 733 MPI_Comm comm,subcomm; 734 Mat Bperm,Bred,B,P; 735 PetscInt nr,nc; 736 IS isrow,iscol; 737 Mat Blocal,*_Blocal; 738 739 PetscFunctionBegin; 740 ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (DMDA)\n");CHKERRQ(ierr); 741 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 742 subcomm = PetscSubcommChild(sred->psubcomm); 743 ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 744 745 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 746 ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 747 748 P = ctx->permutation; 749 ierr = MatPtAP(B,P,MAT_INITIAL_MATRIX,1.1,&Bperm);CHKERRQ(ierr); 750 751 /* Get submatrices */ 752 isrow = sred->isin; 753 ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 754 755 ierr = MatCreateSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 756 Blocal = *_Blocal; 757 Bred = NULL; 758 if (PCTelescope_isActiveRank(sred)) { 759 PetscInt mm; 760 761 if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;} 762 ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 763 /* ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred);CHKERRQ(ierr); */ 764 ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 765 } 766 *A = Bred; 767 768 ierr = ISDestroy(&iscol);CHKERRQ(ierr); 769 ierr = MatDestroy(&Bperm);CHKERRQ(ierr); 770 ierr = MatDestroyMatrices(1,&_Blocal);CHKERRQ(ierr); 771 PetscFunctionReturn(0); 772 } 773 774 PetscErrorCode PCTelescopeMatCreate_dmda(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 775 { 776 PetscErrorCode ierr; 777 DM dm; 778 PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*); 779 void *dmksp_ctx; 780 781 PetscFunctionBegin; 782 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 783 ierr = DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);CHKERRQ(ierr); 784 /* We assume that dmksp_func = NULL, is equivalent to dmActive = PETSC_FALSE */ 785 if (dmksp_func && !sred->ignore_kspcomputeoperators) { 786 DM dmrepart; 787 Mat Ak; 788 789 *A = NULL; 790 if (PCTelescope_isActiveRank(sred)) { 791 ierr = KSPGetDM(sred->ksp,&dmrepart);CHKERRQ(ierr); 792 if (reuse == MAT_INITIAL_MATRIX) { 793 ierr = DMCreateMatrix(dmrepart,&Ak);CHKERRQ(ierr); 794 *A = Ak; 795 } else if (reuse == MAT_REUSE_MATRIX) { 796 Ak = *A; 797 } 798 /* 799 There is no need to explicitly assemble the operator now, 800 the sub-KSP will call the method provided to KSPSetComputeOperators() during KSPSetUp() 801 */ 802 } 803 } else { 804 ierr = PCTelescopeMatCreate_dmda_dmactivefalse(pc,sred,reuse,A);CHKERRQ(ierr); 805 } 806 PetscFunctionReturn(0); 807 } 808 809 PetscErrorCode PCTelescopeSubNullSpaceCreate_dmda_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 810 { 811 PetscErrorCode ierr; 812 PetscBool has_const; 813 PetscInt i,k,n = 0; 814 const Vec *vecs; 815 Vec *sub_vecs = NULL; 816 MPI_Comm subcomm; 817 PC_Telescope_DMDACtx *ctx; 818 819 PetscFunctionBegin; 820 ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 821 subcomm = PetscSubcommChild(sred->psubcomm); 822 ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 823 824 if (PCTelescope_isActiveRank(sred)) { 825 /* create new vectors */ 826 if (n) { 827 ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 828 } 829 } 830 831 /* copy entries */ 832 for (k=0; k<n; k++) { 833 const PetscScalar *x_array; 834 PetscScalar *LA_sub_vec; 835 PetscInt st,ed; 836 837 /* permute vector into ordering associated with re-partitioned dmda */ 838 ierr = MatMultTranspose(ctx->permutation,vecs[k],ctx->xp);CHKERRQ(ierr); 839 840 /* pull in vector x->xtmp */ 841 ierr = VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 842 ierr = VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 843 844 /* copy vector entries into xred */ 845 ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 846 if (sub_vecs) { 847 if (sub_vecs[k]) { 848 ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 849 ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 850 for (i=0; i<ed-st; i++) { 851 LA_sub_vec[i] = x_array[i]; 852 } 853 ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 854 } 855 } 856 ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 857 } 858 859 if (PCTelescope_isActiveRank(sred)) { 860 /* create new (near) nullspace for redundant object */ 861 ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr); 862 ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 863 if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 864 if (nullspace->rmctx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope"); 865 } 866 PetscFunctionReturn(0); 867 } 868 869 PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat) 870 { 871 PetscErrorCode ierr; 872 Mat B; 873 874 PetscFunctionBegin; 875 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 876 { 877 MatNullSpace nullspace,sub_nullspace; 878 ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 879 if (nullspace) { 880 ierr = PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n");CHKERRQ(ierr); 881 ierr = PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr); 882 if (PCTelescope_isActiveRank(sred)) { 883 ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 884 ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr); 885 } 886 } 887 } 888 { 889 MatNullSpace nearnullspace,sub_nearnullspace; 890 ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr); 891 if (nearnullspace) { 892 ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (DMDA)\n");CHKERRQ(ierr); 893 ierr = PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr); 894 if (PCTelescope_isActiveRank(sred)) { 895 ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr); 896 ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr); 897 } 898 } 899 } 900 PetscFunctionReturn(0); 901 } 902 903 PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y) 904 { 905 PC_Telescope sred = (PC_Telescope)pc->data; 906 PetscErrorCode ierr; 907 Mat perm; 908 Vec xtmp,xp,xred,yred; 909 PetscInt i,st,ed; 910 VecScatter scatter; 911 PetscScalar *array; 912 const PetscScalar *x_array; 913 PC_Telescope_DMDACtx *ctx; 914 915 ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 916 xtmp = sred->xtmp; 917 scatter = sred->scatter; 918 xred = sred->xred; 919 yred = sred->yred; 920 perm = ctx->permutation; 921 xp = ctx->xp; 922 923 PetscFunctionBegin; 924 ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 925 926 /* permute vector into ordering associated with re-partitioned dmda */ 927 ierr = MatMultTranspose(perm,x,xp);CHKERRQ(ierr); 928 929 /* pull in vector x->xtmp */ 930 ierr = VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 931 ierr = VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 932 933 /* copy vector entires into xred */ 934 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 935 if (xred) { 936 PetscScalar *LA_xred; 937 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 938 939 ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 940 for (i=0; i<ed-st; i++) { 941 LA_xred[i] = x_array[i]; 942 } 943 ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 944 } 945 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 946 947 /* solve */ 948 if (PCTelescope_isActiveRank(sred)) { 949 ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 950 ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr); 951 } 952 953 /* return vector */ 954 ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 955 if (yred) { 956 const PetscScalar *LA_yred; 957 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 958 ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 959 for (i=0; i<ed-st; i++) { 960 array[i] = LA_yred[i]; 961 } 962 ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 963 } 964 ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 965 ierr = VecScatterBegin(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 966 ierr = VecScatterEnd(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 967 ierr = MatMult(perm,xp,y);CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 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) 972 { 973 PC_Telescope sred = (PC_Telescope)pc->data; 974 PetscErrorCode ierr; 975 Mat perm; 976 Vec xtmp,xp,yred; 977 PetscInt i,st,ed; 978 VecScatter scatter; 979 const PetscScalar *x_array; 980 PetscBool default_init_guess_value = PETSC_FALSE; 981 PC_Telescope_DMDACtx *ctx; 982 983 PetscFunctionBegin; 984 ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 985 xtmp = sred->xtmp; 986 scatter = sred->scatter; 987 yred = sred->yred; 988 perm = ctx->permutation; 989 xp = ctx->xp; 990 991 if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope_dmda only supports max_it = 1"); 992 *reason = (PCRichardsonConvergedReason)0; 993 994 if (!zeroguess) { 995 ierr = PetscInfo(pc,"PCTelescopeDMDA: Scattering y for non-zero-initial guess\n");CHKERRQ(ierr); 996 /* permute vector into ordering associated with re-partitioned dmda */ 997 ierr = MatMultTranspose(perm,y,xp);CHKERRQ(ierr); 998 999 /* pull in vector x->xtmp */ 1000 ierr = VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1001 ierr = VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1002 1003 /* copy vector entires into xred */ 1004 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 1005 if (yred) { 1006 PetscScalar *LA_yred; 1007 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 1008 ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 1009 for (i=0; i<ed-st; i++) { 1010 LA_yred[i] = x_array[i]; 1011 } 1012 ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 1013 } 1014 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 1015 } 1016 1017 if (PCTelescope_isActiveRank(sred)) { 1018 ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr); 1019 if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 1020 } 1021 1022 ierr = PCApply_Telescope_dmda(pc,x,y);CHKERRQ(ierr); 1023 1024 if (PCTelescope_isActiveRank(sred)) { 1025 ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr); 1026 } 1027 1028 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1029 *outits = 1; 1030 PetscFunctionReturn(0); 1031 } 1032 1033 PetscErrorCode PCReset_Telescope_dmda(PC pc) 1034 { 1035 PetscErrorCode ierr; 1036 PC_Telescope sred = (PC_Telescope)pc->data; 1037 PC_Telescope_DMDACtx *ctx; 1038 1039 PetscFunctionBegin; 1040 ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx; 1041 ierr = VecDestroy(&ctx->xp);CHKERRQ(ierr); 1042 ierr = MatDestroy(&ctx->permutation);CHKERRQ(ierr); 1043 ierr = DMDestroy(&ctx->dmrepart);CHKERRQ(ierr); 1044 ierr = PetscFree3(ctx->range_i_re,ctx->range_j_re,ctx->range_k_re);CHKERRQ(ierr); 1045 ierr = PetscFree3(ctx->start_i_re,ctx->start_j_re,ctx->start_k_re);CHKERRQ(ierr); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 PetscErrorCode DMView_DA_Short_3d(DM dm,PetscViewer v) 1050 { 1051 PetscInt M,N,P,m,n,p,ndof,nsw; 1052 MPI_Comm comm; 1053 PetscMPIInt size; 1054 const char* prefix; 1055 PetscErrorCode ierr; 1056 1057 PetscFunctionBegin; 1058 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1059 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1060 ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 1061 ierr = DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1062 if (prefix) {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 1063 else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 1064 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); 1065 PetscFunctionReturn(0); 1066 } 1067 1068 PetscErrorCode DMView_DA_Short_2d(DM dm,PetscViewer v) 1069 { 1070 PetscInt M,N,m,n,ndof,nsw; 1071 MPI_Comm comm; 1072 PetscMPIInt size; 1073 const char* prefix; 1074 PetscErrorCode ierr; 1075 1076 PetscFunctionBegin; 1077 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1078 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1079 ierr = DMGetOptionsPrefix(dm,&prefix);CHKERRQ(ierr); 1080 ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1081 if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object: (%s) %d MPI processes\n",prefix,size);CHKERRQ(ierr);} 1082 else {ierr = PetscViewerASCIIPrintf(v,"DMDA Object: %d MPI processes\n",size);CHKERRQ(ierr);} 1083 ierr = PetscViewerASCIIPrintf(v," M %D N %D m %D n %D dof %D overlap %D\n",M,N,m,n,ndof,nsw);CHKERRQ(ierr); 1084 PetscFunctionReturn(0); 1085 } 1086 1087 PetscErrorCode DMView_DA_Short(DM dm,PetscViewer v) 1088 { 1089 PetscErrorCode ierr; 1090 PetscInt dim; 1091 1092 PetscFunctionBegin; 1093 ierr = DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1094 switch (dim) { 1095 case 2: ierr = DMView_DA_Short_2d(dm,v);CHKERRQ(ierr); 1096 break; 1097 case 3: ierr = DMView_DA_Short_3d(dm,v);CHKERRQ(ierr); 1098 break; 1099 } 1100 PetscFunctionReturn(0); 1101 } 1102 1103