Actual source code: glvis.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for fdopen() */

  3:  #include <petsc/private/viewerimpl.h>
  4:  #include <petsc/private/petscimpl.h>
  5:  #include <petsc/private/glvisviewerimpl.h>

  7: /* we may eventually make this function public */
  8: static PetscErrorCode PetscViewerASCIISocketOpen(MPI_Comm,const char*,PetscInt,PetscViewer*);

 10: struct _n_PetscViewerGLVis {
 11:   PetscViewerGLVisStatus status;
 12:   PetscViewerGLVisType   type;                                                  /* either PETSC_VIEWER_GLVIS_DUMP or PETSC_VIEWER_GLVIS_SOCKET */
 13:   char                   *name;                                                 /* prefix for filename, or hostname, depending on the type */
 14:   PetscInt               port;                                                  /* used just for the socket case */
 15:   PetscReal              pause;                                                 /* if positive, calls PetscSleep(pause) after each VecView_GLVis call */
 16:   PetscViewer            meshwindow;                                            /* used just by the ASCII dumping */
 17:   PetscObject            dm;                                                    /* DM as passed by PetscViewerGLVisSetDM_Private(): should contain discretization info */
 18:   PetscInt               nwindow;                                               /* number of windows/fields to be visualized */
 19:   PetscViewer            *window;
 20:   char                   **windowtitle;
 21:   char                   **fec_type;                                            /* type of elements to be used for visualization, see FiniteElementCollection::Name() */
 22:   PetscErrorCode         (*g2lfield)(PetscObject,PetscInt,PetscObject[],void*); /* global to local operation for generating dofs to be visualized */
 23:   PetscInt               *spacedim;                                             /* geometrical space dimension (just used to initialize the scene) */
 24:   PetscObject            *Ufield;                                               /* work vectors for visualization */
 25:   PetscInt               snapid;                                                /* snapshot id, use PetscViewerGLVisSetSnapId to change this value*/
 26:   void                   *userctx;                                              /* User context, used by g2lfield */
 27:   PetscErrorCode         (*destroyctx)(void*);                                  /* destroy routine for userctx */
 28:   char*                  fmt;                                                   /* format string for FP values */
 29: };
 30: typedef struct _n_PetscViewerGLVis *PetscViewerGLVis;

 32: /*@
 33:      PetscViewerGLVisSetPrecision - Set the number of digits for floating point values

 35:   Not Collective

 37:   Input Parameters:
 38: +  viewer - the PetscViewer
 39: -  prec   - the number of digits required

 41:   Level: beginner

 43: .seealso: PetscViewerGLVisOpen(), PetscViewerGLVisSetFields(), PetscViewerCreate(), PetscViewerSetType()
 44: @*/
 45: PetscErrorCode PetscViewerGLVisSetPrecision(PetscViewer viewer, PetscInt prec)
 46: {

 51:   PetscTryMethod(viewer,"PetscViewerGLVisSetPrecision_C",(PetscViewer,PetscInt),(viewer,prec));
 52:   return(0);
 53: }

 55: static PetscErrorCode PetscViewerGLVisSetPrecision_GLVis(PetscViewer viewer, PetscInt prec)
 56: {
 57:   PetscErrorCode   ierr;
 58:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

 61:   PetscFree(socket->fmt);
 62:   if (prec > 0) {
 63:     PetscMalloc1(16,&socket->fmt);
 64:     PetscSNPrintf(socket->fmt,16," %%.%De",prec);
 65:   } else {
 66:     PetscStrallocpy(" %g",&socket->fmt);
 67:   }
 68:   return(0);
 69: }

 71: /*@
 72:      PetscViewerGLVisSetSnapId - Set the snapshot id. Only relevant when the viewer is of type PETSC_VIEWER_GLVIS_DUMP

 74:   Logically Collective on PetscViewer

 76:   Input Parameters:
 77: +  viewer - the PetscViewer
 78: -  id     - the current snapshot id in a time-dependent simulation

 80:   Level: beginner

 82: .seealso: PetscViewerGLVisOpen(), PetscViewerGLVisSetFields(), PetscViewerCreate(), PetscViewerSetType()
 83: @*/
 84: PetscErrorCode PetscViewerGLVisSetSnapId(PetscViewer viewer, PetscInt id)
 85: {

 91:   PetscTryMethod(viewer,"PetscViewerGLVisSetSnapId_C",(PetscViewer,PetscInt),(viewer,id));
 92:   return(0);
 93: }

 95: static PetscErrorCode PetscViewerGLVisSetSnapId_GLVis(PetscViewer viewer, PetscInt id)
 96: {
 97:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

100:   socket->snapid = id;
101:   return(0);
102: }

