Actual source code: glvis.c

petsc-3.11.4 2019-09-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:   PetscInt               windowsizes[2];
 22:   char                   **fec_type;                                            /* type of elements to be used for visualization, see FiniteElementCollection::Name() */
 23:   PetscErrorCode         (*g2lfield)(PetscObject,PetscInt,PetscObject[],void*); /* global to local operation for generating dofs to be visualized */
 24:   PetscInt               *spacedim;                                             /* geometrical space dimension (just used to initialize the scene) */
 25:   PetscObject            *Ufield;                                               /* work vectors for visualization */
 26:   PetscInt               snapid;                                                /* snapshot id, use PetscViewerGLVisSetSnapId to change this value*/
 27:   void                   *userctx;                                              /* User context, used by g2lfield */
 28:   PetscErrorCode         (*destroyctx)(void*);                                  /* destroy routine for userctx */
 29:   char*                  fmt;                                                   /* format string for FP values */
 30: };
 31: typedef struct _n_PetscViewerGLVis *PetscViewerGLVis;

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

 36:   Not Collective

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

 42:   Level: beginner

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

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

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

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

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

 75:   Logically Collective on PetscViewer

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

 81:   Level: beginner

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

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

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

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

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

108:   Logically Collective on PetscViewer

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

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

128:   Level: intermediate

130: .seealso: PetscViewerGLVisOpen(), PetscViewerCreate(), PetscViewerSetType(), PetscObjectSetName()
131: @*/
132: 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*))
133: {

139:   if (!fec_type) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"You need to provide the FiniteElementCollection names for the fields");
143:   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));
144:   return(0);
145: }

147: 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*))
148: {
149:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
150:   PetscInt         i;
151:   PetscErrorCode   ierr;

154:   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);
155:   if (!socket->nwindow) {
156:     socket->nwindow = nfields;

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

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

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

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

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

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

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

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

241: PetscErrorCode PetscViewerGLVisPause_Private(PetscViewer viewer)
242: {
243:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
244:   PetscErrorCode   ierr;

247:   if (socket->type == PETSC_VIEWER_GLVIS_SOCKET && socket->pause > 0) {
248:     PetscSleep(socket->pause);
249:   }
250:   return(0);
251: }

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

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

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

274: PetscErrorCode PetscViewerGLVisGetDMWindow_Private(PetscViewer viewer,PetscViewer *view)
275: {
276:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
277:   PetscErrorCode   ierr;

281:   if (!socket->meshwindow) {
282:     if (socket->type == PETSC_VIEWER_GLVIS_SOCKET) {
283:       PetscViewerGLVisGetNewWindow_Private(viewer,&socket->meshwindow);
284:     } else {
285:       size_t    len;
286:       PetscBool isstdout;

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

311: PetscErrorCode PetscViewerGLVisRestoreDMWindow_Private(PetscViewer viewer,PetscViewer *view)
312: {
313:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
314:   PetscErrorCode   ierr;

318:   if (*view && *view != socket->meshwindow) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Viewer was not obtained from PetscViewerGLVisGetDMWindow()");
319:   if (*view) {
320:     PetscViewerFlush(*view);
321:     PetscBarrier((PetscObject)viewer);
322:   }
323:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP) { /* destroy the viewer, as it is associated with a single time step */
324:     PetscViewerDestroy(&socket->meshwindow);
325:   } else if (!*view) { /* something went wrong (SIGPIPE) so we just zero the private pointer */
326:     socket->meshwindow = NULL;
327:   }
328:   *view = NULL;
329:   return(0);
330: }

332: PetscErrorCode PetscViewerGLVisGetType_Private(PetscViewer viewer,PetscViewerGLVisType *type)
333: {
334:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

338:   *type = socket->type;
339:   return(0);
340: }

342: /* 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 */
343: PetscErrorCode PetscViewerGLVisGetStatus_Private(PetscViewer viewer, PetscViewerGLVisStatus *sockstatus)
344: {
345:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

349:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP) {
350:     socket->status = PETSCVIEWERGLVIS_DISCONNECTED;
351:   } else if (socket->status == PETSCVIEWERGLVIS_DISCONNECTED && socket->nwindow) {
352:     PetscInt       i;
353:     PetscBool      lconn,conn;

356:     for (i=0,lconn=PETSC_TRUE;i<socket->nwindow;i++)
357:       if (!socket->window[i])
358:         lconn = PETSC_FALSE;

360:     MPIU_Allreduce(&lconn,&conn,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)viewer));
361:     if (conn) socket->status = PETSCVIEWERGLVIS_CONNECTED;
362:   }
363:   *sockstatus = socket->status;
364:   return(0);
365: }

