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