104: /*@C
105:      PetscViewerGLVisSetFields - Sets the required information to visualize different fields from a vector.

107:   Logically Collective on PetscViewer

109:   Input Parameters:
110: +  viewer     - the PetscViewer
111: .  nf         - number of fields to be visualized
112: .  fec_type   - the type of finite element to be used to visualize the data (see FiniteElementCollection::Name() in MFEM)
113: .  dim        - array of space dimension for field vectors (used to initialize the scene)
114: .  g2lfields  - User routine to compute the local field vectors to be visualized; PetscObject is used in place of Vec on the prototype
115: .  Vfield     - array of work vectors, one for each field
116: .  ctx        - User context to store the relevant data to apply g2lfields
117: -  destroyctx - Destroy function for userctx

119:   Notes:
120:     g2lfields is called on the vector V to be visualized in order to extract the relevant dofs to be put in Vfield[], as
121: .vb
122:   g2lfields((PetscObject)V,nfields,(PetscObject*)Vfield[],ctx).
123: .ve
124:   For vector spaces, the block size of Vfield[i] represents the vector dimension. It misses the Fortran bindings.
125:   The names of the Vfield vectors will be displayed in the window title.

127:   Level: intermediate

129: .seealso: PetscViewerGLVisOpen(), PetscViewerCreate(), PetscViewerSetType(), PetscObjectSetName()
130: @*/
131: PetscErrorCode PetscViewerGLVisSetFields(PetscViewer viewer, PetscInt nf, const char* fec_type[], PetscInt dim[], PetscErrorCode(*g2l)(PetscObject,PetscInt,PetscObject[],void*), PetscObject Vfield[], void* ctx, PetscErrorCode(*destroyctx)(void*))
132: {

138:   if (!fec_type) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"You need to provide the FiniteElementCollection names for the fields");
142:   PetscTryMethod(viewer,"PetscViewerGLVisSetFields_C",(PetscViewer,PetscInt,const char*[],PetscInt[],PetscErrorCode(*)(PetscObject,PetscInt,PetscObject[],void*),PetscObject[],void*,PetscErrorCode(*)(void*)),(viewer,nf,fec_type,dim,g2l,Vfield,ctx,destroyctx));
143:   return(0);
144: }

146: static PetscErrorCode PetscViewerGLVisSetFields_GLVis(PetscViewer viewer, PetscInt nfields, const char* fec_type[], PetscInt dim[], PetscErrorCode(*g2l)(PetscObject,PetscInt,PetscObject[],void*), PetscObject Vfield[], void* ctx, PetscErrorCode(*destroyctx)(void*))
147: {
148:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
149:   PetscInt         i;
150:   PetscErrorCode   ierr;

153:   if (socket->nwindow && socket->nwindow != nfields) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Cannot set number of fields %D with number of windows %D",nfields,socket->nwindow);
154:   if (!socket->nwindow) {
155:     socket->nwindow = nfields;

157:     PetscCalloc5(nfields,&socket->window,nfields,&socket->windowtitle,nfields,&socket->fec_type,nfields,&socket->spacedim,nfields,&socket->Ufield);
158:     for (i=0;i<nfields;i++) {
159:       const char     *name;

161:       PetscObjectGetName(Vfield[i],&name);
162:       PetscStrallocpy(name,&socket->windowtitle[i]);
163:       PetscStrallocpy(fec_type[i],&socket->fec_type[i]);
164:       PetscObjectReference(Vfield[i]);
165:       socket->Ufield[i] = Vfield[i];
166:       socket->spacedim[i] = dim[i];
167:     }
168:   }
169:   /* number of fields are not allowed to vary */
170:   if (nfields != socket->nwindow) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Cannot visualize %D fields using %D socket windows",nfields,socket->nwindow);
171:   socket->g2lfield = g2l;
172:   if (socket->destroyctx && socket->userctx) { (*socket->destroyctx)(socket->userctx); }
173:   socket->userctx = ctx;
174:   socket->destroyctx = destroyctx;
175:   return(0);
176: }

