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