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