178: static PetscErrorCode PetscViewerGLVisInfoDestroy_Private(void *ptr)
179: {
180:   PetscViewerGLVisInfo info = (PetscViewerGLVisInfo)ptr;
181:   PetscErrorCode       ierr;

184:   PetscFree(info->fmt);
185:   PetscFree(info);
186:   return(0);
187: }

189: /* we can decide to prevent specific processes from using the viewer */
190: static PetscErrorCode PetscViewerGLVisAttachInfo_Private(PetscViewer viewer, PetscViewer window)
191: {
192:   PetscViewerGLVis     socket = (PetscViewerGLVis)viewer->data;
193:   PetscErrorCode       ierr;
194:   PetscContainer       container;
195:   PetscViewerGLVisInfo info;

198:   PetscObjectQuery((PetscObject)window,"_glvis_info_container",(PetscObject*)&container);
199:   if (!container) {
200:     PetscNew(&info);
201:     info->enabled = PETSC_TRUE;
202:     info->init    = PETSC_FALSE;
203:     info->pause   = socket->pause;
204:     PetscContainerCreate(PetscObjectComm((PetscObject)window),&container);
205:     PetscContainerSetPointer(container,(void*)info);
206:     PetscContainerSetUserDestroy(container,PetscViewerGLVisInfoDestroy_Private);
207:     PetscObjectCompose((PetscObject)window,"_glvis_info_container",(PetscObject)container);
208:     PetscContainerDestroy(&container);
209:   } else {
210:     PetscContainerGetPointer(container,(void**)&info);
211:   }
212:   PetscFree(info->fmt);
213:   PetscStrallocpy(socket->fmt,&info->fmt);
214:   return(0);
215: }

217: static PetscErrorCode PetscViewerGLVisGetNewWindow_Private(PetscViewer viewer,PetscViewer *view)
218: {
219:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
220:   PetscViewer      window = NULL;
221:   PetscBool        ldis,dis;
222:   PetscErrorCode   ierr;

225:   PetscViewerASCIISocketOpen(PETSC_COMM_SELF,socket->name,socket->port,&window);
226:   /* if we could not estabilish a connection the first time,
227:      we disable the socket viewer */
228:   ldis = ierr ? PETSC_TRUE : PETSC_FALSE;
229:   MPI_Allreduce(&ldis,&dis,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)viewer));
230:   if (dis) {
231:     socket->status = PETSCVIEWERGLVIS_DISABLED;
232:     PetscViewerDestroy(&window);
233:   }
234:   *view = window;
235:   return(0);
236: }

238: PetscErrorCode PetscViewerGLVisPause_Private(PetscViewer viewer)
239: {
240:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
241:   PetscErrorCode   ierr;

244:   if (socket->pause > 0) {
245:     PetscSleep(socket->pause);
246:   }
247:   return(0);
248: }

250: /* DM specific support */
251: PetscErrorCode PetscViewerGLVisSetDM_Private(PetscViewer viewer, PetscObject dm)
252: {
253:   PetscErrorCode   ierr;
254:   PetscViewerGLVis socket  = (PetscViewerGLVis)viewer->data;

257:   if (socket->dm && socket->dm != dm) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Cannot change DM associated with the GLVis viewer");
258:   if (!socket->dm) {
259:     PetscErrorCode (*setupwithdm)(PetscObject,PetscViewer) = NULL;

261:     PetscObjectQueryFunction(dm,"DMSetUpGLVisViewer_C",&setupwithdm);
262:     if (setupwithdm) {
263:       (*setupwithdm)(dm,viewer);
264:     } else SETERRQ1(PetscObjectComm(dm),PETSC_ERR_SUP,"No support for DM type %s",dm->type_name);
265:     PetscObjectReference(dm);
266:     socket->dm = dm;
267:   }
268:   return(0);
269: }