367: PetscErrorCode PetscViewerGLVisGetDM_Private(PetscViewer viewer, PetscObject* dm)
368: {
369:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;

372:   *dm = socket->dm;
373:   return(0);
374: }

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

381:   if (nfield)   *nfield   = socket->nwindow;
382:   if (fec)      *fec      = (const char**)socket->fec_type;
383:   if (spacedim) *spacedim = socket->spacedim;
384:   if (g2lfield) *g2lfield = socket->g2lfield;
385:   if (Ufield)   *Ufield   = socket->Ufield;
386:   if (userctx)  *userctx  = socket->userctx;
387:   return(0);
388: }

390: /* accessor routines for the viewer windows:
391:    PETSC_VIEWER_GLVIS_DUMP   : it returns a new viewer every time
392:    PETSC_VIEWER_GLVIS_SOCKET : it returns the socket, and creates it if not yet done.
393: */
394: PetscErrorCode PetscViewerGLVisGetWindow_Private(PetscViewer viewer,PetscInt wid,PetscViewer* view)
395: {
396:   PetscViewerGLVis       socket = (PetscViewerGLVis)viewer->data;
397:   PetscViewerGLVisStatus status;
398:   PetscErrorCode         ierr;

403:   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);
404:   status = socket->status;
405:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP && socket->window[wid]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Window %D is already in use",wid);
406:   switch (status) {
407:     case PETSCVIEWERGLVIS_DISCONNECTED:
408:       if (socket->window[wid]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"This should not happen");
409:       if (socket->type == PETSC_VIEWER_GLVIS_DUMP) {
410:         size_t    len;
411:         PetscBool isstdout;

413:         PetscStrlen(socket->name,&len);
414:         PetscStrcmp(socket->name,"stdout",&isstdout);
415:         if (!socket->name || !len || isstdout) {
416:           PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&socket->window[wid]);
417:         } else {
418:           PetscMPIInt rank;
419:           char        filename[PETSC_MAX_PATH_LEN];

421:           MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
422:           PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"%s-%s-%d.%06d",socket->name,socket->windowtitle[wid],socket->snapid,rank);
423:           PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&socket->window[wid]);
424:         }
425:       } else {
426:         PetscViewerGLVisGetNewWindow_Private(viewer,&socket->window[wid]);
427:       }
428:       if (socket->window[wid]) {
429:         PetscViewerPushFormat(socket->window[wid],PETSC_VIEWER_ASCII_GLVIS);
430:       }
431:       *view = socket->window[wid];
432:       break;
433:     case PETSCVIEWERGLVIS_CONNECTED:
434:       *view = socket->window[wid];
435:       break;
436:     case PETSCVIEWERGLVIS_DISABLED:
437:       *view = NULL;
438:       break;
439:     default:
440:       SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unhandled socket status %d\n",(int)status);
441:       break;
442:   }
443:   if (*view) {
444:     PetscViewerGLVisAttachInfo_Private(viewer,*view);
445:   }
446:   return(0);
447: }

449: /* Restore the window viewer
450:    PETSC_VIEWER_GLVIS_DUMP  : destroys the temporary created ASCII viewer used for dumping
451:    PETSC_VIEWER_GLVIS_SOCKET: - if the returned window viewer is not NULL, just zeros the pointer.
452:                  - it the returned window viewer is NULL, assumes something went wrong
453:                    with the socket (i.e. SIGPIPE when a user closes the popup window)
454:                    and that the caller already handled it (see VecView_GLVis).
455: */
456: PetscErrorCode PetscViewerGLVisRestoreWindow_Private(PetscViewer viewer,PetscInt wid, PetscViewer* view)
457: {
458:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
459:   PetscErrorCode   ierr;

464:   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);
465:   if (*view && *view != socket->window[wid]) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_USER,"Viewer was not obtained from PetscViewerGLVisGetWindow()");
466:   if (*view) {
467:     PetscViewerFlush(*view);
468:     PetscBarrier((PetscObject)viewer);
469:   }
470:   if (socket->type == PETSC_VIEWER_GLVIS_DUMP) { /* destroy the viewer, as it is associated with a single time step */
471:     PetscViewerDestroy(&socket->window[wid]);
472:   } else if (!*view) { /* something went wrong (SIGPIPE) so we just zero the private pointer */
473:     socket->window[wid] = NULL;
474:   }
475:   *view = NULL;
476:   return(0);
477: }

