Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4OpenGLSceneHandler.cc
Go to the documentation of this file.
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$
28//
29//
30// Andrew Walkden 27th March 1996
31// OpenGL stored scene - creates OpenGL display lists.
32// OpenGL immediate scene - draws immediately to buffer
33// (saving space on server).
34
35#ifdef G4VIS_BUILD_OPENGL_DRIVER
36
37// Included here - problems with HP compiler if not before other includes?
38#include "G4NURBS.hh"
39
40// Here follows a special for Mesa, the OpenGL emulator. Does not affect
41// other OpenGL's, as far as I'm aware. John Allison 18/9/96.
42#define CENTERLINE_CLPP /* CenterLine C++ workaround: */
43// Also seems to be required for HP's CC and AIX xlC, at least.
44
46#include "G4OpenGLViewer.hh"
48#include "G4Point3D.hh"
49#include "G4Normal3D.hh"
50#include "G4Transform3D.hh"
51#include "G4Polyline.hh"
52#include "G4Polymarker.hh"
53#include "G4Text.hh"
54#include "G4Circle.hh"
55#include "G4Square.hh"
56#include "G4VMarker.hh"
57#include "G4Polyhedron.hh"
58#include "G4VisAttributes.hh"
60#include "G4VPhysicalVolume.hh"
61#include "G4LogicalVolume.hh"
62#include "G4VSolid.hh"
63#include "G4Scene.hh"
64#include "G4VisExtent.hh"
65#include "G4AttHolder.hh"
67
68G4OpenGLSceneHandler::G4OpenGLSceneHandler (G4VGraphicsSystem& system,
69 G4int id,
70 const G4String& name):
71G4VSceneHandler (system, id, name),
72fPickName(0),
73// glFlush take about 90% time. Dividing glFlush number by 100 will
74// change the first vis time from 100% to 10+90/100 = 10,9%.
75fEventsDrawInterval(1),
76fEventsWaitingToBeFlushed(0),
77fThreePassCapable(false),
78fSecondPassForTransparencyRequested(false),
79fSecondPassForTransparency(false),
80fThirdPassForNonHiddenMarkersRequested(false),
81fThirdPassForNonHiddenMarkers(false)
82{}
83
84G4OpenGLSceneHandler::~G4OpenGLSceneHandler ()
85{
86 ClearStore ();
87}
88
89const GLubyte G4OpenGLSceneHandler::fStippleMaskHashed [128] = {
90 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
91 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
92 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
93 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
94 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
95 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
96 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
97 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
98 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
99 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
100 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
101 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
102 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
103 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
104 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55,
105 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55
106};
107
108void G4OpenGLSceneHandler::ClearAndDestroyAtts()
109{
110 std::map<GLuint, G4AttHolder*>::iterator i;
111 for (i = fPickMap.begin(); i != fPickMap.end(); ++i) delete i->second;
112 fPickMap.clear();
113}
114
115void G4OpenGLSceneHandler::ScaledFlush()
116{
117 fEventsWaitingToBeFlushed++;
118 if (fEventsWaitingToBeFlushed < fEventsDrawInterval) return;
119 glFlush();
120 fEventsWaitingToBeFlushed = 0;
121}
122
123void G4OpenGLSceneHandler::ProcessScene()
124{
125 fThreePassCapable = true;
126
128
129 // Repeat if required...
130 if (fSecondPassForTransparencyRequested) {
131 fSecondPassForTransparency = true;
133 fSecondPassForTransparency = false;
134 fSecondPassForTransparencyRequested = false;
135 }
136
137 // And again if required...
138 if (fThirdPassForNonHiddenMarkersRequested) {
139 fThirdPassForNonHiddenMarkers = true;
141 fThirdPassForNonHiddenMarkers = false;
142 fThirdPassForNonHiddenMarkersRequested = false;
143 }
144
145 fThreePassCapable = false;
146}
147
148void G4OpenGLSceneHandler::PreAddSolid
149(const G4Transform3D& objectTransformation,
150 const G4VisAttributes& visAttribs)
151{
152 G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
153}
154
155void G4OpenGLSceneHandler::BeginPrimitives
156(const G4Transform3D& objectTransformation)
157{
158 G4VSceneHandler::BeginPrimitives (objectTransformation);
159}
160
161void G4OpenGLSceneHandler::EndPrimitives ()
162{
164}
165
166void G4OpenGLSceneHandler::BeginPrimitives2D
167(const G4Transform3D& objectTransformation)
168{
169 G4VSceneHandler::BeginPrimitives2D (objectTransformation);
170}
171
172void G4OpenGLSceneHandler::EndPrimitives2D ()
173{
175}
176
177G4VSolid* G4OpenGLSceneHandler::CreateSectionSolid ()
178{
179 // Clipping done in G4OpenGLViewer::SetView.
180 // return 0;
181
182 // But...OpenGL no longer seems to reconstruct clipped edges, so,
183 // when the BooleanProcessor is up to it, abandon this and use
184 // generic clipping in G4VSceneHandler::CreateSectionSolid...
186}
187
188G4VSolid* G4OpenGLSceneHandler::CreateCutawaySolid ()
189{
190 // Cutaway done in G4OpenGLViewer::SetView.
191 // return 0;
192
193 // But...if not, when the BooleanProcessor is up to it...
195}
196
197void G4OpenGLSceneHandler::AddPrimitive (const G4Polyline& line)
198{
199 G4int nPoints = line.size ();
200 if (nPoints <= 0) return;
201
202 // Note: colour and depth test treated in sub-class.
203
204 glDisable (GL_LIGHTING);
205
206 // Get vis attributes - pick up defaults if none.
207 const G4VisAttributes* pVA =
208 fpViewer -> GetApplicableVisAttributes (line.GetVisAttributes ());
209
210 G4double lineWidth = GetLineWidth(pVA);
211 // Need access to method in G4OpenGLViewer. static_cast doesn't
212 // work with a virtual base class, so use dynamic_cast. No need to
213 // test the outcome since viewer is guaranteed to be a
214 // G4OpenGLViewer, but test it anyway to keep Coverity happy.
215 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer);
216 if (pGLViewer) pGLViewer->ChangeLineWidth(lineWidth);
217
218 glBegin (GL_LINE_STRIP);
219 for (G4int iPoint = 0; iPoint < nPoints; iPoint++) {
220 G4double x, y, z;
221 x = line[iPoint].x();
222 y = line[iPoint].y();
223 z = line[iPoint].z();
224 glVertex3d (x, y, z);
225 }
226 glEnd ();
227}
228
229void G4OpenGLSceneHandler::AddPrimitive (const G4Polymarker& polymarker)
230{
231 if (polymarker.size() == 0) {
232 return;
233 }
234
235 // Note: colour and depth test treated in sub-class.
236
237 glDisable (GL_LIGHTING);
238
239 // Get vis attributes - pick up defaults if none.
240 const G4VisAttributes* pVA =
241 fpViewer -> GetApplicableVisAttributes (polymarker.GetVisAttributes ());
242
243 G4double lineWidth = GetLineWidth(pVA);
244 // Need access to method in G4OpenGLViewer. static_cast doesn't
245 // work with a virtual base class, so use dynamic_cast. No need to
246 // test the outcome since viewer is guaranteed to be a
247 // G4OpenGLViewer, but test it anyway to keep Coverity happy.
248 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer);
249 if (pGLViewer) pGLViewer->ChangeLineWidth(lineWidth);
250
251 G4VMarker::FillStyle style = polymarker.GetFillStyle();
252
253 // G4bool filled = false; Not actually used - comment out to prevent compiler warnings (JA).
254 static G4bool hashedWarned = false;
255
256 switch (style) {
257 case G4VMarker::noFill:
258 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
259 //filled = false;
260 break;
262 if (!hashedWarned) {
263 G4cout << "Hashed fill style in G4OpenGLSceneHandler."
264 << "\n Not implemented. Using G4VMarker::filled."
265 << G4endl;
266 hashedWarned = true;
267 }
268 // Maybe use
269 //glPolygonStipple (fStippleMaskHashed);
270 // Drop through to filled...
272 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL);
273 //filled = true;
274 break;
275 }
276
277 MarkerSizeType sizeType;
278 G4double size = GetMarkerSize(polymarker, sizeType);
279
280 // Draw...
281 if (sizeType == world) { // Size specified in world coordinates.
282
283 G4int nSides;
284 G4double startPhi;
285 switch (polymarker.GetMarkerType()) {
286 default:
288 size = 1.;
289 // Drop through to circles
291 nSides = GetNoOfSides(pVA);
292 startPhi = 0.;
293 break;
295 nSides = 4;
296 startPhi = -pi / 4.;
297 break;
298 }
299
300 const G4Vector3D& viewpointDirection =
301 fpViewer -> GetViewParameters().GetViewpointDirection();
302 const G4Vector3D& up = fpViewer->GetViewParameters().GetUpVector();
303 const G4double dPhi = twopi / nSides;
304 const G4double radius = size / 2.;
305 G4Vector3D start = radius * (up.cross(viewpointDirection)).unit();
306 G4double phi;
307 G4int i;
308 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
309 glBegin (GL_POLYGON);
310 for (i = 0, phi = startPhi; i < nSides; i++, phi += dPhi) {
311 G4Vector3D r = start; r.rotate(phi, viewpointDirection);
312 G4Vector3D p = polymarker[iPoint] + r;
313 glVertex3d (p.x(), p.y(), p.z());
314 }
315 glEnd ();
316 }
317
318 } else { // Size specified in screen (window) coordinates.
319
320 pGLViewer->ChangePointSize(size);
321
322 //Antialiasing only for circles
323 switch (polymarker.GetMarkerType()) {
324 default:
327 glEnable (GL_POINT_SMOOTH); break;
329 glDisable (GL_POINT_SMOOTH); break;
330 }
331
332 glBegin (GL_POINTS);
333 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
334 G4Point3D centre = polymarker[iPoint];
335 glVertex3d(centre.x(),centre.y(),centre.z());
336 }
337 glEnd();
338 }
339}
340
341void G4OpenGLSceneHandler::AddPrimitive (const G4Text& text) {
342 // Pass to specific viewer via virtual function DrawText.
343 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer);
344 if (pGLViewer) pGLViewer->DrawText(text);
345}
346
347void G4OpenGLSceneHandler::AddPrimitive (const G4Circle& circle) {
348 G4Polymarker oneCircle(circle);
349 oneCircle.push_back(circle.GetPosition());
350 oneCircle.SetMarkerType(G4Polymarker::circles);
351 // Call this AddPrimitive to avoid re-doing sub-class code.
352 G4OpenGLSceneHandler::AddPrimitive(oneCircle);
353}
354
355void G4OpenGLSceneHandler::AddPrimitive (const G4Square& square) {
356 G4Polymarker oneSquare(square);
357 oneSquare.push_back(square.GetPosition());
358 oneSquare.SetMarkerType(G4Polymarker::squares);
359 // Call this AddPrimitive to avoid re-doing sub-class code.
360 G4OpenGLSceneHandler::AddPrimitive(oneSquare);
361}
362
363void G4OpenGLSceneHandler::AddPrimitive (const G4Scale& scale)
364{
366}
367
368//Method for handling G4Polyhedron objects for drawing solids.
369void G4OpenGLSceneHandler::AddPrimitive (const G4Polyhedron& polyhedron) {
370
371 // Assume all facets are planar convex quadrilaterals.
372 // Draw each facet individually
373
374 if (polyhedron.GetNoFacets() == 0) return;
375
376 // Need access to data in G4OpenGLViewer. static_cast doesn't work
377 // with a virtual base class, so use dynamic_cast.
378 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer);
379 if (!pGLViewer) return;
380
381 // Get vis attributes - pick up defaults if none.
382 const G4VisAttributes* pVA =
383 fpViewer -> GetApplicableVisAttributes (polyhedron.GetVisAttributes ());
384
385 // Get view parameters that the user can force through the vis
386 // attributes, thereby over-riding the current view parameter.
387 G4ViewParameters::DrawingStyle drawing_style = GetDrawingStyle (pVA);
388
389 // Note that in stored mode, because this call gets embedded in a display
390 // list, it is the colour _at the time of_ creation of the display list, so
391 // even if the colour is changed, for example, by interaction with a Qt
392 // window, current_colour does not change.
393 GLfloat current_colour [4];
394 glGetFloatv (GL_CURRENT_COLOR, current_colour);
395
396 G4bool isTransparent = false;
397 if (current_colour[3] < 1.) { // This object is transparent
398 isTransparent = true;
399 }
400
401 // This is the colour used to paint surfaces in hlr mode.
402 GLfloat clear_colour[4];
403 glGetFloatv (GL_COLOR_CLEAR_VALUE, clear_colour);
404
405 G4double lineWidth = GetLineWidth(pVA);
406 pGLViewer->ChangeLineWidth(lineWidth);
407
408 G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
409
410 G4bool clipping = pGLViewer->fVP.IsSection() || pGLViewer->fVP.IsCutaway();
411
412 // Lighting disabled unless otherwise requested
413 glDisable (GL_LIGHTING);
414
415 switch (drawing_style) {
417 // Set up as for hidden line removal but paint polygon faces later...
419 glEnable (GL_STENCIL_TEST);
420 // The stencil buffer is cleared in G4OpenGLViewer::ClearView.
421 // The procedure below leaves it clear.
422 glStencilFunc (GL_ALWAYS, 0, 1);
423 glStencilOp (GL_INVERT, GL_INVERT, GL_INVERT);
424 glEnable (GL_DEPTH_TEST);
425 glDepthFunc (GL_LEQUAL);
426 if (isTransparent) {
427 // Transparent...
428 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE);
429 glEnable(GL_COLOR_MATERIAL);
430 glDisable (GL_CULL_FACE);
431 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
432 } else {
433 // Opaque...
434 if (clipping) {
435 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE);
436 glEnable(GL_COLOR_MATERIAL);
437 glDisable (GL_CULL_FACE);
438 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
439 } else {
440 glColorMaterial(GL_FRONT, GL_AMBIENT_AND_DIFFUSE);
441 glEnable(GL_COLOR_MATERIAL);
442 glEnable (GL_CULL_FACE);
443 glCullFace (GL_BACK);
444 glPolygonMode (GL_FRONT, GL_LINE);
445 }
446 }
447 break;
449 glEnable (GL_DEPTH_TEST);
450 glDepthFunc (GL_LEQUAL);
451 if (isTransparent) {
452 // Transparent...
453 glDepthMask (GL_FALSE); // Make depth buffer read-only.
454 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE);
455 glEnable(GL_COLOR_MATERIAL);
456 glDisable (GL_CULL_FACE);
457 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL);
458 } else {
459 // Opaque...
460 glDepthMask (GL_TRUE); // Make depth buffer writable (default).
461 if (clipping) {
462 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE);
463 glEnable(GL_COLOR_MATERIAL);
464 glDisable (GL_CULL_FACE);
465 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL);
466 } else {
467 glColorMaterial(GL_FRONT, GL_AMBIENT_AND_DIFFUSE);
468 glEnable(GL_COLOR_MATERIAL);
469 glEnable (GL_CULL_FACE);
470 glCullFace (GL_BACK);
471 glPolygonMode (GL_FRONT, GL_FILL);
472 }
473 }
474 if (!fProcessing2D) glEnable (GL_LIGHTING);
475 break;
477 default:
478 glEnable (GL_DEPTH_TEST);
479 glDepthFunc (GL_LEQUAL); //??? was GL_ALWAYS
480 glDisable (GL_CULL_FACE);
481 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
482 break;
483 }
484
485 //Loop through all the facets...
486 glBegin (GL_QUADS);
487 G4bool notLastFace;
488 do {
489
490 //First, find vertices, edgeflags and normals and note "not last facet"...
491 G4Point3D vertex[4];
492 G4int edgeFlag[4];
493 G4Normal3D normals[4];
494 G4int nEdges;
495 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals);
496
497 //Loop through the four edges of each G4Facet...
498 for(G4int edgeCount = 0; edgeCount < nEdges; ++edgeCount) {
499 // Check to see if edge is visible or not...
500 if (isAuxEdgeVisible) {
501 edgeFlag[edgeCount] = 1;
502 }
503 if (edgeFlag[edgeCount] > 0) {
504 glEdgeFlag (GL_TRUE);
505 } else {
506 glEdgeFlag (GL_FALSE);
507 }
508 glNormal3d (normals[edgeCount].x(),
509 normals[edgeCount].y(),
510 normals[edgeCount].z());
511 glVertex3d (vertex[edgeCount].x(),
512 vertex[edgeCount].y(),
513 vertex[edgeCount].z());
514 }
515 // HepPolyhedron produces triangles too; in that case add an extra
516 // vertex identical to first...
517 if (nEdges == 3) {
518 G4int edgeCount = 3;
519 normals[edgeCount] = normals[0];
520 vertex[edgeCount] = vertex[0];
521 edgeFlag[edgeCount] = -1;
522 glEdgeFlag (GL_FALSE);
523 glNormal3d (normals[edgeCount].x(),
524 normals[edgeCount].y(),
525 normals[edgeCount].z());
526 glVertex3d (vertex[edgeCount].x(),
527 vertex[edgeCount].y(),
528 vertex[edgeCount].z());
529 }
530 // Trap situation where number of edges is > 4...
531 if (nEdges > 4) {
532 G4cerr <<
533 "G4OpenGLSceneHandler::AddPrimitive(G4Polyhedron): WARNING"
534 "\n G4Polyhedron facet with " << nEdges << " edges" << G4endl;
535 }
536
537 glDisable(GL_COLOR_MATERIAL); // Revert to glMaterial for hlr/sr.
538
539 // Do it all over again (twice) for hlr...
540 if (drawing_style == G4ViewParameters::hlr ||
541 drawing_style == G4ViewParameters::hlhsr) {
542
543 glEnd (); // Placed here to balance glBegin above, allowing GL
544 // state changes below, then glBegin again. Avoids
545 // having glBegin/End pairs *inside* loop in the more
546 // usual case of no hidden line removal.
547
548 // Lighting disabled unless otherwise requested
549 glDisable (GL_LIGHTING);
550
551 // Draw through stencil...
552 glStencilFunc (GL_EQUAL, 0, 1);
553 glStencilOp (GL_KEEP, GL_KEEP, GL_KEEP);
554 if (drawing_style == G4ViewParameters::hlhsr) {
555 if (!fProcessing2D) glEnable (GL_LIGHTING);
556 }
557 glEnable (GL_DEPTH_TEST);
558 glDepthFunc (GL_LEQUAL);
559 if (isTransparent) {
560 // Transparent...
561 glDepthMask (GL_FALSE); // Make depth buffer read-only.
562 glDisable (GL_CULL_FACE);
563 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL);
564 } else {
565 // Opaque...
566 glDepthMask (GL_TRUE); // Make depth buffer writable (default).
567 if (clipping) {
568 glDisable (GL_CULL_FACE);
569 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL);
570 } else {
571 glEnable (GL_CULL_FACE);
572 glCullFace (GL_BACK);
573 glPolygonMode (GL_FRONT, GL_FILL);
574 }
575 }
576 GLfloat* painting_colour;
577 if (drawing_style == G4ViewParameters::hlr) {
578 if (isTransparent) {
579 // Transparent - don't paint...
580 goto end_of_drawing_through_stencil;
581 }
582 painting_colour = clear_colour;
583 } else { // drawing_style == G4ViewParameters::hlhsr
584 painting_colour = current_colour;
585 }
586 if (isTransparent) {
587 // Transparent...
588 glMaterialfv (GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE, painting_colour);
589 } else {
590 // Opaque...
591 glMaterialfv (GL_FRONT, GL_AMBIENT_AND_DIFFUSE, painting_colour);
592 }
593 glColor4fv (painting_colour);
594 glBegin (GL_QUADS);
595 for (int edgeCount = 0; edgeCount < 4; ++edgeCount) {
596 if (edgeFlag[edgeCount] > 0) {
597 glEdgeFlag (GL_TRUE);
598 } else {
599 glEdgeFlag (GL_FALSE);
600 }
601 glNormal3d (normals[edgeCount].x(),
602 normals[edgeCount].y(),
603 normals[edgeCount].z());
604 glVertex3d (vertex[edgeCount].x(),
605 vertex[edgeCount].y(),
606 vertex[edgeCount].z());
607 }
608 glEnd ();
609 end_of_drawing_through_stencil:
610
611 // and once more to reset the stencil bits...
612 glStencilFunc (GL_ALWAYS, 0, 1);
613 glStencilOp (GL_INVERT, GL_INVERT, GL_INVERT);
614 glDepthFunc (GL_LEQUAL); // to make sure line gets drawn.
615 if (isTransparent) {
616 // Transparent...
617 glDisable (GL_CULL_FACE);
618 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
619 } else {
620 // Opaque...
621 if (clipping) {
622 glDisable (GL_CULL_FACE);
623 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE);
624 } else {
625 glEnable (GL_CULL_FACE);
626 glCullFace (GL_BACK);
627 glPolygonMode (GL_FRONT, GL_LINE);
628 }
629 }
630 glDisable (GL_LIGHTING);
631 glColor4fv (current_colour);
632 glBegin (GL_QUADS);
633 for (int edgeCount = 0; edgeCount < 4; ++edgeCount) {
634 if (edgeFlag[edgeCount] > 0) {
635 glEdgeFlag (GL_TRUE);
636 } else {
637 glEdgeFlag (GL_FALSE);
638 }
639 glNormal3d (normals[edgeCount].x(),
640 normals[edgeCount].y(),
641 normals[edgeCount].z());
642 glVertex3d (vertex[edgeCount].x(),
643 vertex[edgeCount].y(),
644 vertex[edgeCount].z());
645 }
646 glEnd ();
647 glDepthFunc (GL_LEQUAL); // Revert for next facet.
648 glBegin (GL_QUADS); // Ready for next facet. GL
649 // says it ignores incomplete
650 // quadrilaterals, so final empty
651 // glBegin/End sequence should be OK.
652 }
653 } while (notLastFace);
654
655 glEnd ();
656 glDisable (GL_STENCIL_TEST); // Revert to default for next primitive.
657 glDepthMask (GL_TRUE); // Revert to default for next primitive.
658 glDisable (GL_LIGHTING); // Revert to default for next primitive.
659}
660
661//Method for handling G4NURBS objects for drawing solids.
662//Knots and Ctrl Pnts MUST be arrays of GLfloats.
663void G4OpenGLSceneHandler::AddPrimitive (const G4NURBS& nurb) {
664
665 GLUnurbsObj *gl_nurb;
666 gl_nurb = gluNewNurbsRenderer ();
667
668 GLfloat *u_knot_array, *u_knot_array_ptr;
669 u_knot_array = u_knot_array_ptr = new GLfloat [nurb.GetnbrKnots(G4NURBS::U)];
670 G4NURBS::KnotsIterator u_iterator (nurb, G4NURBS::U);
671 while (u_iterator.pick (u_knot_array_ptr++)){}
672
673 GLfloat *v_knot_array, *v_knot_array_ptr;
674 v_knot_array = v_knot_array_ptr = new GLfloat [nurb.GetnbrKnots(G4NURBS::V)];
675 G4NURBS::KnotsIterator v_iterator (nurb, G4NURBS::V);
676 while (v_iterator.pick (v_knot_array_ptr++)){}
677
678 GLfloat *ctrl_pnt_array, *ctrl_pnt_array_ptr;
679 ctrl_pnt_array = ctrl_pnt_array_ptr =
680 new GLfloat [nurb.GettotalnbrCtrlPts () * G4NURBS::NofC];
681 G4NURBS::CtrlPtsCoordsIterator c_p_iterator (nurb);
682 while (c_p_iterator.pick (ctrl_pnt_array_ptr++)){}
683
684 // Get vis attributes - pick up defaults if none.
685 const G4VisAttributes* pVA =
686 fpViewer -> GetApplicableVisAttributes (nurb.GetVisAttributes ());
687
688 // Get view parameters that the user can force through the vis
689 // attributes, thereby over-riding the current view parameter.
690 G4ViewParameters::DrawingStyle drawing_style = GetDrawingStyle (pVA);
691 //G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
692
693 //Get colour, etc..
694 const G4Colour& c = pVA -> GetColour ();
695
696 switch (drawing_style) {
697
699 // G4cout << "Hidden line removal not implememented in G4OpenGL.\n"
700 // << "Using hidden surface removal." << G4endl;
702 {
703 if (!fProcessing2D) glEnable (GL_LIGHTING);
704 glEnable (GL_DEPTH_TEST);
705 glEnable (GL_AUTO_NORMAL);
706 glEnable (GL_NORMALIZE);
707 gluNurbsProperty (gl_nurb, GLU_DISPLAY_MODE, GLU_FILL);
708 gluNurbsProperty (gl_nurb, GLU_SAMPLING_TOLERANCE, 50.0);
709 GLfloat materialColour [4];
710 materialColour [0] = c.GetRed ();
711 materialColour [1] = c.GetGreen ();
712 materialColour [2] = c.GetBlue ();
713 materialColour [3] = 1.0; // = c.GetAlpha () for transparency -
714 // but see complication in
715 // AddPrimitive(const G4Polyhedron&).
716 glMaterialfv (GL_FRONT, GL_AMBIENT_AND_DIFFUSE, materialColour);
717 break;
718 }
720 // G4cout << "Hidden line removal not implememented in G4OpenGL.\n"
721 // << "Using wireframe." << G4endl;
723 default:
724 glDisable (GL_LIGHTING);
725// glDisable (GL_DEPTH_TEST);
726 glEnable (GL_DEPTH_TEST);
727 glDisable (GL_AUTO_NORMAL);
728 glDisable (GL_NORMALIZE);
729 gluNurbsProperty (gl_nurb, GLU_DISPLAY_MODE, GLU_OUTLINE_POLYGON);
730 gluNurbsProperty (gl_nurb, GLU_SAMPLING_TOLERANCE, 50.0);
731 glColor4d(c.GetRed(), c.GetGreen(), c.GetBlue(),c.GetAlpha());
732 break;
733 }
734
735 gluBeginSurface (gl_nurb);
736 G4int u_stride = 4;
737 G4int v_stride = nurb.GetnbrCtrlPts(G4NURBS::U) * 4;
738
739 gluNurbsSurface (gl_nurb,
740 nurb.GetnbrKnots (G4NURBS::U), (GLfloat*)u_knot_array,
741 nurb.GetnbrKnots (G4NURBS::V), (GLfloat*)v_knot_array,
742 u_stride,
743 v_stride,
744 ctrl_pnt_array,
745 nurb.GetUorder (),
746 nurb.GetVorder (),
747 GL_MAP2_VERTEX_4);
748
749 gluEndSurface (gl_nurb);
750
751 delete [] u_knot_array; // These should be allocated with smart allocators
752 delete [] v_knot_array; // to avoid memory explosion.
753 delete [] ctrl_pnt_array;
754
755 gluDeleteNurbsRenderer (gl_nurb);
756}
757
758void G4OpenGLSceneHandler::AddCompound(const G4VTrajectory& traj) {
759 G4VSceneHandler::AddCompound(traj); // For now.
760}
761
762void G4OpenGLSceneHandler::AddCompound(const G4VHit& hit) {
763 G4VSceneHandler::AddCompound(hit); // For now.
764}
765
766void G4OpenGLSceneHandler::AddCompound(const G4VDigi& digi) {
767 G4VSceneHandler::AddCompound(digi); // For now.
768}
769
770void G4OpenGLSceneHandler::AddCompound(const G4THitsMap<G4double>& hits) {
771 G4VSceneHandler::AddCompound(hits); // For now.
772}
773
774#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
G4double GetBlue() const
Definition: G4Colour.hh:140
G4double GetAlpha() const
Definition: G4Colour.hh:141
G4double GetRed() const
Definition: G4Colour.hh:138
G4double GetGreen() const
Definition: G4Colour.hh:139
G4int GetnbrKnots(t_direction in_dir) const
Definition: G4NURBS.hh:459
@ NofC
Definition: G4NURBS.hh:106
G4int GetUorder() const
Definition: G4NURBS.hh:431
G4int GetVorder() const
Definition: G4NURBS.hh:432
G4int GettotalnbrCtrlPts() const
Definition: G4NURBS.hh:437
G4int GetnbrCtrlPts(t_direction in_dir) const
Definition: G4NURBS.hh:463
MarkerType GetMarkerType() const
Definition: G4Text.hh:73
Definition: G4VHit.hh:49
FillStyle GetFillStyle() const
G4Point3D GetPosition() const
virtual G4VSolid * CreateSectionSolid()
virtual void EndPrimitives()
virtual void ProcessScene()
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation)
virtual void EndPrimitives2D()
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation)
virtual G4VSolid * CreateCutawaySolid()
virtual void AddCompound(const G4VTrajectory &)
virtual void AddPrimitive(const G4Polyline &)=0
const G4VisAttributes * GetVisAttributes() const
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > & rotate(T a, const BasicVector3D< T > &v)
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=0, G4Normal3D *normals=0) const
const G4double pi