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

Last change on this file since 866 was 858, checked in by garnier, 17 years ago

protection against no scene + Qapp could be now get on external app

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