1 2 /* 3 4 Defines a telescoping preconditioner 5 6 Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 7 creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC). 8 9 During PCSetup, the B operator is scattered onto c'. 10 Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 11 Then KSPSolve() is executed on the c' communicator. 12 13 The preconditioner is deemed telescopint as it only calls KSPSolve() on a single 14 sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 15 This means there will be MPI processes within c, which will be idle during the application of this preconditioner. 16 17 Comments: 18 - The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 19 creation routine. We run the sub KSP on only the ranks within the communicator which have a color equal to zero. 20 21 - The telescoping preconditioner is aware of nullspaces which are attached to the only B operator. 22 In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into 23 a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c') 24 25 - The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D - 26 1D support for 1D DMDAs is not provided), a new DMDA is created on c' (e.g. it is re-partitioned), and this new DM 27 is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support 28 for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 29 30 - By default, B' is defined by simply fusing rows from different MPI processes 31 32 - When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B 33 into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the 34 locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 35 36 Limitations/improvements 37 - VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage. 38 39 - The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P, 40 and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). I opted to use P^T.A.P as it appears 41 VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a 42 consistent manner. 43 44 - Mapping of vectors is performed this way 45 Suppose the parent comm size was 4, and we set a reduction factor of 2, thus would give a comm size on c' of 2. 46 Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2. 47 We perform the scatter to the sub-comm in the following way, 48 [1] Given a vector x defined on comm c 49 50 rank(c) : _________ 0 ______ ________ 1 _______ ________ 2 _____________ ___________ 3 __________ 51 x : [0, 1, 2, 3, 4, 5] [6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17] [18, 19, 20, 21, 22, 23] 52 53 scatter to xtmp defined also on comm c so that we have the following values 54 55 rank(c) : ___________________ 0 ________________ _1_ ______________________ 2 _______________________ __3_ 56 xtmp : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [ 6 ] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] [ 18 ] 57 58 The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') are filled with the first bs values living on that rank, 59 where bs is the block size of the vector. 60 61 [2] Copy the value from rank 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'. 62 Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 63 64 rank(c') : ___________________ 0 _______________ ______________________ 1 _____________________ 65 xred : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] 66 67 */ 68 69 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 70 #include <petscksp.h> /*I "petscksp.h" I*/ 71 #include <petscdm.h> /*I "petscdm.h" I*/ 72 73 #include "telescope.h" 74 75 /* 76 PCTelescopeSetUp_default() 77 PCTelescopeMatCreate_default() 78 79 default 80 81 // scatter in 82 x(comm) -> xtmp(comm) 83 84 xred(subcomm) <- xtmp 85 yred(subcomm) 86 87 yred(subcomm) --> xtmp 88 89 // scatter out 90 xtmp(comm) -> y(comm) 91 */ 92 93 PetscBool isActiveRank(PetscSubcomm scomm) 94 { 95 if (scomm->color == 0) { return PETSC_TRUE; } 96 else { return PETSC_FALSE; } 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "private_PCTelescopeGetSubDM" 101 DM private_PCTelescopeGetSubDM(PC_Telescope sred) 102 { 103 DM subdm = NULL; 104 105 if (!isActiveRank(sred->psubcomm)) { subdm = NULL; } 106 else { 107 switch (sred->sr_type) { 108 case TELESCOPE_DEFAULT: subdm = NULL; 109 break; 110 case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 111 break; 112 case TELESCOPE_DMPLEX: subdm = NULL; 113 break; 114 } 115 } 116 return(subdm); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "PCTelescopeSetUp_default" 121 PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 122 { 123 PetscErrorCode ierr; 124 PetscInt m,M,bs,st,ed; 125 Vec x,xred,yred,xtmp; 126 Mat B; 127 MPI_Comm comm,subcomm; 128 VecScatter scatter; 129 IS isin; 130 131 PetscFunctionBegin; 132 ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr); 133 comm = PetscSubcommParent(sred->psubcomm); 134 subcomm = PetscSubcommChild(sred->psubcomm); 135 136 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 137 ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 138 ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 139 ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 140 141 xred = NULL; 142 m = bs; 143 if (isActiveRank(sred->psubcomm)) { 144 ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr); 145 ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr); 146 ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr); 147 ierr = VecSetFromOptions(xred);CHKERRQ(ierr); 148 ierr = VecGetLocalSize(xred,&m); 149 } 150 151 yred = NULL; 152 if (isActiveRank(sred->psubcomm)) { 153 ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 154 } 155 156 ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 157 ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 158 ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 159 ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 160 161 if (isActiveRank(sred->psubcomm)) { 162 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 163 ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr); 164 } else { 165 ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 166 ierr = ISCreateStride(comm,bs,st,1,&isin);CHKERRQ(ierr); 167 } 168 ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 169 170 ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 171 172 sred->isin = isin; 173 sred->scatter = scatter; 174 sred->xred = xred; 175 sred->yred = yred; 176 sred->xtmp = xtmp; 177 ierr = VecDestroy(&x);CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "PCTelescopeMatCreate_default" 183 PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 184 { 185 PetscErrorCode ierr; 186 MPI_Comm comm,subcomm; 187 Mat Bred,B; 188 PetscInt nr,nc; 189 IS isrow,iscol; 190 Mat Blocal,*_Blocal; 191 192 PetscFunctionBegin; 193 ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr); 194 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 195 subcomm = PetscSubcommChild(sred->psubcomm); 196 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 197 ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 198 isrow = sred->isin; 199 ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 200 ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 201 Blocal = *_Blocal; 202 ierr = PetscFree(_Blocal);CHKERRQ(ierr); 203 Bred = NULL; 204 if (isActiveRank(sred->psubcomm)) { 205 PetscInt mm; 206 207 if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 208 209 ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 210 ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 211 } 212 *A = Bred; 213 ierr = ISDestroy(&iscol);CHKERRQ(ierr); 214 ierr = MatDestroy(&Blocal);CHKERRQ(ierr); 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default" 220 PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 221 { 222 PetscErrorCode ierr; 223 MatNullSpace nullspace,sub_nullspace; 224 Mat A,B; 225 PetscBool has_const; 226 PetscInt i,k,n; 227 const Vec *vecs; 228 Vec *sub_vecs; 229 MPI_Comm subcomm; 230 231 PetscFunctionBegin; 232 ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr); 233 ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 234 if (!nullspace) return(0); 235 236 ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr); 237 subcomm = PetscSubcommChild(sred->psubcomm); 238 ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 239 240 if (isActiveRank(sred->psubcomm)) { 241 if (n) { 242 ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 243 } 244 } 245 246 /* copy entries */ 247 for (k=0; k<n; k++) { 248 const PetscScalar *x_array; 249 PetscScalar *LA_sub_vec; 250 PetscInt st,ed,bs; 251 252 /* pull in vector x->xtmp */ 253 ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 254 ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 255 /* copy vector entires into xred */ 256 ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr); 257 ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 258 if (sub_vecs[k]) { 259 ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 260 ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 261 for (i=0; i<ed-st; i++) { 262 LA_sub_vec[i] = x_array[i]; 263 } 264 ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 265 } 266 ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 267 } 268 269 if (isActiveRank(sred->psubcomm)) { 270 /* create new nullspace for redundant object */ 271 ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr); 272 /* attach redundant nullspace to Bred */ 273 ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 274 ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "PCView_Telescope" 281 static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 282 { 283 PC_Telescope sred = (PC_Telescope)pc->data; 284 PetscErrorCode ierr; 285 PetscBool iascii,isstring; 286 PetscViewer subviewer; 287 288 PetscFunctionBegin; 289 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 290 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 291 if (iascii) { 292 if (!sred->psubcomm) { 293 ierr = PetscViewerASCIIPrintf(viewer," Telescope: preconditioner not yet setup\n");CHKERRQ(ierr); 294 } else { 295 MPI_Comm comm,subcomm; 296 PetscMPIInt comm_size,subcomm_size; 297 DM dm,subdm; 298 299 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 300 subdm = private_PCTelescopeGetSubDM(sred); 301 comm = PetscSubcommParent(sred->psubcomm); 302 subcomm = PetscSubcommChild(sred->psubcomm); 303 ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 304 ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 305 306 ierr = PetscViewerASCIIPrintf(viewer," Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 307 ierr = PetscViewerASCIIPrintf(viewer," Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 308 ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 309 if (isActiveRank(sred->psubcomm)) { 310 ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 311 312 if (dm && sred->ignore_dm) { 313 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: ignoring DM\n");CHKERRQ(ierr); 314 } 315 switch (sred->sr_type) { 316 case TELESCOPE_DEFAULT: 317 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: using default setup\n");CHKERRQ(ierr); 318 break; 319 case TELESCOPE_DMDA: 320 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMDA detected\n");CHKERRQ(ierr); 321 ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr); 322 break; 323 case TELESCOPE_DMPLEX: 324 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMPLEX detected\n");CHKERRQ(ierr); 325 break; 326 } 327 328 ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 329 ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 330 } 331 ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 332 } 333 } 334 PetscFunctionReturn(0); 335 } 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "PCSetUp_Telescope" 339 static PetscErrorCode PCSetUp_Telescope(PC pc) 340 { 341 PC_Telescope sred = (PC_Telescope)pc->data; 342 PetscErrorCode ierr; 343 MPI_Comm comm,subcomm=0; 344 PCTelescopeType sr_type; 345 346 PetscFunctionBegin; 347 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 348 349 /* subcomm definition */ 350 if (!pc->setupcalled) { 351 if (!sred->psubcomm) { 352 ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 353 ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 354 ierr = PetscSubcommSetType(sred->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 355 /* disable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 356 /* ierr = PetscSubcommSetFromOptions(sred->psubcomm);CHKERRQ(ierr); */ 357 ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 358 /* create a new PC that processors in each subcomm have copy of */ 359 } 360 } 361 subcomm = PetscSubcommChild(sred->psubcomm); 362 363 /* internal KSP */ 364 if (!pc->setupcalled) { 365 const char *prefix; 366 367 if (isActiveRank(sred->psubcomm)) { 368 ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 369 ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 370 ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 371 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 372 ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 373 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 374 ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 375 ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 376 } 377 } 378 /* Determine type of setup/update */ 379 if (!pc->setupcalled) { 380 PetscBool has_dm,same; 381 DM dm; 382 383 sr_type = TELESCOPE_DEFAULT; 384 has_dm = PETSC_FALSE; 385 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 386 if (dm) { has_dm = PETSC_TRUE; } 387 if (has_dm) { 388 /* check for dmda */ 389 ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 390 if (same) { 391 ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 392 sr_type = TELESCOPE_DMDA; 393 } 394 /* check for dmplex */ 395 ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 396 if (same) { 397 PetscInfo(pc,"PCTelescope: found DMPLEX\n"); 398 sr_type = TELESCOPE_DMPLEX; 399 } 400 } 401 402 if (sred->ignore_dm) { 403 PetscInfo(pc,"PCTelescope: ignore DM\n"); 404 sr_type = TELESCOPE_DEFAULT; 405 } 406 sred->sr_type = sr_type; 407 } else { 408 sr_type = sred->sr_type; 409 } 410 411 /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */ 412 switch (sr_type) { 413 case TELESCOPE_DEFAULT: 414 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 415 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 416 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 417 sred->pctelescope_reset_type = NULL; 418 break; 419 case TELESCOPE_DMDA: 420 pc->ops->apply = PCApply_Telescope_dmda; 421 sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 422 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 423 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 424 sred->pctelescope_reset_type = PCReset_Telescope_dmda; 425 break; 426 case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available"); 427 break; 428 default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided"); 429 break; 430 } 431 432 /* setup */ 433 if (sred->pctelescope_setup_type) { 434 ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 435 } 436 /* update */ 437 if (!pc->setupcalled) { 438 if (sred->pctelescope_matcreate_type) { 439 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 440 } 441 if (sred->pctelescope_matnullspacecreate_type) { 442 ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 443 } 444 } else { 445 if (sred->pctelescope_matcreate_type) { 446 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 447 } 448 } 449 450 /* common - no construction */ 451 if (isActiveRank(sred->psubcomm)) { 452 ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 453 if (pc->setfromoptionscalled && !pc->setupcalled){ 454 ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 455 } 456 } 457 PetscFunctionReturn(0); 458 } 459 460 #undef __FUNCT__ 461 #define __FUNCT__ "PCApply_Telescope" 462 static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 463 { 464 PC_Telescope sred = (PC_Telescope)pc->data; 465 PetscErrorCode ierr; 466 Vec xtmp,xred,yred; 467 PetscInt i,st,ed,bs; 468 VecScatter scatter; 469 PetscScalar *array; 470 const PetscScalar *x_array; 471 472 PetscFunctionBegin; 473 xtmp = sred->xtmp; 474 scatter = sred->scatter; 475 xred = sred->xred; 476 yred = sred->yred; 477 478 /* pull in vector x->xtmp */ 479 ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 480 ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 481 482 /* copy vector entires into xred */ 483 ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 484 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 485 if (xred) { 486 PetscScalar *LA_xred; 487 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 488 ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 489 for (i=0; i<ed-st; i++) { 490 LA_xred[i] = x_array[i]; 491 } 492 ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 493 } 494 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 495 /* solve */ 496 if (isActiveRank(sred->psubcomm)) { 497 ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 498 } 499 /* return vector */ 500 ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 501 ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 502 if (yred) { 503 const PetscScalar *LA_yred; 504 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 505 ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 506 for (i=0; i<ed-st; i++) { 507 array[i] = LA_yred[i]; 508 } 509 ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 510 } else { 511 for (i=0; i<bs; i++) { 512 array[i] = 0.0; 513 } 514 } 515 ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 516 ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 517 ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 518 PetscFunctionReturn(0); 519 } 520 521 #undef __FUNCT__ 522 #define __FUNCT__ "PCReset_Telescope" 523 static PetscErrorCode PCReset_Telescope(PC pc) 524 { 525 PC_Telescope sred = (PC_Telescope)pc->data; 526 PetscErrorCode ierr; 527 528 ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 529 ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 530 ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); 531 ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); 532 ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); 533 ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); 534 ierr = KSPReset(sred->ksp);CHKERRQ(ierr); 535 if (sred->pctelescope_reset_type) { 536 ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 537 } 538 PetscFunctionReturn(0); 539 } 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "PCDestroy_Telescope" 543 static PetscErrorCode PCDestroy_Telescope(PC pc) 544 { 545 PC_Telescope sred = (PC_Telescope)pc->data; 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 550 ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); 551 ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 552 ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr); 553 ierr = PetscFree(pc->data);CHKERRQ(ierr); 554 PetscFunctionReturn(0); 555 } 556 557 #undef __FUNCT__ 558 #define __FUNCT__ "PCSetFromOptions_Telescope" 559 static PetscErrorCode PCSetFromOptions_Telescope(PetscOptions *PetscOptionsObject,PC pc) 560 { 561 PC_Telescope sred = (PC_Telescope)pc->data; 562 PetscErrorCode ierr; 563 MPI_Comm comm; 564 PetscMPIInt size; 565 566 PetscFunctionBegin; 567 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 568 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 569 ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 570 ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 571 if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 572 ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 573 ierr = PetscOptionsTail();CHKERRQ(ierr); 574 PetscFunctionReturn(0); 575 } 576 577 /* PC simplementation specific API's */ 578 579 static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 580 { 581 PC_Telescope red = (PC_Telescope)pc->data; 582 PetscFunctionBegin; 583 if (ksp) *ksp = red->ksp; 584 PetscFunctionReturn(0); 585 } 586 587 static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 588 { 589 PC_Telescope red = (PC_Telescope)pc->data; 590 PetscFunctionBegin; 591 if (fact) *fact = red->redfactor; 592 PetscFunctionReturn(0); 593 } 594 595 static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 596 { 597 PC_Telescope red = (PC_Telescope)pc->data; 598 PetscMPIInt size; 599 PetscErrorCode ierr; 600 601 PetscFunctionBegin; 602 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 603 if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 604 if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 605 red->redfactor = fact; 606 PetscFunctionReturn(0); 607 } 608 609 static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 610 { 611 PC_Telescope red = (PC_Telescope)pc->data; 612 PetscFunctionBegin; 613 if (v) *v = red->ignore_dm; 614 PetscFunctionReturn(0); 615 } 616 static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 617 { 618 PC_Telescope red = (PC_Telescope)pc->data; 619 PetscFunctionBegin; 620 red->ignore_dm = v; 621 PetscFunctionReturn(0); 622 } 623 624 static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 625 { 626 PC_Telescope red = (PC_Telescope)pc->data; 627 PetscFunctionBegin; 628 *dm = private_PCTelescopeGetSubDM(red); 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 634 635 Not Collective 636 637 Input Parameter: 638 . pc - the preconditioner context 639 640 Output Parameter: 641 . subksp - the KSP defined the smaller set of processes 642 643 Level: advanced 644 645 .keywords: PC, telescopting solve 646 @*/ 647 PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 648 { 649 PetscErrorCode ierr; 650 PetscFunctionBegin; 651 ierr = PetscTryMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 655 /*@ 656 PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 657 658 Not Collective 659 660 Input Parameter: 661 . pc - the preconditioner context 662 663 Output Parameter: 664 . fact - the reduction factor 665 666 Level: advanced 667 668 .keywords: PC, telescoping solve 669 @*/ 670 PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 671 { 672 PetscErrorCode ierr; 673 PetscFunctionBegin; 674 ierr = PetscTryMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 675 PetscFunctionReturn(0); 676 } 677 678 /*@ 679 PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 680 681 Not Collective 682 683 Input Parameter: 684 . pc - the preconditioner context 685 686 Output Parameter: 687 . fact - the reduction factor 688 689 Level: advanced 690 691 .keywords: PC, telescoping solve 692 @*/ 693 PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 694 { 695 PetscErrorCode ierr; 696 PetscFunctionBegin; 697 ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 /*@ 702 PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 703 704 Not Collective 705 706 Input Parameter: 707 . pc - the preconditioner context 708 709 Output Parameter: 710 . v - the flag 711 712 Level: advanced 713 714 .keywords: PC, telescoping solve 715 @*/ 716 PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 717 { 718 PetscErrorCode ierr; 719 PetscFunctionBegin; 720 ierr = PetscTryMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 721 PetscFunctionReturn(0); 722 } 723 724 /*@ 725 PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 726 727 Not Collective 728 729 Input Parameter: 730 . pc - the preconditioner context 731 732 Output Parameter: 733 . v - Use PETSC_TRUE to ignore any DM 734 735 Level: advanced 736 737 .keywords: PC, telescoping solve 738 @*/ 739 PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 740 { 741 PetscErrorCode ierr; 742 PetscFunctionBegin; 743 ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 744 PetscFunctionReturn(0); 745 } 746 747 /*@ 748 PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 749 750 Not Collective 751 752 Input Parameter: 753 . pc - the preconditioner context 754 755 Output Parameter: 756 . subdm - The re-partitioned DM 757 758 Level: advanced 759 760 .keywords: PC, telescoping solve 761 @*/ 762 PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 763 { 764 PetscErrorCode ierr; 765 PetscFunctionBegin; 766 ierr = PetscTryMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 767 PetscFunctionReturn(0); 768 } 769 770 /* -------------------------------------------------------------------------------------*/ 771 /*MC 772 PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 773 774 Options Database: 775 + -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and 776 use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes 777 - -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored 778 779 Level: advanced 780 781 Notes: 782 The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 783 Both the B mat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 784 Currently only support for re-partitioning a DMDA is provided. 785 Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator. 786 KSPSetComputeOperators() is not propagated to the sub KSP. 787 Currently there is no support for the flag -pc_use_amat 788 789 Optimizations: 790 (i) VecPlaceArray() could be used for scatters between the vectors defined on different sized communicators, thereby slightly reducing memory footprint; 791 (ii) Memory re-use and faster set-up would follow if the result of P^T.A.P was not re-allocated each time PCSetUp_Telescope() was called (DMDA). 792 793 Contributed by Dave May 794 795 .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), 796 PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), 797 PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 798 M*/ 799 #undef __FUNCT__ 800 #define __FUNCT__ "PCCreate_Telescope" 801 PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 802 { 803 PetscErrorCode ierr; 804 struct _PC_Telescope *sred; 805 806 PetscFunctionBegin; 807 ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 808 sred->redfactor = 1; 809 sred->ignore_dm = PETSC_FALSE; 810 pc->data = (void*)sred; 811 812 pc->ops->apply = PCApply_Telescope; 813 pc->ops->applytranspose = NULL; 814 pc->ops->setup = PCSetUp_Telescope; 815 pc->ops->destroy = PCDestroy_Telescope; 816 pc->ops->reset = PCReset_Telescope; 817 pc->ops->setfromoptions = PCSetFromOptions_Telescope; 818 pc->ops->view = PCView_Telescope; 819 820 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 821 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 822 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 823 sred->pctelescope_reset_type = NULL; 824 825 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 826 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 827 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 828 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 829 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 830 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 831 PetscFunctionReturn(0); 832 } 833