479: /* default window appearance in the PETSC_VIEWER_GLVIS_SOCKET case */
480: PetscErrorCode PetscViewerGLVisInitWindow_Private(PetscViewer viewer, PetscBool mesh, PetscInt dim, const char *name)
481: {
482:   PetscErrorCode       ierr;
483:   PetscViewerGLVisInfo info;
484:   PetscContainer       container;

487:   PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&container);
488:   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Viewer was not obtained from PetscGLVisViewerGetNewWindow_Private");
489:   PetscContainerGetPointer(container,(void**)&info);
490:   if (info->init) return(0);

492:   /* Configure window */
493:   if (info->size[0] > 0) {
494:     PetscViewerASCIIPrintf(viewer,"window_size %D %D\n",info->size[0],info->size[1]);
495:   }
496:   if (name) {
497:     PetscViewerASCIIPrintf(viewer,"window_title '%s'\n",name);
498:   }

500:   /* Configure default view */
501:   if (mesh) {
502:     switch (dim) {
503:     case 1:
504:       PetscViewerASCIIPrintf(viewer,"keys m\n"); /* show mesh */
505:       break;
506:     case 2:
507:       PetscViewerASCIIPrintf(viewer,"keys m\n"); /* show mesh */
508:       break;
509:     case 3: /* TODO: decide default view in 3D */
510:       break;
511:     }
512:   } else {
513:     PetscViewerASCIIPrintf(viewer,"keys cm\n"); /* show colorbar and mesh */
514:     switch (dim) {
515:     case 1:
516:       PetscViewerASCIIPrintf(viewer,"keys RRjl\n"); /* set to 1D (side view), turn off perspective and light */
517:       break;
518:     case 2:
519:       PetscViewerASCIIPrintf(viewer,"keys Rjl\n"); /* set to 2D (top view), turn off perspective and light */
520:       break;
521:     case 3:
522:       break;
523:     }
524:     PetscViewerASCIIPrintf(viewer,"autoscale value\n"); /* update value-range; keep mesh-extents fixed */
525:   }

527:   { /* Additional keys and commands */
528:     char keys[256] = "", cmds[2*PETSC_MAX_PATH_LEN] = "";
529:     PetscOptions opt = ((PetscObject)viewer)->options;
530:     const char  *pre = ((PetscObject)viewer)->prefix;

532:     PetscOptionsGetString(opt,pre,"-glvis_keys",keys,sizeof(keys),NULL);
533:     PetscOptionsGetString(opt,pre,"-glvis_exec",cmds,sizeof(cmds),NULL);
534:     if (keys[0]) {PetscViewerASCIIPrintf(viewer,"keys %s\n",keys);}
535:     if (cmds[0]) {PetscViewerASCIIPrintf(viewer,"%s\n",cmds);}
536:   }

538:   /* Pause visualization */
539:   if (!mesh && info->pause == -1) {
540:     PetscViewerASCIIPrintf(viewer,"autopause 1\n");
541:   }
542:   if (!mesh && info->pause == 0) {
543:     PetscViewerASCIIPrintf(viewer,"pause\n");
544:   }

546:   info->init = PETSC_TRUE;
547:   return(0);
548: }

550: static PetscErrorCode PetscViewerDestroy_GLVis(PetscViewer viewer)
551: {
552:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
553:   PetscInt         i;
554:   PetscErrorCode   ierr;

557:   for (i=0;i<socket->nwindow;i++) {
558:     PetscViewerDestroy(&socket->window[i]);
559:     PetscFree(socket->windowtitle[i]);
560:     PetscFree(socket->fec_type[i]);
561:     PetscObjectDestroy(&socket->Ufield[i]);
562:   }
563:   PetscFree(socket->name);
564:   PetscFree5(socket->window,socket->windowtitle,socket->fec_type,socket->spacedim,socket->Ufield);
565:   PetscFree(socket->fmt);
566:   PetscViewerDestroy(&socket->meshwindow);
567:   PetscObjectDestroy(&socket->dm);
568:   if (socket->destroyctx && socket->userctx) { (*socket->destroyctx)(socket->userctx); }

570:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetPrecision_C",NULL);
571:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetSnapId_C",NULL);
572:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetFields_C",NULL);
573:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
574:   PetscFree(socket);
575:   viewer->data = NULL;
576:   return(0);
577: }

