1 2 /* 3 This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner. 4 */ 5 #include <petsc/private/pcimpl.h> 6 #include <petscksp.h> /*I "petscksp.h" I*/ 7 8 typedef struct { 9 KSP ksp; 10 PC pc; /* actual preconditioner used on each processor */ 11 Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of PetscObjectComm((PetscObject)pc) */ 12 Vec xdup,ydup; /* parallel vector that congregates xsub or ysub facilitating vector scattering */ 13 Mat pmats; /* matrix and optional preconditioner matrix belong to a subcommunicator */ 14 VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */ 15 PetscBool useparallelmat; 16 PetscSubcomm psubcomm; 17 PetscInt nsubcomm; /* num of data structure PetscSubcomm */ 18 PetscBool shifttypeset; 19 MatFactorShiftType shifttype; 20 } PC_Redundant; 21 22 PetscErrorCode PCFactorSetShiftType_Redundant(PC pc,MatFactorShiftType shifttype) 23 { 24 PC_Redundant *red = (PC_Redundant*)pc->data; 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 if (red->ksp) { 29 PC pc; 30 ierr = KSPGetPC(red->ksp,&pc);CHKERRQ(ierr); 31 ierr = PCFactorSetShiftType(pc,shifttype);CHKERRQ(ierr); 32 } else { 33 red->shifttypeset = PETSC_TRUE; 34 red->shifttype = shifttype; 35 } 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer) 40 { 41 PC_Redundant *red = (PC_Redundant*)pc->data; 42 PetscErrorCode ierr; 43 PetscBool iascii,isstring; 44 PetscViewer subviewer; 45 46 PetscFunctionBegin; 47 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 48 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 49 if (iascii) { 50 if (!red->psubcomm) { 51 ierr = PetscViewerASCIIPrintf(viewer," Not yet setup\n");CHKERRQ(ierr); 52 } else { 53 ierr = PetscViewerASCIIPrintf(viewer," First (color=0) of %D PCs follows\n",red->nsubcomm);CHKERRQ(ierr); 54 ierr = PetscViewerGetSubViewer(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 55 if (!red->psubcomm->color) { /* only view first redundant pc */ 56 ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr); 57 ierr = KSPView(red->ksp,subviewer);CHKERRQ(ierr); 58 ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr); 59 } 60 ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 61 } 62 } else if (isstring) { 63 ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr); 64 } 65 PetscFunctionReturn(0); 66 } 67 68 #include <../src/mat/impls/aij/mpi/mpiaij.h> 69 static PetscErrorCode PCSetUp_Redundant(PC pc) 70 { 71 PC_Redundant *red = (PC_Redundant*)pc->data; 72 PetscErrorCode ierr; 73 PetscInt mstart,mend,mlocal,M; 74 PetscMPIInt size; 75 MPI_Comm comm,subcomm; 76 Vec x; 77 78 PetscFunctionBegin; 79 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 80 81 /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 82 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 83 if (size == 1) red->useparallelmat = PETSC_FALSE; 84 85 if (!pc->setupcalled) { 86 PetscInt mloc_sub; 87 if (!red->psubcomm) { /* create red->psubcomm, new ksp and pc over subcomm */ 88 KSP ksp; 89 ierr = PCRedundantGetKSP(pc,&ksp);CHKERRQ(ierr); 90 } 91 subcomm = PetscSubcommChild(red->psubcomm); 92 93 if (red->useparallelmat) { 94 /* grab the parallel matrix and put it into processors of a subcomminicator */ 95 ierr = MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,subcomm,MAT_INITIAL_MATRIX,&red->pmats);CHKERRQ(ierr); 96 97 ierr = MPI_Comm_size(subcomm,&size);CHKERRQ(ierr); 98 if (size > 1) { 99 PetscBool foundpack; 100 101 PetscBool ismpiaij; 102 ierr = PetscObjectTypeCompare((PetscObject)(red->pmats),MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 103 if (ismpiaij) { 104 Mat_MPIAIJ *a=(Mat_MPIAIJ*)(red->pmats)->data; 105 if (!a->Mvctx_mpi1) { /* create a->Mvctx_mpi1 to be used for MatMult() */ 106 a->Mvctx_mpi1_flg = PETSC_TRUE; 107 ierr = MatSetUpMultiply_MPIAIJ(red->pmats);CHKERRQ(ierr); 108 } 109 } 110 111 ierr = MatGetFactorAvailable(red->pmats,NULL,MAT_FACTOR_LU,&foundpack);CHKERRQ(ierr); 112 if (!foundpack) { /* reset default ksp and pc */ 113 ierr = KSPSetType(red->ksp,KSPGMRES);CHKERRQ(ierr); 114 ierr = PCSetType(red->pc,PCBJACOBI);CHKERRQ(ierr); 115 } else { 116 ierr = PCFactorSetMatSolverPackage(red->pc,NULL);CHKERRQ(ierr); 117 } 118 } 119 120 ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats);CHKERRQ(ierr); 121 122 /* get working vectors xsub and ysub */ 123 ierr = MatCreateVecs(red->pmats,&red->xsub,&red->ysub);CHKERRQ(ierr); 124 125 /* create working vectors xdup and ydup. 126 xdup concatenates all xsub's contigously to form a mpi vector over dupcomm (see PetscSubcommCreate_interlaced()) 127 ydup concatenates all ysub and has empty local arrays because ysub's arrays will be place into it. 128 Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */ 129 ierr = MatGetLocalSize(red->pmats,&mloc_sub,NULL);CHKERRQ(ierr); 130 ierr = VecCreateMPI(PetscSubcommContiguousParent(red->psubcomm),mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 131 ierr = VecCreateMPIWithArray(PetscSubcommContiguousParent(red->psubcomm),1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr); 132 133 /* create vecscatters */ 134 if (!red->scatterin) { /* efficiency of scatterin is independent from psubcomm_type! */ 135 IS is1,is2; 136 PetscInt *idx1,*idx2,i,j,k; 137 138 ierr = MatCreateVecs(pc->pmat,&x,0);CHKERRQ(ierr); 139 ierr = VecGetSize(x,&M);CHKERRQ(ierr); 140 ierr = VecGetOwnershipRange(x,&mstart,&mend);CHKERRQ(ierr); 141 mlocal = mend - mstart; 142 ierr = PetscMalloc2(red->psubcomm->n*mlocal,&idx1,red->psubcomm->n*mlocal,&idx2);CHKERRQ(ierr); 143 j = 0; 144 for (k=0; k<red->psubcomm->n; k++) { 145 for (i=mstart; i<mend; i++) { 146 idx1[j] = i; 147 idx2[j++] = i + M*k; 148 } 149 } 150 ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 151 ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 152 ierr = VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 153 ierr = ISDestroy(&is1);CHKERRQ(ierr); 154 ierr = ISDestroy(&is2);CHKERRQ(ierr); 155 156 /* Impl below is good for PETSC_SUBCOMM_INTERLACED (no inter-process communication) and PETSC_SUBCOMM_CONTIGUOUS (communication within subcomm) */ 157 ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1);CHKERRQ(ierr); 158 ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 159 ierr = VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout);CHKERRQ(ierr); 160 ierr = ISDestroy(&is1);CHKERRQ(ierr); 161 ierr = ISDestroy(&is2);CHKERRQ(ierr); 162 ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 163 ierr = VecDestroy(&x);CHKERRQ(ierr); 164 } 165 } else { /* !red->useparallelmat */ 166 ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat);CHKERRQ(ierr); 167 } 168 } else { /* pc->setupcalled */ 169 if (red->useparallelmat) { 170 MatReuse reuse; 171 /* grab the parallel matrix and put it into processors of a subcomminicator */ 172 /*--------------------------------------------------------------------------*/ 173 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 174 /* destroy old matrices */ 175 ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 176 reuse = MAT_INITIAL_MATRIX; 177 } else { 178 reuse = MAT_REUSE_MATRIX; 179 } 180 ierr = MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,PetscSubcommChild(red->psubcomm),reuse,&red->pmats);CHKERRQ(ierr); 181 ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats);CHKERRQ(ierr); 182 } else { /* !red->useparallelmat */ 183 ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat);CHKERRQ(ierr); 184 } 185 } 186 187 if (pc->setfromoptionscalled) { 188 ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 189 } 190 ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 194 static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 195 { 196 PC_Redundant *red = (PC_Redundant*)pc->data; 197 PetscErrorCode ierr; 198 PetscScalar *array; 199 200 PetscFunctionBegin; 201 if (!red->useparallelmat) { 202 ierr = KSPSolve(red->ksp,x,y);CHKERRQ(ierr); 203 PetscFunctionReturn(0); 204 } 205 206 /* scatter x to xdup */ 207 ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 208 ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 209 210 /* place xdup's local array into xsub */ 211 ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 212 ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 213 214 /* apply preconditioner on each processor */ 215 ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 216 ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 217 ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 218 219 /* place ysub's local array into ydup */ 220 ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 221 ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 222 223 /* scatter ydup to y */ 224 ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 225 ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 226 ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 227 ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 static PetscErrorCode PCApplyTranspose_Redundant(PC pc,Vec x,Vec y) 232 { 233 PC_Redundant *red = (PC_Redundant*)pc->data; 234 PetscErrorCode ierr; 235 PetscScalar *array; 236 237 PetscFunctionBegin; 238 if (!red->useparallelmat) { 239 ierr = KSPSolveTranspose(red->ksp,x,y);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 /* scatter x to xdup */ 244 ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 245 ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 246 247 /* place xdup's local array into xsub */ 248 ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 249 ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 250 251 /* apply preconditioner on each processor */ 252 ierr = KSPSolveTranspose(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 253 ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 254 ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 255 256 /* place ysub's local array into ydup */ 257 ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 258 ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 259 260 /* scatter ydup to y */ 261 ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 262 ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 263 ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 264 ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 265 PetscFunctionReturn(0); 266 } 267 268 static PetscErrorCode PCReset_Redundant(PC pc) 269 { 270 PC_Redundant *red = (PC_Redundant*)pc->data; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 if (red->useparallelmat) { 275 ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 276 ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 277 ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 278 ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 279 ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 280 ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 281 } 282 ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 283 ierr = KSPReset(red->ksp);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 static PetscErrorCode PCDestroy_Redundant(PC pc) 288 { 289 PC_Redundant *red = (PC_Redundant*)pc->data; 290 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 294 ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 295 ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 296 ierr = PetscFree(pc->data);CHKERRQ(ierr); 297 PetscFunctionReturn(0); 298 } 299 300 static PetscErrorCode PCSetFromOptions_Redundant(PetscOptionItems *PetscOptionsObject,PC pc) 301 { 302 PetscErrorCode ierr; 303 PC_Redundant *red = (PC_Redundant*)pc->data; 304 305 PetscFunctionBegin; 306 ierr = PetscOptionsHead(PetscOptionsObject,"Redundant options");CHKERRQ(ierr); 307 ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 308 ierr = PetscOptionsTail();CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 313 { 314 PC_Redundant *red = (PC_Redundant*)pc->data; 315 316 PetscFunctionBegin; 317 red->nsubcomm = nreds; 318 PetscFunctionReturn(0); 319 } 320 321 /*@ 322 PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 323 324 Logically Collective on PC 325 326 Input Parameters: 327 + pc - the preconditioner context 328 - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 329 use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 330 331 Level: advanced 332 333 .keywords: PC, redundant solve 334 @*/ 335 PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 336 { 337 PetscErrorCode ierr; 338 339 PetscFunctionBegin; 340 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 341 if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 342 ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 343 PetscFunctionReturn(0); 344 } 345 346 static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 347 { 348 PC_Redundant *red = (PC_Redundant*)pc->data; 349 PetscErrorCode ierr; 350 351 PetscFunctionBegin; 352 ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 353 ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 354 355 red->scatterin = in; 356 357 ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 358 ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 359 red->scatterout = out; 360 PetscFunctionReturn(0); 361 } 362 363 /*@ 364 PCRedundantSetScatter - Sets the scatter used to copy values into the 365 redundant local solve and the scatter to move them back into the global 366 vector. 367 368 Logically Collective on PC 369 370 Input Parameters: 371 + pc - the preconditioner context 372 . in - the scatter to move the values in 373 - out - the scatter to move them out 374 375 Level: advanced 376 377 .keywords: PC, redundant solve 378 @*/ 379 PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 380 { 381 PetscErrorCode ierr; 382 383 PetscFunctionBegin; 384 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 385 PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 386 PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 387 ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 388 PetscFunctionReturn(0); 389 } 390 391 static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 392 { 393 PetscErrorCode ierr; 394 PC_Redundant *red = (PC_Redundant*)pc->data; 395 MPI_Comm comm,subcomm; 396 const char *prefix; 397 398 PetscFunctionBegin; 399 if (!red->psubcomm) { 400 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 401 402 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 403 ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 404 ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 405 ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr); 406 407 ierr = PetscSubcommSetOptionsPrefix(red->psubcomm,prefix);CHKERRQ(ierr); 408 ierr = PetscSubcommSetFromOptions(red->psubcomm);CHKERRQ(ierr); 409 ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 410 411 /* create a new PC that processors in each subcomm have copy of */ 412 subcomm = PetscSubcommChild(red->psubcomm); 413 414 ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 415 ierr = KSPSetErrorIfNotConverged(red->ksp,pc->erroriffailure);CHKERRQ(ierr); 416 ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 417 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)red->ksp);CHKERRQ(ierr); 418 ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 419 ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 420 ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 421 if (red->shifttypeset) { 422 ierr = PCFactorSetShiftType(red->pc,red->shifttype);CHKERRQ(ierr); 423 red->shifttypeset = PETSC_FALSE; 424 } 425 ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 426 ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 427 } 428 *innerksp = red->ksp; 429 PetscFunctionReturn(0); 430 } 431 432 /*@ 433 PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC. 434 435 Not Collective 436 437 Input Parameter: 438 . pc - the preconditioner context 439 440 Output Parameter: 441 . innerksp - the KSP on the smaller set of processes 442 443 Level: advanced 444 445 .keywords: PC, redundant solve 446 @*/ 447 PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 448 { 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 453 PetscValidPointer(innerksp,2); 454 ierr = PetscUseMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 458 static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 459 { 460 PC_Redundant *red = (PC_Redundant*)pc->data; 461 462 PetscFunctionBegin; 463 if (mat) *mat = red->pmats; 464 if (pmat) *pmat = red->pmats; 465 PetscFunctionReturn(0); 466 } 467 468 /*@ 469 PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 470 471 Not Collective 472 473 Input Parameter: 474 . pc - the preconditioner context 475 476 Output Parameters: 477 + mat - the matrix 478 - pmat - the (possibly different) preconditioner matrix 479 480 Level: advanced 481 482 .keywords: PC, redundant solve 483 @*/ 484 PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 485 { 486 PetscErrorCode ierr; 487 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 490 if (mat) PetscValidPointer(mat,2); 491 if (pmat) PetscValidPointer(pmat,3); 492 ierr = PetscUseMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 493 PetscFunctionReturn(0); 494 } 495 496 /* -------------------------------------------------------------------------------------*/ 497 /*MC 498 PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 499 500 Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 501 502 Options Database: 503 . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 504 use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 505 506 Level: intermediate 507 508 Notes: The default KSP is preonly and the default PC is LU. 509 510 PCFactorSetShiftType() applied to this PC will convey they shift type into the inner PC if it is factorization based. 511 512 Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 513 514 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 515 PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 516 M*/ 517 518 PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc) 519 { 520 PetscErrorCode ierr; 521 PC_Redundant *red; 522 PetscMPIInt size; 523 524 PetscFunctionBegin; 525 ierr = PetscNewLog(pc,&red);CHKERRQ(ierr); 526 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 527 528 red->nsubcomm = size; 529 red->useparallelmat = PETSC_TRUE; 530 pc->data = (void*)red; 531 532 pc->ops->apply = PCApply_Redundant; 533 pc->ops->applytranspose = PCApplyTranspose_Redundant; 534 pc->ops->setup = PCSetUp_Redundant; 535 pc->ops->destroy = PCDestroy_Redundant; 536 pc->ops->reset = PCReset_Redundant; 537 pc->ops->setfromoptions = PCSetFromOptions_Redundant; 538 pc->ops->view = PCView_Redundant; 539 540 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 541 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 542 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 543 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 544 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Redundant);CHKERRQ(ierr); 545 PetscFunctionReturn(0); 546 } 547 548