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