579: static PetscErrorCode PetscViewerSetFromOptions_GLVis(PetscOptionItems *PetscOptionsObject,PetscViewer v)
580: {
581:   PetscErrorCode   ierr;
582:   PetscViewerGLVis socket = (PetscViewerGLVis)v->data;
583:   PetscInt         nsizes = 2, prec = PETSC_DECIDE;
584:   PetscBool        set;

587:   PetscOptionsHead(PetscOptionsObject,"GLVis PetscViewer Options");
588:   PetscOptionsInt("-glvis_precision","Number of digits for floating point values","PetscViewerGLVisSetPrecision",prec,&prec,&set);
589:   if (set) {PetscViewerGLVisSetPrecision(v,prec);}
590:   PetscOptionsIntArray("-glvis_size","Window sizes",NULL,socket->windowsizes,&nsizes,&set);
591:   if (set && (nsizes == 1 || socket->windowsizes[1] < 0)) socket->windowsizes[1] = socket->windowsizes[0];
592:   PetscOptionsReal("-glvis_pause","-1 to pause after each visualization, otherwise sleeps for given seconds",NULL,socket->pause,&socket->pause,NULL);
593:   PetscOptionsName("-glvis_keys","Additional keys to configure visualization",NULL,NULL);
594:   PetscOptionsName("-glvis_exec","Additional commands to configure visualization",NULL,NULL);
595:   PetscOptionsTail();
596:   return(0);
597: }

599: static PetscErrorCode PetscViewerSetFileName_GLVis(PetscViewer viewer, const char name[])
600: {
601:   char             *sport;
602:   PetscViewerGLVis socket = (PetscViewerGLVis)viewer->data;
603:   PetscErrorCode   ierr;

606:   socket->type = PETSC_VIEWER_GLVIS_DUMP;
607:   /* we accept localhost^port */
608:   PetscFree(socket->name);
609:   PetscStrallocpy(name,&socket->name);
610:   PetscStrchr(socket->name,'^',&sport);
611:   if (sport) {
612:     PetscInt port = 19916;
613:     size_t   len;

615:     *sport++ = 0;
616:     PetscStrlen(sport,&len);
617:     if (len) PetscOptionsStringToInt(sport,&port);
618:     if (!ierr) {
619:       socket->port = (port != PETSC_DECIDE && port != PETSC_DEFAULT) ? port : 19916;
620:     } else {
621:       socket->port = 19916;
622:     }
623:     socket->type = PETSC_VIEWER_GLVIS_SOCKET;
624:   }
625:   return(0);
626: }

628: /*@C
629:   PetscViewerGLVisOpen - Opens a GLVis type viewer

631:   Collective on MPI_Comm

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

639:   Output Parameters:
640: -  viewer    - the PetscViewer object

642:   Options Database Keys:
643: +  -glvis_precision <precision> - Sets number of digits for floating point values
644: .  -glvis_size <width,height> - Sets the window size (in pixels)
645: .  -glvis_pause <pause> - Sets time (in seconds) that the program pauses after each visualization
646:        (0 is default, -1 implies every visualization)
647: .  -glvis_keys - Additional keys to configure visualization
648: -  -glvis_exec - Additional commands to configure visualization

650:   Notes:
651:     misses Fortran binding

653:   Level: beginner

655: .seealso: PetscViewerCreate(), PetscViewerSetType(), PetscViewerGLVisType
656: @*/
657: PetscErrorCode PetscViewerGLVisOpen(MPI_Comm comm, PetscViewerGLVisType type, const char name[], PetscInt port, PetscViewer *viewer)
658: {
659:   PetscViewerGLVis socket;
660:   PetscErrorCode   ierr;

663:   PetscViewerCreate(comm,viewer);
664:   PetscViewerSetType(*viewer,PETSCVIEWERGLVIS);

666:   socket       = (PetscViewerGLVis)((*viewer)->data);
667:   socket->type = type;
668:   if (type == PETSC_VIEWER_GLVIS_DUMP || name) {
669:     PetscFree(socket->name);
670:     PetscStrallocpy(name,&socket->name);
671:   }
672:   socket->port = (!port || port == PETSC_DETERMINE || port == PETSC_DECIDE) ? 19916 : port;

674:   PetscViewerSetFromOptions(*viewer);
675:   return(0);
676: }

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

