gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
PViewX3D.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 //
6 // PViewX3D is a extension for Post-processing outputs :
7 // creates a file in X3D format with the same features as
8 // what is visible in post-processing screen.
9 // contact : gilles.marckmann@ec-nantes.fr
10 
11 #include <iostream>
12 #include <limits>
13 #include <ctime>
14 #include <math.h>
15 #include "GmshConfig.h"
16 #include "GmshMessage.h"
17 #include "PView.h"
18 #include "PViewData.h"
19 #include "PViewOptions.h"
20 #include "PViewDataList.h"
21 #include "PViewDataGModel.h"
22 #include "VertexArray.h"
23 #include "StringUtils.h"
24 #include "Context.h"
25 #include "OS.h"
26 #include "SBoundingBox3d.h"
27 #include "PViewX3D.h"
28 
29 using namespace std;
30 
31 static bool almostEqual(double x, double y)
32 {
33  return std::abs(x - y) < CTX::instance()->print.x3dPrecision;
34 }
35 
37  const TriangleToSort *second)
38 {
39  return (first->xmin < second->xmin);
40 }
41 
43  const TriangleToSort *second)
44 {
45  return (first->ymin < second->ymin);
46 }
47 
49  const TriangleToSort *second)
50 {
51  return (first->zmin < second->zmin);
52 }
53 
55  const TriangleToSort *second)
56 {
57  return (first->xmax < second->xmax);
58 }
59 
61  const TriangleToSort *second)
62 {
63  return (first->ymax < second->ymax);
64 }
65 
67  const TriangleToSort *second)
68 {
69  return (first->zmax < second->zmax);
70 }
71 
72 bool PView::writeX3D(const std::string &fileName)
73 {
74  // tags duplicated triangles
75  int _size = 1;
76  if(!CTX::instance()->print.x3dRemoveInnerBorders) {
77  for(std::size_t i = 0; i < PView::list.size(); i++) {
78  VertexArray *va = PView::list[i]->va_triangles;
79  _size += va->getNumVertices() / 3;
80  }
81  }
82  int _count = 0;
83  std::vector<bool> visible(_size);
84  if(!CTX::instance()->print.x3dRemoveInnerBorders) {
85  // evaluate bbox of each triangle
86  std::list<TriangleToSort *> tlist;
87  tlist.clear();
88  for(std::size_t ivp = 0; ivp < PView::list.size(); ivp++) {
89  VertexArray *va = PView::list[ivp]->va_triangles;
90  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 3) {
91  float *p0 = va->getVertexArray(3 * ipt);
92  float *p1 = va->getVertexArray(3 * (ipt + 1));
93  float *p2 = va->getVertexArray(3 * (ipt + 2));
94  TriangleToSort *_current = new TriangleToSort;
95  _current->_index = ipt;
96  _current->_globalIndex = _count;
97  visible[_count] = true;
98  _count++;
99  _current->_ppv = PView::list[ivp];
100  _current->xmin = min(p0[0], min(p1[0], p2[0]));
101  _current->ymin = min(p0[1], min(p1[1], p2[1]));
102  _current->zmin = min(p0[2], min(p1[2], p2[2]));
103  _current->xmax = max(p0[0], max(p1[0], p2[0]));
104  _current->ymax = max(p0[1], max(p1[1], p2[1]));
105  _current->zmax = max(p0[2], max(p1[2], p2[2]));
106 
107  tlist.push_back(_current);
108  }
109  }
110  // sort triangles upon the position of bbbox
111  tlist.sort(compare_zmax_triangle);
112  tlist.sort(compare_ymax_triangle);
113  tlist.sort(compare_xmax_triangle);
114  tlist.sort(compare_zmin_triangle);
115  tlist.sort(compare_ymin_triangle);
116  tlist.sort(compare_xmin_triangle);
117 
118  // estimate and tags triangles which are identical
119  std::list<TriangleToSort *>::iterator pt, nt;
120  for(pt = tlist.begin(); pt != tlist.end(); pt++) {
121  nt = pt;
122  nt++;
123  bool found = false;
124  VertexArray *vap = ((*pt)->_ppv)->va_triangles;
125  int ip = (*pt)->_index;
126  float *p0 = vap->getVertexArray(3 * ip);
127  float *p1 = vap->getVertexArray(3 * (ip + 1));
128  float *p2 = vap->getVertexArray(3 * (ip + 2));
129  int gip = (*pt)->_globalIndex;
130  while(nt != tlist.end() && !found) {
131  int gin = (*nt)->_globalIndex;
132  if((((abs((*pt)->xmin - (*nt)->xmin) < 1.e-9) &&
133  (abs((*pt)->ymin - (*nt)->ymin) < 1.e-9)) &&
134  (abs((*pt)->zmin - (*nt)->zmin) < 1.e-9)) &&
135  (((abs((*pt)->xmax - (*nt)->xmax) < 1.e-9) &&
136  (abs((*pt)->ymax - (*nt)->ymax) < 1.e-9)) &&
137  (abs((*pt)->zmax - (*nt)->zmax) < 1.e-9))) {
138  VertexArray *van = ((*nt)->_ppv)->va_triangles;
139  int in = (*nt)->_index;
140  float *n0 = van->getVertexArray(3 * in);
141  float *n1 = van->getVertexArray(3 * (in + 1));
142  float *n2 = van->getVertexArray(3 * (in + 2));
143 
144  if(almostEqual(p0[0], n0[0]) && almostEqual(p0[1], n0[1]) &&
145  almostEqual(p0[2], n0[2])) {
146  if(almostEqual(p1[0], n1[0]) && almostEqual(p1[1], n1[1]) &&
147  almostEqual(p1[2], n1[2])) {
148  if(almostEqual(p2[0], n2[0]) && almostEqual(p2[1], n2[1]) &&
149  almostEqual(p2[2], n2[2])) {
150  found = true;
151  }
152  }
153  else if(almostEqual(p1[0], n2[0]) && almostEqual(p1[1], n2[1]) &&
154  almostEqual(p1[2], n2[2])) {
155  if(almostEqual(p2[0], n1[0]) && almostEqual(p2[1], n1[1]) &&
156  almostEqual(p2[2], n1[2])) {
157  found = true;
158  }
159  }
160  }
161  else if(almostEqual(p0[0], n1[0]) && almostEqual(p0[1], n1[1]) &&
162  almostEqual(p0[2], n1[2])) {
163  if(almostEqual(p1[0], n0[0]) && almostEqual(p1[1], n0[1]) &&
164  almostEqual(p1[2], n0[2])) {
165  if(almostEqual(p2[0], n2[0]) && almostEqual(p2[1], n2[1]) &&
166  almostEqual(p2[2], n2[2])) {
167  found = true;
168  }
169  }
170  else if(almostEqual(p1[0], n2[0]) && almostEqual(p1[1], n2[1]) &&
171  almostEqual(p1[2], n2[2])) {
172  if(almostEqual(p2[0], n0[0]) && almostEqual(p2[1], n0[1]) &&
173  almostEqual(p2[2], n0[2])) {
174  found = true;
175  }
176  }
177  }
178  else if(almostEqual(p0[0], n2[0]) && almostEqual(p0[1], n2[1]) &&
179  almostEqual(p0[2], n2[2])) {
180  if(almostEqual(p1[0], n0[0]) && almostEqual(p1[1], n0[1]) &&
181  almostEqual(p1[2], n0[2])) {
182  if(almostEqual(p2[0], n1[0]) && almostEqual(p2[1], n1[1]) &&
183  almostEqual(p2[2], n1[2])) {
184  found = true;
185  }
186  }
187  else if(almostEqual(p1[0], n1[0]) && almostEqual(p1[1], n1[1]) &&
188  almostEqual(p1[2], n1[2])) {
189  if(almostEqual(p2[0], n0[0]) && almostEqual(p2[1], n0[1]) &&
190  almostEqual(p2[2], n0[2])) {
191  found = true;
192  }
193  }
194  }
195 
196  if(found) {
197  visible[gip] = false;
198  visible[gin] = false;
199  if(pt != tlist.end()) pt++;
200  }
201  else {
202  nt++;
203  }
204  }
205  else {
206  nt = tlist.end();
207  }
208  }
209  }
210  for(pt = tlist.begin(); pt != tlist.end(); pt++) {
211  // delete (*pt);
212  }
213  }
214  // end tags duplicated triangles
215 
216  // beginning writing x3d file
217  time_t rawtime;
218  struct tm *timeinfo;
219  time(&rawtime);
220  timeinfo = localtime(&rawtime);
221  FILE *fp = Fopen(fileName.c_str(), "w");
222  if(!fp) {
223  Msg::Error("Unable to open file '%s'", fileName.c_str());
224  return false;
225  }
226 
227  // x3 Header
228  fprintf(fp, "<?xml version=\"1.0\" encoding=\"UTF-8\"?>\n");
229  fprintf(fp, "<!DOCTYPE X3D PUBLIC \"ISO//Web3D//DTD X3D 3.3//EN\" "
230  "\"http://www.web3d.org/specifications/x3d-3.3.dtd\">\n");
231  fprintf(fp, "<X3D profile='Interchange' version='3.3' "
232  "xmlns:xsd='http://www.w3.org/2001/XMLSchema-instance' >\n");
233  fprintf(fp, " <head>\n");
234  fprintf(fp, " <meta name='title' content='PView'/> \n");
235  fprintf(fp, " <meta name='description' content='%s'/>\n",
236  fileName.c_str());
237  fprintf(fp, " <meta name='creator' content='gmsh'/> \n");
238  fprintf(fp, " <meta name='created' content=' %s '/>\n",
239  asctime(timeinfo));
240  fprintf(fp, " </head>\n");
241  // Viewport
242  SBoundingBox3d bb(0., 0., 0., 0., 0., 0.);
243  for(auto pvit = PView::list.begin(); pvit < PView::list.end(); pvit++) {
244  PViewData *data = (*pvit)->getData(true);
245  PViewOptions *opt = (*pvit)->getOptions();
246  if(!data->getDirty() && opt->visible) {
247  bb += data->getBoundingBox(opt->timeStep);
248  }
249  }
250  SPoint3 _center = bb.center();
251  double _diagonal = bb.diag();
252  fprintf(fp, " <Scene>\n");
253  fprintf(fp,
254  " <Viewpoint description='Book View' orientation='0 0. 1. 0.' "
255  "position='%g %g %g'/> \n",
256  _center.x(), _center.y(), _center.z() + _diagonal * 1.2);
257  fprintf(fp, " <Background skyColor='.7 .7 1'/> \n");
258 
259  // HUD : Head Up Display
260  // here contour/scalebar legends in frame (-.45,-.28, 0.) and (.45, .28,0.) :
261  // viewport .9 x .56
262  double viewportWidth = .9;
263  double viewportHeight = .56;
264  double font_size = 0.02;
265  if(!CTX::instance()->print.x3dCompatibility) {
266  std::vector<PView *> scales;
267  for(std::size_t i = 0; i < PView::list.size(); i++) {
268  PViewData *data = PView::list[i]->getData();
269  PViewOptions *opt = PView::list[i]->getOptions();
270  if(!data->getDirty() && opt->visible && opt->showScale &&
271  opt->type == PViewOptions::Plot3D && data->getNumElements())
272  scales.push_back(PView::list[i]);
273  }
274 
275  if(!scales.empty()) {
276  fprintf(fp, " <ProtoDeclare appinfo='Heads-up display (HUD)' "
277  "name='HeadsUpDisplay'> \n");
278  fprintf(fp, " <ProtoInterface> \n");
279  fprintf(fp,
280  " <field accessType='inputOutput' appinfo='offset "
281  "position for HUD' "
282  "name='screenOffset' type='SFVec3f' value='%g %g %g'/> \n",
283  _center.x(), _center.y(), 5 * _center.z() + _diagonal * 1.2);
284  fprintf(fp, " <field accessType='inputOutput' appinfo='X3D "
285  "content positioned at HUD offset' "
286  "name='children' type='MFNode'> \n");
287  fprintf(fp, " </field> \n");
288  fprintf(fp, " <field accessType='outputOnly' appinfo='HUD "
289  "position update (in world "
290  "coordinates) relative to original location' "
291  "name='position_changed' type='SFVec3f'/> \n");
292  fprintf(
293  fp,
294  " <field accessType='outputOnly' appinfo='HUD orientation "
295  "update relative to "
296  "original location' name='orientation_changed' type='SFRotation'/> \n");
297  fprintf(fp, " </ProtoInterface> \n");
298  fprintf(fp, " <ProtoBody> \n");
299  fprintf(fp, " <Group bboxCenter=\"%g %g %g\"> \n", _center.x(),
300  _center.y(), _center.z());
301  fprintf(fp, " <Transform DEF='HudContainer'> \n");
302  fprintf(fp, " <Transform> \n");
303  fprintf(fp, " <IS> \n");
304  fprintf(fp, " <connect nodeField='translation' "
305  "protoField='screenOffset'/> \n");
306  fprintf(fp, " </IS> \n");
307  fprintf(fp, " <Group> \n");
308  fprintf(fp, " <IS> \n");
309  fprintf(fp, " <connect nodeField='children' "
310  "protoField='children'/> \n");
311  fprintf(fp, " </IS> \n");
312  fprintf(fp, " </Group> \n");
313  fprintf(fp, " </Transform> \n");
314  fprintf(fp, " </Transform> \n");
315  fprintf(fp, " <ProximitySensor DEF='HereIAm' size='10000000 "
316  "10000000 10000000'> \n");
317  fprintf(fp, " <IS> \n");
318  fprintf(fp, " <connect nodeField='position_changed' "
319  "protoField='position_changed'/> \n");
320  fprintf(fp, " <connect nodeField='orientation_changed' "
321  "protoField='orientation_changed'/> \n");
322  fprintf(fp, " </IS> \n");
323  fprintf(fp, " </ProximitySensor> \n");
324  fprintf(fp, " <ROUTE fromField='orientation_changed' "
325  "fromNode='HereIAm' toField='rotation' "
326  "toNode='HudContainer'/> \n");
327  fprintf(fp, " <ROUTE fromField='position_changed' "
328  "fromNode='HereIAm' toField='translation' "
329  "toNode='HudContainer'/> \n");
330  fprintf(fp, " </Group> \n");
331  fprintf(fp, " </ProtoBody> \n");
332  fprintf(fp, " </ProtoDeclare> \n");
333  // fprintf(fp," <Background skyColor='.7 .7 1'/> \n");
334  fprintf(fp,
335  " <Viewpoint description='Book View' orientation='0 0. 1. 0.' "
336  "position='%g %g %g'/> \n",
337  _center.x(), _center.y(), _center.z() + _diagonal * 1.2);
338  fprintf(fp, " <!-- ProtoDeclare is the \"cookie cutter\" template, "
339  "ProtoInstance creates an actual "
340  "occurrence --> \n");
341  fprintf(
342  fp,
343  " <ProtoInstance DEF='HeadsUpDisplay' name='HeadsUpDisplay'> \n");
344  fprintf(fp, " <!-- example: upper left-hand corner of screen (x=-2, "
345  "y=1) and set back z=-5 from "
346  "user view --> \n");
347  fprintf(fp, " <fieldValue name='screenOffset' value='%g %g %g'/> \n",
348  _center.x(), _center.y(), _center.z() - .2 * _diagonal);
349  fprintf(fp, " <fieldValue name='children'> \n");
350 
351  char label[1024];
352  double maxw = 10. * font_size * 3. / 4.;
353  const double tic = viewportWidth / 100;
354  const double bar_size = tic * 1.6;
355  double width = 0., width_prev = 0., width_total = 0.;
356 
357  for(std::size_t i = 0; i < scales.size(); i++) {
358  PView *p = scales[i];
359  PViewData *data = p->getData();
360  PViewOptions *opt = p->getOptions();
361 
362  if(!opt->autoPosition) {
363  double w = viewportWidth / 3;
364  double h = viewportHeight / 11;
365  double x = 0.;
366  double y = -viewportHeight;
367  writeX3DScale(fp, p, x, y, w, h, tic,
368  CTX::instance()->post.horizontalScales, font_size);
369  }
370  else if(CTX::instance()->post.horizontalScales) {
371  double ysep = viewportHeight / 40;
372  double xc = 0.;
373  if(scales.size() == 1) {
374  double w = viewportWidth / 2., h = bar_size;
375  double x = xc - w / 2., y = -viewportHeight / 2 + ysep;
376  writeX3DScale(fp, p, x, y, w, h, tic, 1, font_size);
377  }
378  else {
379  double xsep = maxw / 4. + viewportWidth / 10.;
380  double w = (viewportWidth - 4. * xsep) / 2.;
381  if(w < 30. * viewportWidth / 1000) w = 30. * viewportWidth / 1000;
382  double h = bar_size;
383  double x = xc - (i % 2 ? -xsep / 1.5 : w + xsep / 1.5);
384  double y = -viewportHeight / 2 + ysep +
385  (i / 2) * (bar_size + tic + 2 * font_size + ysep);
386  writeX3DScale(fp, p, x, y, w, h, tic, 1, font_size);
387  }
388  }
389  else {
390  double xsep = viewportWidth / 50;
391  double dy = 2. * font_size;
392  if(scales.size() == 1) {
393  double ysep = (viewportHeight) / 6.;
394  double w = bar_size, h = viewportHeight - 2 * ysep - dy;
395  double x = -viewportWidth / 2 + xsep,
396  y = -viewportHeight / 2 + ysep + dy;
397  writeX3DScale(fp, p, x, y, w, h, tic, 1, font_size);
398  }
399  else {
400  double ysep = viewportHeight / 30.;
401  double w = bar_size;
402  double h = (viewportHeight - 3 * ysep - 2.5 * dy) / 2.;
403  double x = -viewportWidth / 2 + xsep + width_total + (i / 2) * xsep;
404  double y = -viewportHeight / 2 + ysep + dy +
405  (1 - i % 2) * (h + 1.5 * dy + ysep);
406  writeX3DScale(fp, p, x, y, w, h, tic, 1, font_size);
407  }
408  // compute width
409  width_prev = width;
410  width = bar_size + tic + 10. * font_size * 3 / 4;
411  if(opt->showTime) {
412  char tmp[256];
413  sprintf(tmp, opt->format.c_str(), data->getTime(opt->timeStep));
414  sprintf(label, "%s (%s)", data->getName().c_str(), tmp);
415  }
416  else {
417  sprintf(label, "%s", data->getName().c_str());
418  }
419  width = max(width, strlen(label) * font_size * 3 / 4);
420  if(i % 2)
421  width_total += std::max(bar_size + width, bar_size + width_prev);
422  }
423  }
424  fprintf(fp, " </fieldValue> \n");
425  fprintf(fp, " </ProtoInstance> \n");
426  }
427  }
428 
429  // geometrical objects
430  VertexArray *va;
431  PViewData *data;
432  PViewOptions *opt;
433 
434  // points - NOT TREATED YET
435  /*
436  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++){
437  data = PView::list[ipv]->getData(true);
438  opt = PView::list[ipv]->getOptions();
439  if( !data->getDirty() && opt->visible ) {
440  va=PView::list[ipv]->va_points;
441  for(int ipt = 0; ipt < va->getNumVertices(); ipt++){
442  float *p = va->getVertexArray(3 * ipt);
443  double f = 1.;
444  if(opt->pointType > 1){
445  char *n = va->getNormalArray(3 * ipt);
446  f = char2float(*n);
447  }
448  if(opt->pointType == 2){
449  int s = (int)(opt->pointSize * f);
450  if(s){
451  fprintf(fp,"points : %g %g %g\n", p[0], p[1], p[2]);
452  }
453  }
454  else
455  fprintf(fp,"sphere : %g %g %g \n", p[0], p[1], p[2] );
456  }
457  } // enf if dirty
458 
459  }// end loop on PView::list
460  */
461 
462  // lines
463  int _ind = 0;
464  fprintf(fp, " <Shape> \n");
465  fprintf(fp, " <IndexedLineSet coordIndex=' ");
466  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
467  PViewData *data = PView::list[ipv]->getData(true);
468  PViewOptions *opt = PView::list[ipv]->getOptions();
469  if(!data->getDirty() && opt->visible) {
470  va = PView::list[ipv]->va_lines;
471  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 2) {
472  if(opt->lineType != 2 && opt->lineType != 1) {
473  fprintf(fp, "%i %i %i ", _ind, _ind + 1, -1);
474  }
475  _ind += 2;
476  }
477  } // end if dirty
478  } // end for loop on PView::list
479  fprintf(fp, "'>\n");
480  fprintf(fp, " <Coordinate DEF='TurnPoints' point=' ");
481  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
482  PViewData *data = PView::list[ipv]->getData(true);
483  PViewOptions *opt = PView::list[ipv]->getOptions();
484  if(!data->getDirty() && opt->visible) {
485  va = PView::list[ipv]->va_lines;
486  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 2) {
487  if(opt->lineType != 2 && opt->lineType != 1) {
488  float *p0 = va->getVertexArray(3 * ipt);
489  float *p1 = va->getVertexArray(3 * (ipt + 1));
490  fprintf(fp, "%g %g %g %g %g %g ", p0[0], p0[1], p0[2], p1[0], p1[1],
491  p1[2]);
492  }
493  } // end for
494  } // end if dirty
495  } // end for loop on PView::list
496  fprintf(fp, "'/> \n");
497  fprintf(fp, " </IndexedLineSet> \n");
498  fprintf(fp, " <Appearance> \n");
499  fprintf(fp, " <Material emissiveColor='0 0 0'/>\n");
500  fprintf(fp, " <LineProperties containerField='lineProperties'> \n");
501  fprintf(fp, " </LineProperties> \n");
502  fprintf(fp, " </Appearance> \n");
503  fprintf(fp, " </Shape>\n");
504 
505  // vectors - colored arrow replaced by colored cones
506  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
507  data = PView::list[ipv]->getData(true);
508  opt = PView::list[ipv]->getOptions();
509  if(!data->getDirty() && opt->visible) {
510  va = PView::list[ipv]->va_vectors;
511  for(int iv = 0; iv < va->getNumVertices(); iv += 2) {
512  float *s = va->getVertexArray(3 * iv);
513  float *v = va->getVertexArray(3 * (iv + 1));
514  unsigned char *c = va->getColorArray(4 * iv);
515  double rgba[4];
516  UnsignedChar2rgba(c, rgba);
517  double l = sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
518  double lmax = opt->tmpMax;
519  if((l || opt->vectorType == 6) && lmax) {
520  double scale = .5 / _diagonal;
521  double theta = acos(v[1] / l);
522  fprintf(
523  fp,
524  " <Transform rotation='%g %g %g %g' translation='%e %e %e' >\n",
525  v[2], 0., -v[0], theta, s[0] + .5 * scale * v[0],
526  s[1] + .5 * scale * v[1], s[2] + .5 * scale * v[2]);
527  fprintf(fp, " <Shape>\n");
528  fprintf(fp, " <Cone bottomRadius='%g' height='%g'/>\n",
529  .05 * l * scale, l * scale);
530  fprintf(fp, " <Appearance>\n");
531  fprintf(fp, " <Material diffuseColor='%g %g %g'/>\n",
532  rgba[0], rgba[1], rgba[2]);
533  fprintf(fp, " </Appearance>\n");
534  fprintf(fp, " </Shape>\n");
535  fprintf(fp, " </Transform>\n");
536  } // end if l
537  } // end for iv
538  } // end if dirty
539  } // end for loop on PView::list
540 
541  // triangles - colored triangles
542  // count all visible triangles of previous visited PView
543  _count = 0;
544  _ind = 0.;
545  fprintf(fp, " <Transform> \n");
546  fprintf(fp, " <Shape> \n");
547  fprintf(fp, " <Appearance> \n");
548  fprintf(fp, " <Material transparency='%g' ",
549  CTX::instance()->print.x3dTransparency);
550  fprintf(fp, " ambientIntensity='0.15' diffuseColor='.5 .5 .5' "
551  "emissiveColor='0 0 0' ");
552  fprintf(fp, " shininess='0.1' specularColor='.1 .1 .1' /> \n");
553  fprintf(fp, " </Appearance> \n");
554  fprintf(fp, " <IndexedTriangleSet solid='true' ccw='true' "
555  "colorPerVertex='true' \n ");
556  fprintf(
557  fp, " normalPerVertex='true' containerField='geometry' index=' ");
558  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
559  data = PView::list[ipv]->getData(true);
560  opt = PView::list[ipv]->getOptions();
561  if(!data->getDirty() && opt->visible) {
562  va = PView::list[ipv]->va_triangles;
563  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 3) {
564  if((!CTX::instance()->print.x3dRemoveInnerBorders && visible[_count]) ||
566  fprintf(fp, "%i %i %i ", _ind, _ind + 1, _ind + 2);
567  _ind += 3;
568  }
569  _count++;
570  }
571  } // enf if dirty
572  } // end loop on PView::list
573 
574  fprintf(fp, " ' > \n");
575  fprintf(fp, " <Coordinate point='");
576  _count = 0;
577  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
578  data = PView::list[ipv]->getData(true);
579  opt = PView::list[ipv]->getOptions();
580  if(!data->getDirty() && opt->visible) {
581  va = PView::list[ipv]->va_triangles;
582  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 3) {
583  if((!CTX::instance()->print.x3dRemoveInnerBorders && visible[_count]) ||
585  float *p0 = va->getVertexArray(3 * ipt);
586  float *p1 = va->getVertexArray(3 * (ipt + 1));
587  float *p2 = va->getVertexArray(3 * (ipt + 2));
588 
589  fprintf(fp, "%e %e %e %e %e %e %e %e %e ", p0[0], p0[1], p0[2], p1[0],
590  p1[1], p1[2], p2[0], p2[1], p2[2]);
591  }
592  _count++;
593  }
594  } // enf if dirty
595  } // end loop on PView::list
596  fprintf(fp, " '/> \n");
597 
598  fprintf(fp, " <Color color='");
599  _count = 0;
600  for(std::size_t ipv = 0; ipv < PView::list.size(); ipv++) {
601  data = PView::list[ipv]->getData(true);
602  opt = PView::list[ipv]->getOptions();
603  if(!data->getDirty() && opt->visible) {
604  va = PView::list[ipv]->va_triangles;
605  for(int ipt = 0; ipt < va->getNumVertices(); ipt += 3) {
606  if((!CTX::instance()->print.x3dRemoveInnerBorders && visible[_count]) ||
608  unsigned char *c0 = va->getColorArray(4 * ipt);
609  unsigned char *c1 = va->getColorArray(4 * (ipt + 1));
610  unsigned char *c2 = va->getColorArray(4 * (ipt + 2));
611  double rgba0[4], rgba1[4], rgba2[4];
612  UnsignedChar2rgba(c0, rgba0);
613  UnsignedChar2rgba(c1, rgba1);
614  UnsignedChar2rgba(c2, rgba2);
615 
616  fprintf(fp, "%g %g %g %g %g %g %g %g %g ", rgba0[0], rgba0[1],
617  rgba0[2], rgba1[0], rgba1[1], rgba1[2], rgba2[0], rgba2[1],
618  rgba2[2]);
619  }
620  _count++;
621  }
622  } // enf if dirty
623  } // end loop on PView::list
624  fprintf(fp, " '/>\n");
625  fprintf(fp, " </IndexedTriangleSet> \n");
626  fprintf(fp, " </Shape> \n");
627  fprintf(fp, " </Transform> \n");
628 
629  fprintf(fp, " </Scene>\n");
630  fprintf(fp, "</X3D>\n");
631  fclose(fp);
632  return true;
633 }
634 
635 static void writeX3DScale(FILE *fp, PView *p, double xmin, double ymin,
636  double width, double height, double tic,
637  int horizontal, double font_size)
638 {
639  // use adaptive data if available
640  PViewData *data = p->getData(true);
641  PViewOptions *opt = p->getOptions();
642 
643  if(opt->externalViewIndex >= 0) {
644  opt->tmpMin = opt->externalMin;
645  opt->tmpMax = opt->externalMax;
646  }
647  else if(opt->rangeType == PViewOptions::Custom) {
648  opt->tmpMin = opt->customMin;
649  opt->tmpMax = opt->customMax;
650  }
651  else if(opt->rangeType == PViewOptions::PerTimeStep) {
652  opt->tmpMin = data->getMin(opt->timeStep);
653  opt->tmpMax = data->getMax(opt->timeStep);
654  }
655  else {
656  opt->tmpMin = data->getMin();
657  opt->tmpMax = data->getMax();
658  }
659 
660  writeX3DScaleBar(fp, p, xmin, ymin, width, height, tic, horizontal);
661  writeX3DScaleValues(fp, p, xmin, ymin, width, height, tic, horizontal,
662  font_size);
663  writeX3DScaleLabel(fp, p, xmin, ymin, width, height, tic, horizontal,
664  font_size);
665 }
666 
667 static void writeX3DScaleBar(FILE *fp, PView *p, double xmin, double ymin,
668  double width, double height, double tic,
669  int horizontal)
670 {
671  PViewOptions *opt = p->getOptions();
672  double box = (horizontal ? width : height) / (opt->nbIso ? opt->nbIso : 1);
673 
674  for(int i = 0; i < opt->nbIso; i++) {
678  fprintf(fp, " <Shape> \n");
679  fprintf(
680  fp,
681  " <IndexedFaceSet colorPerVertex='true' normalPerVertex='true' "
682  "coordIndex='0 1 2 3 -1' solid='false' ccw='true' > \n");
683  fprintf(fp, " <Coordinate point='");
684  if(horizontal) {
685  fprintf(fp, "%e %e %e %e %e %e %e %e %e %e %e %e ", xmin + i * box,
686  ymin, 0., xmin + (i + 1) * box, ymin, 0., xmin + (i + 1) * box,
687  ymin + height, 0., xmin + i * box, ymin + height, 0.);
688  }
689  else {
690  fprintf(fp, "%e %e %e %e %e %e %e %e %e %e %e %e ", xmin,
691  ymin + i * box, 0., xmin + width, ymin + i * box, 0.,
692  xmin + width, ymin + (i + 1) * box, 0., xmin,
693  ymin + (i + 1) * box, 0.);
694  }
695  fprintf(fp, " '/> \n");
696 
699  double rgba[4] = {.5, .5, .5, 1.};
700  unsigned int col = opt->getColor(i, opt->nbIso);
701  unsignedInt2RGBA(col, rgba[0], rgba[1], rgba[2], rgba[3]);
702  fprintf(fp,
703  " <Color color=' %g %g %g %g %g %g %g %g %g %g %g "
704  "%g '/>\n",
705  rgba[0], rgba[1], rgba[2], rgba[0], rgba[1], rgba[2], rgba[0],
706  rgba[1], rgba[2], rgba[0], rgba[1], rgba[2]);
707  }
708  else if(opt->intervalsType == PViewOptions::Continuous) {
709  double dv = (opt->tmpMax - opt->tmpMin) / (opt->nbIso ? opt->nbIso : 1);
710  double v1 = opt->tmpMin + i * dv;
711  double v2 = opt->tmpMin + (i + 1) * dv;
712  unsigned int col1 = opt->getColor(v1, opt->tmpMin, opt->tmpMax, true);
713  unsigned int col2 = opt->getColor(v2, opt->tmpMin, opt->tmpMax, true);
714  double rgba1[4] = {.5, .5, .5, 1.};
715  double rgba2[4] = {.5, .5, .5, 1.};
716  unsignedInt2RGBA(col1, rgba1[0], rgba1[1], rgba1[2], rgba1[3]);
717  unsignedInt2RGBA(col2, rgba2[0], rgba2[1], rgba2[2], rgba2[3]);
718  fprintf(fp,
719  " <Color color=' %g %g %g %g %g %g %g %g %g %g %g "
720  "%g '/>\n",
721  rgba1[0], rgba1[1], rgba1[2], rgba2[0], rgba2[1], rgba2[2],
722  rgba2[0], rgba2[1], rgba2[2], rgba1[0], rgba1[1], rgba1[2]);
723  }
724  fprintf(fp, " </IndexedFaceSet> \n");
725  }
726  else {
727  fprintf(fp, " <Shape> \n");
728  fprintf(fp, " <IndexedLineSet colorPerVertex='true' "
729  "coordIndex='0 1 -1' > \n");
730  fprintf(fp, " <Coordinate point='");
731  if(horizontal) {
732  fprintf(fp, "%e %e %e %e %e %e ", xmin + box / 2. + i * box, ymin, 0.,
733  xmin + box / 2. + i * box, ymin + height, 0.);
734  }
735  else {
736  fprintf(fp, "%e %e %e %e %e %e ", xmin, ymin + box / 2. + i * box, 0.,
737  xmin + width, ymin + box / 2. + i * box, 0.);
738  }
739  fprintf(fp, " '/> \n");
740  double rgba[4] = {.5, .5, .5, 1.};
741  unsigned int col = opt->getColor(i, opt->nbIso);
742  unsignedInt2RGBA(col, rgba[0], rgba[1], rgba[2], rgba[3]);
743  fprintf(fp, " <Color color=' %g %g %g %g %g %g '/>\n", rgba[0],
744  rgba[1], rgba[2], rgba[0], rgba[1], rgba[2]);
745  fprintf(fp, " </IndexedLineSet> \n");
746  }
747  fprintf(fp, " </Shape> \n");
748  }
749 }
750 
751 static void writeX3DScaleValues(FILE *fp, PView *p, double xmin, double ymin,
752  double width, double height, double tic,
753  int horizontal, double font_size)
754 {
755  PViewOptions *opt = p->getOptions();
756  if(!opt->nbIso) return;
757  double font_h = font_size; // total font height
758  double font_a = font_size * 3. / 4.; // height above ref pt
759  char label[1024];
760  int nbv = opt->nbIso;
761  double maxw = 0.;
762  for(int i = 0; i < nbv + 1; i++) {
763  double v = opt->getScaleValue(i, nbv + 1, opt->tmpMin, opt->tmpMax);
764  sprintf(label, opt->format.c_str(), v);
765  maxw = max(maxw, strlen(label) * font_size * 3. / 4.);
766  }
767  double f = (opt->intervalsType == PViewOptions::Discrete ||
770  2 :
771  2.5;
772 
773  if(horizontal && width < nbv * maxw) {
774  if(width < f * maxw)
775  nbv = 1;
776  else
777  nbv = 2;
778  }
779  else if(!horizontal && height < nbv * font_h) {
780  if(height < f * font_h)
781  nbv = 1;
782  else
783  nbv = 2;
784  }
785  double box = (horizontal ? width : height) / opt->nbIso;
786  double vbox = (horizontal ? width : height) / nbv;
787 
788  // glColor4ubv((GLubyte *) & CTX::instance()->color.text);
789 
793  for(int i = 0; i < nbv + 1; i++) {
794  double v = opt->getScaleValue(i, nbv + 1, opt->tmpMin, opt->tmpMax);
795  sprintf(label, opt->format.c_str(), v);
796  if(horizontal) {
797  writeX3DStringCenter(fp, label, xmin + i * vbox, ymin + height + tic,
798  0., font_h);
799  }
800  else {
801  writeX3DStringCenter(fp, label, xmin + width + tic,
802  ymin + i * vbox - font_a / 3., 0., font_h);
803  }
804  }
805  }
806  else {
807  if(opt->nbIso > 2 && (nbv == 1 || nbv == 2)) {
808  vbox = (vbox * nbv - box) / nbv;
809  nbv++;
810  }
811  for(int i = 0; i < nbv; i++) {
812  double v = opt->getScaleValue(i, nbv, opt->tmpMin, opt->tmpMax);
813  sprintf(label, opt->format.c_str(), v);
814  if(horizontal) {
815  writeX3DStringCenter(fp, label, xmin + box / 2. + i * vbox,
816  ymin + height + tic, 0., font_h);
817  }
818  else {
819  writeX3DStringCenter(fp, label, xmin + width + tic,
820  ymin + box / 2. + i * vbox - font_a / 3., 0.,
821  font_h);
822  }
823  }
824  }
825 }
826 
827 static void writeX3DScaleLabel(FILE *fp, PView *p, double xmin, double ymin,
828  double width, double height, double tic,
829  int horizontal, double font_size)
830 {
831  PViewOptions *opt = p->getOptions();
832  PViewData *data;
833 
834  // requested by Laurent: but is this really what we should be doing?
835  if(opt->externalViewIndex >= 0 &&
836  opt->externalViewIndex < (int)PView::list.size())
837  data = PView::list[opt->externalViewIndex]->getData();
838  else
839  data = p->getData();
840  double font_h = font_size;
841  font_h *= 1.2;
842  char label[1024];
843  int nt = data->getNumTimeSteps();
844  if((opt->showTime == 1 && nt > 1) || opt->showTime == 2) {
845  char tmp[256];
846  sprintf(tmp, opt->format.c_str(), data->getTime(opt->timeStep));
847  sprintf(label, "%s (%s)", data->getName().c_str(), tmp);
848  }
849  else if((opt->showTime == 3 && nt > 1) || opt->showTime == 4) {
850  sprintf(label, "%s (%d/%d)", data->getName().c_str(), opt->timeStep,
851  data->getNumTimeSteps() - 1);
852  }
853  else
854  sprintf(label, "%s", data->getName().c_str());
855  if(horizontal) {
856  writeX3DStringCenter(fp, label, xmin + width / 2.,
857  ymin + height + tic + .9 * font_h, 0., font_h);
858  }
859  else {
860  writeX3DStringCenter(fp, label, xmin, ymin - 2 * font_h, 0., font_h);
861  }
862 }
863 
864 static void writeX3DStringCenter(FILE *fp, char *label, double x, double y,
865  double z, double font_size)
866 {
867  fprintf(fp, " <Transform translation='%g %g %g'> \n", x, y, 0.);
868  fprintf(fp, " <Shape> \n");
869  fprintf(fp, " <Text string='\"%s\"'>\n", label);
870  fprintf(
871  fp,
872  " <FontStyle justify='\"MIDDLE\" \"MIDDLE\"' size=' %d '/> \n",
873  (int)font_size);
874  fprintf(fp, " </Text>\n");
875  fprintf(fp, " <Appearance>\n");
876  fprintf(fp, " <Material diffuseColor='0. 0. 0. '/>\n");
877  fprintf(fp, " </Appearance>\n");
878  fprintf(fp, " </Shape>\n");
879  fprintf(fp, " </Transform> \n");
880 }
writeX3DScaleValues
static void writeX3DScaleValues(FILE *fp, PView *p, double xmin, double ymin, double width, double height, double tic, int horizontal, double font_size)
Definition: PViewX3D.cpp:751
PViewOptions::timeStep
int timeStep
Definition: PViewOptions.h:61
PViewOptions::externalMin
double externalMin
Definition: PViewOptions.h:47
SBoundingBox3d::diag
double diag() const
Definition: SBoundingBox3d.h:93
PView
Definition: PView.h:27
CTX::x3dCompatibility
int x3dCompatibility
Definition: Context.h:352
PViewData::getNumTimeSteps
virtual int getNumTimeSteps()=0
CTX::horizontalScales
int horizontalScales
Definition: Context.h:316
OS.h
UnsignedChar2rgba
static void UnsignedChar2rgba(unsigned char *glc, double *rgba)
Definition: PViewX3D.h:14
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
CTX::x3dRemoveInnerBorders
int x3dRemoveInnerBorders
Definition: Context.h:352
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
PViewOptions::tmpMax
double tmpMax
Definition: PViewOptions.h:47
VertexArray.h
box
Definition: gl2gif.cpp:311
SPoint3
Definition: SPoint3.h:14
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
VertexArray
Definition: VertexArray.h:151
PView.h
TriangleToSort::_ppv
PView * _ppv
Definition: PViewX3D.h:52
PViewOptions::showScale
int showScale
Definition: PViewOptions.h:58
CTX::x3dTransparency
double x3dTransparency
Definition: Context.h:353
TriangleToSort::zmax
float zmax
Definition: PViewX3D.h:56
VertexArray::getNumVertices
int getNumVertices()
Definition: VertexArray.h:174
GmshMessage.h
PViewData.h
PViewOptions::nbIso
int nbIso
Definition: PViewOptions.h:54
compare_xmax_triangle
bool compare_xmax_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:54
TriangleToSort::_globalIndex
int _globalIndex
Definition: PViewX3D.h:54
Fopen
FILE * Fopen(const char *f, const char *mode)
Definition: OS.cpp:273
TriangleToSort::ymin
float ymin
Definition: PViewX3D.h:55
PViewOptions::customMax
double customMax
Definition: PViewOptions.h:47
PViewOptions::externalViewIndex
int externalViewIndex
Definition: PViewOptions.h:73
PViewData::getDirty
virtual bool getDirty()
Definition: PViewData.h:62
TriangleToSort
Definition: PViewX3D.h:50
SBoundingBox3d::center
SPoint3 center() const
Definition: SBoundingBox3d.h:92
CTX::instance
static CTX * instance()
Definition: Context.cpp:122
PViewOptions::Numeric
@ Numeric
Definition: PViewOptions.h:19
SPoint3::x
double x(void) const
Definition: SPoint3.h:125
PViewData::getTime
virtual double getTime(int step)
Definition: PViewData.h:94
PViewData::getMax
virtual double getMax(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PViewData::getBoundingBox
virtual SBoundingBox3d getBoundingBox(int step=-1)=0
writeX3DScaleLabel
static void writeX3DScaleLabel(FILE *fp, PView *p, double xmin, double ymin, double width, double height, double tic, int horizontal, double font_size)
Definition: PViewX3D.cpp:827
SBoundingBox3d.h
VertexArray::getColorArray
unsigned char * getColorArray(int i=0)
Definition: VertexArray.h:192
writeX3DScaleBar
static void writeX3DScaleBar(FILE *fp, PView *p, double xmin, double ymin, double width, double height, double tic, int horizontal)
Definition: PViewX3D.cpp:667
PView::getData
PViewData * getData(bool useAdaptiveIfAvailable=false)
Definition: PView.cpp:233
PViewOptions::rangeType
int rangeType
Definition: PViewOptions.h:59
PViewOptions::intervalsType
int intervalsType
Definition: PViewOptions.h:54
PViewOptions.h
PViewDataList.h
PViewOptions::format
std::string format
Definition: PViewOptions.h:42
compare_zmax_triangle
bool compare_zmax_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:66
compare_ymax_triangle
bool compare_ymax_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:60
PView::writeX3D
static bool writeX3D(const std::string &fileName)
Definition: PViewX3D.cpp:72
SPoint3::y
double y(void) const
Definition: SPoint3.h:127
writeX3DScale
static void writeX3DScale(FILE *fp, PView *p, double xmin, double ymin, double width, double height, double tic, int horizontal, double font_size)
Definition: PViewX3D.cpp:635
PViewDataGModel.h
compare_ymin_triangle
bool compare_ymin_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:42
almostEqual
static bool almostEqual(double x, double y)
Definition: PViewX3D.cpp:31
PViewOptions::customMin
double customMin
Definition: PViewOptions.h:47
PViewOptions::showTime
int showTime
Definition: PViewOptions.h:58
CTX::print
struct CTX::@2 print
PViewData
Definition: PViewData.h:29
compare_zmin_triangle
bool compare_zmin_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:48
StringUtils.h
TriangleToSort::xmax
float xmax
Definition: PViewX3D.h:56
VertexArray::getVertexArray
float * getVertexArray(int i=0)
Definition: VertexArray.h:182
PViewOptions::type
int type
Definition: PViewOptions.h:40
PViewOptions::PerTimeStep
@ PerTimeStep
Definition: PViewOptions.h:37
Context.h
PViewOptions::autoPosition
int autoPosition
Definition: PViewOptions.h:40
PViewOptions::externalMax
double externalMax
Definition: PViewOptions.h:47
PViewOptions::Discrete
@ Discrete
Definition: PViewOptions.h:19
unsignedInt2RGBA
static void unsignedInt2RGBA(unsigned int &color, double &r, double &g, double &b, double &a)
Definition: PViewX3D.h:22
compare_xmin_triangle
bool compare_xmin_triangle(const TriangleToSort *first, const TriangleToSort *second)
Definition: PViewX3D.cpp:36
std
Definition: picojson.h:1135
z
const double z
Definition: GaussQuadratureQuad.cpp:56
TriangleToSort::ymax
float ymax
Definition: PViewX3D.h:56
PViewOptions::lineType
int lineType
Definition: PViewOptions.h:68
PViewOptions
Definition: PViewOptions.h:16
PView::getOptions
PViewOptions * getOptions()
Definition: PView.h:81
PViewOptions::vectorType
int vectorType
Definition: PViewOptions.h:60
PViewOptions::Continuous
@ Continuous
Definition: PViewOptions.h:19
SPoint3::z
double z(void) const
Definition: SPoint3.h:129
PViewOptions::getScaleValue
double getScaleValue(int iso, int numIso, double min, double max)
Definition: PViewOptions.cpp:33
PViewData::getMin
virtual double getMin(int step=-1, bool onlyVisible=false, int tensorRep=0, int forceNumComponents=0, int componentMap[9]=nullptr)=0
PViewData::getNumElements
virtual int getNumElements(int step=-1, int ent=-1)
Definition: PViewData.h:131
writeX3DStringCenter
static void writeX3DStringCenter(FILE *fp, char *label, double x, double y, double z, double font_size)
Definition: PViewX3D.cpp:864
PViewData::getName
virtual std::string getName()
Definition: PViewData.h:70
PViewOptions::tmpMin
double tmpMin
Definition: PViewOptions.h:47
PViewOptions::Custom
@ Custom
Definition: PViewOptions.h:37
PView::list
static std::vector< PView * > list
Definition: PView.h:112
CTX::x3dPrecision
double x3dPrecision
Definition: Context.h:353
PViewOptions::getColor
unsigned int getColor(double val, double min, double max, bool forceLinear=false, int numColors=-1)
Definition: PViewOptions.cpp:86
c1
const double c1
Definition: GaussQuadratureHex.cpp:16
TriangleToSort::_index
int _index
Definition: PViewX3D.h:53
TriangleToSort::zmin
float zmin
Definition: PViewX3D.h:55
TriangleToSort::xmin
float xmin
Definition: PViewX3D.h:55
SBoundingBox3d
Definition: SBoundingBox3d.h:21
PViewX3D.h
PViewOptions::visible
int visible
Definition: PViewOptions.h:54
PViewOptions::Plot3D
@ Plot3D
Definition: PViewOptions.h:18
scale
static void scale(std::vector< double > &x, double s)
Definition: ConjugateGradients.cpp:21