Thanks to visit codestin.com
Credit goes to glvis.github.io

GLVis  v4.2
Accurate and flexible finite element visualization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
vsdata.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-443271.
4 //
5 // This file is part of the GLVis visualization tool and library. For more
6 // information and source code availability see https://glvis.org.
7 //
8 // GLVis is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include <cstdlib>
13 #include <cmath>
14 
15 #include <iomanip>
16 #include <sstream>
17 #include <limits>
18 using namespace std;
19 
20 #include "vsdata.hpp"
21 #include "aux_vis.hpp"
22 #include "material.hpp"
23 #include "palettes.hpp"
24 
25 #include "gl2ps.h"
26 
27 const char *strings_off_on[] = { "off", "on" };
28 
30 {
31  double am = fabs(minv);
32  if (am < fabs(maxv)) { am = fabs(maxv); }
33  if (float(am) < 100*numeric_limits<float>::min()) { am = 1e-3; }
34  if ((maxv-minv) < 1e-5*am)
35  {
36  // Shading quality may be bad since OpenGL uses single precision. We
37  // should probably pre-scale the solution before feeding it to OpenGL
38  int old_prec = cout.precision(12);
39  cout << "[minv,maxv] = " << "[" << minv << "," << maxv
40  << "] (maxv-minv = " << maxv-minv << ")\n --> ";
41  minv -= 0.49999e-5*am;
42  maxv += 0.50001e-5*am;
43  cout << "[" << minv << "," << maxv << "]" << endl;
44  cout.precision(old_prec);
45  }
46 }
47 
49 {
50  if (autoscale == 1)
51  {
52  FindNewBoxAndValueRange(prepare);
53  }
54  else if (autoscale == 2)
55  {
56  FindNewValueRange(prepare);
57  }
58  else if (autoscale == 3)
59  {
60  FindMeshBox(prepare);
61  }
62 }
63 
65 {
66  if (autoscale == 1 || autoscale == 3)
67  {
68  FindNewBoxAndValueRange(prepare);
69  }
70  else
71  {
72  FindNewValueRange(prepare);
73  }
74 }
75 
77  glm::mat4 xfrm)
78 {
79  const int n = 8;
80  const double step = 2*M_PI/n;
81  const double nz = (1.0/4.0);
82  double point = step;
83  int i;
84 
85  glm::mat3 normXfrm = glm::inverseTranspose(glm::mat3(xfrm));
86 
87  glm::vec3 start1vtx = glm::vec3(xfrm * glm::vec4(0, 0, 0, 1));
88  glm::vec3 start1norm = glm::vec3(normXfrm * glm::vec3(0, 0, 1));
89  glm::vec3 start2vtx = glm::vec3(xfrm * glm::vec4(1, 0, -4, 1));
90  glm::vec3 start2norm = glm::vec3(normXfrm * glm::vec3(1, 0, nz));
91 
92  builder.glBegin(GL_TRIANGLE_FAN);
93  builder.glNormal3d(start1norm[0], start1norm[1], start1norm[2]);
94  builder.glVertex3d(start1vtx[0], start1vtx[1], start1vtx[2]);
95  builder.glNormal3d(start2norm[0], start2norm[1], start2norm[2]);
96  builder.glVertex3d(start2vtx[0], start2vtx[1], start2vtx[2]);
97 
98  for (i = 1; i < n; i++)
99  {
100  glm::vec3 baseVtx = glm::vec3(xfrm * glm::vec4(cos(point), sin(point), -4, 1));
101  glm::vec3 baseNorm = glm::vec3(normXfrm * glm::vec3(cos(point), sin(point),
102  nz));
103  builder.glNormal3d(baseNorm[0], baseNorm[1], baseNorm[2]);
104  builder.glVertex3d(baseVtx[0], baseVtx[1], baseVtx[2]);
105  point += step;
106  }
107  builder.glNormal3d(start2norm[0], start2norm[1], start2norm[2]);
108  builder.glVertex3d(start2vtx[0], start2vtx[1], start2vtx[2]);
109  builder.glEnd();
110 }
111 
112 // Draw an arrow starting at point (px, py, pz) with orientation (vx, vy, vz)
113 // and length "length".
115  double px, double py, double pz,
116  double vx, double vy, double vz,
117  double length,
118  double cone_scale)
119 {
120  double xc = 0.5*(bb.x[0]+bb.x[1]);
121  double yc = 0.5*(bb.y[0]+bb.y[1]);
122  double zc = 0.5*(bb.z[0]+bb.z[1]);
123  glm::mat4 xfrm(1.0);
124  xfrm = glm::translate(xfrm, glm::vec3(xc, yc, zc));
125  xfrm = glm::scale(xfrm, glm::vec3(1.0/xscale, 1.0/yscale, 1.0/zscale));
126 
127  double rlen = length/sqrt(vx*vx+vy*vy+vz*vz);
128  px = (px-xc)*xscale;
129  py = (py-yc)*yscale;
130  pz = (pz-zc)*zscale;
131  vx *= rlen*xscale;
132  vy *= rlen*yscale;
133  vz *= rlen*zscale;
134 
135  if (arrow_scaling_type == 0)
136  {
137  length = sqrt(vx*vx+vy*vy+vz*vz);
138  }
139 
140  xfrm = glm::translate(xfrm, glm::vec3(px, py, pz));
141 
142  double rhos = sqrt (vx*vx+vy*vy+vz*vz);
143  float phi = acos(vz/rhos);
144  float theta;
145  theta = atan2 (vy, vx);
146 
147  xfrm = glm::rotate(xfrm, theta, glm::vec3(0, 0, 1));
148  xfrm = glm::rotate(xfrm, phi, glm::vec3(0, 1, 0));
149 
150  xfrm = glm::scale(xfrm, glm::vec3(length));
151 
152 
153  if (arrow_type == 1)
154  {
155  xfrm = glm::translate(xfrm, glm::vec3(0, 0, -0.5));
156  }
157 
158  glm::vec4 pt1 = xfrm * glm::vec4(0, 0, 0, 1);
159  glm::vec4 pt2 = xfrm * glm::vec4(0, 0, 1, 1);
160 
161  builder.glBegin(GL_LINES);
162  builder.glVertex3d(pt1[0], pt1[1], pt1[2]);
163  builder.glVertex3d(pt2[0], pt2[1], pt2[2]);
164  builder.glEnd();
165 
166  xfrm = glm::translate(xfrm, glm::vec3(0, 0, 1));
167  xfrm = glm::scale(xfrm, glm::vec3(cone_scale));
168 
169  Cone(builder, xfrm);
170 }
171 
173  double px, double py, double pz,
174  double vx, double vy, double vz,
175  double length,
176  double cone_scale)
177 {
178  glm::mat4 xfrm(1.0);
179  xfrm = glm::translate(xfrm, glm::vec3(px, py, pz));
180 
181  double rhos = sqrt (vx*vx+vy*vy+vz*vz);
182  float phi = acos(vz/rhos);
183  float theta;
184  theta = atan2 (vy, vx);
185 
186  xfrm = glm::rotate(xfrm, theta, glm::vec3(0, 0, 1));
187  xfrm = glm::rotate(xfrm, phi, glm::vec3(0, 1, 0));
188 
189  xfrm = glm::scale(xfrm, glm::vec3(length));
190 
191  glm::vec4 pt1 = xfrm * glm::vec4(0, 0, 0, 1);
192  glm::vec4 pt2 = xfrm * glm::vec4(0, 0, 1, 1);
193 
194  builder.glBegin(GL_LINES);
195  builder.glVertex3d(pt1[0], pt1[1], pt1[2]);
196  builder.glVertex3d(pt2[0], pt2[1], pt2[2]);
197  builder.glEnd();
198 
199  xfrm = glm::translate(xfrm, glm::vec3(0, 0, 1));
200  xfrm = glm::scale(xfrm, glm::vec3(cone_scale));
201 
202  Cone(builder, xfrm);
203 }
204 
206  double px, double py, double pz,
207  double vx, double vy, double vz,
208  double length,
209  double cone_scale)
210 {
211  double rhos = sqrt (vx*vx+vy*vy+vz*vz);
212  if (rhos == 0.0)
213  {
214  return;
215  }
216  double phi = acos(vz/rhos), theta = atan2(vy, vx);
217  const int n = 8;
218  const double step = 2*M_PI/n, nz = (1.0/4.0);
219  double point = step, cone[n+4][3], normal[n+2][3];
220  int i, j, k;
221 
222  cone[0][0] = 0; cone[0][1] = 0; cone[0][2] = 1;
223  cone[1][0] = cone_scale; cone[1][1] = 0; cone[1][2] = -4*cone_scale + 1;
224  normal[0][0] = 0.0/cone_scale;
225  normal[0][1] = 0.0/cone_scale;
226  normal[0][2] = 1.0/cone_scale;
227  normal[1][0] = 1.0/cone_scale;
228  normal[1][1] = 0.0/cone_scale;
229  normal[1][2] = nz/cone_scale;
230 
231  for (i=2; i<n+1; i++)
232  {
233  normal[i][0] = cos(point)/cone_scale;
234  normal[i][1] = sin(point)/cone_scale;
235  normal[i][2] = nz/cone_scale;
236 
237  cone[i][0] = cos(point)*cone_scale;
238  cone[i][1] = sin(point)*cone_scale;
239  cone[i][2] = -4*cone_scale + 1;
240  point += step;
241  }
242  cone[n+1][0] = cone_scale; cone[n+1][1] = 0; cone[n+1][2] =-4*cone_scale + 1;
243  normal[n+1][0] = 1.0/cone_scale;
244  normal[n+1][1] = 0.0/cone_scale;
245  normal[n+1][2] = nz/cone_scale;
246 
247  cone[n+2][0] = 0; cone[n+2][1] = 0; cone[n+2][2] = 0;
248  cone[n+3][0] = 0; cone[n+3][1] = 0; cone[n+3][2] = 1;
249 
250  if (arrow_scaling_type == 0)
251  {
252  length = rhos;
253  }
254 
255  // double xc = 0.5*(x[0]+x[1]), yc = 0.5*(y[0]+y[1]), zc = 0.5*(z[0]+z[1]);
256  double coord[3];
257  // double rlen = length/rhos;
258 
259  // px = (px-xc)*xscale; py = (py-yc)*yscale; pz = (pz-zc)*zscale;
260  // vx *= rlen*xscale; vy *= rlen*yscale; vz *= rlen*zscale;
261 
262  if (arrow_type == 1)
263  for (i=0; i<n+4; i++)
264  {
265  cone[i][2] -= 0.5;
266  }
267 
268  double M[3][3]= {{cos(theta)*cos(phi), -sin(theta), cos(theta)*sin(phi)},
269  {sin(theta)*cos(phi), cos(theta), sin(theta)*sin(phi)},
270  { -sin(phi), 0., cos(phi)}
271  };
272  double v[3] = { M[0][2]/xscale, M[1][2]/yscale, M[2][2]/zscale };
273  length /= sqrt(v[0]*v[0]+v[1]*v[1]+v[2]*v[2]);
274 
275  for (i=0; i<n+4; i++)
276  {
277  for (j=0; j<3; j++)
278  {
279  coord[j] = cone[i][j] * length;
280  }
281 
282  for (k=0; k<3; k++)
283  {
284  cone[i][k] = 0.;
285  for (j=0; j<3; j++)
286  {
287  cone[i][k] += M[k][j] * coord[j];
288  }
289  }
290  // cone[i][0] = (cone[i][0] + px)/xscale + xc;
291  // cone[i][1] = (cone[i][1] + py)/yscale + yc;
292  // cone[i][2] = (cone[i][2] + pz)/zscale + zc;
293  cone[i][0] = cone[i][0]/xscale + px;
294  cone[i][1] = cone[i][1]/yscale + py;
295  cone[i][2] = cone[i][2]/zscale + pz;
296  }
297 
298  for (i=0; i<=n+1; i++)
299  {
300  for (j=0; j<3; j++)
301  {
302  coord[j] = normal[i][j];
303  }
304 
305  for (k=0; k<3; k++)
306  {
307  normal[i][k] = 0.;
308  for (j=0; j<3; j++)
309  {
310  normal[i][k] += M[k][j] * coord[j];
311  }
312  }
313  normal[i][0] *= xscale;
314  normal[i][1] *= yscale;
315  normal[i][2] *= zscale;
316  }
317 
318  builder.glBegin(GL_TRIANGLE_FAN);
319  for (i=0; i<=n+1; i++)
320  {
321  Normalize(normal[i]);
322  builder.glNormal3dv(normal[i]);
323  builder.glVertex3dv(cone[i]);
324  }
325  builder.glEnd();
326 
327  builder.glBegin(GL_LINES);
328  builder.glVertex3dv(cone[n+2]);
329  builder.glVertex3dv(cone[n+3]);
330  builder.glEnd();
331 }
332 
334  double maxval,
335  Array<double> *mesh_level,
336  Array<double> *lsurf_levels)
337 {
338 
339  int i;
340 
341  float miny;
342  float maxy;
343  float minx;
344  float maxx;
345  float posz = -4.0;
346 
347  if (OrthogonalProjection)
348  {
349  miny = -.65;
350  maxy = .65;
351  minx = 0.73;
352  maxx = 0.80;
353  }
354  else
355  {
356  miny = -1.;
357  maxy = 1.;
358  minx = 1.2;
359  maxx = 1.3;
360  }
361  color_bar.clear();
362  color_bar.addQuad<gl3::VertexTex>(
363  {{minx, miny, posz},{0.f, 0.f}},
364  {{maxx, miny, posz},{0.f, 0.f}},
365  {{maxx, maxy, posz},{1.f, 0.f}},
366  {{minx, maxy, posz},{1.f, 0.f}}
367  );
368 
369  static const int border = 2;
370 
371  if (border == 1)
372  {
373  color_bar.addLines<gl3::Vertex>(
374  {
375  {minx, miny, posz}, {maxx, miny, posz},
376  {maxx, miny, posz}, {maxx, maxy, posz},
377  {maxx, maxy, posz}, {minx, maxy, posz},
378  {minx, maxy, posz}, {minx, miny, posz}
379  });
380  }
381  else if (border == 2)
382  {
383  color_bar.addLine<gl3::Vertex>({minx, miny, posz}, {minx, maxy, posz});
384  color_bar.addLine<gl3::Vertex>({maxx, miny, posz}, {maxx, maxy, posz});
385  }
386 
387  if (lsurf_levels)
388  {
389  for (i = 0; i < lsurf_levels->Size(); i++)
390  {
391  float Y = miny + (maxy - miny) * LogUVal((*lsurf_levels)[i]);
392  color_bar.addLine<gl3::Vertex>({minx, Y, posz}, {maxx, Y, posz});
393  }
394  }
395  if (mesh_level)
396  {
397  for (i = 0; i < mesh_level->Size(); i++)
398  {
399  float Y = miny + (maxy - miny) * LogUVal((*mesh_level)[i]);
400  color_bar.addLine<gl3::Vertex>({minx, Y, posz}, {maxx, Y, posz});
401  }
402  }
403 
404  const double text_x = maxx + 0.4*(maxx-minx);
405  double val;
406  double Y;
407  if (!mesh_level)
408  {
409  for (i = 0; i <= 4; i++)
410  {
411  Y = miny + i * (maxy-miny) / 4;
412 
413  val = ULogVal(i / 4.0);
414 
415  ostringstream buf;
416  buf << setprecision(4) << val;
417  color_bar.addText(text_x,Y,posz, buf.str());
418  }
419  }
420  else
421  {
422  for (i = 0; i < mesh_level->Size(); i++)
423  {
424  val = (*mesh_level)[i];
425  Y = miny + (maxy - miny) * LogUVal(val);
426 
427  ostringstream buf;
428  buf << setprecision(4) << val;
429  color_bar.addText(text_x,Y,posz, buf.str());
430  }
431  }
432 
433  if (lsurf_levels)
434  {
435  for (i = 0; i < lsurf_levels->Size(); i++)
436  {
437  val = (*lsurf_levels)[i];
438  Y = miny + (maxy - miny) * LogUVal(val);
439 
440  ostringstream buf;
441  buf << setprecision(4) << val;
442  color_bar.addText(text_x,Y,posz, buf.str());
443  }
444  }
445  updated_bufs.emplace_back(&color_bar);
446 }
447 
448 // Draw a centered caption at the top (visible with the colorbar)
450 {
451  bool empty = plot_caption.empty();
452  colorbar = (colorbar ? empty+1 : !empty);
453 
454  string caption(plot_caption);
455  if (!extra_caption.empty())
456  {
457  caption += " (" + extra_caption + ")";
458  }
459 
460  caption_buf.clear();
461  caption_buf.addText(0, 0, 0, caption);
462  updated_bufs.emplace_back(&caption_buf);
463  GetFont()->getObjectSize(caption, caption_w, caption_h);
464 }
465 
467 extern thread_local VisualizationScene * locscene;
468 
470 {
472  SendExposeEvent();
473 }
474 
476 {
477  cout << "Enter new caption: " << flush;
478  std::getline(cin, plot_caption);
479  vsdata->PrepareCaption(); // turn on or off the caption
480  SendExposeEvent();
481 }
482 
484 {
485  vsdata -> ToggleScaling();
486  SendExposeEvent();
487 }
488 
490 {
491  vsdata -> ToggleDrawAxes();
492  SendExposeEvent();
493 }
494 
495 void Key_Mod_a_Pressed(GLenum state)
496 {
497  if (state & KMOD_CTRL)
498  {
499  static const char *autoscale_modes[] = { "off", "on", "value", "mesh" };
500  int autoscale = vsdata->GetAutoscale();
501  autoscale = (autoscale + 1)%4;
502  cout << "Autoscale: " << flush;
503  vsdata->SetAutoscale(autoscale);
504  cout << autoscale_modes[autoscale] << endl;
505  SendExposeEvent();
506  }
507  else
508  {
510  SendExposeEvent();
511  }
512 }
513 
515 {
516  cout << vsdata->GetHelpString() << flush;
517 }
518 
520 {
521  vsdata -> ToggleLight();
522  SendExposeEvent();
523 }
524 
526 {
527  vsdata->ToggleLogscale(true);
528  SendExposeEvent();
529 }
530 
532 {
533  locscene -> spinning = 0;
535  vsdata -> CenterObject();
536 
537  locscene -> ViewAngle = 45.0;
538  locscene -> ViewScale = 1.0;
539  locscene -> ViewCenterX = 0.0;
540  locscene -> ViewCenterY = 0.0;
541  locscene->cam.Reset();
542  vsdata -> key_r_state = 0;
543  SendExposeEvent();
544 }
545 
547 {
548  locscene->spinning = 0;
550  vsdata->Toggle2DView();
551  SendExposeEvent();
552 }
553 
554 void KeypPressed(GLenum state)
555 {
556  if (state & KMOD_CTRL)
557  {
558  KeyCtrlP();
559  }
560  else
561  {
563  SendExposeEvent();
564  }
565 }
566 
568 {
570  SendExposeEvent();
571 }
572 
573 static void KeyF5Pressed()
574 {
575  int n;
576  double min, max;
577 
578  cout << "Enter min : " << flush;
579  cin >> min;
580  cout << "Enter max : " << flush;
581  cin >> max;
582  cout << "Enter n : " << flush;
583  cin >> n;
584 
585  vsdata -> SetLevelLines (min, max, n, 0);
586 
587  vsdata -> UpdateLevelLines();
588  SendExposeEvent();
589 }
590 
592 {
593  int RepeatPaletteTimes = vsdata->palette.GetRepeatTimes();
594  cout << "Palette is repeated " << RepeatPaletteTimes << " times.\n"
595  << "(Negative value means the palette is flipped.)\n"
596  << "Enter new value: " << flush;
597  cin >> RepeatPaletteTimes;
598  if (RepeatPaletteTimes == 0)
599  {
600  RepeatPaletteTimes = 1;
601  }
602  cout << "Palette will be repeated " << RepeatPaletteTimes
603  << " times now.\n\n";
604  vsdata->palette.SetRepeatTimes(RepeatPaletteTimes);
605 
606  int pal = vsdata->palette.ChoosePalette();
607 
608  int colors_used = vsdata->palette.GetNumColors(pal);
609  int palette_size = vsdata->palette.GetSize(pal);
610  cout << "\nPalette is using " << colors_used << " colors.\n"
611  << "Enter new value (0 = use original " << palette_size
612  << " colors): " << flush;
613  cin >> colors_used;
614  if (colors_used == 1) { colors_used = 0; }
615  vsdata->palette.SetNumColors(colors_used);
616 
617  vsdata->palette.Init();
618  vsdata->palette.SetIndex(pal);
619 
620  colors_used = vsdata->palette.GetNumColors();
621  cout << "Palette will be using " << colors_used << " colors now.\n";
622 
624  SendExposeEvent();
625 }
626 
627 void KeyF7Pressed(GLenum state)
628 {
629  if (state & KMOD_SHIFT)
630  {
631  cout << "Current bounding box:\n"
632  << " min: (" << vsdata->bb.x[0] << ',' << vsdata->bb.y[0] << ','
633  << vsdata->bb.z[0] << ")\n"
634  << " max: (" << vsdata->bb.x[1] << ',' << vsdata->bb.y[1] << ','
635  << vsdata->bb.z[1] << ")\n"
636  << "Enter new bounding box:\n"
637  << "x_min = " << flush;
638  cin >> vsdata->bb.x[0];
639  cout << "y_min = " << flush;
640  cin >> vsdata->bb.y[0];
641  cout << "z_min = " << flush;
642  cin >> vsdata->bb.z[0];
643  cout << "x_max = " << flush;
644  cin >> vsdata->bb.x[1];
645  cout << "y_max = " << flush;
646  cin >> vsdata->bb.y[1];
647  cout << "z_max = " << flush;
648  cin >> vsdata->bb.z[1];
649  cout << "New bounding box:\n"
650  << " min: (" << vsdata->bb.x[0] << ',' << vsdata->bb.y[0] << ','
651  << vsdata->bb.z[0] << ")\n"
652  << " max: (" << vsdata->bb.x[1] << ',' << vsdata->bb.y[1] << ','
653  << vsdata->bb.z[1] << ")\n" << flush;
655  SendExposeEvent();
656  }
657  else
658  {
659  cout << "[minv,maxv] = [" << vsdata->GetMinV() << "," << vsdata->GetMaxV()
660  << "] maxv-minv = " << vsdata->GetMaxV()-vsdata->GetMinV() << "\n"
661  << "New value for minv: " << flush;
662  cin >> vsdata->GetMinV();
663  cout << "New value for maxv: " << flush;
664  cin >> vsdata->GetMaxV();
665  vsdata->UpdateValueRange(true);
666  SendExposeEvent();
667  }
668 }
669 
671 {
672  float x, y, z, w;
673 
674  cout << "Enter light source position\n(0,0,1,w) - from camera\n"
675  "(0,1,0,w) - from above\n(1,0,0,w) - from the right\n"
676  "w = 0/1 defines directional/spot light\n";
677  cout << "x = " << flush;
678  cin >> x;
679  cout << "y = " << flush;
680  cin >> y;
681  cout << "z = " << flush;
682  cin >> z;
683  cout << "w = " << flush;
684  cin >> w;
685 
686  vsdata->SetLight0CustomPos({x, y, z, w});
687  SendExposeEvent();
688 }
689 
691 {
692  int ml;
693 
694  ml = (vsdata->GetLightMatIdx() + 1) % 5;
695  vsdata->SetLightMatIdx(ml);
696  SendExposeEvent();
697  cout << "New material/light : " << ml << endl;
698 }
699 
701 {
703  vsdata->PrepareAxes();
705  SendExposeEvent();
706 }
707 
709 {
710  vsdata->glTF_Export();
711 }
712 
714 {
715  vsdata->PrintState();
716 }
717 
719 {
720  vsdata -> EventUpdateColors();
721  vsdata -> PrepareLines();
722  // vsdata->CPPrepare();
723  SendExposeEvent();
724 }
725 
727 {
728  locscene->matAlpha -= 0.05;
729  if (locscene->matAlpha < 0.0)
730  {
731  locscene->matAlpha = 0.0;
732  }
734  SendExposeEvent();
735 }
736 
738 {
739  locscene->matAlpha += 0.05;
740  if (locscene->matAlpha > 1.0)
741  {
742  locscene->matAlpha = 1.0;
743  }
745  SendExposeEvent();
746 }
747 
749 {
750  bool curr_aa = GetAppWindow()->getRenderer().getAntialiasing();
752 
753  cout << "Multisampling/Antialiasing: "
754  << strings_off_on[!curr_aa ? 1 : 0] << endl;
755 
756  // vsdata -> EventUpdateColors();
757  SendExposeEvent();
758 }
759 
761 {
762  locscene->matAlphaCenter -= 0.25;
763  // vsdata -> EventUpdateColors();
765  SendExposeEvent();
766 #ifdef GLVIS_DEBUG
767  cout << "MatAlphaCenter = " << locscene->matAlphaCenter << endl;
768 #endif
769 }
770 
772 {
773  locscene->matAlphaCenter += 0.25;
774  // vsdata -> EventUpdateColors();
776  SendExposeEvent();
777 #ifdef GLVIS_DEBUG
778  cout << "MatAlphaCenter = " << locscene->matAlphaCenter << endl;
779 #endif
780 }
781 
783 {
784  vsdata->ToggleRuler();
785  SendExposeEvent();
786 }
787 
789 {
791  SendExposeEvent();
792 }
793 
795 {
797  SendExposeEvent();
798 }
799 
801 {
802  if (warn)
803  {
804  cout << "The range [" << minv << ',' << maxv
805  << "] is not appropriate for logarithmic scale!" << endl;
806  }
807  cout << "Logarithmic scale: " << strings_off_on[logscale ? 1 : 0]
808  << endl;
809 }
810 
812 {
813  if (logscale || LogscaleRange())
814  {
815  logscale = !logscale;
816  palette.SetUseLogscale(logscale);
817  SetLogA();
818  SetLevelLines(minv, maxv, nl);
819  UpdateLevelLines();
820  EventUpdateColors();
821  if (print)
822  {
823  PrintLogscale(false);
824  }
825  }
826  else if (print)
827  {
828  PrintLogscale(true);
829  }
830  PrepareRuler();
831 }
832 
834 {
835  ruler_on = (ruler_on + 1) % 3;
836  PrepareRuler();
837 }
838 
840 {
841  cout << "Current ruler position: (" << ruler_x << ','
842  << ruler_y << ',' << ruler_z << ")\n";
843  cout << "x = " << flush; cin >> ruler_x;
844  cout << "y = " << flush; cin >> ruler_y;
845  cout << "z = " << flush; cin >> ruler_z;
846  if (ruler_x < bb.x[0])
847  {
848  ruler_x = bb.x[0];
849  }
850  else if (ruler_x > bb.x[1])
851  {
852  ruler_x = bb.x[1];
853  }
854  if (ruler_y < bb.y[0])
855  {
856  ruler_y = bb.y[0];
857  }
858  else if (ruler_y > bb.y[1])
859  {
860  ruler_y = bb.y[1];
861  }
862  if (ruler_z < bb.z[0])
863  {
864  ruler_z = bb.z[0];
865  }
866  else if (ruler_z > bb.z[1])
867  {
868  ruler_z = bb.z[1];
869  }
870  cout << "New ruler position: (" << ruler_x << ','
871  << ruler_y << ',' << ruler_z << ")" << endl;
872  PrepareRuler();
873 }
874 
876 {
877  float pos_z = LogVal(ruler_z, log_z);
878  float x_f[2] = {(float) bb.x[0], (float) bb.x[1]};
879  float y_f[2] = {(float) bb.y[0], (float) bb.y[1]};
880  float z_f[2] = {(float) bb.z[0], (float) bb.z[1]};
881  float ruler_x_f = (float) ruler_x,
882  ruler_y_f = (float) ruler_y;
883  ruler_buf.clear();
884  if (ruler_on == 2)
885  {
886  std::array<uint8_t, 4> color = gl3::ColorU8(0.8, 0.8, 0.8, 1.0);
887  std::array<float, 3> norm = { 0, 0, 1 };
888  ruler_buf.addQuad<gl3::VertexNormColor>(
889  {{x_f[0], y_f[0], pos_z},norm,color},
890  {{x_f[1], y_f[0], pos_z},norm,color},
891  {{x_f[1], y_f[1], pos_z},norm,color},
892  {{x_f[0], y_f[1], pos_z},norm,color}
893  );
894 
895  std::array<float, 3> norm_2 = { 0, 1, 0 };
896  ruler_buf.addQuad<gl3::VertexNormColor>(
897  {{x_f[0], ruler_y_f, z_f[0]},norm_2,color},
898  {{x_f[0], ruler_y_f, z_f[1]},norm_2,color},
899  {{x_f[1], ruler_y_f, z_f[1]},norm_2,color},
900  {{x_f[1], ruler_y_f, z_f[0]},norm_2,color}
901  );
902 
903  std::array<float, 3> norm_3 = { 1, 0, 0 };
904  ruler_buf.addQuad<gl3::VertexNormColor>(
905  {{ruler_x_f, y_f[0], z_f[0]},norm_3,color},
906  {{ruler_x_f, y_f[1], z_f[0]},norm_3,color},
907  {{ruler_x_f, y_f[1], z_f[1]},norm_3,color},
908  {{ruler_x_f, y_f[0], z_f[1]},norm_3,color}
909  );
910 
911  ruler_buf.addLines<gl3::Vertex>(
912  {
913  {x_f[0], y_f[0], pos_z}, {x_f[1], y_f[0], pos_z},
914  {x_f[1], y_f[0], pos_z}, {x_f[1], y_f[1], pos_z},
915  {x_f[1], y_f[1], pos_z}, {x_f[0], y_f[1], pos_z},
916  {x_f[0], y_f[1], pos_z}, {x_f[0], y_f[0], pos_z},
917 
918  {x_f[0], ruler_y_f, z_f[0]}, {x_f[1], ruler_y_f, z_f[0]},
919  {x_f[1], ruler_y_f, z_f[0]}, {x_f[1], ruler_y_f, z_f[1]},
920  {x_f[1], ruler_y_f, z_f[1]}, {x_f[0], ruler_y_f, z_f[1]},
921  {x_f[0], ruler_y_f, z_f[1]}, {x_f[0], ruler_y_f, z_f[0]},
922 
923  {ruler_x_f, y_f[0], z_f[0]}, {ruler_x_f, y_f[1], z_f[0]},
924  {ruler_x_f, y_f[1], z_f[0]}, {ruler_x_f, y_f[1], z_f[1]},
925  {ruler_x_f, y_f[1], z_f[1]}, {ruler_x_f, y_f[0], z_f[1]},
926  {ruler_x_f, y_f[0], z_f[1]}, {ruler_x_f, y_f[0], z_f[0]}
927  });
928  }
929 
930  ruler_buf.addLines<gl3::Vertex>(
931  {
932  {x_f[0], ruler_y_f, pos_z},
933  {x_f[1], ruler_y_f, pos_z},
934  {ruler_x_f, y_f[0], pos_z},
935  {ruler_x_f, y_f[1], pos_z},
936  {ruler_x_f, ruler_y_f, z_f[0]},
937  {ruler_x_f, ruler_y_f, z_f[1]}
938  });
939 
940  updated_bufs.emplace_back(&ruler_buf);
941 }
942 
944 {
945  gl3::GlMatrix newrot;
946  newrot.identity();
947  translmat = newrot.mtx;
948 
949  switch (key_r_state)
950  {
951  case 0:
952  break;
953 
954  case 1:
955  newrot.rotate(-90.0, 1.0f, 0.0f, 0.0f);
956  break;
957 
958  case 2:
959  newrot.rotate(-90.0, 1.0f, 0.0f, 0.0f);
960  newrot.rotate(-90.0, 0.0f, 0.0f, 1.0f);
961  break;
962 
963  case 3:
964  newrot.rotate(-90.0, 1.0f, 0.0f, 0.0f);
965  newrot.rotate(-180.0, 0.0f, 0.0f, 1.0f);
966  break;
967 
968  case 4:
969  newrot.rotate(-90.0, 1.0f, 0.0f, 0.0f);
970  newrot.rotate(-270.0, 0.0f, 0.0f, 1.0f);
971  break;
972 
973  case 5:
974  newrot.rotate(180.0, 1.0f, 0.0f, 0.0f);
975  break;
976  }
977 
978  // if (locscene -> view != 2) // make 'R' work the same in 2D and 3D
979  key_r_state = (key_r_state+1)%6;
980 
981  rotmat = newrot.mtx;
982 }
983 
985 {
986  int w, h;
987  wnd->getWindowSize(w, h);
988  gl3::SceneInfo scene {};
989  scene.needs_buffering = std::move(updated_bufs);
990  updated_bufs.clear();
991 
992  gl3::RenderParams params {};
993  params.model_view.identity();
994  params.projection.identity();
995  params.mesh_material = BLK_MAT;
996  params.num_pt_lights = 0;
997  params.static_color = this->GetLineColor();
998  params.use_clip_plane = false;
999  params.contains_translucent = false;
1000 
1001  if (colorbar)
1002  {
1003  // add color bar to draw list
1004  params.projection.mtx = proj_mtx;
1005  scene.queue.emplace_back(params, &color_bar);
1006  params.projection.identity();
1007  }
1008  if (colorbar == 1)
1009  {
1010  // caption size is in screen pixels and needs to be centered with
1011  // GL pixel size
1012  int gl_w, gl_h;
1013  wnd->getGLDrawSize(gl_w, gl_h);
1014  // add caption to draw list
1015  double v_pos = 2.;
1016  double line_h = GetFont()->getFontLineSpacing();
1017  params.model_view.translate(-(double)caption_w / gl_w,
1018  1.0 - 2 * v_pos * line_h / gl_h, 0.0);
1019  scene.queue.emplace_back(params, &caption_buf);
1020  }
1021  params.contains_translucent = true;
1022  if (drawaxes && drawaxes != 3)
1023  {
1024  // add coordinate cross to draw list
1025  params.projection.ortho(-1.,1.,-1.,1.,-2.,2.);
1026  params.model_view.identity();
1027  params.model_view.translate(-1, -1, 0.0);
1028  params.model_view.scale(40.0 / w, 40.0 / h, 1);
1029  params.model_view.translate(2.0, 2.0, 0.0);
1030  params.model_view.mult(cam.RotMatrix());
1031  params.model_view.mult(rotmat);
1032  scene.queue.emplace_back(params, &coord_cross_buf);
1033  }
1034  params.projection.mtx = proj_mtx;
1035  params.model_view.mtx = GetModelViewMtx();
1036  if (drawaxes)
1037  {
1038  // add axes to draw list
1039  scene.queue.emplace_back(params, &axes_buf);
1040  }
1041  params.contains_translucent = false;
1042  if (ruler_on)
1043  {
1044  // add ruler to draw list
1045  params = GetMeshDrawParams();
1046  params.use_clip_plane = false;
1047  params.contains_translucent = false;
1048  scene.queue.emplace_back(params, &ruler_buf);
1049  }
1050  return scene;
1051 }
1052 
1054  glTF_Builder &bld,
1055  glTF_Builder::buffer_id buffer,
1056  glTF_Builder::material_id black_mat)
1057 {
1059 
1060  auto box_node = AddModelNode(bld, "Box");
1061  auto box_mesh = bld.addMesh("Box Mesh");
1062  bld.addNodeMesh(box_node, box_mesh);
1063 
1064  int nlines = AddLines(
1065  bld,
1066  box_mesh,
1067  buffer,
1068  (drawaxes != 3) ? black_mat : white_mat,
1069  axes_buf);
1070  if (nlines == 0)
1071  {
1072  cout << "glTF export: no box found to export!" << endl;
1073  }
1074 }
1075 
1077  glTF_Builder &bld,
1078  glTF_Builder::buffer_id buffer,
1079  glTF_Builder::material_id palette_mat,
1080  const gl3::GlDrawable &gl_drawable)
1081 {
1082  auto elements_node = AddModelNode(bld, "Elements");
1083  auto elements_mesh = bld.addMesh("Elements Mesh");
1084  bld.addNodeMesh(elements_node, elements_mesh);
1085 
1086  int ntria = AddTriangles(
1087  bld,
1088  elements_mesh,
1089  buffer,
1090  palette_mat,
1091  gl_drawable);
1092  if (ntria == 0)
1093  {
1094  cout << "glTF export: no elements found to export!" << endl;
1095  }
1096 }
1097 
1099  glTF_Builder &bld,
1100  glTF_Builder::buffer_id buffer,
1101  glTF_Builder::material_id black_mat,
1102  const gl3::GlDrawable &gl_drawable)
1103 {
1104  auto lines_node = AddModelNode(bld, "Lines");
1105  auto lines_mesh = bld.addMesh("Lines Mesh");
1106  bld.addNodeMesh(lines_node, lines_mesh);
1107 
1108  int nlines = AddLines(
1109  bld,
1110  lines_mesh,
1111  buffer,
1112  black_mat,
1113  gl_drawable);
1114  if (nlines == 0)
1115  {
1116  cout << "glTF export: no mesh/level lines found to export!" << endl;
1117  }
1118 }
1119 
1121 {
1122  cout << "glTF export is not yet implemented for this visualization mode."
1123  << endl;
1124 }
1125 
1127 {
1128  int isSmooth = palette.GetSmoothSetting();
1129  if (isSmooth)
1130  {
1131  palette.UseDiscrete();
1132  cout << "Texture type : discrete" << endl;
1133  }
1134  else
1135  {
1136  palette.UseSmooth();
1137  cout << "Texture type : smooth" << endl;
1138  }
1139 }
1140 
1142 {
1143  if (autoscale != _autoscale)
1144  {
1145  autoscale = _autoscale;
1146  DoAutoscale(true);
1147  }
1148 }
1149 
1151  Mesh & m, Vector & s)
1152  : a_label_x("x"), a_label_y("y"), a_label_z("z")
1153 {
1154  mesh = &m;
1155  sol = &s;
1156 
1157  Init();
1158 }
1159 
1161 {
1162  vsdata = this;
1163  wnd = GetAppWindow();
1164 
1166  scaling = 0;
1167  drawaxes = colorbar = 0;
1168  auto_ref_max = 16;
1169  auto_ref_max_surf_elem = 20000;
1170  minv = 0.0;
1171  maxv = 1.0;
1172  logscale = false;
1173  SetLogA();
1174  PrepareCaption(); // turn on or off the caption
1175 
1176  CuttingPlane = NULL;
1177 
1178  key_r_state = 0;
1179 
1180  // static int init = 0;
1181  // if (!init)
1182  {
1183  // init = 1;
1184 
1185  wnd->setOnKeyDown('l', KeylPressed);
1186  wnd->setOnKeyDown('L', KeyLPressed);
1187 
1188  wnd->setOnKeyDown('s', KeySPressed);
1189 
1190  // wnd->setOnKeyDown('a', KeyaPressed);
1192  wnd->setOnKeyDown('A', KeyAPressed);
1193 
1194  wnd->setOnKeyDown('r', KeyrPressed);
1195  wnd->setOnKeyDown('R', KeyRPressed);
1196 
1197  wnd->setOnKeyDown('p', KeypPressed);
1198  wnd->setOnKeyDown('P', KeyPPressed);
1199 
1200  wnd->setOnKeyDown('h', KeyHPressed);
1201  wnd->setOnKeyDown('H', KeyHPressed);
1202 
1203  wnd->setOnKeyDown(SDLK_F5, KeyF5Pressed);
1204  wnd->setOnKeyDown(SDLK_F6, KeyF6Pressed);
1205  wnd->setOnKeyDown(SDLK_F7, KeyF7Pressed);
1206 
1207  wnd->setOnKeyDown(SDLK_BACKSLASH, KeyBackslashPressed);
1208  wnd->setOnKeyDown('t', KeyTPressed);
1209  wnd->setOnKeyDown('T', KeyTPressed);
1210 
1211  wnd->setOnKeyDown('g', KeygPressed);
1212  wnd->setOnKeyDown('G', KeyGPressed);
1213 
1214  wnd->setOnKeyDown('c', KeycPressed);
1215  wnd->setOnKeyDown('C', KeyCPressed);
1216 
1217  wnd->setOnKeyDown('k', KeykPressed);
1218  wnd->setOnKeyDown('K', KeyKPressed);
1219 
1220  wnd->setOnKeyDown(SDLK_F1, KeyF1Pressed);
1221  wnd->setOnKeyDown(SDLK_F2, KeyF2Pressed);
1222 
1223  wnd->setOnKeyDown(SDLK_COMMA, KeyCommaPressed);
1224  wnd->setOnKeyDown(SDLK_LESS, KeyLessPressed);
1227 
1228  wnd->setOnKeyDown(SDLK_EXCLAIM, KeyToggleTexture);
1229  }
1230 
1231  //Set_Light();
1232 
1233  //glEnable (GL_COLOR_MATERIAL);
1234  //glShadeModel (GL_SMOOTH);
1235 
1236  //gl->enableLight();
1237  //gl->enableDepthTest();
1238  //glEnable(GL_AUTO_NORMAL);
1239  //glEnable(GL_NORMALIZE);
1240 
1241  if (GetMultisample() > 0)
1242  {
1243  glDisable(GL_MULTISAMPLE);
1244  }
1245  // add black fog
1246  // glEnable(GL_FOG);
1247  // GLfloat fogcol[4] = {0,0,0,1};
1248  // glFogfv(GL_FOG_COLOR, fogcol);
1249  // glFogf(GL_FOG_DENSITY,1.0f);
1250 
1251  SetLevelLines(minv, maxv, 15);
1252 
1253  FindNewBox(false);
1254  ruler_on = 0;
1255  ruler_x = 0.5 * (bb.x[0] + bb.x[1]);
1256  ruler_y = 0.5 * (bb.y[0] + bb.y[1]);
1257  ruler_z = 0.5 * (bb.z[0] + bb.z[1]);
1258 
1259  PrepareRuler();
1260 
1261  autoscale = 1;
1262 }
1263 
1265 {
1266  delete CuttingPlane;
1267 }
1268 
1270 {
1271  // double eps = 1e-12;
1272  double eps = 0.0;
1273 
1274  // Find the new scaling
1275  if (scaling)
1276  {
1277  // Scale all sides of the box to 1.
1278  xscale = yscale = zscale = 1.;
1279  if ((bb.x[1]-bb.x[0])>eps) { xscale /= (bb.x[1]-bb.x[0]); }
1280  if ((bb.y[1]-bb.y[0])>eps) { yscale /= (bb.y[1]-bb.y[0]); }
1281  if ((bb.z[1]-bb.z[0])>eps) { zscale /= (bb.z[1]-bb.z[0]); }
1282  }
1283  else
1284  {
1285  // Find the largest side of the box in xscale
1286  xscale = bb.x[1]-bb.x[0];
1287  yscale = bb.y[1]-bb.y[0];
1288  zscale = bb.z[1]-bb.z[0];
1289  if (xscale < yscale) { xscale = yscale; }
1290  if (xscale < zscale) { xscale = zscale; }
1291  // Set proportional scaling so that the largest side of the box is 1.
1292  if (xscale > eps)
1293  {
1294  xscale = ( 1.0 / xscale );
1295  }
1296  else
1297  {
1298  xscale = 1.0;
1299  }
1300  zscale = yscale = xscale;
1301  }
1302 }
1303 
1304 
1306 {
1307  minv = min;
1308  maxv = max;
1309 
1310  UpdateValueRange(true);
1311 }
1312 
1314  const char * a_y,
1315  const char * a_z)
1316 {
1317  a_label_x = a_x;
1318  a_label_y = a_y;
1319  a_label_z = a_z;
1320  PrepareAxes();
1321 }
1322 
1324 {
1325  // std::array<float, 4> blk = GetLineColor();
1326  axes_buf.clear();
1327 
1329  if (drawaxes == 3)
1330  {
1331  // glLineStipple(1, 255);
1332  bld.glBegin(GL_LINES);
1333  bld.glColor3f(1., 0., 0.);
1334  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[0]);
1335  bld.glColor3f(0.75, 0.75, 0.75); // bld.glColor4fv(blk.data());
1336  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[0]);
1337  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[0]);
1338  bld.glColor3f(0., 1., 0.);
1339  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[0]);
1340  bld.glEnd();
1341  // bld.glColor4fv(blk.data());
1342  // bld.setUseColor(false);
1343  // bld.glEnable(GL_LINE_STIPPLE);
1344  bld.glBegin(GL_LINE_STRIP);
1345  bld.glColor3f(0.75, 0.75, 0.75);
1346  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[0]);
1347  bld.glVertex3d(bb.x[1], bb.y[1], bb.z[0]);
1348  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[0]);
1349  bld.glEnd();
1350  }
1351  else
1352  {
1353  // bld.setUseColor(false);
1354  bld.glBegin(GL_LINE_LOOP);
1355  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[0]);
1356  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[0]);
1357  bld.glVertex3d(bb.x[1], bb.y[1], bb.z[0]);
1358  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[0]);
1359  bld.glEnd();
1360  }
1361  bld.glBegin(GL_LINE_LOOP);
1362  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[1]);
1363  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[1]);
1364  bld.glVertex3d(bb.x[1], bb.y[1], bb.z[1]);
1365  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[1]);
1366  bld.glEnd();
1367 
1368  if (drawaxes == 3)
1369  {
1370  // bld.setUseColor(true);
1371  // bld.glDisable(GL_LINE_STIPPLE);
1372  bld.glBegin(GL_LINES);
1373  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[1]);
1374  bld.glColor3f(0., 0., 1.);
1375  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[0]);
1376  bld.glEnd();
1377  // bld.setUseColor(false);
1378  // bld.glEnable(GL_LINE_STIPPLE);
1379  bld.glColor3f(0.75, 0.75, 0.75);
1380  bld.glBegin(GL_LINES);
1381  }
1382  else
1383  {
1384  bld.glBegin(GL_LINES);
1385  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[0]);
1386  bld.glVertex3d(bb.x[0], bb.y[0], bb.z[1]);
1387  }
1388  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[0]);
1389  bld.glVertex3d(bb.x[1], bb.y[0], bb.z[1]);
1390  bld.glVertex3d(bb.x[1], bb.y[1], bb.z[0]);
1391  bld.glVertex3d(bb.x[1], bb.y[1], bb.z[1]);
1392  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[0]);
1393  bld.glVertex3d(bb.x[0], bb.y[1], bb.z[1]);
1394  bld.glEnd();
1395  if (drawaxes == 3)
1396  {
1397  // bld.glDisable(GL_LINE_STIPPLE);
1398  }
1399 
1400  // Write the coordinates of the lower left and upper right corner.
1401  // glEnable (GL_COLOR_MATERIAL);
1402  // GLfloat textcol[3] = { 0, 0, 0 };
1403  // glColor3fv (textcol);
1404 
1405  if (drawaxes == 1)
1406  {
1407  int desc = GetFont()->getFontDescender();
1408  int ox = -desc/2;
1409  int oy = -3*desc/2;
1410  ostringstream buf;
1411  buf << setprecision(4)
1412  << "(" << bb.x[0] << "," << bb.y[0] << "," << bb.z[0] << ")" ;
1413  axes_buf.addText(bb.x[0], bb.y[0], bb.z[0], ox, oy, buf.str());
1414 
1415  ostringstream buf1;
1416  buf1 << setprecision(4)
1417  << "(" << bb.x[1] << "," << bb.y[1] << "," << bb.z[1] << ")" ;
1418  axes_buf.addText(bb.x[1], bb.y[1], bb.z[1], ox, oy, buf1.str());
1419  }
1420  updated_bufs.emplace_back(&axes_buf);
1421 
1422  constexpr float len = 1.2f;
1423  constexpr float l = .9f, cl = .27f, cb = l-cl;
1426  {
1427  {0, 0, 0}, {0, 0, l},
1428  {0, 0, 0}, {0, l, 0},
1429  {0, 0, 0}, {l, 0, 0}
1430  });
1431  coord_cross_buf.addCone(0,0,cb,0,0,1, cl);
1432  coord_cross_buf.addCone(0,cb,0,0,1,0, cl);
1433  coord_cross_buf.addCone(cb,0,0,1,0,0, cl);
1434  coord_cross_buf.addText(len, 0.0f, 0.0f, a_label_x);
1435  coord_cross_buf.addText(0.0f, len, 0.0f, a_label_y);
1436  coord_cross_buf.addText(0.0f, 0.0f, len, a_label_z);
1437  updated_bufs.emplace_back(&coord_cross_buf);
1438 }
1439 
1441  gl3::GlBuilder& builder, double * point, int n, Array<double> &mesh_level,
1442  bool log_vals)
1443 {
1444  int l, k, k1;
1445  double curve, t;
1446  double p[3];
1447 
1448  for (l = 0; l < mesh_level.Size(); l++)
1449  {
1450  // Using GL_LINE_STRIP (explicitly closed for more than 2 points)
1451  // should produce the same result, however visually the level lines
1452  // have discontinuities. Using GL_LINE_LOOP does not have that problem.
1453  builder.glBegin(GL_LINE_LOOP);
1454  curve = LogVal(mesh_level[l], log_vals);
1455  for (k = 0; k < n; k++)
1456  {
1457  k1 = (k+1)%n;
1458  if ( (curve <=point[4*k+3] && curve >= point[4*k1+3]) ||
1459  (curve >=point[4*k+3] && curve <= point[4*k1+3]) )
1460  {
1461  if ((curve - point[4*k1+3]) == 0.)
1462  {
1463  t = 1.;
1464  }
1465  else if ((curve - point[4*k+3]) == 0.)
1466  {
1467  t = 0.;
1468  }
1469  else
1470  {
1471  t = (curve - point[4*k+3]) / (point[4*k1+3]-point[4*k+3]);
1472  }
1473  p[0] = (1.0-t)*point[4*k+0]+t*point[4*k1+0];
1474  p[1] = (1.0-t)*point[4*k+1]+t*point[4*k1+1];
1475  p[2] = (1.0-t)*point[4*k+2]+t*point[4*k1+2];
1476  builder.glVertex3d(p[0], p[1], p[2]);
1477  }
1478  }
1479  builder.glEnd();
1480  }
1481 }
1482 
1484  double min, double max, int n, int adj)
1485 {
1486  int i;
1487  double t, eps;
1488 
1489  if (min < minv)
1490  {
1491  min = minv;
1492  cout << "min set to minv : " << min << endl;
1493  }
1494  if (max > maxv)
1495  {
1496  max = maxv;
1497  cout << "max set to maxv : " << max << endl;
1498  }
1499 
1500  nl = n;
1501  level.SetSize(nl+1);
1502  for (i = 0; i <= nl; i++)
1503  {
1504  t = (double) i / nl;
1505  level[i] = min * (1.0 - t) + t * max;
1506  }
1507 
1508  if (adj)
1509  {
1510  eps = 1.0E-5;
1511  level[0] = level[0] * (1.0 - eps) + level[1] * eps;
1512  level[nl] = level[nl] * (1.0 - eps) + level[nl-1] * eps;
1513  }
1514 
1515  if (logscale)
1516  {
1517  for (i = 0; i <= nl; i++)
1518  {
1519  level[i] = _ULogVal((level[i] - minv) / (maxv - minv));
1520  }
1521  }
1522 }
1523 
1525 {
1526  cout << "\nkeys: " << GetAppWindow()->getSavedKeys() << "\n"
1527  << "\nlight " << strings_off_on[use_light ? 1 : 0]
1528  << "\nperspective " << strings_off_on[OrthogonalProjection ? 0 : 1]
1529  << "\nviewcenter " << ViewCenterX << ' ' << ViewCenterY
1530  << "\nzoom " << (OrthogonalProjection ? ViewScale :
1531  tan(M_PI / 8.) / tan(ViewAngle * (M_PI / 360.0)))
1532  << "\nvaluerange " << minv << ' ' << maxv;
1533  const float *r = glm::value_ptr(rotmat);
1534  ios::fmtflags fmt = cout.flags();
1535  cout << fixed << showpos
1536  << "\nrotmat " << r[ 0] << ' ' << r[ 1] << ' ' << r[ 2] << ' ' << r[ 3]
1537  << "\n " << r[ 4] << ' ' << r[ 5] << ' ' << r[ 6] << ' ' << r[ 7]
1538  << "\n " << r[ 8] << ' ' << r[ 9] << ' ' << r[10] << ' ' << r[11]
1539  << "\n " << r[12] << ' ' << r[13] << ' ' << r[14] << ' ' << r[15]
1540  << '\n' << endl;
1541  cam.Print();
1542  cout.flags(fmt);
1543  mesh->PrintInfo(cout);
1544 }
1545 
1547  int i, int fn, int di)
1548 {
1549  int dim = mesh->Dimension();
1550  int sdim = mesh->SpaceDimension();
1551 
1552  if (shrink != 1.0)
1553  {
1554  if (dim == 2)
1555  {
1556  int d, k;
1557 
1558  for (d = 0; d < sdim; d++)
1559  {
1560  double cd = 0.0;
1561  for (k = 0; k < pointmat.Width(); k++)
1562  {
1563  cd += pointmat(d,k);
1564  }
1565  cd /= pointmat.Width();
1566 
1567  for (k = 0; k < pointmat.Width(); k++)
1568  {
1569  pointmat(d,k) = shrink*pointmat(d,k) + (1-shrink)*cd;
1570  }
1571  }
1572  }
1573  else
1574  {
1575  int attr = mesh->GetBdrAttribute(i);
1576  for (int k = 0; k < pointmat.Width(); k++)
1577  for (int d = 0; d < sdim; d++)
1578  {
1579  pointmat(d,k) = shrink*pointmat(d,k) + (1-shrink)*bdrc(d,attr-1);
1580  }
1581  }
1582  }
1583 
1584  if (shrinkmat != 1.0)
1585  {
1586  int attr, elem1, elem2;
1587  if (dim == 2 || sdim == 2)
1588  {
1589  attr = mesh->GetAttribute(i);
1590  }
1591  else
1592  {
1593  mesh->GetFaceElements(fn, &elem1, &elem2);
1594  if (di == 0)
1595  {
1596  attr = mesh->GetAttribute(elem1);
1597  }
1598  else
1599  {
1600  attr = mesh->GetAttribute(elem2);
1601  }
1602  }
1603 
1604  for (int k = 0; k < pointmat.Width(); k++)
1605  for (int d = 0; d < pointmat.Height(); d++)
1606  {
1607  pointmat(d,k) = shrinkmat*pointmat(d,k) + (1-shrinkmat)*matc(d,attr-1);
1608  }
1609  }
1610 }
1611 
1613 {
1614  DenseMatrix pointmat;
1615  Vector nbdrc(mesh->bdr_attributes.Max());
1616  int sdim = mesh->SpaceDimension();
1617 
1618  bdrc.SetSize(sdim,mesh->bdr_attributes.Max());
1619  bdrc = 0.0;
1620  nbdrc = 0.0;
1621 
1622  for (int i = 0; i < mesh -> GetNBE(); i++)
1623  {
1624  mesh->GetBdrPointMatrix(i, pointmat);
1625  nbdrc(mesh->GetBdrAttribute(i)-1) += pointmat.Width();
1626  for (int k = 0; k < pointmat.Width(); k++)
1627  for (int d = 0; d < sdim; d++)
1628  {
1629  bdrc(d,mesh->GetBdrAttribute(i)-1) += pointmat(d,k);
1630  }
1631  }
1632 
1633  for (int i = 0; i < mesh->bdr_attributes.Max(); i++)
1634  if (nbdrc(i) != 0)
1635  for (int d = 0; d < sdim; d++)
1636  {
1637  bdrc(d,i) /= nbdrc(i);
1638  }
1639 }
1640 
1642 {
1643  DenseMatrix pointmat;
1644  Vector nmatc(mesh->attributes.Max());
1645  int sdim = mesh->SpaceDimension();
1646 
1647  matc.SetSize(sdim,mesh->attributes.Max());
1648  matc = 0.0;
1649  nmatc = 0.0;
1650 
1651  for (int i = 0; i < mesh -> GetNE(); i++)
1652  {
1653  mesh->GetPointMatrix(i, pointmat);
1654  nmatc(mesh->GetAttribute(i)-1) += pointmat.Width();
1655  for (int k = 0; k < pointmat.Width(); k++)
1656  for (int d = 0; d < sdim; d++)
1657  {
1658  matc(d,mesh->GetAttribute(i)-1) += pointmat(d,k);
1659  }
1660  }
1661 
1662  for (int i = 0; i < mesh->attributes.Max(); i++)
1663  if (nmatc(i) != 0)
1664  for (int d = 0; d < sdim; d++)
1665  {
1666  matc(d,i) /= nmatc(i);
1667  }
1668 }
1669 
1670 
1671 Plane::Plane(double A,double B,double C,double D)
1672 {
1673  eqn[0] = A;
1674  eqn[1] = B;
1675  eqn[2] = C;
1676  eqn[3] = D;
1677 
1678  CartesianToSpherical();
1679 
1680  double x[2] = {vsdata -> bb.x[0], vsdata -> bb.x[1]};
1681  double y[2] = {vsdata -> bb.y[0], vsdata -> bb.y[1]};
1682  double z[2] = {vsdata -> bb.z[0], vsdata -> bb.z[1]};
1683  bbox_diam = sqrt ( (x[1]-x[0])*(x[1]-x[0]) +
1684  (y[1]-y[0])*(y[1]-y[0]) +
1685  (z[1]-z[0])*(z[1]-z[0]) );
1686 
1687  x0 = (x[0]+x[1])/2.0;
1688  y0 = (y[0]+y[1])/2.0;
1689  z0 = (z[0]+z[1])/2.0;
1690 
1691  phi_step = M_PI / 36;
1692  theta_step = M_PI / 36;
1693  rho_step = bbox_diam / 200;
1694 }
1695 
1696 void Plane::CartesianToSpherical()
1697 {
1698  rho = sqrt(eqn[0]*eqn[0]+eqn[1]*eqn[1]+eqn[2]*eqn[2]);
1699  phi = asin(eqn[2]/rho);
1700  theta = atan2(eqn[1], eqn[0]);
1701 }
1702 
1703 void Plane::SphericalToCartesian()
1704 {
1705  eqn[0] = rho * cos(phi) * cos(theta);
1706  eqn[1] = rho * cos(phi) * sin(theta);
1707  eqn[2] = rho * sin(phi);
1708  eqn[3] = - (eqn[0] * x0 + eqn[1] * y0 + eqn[2] * z0);
1709 }
1710 
1712 {
1713  phi += phi_step;
1714  SphericalToCartesian();
1715 }
1716 
1718 {
1719  phi -= phi_step;
1720  SphericalToCartesian();
1721 }
1722 
1724 {
1725  theta += theta_step;
1726  SphericalToCartesian();
1727 }
1728 
1730 {
1731  theta -= theta_step;
1732  SphericalToCartesian();
1733 }
1734 
1736 {
1737  double k = (rho_step) / (rho*rho);
1738  x0 -= eqn[0] * k;
1739  y0 -= eqn[1] * k;
1740  z0 -= eqn[2] * k;
1741  eqn[3] += rho_step;
1742  CartesianToSpherical();
1743 }
1744 
1746 {
1747  double k = (rho_step) / (rho*rho);
1748  x0 += eqn[0] * k;
1749  y0 += eqn[1] * k;
1750  z0 += eqn[2] * k;
1751  eqn[3] -= rho_step;
1752  CartesianToSpherical();
1753 }
void SendExposeEvent()
Send expose event. In our case MyReshape is executed and Draw after it.
Definition: aux_vis.cpp:346
void KeyCommaPressed()
Definition: vsdata.cpp:760
void DecreasePhi()
Definition: vsdata.cpp:1717
void KeyF6Pressed()
Definition: vsdata.cpp:591
void identity()
Sets the matrix to the identity matrix.
Definition: types.hpp:164
void addLines(const std::vector< Vert > &verts)
Definition: types.hpp:679
void KeyHPressed()
Definition: vsdata.cpp:514
void SetValueRange(double, double)
Definition: vsdata.cpp:1305
void ShrinkPoints(DenseMatrix &pointmat, int i, int fn, int di)
Shrink the set of points towards attributes centers of gravity.
Definition: vsdata.cpp:1546
struct VisualizationScene::@3 bb
Bounding box.
void KeyTPressed()
Definition: vsdata.cpp:690
thread_local VisualizationSceneScalarData * vsdata
Definition: vsdata.cpp:466
void glVertex3d(double x, double y, double z)
Definition: types.hpp:332
int ChoosePalette()
Definition: palettes.cpp:7503
int Normalize(double v[])
Definition: geom_utils.hpp:38
void getGLDrawSize(int &w, int &h)
Definition: sdl.cpp:563
float getFontLineSpacing() const
Definition: font.hpp:99
int GetSize(int pal=-1) const
Gets the total number of colors in the current palette color array.
Definition: palettes.cpp:7815
void KeyrPressed()
Definition: vsdata.cpp:531
int GetNumColors(int pal=-1) const
Gets the number of colors used in the current palette color array.
Definition: palettes.hpp:42
thread_local SdlWindow * wnd
Definition: aux_vis.cpp:50
void PrepareColorBar(double minval, double maxval, Array< double > *level=NULL, Array< double > *levels=NULL)
Definition: vsdata.cpp:333
void KeyAPressed()
Definition: vsdata.cpp:748
void KeykPressed()
Definition: vsdata.cpp:726
void glColor3f(float r, float g, float b)
Definition: types.hpp:425
gl3::MeshRenderer & getRenderer()
Definition: sdl.hpp:226
virtual void ToggleLogscale(bool print)
Definition: vsdata.cpp:811
void glEnd()
Definition: types.hpp:309
virtual gl3::SceneInfo GetSceneObjs()
Definition: vsdata.cpp:984
gl3::GlDrawable coord_cross_buf
Definition: vsdata.hpp:70
PaletteState palette
Definition: openglvis.hpp:178
void KeyGravePressed()
Definition: vsdata.cpp:782
double LogVal(const double &z, const bool &log_val)
Definition: vsdata.hpp:121
virtual std::string GetHelpString() const
Definition: vsdata.hpp:144
Plane(double A, double B, double C, double D)
Definition: vsdata.cpp:1671
void ComputeBdrAttrCenter()
Compute the center of gravity for each boundary attribute.
Definition: vsdata.cpp:1612
void glNormal3dv(const double *d)
Definition: types.hpp:407
void IncreasePhi()
Definition: vsdata.cpp:1711
float getFontDescender() const
Definition: font.hpp:101
std::array< uint8_t, 4 > ColorU8(float rgba[])
Definition: types.hpp:181
void KeyF2Pressed()
Definition: vsdata.cpp:718
void glBegin(GLenum e)
Definition: types.hpp:295
void DrawPolygonLevelLines(gl3::GlBuilder &builder, double *point, int n, Array< double > &level, bool log_vals)
Definition: vsdata.cpp:1440
void addText(float x, float y, float z, const std::string &text)
Adds a string at the given position in object coordinates.
Definition: types.hpp:658
void KeyKPressed()
Definition: vsdata.cpp:737
void SetRepeatTimes(int rpt)
Definition: palettes.hpp:47
void PrintLogscale(bool warn)
Definition: vsdata.cpp:800
virtual void EventUpdateColors()
Definition: vsdata.hpp:180
void Arrow(gl3::GlBuilder &builder, double px, double py, double pz, double vx, double vy, double vz, double length, double cone_scale=0.075)
Definition: vsdata.cpp:205
thread_local string plot_caption
Definition: glvis.cpp:60
int GetRepeatTimes() const
Definition: palettes.hpp:46
int GetMultisample()
Definition: aux_vis.cpp:1573
virtual void FindNewBox(bool prepare)=0
void PrevIndex()
Definition: palettes.cpp:7863
void KeyaPressed()
Definition: vsdata.cpp:489
void NextIndex()
Definition: palettes.cpp:7858
void getWindowSize(int &w, int &h)
Definition: sdl.cpp:534
double shrinkmat
Shrink factor with respect to the element (material) attributes centers.
Definition: vsdata.hpp:136
double _ULogVal(const double &u)
Definition: vsdata.hpp:110
void SetLevelLines(double min, double max, int n, int adj=1)
Definition: vsdata.cpp:1483
void KeyToggleTexture()
Definition: vsdata.cpp:794
void addCone(float x, float y, float z, float vx, float vy, float vz, float cone_scale=0.075)
Definition: types.cpp:19
void glNormal3d(double nx, double ny, double nz)
Definition: types.hpp:401
virtual void EventUpdateBackground()
Definition: vsdata.hpp:179
void DecreaseTheta()
Definition: vsdata.cpp:1729
void SetLight0CustomPos(std::array< float, 4 > pos)
Definition: openglvis.cpp:1029
void Arrow2(gl3::GlBuilder &builder, double px, double py, double pz, double vx, double vy, double vz, double length, double cone_scale=0.075)
Definition: vsdata.cpp:172
void KeyRPressed()
Definition: vsdata.cpp:546
bool getAntialiasing()
Definition: renderer.hpp:238
GlBuilder createBuilder()
Definition: types.hpp:731
virtual void glTF_Export()
Definition: vsdata.cpp:1120
GlVisFont * GetFont()
Definition: aux_vis.cpp:1678
glm::mat4 mtx
Definition: types.hpp:121
gl3::GlDrawable axes_buf
Definition: vsdata.hpp:69
void Print()
Definition: openglvis.cpp:106
void SetIndex(int num)
Sets the palette texture to bind.
Definition: palettes.hpp:31
void RemoveIdleFunc(void(*Func)(void))
Definition: aux_vis.cpp:497
mesh_id addMesh(const std::string &meshName)
Definition: gltf.cpp:322
virtual void PrintState()
Definition: vsdata.cpp:1524
void KeycPressed()
Definition: vsdata.cpp:469
std::string getSavedKeys() const
Returns the keyboard events that have been logged by the window.
Definition: sdl.hpp:238
SdlWindow * wnd
Definition: openglvis.hpp:68
void glTF_ExportElements(glTF_Builder &bld, glTF_Builder::buffer_id buffer, glTF_Builder::material_id palette_mat, const gl3::GlDrawable &gl_drawable)
Definition: vsdata.cpp:1076
void rotate(float angle, double x, double y, double z)
Applies a rotation transform to the matrix.
Definition: types.hpp:124
vector< gl3::GlDrawable * > updated_bufs
Definition: vsdata.hpp:68
void Arrow3(gl3::GlBuilder &builder, double px, double py, double pz, double vx, double vy, double vz, double length, double cone_scale=0.075)
Definition: vsdata.cpp:114
void Key_Mod_a_Pressed(GLenum state)
Definition: vsdata.cpp:495
void KeyTildePressed()
Definition: vsdata.cpp:788
void SetLightMatIdx(unsigned i)
Definition: openglvis.cpp:1020
void KeyF1Pressed()
Definition: vsdata.cpp:713
SdlWindow * GetAppWindow()
Definition: aux_vis.cpp:58
Crude fixed-function OpenGL emulation helper.
Definition: types.hpp:261
void Init()
Initializes the palette textures.
Definition: palettes.cpp:7694
void GenerateAlphaTexture()
Definition: openglvis.hpp:220
thread_local string extra_caption
Definition: glvis.cpp:61
void MainLoop()
Definition: aux_vis.cpp:515
void KeySPressed()
Definition: vsdata.cpp:483
void KeylPressed()
Definition: vsdata.cpp:519
void KeyGPressed()
Definition: vsdata.cpp:708
void KeygPressed()
Definition: vsdata.cpp:700
static constexpr unsigned INVALID_ID
Definition: gltf.hpp:218
void SetNumColors(int numColors)
Sets the number of colors to use in the current palette color array.
Definition: palettes.hpp:45
void addNodeMesh(node_id node, mesh_id mesh)
Definition: gltf.cpp:433
void getObjectSize(const std::string &text, int &w, int &h)
Get the width and height of the bounding box containing the rendered text.
Definition: font.cpp:148
void glTF_ExportBox(glTF_Builder &bld, glTF_Builder::buffer_id buffer, glTF_Builder::material_id black_mat)
Definition: vsdata.cpp:1053
void DecreaseDistance()
Definition: vsdata.cpp:1745
void KeyF7Pressed(GLenum state)
Definition: vsdata.cpp:627
void KeyLPressed()
Definition: vsdata.cpp:525
void KeyCPressed()
Definition: vsdata.cpp:475
void SetAutoscale(int _autoscale)
Definition: vsdata.cpp:1141
void glTF_ExportMesh(glTF_Builder &bld, glTF_Builder::buffer_id buffer, glTF_Builder::material_id black_mat, const gl3::GlDrawable &gl_drawable)
Definition: vsdata.cpp:1098
void KeyPPressed()
Definition: vsdata.cpp:567
void KeyLessPressed()
Definition: vsdata.cpp:771
void IncreaseDistance()
Definition: vsdata.cpp:1735
void KeyCtrlP()
Definition: aux_vis.cpp:1190
void DoAutoscaleValue(bool prepare)
Definition: vsdata.cpp:64
void glVertex3dv(const double *d)
Definition: types.hpp:396
vector< GlDrawable * > needs_buffering
Definition: renderer.hpp:62
virtual void SetNewScalingFromBox()
Definition: vsdata.cpp:1269
const char * strings_off_on[]
Definition: vsdata.cpp:27
void clear()
Clears the drawable object.
Definition: types.hpp:737
virtual void UpdateValueRange(bool prepare)=0
void setOnKeyDown(int key, Delegate func)
Definition: sdl.hpp:188
void SetAxisLabels(const char *a_x, const char *a_y, const char *a_z)
Definition: vsdata.cpp:1313
virtual void PrepareRuler()
Definition: vsdata.hpp:274
void ComputeElemAttrCenter()
Compute the center of gravity for each element attribute.
Definition: vsdata.cpp:1641
virtual ~VisualizationSceneScalarData()
Definition: vsdata.cpp:1264
void KeyBackslashPressed()
Definition: vsdata.cpp:670
void IncreaseTheta()
Definition: vsdata.cpp:1723
void setAntialiasing(bool aa_status)
Definition: renderer.cpp:45
GlMatrix model_view
Definition: renderer.hpp:40
thread_local VisualizationScene * locscene
Definition: aux_vis.cpp:38
void KeypPressed(GLenum state)
Definition: vsdata.cpp:554
Array< double > level
Definition: vsdata.hpp:81
void DoAutoscale(bool prepare)
Definition: vsdata.cpp:48
void Cone(gl3::GlBuilder &builder, glm::mat4 transform)
Definition: vsdata.cpp:76
void Reset()
Definition: openglvis.cpp:25