Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpenGLStoredViewer.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//
28//
29// Andrew Walkden 7th February 1997
30// Class G4OpenGLStoredViewer : Encapsulates the `storedness' of
31// an OpenGL view, for inheritance by
32// derived (X, Xm...) classes.
33
35
38#include "G4Text.hh"
39#include "G4Circle.hh"
40#include "G4UnitsTable.hh"
41#include "G4Scene.hh"
43
45(G4OpenGLStoredSceneHandler& sceneHandler):
46G4VViewer (sceneHandler, -1),
47G4OpenGLViewer (sceneHandler),
48fG4OpenGLStoredSceneHandler (sceneHandler),
49fDepthTestEnable(true)
50{
51 fLastVP = fDefaultVP; // Update in sub-class after KernelVisitDecision
52}
53
55
57
58 // If there's a significant difference with the last view parameters
59 // of either the scene handler or this viewer, trigger a rebuild.
60
64 }
65}
66
68
69 if (
70 (lastVP.GetDrawingStyle () != fVP.GetDrawingStyle ()) ||
72 (lastVP.IsAuxEdgeVisible () != fVP.IsAuxEdgeVisible ()) ||
73 (lastVP.IsCulling () != fVP.IsCulling ()) ||
74 (lastVP.IsCullingInvisible () != fVP.IsCullingInvisible ()) ||
75 (lastVP.IsDensityCulling () != fVP.IsDensityCulling ()) ||
76 (lastVP.IsCullingCovered () != fVP.IsCullingCovered ()) ||
77 (lastVP.GetCBDAlgorithmNumber() !=
79 // Note: Section and Cutaway can reveal back-facing faces. If
80 // backface culling is implemented, the image can look strange because
81 // the back-facing faces are not there. For the moment, we have disabled
82 // (commented out) backface culling (it seems not to affect performance -
83 // in fact, performance seems to improve), so there is no problem.
84 (lastVP.IsSection () != fVP.IsSection ()) ||
85 // Section (DCUT) is NOT implemented locally so we need to visit the kernel.
86 // (lastVP.IsCutaway () != fVP.IsCutaway ()) ||
87 // Cutaways are implemented locally so we do not need to visit the kernel.
88 (lastVP.IsExplode () != fVP.IsExplode ()) ||
89 (lastVP.GetNoOfSides () != fVP.GetNoOfSides ()) ||
92 (lastVP.IsMarkerNotHidden () != fVP.IsMarkerNotHidden ()) ||
98 (lastVP.IsPicking () != fVP.IsPicking ()) ||
99 (lastVP.GetVisAttributesModifiers() !=
101 (lastVP.IsSpecialMeshRendering() !=
105 )
106 return true;
107
108 if (lastVP.IsDensityCulling () &&
109 (lastVP.GetVisibleDensity () != fVP.GetVisibleDensity ()))
110 return true;
111
112// /**************************************************************
113// If section (DCUT) is implemented locally, comment this out.
114 if (lastVP.IsSection () &&
115 (lastVP.GetSectionPlane () != fVP.GetSectionPlane ()))
116 return true;
117// ***************************************************************/
118
119 /**************************************************************
120 If cutaways are implemented locally, comment this out.
121 if (lastVP.IsCutaway ()) {
122 if (vp.GetCutawayMode() != fVP.GetCutawayMode()) return true;
123 if (lastVP.GetCutawayPlanes ().size () !=
124 fVP.GetCutawayPlanes ().size ()) return true;
125 for (size_t i = 0; i < lastVP.GetCutawayPlanes().size(); ++i)
126 if (lastVP.GetCutawayPlanes()[i] != fVP.GetCutawayPlanes()[i])
127 return true;
128 }
129 ***************************************************************/
130
131 if (lastVP.GetCBDAlgorithmNumber() > 0) {
132 if (lastVP.GetCBDParameters().size() != fVP.GetCBDParameters().size()) return true;
133 else if (lastVP.GetCBDParameters() != fVP.GetCBDParameters()) return true;
134 }
135
136 if (lastVP.IsExplode () &&
137 (lastVP.GetExplodeFactor () != fVP.GetExplodeFactor ()))
138 return true;
139
140 if (lastVP.IsSpecialMeshRendering() &&
142 return true;
143
144 // Time window parameters operate on the existing database so no need
145 // to rebuild even if they change.
146
147 return false;
148}
149
151
152 // We moved these from G4OpenGLViewer to G4ViewParamaters. To avoid
153 // editing many lines below we introduce these convenient aliases.
154#define CONVENIENT_DOUBLE_ALIAS(q) const G4double& f##q = fVP.Get##q();
155#define CONVENIENT_BOOL_ALIAS(q) const G4bool& f##q = fVP.Is##q();
156 CONVENIENT_DOUBLE_ALIAS(StartTime)
158 CONVENIENT_DOUBLE_ALIAS(FadeFactor)
159 CONVENIENT_BOOL_ALIAS(DisplayHeadTime)
160 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeX)
161 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeY)
162 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeSize)
163 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeRed)
164 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeGreen)
165 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeBlue)
166 CONVENIENT_BOOL_ALIAS(DisplayLightFront)
167 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontX)
168 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontY)
169 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontZ)
170 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontT)
171 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontRed)
172 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontGreen)
173 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontBlue)
174
175 const G4Planes& cutaways = fVP.GetCutawayPlanes();
176 G4bool cutawayUnion = fVP.IsCutaway() &&
178 const size_t nCutaways = cutawayUnion? cutaways.size(): 1;
179 G4int iPass = 1;
180 G4bool secondPassForTransparencyRequested = false;
181 G4bool thirdPassForNonHiddenMarkersRequested = false;
182 fDepthTestEnable = true;
183 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
184 do {
185 for (size_t iCutaway = 0; iCutaway < nCutaways; ++iCutaway) {
186
187 if (cutawayUnion) {
188 double a[4];
189 a[0] = cutaways[iCutaway].a();
190 a[1] = cutaways[iCutaway].b();
191 a[2] = cutaways[iCutaway].c();
192 a[3] = cutaways[iCutaway].d();
193 glClipPlane (GL_CLIP_PLANE2, a);
194 glEnable (GL_CLIP_PLANE2);
195 }
196
197 G4bool isPicking = fVP.IsPicking();
198
199 for (size_t iPO = 0;
200 iPO < fG4OpenGLStoredSceneHandler.fPOList.size(); ++iPO) {
201 if (POSelected(iPO)) {
204 G4Colour c = po.fColour;
206 const G4bool isTransparent = c.GetAlpha() < 1.;
207 if ( iPass == 1) {
208 if (isTransparent && transparency_enabled) {
209 secondPassForTransparencyRequested = true;
210 continue;
211 }
213 thirdPassForNonHiddenMarkersRequested = true;
214 continue;
215 }
216 } else if (iPass == 2) { // Second pass for transparency.
217 if (!isTransparent) {
218 continue;
219 }
220 } else { // Third pass for non-hidden markers
221 if (!po.fMarkerOrPolyline) {
222 continue;
223 }
224 }
225 if (isPicking) glLoadName(po.fPickName);
227 glColor4d(c.GetRed(),c.GetGreen(),c.GetBlue(),c.GetAlpha());
228 } else {
229 glColor3d(c.GetRed(),c.GetGreen(),c.GetBlue());
230 }
232 if (fDepthTestEnable !=false) {
233 glDisable (GL_DEPTH_TEST);
234 fDepthTestEnable = false;
235 }
236 } else {
237 if (fDepthTestEnable !=true) {
238 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
239 fDepthTestEnable = true;
240 }
241 }
242 if (po.fpG4TextPlus) {
243 if (po.fpG4TextPlus->fProcessing2D) {
244 glMatrixMode (GL_PROJECTION);
245 glPushMatrix();
246 glLoadIdentity();
247 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
248 glMatrixMode (GL_MODELVIEW);
249 glPushMatrix();
250 glLoadIdentity();
252 glMultMatrixd (oglt.GetGLMatrix ());
253 // This text is from a PODL. We don't want to create a new PODL.
255 } else {
256 glPushMatrix();
258 glMultMatrixd (oglt.GetGLMatrix ());
259 // This text is from a PODL. We don't want to create a new PODL.
261 glPopMatrix();
262 }
263
264 if (po.fpG4TextPlus->fProcessing2D) {
265 glMatrixMode (GL_PROJECTION);
266 glPopMatrix();
267 glMatrixMode (GL_MODELVIEW);
268 glPopMatrix();
269 }
270 } else {
271 glPushMatrix();
273 glMultMatrixd (oglt.GetGLMatrix ());
274 glCallList (po.fDisplayListId);
275 glPopMatrix();
276 }
277 }
278 }
279
280 G4Transform3D lastMatrixTransform;
281 G4bool first = true;
282
283 for (size_t iTO = 0;
284 iTO < fG4OpenGLStoredSceneHandler.fTOList.size(); ++iTO) {
285 if (TOSelected(iTO)) {
288 const G4Colour& c = to.fColour;
289 const G4bool isTransparent = c.GetAlpha() < 1.;
290 if ( iPass == 1) {
291 if (isTransparent && transparency_enabled) {
292 secondPassForTransparencyRequested = true;
293 continue;
294 }
296 thirdPassForNonHiddenMarkersRequested = true;
297 continue;
298 }
299 } else if (iPass == 2) { // Second pass for transparency.
300 if (!isTransparent) {
301 continue;
302 }
303 } else { // Third pass for non-hidden markers
304 if (!to.fMarkerOrPolyline) {
305 continue;
306 }
307 }
309 if (fDepthTestEnable !=false) {
310 glDisable (GL_DEPTH_TEST);
311 fDepthTestEnable = false;
312 }
313 } else {
314 if (fDepthTestEnable !=true) {
315 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
316 fDepthTestEnable = true;
317 }
318 }
319 if (to.fEndTime >= fStartTime && to.fStartTime <= fEndTime) {
320 if (fVP.IsPicking()) glLoadName(to.fPickName);
321 if (to.fpG4TextPlus) {
322 if (to.fpG4TextPlus->fProcessing2D) {
323 glMatrixMode (GL_PROJECTION);
324 glPushMatrix();
325 glLoadIdentity();
326 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
327 glMatrixMode (GL_MODELVIEW);
328 glPushMatrix();
329 glLoadIdentity();
330 }
332 glMultMatrixd (oglt.GetGLMatrix ());
333 // This text is from a TODL. We don't want to create a new TODL.
335 if (to.fpG4TextPlus->fProcessing2D) {
336 glMatrixMode (GL_PROJECTION);
337 glPopMatrix();
338 glMatrixMode (GL_MODELVIEW);
339 glPopMatrix();
340 }
341 } else {
342 if (to.fTransform != lastMatrixTransform) {
343 if (! first) {
344 glPopMatrix();
345 }
346 first = false;
347 glPushMatrix();
349 glMultMatrixd (oglt.GetGLMatrix ());
350 }
351 const G4Colour& cc = to.fColour;
352 if (fFadeFactor > 0. && to.fEndTime < fEndTime) {
353 // Brightness scaling factor
354 G4double bsf = 1. - fFadeFactor *
355 ((fEndTime - to.fEndTime) / (fEndTime - fStartTime));
356 const G4Colour& bg = fVP.GetBackgroundColour();
358 glColor4d
359 (bsf * cc.GetRed() + (1. - bsf) * bg.GetRed(),
360 bsf * cc.GetGreen() + (1. - bsf) * bg.GetGreen(),
361 bsf * cc.GetBlue() + (1. - bsf) * bg.GetBlue(),
362 bsf * cc.GetAlpha() + (1. - bsf) * bg.GetAlpha());
363 } else {
364 glColor3d
365 (bsf * cc.GetRed() + (1. - bsf) * bg.GetRed(),
366 bsf * cc.GetGreen() + (1. - bsf) * bg.GetGreen(),
367 bsf * cc.GetBlue() + (1. - bsf) * bg.GetBlue());
368 }
369 } else {
371 glColor4d(cc.GetRed(),cc.GetGreen(),cc.GetBlue(),cc.GetAlpha());
372 } else {
373 glColor3d(cc.GetRed(),cc.GetGreen(),cc.GetBlue());
374 }
375 }
376 glCallList (to.fDisplayListId);
377 }
378 if (to.fTransform != lastMatrixTransform) {
379 lastMatrixTransform = to.fTransform;
380 }
381 }
382 }
383 }
384 if (! first) {
385 glPopMatrix();
386 }
387
388 if (cutawayUnion) glDisable (GL_CLIP_PLANE2);
389 } // iCutaway
390
391 if (iPass == 2) secondPassForTransparencyRequested = false; // Done.
392 if (iPass == 3) thirdPassForNonHiddenMarkersRequested = false; // Done.
393
394 if (secondPassForTransparencyRequested) iPass = 2;
395 else if (thirdPassForNonHiddenMarkersRequested) iPass = 3;
396 else break;
397
398 } while (true);
399
400 // Display time at "head" of time range, which is fEndTime...
401 if (fDisplayHeadTime && fEndTime < G4VisAttributes::fVeryLongTime) {
402 glMatrixMode (GL_PROJECTION);
403 glPushMatrix();
404 glLoadIdentity();
405 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
406 glMatrixMode (GL_MODELVIEW);
407 glPushMatrix();
408 glLoadIdentity();
409 G4Text headTimeText(G4BestUnit(fEndTime,"Time"),
410 G4Point3D(fDisplayHeadTimeX, fDisplayHeadTimeY, 0.));
411 headTimeText.SetScreenSize(fDisplayHeadTimeSize);
413 (fDisplayHeadTimeRed,
414 fDisplayHeadTimeGreen,
415 fDisplayHeadTimeBlue));
416 headTimeText.SetVisAttributes(&visAtts);
417 AddPrimitiveForASingleFrame(headTimeText);
418 glMatrixMode (GL_PROJECTION);
419 glPopMatrix();
420 glMatrixMode (GL_MODELVIEW);
421 glPopMatrix();
422 }
423
424 // Display light front...
425 if (fDisplayLightFront && fEndTime < G4VisAttributes::fVeryLongTime) {
426 G4double lightFrontRadius = (fEndTime - fDisplayLightFrontT) * c_light;
427 if (lightFrontRadius > 0.) {
428 G4Point3D lightFrontCentre(fDisplayLightFrontX, fDisplayLightFrontY, fDisplayLightFrontZ);
429 G4Point3D circleCentre = lightFrontCentre;
430 G4double circleRadius = lightFrontRadius;
431 if (fVP.GetFieldHalfAngle() > 0.) {
432 // Perspective view. Find horizon centre and radius...
436 if(sceneRadius <= 0.) sceneRadius = 1.;
437 G4double cameraDistance = fVP.GetCameraDistance(sceneRadius);
438 G4Point3D cameraPosition = targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
439 G4Vector3D lightFrontToCameraDirection = cameraPosition - lightFrontCentre;
440 G4double lightFrontCentreDistance = lightFrontToCameraDirection.mag();
441 /*
442 G4cout << "cameraPosition: " << cameraPosition
443 << ", lightFrontCentre: " << lightFrontCentre
444 << ", lightFrontRadius: " << lightFrontRadius
445 << ", lightFrontCentreDistance: " << lightFrontCentreDistance
446 << ", dot: " << lightFrontToCameraDirection * fVP.GetViewpointDirection()
447 << G4endl;
448 */
449 if (lightFrontToCameraDirection * fVP.GetViewpointDirection() > 0. && lightFrontRadius < lightFrontCentreDistance) {
450 // Light front in front of camera...
451 G4double sineHorizonAngle = lightFrontRadius / lightFrontCentreDistance;
452 circleCentre = lightFrontCentre + (lightFrontRadius * sineHorizonAngle) * lightFrontToCameraDirection.unit();
453 circleRadius = lightFrontRadius * std::sqrt(1. - std::pow(sineHorizonAngle, 2));
454 /*
455 G4cout << "sineHorizonAngle: " << sineHorizonAngle
456 << ", circleCentre: " << circleCentre
457 << ", circleRadius: " << circleRadius
458 << G4endl;
459 */
460 } else {
461 circleRadius = -1.;
462 }
463 }
464 if (circleRadius > 0.) {
465 G4Circle lightFront(circleCentre);
466 lightFront.SetWorldRadius(circleRadius);
468 (fDisplayLightFrontRed,
469 fDisplayLightFrontGreen,
470 fDisplayLightFrontBlue));
471 lightFront.SetVisAttributes(visAtts);
472 AddPrimitiveForASingleFrame(lightFront);
473 }
474 }
475 }
476}
477
479{
480 // We don't want this to get into a display list or a TODL or a PODL so
481 // use the fMemoryForDisplayLists flag.
483 fG4OpenGLStoredSceneHandler.G4OpenGLStoredSceneHandler::AddPrimitive(text);
485}
486
488{
489 // We don't want this to get into a display list or a TODL or a PODL so
490 // use the fMemoryForDisplayLists flag.
492 fG4OpenGLStoredSceneHandler.G4OpenGLStoredSceneHandler::AddPrimitive(circle);
494}
#define CONVENIENT_BOOL_ALIAS(q)
#define CONVENIENT_DOUBLE_ALIAS(q)
#define G4OPENGL_FLT_BIG
Definition G4OpenGL.hh:81
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
#define G4BestUnit(a, b)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
std::vector< G4Plane3D > G4Planes
G4double GetBlue() const
Definition G4Colour.hh:154
G4double GetAlpha() const
Definition G4Colour.hh:155
G4double GetRed() const
Definition G4Colour.hh:152
G4double GetGreen() const
Definition G4Colour.hh:153
virtual G4bool CompareForKernelVisit(G4ViewParameters &)
virtual G4bool TOSelected(size_t)
G4OpenGLStoredViewer(G4OpenGLStoredSceneHandler &scene)
virtual G4bool POSelected(size_t)
virtual void DisplayTimePOColourModification(G4Colour &, size_t)
void AddPrimitiveForASingleFrame(const G4Text &text)
G4OpenGLStoredSceneHandler & fG4OpenGLStoredSceneHandler
const GLdouble * GetGLMatrix()
G4bool transparency_enabled
void g4GlOrtho(GLdouble left, GLdouble right, GLdouble bottom, GLdouble top, GLdouble near, GLdouble far)
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
void SetWorldRadius(G4double)
void SetScreenSize(G4double)
G4Scene * GetScene() const
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:253
void NeedKernelVisit()
Definition G4VViewer.cc:81
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:258
G4ViewParameters fVP
Definition G4VViewer.hh:257
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
CutawayMode GetCutawayMode() const
G4double GetCameraDistance(G4double radius) const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsMarkerNotHidden() const
G4double GetGlobalLineWidthScale() const
G4bool IsCutaway() const
const G4Colour & GetBackgroundColour() const
const G4Vector3D & GetViewpointDirection() const
G4bool IsSection() const
const G4Point3D & GetCurrentTargetPoint() const
G4bool IsPicking() const
G4double GetFieldHalfAngle() const
G4bool IsCulling() const
const G4VisAttributes * GetDefaultTextVisAttributes() const
G4bool IsExplode() const
const std::vector< G4double > & GetCBDParameters() const
G4int GetCBDAlgorithmNumber() const
const std::vector< G4ModelingParameters::PVNameCopyNo > & GetSpecialMeshVolumes() const
G4double GetGlobalMarkerScale() const
G4bool IsCullingInvisible() const
const G4VisAttributes * GetDefaultVisAttributes() const
const G4Planes & GetCutawayPlanes() const
G4bool IsDensityCulling() const
G4double GetVisibleDensity() const
SMROption GetSpecialMeshRenderingOption() const
G4bool IsCullingCovered() const
const G4Plane3D & GetSectionPlane() const
DrawingStyle GetDrawingStyle() const
G4bool IsAuxEdgeVisible() const
const G4Colour & GetColour() const
static constexpr G4double fVeryLongTime
G4double GetExtentRadius() const
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
BasicVector3D< T > unit() const