source: trunk/geant4/visualization/OpenGL/src/G4OpenGLViewer.cc @ 715

Last change on this file since 715 was 712, checked in by garnier, 16 years ago

correction du ticket #119

  • Property svn:mime-type set to text/cpp
File size: 26.9 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4OpenGLViewer.cc,v 1.34 2007/05/24 18:27:13 allison Exp $
28// GEANT4 tag $Name: HEAD $
29//
30//
31// Andrew Walkden  27th March 1996
32// OpenGL view - opens window, hard copy, etc.
33
34#ifdef G4VIS_BUILD_OPENGL_DRIVER
35
36#include "G4ios.hh"
37#include "G4OpenGLViewer.hh"
38#include "G4OpenGLSceneHandler.hh"
39#include "G4OpenGLTransform3D.hh"
40
41#include "G4Scene.hh"
42#include "G4VisExtent.hh"
43#include "G4LogicalVolume.hh"
44#include "G4VSolid.hh"
45#include "G4Point3D.hh"
46#include "G4Normal3D.hh"
47#include "G4Plane3D.hh"
48#include "G4AttHolder.hh"
49#include "G4AttCheck.hh"
50#include <sstream>
51
52static const char* gouraudtriangleEPS[] =
53{
54  "/bd{bind def}bind def /triangle { aload pop   setrgbcolor  aload pop 5 3",
55  "roll 4 2 roll 3 2 roll exch moveto lineto lineto closepath fill } bd",
56  "/computediff1 { 2 copy sub abs threshold ge {pop pop pop true} { exch 2",
57  "index sub abs threshold ge { pop pop true} { sub abs threshold ge } ifelse",
58  "} ifelse } bd /computediff3 { 3 copy 0 get 3 1 roll 0 get 3 1 roll 0 get",
59  "computediff1 {true} { 3 copy 1 get 3 1 roll 1 get 3 1 roll 1 get",
60  "computediff1 {true} { 3 copy 2 get 3 1 roll  2 get 3 1 roll 2 get",
61  "computediff1 } ifelse } ifelse } bd /middlecolor { aload pop 4 -1 roll",
62  "aload pop 4 -1 roll add 2 div 5 1 roll 3 -1 roll add 2 div 3 1 roll add 2",
63  "div 3 1 roll exch 3 array astore } bd /gouraudtriangle { computediff3 { 4",
64  "-1 roll aload 7 1 roll 6 -1 roll pop 3 -1 roll pop add 2 div 3 1 roll add",
65  "2 div exch 3 -1 roll aload 7 1 roll exch pop 4 -1 roll pop add 2 div 3 1",
66  "roll add 2 div exch 3 -1 roll aload 7 1 roll pop 3 -1 roll pop add 2 div 3",
67  "1 roll add 2 div exch 7 3 roll 10 -3 roll dup 3 index middlecolor 4 1 roll",
68  "2 copy middlecolor 4 1 roll 3 copy pop middlecolor 4 1 roll 13 -1 roll",
69  "aload pop 17 index 6 index 15 index 19 index 6 index 17 index 6 array",
70  "astore 10 index 10 index 14 index gouraudtriangle 17 index 5 index 17",
71  "index 19 index 5 index 19 index 6 array astore 10 index 9 index 13 index",
72  "gouraudtriangle 13 index 16 index 5 index 15 index 18 index 5 index 6",
73  "array astore 12 index 12 index 9 index gouraudtriangle 17 index 16 index",
74  "15 index 19 index 18 index 17 index 6 array astore 10 index 12 index 14",
75  "index gouraudtriangle 18 {pop} repeat } { aload pop 5 3 roll aload pop 7 3",
76  "roll aload pop 9 3 roll 4 index 6 index 4 index add add 3 div 10 1 roll 7",
77  "index 5 index 3 index add add 3 div 10 1 roll 6 index 4 index 2 index add",
78  "add 3 div 10 1 roll 9 {pop} repeat 3 array astore triangle } ifelse } bd",
79  NULL
80};
81
82G4OpenGLViewer::G4OpenGLViewer (G4OpenGLSceneHandler& scene):
83G4VViewer (scene, -1),
84pointSize (0),
85print_colour (true),
86vectored_ps (true),
87fOpenGLSceneHandler(scene),
88background (G4Colour(0.,0.,0.)),
89transparency_enabled (true),
90antialiasing_enabled (false),
91haloing_enabled (false),
92fStartTime(-G4OPENGL_DBL_MAX),
93fEndTime(G4OPENGL_DBL_MAX),
94fFadeFactor(0.),
95fDisplayHeadTime(false),
96fDisplayHeadTimeX(-0.9),
97fDisplayHeadTimeY(-0.9),
98fDisplayHeadTimeSize(24.),
99fDisplayHeadTimeRed(0.),
100fDisplayHeadTimeGreen(1.),
101fDisplayHeadTimeBlue(1.),
102fDisplayLightFront(false),
103fDisplayLightFrontX(0.),
104fDisplayLightFrontY(0.),
105fDisplayLightFrontZ(0.),
106fDisplayLightFrontT(0.),
107fDisplayLightFrontRed(0.),
108fDisplayLightFrontGreen(1.),
109fDisplayLightFrontBlue(0.)
110{
111  // Make changes to view parameters for OpenGL...
112  fVP.SetAutoRefresh(true);
113  fDefaultVP.SetAutoRefresh(true);
114
115  //  glClearColor (0.0, 0.0, 0.0, 0.0);
116  //  glClearDepth (1.0);
117  //  glDisable (GL_BLEND);
118  //  glDisable (GL_LINE_SMOOTH);
119  //  glDisable (GL_POLYGON_SMOOTH);
120
121  strcpy (print_string, "G4OpenGL.eps");
122}
123
124G4OpenGLViewer::~G4OpenGLViewer () {}
125
126void G4OpenGLViewer::InitializeGLView ()
127{
128  glClearColor (0.0, 0.0, 0.0, 0.0);
129  glClearDepth (1.0);
130  glDisable (GL_BLEND);
131  glDisable (GL_LINE_SMOOTH);
132  glDisable (GL_POLYGON_SMOOTH);
133
134
135void G4OpenGLViewer::ClearView () {
136  glClearColor (background.GetRed(),
137                background.GetGreen(),
138                background.GetBlue(),
139                1.);
140  glClearDepth (1.0);
141  //Below line does not compile with Mesa includes.
142  //glClear (GL_COLOR_BUFFER_BIT | GL_DEPTH_BUFFER_BIT);
143  glClear (GL_COLOR_BUFFER_BIT);
144  glClear (GL_DEPTH_BUFFER_BIT);
145  glClear (GL_STENCIL_BUFFER_BIT);
146  glFlush ();
147}
148
149void G4OpenGLViewer::SetView () {
150 
151  // Calculates view representation based on extent of object being
152  // viewed and (initial) viewpoint.  (Note: it can change later due
153  // to user interaction via visualization system's GUI.)
154 
155  // Lighting.
156  GLfloat lightPosition [4];
157  lightPosition [0] = fVP.GetActualLightpointDirection().x();
158  lightPosition [1] = fVP.GetActualLightpointDirection().y();
159  lightPosition [2] = fVP.GetActualLightpointDirection().z();
160  lightPosition [3] = 0.;
161  // Light position is "true" light direction, so must come after gluLookAt.
162  GLfloat ambient [] = { 0.2, 0.2, 0.2, 1.};
163  GLfloat diffuse [] = { 0.8, 0.8, 0.8, 1.};
164  glEnable (GL_LIGHT0);
165  glLightfv (GL_LIGHT0, GL_AMBIENT, ambient);
166  glLightfv (GL_LIGHT0, GL_DIFFUSE, diffuse);
167 
168  // Get radius of scene, etc.
169  // Note that this procedure properly takes into account zoom, dolly and pan.
170  const G4Point3D targetPoint
171    = fSceneHandler.GetScene()->GetStandardTargetPoint()
172    + fVP.GetCurrentTargetPoint ();
173  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
174  if(radius<=0.) radius = 1.;
175  const G4double cameraDistance = fVP.GetCameraDistance (radius);
176  const G4Point3D cameraPosition =
177    targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
178  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
179  const GLdouble pfar    = fVP.GetFarDistance  (cameraDistance, pnear, radius);
180  const GLdouble right  = fVP.GetFrontHalfHeight (pnear, radius);
181  const GLdouble left   = -right;
182  const GLdouble bottom = left;
183  const GLdouble top    = right;
184 
185  glMatrixMode (GL_PROJECTION); // set up Frustum.
186  glLoadIdentity();
187
188  const G4Vector3D scale = fVP.GetScaleFactor();
189  glScaled(scale.x(),scale.y(),scale.z());
190 
191  if (fVP.GetFieldHalfAngle() == 0.) {
192    glOrtho (left, right, bottom, top, pnear, pfar);
193  }
194  else {
195    glFrustum (left, right, bottom, top, pnear, pfar);
196  }
197 
198  printf("Set view to : %f %f %f %f %f %f \n",left, right, bottom, top, pnear, pfar);
199  glMatrixMode (GL_MODELVIEW); // apply further transformations to scene.
200  glLoadIdentity();
201 
202  const G4Normal3D& upVector = fVP.GetUpVector (); 
203  G4Point3D gltarget;
204  if (cameraDistance > 1.e-6 * radius) {
205    gltarget = targetPoint;
206  }
207  else {
208    gltarget = targetPoint - radius * fVP.GetViewpointDirection().unit();
209  }
210
211  const G4Point3D& pCamera = cameraPosition;  // An alias for brevity.
212  gluLookAt (pCamera.x(),  pCamera.y(),  pCamera.z(),       // Viewpoint.
213             gltarget.x(), gltarget.y(), gltarget.z(),      // Target point.
214             upVector.x(), upVector.y(), upVector.z());     // Up vector.
215 
216  // Light position is "true" light direction, so must come after gluLookAt.
217  glLightfv (GL_LIGHT0, GL_POSITION, lightPosition);
218
219  // OpenGL no longer seems to reconstruct clipped edges, so, when the
220  // BooleanProcessor is up to it, abandon this and use generic
221  // clipping in G4OpenGLSceneHandler::CreateSectionPolyhedron.  Also,
222  // force kernel visit on change of clipping plane in
223  // G4OpenGLStoredViewer::CompareForKernelVisit.
224  if (fVP.IsSection () ) {  // pair of back to back clip planes.
225    const G4Plane3D& s = fVP.GetSectionPlane ();
226    double sArray[4];
227    sArray[0] = s.a();
228    sArray[1] = s.b();
229    sArray[2] = s.c();
230    sArray[3] = s.d() + radius * 1.e-05;
231    glClipPlane (GL_CLIP_PLANE0, sArray);
232    glEnable (GL_CLIP_PLANE0);
233    sArray[0] = -s.a();
234    sArray[1] = -s.b();
235    sArray[2] = -s.c();
236    sArray[3] = -s.d() + radius * 1.e-05;
237    glClipPlane (GL_CLIP_PLANE1, sArray);
238    glEnable (GL_CLIP_PLANE1);
239  } else {
240    glDisable (GL_CLIP_PLANE0);
241    glDisable (GL_CLIP_PLANE1);
242  }
243
244  const G4Planes& cutaways = fVP.GetCutawayPlanes();
245  size_t nPlanes = cutaways.size();
246  if (fVP.IsCutaway() &&
247      fVP.GetCutawayMode() == G4ViewParameters::cutawayIntersection &&
248      nPlanes > 0) {
249    double a[4];
250    a[0] = cutaways[0].a();
251    a[1] = cutaways[0].b();
252    a[2] = cutaways[0].c();
253    a[3] = cutaways[0].d();
254    glClipPlane (GL_CLIP_PLANE2, a);
255    glEnable (GL_CLIP_PLANE2);
256    if (nPlanes > 1) {
257      a[0] = cutaways[1].a();
258      a[1] = cutaways[1].b();
259      a[2] = cutaways[1].c();
260      a[3] = cutaways[1].d();
261      glClipPlane (GL_CLIP_PLANE3, a);
262      glEnable (GL_CLIP_PLANE3);
263    }
264    if (nPlanes > 2) {
265      a[0] = cutaways[2].a();
266      a[1] = cutaways[2].b();
267      a[2] = cutaways[2].c();
268      a[3] = cutaways[2].d();
269      glClipPlane (GL_CLIP_PLANE4, a);
270      glEnable (GL_CLIP_PLANE4);
271    }
272  } else {
273    glDisable (GL_CLIP_PLANE2);
274    glDisable (GL_CLIP_PLANE3);
275    glDisable (GL_CLIP_PLANE4);
276  }
277
278  // Background.
279  background = fVP.GetBackgroundColour ();
280
281}
282
283void G4OpenGLViewer::HaloingFirstPass () {
284 
285  //To perform haloing, first Draw all information to the depth buffer
286  //alone, using a chunky line width, and then Draw all info again, to
287  //the colour buffer, setting a thinner line width an the depth testing
288  //function to less than or equal, so if two lines cross, the one
289  //passing behind the other will not pass the depth test, and so not
290  //get rendered either side of the infront line for a short distance.
291
292  //First, disable writing to the colo(u)r buffer...
293  glColorMask (GL_FALSE, GL_FALSE, GL_FALSE, GL_FALSE);
294
295  //Now enable writing to the depth buffer...
296  glDepthMask (GL_TRUE);
297  glDepthFunc (GL_LESS);
298  glClearDepth (1.0);
299
300  //Finally, set the line width to something wide...
301  glLineWidth (3.0);
302
303}
304
305void G4OpenGLViewer::HaloingSecondPass () {
306
307  //And finally, turn the colour buffer back on with a sesible line width...
308  glColorMask (GL_TRUE, GL_TRUE, GL_TRUE, GL_TRUE);
309  glDepthFunc (GL_LEQUAL);
310  glLineWidth (1.0);
311
312}
313
314void G4OpenGLViewer::Pick(GLdouble x, GLdouble y)
315{
316  //G4cout << "X: " << x << ", Y: " << y << G4endl;
317  const G4int BUFSIZE = 512;
318  GLuint selectBuffer[BUFSIZE];
319  glSelectBuffer(BUFSIZE, selectBuffer);
320  glRenderMode(GL_SELECT);
321  glInitNames();
322  glPushName(0);
323  glMatrixMode(GL_PROJECTION);
324  G4double currentProjectionMatrix[16];
325  glGetDoublev(GL_PROJECTION_MATRIX, currentProjectionMatrix);
326  glPushMatrix();
327  glLoadIdentity();
328  GLint viewport[4];
329  glGetIntegerv(GL_VIEWPORT, viewport);
330  // Define 5x5 pixel pick area
331  gluPickMatrix(x, viewport[3] - y, 5., 5., viewport);
332  glMultMatrixd(currentProjectionMatrix);
333  glMatrixMode(GL_MODELVIEW);
334  DrawView();
335  GLint hits = glRenderMode(GL_RENDER);
336  if (hits < 0)
337    G4cout << "Too many hits.  Zoom in to reduce overlaps." << G4cout;
338  else if (hits > 0) {
339    //G4cout << hits << " hit(s)" << G4endl;
340    GLuint* p = selectBuffer;
341    for (GLint i = 0; i < hits; ++i) {
342      GLuint nnames = *p++;
343      *p++; //OR GLuint zmin = *p++;
344      *p++; //OR GLuint zmax = *p++;
345      //G4cout << "Hit " << i << ": " << nnames << " names"
346      //     << "\nzmin: " << zmin << ", zmax: " << zmax << G4endl;
347      for (GLuint j = 0; j < nnames; ++j) {
348        GLuint name = *p++;
349        //G4cout << "Name " << j << ": PickName: " << name << G4endl;
350        std::map<GLuint, G4AttHolder*>::iterator iter =
351          fOpenGLSceneHandler.fPickMap.find(name);
352        if (iter != fOpenGLSceneHandler.fPickMap.end()) {
353          G4AttHolder* attHolder = iter->second;
354          if(attHolder && attHolder->GetAttDefs().size()) {
355            for (size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
356              G4cout << G4AttCheck(attHolder->GetAttValues()[i],
357                                   attHolder->GetAttDefs()[i]);
358            }
359          }
360        }
361      }
362      G4cout << G4endl;
363    }
364  }
365  glMatrixMode(GL_PROJECTION);
366  glPopMatrix();
367  glMatrixMode(GL_MODELVIEW);
368}
369
370void G4OpenGLViewer::print() {
371
372  // Print vectored PostScript
373 
374  G4int size = 5000000;
375
376  GLfloat* feedback_buffer;
377  GLint returned;
378  FILE* file;
379 
380  feedback_buffer = new GLfloat[size];
381  glFeedbackBuffer (size, GL_3D_COLOR, feedback_buffer);
382  glRenderMode (GL_FEEDBACK);
383 
384  DrawView();
385  returned = glRenderMode (GL_RENDER);
386 
387  if (print_string) {
388    file = fopen (print_string, "w");
389    if (file) {
390      spewWireframeEPS (file, returned, feedback_buffer, "rendereps");
391    } else {
392      printf("Could not open %s\n", print_string);
393    }
394  } else {
395    printBuffer (returned, feedback_buffer);
396  }
397
398  delete[] feedback_buffer;
399}
400
401void G4OpenGLViewer::print3DcolorVertex(GLint size, GLint * count, GLfloat * buffer)
402{
403  G4int i;
404
405  printf("  ");
406  for (i = 0; i < 7; i++) {
407    printf("%4.2f ", buffer[size - (*count)]);
408    *count = *count - 1;
409  }
410  printf("\n");
411}
412
413void G4OpenGLViewer::spewWireframeEPS (FILE* file, GLint size, GLfloat* buffer, const char* cr) {
414
415  GLfloat EPS_GOURAUD_THRESHOLD=0.1;
416
417  GLfloat clearColor[4], viewport[4];
418  GLfloat lineWidth;
419  G4int i;
420
421  glGetFloatv (GL_VIEWPORT, viewport);
422  glGetFloatv (GL_COLOR_CLEAR_VALUE, clearColor);
423  glGetFloatv (GL_LINE_WIDTH, &lineWidth);
424  glGetFloatv (GL_POINT_SIZE, &pointSize);
425
426  fputs ("%!PS-Adobe-2.0 EPSF-2.0\n", file);
427  fprintf (file, "%%%%Creator: %s (using OpenGL feedback)\n", cr);
428  fprintf (file, "%%%%BoundingBox: %g %g %g %g\n", viewport[0], viewport[1], viewport[2], viewport[3]);
429  fputs ("%%EndComments\n", file);
430  fputs ("\n", file);
431  fputs ("gsave\n", file);
432  fputs ("\n", file);
433
434  fputs ("% the gouraudtriangle PostScript fragment below is free\n", file);
435  fputs ("% written by Frederic Delhoume (delhoume@ilog.fr)\n", file);
436  fprintf (file, "/threshold %g def\n", EPS_GOURAUD_THRESHOLD);
437  for (i=0; gouraudtriangleEPS[i]; i++) {
438    fprintf (file, "%s\n", gouraudtriangleEPS[i]);
439  }
440
441  fprintf(file, "\n%g setlinewidth\n", lineWidth);
442 
443  fprintf (file, "%g %g %g setrgbcolor\n", clearColor[0], clearColor[1], clearColor[2]);
444  fprintf (file, "%g %g %g %g rectfill\n\n", viewport[0], viewport[1], viewport[2], viewport[3]);
445
446  spewSortedFeedback (file, size, buffer);
447
448  fputs ("grestore\n\n", file);
449  fputs ("showpage\n", file);
450
451  fclose(file);
452}
453
454void G4OpenGLViewer::printBuffer (GLint size, GLfloat* buffer) {
455
456  GLint count;
457  G4int token, nvertices;
458
459  count=size;
460  while(count) {
461    token=G4int (buffer[size-count]);
462    count--;
463    switch (token) {
464
465    case GL_PASS_THROUGH_TOKEN:
466      printf ("GL_PASS_THROUGH_TOKEN\n");
467      printf ("  %4.2f\n", buffer[size-count]);
468      count--;
469      break;
470
471    case GL_POINT_TOKEN:
472      printf ("GL_POINT_TOKEN\n");
473      print3DcolorVertex (size, &count, buffer);
474      break;
475
476    case GL_LINE_TOKEN:
477      printf ("GL_LINE_TOKEN\n");
478      print3DcolorVertex (size, &count, buffer);
479      print3DcolorVertex (size, &count, buffer);
480      break;
481     
482    case GL_LINE_RESET_TOKEN:
483      printf ("GL_LINE_RESET_TOKEN\n");
484      print3DcolorVertex (size, &count, buffer);
485      print3DcolorVertex (size, &count, buffer);
486      break;
487
488    case GL_POLYGON_TOKEN:
489      printf ("GL_POLYGON_TOKEN\n");
490      nvertices=G4int (buffer[size-count]);
491      count--;
492      for (; nvertices>0; nvertices--) {
493        print3DcolorVertex (size, &count, buffer);
494      }
495    }
496  }
497}
498
499G4float* G4OpenGLViewer::spewPrimitiveEPS (FILE* file, GLfloat* loc) {
500 
501  G4int token;
502  G4int nvertices, i;
503  GLfloat red, green, blue, intensity;
504  G4int smooth;
505  GLfloat dx, dy, dr, dg, db, absR, absG, absB, colormax;
506  G4int steps;
507  Feedback3Dcolor *vertex;
508  GLfloat xstep(0.), ystep(0.), rstep(0.), gstep(0.), bstep(0.);
509  GLfloat xnext(0.), ynext(0.), rnext(0.), gnext(0.), bnext(0.), distance(0.);
510
511  token=G4int (*loc);
512  loc++;
513  switch (token) {
514  case GL_LINE_RESET_TOKEN:
515  case GL_LINE_TOKEN:
516    vertex=(Feedback3Dcolor*)loc;
517    dr=vertex[1].red - vertex[0].red;
518    dg=vertex[1].green - vertex[0].green;
519    db=vertex[1].blue - vertex[0].blue;
520
521    if (!print_colour) {
522      dr+=(dg+db);
523      dr/=3.0;
524      dg=dr;
525      db=dr;
526    }
527
528    if (dr!=0 || dg!=0 || db!=0) {
529      dx=vertex[1].x - vertex[0].x;
530      dy=vertex[1].y - vertex[0].y;
531      distance=std::sqrt(dx*dx + dy*dy);
532
533      absR=std::fabs(dr);
534      absG=std::fabs(dg);
535      absB=std::fabs(db);
536
537      #define Max(a, b) (((a)>(b))?(a):(b))
538
539      #define EPS_SMOOTH_LINE_FACTOR 0.06
540
541      colormax=Max(absR, Max(absG, absB));
542      steps=Max(1, G4int (colormax*distance*EPS_SMOOTH_LINE_FACTOR));
543     
544      xstep=dx/steps;
545      ystep=dy/steps;
546
547      rstep=dr/steps;
548      gstep=dg/steps;
549      bstep=db/steps;
550
551      xnext=vertex[0].x;
552      ynext=vertex[0].y;
553      rnext=vertex[0].red;
554      gnext=vertex[0].green;
555      bnext=vertex[0].blue;
556
557      if (!print_colour) {
558        rnext+=(gnext+bnext);
559        rnext/=3.0;
560        gnext=rnext;
561        bnext=rnext;
562      }
563
564      xnext -= xstep/2.0;
565      ynext -= ystep/2.0;
566      rnext -= rstep/2.0;
567      gnext -= gstep/2.0;
568      bnext -= bstep/2.0;
569    } else {
570      steps=0;
571    }
572    if (print_colour) {
573      fprintf (file, "%g %g %g setrgbcolor\n",
574               vertex[0].red, vertex[0].green, vertex[0].blue);
575    } else {
576      intensity = (vertex[0].red + vertex[0].green + vertex[0].blue) / 3.0;
577      fprintf (file, "%g %g %g setrgbcolor\n",
578               intensity, intensity, intensity);
579    }     
580    fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
581
582    for (i=0; i<steps; i++) {
583
584      xnext += xstep;
585      ynext += ystep;
586      rnext += rstep;
587      gnext += gstep;
588      bnext += bstep;
589
590      fprintf (file, "%g %g lineto stroke\n", xnext, ynext);
591      fprintf (file, "%g %g %g setrgbcolor\n", rnext, gnext, bnext);
592      fprintf (file, "%g %g moveto\n", xnext, ynext);
593    }
594    fprintf (file, "%g %g lineto stroke\n", vertex[1].x, vertex[1].y);
595
596    loc += 14;
597    break;
598
599  case GL_POLYGON_TOKEN:
600    nvertices = G4int (*loc);
601    loc++;
602    vertex=(Feedback3Dcolor*)loc;
603    if (nvertices>0) {
604      red=vertex[0].red;
605      green=vertex[0].green;
606      blue=vertex[0].blue;
607      smooth=0;
608     
609      if (!print_colour) {
610        red+=(green+blue);
611        red/=3.0;
612        green=red;
613        blue=red;
614      }
615     
616      if (print_colour) {
617        for (i=1; i<nvertices; i++) {
618          if (red!=vertex[i].red || green!=vertex[i].green || blue!=vertex[i].blue) {
619            smooth=1;
620            break;
621          }
622        }
623      } else {
624        for (i=1; i<nvertices; i++) {
625          intensity = vertex[i].red + vertex[i].green + vertex[i].blue;
626          intensity/=3.0;
627          if (red!=intensity) {
628            smooth=1;
629            break;
630          }
631        }
632      }
633
634      if (smooth) {
635        G4int triOffset;
636        for (i=0; i<nvertices-2; i++) {
637          triOffset = i*7;
638          fprintf (file, "[%g %g %g %g %g %g]",
639                   vertex[0].x, vertex[i+1].x, vertex[i+2].x,
640                   vertex[0].y, vertex[i+1].y, vertex[i+2].y);
641          if (print_colour) {
642            fprintf (file, " [%g %g %g] [%g %g %g] [%g %g %g] gouraudtriangle\n",
643                     vertex[0].red, vertex[0].green, vertex[0].blue,
644                     vertex[i+1].red, vertex[i+1].green, vertex[i+1].blue,
645                     vertex[i+2].red, vertex[i+2].green, vertex[i+2].blue);
646          } else {
647
648            intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
649            intensity/=3.0;
650            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
651
652            intensity = vertex[1].red + vertex[1].green + vertex[1].blue;
653            intensity/=3.0;
654            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
655
656            intensity = vertex[2].red + vertex[2].green + vertex[2].blue;
657            intensity/=3.0;
658            fprintf (file, " [%g %g %g] gouraudtriangle\n", intensity, intensity, intensity);
659          }
660        }
661      } else {
662        fprintf (file, "newpath\n");
663        fprintf (file, "%g %g %g setrgbcolor\n", red, green, blue);
664        fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
665        for (i=1; i<nvertices; i++) {
666          fprintf (file, "%g %g lineto\n", vertex[i].x, vertex[i].y);
667        }
668        fprintf (file, "closepath fill\n\n");
669      }
670    }
671    loc += nvertices*7;
672    break;
673
674  case GL_POINT_TOKEN:
675    vertex=(Feedback3Dcolor*)loc;
676    if (print_colour) {
677      fprintf (file, "%g %g %g setrgbcolor\n", vertex[0].red, vertex[0].green, vertex[0].blue);
678    } else {
679      intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
680      intensity/=3.0;
681      fprintf (file, "%g %g %g setrgbcolor\n", intensity, intensity, intensity);
682    }     
683    fprintf(file, "%g %g %g 0 360 arc fill\n\n", vertex[0].x, vertex[0].y, pointSize / 2.0);
684    loc += 7;           /* Each vertex element in the feedback
685                           buffer is 7 GLfloats. */
686    break;
687  default:
688    /* XXX Left as an excersie to the reader. */
689    static G4bool spewPrimitiveEPSWarned = false;
690    if (!spewPrimitiveEPSWarned) {
691      std::ostringstream oss;
692      oss <<
693        "Incomplete implementation.  Unexpected token (" << token << ")."
694        "\n  (Seems to be caused by text.)";
695      G4Exception("G4OpenGLViewer::spewPrimitiveEPS",
696                  "Unexpected token",
697                  JustWarning,
698                  oss.str().c_str());
699      spewPrimitiveEPSWarned = true;
700    }
701  }
702  return loc;
703}
704
705typedef struct G4OpenGLViewerDepthIndex {
706  GLfloat *ptr;
707  GLfloat depth;
708} DepthIndex;
709
710extern "C" {
711  int G4OpenGLViewercompare(const void *a, const void *b)
712  {
713    const DepthIndex *p1 = (DepthIndex *) a;
714    const DepthIndex *p2 = (DepthIndex *) b;
715    GLfloat diff = p2->depth - p1->depth;
716   
717    if (diff > 0.0) {
718      return 1;
719    } else if (diff < 0.0) {
720      return -1;
721    } else {
722      return 0;
723    }
724  }
725}
726
727GLdouble G4OpenGLViewer::getSceneNearWidth()
728{
729  const G4Point3D targetPoint
730    = fSceneHandler.GetScene()->GetStandardTargetPoint()
731    + fVP.GetCurrentTargetPoint ();
732  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
733  if(radius<=0.) radius = 1.;
734  const G4double cameraDistance = fVP.GetCameraDistance (radius);
735  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
736  return 2 * fVP.GetFrontHalfHeight (pnear, radius);
737}
738
739GLdouble G4OpenGLViewer::getSceneFarWidth()
740{
741  const G4Point3D targetPoint
742    = fSceneHandler.GetScene()->GetStandardTargetPoint()
743    + fVP.GetCurrentTargetPoint ();
744  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
745  if(radius<=0.) radius = 1.;
746  const G4double cameraDistance = fVP.GetCameraDistance (radius);
747  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
748  const GLdouble pfar    = fVP.GetFarDistance  (cameraDistance, pnear, radius);
749  return 2 * fVP.GetFrontHalfHeight (pfar, radius);
750}
751
752
753GLdouble G4OpenGLViewer::getSceneDepth()
754{
755  const G4Point3D targetPoint
756    = fSceneHandler.GetScene()->GetStandardTargetPoint()
757    + fVP.GetCurrentTargetPoint ();
758  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
759  if(radius<=0.) radius = 1.;
760  const G4double cameraDistance = fVP.GetCameraDistance (radius);
761  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
762  return fVP.GetFarDistance  (cameraDistance, pnear, radius)- pnear;
763}
764
765
766void G4OpenGLViewer::spewSortedFeedback(FILE * file, GLint size, GLfloat * buffer)
767{
768  int token;
769  GLfloat *loc, *end;
770  Feedback3Dcolor *vertex;
771  GLfloat depthSum;
772  int nprimitives, item;
773  DepthIndex *prims;
774  int nvertices, i;
775
776  end = buffer + size;
777
778  /* Count how many primitives there are. */
779  nprimitives = 0;
780  loc = buffer;
781  while (loc < end) {
782    token = int (*loc);
783    loc++;
784    switch (token) {
785    case GL_LINE_TOKEN:
786    case GL_LINE_RESET_TOKEN:
787      loc += 14;
788      nprimitives++;
789      break;
790    case GL_POLYGON_TOKEN:
791      nvertices = int (*loc);
792      loc++;
793      loc += (7 * nvertices);
794      nprimitives++;
795      break;
796    case GL_POINT_TOKEN:
797      loc += 7;
798      nprimitives++;
799      break;
800    default:
801      /* XXX Left as an excersie to the reader. */
802      static G4bool spewSortedFeedbackWarned = false;
803      if (!spewSortedFeedbackWarned) {
804        std::ostringstream oss;
805        oss <<
806          "Incomplete implementation.  Unexpected token (" << token << ")."
807          "\n  (Seems to be caused by text.)";
808        G4Exception("G4OpenGLViewer::spewSortedFeedback",
809                    "Unexpected token",
810                    JustWarning,
811                    oss.str().c_str());
812        spewSortedFeedbackWarned = true;
813      }
814      nprimitives++;
815    }
816  }
817
818  /* Allocate an array of pointers that will point back at
819     primitives in the feedback buffer.  There will be one
820     entry per primitive.  This array is also where we keep the
821     primitive's average depth.  There is one entry per
822     primitive  in the feedback buffer. */
823  prims = (DepthIndex *) malloc(sizeof(DepthIndex) * nprimitives);
824
825  item = 0;
826  loc = buffer;
827  while (loc < end) {
828    prims[item].ptr = loc;  /* Save this primitive's location. */
829    token = int (*loc);
830    loc++;
831    switch (token) {
832    case GL_LINE_TOKEN:
833    case GL_LINE_RESET_TOKEN:
834      vertex = (Feedback3Dcolor *) loc;
835      depthSum = vertex[0].z + vertex[1].z;
836      prims[item].depth = depthSum / 2.0;
837      loc += 14;
838      break;
839    case GL_POLYGON_TOKEN:
840      nvertices = int (*loc);
841      loc++;
842      vertex = (Feedback3Dcolor *) loc;
843      depthSum = vertex[0].z;
844      for (i = 1; i < nvertices; i++) {
845        depthSum += vertex[i].z;
846      }
847      prims[item].depth = depthSum / nvertices;
848      loc += (7 * nvertices);
849      break;
850    case GL_POINT_TOKEN:
851      vertex = (Feedback3Dcolor *) loc;
852      prims[item].depth = vertex[0].z;
853      loc += 7;
854      break;
855    default:
856      /* XXX Left as an excersie to the reader. */
857      assert(1);
858    }
859    item++;
860  }
861  assert(item == nprimitives);
862
863  /* Sort the primitives back to front. */
864  qsort(prims, nprimitives, sizeof(DepthIndex), G4OpenGLViewercompare);
865
866  /* Understand that sorting by a primitives average depth
867     doesn't allow us to disambiguate some cases like self
868     intersecting polygons.  Handling these cases would require
869     breaking up the primitives.  That's too involved for this
870     example.  Sorting by depth is good enough for lots of
871     applications. */
872
873  /* Emit the Encapsulated PostScript for the primitives in
874     back to front order. */
875  for (item = 0; item < nprimitives; item++) {
876    (void) spewPrimitiveEPS(file, prims[item].ptr);
877  }
878
879  free(prims);
880}
881
882#endif
Note: See TracBrowser for help on using the repository browser.