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

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

test print PS with gl2ps

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