Actual source code: viewreg.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <petsc/private/viewerimpl.h>
  3:  #include <petsc/private/hashtable.h>
  4: #if defined(PETSC_HAVE_SAWS)
  5:  #include <petscviewersaws.h>
  6: #endif

  8: PetscFunctionList PetscViewerList = 0;


 11: PetscOptionsHelpPrinted PetscOptionsHelpPrintedSingleton = NULL;
 12: KHASH_SET_INIT_STR(HTPrinted)
 13: struct  _n_PetscOptionsHelpPrinted{
 14:   khash_t(HTPrinted) *printed;
 15:   PetscSegBuffer     strings;
 16: };

 18: PetscErrorCode PetscOptionsHelpPrintedDestroy(PetscOptionsHelpPrinted *hp)
 19: {

 23:   if (!*hp) return(0);
 24:   kh_destroy(HTPrinted,(*hp)->printed);
 25:   PetscSegBufferDestroy(&(*hp)->strings);
 26:   PetscFree(*hp);
 27:   return(0);
 28: }

 30: /*@C
 31:       PetscOptionsHelpPrintedCreate - Creates an object used to manage tracking which help messages have
 32:          been printed so they will not be printed again.

 34:      Not collective

 36:     Level: developer

 38: .seealso: PetscOptionsHelpPrintedCheck(), PetscOptionsHelpPrintChecked()
 39: @*/
 40: PetscErrorCode PetscOptionsHelpPrintedCreate(PetscOptionsHelpPrinted *hp)
 41: {
 42:   PetscErrorCode             ierr;

 45:   PetscNew(hp);
 46:   (*hp)->printed = kh_init(HTPrinted);
 47:   PetscSegBufferCreate(sizeof(char),10000,&(*hp)->strings);
 48:   return(0);
 49: }

 51: /*@C
 52:       PetscOptionsHelpPrintedCheck - Checks if a particular pre, name pair has previous been entered (meaning the help message was printed)

 54:      Not collective

 56:     Input Parameters:
 57: +     hp - the object used to manage tracking what help messages have been printed
 58: .     pre - the prefix part of the string, many be NULL
 59: -     name - the string to look for (cannot be NULL)

 61:     Output Parameter:
 62: .     found - PETSC_TRUE if the string was already set

 64:     Level: intermediate


 67: .seealso: PetscOptionsHelpPrintedCreate()
 68: @*/
 69: PetscErrorCode PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrinted hp,const char *pre,const char* name,PetscBool *found)
 70: {
 71:   size_t          l1,l2;
 72: #if !defined(PETSC_HAVE_THREADSAFETY)
 73:   char            *both;
 74:   int             newitem;
 75: #endif
 76:   PetscErrorCode  ierr;

 79:   PetscStrlen(pre,&l1);
 80:   PetscStrlen(name,&l2);
 81:   if (l1+l2 == 0) {
 82:     *found = PETSC_FALSE;
 83:     return(0);
 84:   }
 85: #if !defined(PETSC_HAVE_THREADSAFETY)
 86:   PetscSegBufferGet(hp->strings,l1+l2+1,&both);
 87:   PetscStrcpy(both,pre);
 88:   PetscStrcat(both,name);
 89:   kh_put(HTPrinted,hp->printed,both,&newitem);
 90:   if (!newitem) {
 91:     PetscSegBufferUnuse(hp->strings,l1+l2+1);
 92:   }
 93:   *found = newitem ? PETSC_FALSE : PETSC_TRUE;
 94: #else
 95:   *found = PETSC_FALSE;
 96: #endif
 97:   return(0);
 98: }

100: static PetscBool noviewer = PETSC_FALSE;
101: static PetscBool noviewers[PETSCVIEWERGETVIEWEROFFPUSHESMAX];
102: static PetscInt  inoviewers = 0;

