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; 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 PetscInt mlocal_sub; 71 PetscMPIInt subsize,subrank; 72 PetscInt rstart_sub,rend_sub,mloc_sub; 73 const char *prefix; 74 75 PetscFunctionBegin; 76 ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 77 ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 78 79 if (!pc->setupcalled) { 80 ierr = PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);CHKERRQ(ierr); 81 ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 82 83 /* create a new PC that processors in each subcomm have copy of */ 84 subcomm = red->psubcomm->comm; 85 ierr = PCCreate(subcomm,&red->pc);CHKERRQ(ierr); 86 ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 87 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 88 ierr = PCSetOptionsPrefix(red->pc,prefix);CHKERRQ(ierr); 89 ierr = PCAppendOptionsPrefix(red->pc,"redundant_");CHKERRQ(ierr); 90 ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr); 91 92 /* create working vectors xsub/ysub and xdup/ydup */ 93 ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 94 ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr); 95 96 /* get local size of xsub/ysub */ 97 ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 98 ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 99 rstart_sub = pc->pmat->rmap.range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */ 100 if (subrank+1 < subsize){ 101 rend_sub = pc->pmat->rmap.range[red->psubcomm->n*(subrank+1)]; 102 } else { 103 rend_sub = m; 104 } 105 mloc_sub = rend_sub - rstart_sub; 106 ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 107 /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 108 ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr); 109 110 /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 111 Note: we use communicator dupcomm, not pc->comm! */ 112 ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 113 ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr); 114 115 /* create vec scatters */ 116 if (!red->scatterin){ 117 IS is1,is2; 118 PetscInt *idx1,*idx2,i,j,k; 119 ierr = PetscMalloc(2*red->psubcomm->n*mlocal*sizeof(PetscInt),&idx1);CHKERRQ(ierr); 120 idx2 = idx1 + red->psubcomm->n*mlocal; 121 j = 0; 122 for (k=0; k<red->psubcomm->n; k++){ 123 for (i=mstart; i<mend; i++){ 124 idx1[j] = i; 125 idx2[j++] = i + m*k; 126 } 127 } 128 ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);CHKERRQ(ierr); 129 ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);CHKERRQ(ierr); 130 ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 131 ierr = ISDestroy(is1);CHKERRQ(ierr); 132 ierr = ISDestroy(is2);CHKERRQ(ierr); 133 134 ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 135 ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 136 ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 137 ierr = ISDestroy(is1);CHKERRQ(ierr); 138 ierr = ISDestroy(is2);CHKERRQ(ierr); 139 ierr = PetscFree(idx1);CHKERRQ(ierr); 140 } 141 } 142 ierr = VecDestroy(vec);CHKERRQ(ierr); 143 144 /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 145 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 146 if (size == 1) { 147 red->useparallelmat = PETSC_FALSE; 148 } 149 150 if (red->useparallelmat) { 151 if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 152 /* destroy old matrices */ 153 if (red->pmats) { 154 ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 155 } 156 } else if (pc->setupcalled == 1) { 157 reuse = MAT_REUSE_MATRIX; 158 str = SAME_NONZERO_PATTERN; 159 } 160 161 /* grab the parallel matrix and put it into processors of a subcomminicator */ 162 /*--------------------------------------------------------------------------*/ 163 ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 164 ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr); 165 166 /* tell PC of the subcommunicator its operators */ 167 ierr = PCSetOperators(red->pc,red->pmats,red->pmats,str);CHKERRQ(ierr); 168 } else { 169 ierr = PCSetOperators(red->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 170 } 171 ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr); 172 ierr = PCSetUp(red->pc);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(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr); 187 ierr = VecScatterEnd(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);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 = PCApply(red->pc,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->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr); 204 ierr = VecScatterEnd(red->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);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__ "PCDestroy_Redundant" 212 static PetscErrorCode PCDestroy_Redundant(PC pc) 213 { 214 PC_Redundant *red = (PC_Redundant*)pc->data; 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 if (red->scatterin) {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);} 219 if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);} 220 if (red->ysub) {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);} 221 if (red->xsub) {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);} 222 if (red->xdup) {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);} 223 if (red->ydup) {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);} 224 if (red->pmats) { 225 ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 226 } 227 if (red->psubcomm) {ierr = PetscSubcommDestroy(red->psubcomm);CHKERRQ(ierr);} 228 if (red->pc) {ierr = PCDestroy(red->pc);CHKERRQ(ierr);} 229 ierr = PetscFree(red);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "PCSetFromOptions_Redundant" 235 static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 236 { 237 PetscErrorCode ierr; 238 PC_Redundant *red = (PC_Redundant*)pc->data; 239 240 PetscFunctionBegin; 241 ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 242 ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 243 ierr = PetscOptionsTail();CHKERRQ(ierr); 244 PetscFunctionReturn(0); 245 } 246 247 EXTERN_C_BEGIN 248 #undef __FUNCT__ 249 #define __FUNCT__ "PCRedundantSetNumber_Redundant" 250 PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 251 { 252 PC_Redundant *red = (PC_Redundant*)pc->data; 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 red->nsubcomm = nreds; 257 PetscFunctionReturn(0); 258 } 259 EXTERN_C_END 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "PCRedundantSetNumber" 263 /*@ 264 PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 265 266 Collective on PC 267 268 Input Parameters: 269 + pc - the preconditioner context 270 - nredundant - number of redundant preconditioner contexts 271 272 Level: advanced 273 274 .keywords: PC, redundant solve 275 @*/ 276 PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC pc,PetscInt nredundant) 277 { 278 PetscErrorCode ierr,(*f)(PC,PetscInt); 279 280 PetscFunctionBegin; 281 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 282 if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 283 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);CHKERRQ(ierr); 284 if (f) { 285 ierr = (*f)(pc,nredundant);CHKERRQ(ierr); 286 } 287 PetscFunctionReturn(0); 288 } 289 290 EXTERN_C_BEGIN 291 #undef __FUNCT__ 292 #define __FUNCT__ "PCRedundantSetScatter_Redundant" 293 PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 294 { 295 PC_Redundant *red = (PC_Redundant*)pc->data; 296 PetscErrorCode ierr; 297 298 PetscFunctionBegin; 299 red->scatterin = in; 300 red->scatterout = out; 301 ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 302 ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 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_comm - number of sub communicators to use 439 440 Level: intermediate 441 442 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 443 PCRedundantGetPC(), PCRedundantGetOperators() 444 M*/ 445 446 EXTERN_C_BEGIN 447 #undef __FUNCT__ 448 #define __FUNCT__ "PCCreate_Redundant" 449 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc) 450 { 451 PetscErrorCode ierr; 452 PC_Redundant *red; 453 PetscMPIInt size; 454 455 PetscFunctionBegin; 456 ierr = PetscNew(PC_Redundant,&red);CHKERRQ(ierr); 457 ierr = PetscLogObjectMemory(pc,sizeof(PC_Redundant));CHKERRQ(ierr); 458 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 459 red->nsubcomm = size; 460 red->useparallelmat = PETSC_TRUE; 461 pc->data = (void*)red; 462 463 pc->ops->apply = PCApply_Redundant; 464 pc->ops->applytranspose = 0; 465 pc->ops->setup = PCSetUp_Redundant; 466 pc->ops->destroy = PCDestroy_Redundant; 467 pc->ops->setfromoptions = PCSetFromOptions_Redundant; 468 pc->ops->view = PCView_Redundant; 469 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant", 470 PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 471 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant", 472 PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 473 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant", 474 PCRedundantGetPC_Redundant);CHKERRQ(ierr); 475 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant", 476 PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 477 PetscFunctionReturn(0); 478 } 479 EXTERN_C_END 480