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