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

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

corrections pour Qt3 mises sous cvs

  • Property svn:mime-type set to text/cpp
File size: 26.8 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.36 2008/03/10 16:57:04 lgarnier Exp $
28// GEANT4 tag $Name:  $
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  glMatrixMode (GL_MODELVIEW); // apply further transformations to scene.
199  glLoadIdentity();
200 
201  const G4Normal3D& upVector = fVP.GetUpVector (); 
202  G4Point3D gltarget;
203  if (cameraDistance > 1.e-6 * radius) {
204    gltarget = targetPoint;
205  }
206  else {
207    gltarget = targetPoint - radius * fVP.GetViewpointDirection().unit();
208  }
209
210  const G4Point3D& pCamera = cameraPosition;  // An alias for brevity.
211  gluLookAt (pCamera.x(),  pCamera.y(),  pCamera.z(),       // Viewpoint.
212             gltarget.x(), gltarget.y(), gltarget.z(),      // Target point.
213             upVector.x(), upVector.y(), upVector.z());     // Up vector.
214 
215  // Light position is "true" light direction, so must come after gluLookAt.
216  glLightfv (GL_LIGHT0, GL_POSITION, lightPosition);
217
218  // OpenGL no longer seems to reconstruct clipped edges, so, when the
219  // BooleanProcessor is up to it, abandon this and use generic
220  // clipping in G4OpenGLSceneHandler::CreateSectionPolyhedron.  Also,
221  // force kernel visit on change of clipping plane in
222  // G4OpenGLStoredViewer::CompareForKernelVisit.
223  if (fVP.IsSection () ) {  // pair of back to back clip planes.
224    const G4Plane3D& s = fVP.GetSectionPlane ();
225    double sArray[4];
226    sArray[0] = s.a();
227    sArray[1] = s.b();
228    sArray[2] = s.c();
229    sArray[3] = s.d() + radius * 1.e-05;
230    glClipPlane (GL_CLIP_PLANE0, sArray);
231    glEnable (GL_CLIP_PLANE0);
232    sArray[0] = -s.a();
233    sArray[1] = -s.b();
234    sArray[2] = -s.c();
235    sArray[3] = -s.d() + radius * 1.e-05;
236    glClipPlane (GL_CLIP_PLANE1, sArray);
237    glEnable (GL_CLIP_PLANE1);
238  } else {
239    glDisable (GL_CLIP_PLANE0);
240    glDisable (GL_CLIP_PLANE1);
241  }
242
243  const G4Planes& cutaways = fVP.GetCutawayPlanes();
244  size_t nPlanes = cutaways.size();
245  if (fVP.IsCutaway() &&
246      fVP.GetCutawayMode() == G4ViewParameters::cutawayIntersection &&
247      nPlanes > 0) {
248    double a[4];
249    a[0] = cutaways[0].a();
250    a[1] = cutaways[0].b();
251    a[2] = cutaways[0].c();
252    a[3] = cutaways[0].d();
253    glClipPlane (GL_CLIP_PLANE2, a);
254    glEnable (GL_CLIP_PLANE2);
255    if (nPlanes > 1) {
256      a[0] = cutaways[1].a();
257      a[1] = cutaways[1].b();
258      a[2] = cutaways[1].c();
259      a[3] = cutaways[1].d();
260      glClipPlane (GL_CLIP_PLANE3, a);
261      glEnable (GL_CLIP_PLANE3);
262    }
263    if (nPlanes > 2) {
264      a[0] = cutaways[2].a();
265      a[1] = cutaways[2].b();
266      a[2] = cutaways[2].c();
267      a[3] = cutaways[2].d();
268      glClipPlane (GL_CLIP_PLANE4, a);
269      glEnable (GL_CLIP_PLANE4);
270    }
271  } else {
272    glDisable (GL_CLIP_PLANE2);
273    glDisable (GL_CLIP_PLANE3);
274    glDisable (GL_CLIP_PLANE4);
275  }
276
277  // Background.
278  background = fVP.GetBackgroundColour ();
279
280}
281
282void G4OpenGLViewer::HaloingFirstPass () {
283 
284  //To perform haloing, first Draw all information to the depth buffer
285  //alone, using a chunky line width, and then Draw all info again, to
286  //the colour buffer, setting a thinner line width an the depth testing
287  //function to less than or equal, so if two lines cross, the one
288  //passing behind the other will not pass the depth test, and so not
289  //get rendered either side of the infront line for a short distance.
290
291  //First, disable writing to the colo(u)r buffer...
292  glColorMask (GL_FALSE, GL_FALSE, GL_FALSE, GL_FALSE);
293
294  //Now enable writing to the depth buffer...
295  glDepthMask (GL_TRUE);
296  glDepthFunc (GL_LESS);
297  glClearDepth (1.0);
298
299  //Finally, set the line width to something wide...
300  glLineWidth (3.0);
301
302}
303
304void G4OpenGLViewer::HaloingSecondPass () {
305
306  //And finally, turn the colour buffer back on with a sesible line width...
307  glColorMask (GL_TRUE, GL_TRUE, GL_TRUE, GL_TRUE);
308  glDepthFunc (GL_LEQUAL);
309  glLineWidth (1.0);
310
311}
312
313void G4OpenGLViewer::Pick(GLdouble x, GLdouble y)
314{
315  //G4cout << "X: " << x << ", Y: " << y << G4endl;
316  const G4int BUFSIZE = 512;
317  GLuint selectBuffer[BUFSIZE];
318  glSelectBuffer(BUFSIZE, selectBuffer);
319  glRenderMode(GL_SELECT);
320  glInitNames();
321  glPushName(0);
322  glMatrixMode(GL_PROJECTION);
323  G4double currentProjectionMatrix[16];
324  glGetDoublev(GL_PROJECTION_MATRIX, currentProjectionMatrix);
325  glPushMatrix();
326  glLoadIdentity();
327  GLint viewport[4];
328  glGetIntegerv(GL_VIEWPORT, viewport);
329  // Define 5x5 pixel pick area
330  gluPickMatrix(x, viewport[3] - y, 5., 5., viewport);
331  glMultMatrixd(currentProjectionMatrix);
332  glMatrixMode(GL_MODELVIEW);
333  DrawView();
334  GLint hits = glRenderMode(GL_RENDER);
335  if (hits < 0)
336    G4cout << "Too many hits.  Zoom in to reduce overlaps." << G4cout;
337  else if (hits > 0) {
338    //G4cout << hits << " hit(s)" << G4endl;
339    GLuint* p = selectBuffer;
340    for (GLint i = 0; i < hits; ++i) {
341      GLuint nnames = *p++;
342      *p++; //OR GLuint zmin = *p++;
343      *p++; //OR GLuint zmax = *p++;
344      //G4cout << "Hit " << i << ": " << nnames << " names"
345      //     << "\nzmin: " << zmin << ", zmax: " << zmax << G4endl;
346      for (GLuint j = 0; j < nnames; ++j) {
347        GLuint name = *p++;
348        //G4cout << "Name " << j << ": PickName: " << name << G4endl;
349        std::map<GLuint, G4AttHolder*>::iterator iter =
350          fOpenGLSceneHandler.fPickMap.find(name);
351        if (iter != fOpenGLSceneHandler.fPickMap.end()) {
352          G4AttHolder* attHolder = iter->second;
353          if(attHolder && attHolder->GetAttDefs().size()) {
354            for (size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
355              G4cout << G4AttCheck(attHolder->GetAttValues()[i],
356                                   attHolder->GetAttDefs()[i]);
357            }
358          }
359        }
360      }
361      G4cout << G4endl;
362    }
363  }
364  glMatrixMode(GL_PROJECTION);
365  glPopMatrix();
366  glMatrixMode(GL_MODELVIEW);
367}
368
369void G4OpenGLViewer::print() {
370
371  // Print vectored PostScript
372 
373  G4int size = 5000000;
374  GLfloat* feedback_buffer = new GLfloat[size];
375  glFeedbackBuffer (size, GL_3D_COLOR, feedback_buffer);
376  glRenderMode (GL_FEEDBACK);
377 
378  DrawView();
379
380  GLint returned;
381  returned = glRenderMode (GL_RENDER);
382 
383  FILE* file;
384  if (print_string) {
385    file = fopen (print_string, "w");
386    if (file) {
387      spewWireframeEPS (file, returned, feedback_buffer, "rendereps");
388    } else {
389      printf("Could not open %s\n", print_string);
390    }
391  } else {
392    printBuffer (returned, feedback_buffer);
393  }
394
395  delete[] feedback_buffer;
396}
397
398void G4OpenGLViewer::print3DcolorVertex(GLint size, GLint * count, GLfloat * buffer)
399{
400  G4int i;
401
402  printf("  ");
403  for (i = 0; i < 7; i++) {
404    printf("%4.2f ", buffer[size - (*count)]);
405    *count = *count - 1;
406  }
407  printf("\n");
408}
409
410void G4OpenGLViewer::spewWireframeEPS (FILE* file, GLint size, GLfloat* buffer, const char* cr) {
411
412  GLfloat EPS_GOURAUD_THRESHOLD=0.1;
413
414  GLfloat clearColor[4], viewport[4];
415  GLfloat lineWidth;
416  G4int i;
417
418  glGetFloatv (GL_VIEWPORT, viewport);
419  glGetFloatv (GL_COLOR_CLEAR_VALUE, clearColor);
420  glGetFloatv (GL_LINE_WIDTH, &lineWidth);
421  glGetFloatv (GL_POINT_SIZE, &pointSize);
422
423  fputs ("%!PS-Adobe-2.0 EPSF-2.0\n", file);
424  fprintf (file, "%%%%Creator: %s (using OpenGL feedback)\n", cr);
425  fprintf (file, "%%%%BoundingBox: %g %g %g %g\n", viewport[0], viewport[1], viewport[2], viewport[3]);
426  fputs ("%%EndComments\n", file);
427  fputs ("\n", file);
428  fputs ("gsave\n", file);
429  fputs ("\n", file);
430
431  fputs ("% the gouraudtriangle PostScript fragment below is free\n", file);
432  fputs ("% written by Frederic Delhoume (delhoume@ilog.fr)\n", file);
433  fprintf (file, "/threshold %g def\n", EPS_GOURAUD_THRESHOLD);
434  for (i=0; gouraudtriangleEPS[i]; i++) {
435    fprintf (file, "%s\n", gouraudtriangleEPS[i]);
436  }
437
438  fprintf(file, "\n%g setlinewidth\n", lineWidth);
439 
440  fprintf (file, "%g %g %g setrgbcolor\n", clearColor[0], clearColor[1], clearColor[2]);
441  fprintf (file, "%g %g %g %g rectfill\n\n", viewport[0], viewport[1], viewport[2], viewport[3]);
442
443  spewSortedFeedback (file, size, buffer);
444
445  fputs ("grestore\n\n", file);
446  fputs ("showpage\n", file);
447
448  fclose(file);
449}
450
451void G4OpenGLViewer::printBuffer (GLint size, GLfloat* buffer) {
452
453  GLint count;
454  G4int token, nvertices;
455
456  count=size;
457  while(count) {
458    token=G4int (buffer[size-count]);
459    count--;
460    switch (token) {
461
462    case GL_PASS_THROUGH_TOKEN:
463      printf ("GL_PASS_THROUGH_TOKEN\n");
464      printf ("  %4.2f\n", buffer[size-count]);
465      count--;
466      break;
467
468    case GL_POINT_TOKEN:
469      printf ("GL_POINT_TOKEN\n");
470      print3DcolorVertex (size, &count, buffer);
471      break;
472
473    case GL_LINE_TOKEN:
474      printf ("GL_LINE_TOKEN\n");
475      print3DcolorVertex (size, &count, buffer);
476      print3DcolorVertex (size, &count, buffer);
477      break;
478     
479    case GL_LINE_RESET_TOKEN:
480      printf ("GL_LINE_RESET_TOKEN\n");
481      print3DcolorVertex (size, &count, buffer);
482      print3DcolorVertex (size, &count, buffer);
483      break;
484
485    case GL_POLYGON_TOKEN:
486      printf ("GL_POLYGON_TOKEN\n");
487      nvertices=G4int (buffer[size-count]);
488      count--;
489      for (; nvertices>0; nvertices--) {
490        print3DcolorVertex (size, &count, buffer);
491      }
492    }
493  }
494}
495
496G4float* G4OpenGLViewer::spewPrimitiveEPS (FILE* file, GLfloat* loc) {
497 
498  G4int token;
499  G4int nvertices, i;
500  GLfloat red, green, blue, intensity;
501  G4int smooth;
502  GLfloat dx, dy, dr, dg, db, absR, absG, absB, colormax;
503  G4int steps;
504  Feedback3Dcolor *vertex;
505  GLfloat xstep(0.), ystep(0.), rstep(0.), gstep(0.), bstep(0.);
506  GLfloat xnext(0.), ynext(0.), rnext(0.), gnext(0.), bnext(0.), distance(0.);
507
508  token=G4int (*loc);
509  loc++;
510  switch (token) {
511  case GL_LINE_RESET_TOKEN:
512  case GL_LINE_TOKEN:
513    vertex=(Feedback3Dcolor*)loc;
514    dr=vertex[1].red - vertex[0].red;
515    dg=vertex[1].green - vertex[0].green;
516    db=vertex[1].blue - vertex[0].blue;
517
518    if (!print_colour) {
519      dr+=(dg+db);
520      dr/=3.0;
521      dg=dr;
522      db=dr;
523    }
524
525    if (dr!=0 || dg!=0 || db!=0) {
526      dx=vertex[1].x - vertex[0].x;
527      dy=vertex[1].y - vertex[0].y;
528      distance=std::sqrt(dx*dx + dy*dy);
529
530      absR=std::fabs(dr);
531      absG=std::fabs(dg);
532      absB=std::fabs(db);
533
534      #define Max(a, b) (((a)>(b))?(a):(b))
535
536      #define EPS_SMOOTH_LINE_FACTOR 0.06
537
538      colormax=Max(absR, Max(absG, absB));
539      steps=Max(1, G4int (colormax*distance*EPS_SMOOTH_LINE_FACTOR));
540     
541      xstep=dx/steps;
542      ystep=dy/steps;
543
544      rstep=dr/steps;
545      gstep=dg/steps;
546      bstep=db/steps;
547
548      xnext=vertex[0].x;
549      ynext=vertex[0].y;
550      rnext=vertex[0].red;
551      gnext=vertex[0].green;
552      bnext=vertex[0].blue;
553
554      if (!print_colour) {
555        rnext+=(gnext+bnext);
556        rnext/=3.0;
557        gnext=rnext;
558        bnext=rnext;
559      }
560
561      xnext -= xstep/2.0;
562      ynext -= ystep/2.0;
563      rnext -= rstep/2.0;
564      gnext -= gstep/2.0;
565      bnext -= bstep/2.0;
566    } else {
567      steps=0;
568    }
569    if (print_colour) {
570      fprintf (file, "%g %g %g setrgbcolor\n",
571               vertex[0].red, vertex[0].green, vertex[0].blue);
572    } else {
573      intensity = (vertex[0].red + vertex[0].green + vertex[0].blue) / 3.0;
574      fprintf (file, "%g %g %g setrgbcolor\n",
575               intensity, intensity, intensity);
576    }     
577    fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
578
579    for (i=0; i<steps; i++) {
580
581      xnext += xstep;
582      ynext += ystep;
583      rnext += rstep;
584      gnext += gstep;
585      bnext += bstep;
586
587      fprintf (file, "%g %g lineto stroke\n", xnext, ynext);
588      fprintf (file, "%g %g %g setrgbcolor\n", rnext, gnext, bnext);
589      fprintf (file, "%g %g moveto\n", xnext, ynext);
590    }
591    fprintf (file, "%g %g lineto stroke\n", vertex[1].x, vertex[1].y);
592
593    loc += 14;
594    break;
595
596  case GL_POLYGON_TOKEN:
597    nvertices = G4int (*loc);
598    loc++;
599    vertex=(Feedback3Dcolor*)loc;
600    if (nvertices>0) {
601      red=vertex[0].red;
602      green=vertex[0].green;
603      blue=vertex[0].blue;
604      smooth=0;
605     
606      if (!print_colour) {
607        red+=(green+blue);
608        red/=3.0;
609        green=red;
610        blue=red;
611      }
612     
613      if (print_colour) {
614        for (i=1; i<nvertices; i++) {
615          if (red!=vertex[i].red || green!=vertex[i].green || blue!=vertex[i].blue) {
616            smooth=1;
617            break;
618          }
619        }
620      } else {
621        for (i=1; i<nvertices; i++) {
622          intensity = vertex[i].red + vertex[i].green + vertex[i].blue;
623          intensity/=3.0;
624          if (red!=intensity) {
625            smooth=1;
626            break;
627          }
628        }
629      }
630
631      if (smooth) {
632        G4int triOffset;
633        for (i=0; i<nvertices-2; i++) {
634          triOffset = i*7;
635          fprintf (file, "[%g %g %g %g %g %g]",
636                   vertex[0].x, vertex[i+1].x, vertex[i+2].x,
637                   vertex[0].y, vertex[i+1].y, vertex[i+2].y);
638          if (print_colour) {
639            fprintf (file, " [%g %g %g] [%g %g %g] [%g %g %g] gouraudtriangle\n",
640                     vertex[0].red, vertex[0].green, vertex[0].blue,
641                     vertex[i+1].red, vertex[i+1].green, vertex[i+1].blue,
642                     vertex[i+2].red, vertex[i+2].green, vertex[i+2].blue);
643          } else {
644
645            intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
646            intensity/=3.0;
647            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
648
649            intensity = vertex[1].red + vertex[1].green + vertex[1].blue;
650            intensity/=3.0;
651            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
652
653            intensity = vertex[2].red + vertex[2].green + vertex[2].blue;
654            intensity/=3.0;
655            fprintf (file, " [%g %g %g] gouraudtriangle\n", intensity, intensity, intensity);
656          }
657        }
658      } else {
659        fprintf (file, "newpath\n");
660        fprintf (file, "%g %g %g setrgbcolor\n", red, green, blue);
661        fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
662        for (i=1; i<nvertices; i++) {
663          fprintf (file, "%g %g lineto\n", vertex[i].x, vertex[i].y);
664        }
665        fprintf (file, "closepath fill\n\n");
666      }
667    }
668    loc += nvertices*7;
669    break;
670
671  case GL_POINT_TOKEN:
672    vertex=(Feedback3Dcolor*)loc;
673    if (print_colour) {
674      fprintf (file, "%g %g %g setrgbcolor\n", vertex[0].red, vertex[0].green, vertex[0].blue);
675    } else {
676      intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
677      intensity/=3.0;
678      fprintf (file, "%g %g %g setrgbcolor\n", intensity, intensity, intensity);
679    }     
680    fprintf(file, "%g %g %g 0 360 arc fill\n\n", vertex[0].x, vertex[0].y, pointSize / 2.0);
681    loc += 7;           /* Each vertex element in the feedback
682                           buffer is 7 GLfloats. */
683    break;
684  default:
685    /* XXX Left as an excersie to the reader. */
686    static G4bool spewPrimitiveEPSWarned = false;
687    if (!spewPrimitiveEPSWarned) {
688      std::ostringstream oss;
689      oss <<
690        "Incomplete implementation.  Unexpected token (" << token << ")."
691        "\n  (Seems to be caused by text.)";
692      G4Exception("G4OpenGLViewer::spewPrimitiveEPS",
693                  "Unexpected token",
694                  JustWarning,
695                  oss.str().c_str());
696      spewPrimitiveEPSWarned = true;
697    }
698  }
699  return loc;
700}
701
702typedef struct G4OpenGLViewerDepthIndex {
703  GLfloat *ptr;
704  GLfloat depth;
705} DepthIndex;
706
707extern "C" {
708  int G4OpenGLViewercompare(const void *a, const void *b)
709  {
710    const DepthIndex *p1 = (DepthIndex *) a;
711    const DepthIndex *p2 = (DepthIndex *) b;
712    GLfloat diff = p2->depth - p1->depth;
713   
714    if (diff > 0.0) {
715      return 1;
716    } else if (diff < 0.0) {
717      return -1;
718    } else {
719      return 0;
720    }
721  }
722}
723
724GLdouble G4OpenGLViewer::getSceneNearWidth()
725{
726  const G4Point3D targetPoint
727    = fSceneHandler.GetScene()->GetStandardTargetPoint()
728    + fVP.GetCurrentTargetPoint ();
729  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
730  if(radius<=0.) radius = 1.;
731  const G4double cameraDistance = fVP.GetCameraDistance (radius);
732  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
733  return 2 * fVP.GetFrontHalfHeight (pnear, radius);
734}
735
736GLdouble G4OpenGLViewer::getSceneFarWidth()
737{
738  const G4Point3D targetPoint
739    = fSceneHandler.GetScene()->GetStandardTargetPoint()
740    + fVP.GetCurrentTargetPoint ();
741  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
742  if(radius<=0.) radius = 1.;
743  const G4double cameraDistance = fVP.GetCameraDistance (radius);
744  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
745  const GLdouble pfar    = fVP.GetFarDistance  (cameraDistance, pnear, radius);
746  return 2 * fVP.GetFrontHalfHeight (pfar, radius);
747}
748
749
750GLdouble G4OpenGLViewer::getSceneDepth()
751{
752  const G4Point3D targetPoint
753    = fSceneHandler.GetScene()->GetStandardTargetPoint()
754    + fVP.GetCurrentTargetPoint ();
755  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
756  if(radius<=0.) radius = 1.;
757  const G4double cameraDistance = fVP.GetCameraDistance (radius);
758  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
759  return fVP.GetFarDistance  (cameraDistance, pnear, radius)- pnear;
760}
761
762
763void G4OpenGLViewer::spewSortedFeedback(FILE * file, GLint size, GLfloat * buffer)
764{
765  int token;
766  GLfloat *loc, *end;
767  Feedback3Dcolor *vertex;
768  GLfloat depthSum;
769  int nprimitives, item;
770  DepthIndex *prims;
771  int nvertices, i;
772
773  end = buffer + size;
774
775  /* Count how many primitives there are. */
776  nprimitives = 0;
777  loc = buffer;
778  while (loc < end) {
779    token = int (*loc);
780    loc++;
781    switch (token) {
782    case GL_LINE_TOKEN:
783    case GL_LINE_RESET_TOKEN:
784      loc += 14;
785      nprimitives++;
786      break;
787    case GL_POLYGON_TOKEN:
788      nvertices = int (*loc);
789      loc++;
790      loc += (7 * nvertices);
791      nprimitives++;
792      break;
793    case GL_POINT_TOKEN:
794      loc += 7;
795      nprimitives++;
796      break;
797    default:
798      /* XXX Left as an excersie to the reader. */
799      static G4bool spewSortedFeedbackWarned = false;
800      if (!spewSortedFeedbackWarned) {
801        std::ostringstream oss;
802        oss <<
803          "Incomplete implementation.  Unexpected token (" << token << ")."
804          "\n  (Seems to be caused by text.)";
805        G4Exception("G4OpenGLViewer::spewSortedFeedback",
806                    "Unexpected token",
807                    JustWarning,
808                    oss.str().c_str());
809        spewSortedFeedbackWarned = true;
810      }
811      nprimitives++;
812    }
813  }
814
815  /* Allocate an array of pointers that will point back at
816     primitives in the feedback buffer.  There will be one
817     entry per primitive.  This array is also where we keep the
818     primitive's average depth.  There is one entry per
819     primitive  in the feedback buffer. */
820  prims = (DepthIndex *) malloc(sizeof(DepthIndex) * nprimitives);
821
822  item = 0;
823  loc = buffer;
824  while (loc < end) {
825    prims[item].ptr = loc;  /* Save this primitive's location. */
826    token = int (*loc);
827    loc++;
828    switch (token) {
829    case GL_LINE_TOKEN:
830    case GL_LINE_RESET_TOKEN:
831      vertex = (Feedback3Dcolor *) loc;
832      depthSum = vertex[0].z + vertex[1].z;
833      prims[item].depth = depthSum / 2.0;
834      loc += 14;
835      break;
836    case GL_POLYGON_TOKEN:
837      nvertices = int (*loc);
838      loc++;
839      vertex = (Feedback3Dcolor *) loc;
840      depthSum = vertex[0].z;
841      for (i = 1; i < nvertices; i++) {
842        depthSum += vertex[i].z;
843      }
844      prims[item].depth = depthSum / nvertices;
845      loc += (7 * nvertices);
846      break;
847    case GL_POINT_TOKEN:
848      vertex = (Feedback3Dcolor *) loc;
849      prims[item].depth = vertex[0].z;
850      loc += 7;
851      break;
852    default:
853      /* XXX Left as an excersie to the reader. */
854      assert(1);
855    }
856    item++;
857  }
858  assert(item == nprimitives);
859
860  /* Sort the primitives back to front. */
861  qsort(prims, nprimitives, sizeof(DepthIndex), G4OpenGLViewercompare);
862
863  /* Understand that sorting by a primitives average depth
864     doesn't allow us to disambiguate some cases like self
865     intersecting polygons.  Handling these cases would require
866     breaking up the primitives.  That's too involved for this
867     example.  Sorting by depth is good enough for lots of
868     applications. */
869
870  /* Emit the Encapsulated PostScript for the primitives in
871     back to front order. */
872  for (item = 0; item < nprimitives; item++) {
873    (void) spewPrimitiveEPS(file, prims[item].ptr);
874  }
875
876  free(prims);
877}
878
879#endif
Note: See TracBrowser for help on using the repository browser.