271: PetscErrorCode PetscViewerGLVisGetDMWindow_Private(PetscViewer viewer,PetscViewer* view)
272: {
273:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
274:   PetscErrorCode   ierr;

277:   if (!socket->meshwindow) {
278:     if (socket->type == PETSC_VIEWER_GLVIS_SOCKET) {
279:       PetscViewerGLVisGetNewWindow_Private(viewer,&socket->meshwindow);
280:     } else {
281:       size_t    len;
282:       PetscBool isstdout;

284:       PetscStrlen(socket->name,&len);
285:       PetscStrcmp(socket->name,"stdout",&isstdout);
286:       if (!socket->name || !len || isstdout) {
287:         PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&socket->meshwindow);
288:       } else {
289:         PetscMPIInt rank;
290:         char        filename[PETSC_MAX_PATH_LEN];
291:         MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
292:         PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"%s-mesh.%06d",socket->name,rank);
293:         PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&socket->meshwindow);
294:       }
295:     }
296:     if (socket->meshwindow) {
297:       PetscViewerPushFormat(socket->meshwindow,PETSC_VIEWER_ASCII_GLVIS);
298:     }
299:   }
300:   if (socket->meshwindow) {
301:     PetscViewerGLVisAttachInfo_Private(viewer,socket->meshwindow);
302:   }
303:   *view = socket->meshwindow;
304:   return(0);
305: }

307: PetscErrorCode PetscViewerGLVisGetType_Private(PetscViewer viewer,PetscViewerGLVisType *type)
308: {
309:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

313:   *type = socket->type;
314:   return(0);
315: }

317: /* This function is only relevant in the SOCKET_GLIVS case. The status is computed the first time it is requested, as GLVis currently has issues when connecting the first time through the socket */
318: PetscErrorCode PetscViewerGLVisGetStatus_Private(PetscViewer viewer, PetscViewerGLVisStatus *sockstatus)
319: {
320:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

324:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP) {
325:     socket->status = PETSCVIEWERGLVIS_DISCONNECTED;
326:   } else if (socket->status == PETSCVIEWERGLVIS_DISCONNECTED && socket->nwindow) {
327:     PetscInt       i;
328:     PetscBool      lconn,conn;

331:     for (i=0,lconn=PETSC_TRUE;i<socket->nwindow;i++)
332:       if (!socket->window[i])
333:         lconn = PETSC_FALSE;

335:     MPI_Allreduce(&lconn,&conn,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)viewer));
336:     if (conn) socket->status = PETSCVIEWERGLVIS_CONNECTED;
337:   }
338:   *sockstatus = socket->status;
339:   return(0);
340: }

342: PetscErrorCode PetscViewerGLVisGetDM_Private(PetscViewer viewer, PetscObject* dm)
343: {
344:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

347:   *dm = socket->dm;
348:   return(0);
349: }

351: PetscErrorCode PetscViewerGLVisGetFields_Private(PetscViewer viewer, PetscInt* nfield, const char **fec[], PetscInt *spacedim[], PetscErrorCode(**g2lfield)(PetscObject,PetscInt,PetscObject[],void*), PetscObject *Ufield[], void **userctx)
352: {
353:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

356:   if (nfield)   *nfield   = socket->nwindow;
357:   if (fec)      *fec      = (const char**)socket->fec_type;
358:   if (spacedim) *spacedim = socket->spacedim;
359:   if (g2lfield) *g2lfield = socket->g2lfield;
360:   if (Ufield)   *Ufield   = socket->Ufield;
361:   if (userctx)  *userctx  = socket->userctx;
362:   return(0);
363: }

