Geant4 Cross Reference |
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 // 28 // 29 // Andrew Walkden 27th March 1996 30 // OpenGL stored scene - creates OpenGL display lists. 31 // OpenGL immediate scene - draws immediately to buffer 32 // (saving space on server). 33 34 35 # include "G4OpenGLSceneHandler.hh" 36 # include "G4OpenGLViewer.hh" 37 # include "G4OpenGLTransform3D.hh" 38 # include "G4Point3D.hh" 39 # include "G4Normal3D.hh" 40 # include "G4Transform3D.hh" 41 # include "G4Polyline.hh" 42 # include "G4Polymarker.hh" 43 # include "G4Text.hh" 44 # include "G4Circle.hh" 45 # include "G4Square.hh" 46 # include "G4VMarker.hh" 47 # include "G4Polyhedron.hh" 48 # include "G4VisAttributes.hh" 49 # include "G4PhysicalVolumeModel.hh" 50 # include "G4VPhysicalVolume.hh" 51 # include "G4LogicalVolume.hh" 52 # include "G4VSolid.hh" 53 # include "G4Scene.hh" 54 # include "G4VisExtent.hh" 55 # include "G4AttHolder.hh" 56 # include "G4PhysicalConstants.hh" 57 # include "G4RunManager.hh" 58 # include "G4Run.hh" 59 # include "G4RunManagerFactory.hh" 60 # include "G4Mesh.hh" 61 # include "G4PseudoScene.hh" 62 # include "G4VisManager.hh" 63 64 const GLubyte G4OpenGLSceneHandler::fStippleMaskHashed [128] = { 65 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 66 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 67 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 68 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 69 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 70 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 71 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 72 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 73 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 74 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 75 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 76 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 77 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 78 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 79 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55, 80 0x55,0x55,0x55,0x55,0x55,0x55,0x55,0x55 81 }; 82 83 G4OpenGLSceneHandler::G4OpenGLSceneHandler (G4VGraphicsSystem& system, 84 G4int id, 85 const G4String& name): 86 G4VSceneHandler (system, id, name), 87 fPickName(0), 88 fThreePassCapable(false), 89 fSecondPassForTransparencyRequested(false), 90 fSecondPassForTransparency(false), 91 fThirdPassForNonHiddenMarkersRequested(false), 92 fThirdPassForNonHiddenMarkers(false), 93 fEdgeFlag(true) 94 { 95 } 96 97 G4OpenGLSceneHandler::~G4OpenGLSceneHandler () 98 { 99 ClearStore (); 100 } 101 102 void G4OpenGLSceneHandler::ClearAndDestroyAtts() 103 { 104 std::map<GLuint, G4AttHolder*>::iterator i; 105 for (i = fPickMap.begin(); i != fPickMap.end(); ++i) delete i->second; 106 fPickMap.clear(); 107 } 108 109 G4int G4OpenGLSceneHandler::fEntitiesFlushInterval = 100; 110 G4OpenGLSceneHandler::FlushAction 111 G4OpenGLSceneHandler::fFlushAction = G4OpenGLSceneHandler::NthEvent; 112 113 void G4OpenGLSceneHandler::ScaledFlush() 114 { 115 if (fReadyForTransients) { 116 117 // Drawing transients, e.g., trajectories. 118 119 if (!fpScene) { 120 // No scene - shouldn't happen 121 glFlush(); 122 return; 123 } 124 // Get event from modeling parameters 125 if (!fpModel) { 126 // No model - shouldn't happen 127 glFlush(); 128 return; 129 } 130 const G4ModelingParameters* modelingParameters = 131 fpModel->GetModelingParameters(); 132 if (!modelingParameters) { 133 // No modeling parameters - shouldn't happen 134 glFlush(); 135 return; 136 } 137 const G4Event* thisEvent = modelingParameters->GetEvent(); 138 if (!thisEvent) { 139 // No event, so not in event loop. 140 if (fFlushAction == endOfEvent) { 141 fFlushAction = endOfRun; 142 } else if (fFlushAction == NthEvent) { 143 fFlushAction = NthPrimitive; 144 } 145 } 146 G4RunManager* runMan = G4RunManagerFactory::GetMasterRunManager(); 147 if (!runMan) { 148 // No run manager - shouldn't happen 149 glFlush(); 150 return; 151 } 152 const G4Run* thisRun = runMan->GetCurrentRun(); 153 if (!thisRun) { 154 // No run, so not in event loop. 155 if (fFlushAction == endOfRun) { 156 fFlushAction = NthPrimitive; 157 } else if (fFlushAction == NthEvent) { 158 fFlushAction = NthPrimitive; 159 } 160 } 161 162 switch (fFlushAction) { 163 case endOfEvent: 164 // If "/vis/scene/endOfEventAction refresh", primitives are flushed at 165 // end of run anyway, so only scale if false. 166 if (!fpScene->GetRefreshAtEndOfEvent()) { 167 // But if "/vis/scene/endOfEventAction accumulate", ShowView is not 168 // called until end of run, so we have to watch for a new event. 169 // Get event from modeling parameters 170 G4int thisEventID = thisEvent->GetEventID(); 171 static G4int lastEventID = 0; 172 if (thisEventID != lastEventID) { 173 glFlush(); 174 lastEventID = thisEventID; 175 } 176 } 177 break; 178 case endOfRun: 179 // If "/vis/scene/endOfRunAction refresh", primitives are flushed at 180 // end of run anyway, so only scale if false. 181 if (!fpScene->GetRefreshAtEndOfRun()) { 182 // If "/vis/scene/endOfRunAction accumulate", ShowView is never called 183 // so we have to watch for a new run. 184 G4int thisRunID = thisRun->GetRunID(); 185 static G4int lastRunID = 0; 186 if (thisRunID != lastRunID) { 187 glFlush(); 188 lastRunID = thisRunID; 189 } 190 } 191 break; 192 case eachPrimitive: 193 // This is equivalent to numeric with fEntitiesFlushInterval == 1. 194 fEntitiesFlushInterval = 1; 195 [[fallthrough]]; // Fall through to NthPrimitive. 196 case NthPrimitive: 197 { // Encapsulate in scope {} brackets to satisfy Windows. 198 static G4int primitivesWaitingToBeFlushed = 0; 199 primitivesWaitingToBeFlushed++; 200 if (primitivesWaitingToBeFlushed < fEntitiesFlushInterval) return; 201 glFlush(); 202 primitivesWaitingToBeFlushed = 0; 203 break; 204 } 205 case NthEvent: 206 // If "/vis/scene/endOfEventAction refresh", primitives are flushed at 207 // end of event anyway, so only scale if false. 208 if (!fpScene->GetRefreshAtEndOfEvent()) { 209 G4int thisEventID = thisEvent->GetEventID(); 210 static G4int lastEventID = 0; 211 if (thisEventID != lastEventID) { 212 static G4int eventsWaitingToBeFlushed = 0; 213 eventsWaitingToBeFlushed++; 214 if (eventsWaitingToBeFlushed < fEntitiesFlushInterval) return; 215 glFlush(); 216 eventsWaitingToBeFlushed = 0; 217 lastEventID = thisEventID; 218 } 219 } 220 break; 221 case never: 222 break; 223 default: 224 break; 225 } 226 227 } 228 229 else 230 231 { 232 233 // For run duration model drawing (detector drawing): 234 // Immediate mode: a huge speed up is obtained if flushes are scaled. 235 // Stored mode: no discernable difference since drawing is done to the 236 // back buffer and then swapped. 237 // So eachPrimitive and NthPrimitive make sense. But endOfEvent and 238 // endOfRun are treated as "no action", i.e., a flush will only be issued, 239 // as happens anyway, when drawing is complete. 240 241 switch (fFlushAction) { 242 case endOfEvent: 243 break; 244 case endOfRun: 245 break; 246 case eachPrimitive: 247 // This is equivalent to NthPrimitive with fEntitiesFlushInterval == 1. 248 fEntitiesFlushInterval = 1; 249 [[fallthrough]]; // Fall through to NthPrimitive. 250 case NthPrimitive: 251 { // Encapsulate in scope {} brackets to satisfy Windows. 252 static G4int primitivesWaitingToBeFlushed = 0; 253 primitivesWaitingToBeFlushed++; 254 if (primitivesWaitingToBeFlushed < fEntitiesFlushInterval) return; 255 glFlush(); 256 primitivesWaitingToBeFlushed = 0; 257 break; 258 } 259 case NthEvent: 260 break; 261 case never: 262 break; 263 default: 264 break; 265 } 266 267 } 268 } 269 270 void G4OpenGLSceneHandler::ProcessScene() 271 { 272 fThreePassCapable = true; 273 274 G4VSceneHandler::ProcessScene(); 275 276 // Repeat if required... 277 if (fSecondPassForTransparencyRequested) { 278 fSecondPassForTransparency = true; 279 G4VSceneHandler::ProcessScene(); 280 fSecondPassForTransparency = false; 281 fSecondPassForTransparencyRequested = false; 282 } 283 284 // And again if required... 285 if (fThirdPassForNonHiddenMarkersRequested) { 286 fThirdPassForNonHiddenMarkers = true; 287 G4VSceneHandler::ProcessScene(); 288 fThirdPassForNonHiddenMarkers = false; 289 fThirdPassForNonHiddenMarkersRequested = false; 290 } 291 292 fThreePassCapable = false; 293 } 294 295 void G4OpenGLSceneHandler::PreAddSolid 296 (const G4Transform3D& objectTransformation, 297 const G4VisAttributes& visAttribs) 298 { 299 G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs); 300 } 301 302 void G4OpenGLSceneHandler::BeginPrimitives 303 (const G4Transform3D& objectTransformation) 304 { 305 G4VSceneHandler::BeginPrimitives (objectTransformation); 306 } 307 308 void G4OpenGLSceneHandler::EndPrimitives () 309 { 310 G4VSceneHandler::EndPrimitives (); 311 } 312 313 void G4OpenGLSceneHandler::BeginPrimitives2D 314 (const G4Transform3D& objectTransformation) 315 { 316 G4VSceneHandler::BeginPrimitives2D (objectTransformation); 317 } 318 319 void G4OpenGLSceneHandler::EndPrimitives2D () 320 { 321 G4VSceneHandler::EndPrimitives2D (); 322 } 323 324 G4DisplacedSolid* G4OpenGLSceneHandler::CreateSectionSolid () 325 { 326 return G4VSceneHandler::CreateSectionSolid(); 327 // If clipping done in G4OpenGLViewer::SetView 328 // return 0; 329 // Note: if you change this, you must also change 330 // G4OpenGLStoredViewer::CompareForKernelVisit 331 } 332 333 G4DisplacedSolid* G4OpenGLSceneHandler::CreateCutawaySolid () 334 { 335 // return G4VSceneHandler::CreateCutawaySolid(); 336 // If cutaway done in G4OpenGLViewer::SetView. 337 return 0; 338 // Note: if you change this, you must also change 339 // G4OpenGLStoredViewer::CompareForKernelVisit 340 } 341 342 void G4OpenGLSceneHandler::AddPrimitive (const G4Polyline& line) 343 { 344 std::size_t nPoints = line.size (); 345 if (nPoints <= 0) return; 346 347 // Note: colour and depth test treated in sub-class. 348 349 glDisable (GL_LIGHTING); 350 351 G4double lineWidth = GetLineWidth(fpVisAttribs); 352 // Need access to method in G4OpenGLViewer. static_cast doesn't 353 // work with a virtual base class, so use dynamic_cast. No need to 354 // test the outcome since viewer is guaranteed to be a 355 // G4OpenGLViewer, but test it anyway to keep Coverity happy. 356 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer); 357 if (pGLViewer) pGLViewer->ChangeLineWidth(lineWidth); 358 359 fEdgeFlag = true; 360 glBegin (GL_LINE_STRIP); 361 // No ned glEdgeFlag for lines : 362 // Boundary and nonboundary edge flags on vertices are significant only if GL_POLYGON_MODE is set to GL_POINT or GL_LINE. See glPolygonMode. 363 364 // glEdgeFlag (GL_TRUE); 365 for (std::size_t iPoint = 0; iPoint < nPoints; ++iPoint) { 366 G4double x, y, z; 367 x = line[iPoint].x(); 368 y = line[iPoint].y(); 369 z = line[iPoint].z(); 370 glVertex3d (x, y, z); 371 } 372 glEnd (); 373 } 374 375 void G4OpenGLSceneHandler::AddPrimitive (const G4Polymarker& polymarker) 376 { 377 if (polymarker.size() == 0) { 378 return; 379 } 380 381 // Note: colour and depth test treated in sub-class. 382 383 glDisable (GL_LIGHTING); 384 385 MarkerSizeType sizeType; 386 G4double size = GetMarkerSize(polymarker, sizeType); 387 388 // Need access to method in G4OpenGLViewer. static_cast doesn't 389 // work with a virtual base class, so use dynamic_cast. No need to 390 // test the outcome since viewer is guaranteed to be a 391 // G4OpenGLViewer, but test it anyway to keep Coverity happy. 392 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer); 393 if (!pGLViewer) return; 394 395 if (sizeType == world) { // Size specified in world coordinates. 396 G4double lineWidth = GetLineWidth(fpVisAttribs); 397 pGLViewer->ChangeLineWidth(lineWidth); 398 399 G4VMarker::FillStyle style = polymarker.GetFillStyle(); 400 401 // G4bool filled = false; Not actually used - comment out to prevent compiler warnings (JA). 402 static G4bool hashedWarned = false; 403 404 switch (style) { 405 case G4VMarker::noFill: 406 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 407 glEdgeFlag (GL_TRUE); 408 //filled = false; 409 break; 410 case G4VMarker::hashed: 411 if (!hashedWarned) { 412 G4cout << "Hashed fill style in G4OpenGLSceneHandler." 413 << "\n Not implemented. Using G4VMarker::filled." 414 << G4endl; 415 hashedWarned = true; 416 } 417 // Maybe use 418 //glPolygonStipple (fStippleMaskHashed); 419 [[fallthrough]]; // Drop through to filled... 420 case G4VMarker::filled: 421 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL); 422 //filled = true; 423 break; 424 } 425 } 426 427 // Draw... 428 if (sizeType == world) { // Size specified in world coordinates. 429 430 G4int nSides; 431 G4double startPhi; 432 switch (polymarker.GetMarkerType()) { 433 default: 434 case G4Polymarker::dots: 435 size = 1.; 436 [[fallthrough]]; // Fall through to circles 437 case G4Polymarker::circles: 438 nSides = GetNoOfSides(fpVisAttribs); 439 startPhi = 0.; 440 break; 441 case G4Polymarker::squares: 442 nSides = 4; 443 startPhi = -pi / 4.; 444 break; 445 } 446 447 const G4Vector3D& viewpointDirection = 448 fpViewer -> GetViewParameters().GetViewpointDirection(); 449 const G4Vector3D& up = fpViewer->GetViewParameters().GetUpVector(); 450 const G4double dPhi = twopi / nSides; 451 const G4double radius = size / 2.; 452 G4Vector3D start = radius * (up.cross(viewpointDirection)).unit(); 453 G4double phi; 454 G4int i; 455 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) { 456 fEdgeFlag = true; 457 glBegin (GL_POLYGON); 458 for (i = 0, phi = startPhi; i < nSides; i++, phi += dPhi) { 459 G4Vector3D r = start; r.rotate(phi, viewpointDirection); 460 G4Vector3D p = polymarker[iPoint] + r; 461 glVertex3d (p.x(), p.y(), p.z()); 462 } 463 glEnd (); 464 } 465 466 } else { // Size specified in screen (window) coordinates. 467 468 pGLViewer->ChangePointSize(size); 469 470 //Antialiasing only for circles 471 switch (polymarker.GetMarkerType()) { 472 default: 473 case G4Polymarker::dots: 474 case G4Polymarker::circles: 475 glEnable (GL_POINT_SMOOTH); break; 476 case G4Polymarker::squares: 477 glDisable (GL_POINT_SMOOTH); break; 478 } 479 glBegin (GL_POINTS); 480 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) { 481 G4Point3D centre = polymarker[iPoint]; 482 glVertex3d(centre.x(),centre.y(),centre.z()); 483 } 484 glEnd(); 485 } 486 } 487 488 void G4OpenGLSceneHandler::AddPrimitive (const G4Text& text) { 489 // Pass to specific viewer via virtual function DrawText. 490 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer); 491 if (pGLViewer) pGLViewer->DrawText(text); 492 } 493 494 void G4OpenGLSceneHandler::AddPrimitive (const G4Circle& circle) { 495 G4Polymarker oneCircle(circle); 496 oneCircle.push_back(circle.GetPosition()); 497 oneCircle.SetMarkerType(G4Polymarker::circles); 498 // Call this AddPrimitive to avoid re-doing sub-class code. 499 G4OpenGLSceneHandler::AddPrimitive(oneCircle); 500 } 501 502 void G4OpenGLSceneHandler::AddPrimitive (const G4Square& square) { 503 G4Polymarker oneSquare(square); 504 oneSquare.push_back(square.GetPosition()); 505 oneSquare.SetMarkerType(G4Polymarker::squares); 506 // Call this AddPrimitive to avoid re-doing sub-class code. 507 G4OpenGLSceneHandler::AddPrimitive(oneSquare); 508 } 509 510 //Method for handling G4Polyhedron objects for drawing solids. 511 void G4OpenGLSceneHandler::AddPrimitive (const G4Polyhedron& polyhedron) { 512 513 // Assume all facets are planar convex quadrilaterals. 514 // Draw each facet individually 515 516 if (polyhedron.GetNoFacets() == 0) return; 517 518 // Need access to data in G4OpenGLViewer. static_cast doesn't work 519 // with a virtual base class, so use dynamic_cast. 520 G4OpenGLViewer* pGLViewer = dynamic_cast<G4OpenGLViewer*>(fpViewer); 521 if (!pGLViewer) return; 522 523 // Get view parameters that the user can force through the vis 524 // attributes, thereby over-riding the current view parameter. 525 G4ViewParameters::DrawingStyle drawing_style = GetDrawingStyle (fpVisAttribs); 526 527 // Note that in stored mode, because this call gets embedded in a display 528 // list, it is the colour _at the time of_ creation of the display list, so 529 // even if the colour is changed, for example, by interaction with a Qt 530 // window, current_colour does not change. 531 GLfloat* painting_colour; 532 GLfloat clear_colour[4]; 533 GLfloat current_colour[4]; 534 glGetFloatv (GL_CURRENT_COLOR, current_colour); 535 536 G4bool isTransparent = false; 537 if (current_colour[3] < 1.) { // This object is transparent 538 isTransparent = true; 539 } 540 541 if (drawing_style == G4ViewParameters::hlr) { 542 // This is the colour used to paint surfaces in hlr mode. 543 glGetFloatv (GL_COLOR_CLEAR_VALUE, clear_colour); 544 painting_colour = clear_colour; 545 } else { // drawing_style == G4ViewParameters::hlhsr 546 painting_colour = current_colour; 547 } 548 549 G4double lineWidth = GetLineWidth(fpVisAttribs); 550 pGLViewer->ChangeLineWidth(lineWidth); 551 552 G4bool isAuxEdgeVisible = GetAuxEdgeVisible (fpVisAttribs); 553 554 G4bool clipping = pGLViewer->fVP.IsSection() || pGLViewer->fVP.IsCutaway(); 555 556 // Lighting disabled unless otherwise requested 557 glDisable (GL_LIGHTING); 558 559 switch (drawing_style) { 560 case (G4ViewParameters::hlhsr): 561 // Set up as for hidden line removal but paint polygon faces later... 562 case (G4ViewParameters::hlr): 563 glEnable (GL_STENCIL_TEST); 564 // The stencil buffer is cleared in G4OpenGLViewer::ClearView. 565 // The procedure below leaves it clear. 566 glStencilFunc (GL_ALWAYS, 0, 1); 567 glStencilOp (GL_INVERT, GL_INVERT, GL_INVERT); 568 glEnable (GL_DEPTH_TEST); 569 glDepthFunc (GL_LEQUAL); 570 if (isTransparent) { 571 // Transparent... 572 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE); 573 glEnable(GL_COLOR_MATERIAL); 574 //glDisable (GL_CULL_FACE); 575 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 576 } else { 577 // Opaque... 578 if (clipping) { 579 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE); 580 glEnable(GL_COLOR_MATERIAL); 581 //glDisable (GL_CULL_FACE); 582 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 583 } else { 584 glColorMaterial(GL_FRONT, GL_AMBIENT_AND_DIFFUSE); 585 glEnable(GL_COLOR_MATERIAL); 586 //glEnable (GL_CULL_FACE); 587 //glCullFace (GL_BACK); 588 glPolygonMode (GL_FRONT, GL_LINE); 589 } 590 } 591 break; 592 case (G4ViewParameters::hsr): 593 glEnable (GL_DEPTH_TEST); 594 glDepthFunc (GL_LEQUAL); 595 if (isTransparent) { 596 // Transparent... 597 glDepthMask (GL_FALSE); // Make depth buffer read-only. 598 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE); 599 glEnable(GL_COLOR_MATERIAL); 600 //glDisable (GL_CULL_FACE); 601 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL); 602 } else { 603 // Opaque... 604 glDepthMask (GL_TRUE); // Make depth buffer writable (default). 605 if (clipping) { 606 glColorMaterial(GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE); 607 glEnable(GL_COLOR_MATERIAL); 608 //glDisable (GL_CULL_FACE); 609 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL); 610 } else { 611 glColorMaterial(GL_FRONT, GL_AMBIENT_AND_DIFFUSE); 612 glEnable(GL_COLOR_MATERIAL); 613 //glEnable (GL_CULL_FACE); 614 //glCullFace (GL_BACK); 615 glPolygonMode (GL_FRONT, GL_FILL); 616 } 617 } 618 if (!fProcessing2D) glEnable (GL_LIGHTING); 619 break; 620 case (G4ViewParameters::wireframe): 621 default: 622 glEnable (GL_DEPTH_TEST); 623 glDepthFunc (GL_LEQUAL); //??? was GL_ALWAYS 624 //glDisable (GL_CULL_FACE); 625 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 626 break; 627 } 628 629 //Loop through all the facets... 630 fEdgeFlag = true; 631 glBegin (GL_QUADS); 632 glEdgeFlag (GL_TRUE); 633 G4bool notLastFace; 634 do { 635 636 //First, find vertices, edgeflags and normals and note "not last facet"... 637 G4Point3D vertex[4]; 638 G4int edgeFlag[4]; 639 G4Normal3D normals[4]; 640 G4int nEdges; 641 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals); 642 643 //Loop through the four edges of each G4Facet... 644 for(G4int edgeCount = 0; edgeCount < nEdges; ++edgeCount) { 645 // Check to see if edge is visible or not... 646 if (isAuxEdgeVisible) { 647 edgeFlag[edgeCount] = 1; 648 } 649 if (edgeFlag[edgeCount] > 0) { 650 if (fEdgeFlag != true) { 651 glEdgeFlag (GL_TRUE); 652 fEdgeFlag = true; 653 } 654 } else { 655 if (fEdgeFlag != false) { 656 glEdgeFlag (GL_FALSE); 657 fEdgeFlag = false; 658 } 659 } 660 glNormal3d (normals[edgeCount].x(), 661 normals[edgeCount].y(), 662 normals[edgeCount].z()); 663 glVertex3d (vertex[edgeCount].x(), 664 vertex[edgeCount].y(), 665 vertex[edgeCount].z()); 666 } 667 668 // HepPolyhedron produces triangles too; in that case add an extra 669 // vertex identical to first... 670 if (nEdges == 3) { 671 G4int edgeCount = 3; 672 normals[edgeCount] = normals[0]; 673 vertex[edgeCount] = vertex[0]; 674 edgeFlag[edgeCount] = -1; 675 if (fEdgeFlag != false) { 676 glEdgeFlag (GL_FALSE); 677 fEdgeFlag = false; 678 } 679 680 glNormal3d (normals[edgeCount].x(), 681 normals[edgeCount].y(), 682 normals[edgeCount].z()); 683 glVertex3d (vertex[edgeCount].x(), 684 vertex[edgeCount].y(), 685 vertex[edgeCount].z()); 686 } 687 // Trap situation where number of edges is > 4... 688 if (nEdges > 4) { 689 G4cerr << 690 "G4OpenGLSceneHandler::AddPrimitive(G4Polyhedron): WARNING" 691 "\n G4Polyhedron facet with " << nEdges << " edges" << G4endl; 692 } 693 694 695 // Do it all over again (twice) for hlr... 696 if (drawing_style == G4ViewParameters::hlr || 697 drawing_style == G4ViewParameters::hlhsr) { 698 699 glDisable(GL_COLOR_MATERIAL); // Revert to glMaterial for hlr/sr. 700 glEnd (); // Placed here to balance glBegin above, allowing GL 701 702 // state changes below, then glBegin again. Avoids 703 // having glBegin/End pairs *inside* loop in the more 704 // usual case of no hidden line removal. 705 706 // Lighting disabled unless otherwise requested 707 glDisable (GL_LIGHTING); 708 709 // Draw through stencil... 710 glStencilFunc (GL_EQUAL, 0, 1); 711 glStencilOp (GL_KEEP, GL_KEEP, GL_KEEP); 712 if (drawing_style == G4ViewParameters::hlhsr) { 713 if (!fProcessing2D) glEnable (GL_LIGHTING); 714 } 715 glEnable (GL_DEPTH_TEST); 716 glDepthFunc (GL_LEQUAL); 717 if (isTransparent) { 718 // Transparent... 719 glDepthMask (GL_FALSE); // Make depth buffer read-only. 720 //glDisable (GL_CULL_FACE); 721 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL); 722 } else { 723 // Opaque... 724 glDepthMask (GL_TRUE); // Make depth buffer writable (default). 725 if (clipping) { 726 //glDisable (GL_CULL_FACE); 727 glPolygonMode (GL_FRONT_AND_BACK, GL_FILL); 728 } else { 729 //glEnable (GL_CULL_FACE); 730 //glCullFace (GL_BACK); 731 glPolygonMode (GL_FRONT, GL_FILL); 732 } 733 } 734 if (drawing_style == G4ViewParameters::hlr) { 735 if (isTransparent) { 736 // Transparent - don't paint... 737 goto end_of_drawing_through_stencil; 738 } 739 } 740 if (isTransparent) { 741 // Transparent... 742 glMaterialfv (GL_FRONT_AND_BACK, GL_AMBIENT_AND_DIFFUSE, painting_colour); 743 } else { 744 // Opaque... 745 glMaterialfv (GL_FRONT, GL_AMBIENT_AND_DIFFUSE, painting_colour); 746 } 747 glColor4fv (painting_colour); 748 glBegin (GL_QUADS); 749 glEdgeFlag (GL_TRUE); 750 fEdgeFlag = true; 751 752 for (int edgeCount = 0; edgeCount < 4; ++edgeCount) { 753 if (edgeFlag[edgeCount] > 0) { 754 if (fEdgeFlag != true) { 755 glEdgeFlag (GL_TRUE); 756 fEdgeFlag = true; 757 } 758 } else { 759 if (fEdgeFlag != false) { 760 glEdgeFlag (GL_FALSE); 761 fEdgeFlag = false; 762 } 763 } 764 glNormal3d (normals[edgeCount].x(), 765 normals[edgeCount].y(), 766 normals[edgeCount].z()); 767 glVertex3d (vertex[edgeCount].x(), 768 vertex[edgeCount].y(), 769 vertex[edgeCount].z()); 770 } 771 glEnd (); 772 end_of_drawing_through_stencil: 773 774 // and once more to reset the stencil bits... 775 glStencilFunc (GL_ALWAYS, 0, 1); 776 glStencilOp (GL_INVERT, GL_INVERT, GL_INVERT); 777 glDepthFunc (GL_LEQUAL); // to make sure line gets drawn. 778 if (isTransparent) { 779 // Transparent... 780 //glDisable (GL_CULL_FACE); 781 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 782 } else { 783 // Opaque... 784 if (clipping) { 785 //glDisable (GL_CULL_FACE); 786 glPolygonMode (GL_FRONT_AND_BACK, GL_LINE); 787 } else { 788 //glEnable (GL_CULL_FACE); 789 //glCullFace (GL_BACK); 790 glPolygonMode (GL_FRONT, GL_LINE); 791 } 792 } 793 glDisable (GL_LIGHTING); 794 glColor4fv (current_colour); 795 fEdgeFlag = true; 796 glBegin (GL_QUADS); 797 glEdgeFlag (GL_TRUE); 798 fEdgeFlag = true; 799 for (int edgeCount = 0; edgeCount < 4; ++edgeCount) { 800 if (edgeFlag[edgeCount] > 0) { 801 if (fEdgeFlag != true) { 802 glEdgeFlag (GL_TRUE); 803 fEdgeFlag = true; 804 } 805 } else { 806 if (fEdgeFlag != false) { 807 glEdgeFlag (GL_FALSE); 808 fEdgeFlag = false; 809 } 810 } 811 glNormal3d (normals[edgeCount].x(), 812 normals[edgeCount].y(), 813 normals[edgeCount].z()); 814 glVertex3d (vertex[edgeCount].x(), 815 vertex[edgeCount].y(), 816 vertex[edgeCount].z()); 817 } 818 glEnd (); 819 820 glDepthFunc (GL_LEQUAL); // Revert for next facet. 821 fEdgeFlag = true; 822 glBegin (GL_QUADS); // Ready for next facet. GL 823 glEdgeFlag (GL_TRUE); 824 fEdgeFlag = true; 825 // says it ignores incomplete 826 // quadrilaterals, so final empty 827 // glBegin/End sequence should be OK. 828 } 829 } while (notLastFace); 830 831 glEnd (); 832 glDisable (GL_STENCIL_TEST); // Revert to default for next primitive. 833 glDepthMask (GL_TRUE); // Revert to default for next primitive. 834 glDisable (GL_LIGHTING); // Revert to default for next primitive. 835 } 836 837 void G4OpenGLSceneHandler::AddCompound(const G4VTrajectory& traj) { 838 G4VSceneHandler::AddCompound(traj); // For now. 839 } 840 841 void G4OpenGLSceneHandler::AddCompound(const G4VHit& hit) { 842 G4VSceneHandler::AddCompound(hit); // For now. 843 } 844 845 void G4OpenGLSceneHandler::AddCompound(const G4VDigi& digi) { 846 G4VSceneHandler::AddCompound(digi); // For now. 847 } 848 849 void G4OpenGLSceneHandler::AddCompound(const G4THitsMap<G4double>& hits) { 850 G4VSceneHandler::AddCompound(hits); // For now. 851 } 852 853 void G4OpenGLSceneHandler::AddCompound(const G4THitsMap<G4StatDouble>& hits) { 854 G4VSceneHandler::AddCompound(hits); // For now. 855 } 856 857 void G4OpenGLSceneHandler::AddCompound(const G4Mesh& mesh) { 858 StandardSpecialMeshRendering(mesh); 859 } 860