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