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

Last change on this file since 682 was 678, checked in by garnier, 18 years ago

modifs de test pour rubberband et mise en place du button mouse zoom

  • Property svn:mime-type set to text/cpp
File size: 25.4 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: $
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 printf("G4OpenGLViewer::Pick\n");
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
727void G4OpenGLViewer::spewSortedFeedback(FILE * file, GLint size, GLfloat * buffer)
728{
729 int token;
730 GLfloat *loc, *end;
731 Feedback3Dcolor *vertex;
732 GLfloat depthSum;
733 int nprimitives, item;
734 DepthIndex *prims;
735 int nvertices, i;
736
737 end = buffer + size;
738
739 /* Count how many primitives there are. */
740 nprimitives = 0;
741 loc = buffer;
742 while (loc < end) {
743 token = int (*loc);
744 loc++;
745 switch (token) {
746 case GL_LINE_TOKEN:
747 case GL_LINE_RESET_TOKEN:
748 loc += 14;
749 nprimitives++;
750 break;
751 case GL_POLYGON_TOKEN:
752 nvertices = int (*loc);
753 loc++;
754 loc += (7 * nvertices);
755 nprimitives++;
756 break;
757 case GL_POINT_TOKEN:
758 loc += 7;
759 nprimitives++;
760 break;
761 default:
762 /* XXX Left as an excersie to the reader. */
763 static G4bool spewSortedFeedbackWarned = false;
764 if (!spewSortedFeedbackWarned) {
765 std::ostringstream oss;
766 oss <<
767 "Incomplete implementation. Unexpected token (" << token << ")."
768 "\n (Seems to be caused by text.)";
769 G4Exception("G4OpenGLViewer::spewSortedFeedback",
770 "Unexpected token",
771 JustWarning,
772 oss.str().c_str());
773 spewSortedFeedbackWarned = true;
774 }
775 nprimitives++;
776 }
777 }
778
779 /* Allocate an array of pointers that will point back at
780 primitives in the feedback buffer. There will be one
781 entry per primitive. This array is also where we keep the
782 primitive's average depth. There is one entry per
783 primitive in the feedback buffer. */
784 prims = (DepthIndex *) malloc(sizeof(DepthIndex) * nprimitives);
785
786 item = 0;
787 loc = buffer;
788 while (loc < end) {
789 prims[item].ptr = loc; /* Save this primitive's location. */
790 token = int (*loc);
791 loc++;
792 switch (token) {
793 case GL_LINE_TOKEN:
794 case GL_LINE_RESET_TOKEN:
795 vertex = (Feedback3Dcolor *) loc;
796 depthSum = vertex[0].z + vertex[1].z;
797 prims[item].depth = depthSum / 2.0;
798 loc += 14;
799 break;
800 case GL_POLYGON_TOKEN:
801 nvertices = int (*loc);
802 loc++;
803 vertex = (Feedback3Dcolor *) loc;
804 depthSum = vertex[0].z;
805 for (i = 1; i < nvertices; i++) {
806 depthSum += vertex[i].z;
807 }
808 prims[item].depth = depthSum / nvertices;
809 loc += (7 * nvertices);
810 break;
811 case GL_POINT_TOKEN:
812 vertex = (Feedback3Dcolor *) loc;
813 prims[item].depth = vertex[0].z;
814 loc += 7;
815 break;
816 default:
817 /* XXX Left as an excersie to the reader. */
818 assert(1);
819 }
820 item++;
821 }
822 assert(item == nprimitives);
823
824 /* Sort the primitives back to front. */
825 qsort(prims, nprimitives, sizeof(DepthIndex), G4OpenGLViewercompare);
826
827 /* Understand that sorting by a primitives average depth
828 doesn't allow us to disambiguate some cases like self
829 intersecting polygons. Handling these cases would require
830 breaking up the primitives. That's too involved for this
831 example. Sorting by depth is good enough for lots of
832 applications. */
833
834 /* Emit the Encapsulated PostScript for the primitives in
835 back to front order. */
836 for (item = 0; item < nprimitives; item++) {
837 (void) spewPrimitiveEPS(file, prims[item].ptr);
838 }
839
840 free(prims);
841}
842
843#endif
Note: See TracBrowser for help on using the repository browser.