Actual source code: petscviewer.h

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: /*
  2:      PetscViewers are objects where other objects can be looked at or stored.
  3: */


  8: #include <petscsys.h>
  9: #include <petscviewertypes.h>

 11: PETSC_EXTERN PetscClassId PETSC_VIEWER_CLASSID;

 13: /*J
 14:     PetscViewerType - String with the name of a PETSc PETScViewer

 16:    Level: beginner

 18: .seealso: PetscViewerSetType(), PetscViewer, PetscViewerRegister(), PetscViewerCreate()
 19: J*/
 20: typedef const char* PetscViewerType;
 21: #define PETSCVIEWERSOCKET       "socket"
 22: #define PETSCVIEWERASCII        "ascii"
 23: #define PETSCVIEWERBINARY       "binary"
 24: #define PETSCVIEWERSTRING       "string"
 25: #define PETSCVIEWERDRAW         "draw"
 26: #define PETSCVIEWERVU           "vu"
 27: #define PETSCVIEWERMATHEMATICA  "mathematica"
 28: #define PETSCVIEWERNETCDF       "netcdf"
 29: #define PETSCVIEWERHDF5         "hdf5"
 30: #define PETSCVIEWERVTK          "vtk"
 31: #define PETSCVIEWERMATLAB       "matlab"
 32: #define PETSCVIEWERSAWS          "saws"

 34: PETSC_EXTERN PetscFunctionList PetscViewerList;
 35: PETSC_EXTERN PetscErrorCode PetscViewerInitializePackage(void);

 37: PETSC_EXTERN PetscErrorCode PetscViewerRegister(const char[],PetscErrorCode (*)(PetscViewer));

 39: PETSC_EXTERN PetscErrorCode PetscViewerCreate(MPI_Comm,PetscViewer*);
 40: PETSC_EXTERN PetscErrorCode PetscViewerSetFromOptions(PetscViewer);
 41: PETSC_EXTERN PetscErrorCode PetscViewerASCIIOpenWithFILE(MPI_Comm,FILE*,PetscViewer*);

 43: PETSC_EXTERN PetscErrorCode PetscViewerASCIIOpen(MPI_Comm,const char[],PetscViewer*);
 44: PETSC_EXTERN PetscErrorCode PetscViewerASCIISetFILE(PetscViewer,FILE*);
 45: PETSC_EXTERN PetscErrorCode PetscViewerBinaryOpen(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
 46: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetFlowControl(PetscViewer,PetscInt*);
 47: PETSC_EXTERN PetscErrorCode PetscViewerBinarySetFlowControl(PetscViewer,PetscInt);
 48: PETSC_EXTERN PetscErrorCode PetscViewerBinarySetUseMPIIO(PetscViewer,PetscBool);
 49: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetUseMPIIO(PetscViewer,PetscBool *);
 50: #if defined(PETSC_HAVE_MPIIO)
 51: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetMPIIODescriptor(PetscViewer,MPI_File*);
 52: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetMPIIOOffset(PetscViewer,MPI_Offset*);
 53: PETSC_EXTERN PetscErrorCode PetscViewerBinaryAddMPIIOOffset(PetscViewer,MPI_Offset);
 54: #endif

 56: PETSC_EXTERN PetscErrorCode PetscViewerSocketOpen(MPI_Comm,const char[],int,PetscViewer*);
 57: PETSC_EXTERN PetscErrorCode PetscViewerStringOpen(MPI_Comm,char[],size_t,PetscViewer*);
 58: PETSC_EXTERN PetscErrorCode PetscViewerDrawOpen(MPI_Comm,const char[],const char[],int,int,int,int,PetscViewer*);
 59: #include <petscdrawtypes.h>
 60: PETSC_EXTERN PetscErrorCode PetscViewerDrawSetDrawType(PetscViewer,PetscDrawType);
 61: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaOpen(MPI_Comm, int, const char[], const char[], PetscViewer *);
 62: PETSC_EXTERN PetscErrorCode PetscViewerSiloOpen(MPI_Comm, const char[], PetscViewer *);
 63: PETSC_EXTERN PetscErrorCode PetscViewerMatlabOpen(MPI_Comm,const char[],PetscFileMode,PetscViewer*);

 65: PETSC_EXTERN PetscErrorCode PetscViewerGetType(PetscViewer,PetscViewerType*);
 66: PETSC_EXTERN PetscErrorCode PetscViewerSetType(PetscViewer,PetscViewerType);
 67: PETSC_EXTERN PetscErrorCode PetscViewerDestroy(PetscViewer*);
 68: PETSC_EXTERN PetscErrorCode PetscViewerGetSingleton(PetscViewer,PetscViewer*);
 69: PETSC_EXTERN PetscErrorCode PetscViewerRestoreSingleton(PetscViewer,PetscViewer*);
 70: PETSC_EXTERN PetscErrorCode PetscViewerGetSubcomm(PetscViewer,MPI_Comm,PetscViewer*);
 71: PETSC_EXTERN PetscErrorCode PetscViewerRestoreSubcomm(PetscViewer,MPI_Comm,PetscViewer*);

 73: PETSC_EXTERN PetscErrorCode PetscViewerSetUp(PetscViewer);
 74: PETSC_EXTERN PetscErrorCode PetscViewerView(PetscViewer,PetscViewer);
 75: PETSC_STATIC_INLINE PetscErrorCode PetscViewerViewFromOptions(PetscViewer A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

 77: PETSC_EXTERN PetscErrorCode PetscViewerSetOptionsPrefix(PetscViewer,const char[]);
 78: PETSC_EXTERN PetscErrorCode PetscViewerAppendOptionsPrefix(PetscViewer,const char[]);
 79: PETSC_EXTERN PetscErrorCode PetscViewerGetOptionsPrefix(PetscViewer,const char*[]);

 81: /*E
 82:     PetscViewerFormat - Way a viewer presents the object

 84:    Level: beginner

 86:    The values below are also listed in petsc/finclude/petscviewer.h. If another values is added below it
 87:    must also be added there.

 89: .seealso: PetscViewerSetFormat(), PetscViewer, PetscViewerType, PetscViewerPushFormat(), PetscViewerPopFormat()
 90: E*/
 91: typedef enum {
 92:   PETSC_VIEWER_DEFAULT,
 93:   PETSC_VIEWER_ASCII_MATLAB,
 94:   PETSC_VIEWER_ASCII_MATHEMATICA,
 95:   PETSC_VIEWER_ASCII_IMPL,
 96:   PETSC_VIEWER_ASCII_INFO,
 97:   PETSC_VIEWER_ASCII_INFO_DETAIL,
 98:   PETSC_VIEWER_ASCII_COMMON,
 99:   PETSC_VIEWER_ASCII_SYMMODU,
100:   PETSC_VIEWER_ASCII_INDEX,
101:   PETSC_VIEWER_ASCII_DENSE,
102:   PETSC_VIEWER_ASCII_MATRIXMARKET,
103:   PETSC_VIEWER_ASCII_VTK,
104:   PETSC_VIEWER_ASCII_VTK_CELL,
105:   PETSC_VIEWER_ASCII_VTK_COORDS,
106:   PETSC_VIEWER_ASCII_PCICE,
107:   PETSC_VIEWER_ASCII_PYTHON,
108:   PETSC_VIEWER_ASCII_FACTOR_INFO,
109:   PETSC_VIEWER_ASCII_LATEX,
110:   PETSC_VIEWER_DRAW_BASIC,
111:   PETSC_VIEWER_DRAW_LG,
112:   PETSC_VIEWER_DRAW_CONTOUR,
113:   PETSC_VIEWER_DRAW_PORTS,
114:   PETSC_VIEWER_VTK_VTS,
115:   PETSC_VIEWER_VTK_VTR,
116:   PETSC_VIEWER_VTK_VTU,
117:   PETSC_VIEWER_BINARY_MATLAB,
118:   PETSC_VIEWER_NATIVE,
119:   PETSC_VIEWER_HDF5_VIZ,
120:   PETSC_VIEWER_NOFORMAT
121:   } PetscViewerFormat;
122: PETSC_EXTERN const char *const PetscViewerFormats[];

124: PETSC_EXTERN PetscErrorCode PetscViewerSetFormat(PetscViewer,PetscViewerFormat);
125: PETSC_EXTERN PetscErrorCode PetscViewerPushFormat(PetscViewer,PetscViewerFormat);
126: PETSC_EXTERN PetscErrorCode PetscViewerPopFormat(PetscViewer);
127: PETSC_EXTERN PetscErrorCode PetscViewerGetFormat(PetscViewer,PetscViewerFormat*);
128: PETSC_EXTERN PetscErrorCode PetscViewerFlush(PetscViewer);

130: PETSC_EXTERN PetscErrorCode PetscOptionsGetViewer(MPI_Comm,const char[],const char[],PetscViewer*,PetscViewerFormat*,PetscBool*);
131: #define PetscOptionsViewer(a,b,c,d,e,f) PetscOptionsViewer_Private(PetscOptionsObject,a,b,c,d,e,f);
132: PETSC_EXTERN PetscErrorCode PetscOptionsViewer_Private(PetscOptions*,const char[],const char[],const char[],PetscViewer*,PetscViewerFormat *,PetscBool *);

134: /*
135:    Operations explicit to a particular class of viewers
136: */

138: PETSC_EXTERN PetscErrorCode PetscViewerASCIIGetPointer(PetscViewer,FILE**);
139: PETSC_EXTERN PetscErrorCode PetscViewerFileGetMode(PetscViewer,PetscFileMode*);
140: PETSC_EXTERN PetscErrorCode PetscViewerFileSetMode(PetscViewer,PetscFileMode);
141: PETSC_EXTERN PetscErrorCode PetscViewerRead(PetscViewer,void*,PetscInt,PetscInt*,PetscDataType);
142: PETSC_EXTERN PetscErrorCode PetscViewerASCIIPrintf(PetscViewer,const char[],...);
143: PETSC_EXTERN PetscErrorCode PetscViewerASCIISynchronizedPrintf(PetscViewer,const char[],...);
144: PETSC_EXTERN PetscErrorCode PetscViewerASCIISynchronizedAllow(PetscViewer,PetscBool);
145: PETSC_EXTERN PetscErrorCode PetscViewerASCIIPushTab(PetscViewer);
146: PETSC_EXTERN PetscErrorCode PetscViewerASCIIPopTab(PetscViewer);
147: PETSC_EXTERN PetscErrorCode PetscViewerASCIIUseTabs(PetscViewer,PetscBool );
148: PETSC_EXTERN PetscErrorCode PetscViewerASCIISetTab(PetscViewer,PetscInt);
149: PETSC_EXTERN PetscErrorCode PetscViewerASCIIGetTab(PetscViewer,PetscInt*);
150: PETSC_EXTERN PetscErrorCode PetscViewerASCIIAddTab(PetscViewer,PetscInt);
151: PETSC_EXTERN PetscErrorCode PetscViewerASCIISubtractTab(PetscViewer,PetscInt);
152: PETSC_EXTERN PetscErrorCode PetscViewerASCIIRead(PetscViewer,void *,PetscInt,PetscInt*,PetscDataType);
153: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetDescriptor(PetscViewer,int*);
154: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetInfoPointer(PetscViewer,FILE **);
155: PETSC_EXTERN PetscErrorCode PetscViewerBinaryRead(PetscViewer,void*,PetscInt,PetscInt*,PetscDataType);
156: PETSC_EXTERN PetscErrorCode PetscViewerBinaryWrite(PetscViewer,void*,PetscInt,PetscDataType,PetscBool );
157: PETSC_EXTERN PetscErrorCode PetscViewerStringSPrintf(PetscViewer,const char[],...);
158: PETSC_EXTERN PetscErrorCode PetscViewerStringSetString(PetscViewer,char[],PetscInt);
159: PETSC_EXTERN PetscErrorCode PetscViewerDrawClear(PetscViewer);
160: PETSC_EXTERN PetscErrorCode PetscViewerDrawSetHold(PetscViewer,PetscBool);
161: PETSC_EXTERN PetscErrorCode PetscViewerDrawGetHold(PetscViewer,PetscBool*);
162: PETSC_EXTERN PetscErrorCode PetscViewerDrawSetPause(PetscViewer,PetscReal);
163: PETSC_EXTERN PetscErrorCode PetscViewerDrawGetPause(PetscViewer,PetscReal*);
164: PETSC_EXTERN PetscErrorCode PetscViewerDrawSetInfo(PetscViewer,const char[],const char[],int,int,int,int);
165: PETSC_EXTERN PetscErrorCode PetscViewerDrawResize(PetscViewer,int,int);
166: PETSC_EXTERN PetscErrorCode PetscViewerDrawSetBounds(PetscViewer,PetscInt,const PetscReal*);
167: PETSC_EXTERN PetscErrorCode PetscViewerDrawGetBounds(PetscViewer,PetscInt*,const PetscReal**);
168: PETSC_EXTERN PetscErrorCode PetscViewerSocketSetConnection(PetscViewer,const char[],int);
169: PETSC_EXTERN PetscErrorCode PetscViewerBinarySkipInfo(PetscViewer);
170: PETSC_EXTERN PetscErrorCode PetscViewerBinarySetSkipInfo(PetscViewer,PetscBool);
171: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetSkipInfo(PetscViewer,PetscBool*);
172: PETSC_EXTERN PetscErrorCode PetscViewerBinarySetSkipOptions(PetscViewer,PetscBool );
173: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetSkipOptions(PetscViewer,PetscBool *);
174: PETSC_EXTERN PetscErrorCode PetscViewerBinarySetSkipHeader(PetscViewer,PetscBool);
175: PETSC_EXTERN PetscErrorCode PetscViewerBinaryGetSkipHeader(PetscViewer,PetscBool*);
176: PETSC_EXTERN PetscErrorCode PetscViewerBinaryReadStringArray(PetscViewer,char***);
177: PETSC_EXTERN PetscErrorCode PetscViewerBinaryWriteStringArray(PetscViewer,char**);

179: PETSC_EXTERN PetscErrorCode PetscViewerFileSetName(PetscViewer,const char[]);
180: PETSC_EXTERN PetscErrorCode PetscViewerFileGetName(PetscViewer,const char**);

182: PETSC_EXTERN PetscErrorCode PetscViewerVUGetPointer(PetscViewer, FILE**);
183: PETSC_EXTERN PetscErrorCode PetscViewerVUSetVecSeen(PetscViewer, PetscBool );
184: PETSC_EXTERN PetscErrorCode PetscViewerVUGetVecSeen(PetscViewer, PetscBool  *);
185: PETSC_EXTERN PetscErrorCode PetscViewerVUPrintDeferred(PetscViewer, const char [], ...);
186: PETSC_EXTERN PetscErrorCode PetscViewerVUFlushDeferred(PetscViewer);

188: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaInitializePackage(void);
189: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaFinalizePackage(void);
190: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaGetName(PetscViewer, const char **);
191: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaSetName(PetscViewer, const char []);
192: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaClearName(PetscViewer);
193: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaSkipPackets(PetscViewer, int);

195: PETSC_EXTERN PetscErrorCode PetscViewerSiloGetName(PetscViewer, char **);
196: PETSC_EXTERN PetscErrorCode PetscViewerSiloSetName(PetscViewer, const char []);
197: PETSC_EXTERN PetscErrorCode PetscViewerSiloClearName(PetscViewer);
198: PETSC_EXTERN PetscErrorCode PetscViewerSiloGetMeshName(PetscViewer, char **);
199: PETSC_EXTERN PetscErrorCode PetscViewerSiloSetMeshName(PetscViewer, const char []);
200: PETSC_EXTERN PetscErrorCode PetscViewerSiloClearMeshName(PetscViewer);

202: PETSC_EXTERN PetscErrorCode PetscViewerNetcdfOpen(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
203: PETSC_EXTERN PetscErrorCode PetscViewerNetcdfGetID(PetscViewer, int *);

205: typedef enum {PETSC_VTK_POINT_FIELD, PETSC_VTK_POINT_VECTOR_FIELD, PETSC_VTK_CELL_FIELD, PETSC_VTK_CELL_VECTOR_FIELD} PetscViewerVTKFieldType;
206: PETSC_EXTERN PetscErrorCode PetscViewerVTKAddField(PetscViewer,PetscObject,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType,PetscObject);
207: PETSC_EXTERN PetscErrorCode PetscViewerVTKOpen(MPI_Comm,const char[],PetscFileMode,PetscViewer*);

209: /*
210:      These are all the default viewers that do not have to be explicitly opened
211: */
212: PETSC_EXTERN PetscViewer    PETSC_VIEWER_STDOUT_(MPI_Comm);
213: PETSC_EXTERN PetscErrorCode PetscViewerASCIIGetStdout(MPI_Comm,PetscViewer*);
214: PETSC_EXTERN PetscViewer    PETSC_VIEWER_STDERR_(MPI_Comm);
215: PETSC_EXTERN PetscErrorCode PetscViewerASCIIGetStderr(MPI_Comm,PetscViewer*);
216: PETSC_EXTERN PetscViewer    PETSC_VIEWER_DRAW_(MPI_Comm);
217: PETSC_EXTERN PetscViewer    PETSC_VIEWER_SOCKET_(MPI_Comm);
218: PETSC_EXTERN PetscViewer    PETSC_VIEWER_BINARY_(MPI_Comm);
219: PETSC_EXTERN PetscViewer    PETSC_VIEWER_MATLAB_(MPI_Comm);
220: PETSC_EXTERN PetscViewer    PETSC_VIEWER_HDF5_(MPI_Comm);
221: PETSC_EXTERN PetscViewer   PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE;

223: #define PETSC_VIEWER_STDERR_SELF  PETSC_VIEWER_STDERR_(PETSC_COMM_SELF)
224: #define PETSC_VIEWER_STDERR_WORLD PETSC_VIEWER_STDERR_(PETSC_COMM_WORLD)

226: /*MC
227:   PETSC_VIEWER_STDOUT_WORLD  - same as PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)

229:   Level: beginner
230: M*/
231: #define PETSC_VIEWER_STDOUT_WORLD PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)

233: /*MC
234:   PETSC_VIEWER_STDOUT_SELF  - same as PETSC_VIEWER_STDOUT_(PETSC_COMM_SELF)

236:   Level: beginner
237: M*/
238: #define PETSC_VIEWER_STDOUT_SELF  PETSC_VIEWER_STDOUT_(PETSC_COMM_SELF)

240: /*MC
241:   PETSC_VIEWER_DRAW_WORLD  - same as PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD)

243:   Level: intermediate
244: M*/
245: #define PETSC_VIEWER_DRAW_WORLD   PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD)

247: /*MC
248:   PETSC_VIEWER_DRAW_SELF  - same as PETSC_VIEWER_DRAW_(PETSC_COMM_SELF)

250:   Level: intermediate
251: M*/
252: #define PETSC_VIEWER_DRAW_SELF    PETSC_VIEWER_DRAW_(PETSC_COMM_SELF)

254: /*MC
255:   PETSC_VIEWER_SOCKET_WORLD  - same as PETSC_VIEWER_SOCKET_(PETSC_COMM_WORLD)

257:   Level: intermediate
258: M*/
259: #define PETSC_VIEWER_SOCKET_WORLD PETSC_VIEWER_SOCKET_(PETSC_COMM_WORLD)

261: /*MC
262:   PETSC_VIEWER_SOCKET_SELF  - same as PETSC_VIEWER_SOCKET_(PETSC_COMM_SELF)

264:   Level: intermediate
265: M*/
266: #define PETSC_VIEWER_SOCKET_SELF  PETSC_VIEWER_SOCKET_(PETSC_COMM_SELF)

268: /*MC
269:   PETSC_VIEWER_BINARY_WORLD  - same as PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD)

271:   Level: intermediate
272: M*/
273: #define PETSC_VIEWER_BINARY_WORLD PETSC_VIEWER_BINARY_(PETSC_COMM_WORLD)

275: /*MC
276:   PETSC_VIEWER_BINARY_SELF  - same as PETSC_VIEWER_BINARY_(PETSC_COMM_SELF)

278:   Level: intermediate
279: M*/
280: #define PETSC_VIEWER_BINARY_SELF  PETSC_VIEWER_BINARY_(PETSC_COMM_SELF)

282: /*MC
283:   PETSC_VIEWER_MATLAB_WORLD  - same as PETSC_VIEWER_MATLAB_(PETSC_COMM_WORLD)

285:   Level: intermediate
286: M*/
287: #define PETSC_VIEWER_MATLAB_WORLD PETSC_VIEWER_MATLAB_(PETSC_COMM_WORLD)

289: /*MC
290:   PETSC_VIEWER_MATLAB_SELF  - same as PETSC_VIEWER_MATLAB_(PETSC_COMM_SELF)

292:   Level: intermediate
293: M*/
294: #define PETSC_VIEWER_MATLAB_SELF  PETSC_VIEWER_MATLAB_(PETSC_COMM_SELF)

296: #define PETSC_VIEWER_MATHEMATICA_WORLD (PetscViewerInitializeMathematicaWorld_Private(),PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE)

300: PETSC_STATIC_INLINE PetscErrorCode PetscViewerFlowControlStart(PetscViewer viewer,PetscInt *mcnt,PetscInt *cnt)
301: {
304:   PetscViewerBinaryGetFlowControl(viewer,mcnt);
305:   PetscViewerBinaryGetFlowControl(viewer,cnt);
306:   return(0);
307: }

311: PETSC_STATIC_INLINE PetscErrorCode PetscViewerFlowControlStepMaster(PetscViewer viewer,PetscInt i,PetscInt *mcnt,PetscInt cnt)
312: {
314:   MPI_Comm       comm;

317:   PetscObjectGetComm((PetscObject)viewer,&comm);
318:   if (i >= *mcnt) {
319:     *mcnt += cnt;
320:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
321:   }
322:   return(0);
323: }

327: PETSC_STATIC_INLINE PetscErrorCode PetscViewerFlowControlEndMaster(PetscViewer viewer,PetscInt *mcnt)
328: {
330:   MPI_Comm       comm;
332:   PetscObjectGetComm((PetscObject)viewer,&comm);
333:   *mcnt = 0;
334:   MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
335:   return(0);
336: }

340: PETSC_STATIC_INLINE PetscErrorCode PetscViewerFlowControlStepWorker(PetscViewer viewer,PetscMPIInt rank,PetscInt *mcnt)
341: {
343:   MPI_Comm       comm;
345:   PetscObjectGetComm((PetscObject)viewer,&comm);
346:   while (PETSC_TRUE) {
347:     if (rank < *mcnt) break;
348:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
349:   }
350:   return(0);
351: }

355: PETSC_STATIC_INLINE PetscErrorCode PetscViewerFlowControlEndWorker(PetscViewer viewer,PetscInt *mcnt)
356: {
358:   MPI_Comm       comm;
360:   PetscObjectGetComm((PetscObject)viewer,&comm);
361:   while (PETSC_TRUE) {
362:     MPI_Bcast(mcnt,1,MPIU_INT,0,comm);
363:     if (!*mcnt) break;
364:   }
365:   return(0);
366: }

368: /*
369:    PetscViewer writes to MATLAB .mat file
370: */
371: PETSC_EXTERN PetscErrorCode PetscViewerMatlabPutArray(PetscViewer,int,int,const PetscScalar*,const char*);
372: PETSC_EXTERN PetscErrorCode PetscViewerMatlabGetArray(PetscViewer,int,int,PetscScalar*,const char*);
373: PETSC_EXTERN PetscErrorCode PetscViewerMatlabPutVariable(PetscViewer,const char*,void*);

375: #if defined(PETSC_HAVE_SAWS)
376: PETSC_EXTERN PetscErrorCode PetscObjectViewSAWs(PetscObject,PetscViewer);
377: #endif

379: /*S
380:      PetscViewers - Abstract collection of PetscViewers. It is just an expandable array of viewers.

382:    Level: intermediate

384:   Concepts: viewing

386: .seealso:  PetscViewerCreate(), PetscViewerSetType(), PetscViewerType, PetscViewer, PetscViewersCreate(),
387:            PetscViewersGetViewer()
388: S*/
389: typedef struct _n_PetscViewers* PetscViewers;
390: PETSC_EXTERN PetscErrorCode PetscViewersCreate(MPI_Comm,PetscViewers*);
391: PETSC_EXTERN PetscErrorCode PetscViewersDestroy(PetscViewers*);
392: PETSC_EXTERN PetscErrorCode PetscViewersGetViewer(PetscViewers,PetscInt,PetscViewer*);

394: /* Reset __FUNCT__ in case the user does not define it themselves */

398: #endif