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