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