681:   Collective on MPI_Comm

683:   Input Parameter:
684: . comm - the MPI communicator to share the GLVIS PetscViewer

686:   Level: intermediate

688:   Notes:
689:     misses Fortran bindings

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

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

701: .seealso: PetscViewerGLVISOpen(), PetscViewerGLVisType, PetscViewerCreate(), PetscViewerDestroy()
702: */
703: PetscViewer PETSC_VIEWER_GLVIS_(MPI_Comm comm)
704: {
705:   PetscErrorCode       ierr;
706:   PetscBool            flg;
707:   PetscViewer          viewer;
708:   PetscViewerGLVisType type;
709:   char                 fname[PETSC_MAX_PATH_LEN],sport[16];
710:   PetscInt             port = 19916; /* default for GLVis */

713:   PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
714:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
715:   if (!flg) {
716:     type = PETSC_VIEWER_GLVIS_SOCKET;
717:     PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_HOSTNAME",fname,PETSC_MAX_PATH_LEN,&flg);
718:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
719:     if (!flg) {
720:       PetscStrcpy(fname,"localhost");
721:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
722:     }
723:     PetscOptionsGetenv(comm,"PETSC_VIEWER_GLVIS_PORT",sport,16,&flg);
724:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
725:     if (flg) {
726:       PetscOptionsStringToInt(sport,&port);
727:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
728:     }
729:   } else {
730:     type = PETSC_VIEWER_GLVIS_DUMP;
731:   }
732:   PetscViewerGLVisOpen(comm,type,fname,port,&viewer);
733:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
734:   PetscObjectRegisterDestroy((PetscObject)viewer);
735:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_GLVIS_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
736:   PetscFunctionReturn(viewer);
737: }

739: PETSC_EXTERN PetscErrorCode PetscViewerCreate_GLVis(PetscViewer viewer)
740: {
741:   PetscViewerGLVis socket;
742:   PetscErrorCode   ierr;

745:   PetscNewLog(viewer,&socket);

747:   /* defaults to socket viewer */
748:   PetscStrallocpy("localhost",&socket->name);
749:   socket->port  = 19916; /* GLVis default listening port */
750:   socket->type  = PETSC_VIEWER_GLVIS_SOCKET;
751:   socket->pause = 0; /* just pause the first time */

753:   socket->windowsizes[0] = 600;
754:   socket->windowsizes[1] = 600;

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

759:   viewer->data                = (void*)socket;
760:   viewer->ops->destroy        = PetscViewerDestroy_GLVis;
761:   viewer->ops->setfromoptions = PetscViewerSetFromOptions_GLVis;

763:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetPrecision_C",PetscViewerGLVisSetPrecision_GLVis);
764:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetSnapId_C",PetscViewerGLVisSetSnapId_GLVis);
765:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGLVisSetFields_C",PetscViewerGLVisSetFields_GLVis);
766:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",PetscViewerSetFileName_GLVis);
767:   return(0);
768: }

770: /* this is a private implementation of a SOCKET with ASCII data format
771:    GLVis does not currently handle binary socket streams */
772: #if defined(PETSC_HAVE_UNISTD_H)
773: #include <unistd.h>
774: #endif

776: #if !defined(PETSC_HAVE_WINDOWS_H)
777: static PetscErrorCode (*PetscViewerDestroy_ASCII)(PetscViewer);

