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