1 2 #include <petsc/private/petscimpl.h> 3 #include <petsc/private/matimpl.h> 4 #include <petsc/private/pcimpl.h> 5 #include <petscksp.h> /*I "petscksp.h" I*/ 6 #include <petscdm.h> /*I "petscdm.h" I*/ 7 #include "../src/ksp/pc/impls/telescope/telescope.h" 8 9 static PetscBool cited = PETSC_FALSE; 10 static const char citation[] = 11 "@inproceedings{MaySananRuppKnepleySmith2016,\n" 12 " title = {Extreme-Scale Multigrid Components within PETSc},\n" 13 " author = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n" 14 " booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n" 15 " series = {PASC '16},\n" 16 " isbn = {978-1-4503-4126-4},\n" 17 " location = {Lausanne, Switzerland},\n" 18 " pages = {5:1--5:12},\n" 19 " articleno = {5},\n" 20 " numpages = {12},\n" 21 " url = {http://doi.acm.org/10.1145/2929908.2929913},\n" 22 " doi = {10.1145/2929908.2929913},\n" 23 " acmid = {2929913},\n" 24 " publisher = {ACM},\n" 25 " address = {New York, NY, USA},\n" 26 " keywords = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n" 27 " year = {2016}\n" 28 "}\n"; 29 30 /* 31 PCTelescopeSetUp_default() 32 PCTelescopeMatCreate_default() 33 34 default 35 36 // scatter in 37 x(comm) -> xtmp(comm) 38 39 xred(subcomm) <- xtmp 40 yred(subcomm) 41 42 yred(subcomm) --> xtmp 43 44 // scatter out 45 xtmp(comm) -> y(comm) 46 */ 47 48 PetscBool PetscSubComm_isActiveRank(PetscSubcomm scomm) 49 { 50 if (scomm->color == 0) { return PETSC_TRUE; } 51 else { return PETSC_FALSE; } 52 } 53 54 PetscBool isActiveRank(PC_Telescope sred) 55 { 56 if (sred->psubcomm) { 57 return(PetscSubComm_isActiveRank(sred->psubcomm)); 58 } else { 59 if (sred->subcomm != MPI_COMM_NULL) { return PETSC_TRUE; } 60 else { return PETSC_FALSE; } 61 } 62 } 63 64 /* 65 Collective on MPI_Comm[comm_f] 66 Notes 67 * Using comm_f = MPI_COMM_NULL will result in an error 68 * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid. 69 * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid. 70 */ 71 PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool *isvalid) 72 { 73 int valid = 1; 74 MPI_Group group_f,group_c; 75 PetscErrorCode ierr; 76 int errorcode; 77 PetscMPIInt count,k,size_f = 0,size_c = 0,size_c_sum = 0; 78 int *ranks_f = NULL,*ranks_c = NULL; 79 80 if (comm_f == MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"comm_f cannot be MPI_COMM_NULL"); 81 82 ierr = MPI_Comm_group(comm_f,&group_f);CHKERRQ(ierr); 83 if (comm_c != MPI_COMM_NULL) { 84 ierr = MPI_Comm_group(comm_c,&group_c);CHKERRQ(ierr); 85 } 86 87 ierr = MPI_Comm_size(comm_f,&size_f);CHKERRQ(ierr); 88 if (comm_c != MPI_COMM_NULL) { 89 ierr = MPI_Comm_size(comm_c,&size_c);CHKERRQ(ierr); 90 } 91 92 /* check not all comm_c's are NULL */ 93 size_c_sum = size_c; 94 ierr = MPI_Allreduce(MPI_IN_PLACE,&size_c_sum,1,MPI_INT,MPI_SUM,comm_f);CHKERRQ(ierr); 95 if (size_c_sum == 0) { 96 valid = 0; 97 } 98 99 /* check we can map at least 1 rank in comm_c to comm_f */ 100 ierr = PetscMalloc1(size_f,&ranks_f);CHKERRQ(ierr); 101 ierr = PetscMalloc1(size_c,&ranks_c);CHKERRQ(ierr); 102 for (k=0; k<size_f; k++) { 103 ranks_f[k] = MPI_UNDEFINED; 104 } 105 for (k=0; k<size_c; k++) { 106 ranks_c[k] = (int)k; 107 } 108 109 count = 0; 110 if (comm_c != MPI_COMM_NULL) { 111 errorcode = MPI_Group_translate_ranks(group_c,size_c,ranks_c,group_f,ranks_f); 112 for (k=0; k<size_f; k++) { 113 if (ranks_f[k] == MPI_UNDEFINED) { 114 count++; 115 } 116 } 117 } 118 if (count == size_f) { 119 valid = 0; 120 } 121 122 ierr = MPI_Allreduce(MPI_IN_PLACE,&valid,1,MPI_INT,MPI_MIN,comm_f);CHKERRQ(ierr); 123 if (valid == 1) { *isvalid = PETSC_TRUE; } 124 else { *isvalid = PETSC_FALSE; } 125 126 ierr = PetscFree(ranks_f);CHKERRQ(ierr); 127 ierr = PetscFree(ranks_c);CHKERRQ(ierr); 128 ierr = MPI_Group_free(&group_f);CHKERRQ(ierr); 129 if (comm_c != MPI_COMM_NULL) { 130 ierr = MPI_Group_free(&group_c);CHKERRQ(ierr); 131 } 132 PetscFunctionReturn(0); 133 } 134 135 DM private_PCTelescopeGetSubDM(PC_Telescope sred) 136 { 137 DM subdm = NULL; 138 139 if (!isActiveRank(sred)) { subdm = NULL; } 140 else { 141 switch (sred->sr_type) { 142 case TELESCOPE_DEFAULT: subdm = NULL; 143 break; 144 case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart; 145 break; 146 case TELESCOPE_DMPLEX: subdm = NULL; 147 break; 148 case TELESCOPE_COARSEDM: if (sred->ksp) { KSPGetDM(sred->ksp,&subdm); } 149 break; 150 } 151 } 152 return(subdm); 153 } 154 155 PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred) 156 { 157 PetscErrorCode ierr; 158 PetscInt m,M,bs,st,ed; 159 Vec x,xred,yred,xtmp; 160 Mat B; 161 MPI_Comm comm,subcomm; 162 VecScatter scatter; 163 IS isin; 164 165 PetscFunctionBegin; 166 ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr); 167 comm = PetscSubcommParent(sred->psubcomm); 168 subcomm = PetscSubcommChild(sred->psubcomm); 169 170 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 171 ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 172 ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr); 173 ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr); 174 175 xred = NULL; 176 m = 0; 177 if (isActiveRank(sred)) { 178 ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr); 179 ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr); 180 ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr); 181 ierr = VecSetFromOptions(xred);CHKERRQ(ierr); 182 ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr); 183 } 184 185 yred = NULL; 186 if (isActiveRank(sred)) { 187 ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr); 188 } 189 190 ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr); 191 ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr); 192 ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr); 193 ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr); 194 195 if (isActiveRank(sred)) { 196 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 197 ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr); 198 } else { 199 ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr); 200 ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr); 201 } 202 ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr); 203 204 ierr = VecScatterCreateWithData(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr); 205 206 sred->isin = isin; 207 sred->scatter = scatter; 208 sred->xred = xred; 209 sred->yred = yred; 210 sred->xtmp = xtmp; 211 ierr = VecDestroy(&x);CHKERRQ(ierr); 212 PetscFunctionReturn(0); 213 } 214 215 PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A) 216 { 217 PetscErrorCode ierr; 218 MPI_Comm comm,subcomm; 219 Mat Bred,B; 220 PetscInt nr,nc; 221 IS isrow,iscol; 222 Mat Blocal,*_Blocal; 223 224 PetscFunctionBegin; 225 ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr); 226 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 227 subcomm = PetscSubcommChild(sred->psubcomm); 228 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 229 ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr); 230 isrow = sred->isin; 231 ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr); 232 ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr); 233 Blocal = *_Blocal; 234 ierr = PetscFree(_Blocal);CHKERRQ(ierr); 235 Bred = NULL; 236 if (isActiveRank(sred)) { 237 PetscInt mm; 238 239 if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; } 240 241 ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr); 242 ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr); 243 } 244 *A = Bred; 245 ierr = ISDestroy(&iscol);CHKERRQ(ierr); 246 ierr = MatDestroy(&Blocal);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace) 251 { 252 PetscErrorCode ierr; 253 PetscBool has_const; 254 const Vec *vecs; 255 Vec *sub_vecs = NULL; 256 PetscInt i,k,n = 0; 257 MPI_Comm subcomm; 258 259 PetscFunctionBegin; 260 subcomm = PetscSubcommChild(sred->psubcomm); 261 ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr); 262 263 if (isActiveRank(sred)) { 264 if (n) { 265 ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr); 266 } 267 } 268 269 /* copy entries */ 270 for (k=0; k<n; k++) { 271 const PetscScalar *x_array; 272 PetscScalar *LA_sub_vec; 273 PetscInt st,ed; 274 275 /* pull in vector x->xtmp */ 276 ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 277 ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 278 if (sub_vecs) { 279 /* copy vector entries into xred */ 280 ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 281 if (sub_vecs[k]) { 282 ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr); 283 ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 284 for (i=0; i<ed-st; i++) { 285 LA_sub_vec[i] = x_array[i]; 286 } 287 ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr); 288 } 289 ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr); 290 } 291 } 292 293 if (isActiveRank(sred)) { 294 /* create new (near) nullspace for redundant object */ 295 ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr); 296 ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr); 297 if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope"); 298 if (nullspace->rmctx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope"); 299 } 300 PetscFunctionReturn(0); 301 } 302 303 static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat) 304 { 305 PetscErrorCode ierr; 306 Mat B; 307 308 PetscFunctionBegin; 309 ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr); 310 311 /* Propagate the nullspace if it exists */ 312 { 313 MatNullSpace nullspace,sub_nullspace; 314 ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr); 315 if (nullspace) { 316 ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr); 317 ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr); 318 if (isActiveRank(sred)) { 319 ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr); 320 ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr); 321 } 322 } 323 } 324 325 /* Propagate the near nullspace if it exists */ 326 { 327 MatNullSpace nearnullspace,sub_nearnullspace; 328 ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr); 329 if (nearnullspace) { 330 ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr); 331 ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr); 332 if (isActiveRank(sred)) { 333 ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr); 334 ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr); 335 } 336 } 337 } 338 PetscFunctionReturn(0); 339 } 340 341 static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer) 342 { 343 PC_Telescope sred = (PC_Telescope)pc->data; 344 PetscErrorCode ierr; 345 PetscBool iascii,isstring; 346 PetscViewer subviewer; 347 348 PetscFunctionBegin; 349 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 350 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 351 if (iascii) { 352 { 353 MPI_Comm comm,subcomm; 354 PetscMPIInt comm_size,subcomm_size; 355 DM dm = NULL,subdm = NULL; 356 357 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 358 subdm = private_PCTelescopeGetSubDM(sred); 359 360 if (sred->psubcomm) { 361 comm = PetscSubcommParent(sred->psubcomm); 362 subcomm = PetscSubcommChild(sred->psubcomm); 363 ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr); 364 ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr); 365 366 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 367 ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr); 368 ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr); 369 switch (sred->subcommtype) { 370 case PETSC_SUBCOMM_INTERLACED : 371 ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: type = interlaced\n",sred->subcommtype);CHKERRQ(ierr); 372 break; 373 case PETSC_SUBCOMM_CONTIGUOUS : 374 ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm type = contiguous\n",sred->subcommtype);CHKERRQ(ierr); 375 break; 376 default : 377 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope"); 378 } 379 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 380 } else { 381 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 382 subcomm = sred->subcomm; 383 if (!isActiveRank(sred)) { 384 subcomm = PETSC_COMM_SELF; 385 } 386 387 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 388 ierr = PetscViewerASCIIPrintf(viewer,"subcomm: using user provided sub-communicator\n");CHKERRQ(ierr); 389 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 390 } 391 392 ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 393 if (isActiveRank(sred)) { 394 ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 395 396 if (dm && sred->ignore_dm) { 397 ierr = PetscViewerASCIIPrintf(subviewer,"ignoring DM\n");CHKERRQ(ierr); 398 } 399 if (sred->ignore_kspcomputeoperators) { 400 ierr = PetscViewerASCIIPrintf(subviewer,"ignoring KSPComputeOperators\n");CHKERRQ(ierr); 401 } 402 switch (sred->sr_type) { 403 case TELESCOPE_DEFAULT: 404 ierr = PetscViewerASCIIPrintf(subviewer,"setup type: default\n");CHKERRQ(ierr); 405 break; 406 case TELESCOPE_DMDA: 407 ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMDA auto-repartitioning\n");CHKERRQ(ierr); 408 ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr); 409 break; 410 case TELESCOPE_DMPLEX: 411 ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMPLEX auto-repartitioning\n");CHKERRQ(ierr); 412 break; 413 case TELESCOPE_COARSEDM: 414 ierr = PetscViewerASCIIPrintf(subviewer,"setup type: coarse DM\n");CHKERRQ(ierr); 415 break; 416 } 417 418 if (dm) { 419 PetscObject obj = (PetscObject)dm; 420 ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object:");CHKERRQ(ierr); 421 PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 422 if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 423 if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 424 if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 425 ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr); 426 PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 427 } else { 428 ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object: NULL\n");CHKERRQ(ierr); 429 } 430 if (subdm) { 431 PetscObject obj = (PetscObject)subdm; 432 ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object:");CHKERRQ(ierr); 433 PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE); 434 if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); } 435 if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); } 436 if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); } 437 ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr); 438 PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE); 439 } else { 440 ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object: NULL\n");CHKERRQ(ierr); 441 } 442 443 ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr); 444 ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 445 } 446 ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr); 447 } 448 } 449 PetscFunctionReturn(0); 450 } 451 452 static PetscErrorCode PCSetUp_Telescope(PC pc) 453 { 454 PC_Telescope sred = (PC_Telescope)pc->data; 455 PetscErrorCode ierr; 456 MPI_Comm comm,subcomm=0; 457 PCTelescopeType sr_type; 458 459 PetscFunctionBegin; 460 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 461 462 /* Determine type of setup/update */ 463 if (!pc->setupcalled) { 464 PetscBool has_dm,same; 465 DM dm; 466 467 sr_type = TELESCOPE_DEFAULT; 468 has_dm = PETSC_FALSE; 469 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 470 if (dm) { has_dm = PETSC_TRUE; } 471 if (has_dm) { 472 /* check for dmda */ 473 ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr); 474 if (same) { 475 ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr); 476 sr_type = TELESCOPE_DMDA; 477 } 478 /* check for dmplex */ 479 ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr); 480 if (same) { 481 ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr); 482 sr_type = TELESCOPE_DMPLEX; 483 } 484 485 if (sred->use_coarse_dm) { 486 ierr = PetscInfo(pc,"PCTelescope: using coarse DM\n");CHKERRQ(ierr); 487 sr_type = TELESCOPE_COARSEDM; 488 } 489 490 if (sred->ignore_dm) { 491 ierr = PetscInfo(pc,"PCTelescope: ignoring DM\n");CHKERRQ(ierr); 492 sr_type = TELESCOPE_DEFAULT; 493 } 494 } 495 sred->sr_type = sr_type; 496 } else { 497 sr_type = sred->sr_type; 498 } 499 500 /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */ 501 switch (sr_type) { 502 case TELESCOPE_DEFAULT: 503 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 504 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 505 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 506 sred->pctelescope_reset_type = NULL; 507 break; 508 case TELESCOPE_DMDA: 509 pc->ops->apply = PCApply_Telescope_dmda; 510 pc->ops->applyrichardson = PCApplyRichardson_Telescope_dmda; 511 sred->pctelescope_setup_type = PCTelescopeSetUp_dmda; 512 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_dmda; 513 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda; 514 sred->pctelescope_reset_type = PCReset_Telescope_dmda; 515 break; 516 case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available"); 517 break; 518 case TELESCOPE_COARSEDM: 519 pc->ops->apply = PCApply_Telescope_CoarseDM; 520 pc->ops->applyrichardson = PCApplyRichardson_Telescope_CoarseDM; 521 sred->pctelescope_setup_type = PCTelescopeSetUp_CoarseDM; 522 sred->pctelescope_matcreate_type = NULL; 523 sred->pctelescope_matnullspacecreate_type = NULL;/*PCTelescopeMatNullSpaceCreate_CoarseDM;*/ 524 sred->pctelescope_reset_type = PCReset_Telescope_CoarseDM; 525 break; 526 default: SETERRQ(comm,PETSC_ERR_SUP,"Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM"); 527 break; 528 } 529 530 /* subcomm definition */ 531 if (!pc->setupcalled) { 532 if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) { 533 if (!sred->psubcomm) { 534 ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr); 535 ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr); 536 ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr); 537 ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 538 sred->subcomm = PetscSubcommChild(sred->psubcomm); 539 } 540 } else { /* query PC for DM, check communicators */ 541 DM dm,dm_coarse_partition = NULL; 542 MPI_Comm comm_fine,comm_coarse_partition = MPI_COMM_NULL; 543 PetscMPIInt csize_fine=0,csize_coarse_partition=0,cs[2],csg[2],cnt=0; 544 PetscBool isvalidsubcomm; 545 546 ierr = PCGetDM(pc,&dm);CHKERRQ(ierr); 547 comm_fine = PetscObjectComm((PetscObject)dm); 548 ierr = DMGetCoarseDM(dm,&dm_coarse_partition);CHKERRQ(ierr); 549 if (dm_coarse_partition) { cnt = 1; } 550 ierr = MPI_Allreduce(MPI_IN_PLACE,&cnt,1,MPI_INT,MPI_SUM,comm_fine);CHKERRQ(ierr); 551 if (cnt == 0) SETERRQ(comm_fine,PETSC_ERR_SUP,"Zero instances of a coarse DM were found"); 552 553 ierr = MPI_Comm_size(comm_fine,&csize_fine);CHKERRQ(ierr); 554 if (dm_coarse_partition) { 555 comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition); 556 ierr = MPI_Comm_size(comm_coarse_partition,&csize_coarse_partition);CHKERRQ(ierr); 557 } 558 559 cs[0] = csize_fine; 560 cs[1] = csize_coarse_partition; 561 ierr = MPI_Allreduce(cs,csg,2,MPI_INT,MPI_MAX,comm_fine);CHKERRQ(ierr); 562 if (csg[0] == csg[1]) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM uses the same size communicator as the parent DM attached to the PC"); 563 564 ierr = PCTelescopeTestValidSubcomm(comm_fine,comm_coarse_partition,&isvalidsubcomm);CHKERRQ(ierr); 565 if (!isvalidsubcomm) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM communicator is not a sub-communicator of parentDM->comm"); 566 sred->subcomm = comm_coarse_partition; 567 } 568 } 569 subcomm = sred->subcomm; 570 571 /* internal KSP */ 572 if (!pc->setupcalled) { 573 const char *prefix; 574 575 if (isActiveRank(sred)) { 576 ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr); 577 ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr); 578 ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 579 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr); 580 ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr); 581 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 582 ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr); 583 ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr); 584 } 585 } 586 587 /* setup */ 588 if (sred->pctelescope_setup_type) { 589 ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr); 590 } 591 /* update */ 592 if (!pc->setupcalled) { 593 if (sred->pctelescope_matcreate_type) { 594 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr); 595 } 596 if (sred->pctelescope_matnullspacecreate_type) { 597 ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 598 } 599 } else { 600 if (sred->pctelescope_matcreate_type) { 601 ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr); 602 } 603 } 604 605 /* common - no construction */ 606 if (isActiveRank(sred)) { 607 ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr); 608 if (pc->setfromoptionscalled && !pc->setupcalled){ 609 ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr); 610 } 611 } 612 613 #if 0 614 /* we perform this last as Bred is not available with KSPSetComputeOperators() until KSPSetUp has been called */ 615 if (!pc->setupcalled) { 616 if (isActiveRank(sred)) { 617 ierr = KSPSetUp(sred->ksp);CHKERRQ(ierr); 618 } 619 if (sred->pctelescope_matnullspacecreate_type) { 620 ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr); 621 } 622 } 623 #endif 624 625 PetscFunctionReturn(0); 626 } 627 628 static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y) 629 { 630 PC_Telescope sred = (PC_Telescope)pc->data; 631 PetscErrorCode ierr; 632 Vec xtmp,xred,yred; 633 PetscInt i,st,ed; 634 VecScatter scatter; 635 PetscScalar *array; 636 const PetscScalar *x_array; 637 638 PetscFunctionBegin; 639 ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 640 641 xtmp = sred->xtmp; 642 scatter = sred->scatter; 643 xred = sred->xred; 644 yred = sred->yred; 645 646 /* pull in vector x->xtmp */ 647 ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 648 ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 649 650 /* copy vector entries into xred */ 651 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 652 if (xred) { 653 PetscScalar *LA_xred; 654 ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr); 655 ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr); 656 for (i=0; i<ed-st; i++) { 657 LA_xred[i] = x_array[i]; 658 } 659 ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr); 660 } 661 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 662 /* solve */ 663 if (isActiveRank(sred)) { 664 ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr); 665 ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr); 666 } 667 /* return vector */ 668 ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr); 669 if (yred) { 670 const PetscScalar *LA_yred; 671 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 672 ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr); 673 for (i=0; i<ed-st; i++) { 674 array[i] = LA_yred[i]; 675 } 676 ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr); 677 } 678 ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr); 679 ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 680 ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 681 PetscFunctionReturn(0); 682 } 683 684 static PetscErrorCode PCApplyRichardson_Telescope(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 685 { 686 PC_Telescope sred = (PC_Telescope)pc->data; 687 PetscErrorCode ierr; 688 Vec xtmp,yred; 689 PetscInt i,st,ed; 690 VecScatter scatter; 691 const PetscScalar *x_array; 692 PetscBool default_init_guess_value; 693 694 PetscFunctionBegin; 695 xtmp = sred->xtmp; 696 scatter = sred->scatter; 697 yred = sred->yred; 698 699 if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1"); 700 *reason = (PCRichardsonConvergedReason)0; 701 702 if (!zeroguess) { 703 ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr); 704 /* pull in vector y->xtmp */ 705 ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 706 ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 707 708 /* copy vector entries into xred */ 709 ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr); 710 if (yred) { 711 PetscScalar *LA_yred; 712 ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr); 713 ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr); 714 for (i=0; i<ed-st; i++) { 715 LA_yred[i] = x_array[i]; 716 } 717 ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr); 718 } 719 ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr); 720 } 721 722 if (isActiveRank(sred)) { 723 ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr); 724 if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr); 725 } 726 727 ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr); 728 729 if (isActiveRank(sred)) { 730 ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr); 731 } 732 733 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 734 *outits = 1; 735 PetscFunctionReturn(0); 736 } 737 738 static PetscErrorCode PCReset_Telescope(PC pc) 739 { 740 PC_Telescope sred = (PC_Telescope)pc->data; 741 PetscErrorCode ierr; 742 743 ierr = ISDestroy(&sred->isin);CHKERRQ(ierr); 744 ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr); 745 ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); 746 ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); 747 ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); 748 ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); 749 ierr = KSPReset(sred->ksp);CHKERRQ(ierr); 750 if (sred->pctelescope_reset_type) { 751 ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr); 752 } 753 PetscFunctionReturn(0); 754 } 755 756 static PetscErrorCode PCDestroy_Telescope(PC pc) 757 { 758 PC_Telescope sred = (PC_Telescope)pc->data; 759 PetscErrorCode ierr; 760 761 PetscFunctionBegin; 762 ierr = PCReset_Telescope(pc);CHKERRQ(ierr); 763 ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); 764 ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr); 765 ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr); 766 ierr = PetscFree(pc->data);CHKERRQ(ierr); 767 PetscFunctionReturn(0); 768 } 769 770 static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc) 771 { 772 PC_Telescope sred = (PC_Telescope)pc->data; 773 PetscErrorCode ierr; 774 MPI_Comm comm; 775 PetscMPIInt size; 776 PetscBool flg; 777 PetscSubcommType subcommtype; 778 779 PetscFunctionBegin; 780 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 781 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 782 ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr); 783 ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr); 784 if (flg) { 785 ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr); 786 } 787 ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr); 788 if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size"); 789 ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr); 790 ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr); 791 ierr = PetscOptionsBool("-pc_telescope_use_coarse_dm","Define sub-communicator from the coarse DM","PCTelescopeSetUseCoarseDM",sred->use_coarse_dm,&sred->use_coarse_dm,0);CHKERRQ(ierr); 792 ierr = PetscOptionsTail();CHKERRQ(ierr); 793 PetscFunctionReturn(0); 794 } 795 796 /* PC simplementation specific API's */ 797 798 static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp) 799 { 800 PC_Telescope red = (PC_Telescope)pc->data; 801 PetscFunctionBegin; 802 if (ksp) *ksp = red->ksp; 803 PetscFunctionReturn(0); 804 } 805 806 static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype) 807 { 808 PC_Telescope red = (PC_Telescope)pc->data; 809 PetscFunctionBegin; 810 if (subcommtype) *subcommtype = red->subcommtype; 811 PetscFunctionReturn(0); 812 } 813 814 static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype) 815 { 816 PC_Telescope red = (PC_Telescope)pc->data; 817 818 PetscFunctionBegin; 819 if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"You cannot change the subcommunicator type for PCTelescope after it has been set up."); 820 red->subcommtype = subcommtype; 821 PetscFunctionReturn(0); 822 } 823 824 static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact) 825 { 826 PC_Telescope red = (PC_Telescope)pc->data; 827 PetscFunctionBegin; 828 if (fact) *fact = red->redfactor; 829 PetscFunctionReturn(0); 830 } 831 832 static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact) 833 { 834 PC_Telescope red = (PC_Telescope)pc->data; 835 PetscMPIInt size; 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 840 if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact); 841 if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact); 842 red->redfactor = fact; 843 PetscFunctionReturn(0); 844 } 845 846 static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v) 847 { 848 PC_Telescope red = (PC_Telescope)pc->data; 849 PetscFunctionBegin; 850 if (v) *v = red->ignore_dm; 851 PetscFunctionReturn(0); 852 } 853 854 static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v) 855 { 856 PC_Telescope red = (PC_Telescope)pc->data; 857 PetscFunctionBegin; 858 red->ignore_dm = v; 859 PetscFunctionReturn(0); 860 } 861 862 static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool *v) 863 { 864 PC_Telescope red = (PC_Telescope)pc->data; 865 PetscFunctionBegin; 866 if (v) *v = red->use_coarse_dm; 867 PetscFunctionReturn(0); 868 } 869 870 static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v) 871 { 872 PC_Telescope red = (PC_Telescope)pc->data; 873 PetscFunctionBegin; 874 red->use_coarse_dm = v; 875 PetscFunctionReturn(0); 876 } 877 878 static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v) 879 { 880 PC_Telescope red = (PC_Telescope)pc->data; 881 PetscFunctionBegin; 882 if (v) *v = red->ignore_kspcomputeoperators; 883 PetscFunctionReturn(0); 884 } 885 886 static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v) 887 { 888 PC_Telescope red = (PC_Telescope)pc->data; 889 PetscFunctionBegin; 890 red->ignore_kspcomputeoperators = v; 891 PetscFunctionReturn(0); 892 } 893 894 static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm) 895 { 896 PC_Telescope red = (PC_Telescope)pc->data; 897 PetscFunctionBegin; 898 *dm = private_PCTelescopeGetSubDM(red); 899 PetscFunctionReturn(0); 900 } 901 902 /*@ 903 PCTelescopeGetKSP - Gets the KSP created by the telescoping PC. 904 905 Not Collective 906 907 Input Parameter: 908 . pc - the preconditioner context 909 910 Output Parameter: 911 . subksp - the KSP defined the smaller set of processes 912 913 Level: advanced 914 915 .keywords: PC, telescopting solve 916 @*/ 917 PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp) 918 { 919 PetscErrorCode ierr; 920 PetscFunctionBegin; 921 ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr); 922 PetscFunctionReturn(0); 923 } 924 925 /*@ 926 PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by. 927 928 Not Collective 929 930 Input Parameter: 931 . pc - the preconditioner context 932 933 Output Parameter: 934 . fact - the reduction factor 935 936 Level: advanced 937 938 .keywords: PC, telescoping solve 939 @*/ 940 PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact) 941 { 942 PetscErrorCode ierr; 943 PetscFunctionBegin; 944 ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 948 /*@ 949 PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by. 950 951 Not Collective 952 953 Input Parameter: 954 . pc - the preconditioner context 955 956 Output Parameter: 957 . fact - the reduction factor 958 959 Level: advanced 960 961 .keywords: PC, telescoping solve 962 @*/ 963 PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact) 964 { 965 PetscErrorCode ierr; 966 PetscFunctionBegin; 967 ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 /*@ 972 PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used. 973 974 Not Collective 975 976 Input Parameter: 977 . pc - the preconditioner context 978 979 Output Parameter: 980 . v - the flag 981 982 Level: advanced 983 984 .keywords: PC, telescoping solve 985 @*/ 986 PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v) 987 { 988 PetscErrorCode ierr; 989 PetscFunctionBegin; 990 ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 991 PetscFunctionReturn(0); 992 } 993 994 /*@ 995 PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC. 996 997 Not Collective 998 999 Input Parameter: 1000 . pc - the preconditioner context 1001 1002 Output Parameter: 1003 . v - Use PETSC_TRUE to ignore any DM 1004 1005 Level: advanced 1006 1007 .keywords: PC, telescoping solve 1008 @*/ 1009 PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v) 1010 { 1011 PetscErrorCode ierr; 1012 PetscFunctionBegin; 1013 ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 1014 PetscFunctionReturn(0); 1015 } 1016 1017 /*@ 1018 PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used. 1019 1020 Not Collective 1021 1022 Input Parameter: 1023 . pc - the preconditioner context 1024 1025 Output Parameter: 1026 . v - the flag 1027 1028 Level: advanced 1029 1030 .keywords: PC, telescoping solve 1031 @*/ 1032 PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc,PetscBool *v) 1033 { 1034 PetscErrorCode ierr; 1035 PetscFunctionBegin; 1036 ierr = PetscUseMethod(pc,"PCTelescopeGetUseCoarseDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 1037 PetscFunctionReturn(0); 1038 } 1039 1040 /*@ 1041 PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM 1042 1043 Not Collective 1044 1045 Input Parameter: 1046 . pc - the preconditioner context 1047 1048 Output Parameter: 1049 . v - Use PETSC_TRUE to ignore any DM 1050 1051 Notes: 1052 When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope 1053 will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and 1054 -pc_telescope_subcomm_type will no longer have any meaning. 1055 It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes. 1056 An error will occur of the size of the communicator associated with the coarse DM 1057 is the same as that of the parent DM. 1058 Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent. 1059 This will be checked at the time the preconditioner is setup and an error will occur if 1060 the coarse DM does not define a sub-communicator of that used by the parent DM. 1061 1062 The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of 1063 the DM used (e.g. it supports DMSHELL, DMPLEX, etc). 1064 1065 Support is currently only provided for the case when you are using KSPSetComputeOperators() 1066 1067 The user is required to compose a function with the parent DM to facilitate the transfer of fields (Vec) between the different decompositions defined by the fine and coarse DMs. 1068 In the user code, this is achieved via 1069 .vb 1070 { 1071 DM dm_fine; 1072 PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method); 1073 } 1074 .ve 1075 The signature of the user provided field scatter method is 1076 .vb 1077 PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse); 1078 .ve 1079 The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE. 1080 SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM. 1081 1082 Optionally, the user may also compose a function with the parent DM to facilitate the transfer 1083 of state variables between the fine and coarse DMs. 1084 In the context of a finite element discretization, an example state variable might be 1085 values associated with quadrature points within each element. 1086 A user provided state scatter method is composed via 1087 .vb 1088 { 1089 DM dm_fine; 1090 PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method); 1091 } 1092 .ve 1093 The signature of the user provided state scatter method is 1094 .vb 1095 PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse); 1096 .ve 1097 SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM. 1098 The user is only required to support mode = SCATTER_FORWARD. 1099 No assumption is made about the data type of the state variables. 1100 These must be managed by the user and must be accessible from the DM. 1101 1102 Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be 1103 associated with the sub-KSP residing within PCTelescope. 1104 In general, PCTelescope assumes that the context on the fine and coarse DM used with 1105 KSPSetComputeOperators() should be "similar" in type or origin. 1106 Specifically the following rules are used to infer what context on the sub-KSP should be. 1107 1108 First the contexts from the KSP and the fine and coarse DMs are retrieved. 1109 Note that the special case of a DMSHELL context is queried. 1110 1111 .vb 1112 DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx); 1113 DMGetApplicationContext(dm_fine,&dmfine_appctx); 1114 DMShellGetContext(dm_fine,&dmfine_shellctx); 1115 1116 DMGetApplicationContext(dm_coarse,&dmcoarse_appctx); 1117 DMShellGetContext(dm_coarse,&dmcoarse_shellctx); 1118 .ve 1119 1120 The following rules are then enforced: 1121 1122 1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP: 1123 KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL); 1124 1125 2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx, 1126 check that dmcoarse_appctx is also non-NULL. If this is true, then: 1127 KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx); 1128 1129 3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx, 1130 check that dmcoarse_shellctx is also non-NULL. If this is true, then: 1131 KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx); 1132 1133 If neither of the above three tests passed, then PCTelescope cannot safely determine what 1134 context should be provided to KSPSetComputeOperators() for use with the sub-KSP. 1135 In this case, an additional mechanism is provided via a composed function which will return 1136 the actual context to be used. To use this feature you must compose the "getter" function 1137 with the coarse DM, e.g. 1138 .vb 1139 { 1140 DM dm_coarse; 1141 PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter); 1142 } 1143 .ve 1144 The signature of the user provided method is 1145 .vb 1146 PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext); 1147 .ve 1148 1149 Level: advanced 1150 1151 .keywords: PC, telescoping solve 1152 @*/ 1153 PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc,PetscBool v) 1154 { 1155 PetscErrorCode ierr; 1156 PetscFunctionBegin; 1157 ierr = PetscTryMethod(pc,"PCTelescopeSetUseCoarseDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 1158 PetscFunctionReturn(0); 1159 } 1160 1161 /*@ 1162 PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used. 1163 1164 Not Collective 1165 1166 Input Parameter: 1167 . pc - the preconditioner context 1168 1169 Output Parameter: 1170 . v - the flag 1171 1172 Level: advanced 1173 1174 .keywords: PC, telescoping solve 1175 @*/ 1176 PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v) 1177 { 1178 PetscErrorCode ierr; 1179 PetscFunctionBegin; 1180 ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr); 1181 PetscFunctionReturn(0); 1182 } 1183 1184 /*@ 1185 PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators. 1186 1187 Not Collective 1188 1189 Input Parameter: 1190 . pc - the preconditioner context 1191 1192 Output Parameter: 1193 . v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc 1194 1195 Level: advanced 1196 1197 .keywords: PC, telescoping solve 1198 @*/ 1199 PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v) 1200 { 1201 PetscErrorCode ierr; 1202 PetscFunctionBegin; 1203 ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr); 1204 PetscFunctionReturn(0); 1205 } 1206 1207 /*@ 1208 PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP. 1209 1210 Not Collective 1211 1212 Input Parameter: 1213 . pc - the preconditioner context 1214 1215 Output Parameter: 1216 . subdm - The re-partitioned DM 1217 1218 Level: advanced 1219 1220 .keywords: PC, telescoping solve 1221 @*/ 1222 PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm) 1223 { 1224 PetscErrorCode ierr; 1225 PetscFunctionBegin; 1226 ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr); 1227 PetscFunctionReturn(0); 1228 } 1229 1230 /*@ 1231 PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous) 1232 1233 Logically Collective 1234 1235 Input Parameter: 1236 + pc - the preconditioner context 1237 - subcommtype - the subcommunicator type (see PetscSubcommType) 1238 1239 Level: advanced 1240 1241 .keywords: PC, telescoping solve 1242 1243 .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE 1244 @*/ 1245 PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) 1246 { 1247 PetscErrorCode ierr; 1248 PetscFunctionBegin; 1249 ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr); 1250 PetscFunctionReturn(0); 1251 } 1252 1253 /*@ 1254 PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous) 1255 1256 Not Collective 1257 1258 Input Parameter: 1259 . pc - the preconditioner context 1260 1261 Output Parameter: 1262 . subcommtype - the subcommunicator type (see PetscSubcommType) 1263 1264 Level: advanced 1265 1266 .keywords: PC, telescoping solve 1267 1268 .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE 1269 @*/ 1270 PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) 1271 { 1272 PetscErrorCode ierr; 1273 PetscFunctionBegin; 1274 ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* -------------------------------------------------------------------------------------*/ 1279 /*MC 1280 PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve. 1281 1282 Options Database: 1283 + -pc_telescope_reduction_factor <r> - factor to use communicator size by. e.g. with 64 MPI processes and r=4, the new sub-communicator will have 64/4 = 16 ranks. 1284 - -pc_telescope_ignore_dm - flag to indicate whether an attached DM should be ignored 1285 - -pc_telescope_subcomm_type <interlaced,contiguous> - how to define the reduced communicator. see PetscSubcomm for more. 1286 1287 Level: advanced 1288 1289 Notes: 1290 The preconditioner is deemed telescopic as it only calls KSPSolve() on a single 1291 sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators. 1292 This means there will be MPI processes which will be idle during the application of this preconditioner. 1293 1294 The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator. 1295 Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM. 1296 Currently only support for re-partitioning a DMDA is provided. 1297 Any nullspace attached to the original Bmat operator is extracted, re-partitioned and set on the repartitioned Bmat operator. 1298 KSPSetComputeOperators() is not propagated to the sub KSP. 1299 Currently there is no support for the flag -pc_use_amat 1300 1301 Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation 1302 creates a child sub-communicator (c') containing fewer MPI processes than the original parent preconditioner (PC). 1303 1304 Developer Notes: 1305 During PCSetup, the B operator is scattered onto c'. 1306 Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c'). 1307 Then, KSPSolve() is executed on the c' communicator. 1308 1309 The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED 1310 creation routine by default (this can be changed with -pc_telescope_subcomm_type). We run the sub KSP on only the ranks within the communicator which have a color equal to zero. 1311 1312 The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator. 1313 In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into 1314 a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c') 1315 1316 The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D - 1317 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 1318 is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support 1319 for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs. 1320 1321 By default, B' is defined by simply fusing rows from different MPI processes 1322 1323 When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B 1324 into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the 1325 locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat() 1326 1327 Limitations/improvements include the following. 1328 VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage. 1329 1330 The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P, 1331 and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). We opted to use P^T.A.P as it appears 1332 VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a 1333 consistent manner. 1334 1335 Mapping of vectors is performed in the following way. 1336 Suppose the parent comm size was 4, and we set a reduction factor of 2; this would give a comm size on c' of 2. 1337 Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2. 1338 We perform the scatter to the sub-comm in the following way. 1339 [1] Given a vector x defined on comm c 1340 1341 rank(c) : _________ 0 ______ ________ 1 _______ ________ 2 _____________ ___________ 3 __________ 1342 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] 1343 1344 scatter to xtmp defined also on comm c so that we have the following values 1345 1346 rank(c) : ___________________ 0 ________________ _1_ ______________________ 2 _______________________ __3_ 1347 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] [ ] 1348 1349 The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values 1350 1351 1352 [2] Copy the value from rank 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'. 1353 Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0. 1354 1355 rank(c') : ___________________ 0 _______________ ______________________ 1 _____________________ 1356 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] 1357 1358 1359 Contributed by Dave May 1360 1361 Reference: 1362 Dave A. May, Patrick Sanan, Karl Rupp, Matthew G. Knepley, and Barry F. Smith, "Extreme-Scale Multigrid Components within PETSc". 2016. In Proceedings of the Platform for Advanced Scientific Computing Conference (PASC '16). DOI: 10.1145/2929908.2929913 1363 1364 .seealso: PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT 1365 M*/ 1366 PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) 1367 { 1368 PetscErrorCode ierr; 1369 struct _PC_Telescope *sred; 1370 1371 PetscFunctionBegin; 1372 ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr); 1373 sred->psubcomm = NULL; 1374 sred->subcommtype = PETSC_SUBCOMM_INTERLACED; 1375 sred->subcomm = MPI_COMM_NULL; 1376 sred->redfactor = 1; 1377 sred->ignore_dm = PETSC_FALSE; 1378 sred->ignore_kspcomputeoperators = PETSC_FALSE; 1379 sred->use_coarse_dm = PETSC_FALSE; 1380 pc->data = (void*)sred; 1381 1382 pc->ops->apply = PCApply_Telescope; 1383 pc->ops->applytranspose = NULL; 1384 pc->ops->applyrichardson = PCApplyRichardson_Telescope; 1385 pc->ops->setup = PCSetUp_Telescope; 1386 pc->ops->destroy = PCDestroy_Telescope; 1387 pc->ops->reset = PCReset_Telescope; 1388 pc->ops->setfromoptions = PCSetFromOptions_Telescope; 1389 pc->ops->view = PCView_Telescope; 1390 1391 sred->pctelescope_setup_type = PCTelescopeSetUp_default; 1392 sred->pctelescope_matcreate_type = PCTelescopeMatCreate_default; 1393 sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default; 1394 sred->pctelescope_reset_type = NULL; 1395 1396 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr); 1397 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr); 1398 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr); 1399 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr); 1400 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr); 1401 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr); 1402 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr); 1403 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 1404 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr); 1405 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr); 1406 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetUseCoarseDM_C",PCTelescopeGetUseCoarseDM_Telescope);CHKERRQ(ierr); 1407 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetUseCoarseDM_C",PCTelescopeSetUseCoarseDM_Telescope);CHKERRQ(ierr); 1408 PetscFunctionReturn(0); 1409 } 1410