779: static PetscErrorCode PetscViewerDestroy_ASCII_Socket(PetscViewer viewer)
780: {
781:   FILE *stream;
782:   PetscErrorCode 0;
784:   PetscViewerASCIIGetPointer(viewer,&stream);
785:   if (stream) {
786:     fclose(stream);
787:     if (ierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on stream");
788:   }
789:   PetscViewerDestroy_ASCII(viewer);
790:   return(0);
791: }
792: #endif

794: static PetscErrorCode PetscViewerASCIISocketOpen(MPI_Comm comm,const char* hostname,PetscInt port,PetscViewer* viewer)
795: {
796: #if defined(PETSC_HAVE_WINDOWS_H)
798:   SETERRQ(comm,PETSC_ERR_SUP,"Not implemented for Windows");
799: #else
800:   FILE           *stream = NULL;
801:   int            fd;

807: #if defined(PETSC_USE_SOCKET_VIEWER)
808:   PetscOpenSocket(hostname,port,&fd);
809: #else
810:   SETERRQ(comm,PETSC_ERR_SUP,"Missing Socket viewer");
811: #endif
812:   if (ierr) {
813:     PetscInt sierr;
814:     char     err[1024];

816:     PetscSNPrintf(err,1024,"Cannot connect to socket on %s:%D. Socket visualization is disabled\n",hostname,port);
817:     PetscInfo(NULL,err);
818:     *viewer = NULL;
819:     PetscFunctionReturn(sierr);
820:   } else {
821:     char msg[1024];

823:     PetscSNPrintf(msg,1024,"Successfully connect to socket on %s:%D. Socket visualization is enabled\n",hostname,port);
824:     PetscInfo(NULL,msg);
825:   }
826:   stream = fdopen(fd,"w"); /* Not possible on Windows */
827:   if (!stream) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SYS,"Cannot open stream from socket %s:%d",hostname,port);
828:   PetscViewerASCIIOpenWithFILE(PETSC_COMM_SELF,stream,viewer);
829:   PetscViewerDestroy_ASCII = (*viewer)->ops->destroy;
830:   (*viewer)->ops->destroy = PetscViewerDestroy_ASCII_Socket;
831: #endif
832:   return(0);
833: }

835: #if !defined(PETSC_MISSING_SIGPIPE)

837: #include <signal.h>

839: #if defined(PETSC_HAVE_WINDOWS_H)
840: #define PETSC_DEVNULL "NUL"
841: #else
842: #define PETSC_DEVNULL "/dev/null"
843: #endif

845: static volatile PetscBool PetscGLVisBrokenPipe = PETSC_FALSE;

847: static void (*PetscGLVisSigHandler_save)(int) = NULL;

849: static void PetscGLVisSigHandler_SIGPIPE(PETSC_UNUSED int sig)
850: {
851:   PetscGLVisBrokenPipe = PETSC_TRUE;
852: #if !defined(PETSC_MISSING_SIG_IGN)
853:   signal(SIGPIPE,SIG_IGN);
854: #endif
855: }

857: PetscErrorCode PetscGLVisCollectiveBegin(PETSC_UNUSED MPI_Comm comm,PETSC_UNUSED PetscViewer *win)
858: {
860:   if (PetscGLVisSigHandler_save) SETERRQ1(comm,PETSC_ERR_PLIB,"Nested call to %s()",PETSC_FUNCTION_NAME);
861:   PetscGLVisBrokenPipe = PETSC_FALSE;
862:   PetscGLVisSigHandler_save = signal(SIGPIPE,PetscGLVisSigHandler_SIGPIPE);
863:   return(0);
864: }

866: PetscErrorCode PetscGLVisCollectiveEnd(MPI_Comm comm,PetscViewer *win)
867: {
868:   PetscBool      flag,brokenpipe;

872:   flag = PetscGLVisBrokenPipe;
873:   MPIU_Allreduce(&flag,&brokenpipe,1,MPIU_BOOL,MPI_LOR,comm);
874:   if (brokenpipe) {
875:     FILE *sock, *null = fopen(PETSC_DEVNULL,"w");
876:     PetscViewerASCIIGetPointer(*win,&sock);
877:     PetscViewerASCIISetFILE(*win,null);
878:     PetscViewerDestroy(win);
879:     if (sock) (void)fclose(sock);
880:   }
881:   (void)signal(SIGPIPE,PetscGLVisSigHandler_save);
882:   PetscGLVisSigHandler_save = NULL;
883:   PetscGLVisBrokenPipe = PETSC_FALSE;
884:   return(0);
885: }

887: #else

889: PetscErrorCode PetscGLVisCollectiveBegin(PETSC_UNUSED MPI_Comm comm,PETSC_UNUSED PetscViewer *win)
890: {
892:   return(0);
893: }

895: PetscErrorCode PetscGLVisCollectiveEnd(PETSC_UNUSED MPI_Comm comm,PETSC_UNUSED PetscViewer *win)
896: {
898:   return(0);
899: }

901: #endif