365: /* accessor routines for the viewer windows:
366:    PETSC_VIEWER_GLVIS_DUMP   : it returns a new viewer every time
367:    PETSC_VIEWER_GLVIS_SOCKET : it returns the socket, and creates it if not yet done.
368: */
369: PetscErrorCode PetscViewerGLVisGetWindow_Private(PetscViewer viewer,PetscInt wid,PetscViewer* view)
370: {
371:   PetscViewerGLVis       socket = (PetscViewerGLVis)viewer->data;
372:   PetscViewerGLVisStatus status;
373:   PetscErrorCode         ierr;

378:   if (wid < 0 || wid > socket->nwindow-1) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Cannot get window id %D: allowed range [0,%D)",wid,socket->nwindow-1);
379:   status = socket->status;
380:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP && socket->window[wid]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Window %D is already in use",wid);
381:   switch (status) {
382:     case PETSCVIEWERGLVIS_DISCONNECTED:
383:       if (socket->window[wid]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"This should not happen");
384:       if (socket->type == PETSC_VIEWER_GLVIS_DUMP) {
385:         size_t    len;
386:         PetscBool isstdout;

388:         PetscStrlen(socket->name,&len);
389:         PetscStrcmp(socket->name,"stdout",&isstdout);
390:         if (!socket->name || !len || isstdout) {
391:           PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&socket->window[wid]);
392:         } else {
393:           PetscMPIInt rank;
394:           char        filename[PETSC_MAX_PATH_LEN];

396:           MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
397:           PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"%s-%s-%d.%06d",socket->name,socket->windowtitle[wid],socket->snapid,rank);
398:           PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&socket->window[wid]);
399:         }
400:       } else {
401:         PetscViewerGLVisGetNewWindow_Private(viewer,&socket->window[wid]);
402:       }
403:       if (socket->window[wid]) {
404:         PetscViewerPushFormat(socket->window[wid],PETSC_VIEWER_ASCII_GLVIS);
405:       }
406:       *view = socket->window[wid];
407:       break;
408:     case PETSCVIEWERGLVIS_CONNECTED:
409:       *view = socket->window[wid];
410:       break;
411:     case PETSCVIEWERGLVIS_DISABLED:
412:       *view = NULL;
413:       break;
414:     default:
415:       SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unhandled socket status %d\n",(int)status);
416:       break;
417:   }
418:   if (*view) {
419:     PetscViewerGLVisAttachInfo_Private(viewer,*view);
420:   }
421:   return(0);
422: }

424: /* Restore the window viewer
425:    PETSC_VIEWER_GLVIS_DUMP  : destroys the temporary created ASCII viewer used for dumping
426:    PETSC_VIEWER_GLVIS_SOCKET: - if the returned window viewer is not NULL, just zeros the pointer.
427:                  - it the returned window viewer is NULL, assumes something went wrong
428:                    with the socket (i.e. SIGPIPE when a user closes the popup window)
429:                    and that the caller already handled it (see VecView_GLVis).
430: */
431: PetscErrorCode PetscViewerGLVisRestoreWindow_Private(PetscViewer viewer,PetscInt wid, PetscViewer* view)
432: {
433:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
434:   PetscErrorCode   ierr;

439:   if (wid < 0 || wid > socket->nwindow-1) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Cannot restore window id %D: allowed range [0,%D)",wid,socket->nwindow);
440:   if (*view && *view != socket->window[wid]) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Viewer was not obtained from PetscViewerGLVisGetWindow");
441:   if (*view) {
442:     PetscViewerFlush(*view);
443:     PetscBarrier((PetscObject)viewer);
444:   }
445:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP) { /* destroy the viewer, as it is associated with a single time step */
446:     PetscViewerDestroy(&socket->window[wid]);
447:   } else if (!*view) { /* something went wrong (SIGPIPE) so we just zero the private pointer */
448:     socket->window[wid] = NULL;
449:   }
450:   *view = NULL;
451:   return(0);
452: }

