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