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