Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
display.h
Go to the documentation of this file.
1/** @file display.h
2 */
3/**
4# Server-side display
5
6This file implements the server-side of the interactive display of a
7running Basilisk code. The client-side is typically done by the
8[javascript implementation](jview/README).
9
10This works using two main components:
11
12* [The wsServer WebSocket library](wsServer/README.md)
13* [Vertex buffers](vertexbuffer.h)
14
15The initial state of the display can be controlled using the DISPLAY
16macro, following these rules:
17
18* `-DDISPLAY=1, -DDISPLAY`: play controls, start running immediately.
19* `-DDISPLAY=-1`: play controls, initially paused.
20* `-DDISPLAY=0`: no play controls, start running immediately.
21* `#include "display.h"`: play controls, initially paused.
22
23This can also be changed by setting `display_play = 0;` in `main()`.
24
25## Display port
26
27Multiple servers can run simultaneously on the same machine but must
28use different ports. When connecting, the server looks for a free port
29in the default range (typically 7100 to 7200, see the DISPLAY_RANGE
30macro below). The corresponding connection information (including the
31URL of the [javascript interface](/src/jview/README)) is written in
32the `display.html` file, which can be used to open the connection on
33the client side (i.e. the web browser).
34
35## Display control
36
37The values of integer or double variables can be controlled
38interactively using the `display_control()` macro, for example:
39
40~~~literatec
41display_control (value, 0, 100, "Short Name", "Long tooltip");
42~~~
43
44where the first argument is the name of the variable and the numbers
45indicate the (optional) range (min and max values). The (optional)
46name and tooltip are used to generate the [Basilisk
47View](jview/README) interface. The only compulsory argument is the
48variable name.
49
50## Default display
51
52Default [drawing function calls](draw.h), which will be sent to an empty
53[javascript client](jview/README), can be setup using something like:
54
55~~~literatec
56display ("squares (color = 'u.x', spread = -1);");
57~~~
58
59where the single quotes will be expanded to doubles quotes when
60interpreted. Do not forget the trailing semi-column.
61
62Note that these function calls will be appended to existing ones. To
63overwrite existing calls, use:
64
65~~~literatec
66display ("squares (color = 'u.x', spread = -1);", true);
67~~~
68*/
69
70#ifndef DISPLAY_JS
71# define DISPLAY_JS "http://basilisk.ida.upmc.fr/three.js/editor/index.html"
72#endif
73
74#ifndef DISPLAY_HOST
75# define DISPLAY_HOST "localhost"
76#endif
77
78#ifndef DISPLAY_RANGE
79# define DISPLAY_RANGE "7100:7200"
80#endif
81
82#define debug(...) do { if (Display.debug) fprintf (qerr, __VA_ARGS__), fflush(qerr); } while(0)
83
84#include <netdb.h>
85#include <wsServer/include/ws.h>
86#pragma autolink -L$BASILISK/wsServer -lws
87
88#include "view.h"
89#include "khash.h"
90
91typedef struct {
92 int fd;
93 int iter;
95
97
98static struct {
100 int sock, port;
101 char * error;
103 bool debug;
104} Display = { .sock = -1 };
105
106static void display_display()
107{
108 debug ("**************************\n");
109 for (khiter_t k = kh_begin (Display.objects); k != kh_end (Display.objects);
110 ++k)
111 if (kh_exist (Display.objects, k)) {
112 debug ("%s", kh_key (Display.objects, k));
113 DisplayClient * client = kh_value (Display.objects, k);
114 while (client->fd >= 0) {
115 debug (" %d/%d", client->fd, client->iter);
116 client++;
117 }
118 debug ("\n");
119 }
120 debug ("--------------------------\n");
121}
122
124{
125 if (fseek (fp, 0L, SEEK_END) < 0)
126 return NULL;
127 long bufsize = ftell (fp);
128 if (bufsize <= 0)
129 return NULL;
130 char * buf = malloc (sizeof(char)*(bufsize + 1));
131 if (fseek (fp, 0L, SEEK_SET) < 0) {
132 free (buf);
133 return NULL;
134 }
135 size_t newLen = fread (buf, sizeof(char), bufsize, fp);
136 buf[newLen] = '\0'; /* Just to be safe. */
137 return buf;
138}
139
140static void display_command (const char * command)
141{
142 debug ("display_command (%s)\n", command);
143
145 VertexBuffer.vertex = 0; // fixme
146
147 // Temporarily redirect stderr to catch errors
148 int bak = -1;
149 if (pid() == 0) {
150 fflush (stderr);
151 bak = dup (2);
152 FILE * fp = tmpfile();
153 dup2 (fileno (fp), 2);
154 fclose (fp);
155 }
156
157 char * line = strdup (command);
158 bool status = process_line (line);
159 free (line);
160
161 free (Display.error);
162 Display.error = NULL;
163 if (status) {
164 if (VertexBuffer.type < 0)
165 VertexBuffer.type = 0;
166 }
167 else { // An error occured
168 if (pid() == 0) {
169 fflush (stderr);
170 FILE * fp = fdopen (2, "r");
172 int len = Display.error ? strlen (Display.error) : 0;
173 if (len > 0 && Display.error[len - 1] == '\n')
174 Display.error[len - 1] = '\0';
175 fclose(fp);
176 }
177 else
178 Display.error = strdup ("error (on slave process)");
179 VertexBuffer.type = - 1;
180 }
181
182 if (pid() == 0) {
183 // Restore stderr to its previous value
184 dup2 (bak, 2);
185 close (bak);
186 }
187
188 if (VertexBuffer.type < 0)
189 debug ("error: '%s'\n", Display.error);
190 else {
191 if (pid() > 0) {
192 if (VertexBuffer.normal->len < VertexBuffer.position->len)
193 VertexBuffer.normal->len = 0;
194 if (VertexBuffer.color->len < VertexBuffer.position->len)
195 VertexBuffer.color->len = 0;
196 }
197 debug ("position: %ld, normal: %ld, color: %ld, index: %ld, type: %d\n",
198 VertexBuffer.position->len,
199 VertexBuffer.normal->len,
200 VertexBuffer.color->len,
201 VertexBuffer.index->len, VertexBuffer.type);
202 }
203}
204
205static int ws_send_array (int fd, Array * a, int status, int type,
206 unsigned int * shift)
207{
208#if _MPI
209 if (pid() == 0) {
210 void * p = NULL;
211 long len;
212 for (int pe = 0; pe < npe(); pe++) {
213 if (pe == 0)
214 p = a->p, len = a->len;
215 else {
217 MPI_Recv (&len, 1, MPI_LONG, pe, 22, MPI_COMM_WORLD, &status);
218 if (len > 0) {
219 p = malloc (len);
220 MPI_Recv (p, len, MPI_BYTE, pe, 23, MPI_COMM_WORLD, &status);
221 }
222 else
223 p = NULL;
224 }
225 if (type == 0) // position
226 shift[pe] = (pe > 0 ? shift[pe - 1] : 0) + len/(3*sizeof(float));
227 else if (type == 1 && pe > 0) // index
228 for (unsigned int i = 0; i < len/sizeof(unsigned int); i++)
229 ((unsigned int *) p)[i] += shift[pe - 1];
230 if (status >= 0 && len > 0 && ws_send (fd, p, len) < len)
231 status = -1;
232 if (pe > 0)
233 free (p);
234 }
235 }
236 else { // pid() > 0
237 MPI_Send (&a->len, 1, MPI_LONG, 0, 22, MPI_COMM_WORLD);
238 if (a->len > 0)
239 MPI_Send (a->p, a->len, MPI_BYTE, 0, 23, MPI_COMM_WORLD);
240 }
241#else
242 if (status >= 0 && a->len > 0 && ws_send (fd, a->p, a->len) < a->len)
243 status = -1;
244#endif
245 return status;
246}
247
248static int display_send (const char * command, int fd)
249{
250 int status = 0;
251
252 debug ("sending '%s' to %d\n", command, fd);
253
254 unsigned int commandlen = strlen (command);
255 unsigned int errorlen = Display.error ? strlen (Display.error) : 0;
256
257 int paddedlen = 4*ceil(commandlen/4.);
258 size_t len = 2*sizeof(unsigned int) + paddedlen;
259
260 unsigned int lens[] = {VertexBuffer.position->len,
261 VertexBuffer.normal->len,
262 VertexBuffer.color->len,
263 VertexBuffer.index->len}, glens[4];
264 int type = VertexBuffer.type, gtype;
265#if _MPI
268#else
269 for (int i = 0; i < 4; i++) glens[i] = lens[i];
270 gtype = type;
271#endif
272
273 if (gtype < 0)
274 len += errorlen;
275 else
276 len += 2*sizeof(int) +
277 4*sizeof (unsigned int) + glens[0] + glens[1] + glens[2] + glens[3];
278
279 if (pid() == 0) {
280 if (ws_sendframe_init (fd, len, false, WS_FR_OP_BIN) < 0 ||
281 ws_send (fd, &commandlen, sizeof(unsigned int)) < sizeof(unsigned int) ||
282 ws_send (fd, &errorlen, sizeof(unsigned int)) < sizeof(unsigned int) ||
284 status = -1;
285
286 // padding to a multiple of four
287 for (int i = 0; i < paddedlen - commandlen && status >= 0; i++) {
288 char c = '\0';
289 if (ws_send (fd, &c, 1) < 1)
290 status = -1;
291 }
292 }
293
294 if (gtype < 0) {
295 if (pid() == 0 && status >= 0 &&
296 ws_send (fd, Display.error, errorlen) < errorlen)
297 status = -1;
298 }
299 else {
300 if (pid() == 0 && status >= 0 &&
301 (ws_send (fd, &VertexBuffer.dim, sizeof(int)) < sizeof(int) ||
302 ws_send (fd, &gtype, sizeof(int)) < sizeof(int) ||
303 ws_send (fd, glens, 4*sizeof (unsigned int)) < 4*sizeof (unsigned int)))
304 status = -1;
305 unsigned int * shift = malloc (sizeof(unsigned int)*npe());
306 status = ws_send_array (fd, VertexBuffer.position, status, 0, shift);
307 status = ws_send_array (fd, VertexBuffer.normal, status, -1, shift);
308 status = ws_send_array (fd, VertexBuffer.color, status, -1, shift);
309 status = ws_send_array (fd, VertexBuffer.index, status, 1, shift);
310 free (shift);
311 }
312
313#if _MPI
315#endif
316
317 return status;
318}
319
320static void display_add (const char * command, int fd)
321{
322 debug ("adding '%s'\n", command);
323 khiter_t k = kh_get (strhash, Display.objects, command);
324 if (k == kh_end (Display.objects)) {
325 int ret;
326 k = kh_put (strhash, Display.objects, strdup (command), &ret);
328 client->fd = -1;
329 kh_value (Display.objects, k) = client;
330 }
332 int len = 0;
333 while (client->fd >= 0) {
334 if (client->fd == fd)
335 fd = -1;
336 client++, len++;
337 }
338 if (fd >= 0) { // not already in the array
339 kh_value (Display.objects, k) = clients =
340 realloc (clients, (len + 2)*sizeof (DisplayClient));
341 clients[len].fd = fd;
342 clients[len].iter = -1;
343 clients[len + 1].fd = -1;
344 }
346}
347
348#define JSON_BUILD(...) len += snprintf (build + len, 4096, __VA_ARGS__), \
349 build = realloc (build, len + 4096)
350
351static void array_remove (khiter_t k, int fd)
352{
354 int i = -1, len = 0;
355 while (client->fd >= 0) {
356 if (client->fd == fd) {
357 if (i != -1)
358 debug ("array_remove(): error! found multiple %d in '%s'\n",
359 fd, kh_key (Display.objects, k));
360 i = len;
361 }
362 client++, len++;
363 }
364 if (i < 0)
365 debug ("array_remove(): error! could not find %d in '%s'\n",
366 fd, kh_key (Display.objects, k));
367 else if (len == 1) {
368 free ((void *) kh_key (Display.objects, k));
369 free ((void *) kh_value (Display.objects, k));
370 kh_del (strhash, Display.objects, k);
371 }
372 else
373 for (int j = i; j < len; j++)
374 clients[j] = clients[j + 1];
375}
376
377static void display_remove (const char * command, int fd)
378{
379 debug ("removing '%s'\n", command);
380 khiter_t k = kh_get (strhash, Display.objects, command);
381 if (k == kh_end (Display.objects))
382 debug ("display_remove(): error! could not find '%s' (%d)\n",
383 command, fd);
384 else
385 array_remove (k, fd);
387}
388
389typedef struct {
390 char * name, * tooltip;
391 void * ptr;
392 double min, max;
393 int size;
395
396static char * display_control_json()
397{
398 char * build = malloc (4096);
399 int len = 0;
400 JSON_BUILD ("#{");
401 char sep = ' ';
402 DisplayControl * d = Display.controls->p;
403 for (int i = 0; i < Display.controls->len/sizeof(DisplayControl); i++, d++) {
404 JSON_BUILD ("%c\n \"%s\": { ", sep, d->name); sep = ',';
405 if (d->tooltip)
406 JSON_BUILD ("\"tooltip\": \"%s\", ", d->tooltip);
407 switch (d->size) {
408 case 4:
409 JSON_BUILD ("\"type\": \"int\", \"value\": %d, \"min\": %g, \"max\": %g",
410 *((int *)d->ptr), d->min, d->max);
411 break;
412 case 8:
413 JSON_BUILD ("\"type\": \"double\", \"value\": %g, "
414 "\"min\": %g, \"max\": %g",
415 *((double *)d->ptr), d->min, d->max);
416 break;
417 default:
418 assert (false);
419 }
420 JSON_BUILD (" }");
421 }
422 JSON_BUILD ("}");
423 return build;
424}
425
426static DisplayControl * display_control_lookup (const char * name)
427{
428 DisplayControl * d = Display.controls->p;
429 for (int i = 0; i < Display.controls->len/sizeof(DisplayControl); i++, d++)
430 if (!strcmp (d->name, name))
431 return d;
432 return NULL;
433}
434
436 double min, double max,
437 char * name = NULL, char * tooltip = NULL,
438 char * ptr_name, int size)
439{
441 if (!name)
442 name = ptr_name;
443
444 if (display_control_lookup (name))
445 return;
446
447 d.name = strdup (name);
448 d.tooltip = tooltip ? strdup (tooltip) : NULL;
449 d.ptr = ptr;
450 d.size = size;
451 d.min = min, d.max = max;
453
454 if (pid() == 0) {
456 ws_sendframe_txt (0, controls, true);
457 free (controls);
458 }
459}
460
461#undef display_control
462#define display_control(val, ...) display_control_internal \
463 (&(val), __VA_ARGS__, size = sizeof(val), ptr_name = #val)
464
465static void display_control_update (const char * command, int fd)
466{
467 char * s = strdup (command), * s1 = strchr (command, ':');
468 *s1++ = '\0';
470 if (d == NULL)
471 debug ("display_control_update(): error! could not find '%s' (%d)\n",
472 command, fd);
473 else {
474 debug ("display_control_update (%s) = %s\n", command, s1);
475 double val = atof(s1);
476 if (d->max > d->min)
477 val = clamp (val, d->min, d->max);
478 switch (d->size) {
479 case 4: *((int *)d->ptr) = val; break;
480 case 8: *((double *)d->ptr) = val; break;
481 default: assert (false);
482 }
483
484 if (pid() == 0) {
486 ws_sendframe_txt (- fd, controls, true);
487 free (controls);
488 }
489 }
490 free (s);
491}
492
493static char * bview_interface_json()
494{
495 char * build = malloc (4096);
496 int len = 0;
497
498 JSON_BUILD ("{\n");
499
500 int i = 0;
501 while (bview_interface[i].json) {
502 JSON_BUILD ("%s", i ? ",\n" : "");
503 len += bview_interface[i].json (build + len, 4096);
504 build = realloc (build, len + 4096);
505 JSON_BUILD ("\n");
506 i++;
507 }
508 JSON_BUILD ("}");
509 return build;
510}
511
512void display_onclose (int fd)
513{
514 debug ("closing %d\n", fd);
515 for (khiter_t k = kh_begin (Display.objects); k != kh_end (Display.objects);
516 ++k)
517 if (kh_exist (Display.objects, k))
518 array_remove (k, fd);
520}
521
522void display_onmessage (int fd, const char * msg, size_t size, int type)
523{
524 if (type == WS_FR_OP_TXT) {
525 if (!msg)
526 fprintf (stderr, "error receiving data on websocket\n");
527 else switch (msg[0]) {
528 case '+': display_add (msg + 1, fd); break;
529 case '-': display_remove (msg + 1, fd); break;
530 case '#': display_control_update (msg + 1, fd); break;
531 default: fprintf (stderr,
532 "display_onmessage: error: unknown message type '%s'\n",
533 msg);
534 break;
535 }
536 }
537 else
538 fprintf (stderr, "display_onmessage: error: unexpected message type '%d'\n",
539 type);
540}
541
542void display_onopen (int fd)
543{
544 char * interface = bview_interface_json();
546 int status = 0;
547
548 if (pid() == 0)
549 if (ws_sendframe_txt (fd, interface, false) < 0 ||
550 ws_sendframe_txt (fd, controls, false) < 0 ||
552 status = -1;
553 free (interface);
554 free (controls);
555
556#if _MPI
558#endif
559
560 debug ("open %d status %d\n", fd, status);
561
562 if (status < 0) {
563 display_onclose (fd);
564 if (pid() == 0)
565 close (fd);
566 }
567}
568
569static void display_update (int i)
570{
571 for (khiter_t k = kh_begin (Display.objects); k != kh_end (Display.objects);
572 ++k)
573 if (kh_exist (Display.objects, k)) {
574 DisplayClient * client = kh_value (Display.objects, k);
575 while (client->fd >= 0) {
576 if (client->iter < i)
577 break;
578 client++;
579 }
580 if (client->fd >= 0) { // at least one client needs update
581 const char * command = kh_key (Display.objects, k);
583 client = kh_value (Display.objects, k);
584 while (client->fd >= 0) {
585 if (client->iter < i) {
586 client->iter = i;
587 if (display_send (command, client->fd) < 0) {
588 debug ("error sending '%s' to '%d'\n", command, client->fd);
589 if (pid() == 0)
590 close (client->fd);
592 if (!kh_exist (Display.objects, k))
593 break;
594 }
595 else
596 client++;
597 }
598 else
599 client++;
600 }
602 if (Display.error && kh_exist (Display.objects, k)) {
603 free ((void *) kh_key (Display.objects, k));
604 free ((void *) kh_value (Display.objects, k));
605 kh_del (strhash, Display.objects, k);
606 }
607 }
608 }
609}
610
611#if _MPI
612static Array * pack_messages (struct ws_message * messages)
613{
614 struct ws_message * msg = messages;
615 Array * packed = array_new();
616 while (msg && msg->fd >= 0) {
617 array_append (packed, msg, sizeof(struct ws_message));
618 array_append (packed, msg->msg, msg->size);
619 msg++;
620 }
621 return packed;
622}
623
624static struct ws_message * unpack_messages (Array * packed)
625{
626 Array * array = array_new();
627 char * p = packed->p;
628 while (p - (char *) packed->p < packed->len) {
629 struct ws_message * msg = (struct ws_message *) p;
630 msg->msg = sysmalloc (msg->size + 1);
631 p += sizeof(struct ws_message);
632 memcpy (msg->msg, p, msg->size);
633 msg->msg[msg->size] = '\0';
634 array_append (array, msg, sizeof(struct ws_message));
635 p += msg->size;
636 }
637 struct ws_message msg = {-1};
638 array_append (array, &msg, sizeof(struct ws_message));
639 struct ws_message * messages = array->p;
640 free (array);
641 return messages;
642}
643#endif
644
646{
647 struct ws_message * messages = NULL, * msg;
648 int nmsg = 0;
649 if (pid() == 0) {
651 while (msg && msg->fd >= 0) msg++, nmsg++;
652 }
653
654#if _MPI
655 Array * packed;
656 if (pid() == 0)
658 else
659 packed = calloc (1, sizeof (Array));
660 MPI_Bcast (&packed->len, 1, MPI_LONG, 0, MPI_COMM_WORLD);
661 if (packed->len > 0) {
662 if (pid() != 0)
663 packed->p = malloc (packed->len);
665 if (pid() != 0)
667 }
669#endif
670
671 msg = messages;
672 nmsg = 0;
673 while (msg && msg->fd >= 0) {
674 switch (msg->type) {
675
676 case WS_FR_OP_OPEN:
677 display_onopen (msg->fd); break;
678
679 case WS_FR_OP_TXT: case WS_FR_OP_BIN:
680 display_onmessage (msg->fd, msg->msg, msg->size, msg->type);
681 break;
682
683 case WS_FR_OP_CLSE:
684 display_onclose (msg->fd); break;
685
686 default:
687 assert (false);
688
689 }
690 sysfree (msg->msg);
691 msg++, nmsg++;
692 }
694 return nmsg;
695}
696
698{
699 char hostname[1024];
700 hostname[1023] = '\0';
701 gethostname (hostname, 1023);
702 struct hostent * h = gethostbyname (hostname);
703 if (!h)
705 "src/display.h:%d: warning: gethostbyname(\"%s\") returned NULL\n",
707 fprintf (fp, DISPLAY_JS "?ws://%s:%d", h ? h->h_name : "127.0.0.1",
708 Display.port);
709}
710
711int display_usage = 20; // use 20% of runtime, maximum
712
713#ifdef DISPLAY
714 // negative: pause, zero: play, positive: step
715int display_play = DISPLAY < 0 ? -1 : 0;
716#else
717int display_play = -1;
718#endif
719
720static void display_destroy()
721{
722 for (khiter_t k = kh_begin (Display.objects); k != kh_end (Display.objects);
723 ++k)
724 if (kh_exist (Display.objects, k)) {
725 free ((void *) kh_key (Display.objects, k));
726 free ((void *) kh_value (Display.objects, k));
727 }
728 kh_destroy (strhash, Display.objects);
729
730 DisplayControl * d = Display.controls->p;
731 for (int i = 0; i < Display.controls->len/sizeof(DisplayControl); i++, d++) {
732 free (d->name);
733 free (d->tooltip);
734 }
735 array_free (Display.controls);
736
737 if (pid() == 0) {
738 remove ("display.html");
739
740 // fixme: close connection cleanly
741 close (Display.sock);
742 }
743}
744
747{
748 if (pid() == 0) {
749 const char * port = DISPLAY_RANGE;
750 if (!strchr (port, ':'))
751 Display.sock = ws_socket_open (atoi (port));
752 else {
753 char * s = strdup (port);
754 char * s1 = strchr (s, ':');
755 *s1++ = '\0';
756 int pmax = atoi(s1);
757 Display.port = atoi(s);
758 while ((Display.sock = ws_socket_open (Display.port)) < 0 &&
759 Display.port < pmax)
760 Display.port++;
761 free (s);
762 }
763 if (Display.sock < 0) {
764 char s[80];
765 sprintf (s, "display(): could not open port '%s'", port);
766 perror (s);
767 exit (1);
768 }
769
770 FILE * fp = fopen ("display.html", "w");
771 if (!fp)
772 perror ("display.html");
773 else {
774 fputs ("<head><meta http-equiv=\"refresh\" content=\"0;URL=", fp);
775 display_url (fp);
776 fputs ("\"></head>\n", fp);
777 fclose (fp);
778 }
779 }
780
781 Display.objects = kh_init (strhash);
782 Display.controls = array_new();
783
785
786#ifndef DISPLAY
787 display_control (display_play, -1, 1, "Run/Pause");
788#elif DISPLAY != 0
789 display_control (display_play, -1, 1, "Run/Pause");
790#endif
791#ifndef DISPLAY_NO_CONTROLS
792 display_control (display_usage, 0, 50, "Display %",
793 "maximum % of runtime used by display");
794#endif
795}
796
797/** @brief Event: refresh_display (i++, last) */
799{
800 do {
801 if (display_play)
803 if (display_poll (display_play ? - 1 : 0))
805 } while (display_play < 0);
806
807 static timer global_timer = {0};
808 static double poll_elapsed = 0.;
809 int refresh = (poll_elapsed <=
811#if _MPI
813#endif
814 if (refresh) {
818 }
819}
const vector a
Definition all-mach.h:59
Array * array_new()
Definition array.h:10
void array_free(Array *a)
Definition array.h:18
void * array_append(Array *a, void *elem, size_t size)
Definition array.h:24
scalar h[]
Definition atmosphere.h:6
int npe
Definition balance.h:195
if TRASH define &&NewPid *& val(newpid, 0, 0, 0)) -> pid > 0) @else @ define is_newpid()(((NewPid *)&val(newpid, 0, 0, 0)) ->pid > 0) @endif Array *linear_tree(size_t size, scalar newpid)
Definition balance.h:13
int min
Definition balance.h:192
#define DISPLAY
Definition bview.c:21
define k
free(list1)
int x
Definition common.h:76
void free_solver_func_add(free_solver_func func)
Definition common.h:512
double timer_elapsed(timer t)
Definition common.h:379
static char * display_defaults
Definition common.h:521
static number clamp(number x, number a, number b)
Definition common.h:15
timer timer_start(void)
Definition common.h:368
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line line line type
Definition config.h:571
#define p
Definition tree.h:320
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line calloc(n, s) @ define prealloc(p
#define assert(a)
Definition config.h:107
#define LINENO
Definition config.h:105
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line realloc(p, s) @ define pfree(p
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line line line strdup(s) @ define tracing(...) @ define end_tracing(...) @define tid() 0 @define pid() 0 @define npe() 1 @define mpi_all_reduce(v
static int display_send(const char *command, int fd)
Definition display.h:248
void display_url(FILE *fp)
Definition display.h:697
static void display_command(const char *command)
Definition display.h:140
void event_refresh_display(void)
Event: refresh_display (i++, last)
Definition display.h:798
static void display_display()
Definition display.h:106
kh_strhash_t * objects
Definition display.h:99
init_solver void display_init()
Definition display.h:746
#define JSON_BUILD(...)
Definition display.h:348
char * error
Definition display.h:101
static void display_control_update(const char *command, int fd)
Definition display.h:465
int sock
Definition display.h:100
int port
Definition display.h:100
static DisplayControl * display_control_lookup(const char *name)
Definition display.h:426
static struct @3 Display
void display_onopen(int fd)
Definition display.h:542
#define debug(...)
Definition display.h:82
static void display_destroy()
Definition display.h:720
static char * read_file_into_buffer(FILE *fp)
Definition display.h:123
int display_usage
Definition display.h:711
Array * controls
Definition display.h:102
void display_onmessage(int fd, const char *msg, size_t size, int type)
Definition display.h:522
#define display_control(val,...)
Definition display.h:462
void display_onclose(int fd)
Definition display.h:512
#define DISPLAY_JS
Definition display.h:71
static int ws_send_array(int fd, Array *a, int status, int type, unsigned int *shift)
Definition display.h:205
int display_poll(int timeout)
Definition display.h:645
static void display_update(int i)
Definition display.h:569
static char * display_control_json()
Definition display.h:396
int display_play
Definition display.h:715
static void array_remove(khiter_t k, int fd)
Definition display.h:351
static char * bview_interface_json()
Definition display.h:493
void display_control_internal(void *ptr, double min, double max, char *name=NULL, char *tooltip=NULL, char *ptr_name, int size)
Definition display.h:435
static void display_add(const char *command, int fd)
Definition display.h:320
#define DISPLAY_RANGE
Definition display.h:79
static void display_remove(const char *command, int fd)
Definition display.h:377
struct @2 bview_interface[]
int(* json)(char *s, int len)
Definition draw.h:1815
scalar s
Definition embed-tree.h:56
*cs[i, 0, 0] a *[i -1, 0, 0] j
Definition embed.h:88
double d[2]
Definition embed.h:383
scalar int i
Definition embed.h:74
static int line
Definition errors.c:754
void init_solver()
Definition events.h:298
static int
Definition include.c:978
#define kh_exist(h, x)
Definition khash.h:490
#define kh_del(name, h, k)
Definition khash.h:482
#define KHASH_MAP_INIT_STR(name, khval_t)
Definition khash.h:615
#define kh_init(name)
Definition khash.h:430
khint_t khiter_t
Definition khash.h:154
#define kh_key(h, x)
Definition khash.h:498
#define kh_get(name, h, k)
Definition khash.h:474
#define kh_destroy(name, h)
Definition khash.h:437
#define kh_begin(h)
Definition khash.h:519
#define khash_t(name)
Definition khash.h:423
#define kh_value(h, x)
Definition khash.h:512
#define kh_end(h)
Definition khash.h:526
#define kh_put(name, h, k, r)
Definition khash.h:465
size_t max
Definition mtrace.h:8
FILE * fp
Definition mtrace.h:7
size_t size
size *double * buf
Definition array.h:5
char * name
Definition display.h:390
define n sizeof(Cell))) @define fine(a
void vertex_buffer_free()
void vertex_buffer_setup()
struct @17 VertexBuffer
bool process_line(char *line)
Definition view.h:645
scalar c
Definition vof.h:57
src vof fflush(ferr)