source: trunk/source/visualization/OpenGL/src/G4OpenGLViewer.cc @ 916

Last change on this file since 916 was 916, checked in by garnier, 15 years ago

renommage et suppression de EPS dans Qt

  • Property svn:mime-type set to text/cpp
File size: 33.3 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.47 2009/02/04 16:48:41 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),
84fPrintFilename ("G4OpenGL.eps"),
85fPrintColour (true),
86fVectoredPs (true),
87fOpenGLSceneHandler(scene),
88background (G4Colour(0.,0.,0.)),
89transparency_enabled (true),
90antialiasing_enabled (false),
91haloing_enabled (false),
92fStartTime(-DBL_MAX),
93fEndTime(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.),
110fPointSize (0)
111{
112  // Make changes to view parameters for OpenGL...
113  fVP.SetAutoRefresh(true);
114  fDefaultVP.SetAutoRefresh(true);
115  fWinSize_x = fVP.GetWindowSizeHintX();
116  fWinSize_y = fVP.GetWindowSizeHintY();
117
118  //  glClearColor (0.0, 0.0, 0.0, 0.0);
119  //  glClearDepth (1.0);
120  //  glDisable (GL_BLEND);
121  //  glDisable (GL_LINE_SMOOTH);
122  //  glDisable (GL_POLYGON_SMOOTH);
123
124}
125
126G4OpenGLViewer::~G4OpenGLViewer () {}
127
128void G4OpenGLViewer::InitializeGLView ()
129{
130  glClearColor (0.0, 0.0, 0.0, 0.0);
131  glClearDepth (1.0);
132  glDisable (GL_BLEND);
133  glDisable (GL_LINE_SMOOTH);
134  glDisable (GL_POLYGON_SMOOTH);
135
136
137void G4OpenGLViewer::ClearView () {
138  glClearColor (background.GetRed(),
139                background.GetGreen(),
140                background.GetBlue(),
141                1.);
142  glClearDepth (1.0);
143  //Below line does not compile with Mesa includes.
144  //glClear (GL_COLOR_BUFFER_BIT | GL_DEPTH_BUFFER_BIT);
145  glClear (GL_COLOR_BUFFER_BIT);
146  glClear (GL_DEPTH_BUFFER_BIT);
147  glClear (GL_STENCIL_BUFFER_BIT);
148  glFlush ();
149}
150
151
152/**
153 * Set the viewport of the scene
154 */
155void G4OpenGLViewer::ResizeGLView()
156{
157  int side = fWinSize_x;
158  if (fWinSize_y < fWinSize_x) side = fWinSize_y;
159  glViewport((fWinSize_x - side) / 2, (fWinSize_y - side) / 2, side, side); 
160}
161
162
163void G4OpenGLViewer::SetView () {
164
165  if (!fSceneHandler.GetScene()) {
166    G4cerr << "G4OpenGLStoredViewer: Creating a Viewer without a scene is not allowed. \nPlease use /vis/scene/create before /vis/open/.... "
167           << G4endl;
168    return;
169  }
170  // Calculates view representation based on extent of object being
171  // viewed and (initial) viewpoint.  (Note: it can change later due
172  // to user interaction via visualization system's GUI.)
173 
174  // Lighting.
175  GLfloat lightPosition [4];
176  lightPosition [0] = fVP.GetActualLightpointDirection().x();
177  lightPosition [1] = fVP.GetActualLightpointDirection().y();
178  lightPosition [2] = fVP.GetActualLightpointDirection().z();
179  lightPosition [3] = 0.;
180  // Light position is "true" light direction, so must come after gluLookAt.
181  GLfloat ambient [] = { 0.2, 0.2, 0.2, 1.};
182  GLfloat diffuse [] = { 0.8, 0.8, 0.8, 1.};
183  glEnable (GL_LIGHT0);
184  glLightfv (GL_LIGHT0, GL_AMBIENT, ambient);
185  glLightfv (GL_LIGHT0, GL_DIFFUSE, diffuse);
186 
187  // Get radius of scene, etc.
188  // Note that this procedure properly takes into account zoom, dolly and pan.
189  const G4Point3D targetPoint
190    = fSceneHandler.GetScene()->GetStandardTargetPoint()
191    + fVP.GetCurrentTargetPoint ();
192  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
193  if(radius<=0.) radius = 1.;
194  const G4double cameraDistance = fVP.GetCameraDistance (radius);
195  const G4Point3D cameraPosition =
196    targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
197  const GLdouble pnear  = fVP.GetNearDistance (cameraDistance, radius);
198  const GLdouble pfar   = fVP.GetFarDistance  (cameraDistance, pnear, radius);
199  const GLdouble right  = fVP.GetFrontHalfHeight (pnear, radius);
200  const GLdouble left   = -right;
201  const GLdouble bottom = left;
202  const GLdouble top    = right;
203 
204  // FIXME
205  ResizeGLView();
206  //SHOULD SetWindowsSizeHint()...
207
208  glMatrixMode (GL_PROJECTION); // set up Frustum.
209  glLoadIdentity();
210
211  const G4Vector3D scaleFactor = fVP.GetScaleFactor();
212  glScaled(scaleFactor.x(),scaleFactor.y(),scaleFactor.z());
213 
214  if (fVP.GetFieldHalfAngle() == 0.) {
215    glOrtho (left, right, bottom, top, pnear, pfar);
216  }
217  else {
218    glFrustum (left, right, bottom, top, pnear, pfar);
219  } 
220
221  glMatrixMode (GL_MODELVIEW); // apply further transformations to scene.
222  glLoadIdentity();
223 
224  const G4Normal3D& upVector = fVP.GetUpVector (); 
225  G4Point3D gltarget;
226  if (cameraDistance > 1.e-6 * radius) {
227    gltarget = targetPoint;
228  }
229  else {
230    gltarget = targetPoint - radius * fVP.GetViewpointDirection().unit();
231  }
232
233  const G4Point3D& pCamera = cameraPosition;  // An alias for brevity.
234  gluLookAt (pCamera.x(),  pCamera.y(),  pCamera.z(),       // Viewpoint.
235             gltarget.x(), gltarget.y(), gltarget.z(),      // Target point.
236             upVector.x(), upVector.y(), upVector.z());     // Up vector.
237
238  // Light position is "true" light direction, so must come after gluLookAt.
239  glLightfv (GL_LIGHT0, GL_POSITION, lightPosition);
240
241  // OpenGL no longer seems to reconstruct clipped edges, so, when the
242  // BooleanProcessor is up to it, abandon this and use generic
243  // clipping in G4OpenGLSceneHandler::CreateSectionPolyhedron.  Also,
244  // force kernel visit on change of clipping plane in
245  // G4OpenGLStoredViewer::CompareForKernelVisit.
246  if (fVP.IsSection () ) {  // pair of back to back clip planes.
247    const G4Plane3D& s = fVP.GetSectionPlane ();
248    double sArray[4];
249    sArray[0] = s.a();
250    sArray[1] = s.b();
251    sArray[2] = s.c();
252    sArray[3] = s.d() + radius * 1.e-05;
253    glClipPlane (GL_CLIP_PLANE0, sArray);
254    glEnable (GL_CLIP_PLANE0);
255    sArray[0] = -s.a();
256    sArray[1] = -s.b();
257    sArray[2] = -s.c();
258    sArray[3] = -s.d() + radius * 1.e-05;
259    glClipPlane (GL_CLIP_PLANE1, sArray);
260    glEnable (GL_CLIP_PLANE1);
261  } else {
262    glDisable (GL_CLIP_PLANE0);
263    glDisable (GL_CLIP_PLANE1);
264  }
265
266  const G4Planes& cutaways = fVP.GetCutawayPlanes();
267  size_t nPlanes = cutaways.size();
268  if (fVP.IsCutaway() &&
269      fVP.GetCutawayMode() == G4ViewParameters::cutawayIntersection &&
270      nPlanes > 0) {
271    double a[4];
272    a[0] = cutaways[0].a();
273    a[1] = cutaways[0].b();
274    a[2] = cutaways[0].c();
275    a[3] = cutaways[0].d();
276    glClipPlane (GL_CLIP_PLANE2, a);
277    glEnable (GL_CLIP_PLANE2);
278    if (nPlanes > 1) {
279      a[0] = cutaways[1].a();
280      a[1] = cutaways[1].b();
281      a[2] = cutaways[1].c();
282      a[3] = cutaways[1].d();
283      glClipPlane (GL_CLIP_PLANE3, a);
284      glEnable (GL_CLIP_PLANE3);
285    }
286    if (nPlanes > 2) {
287      a[0] = cutaways[2].a();
288      a[1] = cutaways[2].b();
289      a[2] = cutaways[2].c();
290      a[3] = cutaways[2].d();
291      glClipPlane (GL_CLIP_PLANE4, a);
292      glEnable (GL_CLIP_PLANE4);
293    }
294  } else {
295    glDisable (GL_CLIP_PLANE2);
296    glDisable (GL_CLIP_PLANE3);
297    glDisable (GL_CLIP_PLANE4);
298  }
299
300  // Background.
301  background = fVP.GetBackgroundColour ();
302
303}
304
305void G4OpenGLViewer::HaloingFirstPass () {
306 
307  //To perform haloing, first Draw all information to the depth buffer
308  //alone, using a chunky line width, and then Draw all info again, to
309  //the colour buffer, setting a thinner line width an the depth testing
310  //function to less than or equal, so if two lines cross, the one
311  //passing behind the other will not pass the depth test, and so not
312  //get rendered either side of the infront line for a short distance.
313
314  //First, disable writing to the colo(u)r buffer...
315  glColorMask (GL_FALSE, GL_FALSE, GL_FALSE, GL_FALSE);
316
317  //Now enable writing to the depth buffer...
318  glDepthMask (GL_TRUE);
319  glDepthFunc (GL_LESS);
320  glClearDepth (1.0);
321
322  //Finally, set the line width to something wide...
323  glLineWidth (3.0);
324
325}
326
327void G4OpenGLViewer::HaloingSecondPass () {
328
329  //And finally, turn the colour buffer back on with a sesible line width...
330  glColorMask (GL_TRUE, GL_TRUE, GL_TRUE, GL_TRUE);
331  glDepthFunc (GL_LEQUAL);
332  glLineWidth (1.0);
333
334}
335
336void G4OpenGLViewer::Pick(GLdouble x, GLdouble y)
337{
338  //G4cout << "X: " << x << ", Y: " << y << G4endl;
339  const G4int BUFSIZE = 512;
340  GLuint selectBuffer[BUFSIZE];
341  glSelectBuffer(BUFSIZE, selectBuffer);
342  glRenderMode(GL_SELECT);
343  glInitNames();
344  glPushName(0);
345  glMatrixMode(GL_PROJECTION);
346  G4double currentProjectionMatrix[16];
347  glGetDoublev(GL_PROJECTION_MATRIX, currentProjectionMatrix);
348  glPushMatrix();
349  glLoadIdentity();
350  GLint viewport[4];
351  glGetIntegerv(GL_VIEWPORT, viewport);
352  // Define 5x5 pixel pick area
353  gluPickMatrix(x, viewport[3] - y, 5., 5., viewport);
354  glMultMatrixd(currentProjectionMatrix);
355  glMatrixMode(GL_MODELVIEW);
356  DrawView();
357  GLint hits = glRenderMode(GL_RENDER);
358  if (hits < 0)
359    G4cout << "Too many hits.  Zoom in to reduce overlaps." << G4cout;
360  else if (hits > 0) {
361    //G4cout << hits << " hit(s)" << G4endl;
362    GLuint* p = selectBuffer;
363    for (GLint i = 0; i < hits; ++i) {
364      GLuint nnames = *p++;
365      *p++; //OR GLuint zmin = *p++;
366      *p++; //OR GLuint zmax = *p++;
367      //G4cout << "Hit " << i << ": " << nnames << " names"
368      //     << "\nzmin: " << zmin << ", zmax: " << zmax << G4endl;
369      for (GLuint j = 0; j < nnames; ++j) {
370        GLuint name = *p++;
371        //G4cout << "Name " << j << ": PickName: " << name << G4endl;
372        std::map<GLuint, G4AttHolder*>::iterator iter =
373          fOpenGLSceneHandler.fPickMap.find(name);
374        if (iter != fOpenGLSceneHandler.fPickMap.end()) {
375          G4AttHolder* attHolder = iter->second;
376          if(attHolder && attHolder->GetAttDefs().size()) {
377            for (size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
378              G4cout << G4AttCheck(attHolder->GetAttValues()[i],
379                                   attHolder->GetAttDefs()[i]);
380            }
381          }
382        }
383      }
384      G4cout << G4endl;
385    }
386  }
387  glMatrixMode(GL_PROJECTION);
388  glPopMatrix();
389  glMatrixMode(GL_MODELVIEW);
390}
391
392void G4OpenGLViewer::print() {
393
394  // Print vectored PostScript
395 
396  G4int size = 5000000;
397  GLfloat* feedback_buffer = new GLfloat[size];
398  glFeedbackBuffer (size, GL_3D_COLOR, feedback_buffer);
399  glRenderMode (GL_FEEDBACK);
400 
401  DrawView();
402
403  GLint returned;
404  returned = glRenderMode (GL_RENDER);
405 
406  FILE* file;
407  if (!fPrintFilename.empty()) {
408    file = fopen (fPrintFilename.c_str(), "w");
409    if (file) {
410      spewWireframeEPS (file, returned, feedback_buffer, "rendereps");
411    } else {
412      printf("Could not open %s\n", fPrintFilename.c_str());
413    }
414  } else {
415    printBuffer (returned, feedback_buffer);
416  }
417
418  delete[] feedback_buffer;
419}
420
421void G4OpenGLViewer::print3DcolorVertex(GLint size, GLint * count, GLfloat * buffer)
422{
423  G4int i;
424
425  printf("  ");
426  for (i = 0; i < 7; i++) {
427    printf("%4.2f ", buffer[size - (*count)]);
428    *count = *count - 1;
429  }
430  printf("\n");
431}
432
433void G4OpenGLViewer::spewWireframeEPS (FILE* file, GLint size, GLfloat* buffer, const char* cr) {
434
435  GLfloat EPS_GOURAUD_THRESHOLD=0.1;
436
437  GLfloat clearColor[4], viewport[4];
438  GLfloat lineWidth;
439  G4int i;
440
441  glGetFloatv (GL_VIEWPORT, viewport);
442  glGetFloatv (GL_COLOR_CLEAR_VALUE, clearColor);
443  glGetFloatv (GL_LINE_WIDTH, &lineWidth);
444  glGetFloatv (GL_POINT_SIZE, &fPointSize);
445
446  fputs ("%!PS-Adobe-2.0 EPSF-2.0\n", file);
447  fprintf (file, "%%%%Creator: %s (using OpenGL feedback)\n", cr);
448  fprintf (file, "%%%%BoundingBox: %g %g %g %g\n", viewport[0], viewport[1], viewport[2], viewport[3]);
449  fputs ("%%EndComments\n", file);
450  fputs ("\n", file);
451  fputs ("gsave\n", file);
452  fputs ("\n", file);
453
454  fputs ("% the gouraudtriangle PostScript fragment below is free\n", file);
455  fputs ("% written by Frederic Delhoume (delhoume@ilog.fr)\n", file);
456  fprintf (file, "/threshold %g def\n", EPS_GOURAUD_THRESHOLD);
457  for (i=0; gouraudtriangleEPS[i]; i++) {
458    fprintf (file, "%s\n", gouraudtriangleEPS[i]);
459  }
460
461  fprintf(file, "\n%g setlinewidth\n", lineWidth);
462 
463  fprintf (file, "%g %g %g setrgbcolor\n", clearColor[0], clearColor[1], clearColor[2]);
464  fprintf (file, "%g %g %g %g rectfill\n\n", viewport[0], viewport[1], viewport[2], viewport[3]);
465
466  spewSortedFeedback (file, size, buffer);
467
468  fputs ("grestore\n\n", file);
469  fputs ("showpage\n", file);
470
471  fclose(file);
472}
473
474void G4OpenGLViewer::printBuffer (GLint size, GLfloat* buffer) {
475
476  GLint count;
477  G4int token, nvertices;
478
479  count=size;
480  while(count) {
481    token=G4int (buffer[size-count]);
482    count--;
483    switch (token) {
484
485    case GL_PASS_THROUGH_TOKEN:
486      printf ("GL_PASS_THROUGH_TOKEN\n");
487      printf ("  %4.2f\n", buffer[size-count]);
488      count--;
489      break;
490
491    case GL_POINT_TOKEN:
492      printf ("GL_POINT_TOKEN\n");
493      print3DcolorVertex (size, &count, buffer);
494      break;
495
496    case GL_LINE_TOKEN:
497      printf ("GL_LINE_TOKEN\n");
498      print3DcolorVertex (size, &count, buffer);
499      print3DcolorVertex (size, &count, buffer);
500      break;
501     
502    case GL_LINE_RESET_TOKEN:
503      printf ("GL_LINE_RESET_TOKEN\n");
504      print3DcolorVertex (size, &count, buffer);
505      print3DcolorVertex (size, &count, buffer);
506      break;
507
508    case GL_POLYGON_TOKEN:
509      printf ("GL_POLYGON_TOKEN\n");
510      nvertices=G4int (buffer[size-count]);
511      count--;
512      for (; nvertices>0; nvertices--) {
513        print3DcolorVertex (size, &count, buffer);
514      }
515    }
516  }
517}
518
519G4float* G4OpenGLViewer::spewPrimitiveEPS (FILE* file, GLfloat* loc) {
520 
521  G4int token;
522  G4int nvertices, i;
523  GLfloat red, green, blue, intensity;
524  G4int smooth;
525  GLfloat dx, dy, dr, dg, db, absR, absG, absB, colormax;
526  G4int steps;
527  Feedback3Dcolor *vertex;
528  GLfloat xstep(0.), ystep(0.), rstep(0.), gstep(0.), bstep(0.);
529  GLfloat xnext(0.), ynext(0.), rnext(0.), gnext(0.), bnext(0.), distance(0.);
530
531  token=G4int (*loc);
532  loc++;
533  switch (token) {
534  case GL_LINE_RESET_TOKEN:
535  case GL_LINE_TOKEN:
536    vertex=(Feedback3Dcolor*)loc;
537    dr=vertex[1].red - vertex[0].red;
538    dg=vertex[1].green - vertex[0].green;
539    db=vertex[1].blue - vertex[0].blue;
540
541    if (!fPrintColour) {
542      dr+=(dg+db);
543      dr/=3.0;
544      dg=dr;
545      db=dr;
546    }
547
548    if (dr!=0 || dg!=0 || db!=0) {
549      dx=vertex[1].x - vertex[0].x;
550      dy=vertex[1].y - vertex[0].y;
551      distance=std::sqrt(dx*dx + dy*dy);
552
553      absR=std::fabs(dr);
554      absG=std::fabs(dg);
555      absB=std::fabs(db);
556
557      #define Max(a, b) (((a)>(b))?(a):(b))
558
559      #define EPS_SMOOTH_LINE_FACTOR 0.06
560
561      colormax=Max(absR, Max(absG, absB));
562      steps=Max(1, G4int (colormax*distance*EPS_SMOOTH_LINE_FACTOR));
563     
564      xstep=dx/steps;
565      ystep=dy/steps;
566
567      rstep=dr/steps;
568      gstep=dg/steps;
569      bstep=db/steps;
570
571      xnext=vertex[0].x;
572      ynext=vertex[0].y;
573      rnext=vertex[0].red;
574      gnext=vertex[0].green;
575      bnext=vertex[0].blue;
576
577      if (!fPrintColour) {
578        rnext+=(gnext+bnext);
579        rnext/=3.0;
580        gnext=rnext;
581        bnext=rnext;
582      }
583
584      xnext -= xstep/2.0;
585      ynext -= ystep/2.0;
586      rnext -= rstep/2.0;
587      gnext -= gstep/2.0;
588      bnext -= bstep/2.0;
589    } else {
590      steps=0;
591    }
592    if (fPrintColour) {
593      fprintf (file, "%g %g %g setrgbcolor\n",
594               vertex[0].red, vertex[0].green, vertex[0].blue);
595    } else {
596      intensity = (vertex[0].red + vertex[0].green + vertex[0].blue) / 3.0;
597      fprintf (file, "%g %g %g setrgbcolor\n",
598               intensity, intensity, intensity);
599    }     
600    fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
601
602    for (i=0; i<steps; i++) {
603
604      xnext += xstep;
605      ynext += ystep;
606      rnext += rstep;
607      gnext += gstep;
608      bnext += bstep;
609
610      fprintf (file, "%g %g lineto stroke\n", xnext, ynext);
611      fprintf (file, "%g %g %g setrgbcolor\n", rnext, gnext, bnext);
612      fprintf (file, "%g %g moveto\n", xnext, ynext);
613    }
614    fprintf (file, "%g %g lineto stroke\n", vertex[1].x, vertex[1].y);
615
616    loc += 14;
617    break;
618
619  case GL_POLYGON_TOKEN:
620    nvertices = G4int (*loc);
621    loc++;
622    vertex=(Feedback3Dcolor*)loc;
623    if (nvertices>0) {
624      red=vertex[0].red;
625      green=vertex[0].green;
626      blue=vertex[0].blue;
627      smooth=0;
628     
629      if (!fPrintColour) {
630        red+=(green+blue);
631        red/=3.0;
632        green=red;
633        blue=red;
634      }
635     
636      if (fPrintColour) {
637        for (i=1; i<nvertices; i++) {
638          if (red!=vertex[i].red || green!=vertex[i].green || blue!=vertex[i].blue) {
639            smooth=1;
640            break;
641          }
642        }
643      } else {
644        for (i=1; i<nvertices; i++) {
645          intensity = vertex[i].red + vertex[i].green + vertex[i].blue;
646          intensity/=3.0;
647          if (red!=intensity) {
648            smooth=1;
649            break;
650          }
651        }
652      }
653
654      if (smooth) {
655        G4int triOffset;
656        for (i=0; i<nvertices-2; i++) {
657          triOffset = i*7;
658          fprintf (file, "[%g %g %g %g %g %g]",
659                   vertex[0].x, vertex[i+1].x, vertex[i+2].x,
660                   vertex[0].y, vertex[i+1].y, vertex[i+2].y);
661          if (fPrintColour) {
662            fprintf (file, " [%g %g %g] [%g %g %g] [%g %g %g] gouraudtriangle\n",
663                     vertex[0].red, vertex[0].green, vertex[0].blue,
664                     vertex[i+1].red, vertex[i+1].green, vertex[i+1].blue,
665                     vertex[i+2].red, vertex[i+2].green, vertex[i+2].blue);
666          } else {
667
668            intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
669            intensity/=3.0;
670            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
671
672            intensity = vertex[1].red + vertex[1].green + vertex[1].blue;
673            intensity/=3.0;
674            fprintf (file, " [%g %g %g]", intensity, intensity, intensity);
675
676            intensity = vertex[2].red + vertex[2].green + vertex[2].blue;
677            intensity/=3.0;
678            fprintf (file, " [%g %g %g] gouraudtriangle\n", intensity, intensity, intensity);
679          }
680        }
681      } else {
682        fprintf (file, "newpath\n");
683        fprintf (file, "%g %g %g setrgbcolor\n", red, green, blue);
684        fprintf (file, "%g %g moveto\n", vertex[0].x, vertex[0].y);
685        for (i=1; i<nvertices; i++) {
686          fprintf (file, "%g %g lineto\n", vertex[i].x, vertex[i].y);
687        }
688        fprintf (file, "closepath fill\n\n");
689      }
690    }
691    loc += nvertices*7;
692    break;
693
694  case GL_POINT_TOKEN:
695    vertex=(Feedback3Dcolor*)loc;
696    if (fPrintColour) {
697      fprintf (file, "%g %g %g setrgbcolor\n", vertex[0].red, vertex[0].green, vertex[0].blue);
698    } else {
699      intensity = vertex[0].red + vertex[0].green + vertex[0].blue;
700      intensity/=3.0;
701      fprintf (file, "%g %g %g setrgbcolor\n", intensity, intensity, intensity);
702    }     
703    fprintf(file, "%g %g %g 0 360 arc fill\n\n", vertex[0].x, vertex[0].y, fPointSize / 2.0);
704    loc += 7;           /* Each vertex element in the feedback
705                           buffer is 7 GLfloats. */
706    break;
707  default:
708    /* XXX Left as an excersie to the reader. */
709    static G4bool spewPrimitiveEPSWarned = false;
710    if (!spewPrimitiveEPSWarned) {
711      std::ostringstream oss;
712      oss <<
713        "Incomplete implementation.  Unexpected token (" << token << ")."
714        "\n  (Seems to be caused by text.)";
715      G4Exception("G4OpenGLViewer::spewPrimitiveEPS",
716                  "Unexpected token",
717                  JustWarning,
718                  oss.str().c_str());
719      spewPrimitiveEPSWarned = true;
720    }
721  }
722  return loc;
723}
724
725typedef struct G4OpenGLViewerDepthIndex {
726  GLfloat *ptr;
727  GLfloat depth;
728} DepthIndex;
729
730extern "C" {
731  int G4OpenGLViewercompare(const void *a, const void *b)
732  {
733    const DepthIndex *p1 = (DepthIndex *) a;
734    const DepthIndex *p2 = (DepthIndex *) b;
735    GLfloat diff = p2->depth - p1->depth;
736   
737    if (diff > 0.0) {
738      return 1;
739    } else if (diff < 0.0) {
740      return -1;
741    } else {
742      return 0;
743    }
744  }
745}
746
747GLubyte* G4OpenGLViewer::grabPixels (int inColor, unsigned int width, unsigned int height) {
748 
749  GLubyte* buffer;
750  GLint swapbytes, lsbfirst, rowlength;
751  GLint skiprows, skippixels, alignment;
752  GLenum format;
753  int size;
754
755  if (inColor) {
756    format = GL_RGB;
757    size = width*height*3;
758  } else {
759    format = GL_LUMINANCE;
760    size = width*height*1;
761  }
762
763  buffer = new GLubyte[size];
764  if (buffer == NULL)
765    return NULL;
766
767  glGetIntegerv (GL_UNPACK_SWAP_BYTES, &swapbytes);
768  glGetIntegerv (GL_UNPACK_LSB_FIRST, &lsbfirst);
769  glGetIntegerv (GL_UNPACK_ROW_LENGTH, &rowlength);
770
771  glGetIntegerv (GL_UNPACK_SKIP_ROWS, &skiprows);
772  glGetIntegerv (GL_UNPACK_SKIP_PIXELS, &skippixels);
773  glGetIntegerv (GL_UNPACK_ALIGNMENT, &alignment);
774
775  glPixelStorei (GL_UNPACK_SWAP_BYTES, GL_FALSE);
776  glPixelStorei (GL_UNPACK_LSB_FIRST, GL_FALSE);
777  glPixelStorei (GL_UNPACK_ROW_LENGTH, 0);
778
779  glPixelStorei (GL_UNPACK_SKIP_ROWS, 0);
780  glPixelStorei (GL_UNPACK_SKIP_PIXELS, 0);
781  glPixelStorei (GL_UNPACK_ALIGNMENT, 1);
782
783  glReadPixels (0, 0, (GLsizei)width, (GLsizei)height, format, GL_UNSIGNED_BYTE, (GLvoid*) buffer);
784
785  glPixelStorei (GL_UNPACK_SWAP_BYTES, swapbytes);
786  glPixelStorei (GL_UNPACK_LSB_FIRST, lsbfirst);
787  glPixelStorei (GL_UNPACK_ROW_LENGTH, rowlength);
788 
789  glPixelStorei (GL_UNPACK_SKIP_ROWS, skiprows);
790  glPixelStorei (GL_UNPACK_SKIP_PIXELS, skippixels);
791  glPixelStorei (GL_UNPACK_ALIGNMENT, alignment);
792 
793  return buffer;
794}
795
796int G4OpenGLViewer::generateEPS (const char* filnam,
797                                int inColour,
798                                unsigned int width,
799                                unsigned int height) {
800
801  FILE* fp;
802  GLubyte* pixels;
803  GLubyte* curpix;
804  int components, pos, i;
805
806  pixels = grabPixels (inColour, width, height);
807
808  if (pixels == NULL)
809    return 1;
810  if (inColour) {
811    components = 3;
812  } else {
813    components = 1;
814  }
815 
816  fp = fopen (filnam, "w");
817  if (fp == NULL) {
818    return 2;
819  }
820 
821  fprintf (fp, "%%!PS-Adobe-2.0 EPSF-1.2\n");
822  fprintf (fp, "%%%%Title: %s\n", filnam);
823  fprintf (fp, "%%%%Creator: OpenGL pixmap render output\n");
824  fprintf (fp, "%%%%BoundingBox: 0 0 %d %d\n", width, height);
825  fprintf (fp, "%%%%EndComments\n");
826  fprintf (fp, "gsave\n");
827  fprintf (fp, "/bwproc {\n");
828  fprintf (fp, "    rgbproc\n");
829  fprintf (fp, "    dup length 3 idiv string 0 3 0 \n");
830  fprintf (fp, "    5 -1 roll {\n");
831  fprintf (fp, "    add 2 1 roll 1 sub dup 0 eq\n");
832  fprintf (fp, "    { pop 3 idiv 3 -1 roll dup 4 -1 roll dup\n");
833  fprintf (fp, "       3 1 roll 5 -1 roll } put 1 add 3 0 \n");
834  fprintf (fp, "    { 2 1 roll } ifelse\n");
835  fprintf (fp, "    }forall\n");
836  fprintf (fp, "    pop pop pop\n");
837  fprintf (fp, "} def\n");
838  fprintf (fp, "systemdict /colorimage known not {\n");
839  fprintf (fp, "   /colorimage {\n");
840  fprintf (fp, "       pop\n");
841  fprintf (fp, "       pop\n");
842  fprintf (fp, "       /rgbproc exch def\n");
843  fprintf (fp, "       { bwproc } image\n");
844  fprintf (fp, "   }  def\n");
845  fprintf (fp, "} if\n");
846  fprintf (fp, "/picstr %d string def\n", width * components);
847  fprintf (fp, "%d %d scale\n", width, height);
848  fprintf (fp, "%d %d %d\n", width, height, 8);
849  fprintf (fp, "[%d 0 0 %d 0 0]\n", width, height);
850  fprintf (fp, "{currentfile picstr readhexstring pop}\n");
851  fprintf (fp, "false %d\n", components);
852  fprintf (fp, "colorimage\n");
853 
854  curpix = (GLubyte*) pixels;
855  pos = 0;
856  for (i = width*height*components; i>0; i--) {
857    fprintf (fp, "%02hx ", *(curpix++));
858    if (++pos >= 32) {
859      fprintf (fp, "\n");
860      pos = 0;
861    }
862  }
863  if (pos)
864    fprintf (fp, "\n");
865
866  fprintf (fp, "grestore\n");
867  fprintf (fp, "showpage\n");
868  delete pixels;
869  fclose (fp);
870  return 0;
871}
872
873GLdouble G4OpenGLViewer::getSceneNearWidth()
874{
875  const G4Point3D targetPoint
876    = fSceneHandler.GetScene()->GetStandardTargetPoint()
877    + fVP.GetCurrentTargetPoint ();
878  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
879  if(radius<=0.) radius = 1.;
880  const G4double cameraDistance = fVP.GetCameraDistance (radius);
881  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
882  return 2 * fVP.GetFrontHalfHeight (pnear, radius);
883}
884
885GLdouble G4OpenGLViewer::getSceneFarWidth()
886{
887  const G4Point3D targetPoint
888    = fSceneHandler.GetScene()->GetStandardTargetPoint()
889    + fVP.GetCurrentTargetPoint ();
890  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
891  if(radius<=0.) radius = 1.;
892  const G4double cameraDistance = fVP.GetCameraDistance (radius);
893  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
894  const GLdouble pfar    = fVP.GetFarDistance  (cameraDistance, pnear, radius);
895  return 2 * fVP.GetFrontHalfHeight (pfar, radius);
896}
897
898
899GLdouble G4OpenGLViewer::getSceneDepth()
900{
901  const G4Point3D targetPoint
902    = fSceneHandler.GetScene()->GetStandardTargetPoint()
903    + fVP.GetCurrentTargetPoint ();
904  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
905  if(radius<=0.) radius = 1.;
906  const G4double cameraDistance = fVP.GetCameraDistance (radius);
907  const GLdouble pnear   = fVP.GetNearDistance (cameraDistance, radius);
908  return fVP.GetFarDistance  (cameraDistance, pnear, radius)- pnear;
909}
910
911
912void G4OpenGLViewer::spewSortedFeedback(FILE * file, GLint size, GLfloat * buffer)
913{
914  int token;
915  GLfloat *loc, *end;
916  Feedback3Dcolor *vertex;
917  GLfloat depthSum;
918  int nprimitives, item;
919  DepthIndex *prims;
920  int nvertices, i;
921
922  end = buffer + size;
923
924  /* Count how many primitives there are. */
925  nprimitives = 0;
926  loc = buffer;
927  while (loc < end) {
928    token = int (*loc);
929    loc++;
930    switch (token) {
931    case GL_LINE_TOKEN:
932    case GL_LINE_RESET_TOKEN:
933      loc += 14;
934      nprimitives++;
935      break;
936    case GL_POLYGON_TOKEN:
937      nvertices = int (*loc);
938      loc++;
939      loc += (7 * nvertices);
940      nprimitives++;
941      break;
942    case GL_POINT_TOKEN:
943      loc += 7;
944      nprimitives++;
945      break;
946    default:
947      /* XXX Left as an excersie to the reader. */
948      static G4bool spewSortedFeedbackWarned = false;
949      if (!spewSortedFeedbackWarned) {
950        std::ostringstream oss;
951        oss <<
952          "Incomplete implementation.  Unexpected token (" << token << ")."
953          "\n  (Seems to be caused by text.)";
954        G4Exception("G4OpenGLViewer::spewSortedFeedback",
955                    "Unexpected token",
956                    JustWarning,
957                    oss.str().c_str());
958        spewSortedFeedbackWarned = true;
959      }
960      nprimitives++;
961    }
962  }
963
964  /* Allocate an array of pointers that will point back at
965     primitives in the feedback buffer.  There will be one
966     entry per primitive.  This array is also where we keep the
967     primitive's average depth.  There is one entry per
968     primitive  in the feedback buffer. */
969  prims = (DepthIndex *) malloc(sizeof(DepthIndex) * nprimitives);
970
971  item = 0;
972  loc = buffer;
973  while (loc < end) {
974    prims[item].ptr = loc;  /* Save this primitive's location. */
975    token = int (*loc);
976    loc++;
977    switch (token) {
978    case GL_LINE_TOKEN:
979    case GL_LINE_RESET_TOKEN:
980      vertex = (Feedback3Dcolor *) loc;
981      depthSum = vertex[0].z + vertex[1].z;
982      prims[item].depth = depthSum / 2.0;
983      loc += 14;
984      break;
985    case GL_POLYGON_TOKEN:
986      nvertices = int (*loc);
987      loc++;
988      vertex = (Feedback3Dcolor *) loc;
989      depthSum = vertex[0].z;
990      for (i = 1; i < nvertices; i++) {
991        depthSum += vertex[i].z;
992      }
993      prims[item].depth = depthSum / nvertices;
994      loc += (7 * nvertices);
995      break;
996    case GL_POINT_TOKEN:
997      vertex = (Feedback3Dcolor *) loc;
998      prims[item].depth = vertex[0].z;
999      loc += 7;
1000      break;
1001    default:
1002      /* XXX Left as an excersie to the reader. */
1003      assert(1);
1004    }
1005    item++;
1006  }
1007  assert(item == nprimitives);
1008
1009  /* Sort the primitives back to front. */
1010  qsort(prims, nprimitives, sizeof(DepthIndex), G4OpenGLViewercompare);
1011
1012  /* Understand that sorting by a primitives average depth
1013     doesn't allow us to disambiguate some cases like self
1014     intersecting polygons.  Handling these cases would require
1015     breaking up the primitives.  That's too involved for this
1016     example.  Sorting by depth is good enough for lots of
1017     applications. */
1018
1019  /* Emit the Encapsulated PostScript for the primitives in
1020     back to front order. */
1021  for (item = 0; item < nprimitives; item++) {
1022    (void) spewPrimitiveEPS(file, prims[item].ptr);
1023  }
1024
1025  free(prims);
1026}
1027
1028void G4OpenGLViewer::rotateScene(G4double dx, G4double dy,G4double deltaRotation)
1029{
1030
1031  G4Vector3D vp;
1032  G4Vector3D up;
1033 
1034  G4Vector3D xprime;
1035  G4Vector3D yprime;
1036  G4Vector3D zprime;
1037 
1038  G4double delta_alpha;
1039  G4double delta_theta;
1040 
1041  G4Vector3D new_vp;
1042  G4Vector3D new_up;
1043 
1044  G4double cosalpha;
1045  G4double sinalpha;
1046 
1047  G4Vector3D a1;
1048  G4Vector3D a2;
1049  G4Vector3D delta;
1050  G4Vector3D viewPoint;
1051
1052   
1053  //phi spin stuff here
1054 
1055  vp = fVP.GetViewpointDirection ().unit ();
1056  up = fVP.GetUpVector ().unit ();
1057 
1058  yprime = (up.cross(vp)).unit();
1059  zprime = (vp.cross(yprime)).unit();
1060 
1061  if (fVP.GetLightsMoveWithCamera()) {
1062    delta_alpha = dy * deltaRotation;
1063    delta_theta = -dx * deltaRotation;
1064  } else {
1065    delta_alpha = -dy * deltaRotation;
1066    delta_theta = dx * deltaRotation;
1067  }   
1068 
1069  delta_alpha *= deg;
1070  delta_theta *= deg;
1071 
1072  new_vp = std::cos(delta_alpha) * vp + std::sin(delta_alpha) * zprime;
1073 
1074  // to avoid z rotation flipping
1075  // to allow more than 360° rotation
1076
1077  const G4Point3D targetPoint
1078    = fSceneHandler.GetScene()->GetStandardTargetPoint()
1079    + fVP.GetCurrentTargetPoint ();
1080  G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
1081  if(radius<=0.) radius = 1.;
1082  const G4double cameraDistance = fVP.GetCameraDistance (radius);
1083  const G4Point3D cameraPosition =
1084    targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
1085
1086  if (fVP.GetLightsMoveWithCamera()) {
1087    new_up = (new_vp.cross(yprime)).unit();
1088    if (new_vp.z()*vp.z() <0) {
1089      new_up.set(new_up.x(),-new_up.y(),new_up.z());
1090    }
1091  } else {
1092    new_up = up;
1093    if (new_vp.z()*vp.z() <0) {
1094      new_up.set(new_up.x(),-new_up.y(),new_up.z());
1095    }
1096  }
1097  fVP.SetUpVector(new_up);
1098  ////////////////
1099  // Rotates by fixed azimuthal angle delta_theta.
1100 
1101  cosalpha = new_up.dot (new_vp.unit());
1102  sinalpha = std::sqrt (1. - std::pow (cosalpha, 2));
1103  yprime = (new_up.cross (new_vp.unit())).unit ();
1104  xprime = yprime.cross (new_up);
1105  // Projection of vp on plane perpendicular to up...
1106  a1 = sinalpha * xprime;
1107  // Required new projection...
1108  a2 = sinalpha * (std::cos (delta_theta) * xprime + std::sin (delta_theta) * yprime);
1109  // Required Increment vector...
1110  delta = a2 - a1;
1111  // So new viewpoint is...
1112  viewPoint = new_vp.unit() + delta;
1113 
1114  fVP.SetViewAndLights (viewPoint);
1115}
1116
1117#endif
Note: See TracBrowser for help on using the repository browser.