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 sub_vecs = NULL; 242 /* create new vectors */ 243 if (n != 0) { 244 PetscMalloc(sizeof(Vec)*n,&sub_vecs); 245 for (k=0; k<n; k++) { 246 ierr = VecDuplicate(sred->xred,&sub_vecs[k]);CHKERRQ(ierr); 247 } 248 } 249 } 250 251 /* copy entries */ 252 for (k=0; k<n; k++) { 253 const PetscScalar *x_array; 254 PetscScalar *LA_sub_vec; 255 PetscInt st,ed,bs; 256 257 /* pull in vector x->xtmp */ 258 ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 259 ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 260 /* copy vector entires into xred */ 261 ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr); 262 ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 263 if (sub_vecs[k]) { 264 ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 265 ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 266 for (i=0; i<ed-st; i++) { 267 LA_sub_vec[i] = x_array[i]; 268 } 269 ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 270 } 271 ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 272 } 273 274 if (isActiveRank(sred->psubcomm)) { 275 /* create new nullspace for redundant object */ 276 ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr); 277 /* attach redundant nullspace to Bred */ 278 ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 279 280 for (k=0; k<n; k++) { 281 ierr = VecDestroy(&sub_vecs[k]);CHKERRQ(ierr); 282 } 283 if (sub_vecs) PetscFree(sub_vecs); 284 } 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "PCView_Telescope" 290 static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 291 { 292 PC_Telescope sred = (PC_Telescope)pc->data; 293 PetscErrorCode ierr; 294 PetscBool iascii,isstring; 295 PetscViewer subviewer; 296 297 PetscFunctionBegin; 298 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 299 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 300 if (iascii) { 301 if (!sred->psubcomm) { 302 ierr = PetscViewerASCIIPrintf(viewer," Telescope: preconditioner not yet setup\n");CHKERRQ(ierr); 303 } else { 304 MPI_Comm comm,subcomm; 305 PetscMPIInt comm_size,subcomm_size; 306 DM dm,subdm; 307 308 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 309 subdm = private_PCTelescopeGetSubDM(sred); 310 comm = PetscSubcommParent(sred->psubcomm); 311 subcomm = PetscSubcommChild(sred->psubcomm); 312 ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 313 ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 314 315 ierr = PetscViewerASCIIPrintf(viewer," Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 316 ierr = PetscViewerASCIIPrintf(viewer," Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 317 ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 318 if (isActiveRank(sred->psubcomm)) { 319 ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 320 321 if (dm && sred->ignore_dm) { 322 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: ignoring DM\n");CHKERRQ(ierr); 323 } 324 switch (sred->sr_type) { 325 case TELESCOPE_DEFAULT: 326 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: using default setup\n");CHKERRQ(ierr); 327 break; 328 case TELESCOPE_DMDA: 329 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMDA detected\n");CHKERRQ(ierr); 330 ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr); 331 break; 332 case TELESCOPE_DMPLEX: 333 ierr = PetscViewerASCIIPrintf(subviewer," Telescope: DMPLEX detected\n");CHKERRQ(ierr); 334 break; 335 } 336 337 ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 338 ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 339 } 340 ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 341 } 342 } 343 PetscFunctionReturn(0); 344 } 345 346 #undef __FUNCT__ 347 #define __FUNCT__ "PCSetUp_Telescope" 348 static PetscErrorCode PCSetUp_Telescope(PC pc) 349 { 350 PC_Telescope sred = (PC_Telescope)pc->data; 351 PetscErrorCode ierr; 352 MPI_Comm comm,subcomm=0; 353 PCTelescopeType sr_type; 354 355 PetscFunctionBegin; 356 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 357 358 /* subcomm definition */ 359 if (!pc->setupcalled) { 360 if (!sred->psubcomm) { 361 ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 362 ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 363 ierr = PetscSubcommSetType(sred->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 364 /* disable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 365 /* ierr = PetscSubcommSetFromOptions(sred->psubcomm);CHKERRQ(ierr); */ 366 ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 367 /* create a new PC that processors in each subcomm have copy of */ 368 } 369 } 370 subcomm = PetscSubcommChild(sred->psubcomm); 371 372 /* internal KSP */ 373 if (!pc->setupcalled) { 374 const char *prefix; 375 376 if (isActiveRank(sred->psubcomm)) { 377 ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 378 ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 379 ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 380 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 381 ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 382 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 383 ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 384 ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 385 } 386 } 387 /* Determine type of setup/update */ 388 if (!pc->setupcalled) { 389 PetscBool has_dm,same; 390 DM dm; 391 392 sr_type = TELESCOPE_DEFAULT; 393 has_dm = PETSC_FALSE; 394 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 395 if (dm) { has_dm = PETSC_TRUE; } 396 if (has_dm) { 397 /* check for dmda */ 398 ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 399 if (same) { 400 ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 401 sr_type = TELESCOPE_DMDA; 402 } 403 /* check for dmplex */ 404 ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 405 if (same) { 406 PetscInfo(pc,"PCTelescope: found DMPLEX\n"); 407 sr_type = TELESCOPE_DMPLEX; 408 } 409 } 410 411 if (sred->ignore_dm) { 412 PetscInfo(pc,"PCTelescope: ignore DM\n"); 413 sr_type = TELESCOPE_DEFAULT; 414 } 415 sred->sr_type = sr_type; 416 } else { 417 sr_type = sred->sr_type; 418 } 419 420 /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */ 421 switch (sr_type) { 422 case TELESCOPE_DEFAULT: 423 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 424 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 425 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 426 sred->pctelescope_reset_type = NULL; 427 break; 428 case TELESCOPE_DMDA: 429 pc->ops->apply = PCApply_Telescope_dmda; 430 sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 431 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 432 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 433 sred->pctelescope_reset_type = PCReset_Telescope_dmda; 434 break; 435 case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available"); 436 break; 437 default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided"); 438 break; 439 } 440 441 /* setup */ 442 if (sred->pctelescope_setup_type) { 443 ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 444 } 445 /* update */ 446 if (!pc->setupcalled) { 447 if (sred->pctelescope_matcreate_type) { 448 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 449 } 450 if (sred->pctelescope_matnullspacecreate_type) { 451 ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 452 } 453 } else { 454 if (sred->pctelescope_matcreate_type) { 455 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 456 } 457 } 458 459 /* common - no construction */ 460 if (isActiveRank(sred->psubcomm)) { 461 ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 462 if (pc->setfromoptionscalled && !pc->setupcalled){ 463 ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 464 } 465 } 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "PCApply_Telescope" 471 static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 472 { 473 PC_Telescope sred = (PC_Telescope)pc->data; 474 PetscErrorCode ierr; 475 Vec xtmp,xred,yred; 476 PetscInt i,st,ed,bs; 477 VecScatter scatter; 478 PetscScalar *array; 479 const PetscScalar *x_array; 480 481 PetscFunctionBegin; 482 xtmp = sred->xtmp; 483 scatter = sred->scatter; 484 xred = sred->xred; 485 yred = sred->yred; 486 487 /* pull in vector x->xtmp */ 488 ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 489 ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 490 491 /* copy vector entires into xred */ 492 ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 493 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 494 if (xred) { 495 PetscScalar *LA_xred; 496 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 497 ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 498 for (i=0; i<ed-st; i++) { 499 LA_xred[i] = x_array[i]; 500 } 501 ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 502 } 503 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 504 /* solve */ 505 if (isActiveRank(sred->psubcomm)) { 506 ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 507 } 508 /* return vector */ 509 ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr); 510 ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 511 if (yred) { 512 const PetscScalar *LA_yred; 513 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 514 ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 515 for (i=0; i<ed-st; i++) { 516 array[i] = LA_yred[i]; 517 } 518 ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 519 } else { 520 for (i=0; i<bs; i++) { 521 array[i] = 0.0; 522 } 523 } 524 ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 525 ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 526 ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "PCReset_Telescope" 532 static PetscErrorCode PCReset_Telescope(PC pc) 533 { 534 PC_Telescope sred = (PC_Telescope)pc->data; 535 PetscErrorCode ierr; 536 537 ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 538 ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 539 if (sred->xred) { ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); } 540 if (sred->yred) { ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); } 541 if (sred->xtmp) { ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); } 542 if (sred->Bred) { ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); } 543 if (sred->ksp) { ierr = KSPReset(sred->ksp);CHKERRQ(ierr); } 544 if (sred->pctelescope_reset_type) { 545 ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 546 } 547 PetscFunctionReturn(0); 548 } 549 550 #undef __FUNCT__ 551 #define __FUNCT__ "PCDestroy_Telescope" 552 static PetscErrorCode PCDestroy_Telescope(PC pc) 553 { 554 PC_Telescope sred = (PC_Telescope)pc->data; 555 PetscErrorCode ierr; 556 557 PetscFunctionBegin; 558 ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 559 if (sred->ksp) { ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); } 560 ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 561 if (sred->dm_ctx) PetscFree(sred->dm_ctx); 562 PetscFree(pc->data); 563 PetscFunctionReturn(0); 564 } 565 566 #undef __FUNCT__ 567 #define __FUNCT__ "PCSetFromOptions_Telescope" 568 static PetscErrorCode PCSetFromOptions_Telescope(PetscOptions *PetscOptionsObject,PC pc) 569 { 570 PC_Telescope sred = (PC_Telescope)pc->data; 571 PetscErrorCode ierr; 572 MPI_Comm comm; 573 PetscMPIInt size; 574 575 PetscFunctionBegin; 576 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 577 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 578 ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 579 ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 580 if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 581 ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 582 ierr = PetscOptionsTail();CHKERRQ(ierr); 583 PetscFunctionReturn(0); 584 } 585 586 /* PC simplementation specific API's */ 587 588 static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 589 { 590 PC_Telescope red = (PC_Telescope)pc->data; 591 PetscFunctionBegin; 592 if (ksp) *ksp = red->ksp; 593 PetscFunctionReturn(0); 594 } 595 596 static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 597 { 598 PC_Telescope red = (PC_Telescope)pc->data; 599 PetscFunctionBegin; 600 if (fact) *fact = red->redfactor; 601 PetscFunctionReturn(0); 602 } 603 604 static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 605 { 606 PC_Telescope red = (PC_Telescope)pc->data; 607 PetscMPIInt size; 608 PetscErrorCode ierr; 609 610 PetscFunctionBegin; 611 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 612 if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 613 if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 614 red->redfactor = fact; 615 PetscFunctionReturn(0); 616 } 617 618 static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 619 { 620 PC_Telescope red = (PC_Telescope)pc->data; 621 PetscFunctionBegin; 622 if (v) *v = red->ignore_dm; 623 PetscFunctionReturn(0); 624 } 625 static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 626 { 627 PC_Telescope red = (PC_Telescope)pc->data; 628 PetscFunctionBegin; 629 red->ignore_dm = v; 630 PetscFunctionReturn(0); 631 } 632 633 static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 634 { 635 PC_Telescope red = (PC_Telescope)pc->data; 636 PetscFunctionBegin; 637 *dm = private_PCTelescopeGetSubDM(red); 638 PetscFunctionReturn(0); 639 } 640 641 /*@ 642 PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 643 644 Not Collective 645 646 Input Parameter: 647 . pc - the preconditioner context 648 649 Output Parameter: 650 . subksp - the KSP defined the smaller set of processes 651 652 Level: advanced 653 654 .keywords: PC, telescopting solve 655 @*/ 656 PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 657 { 658 PetscErrorCode ierr; 659 PetscFunctionBegin; 660 ierr = PetscTryMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 661 PetscFunctionReturn(0); 662 } 663 664 /*@ 665 PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 666 667 Not Collective 668 669 Input Parameter: 670 . pc - the preconditioner context 671 672 Output Parameter: 673 . fact - the reduction factor 674 675 Level: advanced 676 677 .keywords: PC, telescoping solve 678 @*/ 679 PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 680 { 681 PetscErrorCode ierr; 682 PetscFunctionBegin; 683 ierr = PetscTryMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 684 PetscFunctionReturn(0); 685 } 686 687 /*@ 688 PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 689 690 Not Collective 691 692 Input Parameter: 693 . pc - the preconditioner context 694 695 Output Parameter: 696 . fact - the reduction factor 697 698 Level: advanced 699 700 .keywords: PC, telescoping solve 701 @*/ 702 PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 703 { 704 PetscErrorCode ierr; 705 PetscFunctionBegin; 706 ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 707 PetscFunctionReturn(0); 708 } 709 710 /*@ 711 PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 712 713 Not Collective 714 715 Input Parameter: 716 . pc - the preconditioner context 717 718 Output Parameter: 719 . v - the flag 720 721 Level: advanced 722 723 .keywords: PC, telescoping solve 724 @*/ 725 PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 726 { 727 PetscErrorCode ierr; 728 PetscFunctionBegin; 729 ierr = PetscTryMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 730 PetscFunctionReturn(0); 731 } 732 733 /*@ 734 PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 735 736 Not Collective 737 738 Input Parameter: 739 . pc - the preconditioner context 740 741 Output Parameter: 742 . v - Use PETSC_TRUE to ignore any DM 743 744 Level: advanced 745 746 .keywords: PC, telescoping solve 747 @*/ 748 PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 749 { 750 PetscErrorCode ierr; 751 PetscFunctionBegin; 752 ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 753 PetscFunctionReturn(0); 754 } 755 756 /*@ 757 PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 758 759 Not Collective 760 761 Input Parameter: 762 . pc - the preconditioner context 763 764 Output Parameter: 765 . subdm - The re-partitioned DM 766 767 Level: advanced 768 769 .keywords: PC, telescoping solve 770 @*/ 771 PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 772 { 773 PetscErrorCode ierr; 774 PetscFunctionBegin; 775 ierr = PetscTryMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 776 PetscFunctionReturn(0); 777 } 778 779 /* -------------------------------------------------------------------------------------*/ 780 /*MC 781 PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 782 783 Options Database: 784 + -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and 785 use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes 786 - -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored 787 788 Level: advanced 789 790 Notes: 791 The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 792 Both the B mat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 793 Currently only support for re-partitioning a DMDA is provided. 794 Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator. 795 KSPSetComputeOperators() is not propagated to the sub KSP. 796 Currently there is no support for the flag -pc_use_amat 797 798 Optimizations: 799 (i) VecPlaceArray() could be used for scatters between the vectors defined on different sized communicators, thereby slightly reducing memory footprint; 800 (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). 801 802 Contributed by Dave May 803 804 .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), 805 PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), 806 PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 807 M*/ 808 #undef __FUNCT__ 809 #define __FUNCT__ "PCCreate_Telescope" 810 PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 811 { 812 PetscErrorCode ierr; 813 struct _PC_Telescope *sred; 814 815 PetscFunctionBegin; 816 ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 817 sred->redfactor = 1; 818 sred->ignore_dm = PETSC_FALSE; 819 pc->data = (void*)sred; 820 821 pc->ops->apply = PCApply_Telescope; 822 pc->ops->applytranspose = NULL; 823 pc->ops->setup = PCSetUp_Telescope; 824 pc->ops->destroy = PCDestroy_Telescope; 825 pc->ops->reset = PCReset_Telescope; 826 pc->ops->setfromoptions = PCSetFromOptions_Telescope; 827 pc->ops->view = PCView_Telescope; 828 829 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 830 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 831 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 832 sred->pctelescope_reset_type = NULL; 833 834 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 835 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 836 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 837 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 838 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 839 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842