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