Actual source code: viewreg.c


  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 = NULL;

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

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

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

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

 33:      Not collective

 35:     Level: developer

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

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

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

 53:      Not collective

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

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

 63:     Level: intermediate

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

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

 98: static PetscBool noviewer = PETSC_FALSE;
 99: static PetscBool noviewers[PETSCVIEWERGETVIEWEROFFPUSHESMAX];
100: static PetscInt  inoviewers = 0;

102: /*@
103:   PetscOptionsPushGetViewerOff - control whether PetscOptionsGetViewer returns a viewer.

105:   Logically Collective

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

110:   Level: developer

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

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

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

128: /*@
129:   PetscOptionsPopGetViewerOff - reset whether PetscOptionsGetViewer returns a viewer.

131:   Logically Collective

133:   Level: developer

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

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

149: /*@
150:   PetscOptionsGetViewerOff - does PetscOptionsGetViewer return a viewer?

152:   Logically Collective

154:   Output Parameter:
155: . flg - whether viewers are returned.

157:   Level: developer

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

163: .seealso: PetscOptionsGetViewer(), PetscOptionsPushGetViewerOff(), PetscOptionsPopGetViewerOff()
164: @*/
165: PetscErrorCode  PetscOptionsGetViewerOff(PetscBool *flg)
166: {
169:   *flg = noviewer;
170:   return(0);
171: }

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

176:    Collective

178:    Input Parameters:
179: +  comm - the communicator to own the viewer
180: .  options - options database, use NULL for default global database
181: .  pre - the string to prepend to the name or NULL
182: -  name - the option one is seeking

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

189:    Level: intermediate

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

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

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

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

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


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

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

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

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

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

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

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

373:    Collective

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

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

381:    Level: advanced

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

385: @*/
386: PetscErrorCode  PetscViewerCreate(MPI_Comm comm,PetscViewer *inviewer)
387: {
388:   PetscViewer    viewer;

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

400: /*@C
401:    PetscViewerSetType - Builds PetscViewer for a particular implementation.

403:    Collective on PetscViewer

405:    Input Parameters:
406: +  viewer      - the PetscViewer context
407: -  type        - for example, PETSCVIEWERASCII

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

413:    Level: advanced

415:    Notes:
416:    See "include/petscviewer.h" for available methods (for instance,
417:    PETSCVIEWERSOCKET)

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

429:   PetscObjectTypeCompare((PetscObject)viewer,type,&match);
430:   if (match) return(0);

432:   /* cleanup any old type that may be there */
433:   if (viewer->data) {
434:     (*viewer->ops->destroy)(viewer);

436:     viewer->ops->destroy = NULL;
437:     viewer->data         = NULL;
438:   }
439:   PetscMemzero(viewer->ops,sizeof(struct _PetscViewerOps));

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

444:   PetscObjectChangeTypeName((PetscObject)viewer,type);
445:   (*r)(viewer);
446:   return(0);
447: }

449: /*@C
450:    PetscViewerRegister - Adds a viewer

452:    Not Collective

454:    Input Parameters:
455: +  name_solver - name of a new user-defined viewer
456: -  routine_create - routine to create method context

458:    Level: developer
459:    Notes:
460:    PetscViewerRegister() may be called multiple times to add several user-defined viewers.

462:    Sample usage:
463: .vb
464:    PetscViewerRegister("my_viewer_type",MyViewerCreate);
465: .ve

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

472: .seealso: PetscViewerRegisterAll()
473:  @*/
474: PetscErrorCode  PetscViewerRegister(const char *sname,PetscErrorCode (*function)(PetscViewer))
475: {

479:   PetscViewerInitializePackage();
480:   PetscFunctionListAdd(&PetscViewerList,sname,function);
481:   return(0);
482: }

484: /*@C
485:    PetscViewerSetFromOptions - Sets the graphics type from the options database.
486:       Defaults to a PETSc X windows graphics.

488:    Collective on PetscViewer

490:    Input Parameter:
491: .     PetscViewer - the graphics context

493:    Level: intermediate

495:    Notes:
496:     Must be called after PetscViewerCreate() before the PetscViewer is used.

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

500: @*/
501: PetscErrorCode  PetscViewerSetFromOptions(PetscViewer viewer)
502: {
503:   PetscErrorCode    ierr;
504:   char              vtype[256];
505:   PetscBool         flg;


510:   if (!PetscViewerList) {
511:     PetscViewerRegisterAll();
512:   }
513:   PetscObjectOptionsBegin((PetscObject)viewer);
514:   PetscOptionsFList("-viewer_type","Type of PetscViewer","None",PetscViewerList,(char*)(((PetscObject)viewer)->type_name ? ((PetscObject)viewer)->type_name : PETSCVIEWERASCII),vtype,256,&flg);
515:   if (flg) {
516:     PetscViewerSetType(viewer,vtype);
517:   }
518:   /* type has not been set? */
519:   if (!((PetscObject)viewer)->type_name) {
520:     PetscViewerSetType(viewer,PETSCVIEWERASCII);
521:   }
522:   if (viewer->ops->setfromoptions) {
523:     (*viewer->ops->setfromoptions)(PetscOptionsObject,viewer);
524:   }

526:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
527:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)viewer);
528:   PetscViewerViewFromOptions(viewer,NULL,"-viewer_view");
529:   PetscOptionsEnd();
530:   return(0);
531: }

533: PetscErrorCode PetscViewerFlowControlStart(PetscViewer viewer,PetscInt *mcnt,PetscInt *cnt)
534: {
537:   PetscViewerBinaryGetFlowControl(viewer,mcnt);
538:   PetscViewerBinaryGetFlowControl(viewer,cnt);
539:   return(0);
540: }

542: PetscErrorCode PetscViewerFlowControlStepMain(PetscViewer viewer,PetscInt i,PetscInt *mcnt,PetscInt cnt)
543: {
545:   MPI_Comm       comm;

548:   PetscObjectGetComm((PetscObject)viewer,&comm);
549:   if (i >= *mcnt) {
550:     *mcnt += cnt;
551:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
552:   }
553:   return(0);
554: }

556: PetscErrorCode PetscViewerFlowControlEndMain(PetscViewer viewer,PetscInt *mcnt)
557: {
559:   MPI_Comm       comm;
561:   PetscObjectGetComm((PetscObject)viewer,&comm);
562:   *mcnt = 0;
563:   MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
564:   return(0);
565: }

567: PetscErrorCode PetscViewerFlowControlStepWorker(PetscViewer viewer,PetscMPIInt rank,PetscInt *mcnt)
568: {
570:   MPI_Comm       comm;
572:   PetscObjectGetComm((PetscObject)viewer,&comm);
573:   while (PETSC_TRUE) {
574:     if (rank < *mcnt) break;
575:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
576:   }
577:   return(0);
578: }

580: PetscErrorCode PetscViewerFlowControlEndWorker(PetscViewer viewer,PetscInt *mcnt)
581: {
583:   MPI_Comm       comm;
585:   PetscObjectGetComm((PetscObject)viewer,&comm);
586:   while (PETSC_TRUE) {
587:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
588:     if (!*mcnt) break;
589:   }
590:   return(0);
591: }