104: /*@
105:   PetscOptionsPushGetViewerOff - control whether PetscOptionsGetViewer returns a viewer.

107:   Logically Collective

109:   Input Parameter:
110: . flg - PETSC_TRUE to turn off viewer creation, PETSC_FALSE to turn it on.

112:   Level: developer

114:   Notes:
115:     Calling XXXViewFromOptions in an inner loop can be very expensive.  This can appear, for example, when using
116:    many small subsolves.  Call this function to control viewer creation in PetscOptionsGetViewer, thus removing the expensive XXXViewFromOptions calls.

118: .seealso: PetscOptionsGetViewer(), PetscOptionsPopGetViewerOff()
119: @*/
120: PetscErrorCode  PetscOptionsPushGetViewerOff(PetscBool flg)
121: {
123:   if (inoviewers > PETSCVIEWERGETVIEWEROFFPUSHESMAX - 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptionsPushGetViewerOff(), perhaps you forgot PetscOptionsPopGetViewerOff()?");

125:   noviewers[inoviewers++] = noviewer;
126:   noviewer = flg;
127:   return(0);
128: }

130: /*@
131:   PetscOptionsPopGetViewerOff - reset whether PetscOptionsGetViewer returns a viewer.

133:   Logically Collective

135:   Level: developer

137:   Notes:
138:     Calling XXXViewFromOptions in an inner loop can be very expensive.  This can appear, for example, when using
139:    many small subsolves.  Call this function to control viewer creation in PetscOptionsGetViewer, thus removing the expensive XXXViewFromOptions calls.

141: .seealso: PetscOptionsGetViewer(), PetscOptionsPushGetViewerOff()
142: @*/
143: PetscErrorCode  PetscOptionsPopGetViewerOff(void)
144: {
146:   if (!inoviewers) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many PetscOptionsPopGetViewerOff(), perhaps you forgot PetscOptionsPushGetViewerOff()?");
147:   noviewer = noviewers[--inoviewers];
148:   return(0);
149: }

151: /*@
152:   PetscOptionsGetViewerOff - does PetscOptionsGetViewer return a viewer?

154:   Logically Collective

156:   Output Parameter:
157: . flg - whether viewers are returned.

159:   Level: developer

161:   Notes:
162:     Calling XXXViewFromOptions in an inner loop can be very expensive.  This can appear, for example, when using
163:    many small subsolves.

165: .seealso: PetscOptionsGetViewer(), PetscOptionsPushGetViewerOff(), PetscOptionsPopGetViewerOff()
166: @*/
167: PetscErrorCode  PetscOptionsGetViewerOff(PetscBool *flg)
168: {
171:   *flg = noviewer;
172:   return(0);
173: }

175: /*@C
176:    PetscOptionsGetViewer - Gets a viewer appropriate for the type indicated by the user

178:    Collective on MPI_Comm

180:    Input Parameters:
181: +  comm - the communicator to own the viewer
182: .  pre - the string to prepend to the name or NULL
183: -  name - the option one is seeking

185:    Output Parameter:
186: +  viewer - the viewer, pass NULL if not needed
187: .  format - the PetscViewerFormat requested by the user, pass NULL if not needed
188: -  set - PETSC_TRUE if found, else PETSC_FALSE

190:    Level: intermediate

192:    Notes:
193:     If no value is provided ascii:stdout is used
194: $       ascii[:[filename][:[format][:append]]]    defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab,
195:                                                   for example ascii::ascii_info prints just the information about the object not all details
196:                                                   unless :append is given filename opens in write mode, overwriting what was already there
197: $       binary[:[filename][:[format][:append]]]   defaults to the file binaryoutput
198: $       draw[:drawtype[:filename]]                for example, draw:tikz, draw:tikz:figure.tex  or draw:x
199: $       socket[:port]                             defaults to the standard output port
200: $       saws[:communicatorname]                    publishes object to the Scientific Application Webserver (SAWs)

202:    Use PetscViewerDestroy() after using the viewer, otherwise a memory leak will occur

204:    You can control whether calls to this function create a viewer (or return early with *set of PETSC_FALSE) with
205:    PetscOptionsPushGetViewerOff.  This is useful if calling many small subsolves, in which case XXXViewFromOptions can take
206:    an appreciable fraction of the runtime.

208:    If PETSc is configured with --with-viewfromoptions=0 this function always returns with *set of PETSC_FALSE

210: .seealso: PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
211:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
212:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
213:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
214:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
215:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
216:           PetscOptionsFList(), PetscOptionsEList(), PetscOptionsPushGetViewerOff(), PetscOptionsPopGetViewerOff(),
217:           PetscOptionsGetViewerOff()
218: @*/
219: PetscErrorCode  PetscOptionsGetViewer(MPI_Comm comm,PetscOptions options,const char pre[],const char name[],PetscViewer *viewer,PetscViewerFormat *format,PetscBool  *set)
220: {
221:   const char                     *value;
222:   PetscErrorCode                 ierr;
223:   PetscBool                      flag,hashelp;


228:   if (viewer) *viewer = NULL;
229:   if (format) *format = PETSC_VIEWER_DEFAULT;
230:   if (set)    *set    = PETSC_FALSE;
231:   PetscOptionsGetViewerOff(&flag);
232:   if (flag) return(0);

234:   PetscOptionsHasHelp(NULL,&hashelp);
235:   if (hashelp) {
236:     PetscBool found;

238:     if (!PetscOptionsHelpPrintedSingleton) {
239:       PetscOptionsHelpPrintedCreate(&PetscOptionsHelpPrintedSingleton);
240:     }
241:     PetscOptionsHelpPrintedCheck(PetscOptionsHelpPrintedSingleton,pre,name,&found);
242:     if (!found) {
243:       if (viewer) {
244:         (*PetscHelpPrintf)(comm,"\n  -%s%s ascii[:[filename][:[format][:append]]]: %s (%s)\n",pre ? pre : "",name+1,"Prints object to stdout or ASCII file","PetscOptionsGetViewer");
245:         (*PetscHelpPrintf)(comm,"    -%s%s binary[:[filename][:[format][:append]]]: %s (%s)\n",pre ? pre : "",name+1,"Saves object to a binary file","PetscOptionsGetViewer");
246:         (*PetscHelpPrintf)(comm,"    -%s%s draw[:[drawtype][:filename|format]] %s (%s)\n",pre ? pre : "",name+1,"Draws object","PetscOptionsGetViewer");
247:         (*PetscHelpPrintf)(comm,"    -%s%s socket[:port]: %s (%s)\n",pre ? pre : "",name+1,"Pushes object to a Unix socket","PetscOptionsGetViewer");
248:         (*PetscHelpPrintf)(comm,"    -%s%s saws[:communicatorname]: %s (%s)\n\n",pre ? pre : "",name+1,"Publishes object to SAWs","PetscOptionsGetViewer");
249:       } else {
250:         (*PetscHelpPrintf)(comm,"    -%s%s\n",pre ? pre : "",name+1);
251:       }
252:     }
253:   }

255:   if (format) *format = PETSC_VIEWER_DEFAULT;
256:   PetscOptionsFindPair(options,pre,name,&value,&flag);
257:   if (flag) {
258:     if (set) *set = PETSC_TRUE;
259:     if (!value) {
260:       if (viewer) {
261:         PetscViewerASCIIGetStdout(comm,viewer);
262:         PetscObjectReference((PetscObject)*viewer);
263:       }
264:     } else {
265:       char       *loc0_vtype,*loc1_fname,*loc2_fmt = NULL,*loc3_fmode = NULL;
266:       PetscInt   cnt;
267:       const char *viewers[] = {PETSCVIEWERASCII,PETSCVIEWERBINARY,PETSCVIEWERDRAW,PETSCVIEWERSOCKET,PETSCVIEWERMATLAB,PETSCVIEWERSAWS,PETSCVIEWERVTK,PETSCVIEWERHDF5,PETSCVIEWERGLVIS,0};

269:       PetscStrallocpy(value,&loc0_vtype);
270:       PetscStrchr(loc0_vtype,':',&loc1_fname);
271:       if (loc1_fname) {
272:         *loc1_fname++ = 0;
273:         PetscStrchr(loc1_fname,':',&loc2_fmt);
274:       }
275:       if (loc2_fmt) {
276:         *loc2_fmt++ = 0;
277:         PetscStrchr(loc2_fmt,':',&loc3_fmode);
278:       }
279:       if (loc3_fmode) *loc3_fmode++ = 0;
280:       PetscStrendswithwhich(*loc0_vtype ? loc0_vtype : "ascii",viewers,&cnt);
281:       if (cnt > (PetscInt) sizeof(viewers)-1) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown viewer type: %s",loc0_vtype);
282:       if (viewer) {
283:         if (!loc1_fname) {
284:           switch (cnt) {
285:           case 0:
286:             PetscViewerASCIIGetStdout(comm,viewer);
287:             break;
288:           case 1:
289:             if (!(*viewer = PETSC_VIEWER_BINARY_(comm))) CHKERRQ(PETSC_ERR_PLIB);
290:             break;
291:           case 2:
292:             if (!(*viewer = PETSC_VIEWER_DRAW_(comm))) CHKERRQ(PETSC_ERR_PLIB);
293:             break;
294: #if defined(PETSC_USE_SOCKET_VIEWER)
295:           case 3:
296:             if (!(*viewer = PETSC_VIEWER_SOCKET_(comm))) CHKERRQ(PETSC_ERR_PLIB);
297:             break;
298: #endif
299: #if defined(PETSC_HAVE_MATLAB_ENGINE)
300:           case 4:
301:             if (!(*viewer = PETSC_VIEWER_MATLAB_(comm))) CHKERRQ(PETSC_ERR_PLIB);
302:             break;
303: #endif
304: #if defined(PETSC_HAVE_SAWS)
305:           case 5:
306:             if (!(*viewer = PETSC_VIEWER_SAWS_(comm))) CHKERRQ(PETSC_ERR_PLIB);
307:             break;
308: #endif
309: #if defined(PETSC_HAVE_HDF5)
310:           case 7:
311:             if (!(*viewer = PETSC_VIEWER_HDF5_(comm))) CHKERRQ(PETSC_ERR_PLIB);
312:             break;
313: #endif
314:           case 8:
315:             if (!(*viewer = PETSC_VIEWER_GLVIS_(comm))) CHKERRQ(PETSC_ERR_PLIB);
316:             break;
317:           default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported viewer %s",loc0_vtype);
318:           }
319:           PetscObjectReference((PetscObject)*viewer);
320:         } else {
321:           if (loc2_fmt && !*loc1_fname && (cnt == 0)) { /* ASCII format without file name */
322:             PetscViewerASCIIGetStdout(comm,viewer);
323:             PetscObjectReference((PetscObject)*viewer);
324:           } else {
325:             PetscFileMode fmode;
326:             PetscViewerCreate(comm,viewer);
327:             PetscViewerSetType(*viewer,*loc0_vtype ? loc0_vtype : "ascii");
328:             fmode = FILE_MODE_WRITE;
329:             if (loc3_fmode && *loc3_fmode) { /* Has non-empty file mode ("write" or "append") */
330:               PetscEnumFind(PetscFileModes,loc3_fmode,(PetscEnum*)&fmode,&flag);
331:               if (!flag) SETERRQ1(comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown file mode: %s",loc3_fmode);
332:             }
333:             if (loc2_fmt) {
334:               PetscBool tk,im;
335:               PetscStrcmp(loc1_fname,"tikz",&tk);
336:               PetscStrcmp(loc1_fname,"image",&im);
337:               if (tk || im) {
338:                 PetscViewerDrawSetInfo(*viewer,NULL,loc2_fmt,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
339:                 *loc2_fmt = 0;
340:               }
341:             }
342:             PetscViewerFileSetMode(*viewer,flag?fmode:FILE_MODE_WRITE);
343:             PetscViewerFileSetName(*viewer,loc1_fname);
344:             if (*loc1_fname) {
345:               PetscViewerDrawSetDrawType(*viewer,loc1_fname);
346:             }
347:             PetscViewerSetFromOptions(*viewer);
348:           }
349:         }
350:       }
351:       if (viewer) {
352:         PetscViewerSetUp(*viewer);
353:       }
354:       if (loc2_fmt && *loc2_fmt) {
355:         PetscViewerFormat tfmt;

357:         PetscEnumFind(PetscViewerFormats,loc2_fmt,(PetscEnum*)&tfmt,&flag);
358:         if (format) *format = tfmt;
359:         if (!flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown viewer format %s",loc2_fmt);
360:       } else if (viewer && (cnt == 6) && format) { /* Get format from VTK viewer */
361:         PetscViewerGetFormat(*viewer,format);
362:       }
363:       PetscFree(loc0_vtype);
364:     }
365:   }
366:   return(0);
367: }

369: /*@
370:    PetscViewerCreate - Creates a viewing context

372:    Collective on MPI_Comm

374:    Input Parameter:
375: .  comm - MPI communicator

377:    Output Parameter:
378: .  inviewer - location to put the PetscViewer context

380:    Level: advanced

382:    Concepts: graphics^creating PetscViewer
383:    Concepts: file input/output^creating PetscViewer
384:    Concepts: sockets^creating PetscViewer

386: .seealso: PetscViewerDestroy(), PetscViewerSetType(), PetscViewerType

388: @*/
389: PetscErrorCode  PetscViewerCreate(MPI_Comm comm,PetscViewer *inviewer)
390: {
391:   PetscViewer    viewer;

395:   *inviewer = 0;
396:   PetscViewerInitializePackage();
397:   PetscHeaderCreate(viewer,PETSC_VIEWER_CLASSID,"PetscViewer","PetscViewer","Viewer",comm,PetscViewerDestroy,NULL);
398:   *inviewer    = viewer;
399:   viewer->data = 0;
400:   return(0);
401: }

403: /*@C
404:    PetscViewerSetType - Builds PetscViewer for a particular implementation.

406:    Collective on PetscViewer

408:    Input Parameter:
409: +  viewer      - the PetscViewer context
410: -  type        - for example, PETSCVIEWERASCII

412:    Options Database Command:
413: .  -viewer_type  <type> - Sets the type; use -help for a list
414:     of available methods (for instance, ascii)

416:    Level: advanced

418:    Notes:
419:    See "include/petscviewer.h" for available methods (for instance,
420:    PETSCVIEWERSOCKET)

422: .seealso: PetscViewerCreate(), PetscViewerGetType(), PetscViewerType, PetscViewerPushFormat()
423: @*/
424: PetscErrorCode  PetscViewerSetType(PetscViewer viewer,PetscViewerType type)
425: {
426:   PetscErrorCode ierr,(*r)(PetscViewer);
427:   PetscBool      match;

432:   PetscObjectTypeCompare((PetscObject)viewer,type,&match);
433:   if (match) return(0);

435:   /* cleanup any old type that may be there */
436:   if (viewer->data) {
437:     (*viewer->ops->destroy)(viewer);

439:     viewer->ops->destroy = NULL;
440:     viewer->data         = 0;
441:   }
442:   PetscMemzero(viewer->ops,sizeof(struct _PetscViewerOps));

444:    PetscFunctionListFind(PetscViewerList,type,&r);
445:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown PetscViewer type given: %s",type);

447:   PetscObjectChangeTypeName((PetscObject)viewer,type);
448:   (*r)(viewer);
449:   return(0);
450: }

452: /*@C
453:    PetscViewerRegister - Adds a viewer

455:    Not Collective

457:    Input Parameters:
458: +  name_solver - name of a new user-defined viewer
459: -  routine_create - routine to create method context

461:    Level: developer
462:    Notes:
463:    PetscViewerRegister() may be called multiple times to add several user-defined viewers.

465:    Sample usage:
466: .vb
467:    PetscViewerRegister("my_viewer_type",MyViewerCreate);
468: .ve

470:    Then, your solver can be chosen with the procedural interface via
471: $     PetscViewerSetType(viewer,"my_viewer_type")
472:    or at runtime via the option
473: $     -viewer_type my_viewer_type

475:   Concepts: registering^Viewers

477: .seealso: PetscViewerRegisterAll(), PetscViewerRegisterDestroy()
478:  @*/
479: PetscErrorCode  PetscViewerRegister(const char *sname,PetscErrorCode (*function)(PetscViewer))
480: {

484:   PetscViewerInitializePackage();
485:   PetscFunctionListAdd(&PetscViewerList,sname,function);
486:   return(0);
487: }

489: /*@C
490:    PetscViewerSetFromOptions - Sets the graphics type from the options database.
491:       Defaults to a PETSc X windows graphics.

493:    Collective on PetscViewer

495:    Input Parameter:
496: .     PetscViewer - the graphics context

498:    Level: intermediate

500:    Notes:
501:     Must be called after PetscViewerCreate() before the PetscViewer is used.

503:   Concepts: PetscViewer^setting options

505: .seealso: PetscViewerCreate(), PetscViewerSetType(), PetscViewerType

507: @*/
508: PetscErrorCode  PetscViewerSetFromOptions(PetscViewer viewer)
509: {
510:   PetscErrorCode    ierr;
511:   char              vtype[256];
512:   PetscBool         flg;


517:   if (!PetscViewerList) {
518:     PetscViewerRegisterAll();
519:   }
520:   PetscObjectOptionsBegin((PetscObject)viewer);
521:   PetscOptionsFList("-viewer_type","Type of PetscViewer","None",PetscViewerList,(char*)(((PetscObject)viewer)->type_name ? ((PetscObject)viewer)->type_name : PETSCVIEWERASCII),vtype,256,&flg);
522:   if (flg) {
523:     PetscViewerSetType(viewer,vtype);
524:   }
525:   /* type has not been set? */
526:   if (!((PetscObject)viewer)->type_name) {
527:     PetscViewerSetType(viewer,PETSCVIEWERASCII);
528:   }
529:   if (viewer->ops->setfromoptions) {
530:     (*viewer->ops->setfromoptions)(PetscOptionsObject,viewer);
531:   }

533:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
534:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)viewer);
535:   PetscViewerViewFromOptions(viewer,NULL,"-viewer_view");
536:   PetscOptionsEnd();
537:   return(0);
538: }

540: PetscErrorCode PetscViewerFlowControlStart(PetscViewer viewer,PetscInt *mcnt,PetscInt *cnt)
541: {
544:   PetscViewerBinaryGetFlowControl(viewer,mcnt);
545:   PetscViewerBinaryGetFlowControl(viewer,cnt);
546:   return(0);
547: }

549: PetscErrorCode PetscViewerFlowControlStepMaster(PetscViewer viewer,PetscInt i,PetscInt *mcnt,PetscInt cnt)
550: {
552:   MPI_Comm       comm;

555:   PetscObjectGetComm((PetscObject)viewer,&comm);
556:   if (i >= *mcnt) {
557:     *mcnt += cnt;
558:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
559:   }
560:   return(0);
561: }

563: PetscErrorCode PetscViewerFlowControlEndMaster(PetscViewer viewer,PetscInt *mcnt)
564: {
566:   MPI_Comm       comm;
568:   PetscObjectGetComm((PetscObject)viewer,&comm);
569:   *mcnt = 0;
570:   MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
571:   return(0);
572: }

574: PetscErrorCode PetscViewerFlowControlStepWorker(PetscViewer viewer,PetscMPIInt rank,PetscInt *mcnt)
575: {
577:   MPI_Comm       comm;
579:   PetscObjectGetComm((PetscObject)viewer,&comm);
580:   while (PETSC_TRUE) {
581:     if (rank < *mcnt) break;
582:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
583:   }
584:   return(0);
585: }

587: PetscErrorCode PetscViewerFlowControlEndWorker(PetscViewer viewer,PetscInt *mcnt)
588: {
590:   MPI_Comm       comm;
592:   PetscObjectGetComm((PetscObject)viewer,&comm);
593:   while (PETSC_TRUE) {
594:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
595:     if (!*mcnt) break;
596:   }
597:   return(0);
598: }