454: /* default window appearance in the PETSC_VIEWER_GLVIS_SOCKET case */
455: PetscErrorCode PetscViewerGLVisInitWindow_Private(PetscViewer viewer, PetscBool mesh, PetscInt dim, const char *name)
456: {
457:   PetscErrorCode       ierr;
458:   PetscViewerGLVisInfo info;
459:   PetscContainer       container;

462:   PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&container);
463:   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Viewer was not obtained from PetscGLVisViewerGetNewWindow_Private");
464:   PetscContainerGetPointer(container,(void**)&info);
465:   if (info->init) {
466:     if (info->pause < 0) {
467:       PetscViewerASCIIPrintf(viewer,"pause\n"); /* pause */
468:     }
469:     return(0);
470:   }
471:   PetscViewerASCIIPrintf(viewer,"window_size 800 800\n");
472:   if (name) {
473:     PetscViewerASCIIPrintf(viewer,"window_title '%s'\n",name);
474:   }
475:   if (mesh) {
476:     switch (dim) {
477:     case 1:
478:       break;
479:     case 2:
480:       PetscViewerASCIIPrintf(viewer,"keys cmeeppppp\n"); /* show colorbar, mesh and ranks */
481:       break;
482:     case 3: /* TODO: decide default view in 3D */
483:       break;
484:     }
485:   } else {
486:     PetscViewerASCIIPrintf(viewer,"keys cm\n"); /* show colorbar and mesh */
487:     switch (dim) {
488:     case 1:
489:       PetscViewerASCIIPrintf(viewer,"keys RRj\n"); /* set to 1D (side view) and turn off perspective */
490:       break;
491:     case 2:
492:       PetscViewerASCIIPrintf(viewer,"keys Rjl\n"); /* set to 2D (top view), turn off perspective and light */
493:       break;
494:     case 3:
495:       break;
496:     }
497:     PetscViewerASCIIPrintf(viewer,"autoscale value\n"); /* update value-range; keep mesh-extents fixed */
498:     if (info->pause == 0) {
499:       PetscViewerASCIIPrintf(viewer,"pause\n"); /* pause */
500:     }
501:   }
502:   info->init = PETSC_TRUE;
503:   return(0);
504: }

506: static PetscErrorCode PetscViewerDestroy_GLVis(PetscViewer viewer)
507: {
508:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
509:   PetscInt         i;
510:   PetscErrorCode   ierr;

513:   for (i=0;i<socket->nwindow;i++) {
514:     PetscViewerDestroy(&socket->window[i]);
515:     PetscFree(socket->windowtitle[i]);
516:     PetscFree(socket->fec_type[i]);
517:     PetscObjectDestroy(&socket->Ufield[i]);
518:   }
519:   PetscFree(socket->name);
520:   PetscFree5(socket->window,socket->windowtitle,socket->fec_type,socket->spacedim,socket->Ufield);
521:   PetscFree(socket->fmt);
522:   PetscViewerDestroy(&socket->meshwindow);
523:   PetscObjectDestroy(&socket->dm);
524:   if (socket->destroyctx && socket->userctx) { (*socket->destroyctx)(socket->userctx); }

526:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetPrecision_C",NULL);
527:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetSnapId_C",NULL);
528:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetFields_C",NULL);
529:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
530:   PetscFree(socket);
531:   viewer->data = NULL;
532:   return(0);
533: }

535: static PetscErrorCode PetscViewerSetFromOptions_GLVis(PetscOptionItems *PetscOptionsObject,PetscViewer v)
536: {
537:   PetscErrorCode   ierr;
538:   PetscViewerGLVis socket = (PetscViewerGLVis)v->data;

541:   PetscOptionsHead(PetscOptionsObject,"GLVis PetscViewer Options");
542:   PetscOptionsReal("-viewer_glvis_pause","-1 to pause after each visualization, otherwise sleeps for given seconds",NULL,socket->pause,&socket->pause,NULL);
543:   PetscOptionsTail();
544:   return(0);
545: }

547: static PetscErrorCode PetscViewerSetFileName_GLVis(PetscViewer viewer, const char name[])
548: {
549:   char             *sport;
550:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
551:   PetscErrorCode   ierr;

554:   socket->type = PETSC_VIEWER_GLVIS_DUMP;
555:   /* we accept localhost^port */
556:   PetscFree(socket->name);
557:   PetscStrallocpy(name,&socket->name);
558:   PetscStrchr(socket->name,'^',&sport);
559:   if (sport) {
560:     PetscInt port = 19916;
561:     size_t   len;

563:     *sport++ = 0;
564:     PetscStrlen(sport,&len);
565:     if (len) PetscOptionsStringToInt(sport,&port);
566:     if (!ierr) {
567:       socket->port = (port != PETSC_DECIDE && port != PETSC_DEFAULT) ? port : 19916;
568:     } else {
569:       socket->port = 19916;
570:     }
571:     socket->type = PETSC_VIEWER_GLVIS_SOCKET;
572:   }
573:   return(0);
574: }

576: /*
577:      PetscViewerGLVisOpen - Opens a GLVis type viewer

579:   Collective on comm

581:   Input Parameters:
582: +  comm      - the MPI communicator
583: .  type      - the viewer type: PETSC_VIEWER_GLVIS_SOCKET for real-time visualization or PETSC_VIEWER_GLVIS_DUMP for dumping to disk
584: .  name      - either the hostname where the GLVis server is running or the base filename for dumping the data for subsequent visualizations
585: -  port      - socket port where the GLVis server is listening. Not referenced when type is PETSC_VIEWER_GLVIS_DUMP

587:   Output Parameters:
588: -  viewer    - the PetscViewer object

590:   Notes:
591:     misses Fortran binding

593:   Level: beginner

595: .seealso: PetscViewerCreate(), PetscViewerSetType(), PetscViewerGLVisType
596: */
597: PETSC_EXTERN PetscErrorCode PetscViewerGLVisOpen(MPI_Comm comm, PetscViewerGLVisType type, const char* name, PetscInt port, PetscViewer* viewer)
598: {
599:   PetscViewerGLVis socket;
600:   PetscErrorCode   ierr;

603:   PetscViewerCreate(comm,viewer);
604:   PetscViewerSetType(*viewer,PETSCVIEWERGLVIS);

606:   socket       = (PetscViewerGLVis)((*viewer)->data);
607:   socket->type = type;
608:   if (type == PETSC_VIEWER_GLVIS_DUMP || name) {
609:     PetscFree(socket->name);
610:     PetscStrallocpy(name,&socket->name);
611:   }
612:   socket->port = (!port || port == PETSC_DETERMINE || port == PETSC_DECIDE) ? 19916 : port;

614:   PetscViewerSetFromOptions(*viewer);
615:   return(0);
616: }

618: /*
619:   PETSC_VIEWER_GLVIS_ - Creates an GLVIS PetscViewer shared by all processors in a communicator.

621:   Collective on MPI_Comm

623:   Input Parameter:
624: . comm - the MPI communicator to share the GLVIS PetscViewer

626:   Level: intermediate

628:   Notes:
629:     misses Fortran bindings

631:   Environmental variables:
632: + PETSC_VIEWER_GLVIS_FILENAME : output filename (if specified dump to disk, and takes precedence on PETSC_VIEWER_GLVIS_HOSTNAME)
633: . PETSC_VIEWER_GLVIS_HOSTNAME : machine where the GLVis server is listening (defaults to localhost)
634: - PETSC_VIEWER_GLVIS_PORT     : port opened by the GLVis server (defaults to 19916)

636:   Notes:
637:   Unlike almost all other PETSc routines, PETSC_VIEWER_GLVIS_ does not return
638:   an error code.  The GLVIS PetscViewer is usually used in the form
639: $       XXXView(XXX object, PETSC_VIEWER_GLVIS_(comm));

641: .seealso: PetscViewerGLVISOpen(), PetscViewerGLVisType, PetscViewerCreate(), PetscViewerDestroy()
642: */
643: PETSC_EXTERN PetscViewer PETSC_VIEWER_GLVIS_(MPI_Comm comm)
644: {
645:   PetscErrorCode       ierr;
646:   PetscBool            flg;
647:   PetscViewer          viewer;
648:   PetscViewerGLVisType type;
649:   char                 fname[PETSC_MAX_PATH_LEN],sport[16];
650:   PetscInt             port = 19916; /* default for GLVis */

653:   PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
654:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
655:   if (!flg) {
656:     type = PETSC_VIEWER_GLVIS_SOCKET;
657:     PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_HOSTNAME",fname,PETSC_MAX_PATH_LEN,&flg);
658:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
659:     if (!flg) {
660:       PetscStrcpy(fname,"localhost");
661:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
662:     }
663:     PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_PORT",sport,16,&flg);
664:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
665:     if (flg) {
666:       PetscOptionsStringToInt(sport,&port);
667:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
668:     }
669:   } else {
670:     type = PETSC_VIEWER_GLVIS_DUMP;
671:   }
672:   PetscViewerGLVisOpen(comm,type,fname,port,&viewer);
673:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
674:   PetscObjectRegisterDestroy((PetscObject)viewer);
675:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
676:   PetscFunctionReturn(viewer);
677: }

679: PETSC_EXTERN PetscErrorCode PetscViewerCreate_GLVis(PetscViewer viewer)
680: {
681:   PetscViewerGLVis socket;
682:   PetscErrorCode   ierr;

685:   PetscNewLog(viewer,&socket);

687:   /* defaults to socket viewer */
688:   PetscStrallocpy("localhost",&socket->name);
689:   socket->port  = 19916; /* GLVis default listening port */
690:   socket->type  = PETSC_VIEWER_GLVIS_SOCKET;
691:   socket->pause = 0; /* just pause the first time */

693:   /* defaults to full precision */
694:   PetscStrallocpy(" %g",&socket->fmt);

696:   viewer->data                = (void*)socket;
697:   viewer->ops->destroy        = PetscViewerDestroy_GLVis;
698:   viewer->ops->setfromoptions = PetscViewerSetFromOptions_GLVis;

700:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetPrecision_C",PetscViewerGLVisSetPrecision_GLVis);
701:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetSnapId_C",PetscViewerGLVisSetSnapId_GLVis);
702:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetFields_C",PetscViewerGLVisSetFields_GLVis);
703:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",PetscViewerSetFileName_GLVis);
704:   return(0);
705: }

707: /* this is a private implementation of a SOCKET with ASCII data format
708:    GLVis does not currently handle binary socket streams */
709: #if defined(PETSC_HAVE_UNISTD_H)
710: #include <unistd.h>
711: #endif

713: #if !defined(PETSC_HAVE_WINDOWS_H)
714: static PetscErrorCode (*PetscViewerDestroy_ASCII)(PetscViewer);

716: static PetscErrorCode PetscViewerDestroy_ASCII_Socket(PetscViewer viewer)
717: {
718:   FILE *stream;
719:   PetscErrorCode 0;
721:   PetscViewerASCIIGetPointer(viewer,&stream);
722:   if (stream) {
723:     fclose(stream);
724:     if (ierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on stream");
725:   }
726:   PetscViewerDestroy_ASCII(viewer);
727:   return(0);
728: }
729: #endif

731: static PetscErrorCode PetscViewerASCIISocketOpen(MPI_Comm comm,const char* hostname,PetscInt port,PetscViewer* viewer)
732: {
733: #if defined(PETSC_HAVE_WINDOWS_H)
735:   SETERRQ(comm,PETSC_ERR_SUP,"Not implemented for Windows");
736: #else
737:   FILE           *stream = NULL;
738:   int            fd;

744: #if defined(PETSC_USE_SOCKET_VIEWER)
745:   PetscOpenSocket(hostname,port,&fd);
746: #else
747:   SETERRQ(comm,PETSC_ERR_SUP,"Missing Socket viewer");
748: #endif
749:   if (ierr) {
750:     PetscInt sierr;
751:     char     err[1024];

753:     PetscSNPrintf(err,1024,"Cannot connect to socket on %s:%D. Socket visualization is disabled\n",hostname,port);
754:     PetscInfo(NULL,err);
755:     *viewer = NULL;
756:     PetscFunctionReturn(sierr);
757:   } else {
758:     char msg[1024];

760:     PetscSNPrintf(msg,1024,"Successfully connect to socket on %s:%D. Socket visualization is enabled\n",hostname,port);
761:     PetscInfo(NULL,msg);
762:   }
763:   stream = fdopen(fd,"w"); /* Not possible on Windows */
764:   if (!stream) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SYS,"Cannot open stream from socket %s:%d",hostname,port);
765:   PetscViewerASCIIOpenWithFILE(PETSC_COMM_SELF,stream,viewer);
766:   PetscViewerDestroy_ASCII = (*viewer)->ops->destroy;
767:   (*viewer)->ops->destroy = PetscViewerDestroy_ASCII_Socket;
768: #endif
769:   return(0);
770: }