Geant4 11.4.0
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),
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
61 if (!fG4OpenGLStoredSceneHandler.fTopPODL ||
64 }
65}
66
68
69 if (
70 (lastVP.GetDrawingStyle () != fVP.GetDrawingStyle ()) ||
71 (lastVP.GetNumberOfCloudPoints() != fVP.GetNumberOfCloudPoints()) ||
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() !=
78 fVP.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 ()) ||
90 (lastVP.GetGlobalMarkerScale() != fVP.GetGlobalMarkerScale()) ||
91 (lastVP.GetGlobalLineWidthScale() != fVP.GetGlobalLineWidthScale()) ||
92 (lastVP.IsMarkerNotHidden () != fVP.IsMarkerNotHidden ()) ||
94 fVP.GetDefaultVisAttributes()->GetColour()) ||
96 fVP.GetDefaultTextVisAttributes()->GetColour()) ||
97 (lastVP.GetBackgroundColour ()!= fVP.GetBackgroundColour ())||
98 (lastVP.IsPicking () != fVP.IsPicking ()) ||
99 (lastVP.GetVisAttributesModifiers() !=
100 fVP.GetVisAttributesModifiers()) ||
101 (lastVP.IsSpecialMeshRendering() !=
102 fVP.IsSpecialMeshRendering()) ||
104 fVP.GetSpecialMeshRenderingOption()) ||
105 (lastVP.GetTransparencyByDepth() != fVP.GetTransparencyByDepth()) ||
106 (lastVP.IsDotsSmooth() != fVP.IsDotsSmooth()) ||
107 (lastVP.GetDotsSize() != fVP.GetDotsSize())
108 )
109 return true;
110
111 if (lastVP.IsDensityCulling () &&
112 (lastVP.GetVisibleDensity () != fVP.GetVisibleDensity ()))
113 return true;
114
115// /**************************************************************
116// If section (DCUT) is implemented locally, comment this out.
117 if (lastVP.IsSection () &&
118 (lastVP.GetSectionPlane () != fVP.GetSectionPlane ()))
119 return true;
120// ***************************************************************/
121
122 /**************************************************************
123 If cutaways are implemented locally, comment this out.
124 if (lastVP.IsCutaway ()) {
125 if (vp.GetCutawayMode() != fVP.GetCutawayMode()) return true;
126 if (lastVP.GetCutawayPlanes ().size () !=
127 fVP.GetCutawayPlanes ().size ()) return true;
128 for (size_t i = 0; i < lastVP.GetCutawayPlanes().size(); ++i)
129 if (lastVP.GetCutawayPlanes()[i] != fVP.GetCutawayPlanes()[i])
130 return true;
131 }
132 ***************************************************************/
133
134 if (lastVP.GetCBDAlgorithmNumber() > 0) {
135 if (lastVP.GetCBDParameters().size() != fVP.GetCBDParameters().size()) return true;
136 else if (lastVP.GetCBDParameters() != fVP.GetCBDParameters()) return true;
137 }
138
139 if (lastVP.IsExplode () &&
140 (lastVP.GetExplodeFactor () != fVP.GetExplodeFactor ()))
141 return true;
142
143 if (lastVP.IsSpecialMeshRendering() &&
144 (lastVP.GetSpecialMeshVolumes() != fVP.GetSpecialMeshVolumes()))
145 return true;
146
147 if (lastVP.GetTransparencyByDepth() > 0. &&
148 lastVP.GetTransparencyByDepthOption() != fVP.GetTransparencyByDepthOption())
149 return true;
150
151 // Time window parameters operate on the existing database so no need
152 // to rebuild even if they change.
153
154 return false;
155}
156
158
159 // We moved these from G4OpenGLViewer to G4ViewParamaters. To avoid
160 // editing many lines below we introduce these convenient aliases.
161#define CONVENIENT_DOUBLE_ALIAS(q) const G4double& f##q = fVP.GetTimeParameters().f##q;
162#define CONVENIENT_BOOL_ALIAS(q) const G4bool& f##q = fVP.GetTimeParameters().f##q;
163 CONVENIENT_DOUBLE_ALIAS(StartTime)
165 CONVENIENT_DOUBLE_ALIAS(FadeFactor)
166 CONVENIENT_BOOL_ALIAS(DisplayHeadTime)
167 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeX)
168 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeY)
169 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeSize)
170 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeRed)
171 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeGreen)
172 CONVENIENT_DOUBLE_ALIAS(DisplayHeadTimeBlue)
173 CONVENIENT_BOOL_ALIAS(DisplayLightFront)
174 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontX)
175 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontY)
176 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontZ)
177 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontT)
178 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontRed)
179 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontGreen)
180 CONVENIENT_DOUBLE_ALIAS(DisplayLightFrontBlue)
181
182 const G4Planes& cutaways = fVP.GetCutawayPlanes();
183 G4bool cutawayUnion = fVP.IsCutaway() &&
184 fVP.GetCutawayMode() == G4ViewParameters::cutawayUnion;
185 const size_t nCutaways = cutawayUnion? cutaways.size(): 1;
186 G4int iPass = 1;
187 G4bool secondPassForTransparencyRequested = false;
188 G4bool thirdPassForNonHiddenMarkersRequested = false;
189 fDepthTestEnable = true;
190 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
191 do {
192 for (size_t iCutaway = 0; iCutaway < nCutaways; ++iCutaway) {
193
194 if (cutawayUnion) {
195 double a[4];
196 a[0] = cutaways[iCutaway].a();
197 a[1] = cutaways[iCutaway].b();
198 a[2] = cutaways[iCutaway].c();
199 a[3] = cutaways[iCutaway].d();
200 glClipPlane (GL_CLIP_PLANE2, a);
201 glEnable (GL_CLIP_PLANE2);
202 }
203
204 G4bool isPicking = fVP.IsPicking();
205
206 for (size_t iPO = 0;
207 iPO < fG4OpenGLStoredSceneHandler.fPOList.size(); ++iPO) {
208 if (POSelected(iPO)) {
210 fG4OpenGLStoredSceneHandler.fPOList[iPO];
211 G4Colour c = po.fColour;
213 const G4bool isTransparent = c.GetAlpha() < 1.;
214 if ( iPass == 1) {
215 if (isTransparent && transparency_enabled) {
216 secondPassForTransparencyRequested = true;
217 continue;
218 }
219 if (po.fMarkerOrPolyline && fVP.IsMarkerNotHidden()) {
220 thirdPassForNonHiddenMarkersRequested = true;
221 continue;
222 }
223 } else if (iPass == 2) { // Second pass for transparency.
224 if (!isTransparent) {
225 continue;
226 }
227 } else { // Third pass for non-hidden markers
228 if (!po.fMarkerOrPolyline) {
229 continue;
230 }
231 }
232 if (isPicking) glLoadName(po.fPickName);
234 glColor4d(c.GetRed(),c.GetGreen(),c.GetBlue(),c.GetAlpha());
235 } else {
236 glColor3d(c.GetRed(),c.GetGreen(),c.GetBlue());
237 }
238 if (po.fMarkerOrPolyline && fVP.IsMarkerNotHidden()) {
239 if (fDepthTestEnable !=false) {
240 glDisable (GL_DEPTH_TEST);
241 fDepthTestEnable = false;
242 }
243 } else {
244 if (fDepthTestEnable !=true) {
245 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
246 fDepthTestEnable = true;
247 }
248 }
249 if (po.fpG4TextPlus) {
250 if (po.fpG4TextPlus->fProcessing2D) {
251 glMatrixMode (GL_PROJECTION);
252 glPushMatrix();
253 glLoadIdentity();
254 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
255 glMatrixMode (GL_MODELVIEW);
256 glPushMatrix();
257 glLoadIdentity();
259 glMultMatrixd (oglt.GetGLMatrix ());
260 // This text is from a PODL. We don't want to create a new PODL.
262 } else {
263 glPushMatrix();
265 glMultMatrixd (oglt.GetGLMatrix ());
266 // This text is from a PODL. We don't want to create a new PODL.
268 glPopMatrix();
269 }
270
271 if (po.fpG4TextPlus->fProcessing2D) {
272 glMatrixMode (GL_PROJECTION);
273 glPopMatrix();
274 glMatrixMode (GL_MODELVIEW);
275 glPopMatrix();
276 }
277 } else {
278 glPushMatrix();
280 glMultMatrixd (oglt.GetGLMatrix ());
281 glCallList (po.fDisplayListId);
282 glPopMatrix();
283 }
284 }
285 }
286
287 G4Transform3D lastMatrixTransform;
288 G4bool first = true;
289
290 for (size_t iTO = 0;
291 iTO < fG4OpenGLStoredSceneHandler.fTOList.size(); ++iTO) {
292 if (TOSelected(iTO)) {
294 fG4OpenGLStoredSceneHandler.fTOList[iTO];
295 const G4Colour& c = to.fColour;
296 const G4bool isTransparent = c.GetAlpha() < 1.;
297 if ( iPass == 1) {
298 if (isTransparent && transparency_enabled) {
299 secondPassForTransparencyRequested = true;
300 continue;
301 }
302 if (to.fMarkerOrPolyline && fVP.IsMarkerNotHidden()) {
303 thirdPassForNonHiddenMarkersRequested = true;
304 continue;
305 }
306 } else if (iPass == 2) { // Second pass for transparency.
307 if (!isTransparent) {
308 continue;
309 }
310 } else { // Third pass for non-hidden markers
311 if (!to.fMarkerOrPolyline) {
312 continue;
313 }
314 }
315 if (to.fMarkerOrPolyline && fVP.IsMarkerNotHidden()) {
316 if (fDepthTestEnable !=false) {
317 glDisable (GL_DEPTH_TEST);
318 fDepthTestEnable = false;
319 }
320 } else {
321 if (fDepthTestEnable !=true) {
322 glEnable (GL_DEPTH_TEST); glDepthFunc (GL_LEQUAL);
323 fDepthTestEnable = true;
324 }
325 }
326 if (to.fEndTime >= fStartTime && to.fStartTime <= fEndTime) {
327 if (fVP.IsPicking()) glLoadName(to.fPickName);
328 if (to.fpG4TextPlus) {
329 if (to.fpG4TextPlus->fProcessing2D) {
330 glMatrixMode (GL_PROJECTION);
331 glPushMatrix();
332 glLoadIdentity();
333 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
334 glMatrixMode (GL_MODELVIEW);
335 glPushMatrix();
336 glLoadIdentity();
337 }
339 glMultMatrixd (oglt.GetGLMatrix ());
340 // This text is from a TODL. We don't want to create a new TODL.
342 if (to.fpG4TextPlus->fProcessing2D) {
343 glMatrixMode (GL_PROJECTION);
344 glPopMatrix();
345 glMatrixMode (GL_MODELVIEW);
346 glPopMatrix();
347 }
348 } else {
349 if (to.fTransform != lastMatrixTransform) {
350 if (! first) {
351 glPopMatrix();
352 }
353 first = false;
354 glPushMatrix();
356 glMultMatrixd (oglt.GetGLMatrix ());
357 }
358 const G4Colour& cc = to.fColour;
359 if (fFadeFactor > 0. && to.fEndTime < fEndTime) {
360 // Brightness scaling factor
361 G4double bsf = 1. - fFadeFactor *
362 ((fEndTime - to.fEndTime) / (fEndTime - fStartTime));
363 const G4Colour& bg = fVP.GetBackgroundColour();
365 glColor4d
366 (bsf * cc.GetRed() + (1. - bsf) * bg.GetRed(),
367 bsf * cc.GetGreen() + (1. - bsf) * bg.GetGreen(),
368 bsf * cc.GetBlue() + (1. - bsf) * bg.GetBlue(),
369 bsf * cc.GetAlpha() + (1. - bsf) * bg.GetAlpha());
370 } else {
371 glColor3d
372 (bsf * cc.GetRed() + (1. - bsf) * bg.GetRed(),
373 bsf * cc.GetGreen() + (1. - bsf) * bg.GetGreen(),
374 bsf * cc.GetBlue() + (1. - bsf) * bg.GetBlue());
375 }
376 } else {
378 glColor4d(cc.GetRed(),cc.GetGreen(),cc.GetBlue(),cc.GetAlpha());
379 } else {
380 glColor3d(cc.GetRed(),cc.GetGreen(),cc.GetBlue());
381 }
382 }
383 glCallList (to.fDisplayListId);
384 }
385 if (to.fTransform != lastMatrixTransform) {
386 lastMatrixTransform = to.fTransform;
387 }
388 }
389 }
390 }
391 if (! first) {
392 glPopMatrix();
393 }
394
395 if (cutawayUnion) glDisable (GL_CLIP_PLANE2);
396 } // iCutaway
397
398 if (iPass == 2) secondPassForTransparencyRequested = false; // Done.
399 if (iPass == 3) thirdPassForNonHiddenMarkersRequested = false; // Done.
400
401 if (secondPassForTransparencyRequested) iPass = 2;
402 else if (thirdPassForNonHiddenMarkersRequested) iPass = 3;
403 else break;
404
405 } while (true);
406
407 // Display time at "head" of time range, which is fEndTime...
408 if (fDisplayHeadTime && fEndTime < G4VisAttributes::fVeryLongTime) {
409 glMatrixMode (GL_PROJECTION);
410 glPushMatrix();
411 glLoadIdentity();
412 g4GlOrtho (-1., 1., -1., 1., -G4OPENGL_FLT_BIG, G4OPENGL_FLT_BIG);
413 glMatrixMode (GL_MODELVIEW);
414 glPushMatrix();
415 glLoadIdentity();
416 G4Text headTimeText(G4BestUnit(fEndTime,"Time"),
417 G4Point3D(fDisplayHeadTimeX, fDisplayHeadTimeY, 0.));
418 headTimeText.SetScreenSize(fDisplayHeadTimeSize);
420 (fDisplayHeadTimeRed,
421 fDisplayHeadTimeGreen,
422 fDisplayHeadTimeBlue));
423 headTimeText.SetVisAttributes(&visAtts);
424 AddPrimitiveForASingleFrame(headTimeText);
425 glMatrixMode (GL_PROJECTION);
426 glPopMatrix();
427 glMatrixMode (GL_MODELVIEW);
428 glPopMatrix();
429 }
430
431 // Display light front...
432 if (fDisplayLightFront && fEndTime < G4VisAttributes::fVeryLongTime) {
433 G4double lightFrontRadius = (fEndTime - fDisplayLightFrontT) * c_light;
434 if (lightFrontRadius > 0.) {
435 G4Point3D lightFrontCentre(fDisplayLightFrontX, fDisplayLightFrontY, fDisplayLightFrontZ);
436 G4Point3D circleCentre = lightFrontCentre;
437 G4double circleRadius = lightFrontRadius;
438 if (fVP.GetFieldHalfAngle() > 0.) {
439 // Perspective view. Find horizon centre and radius...
440 G4Point3D targetPoint = fSceneHandler.GetScene()->GetStandardTargetPoint() +
441 fVP.GetCurrentTargetPoint();
442 G4double sceneRadius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
443 if(sceneRadius <= 0.) sceneRadius = 1.;
444 G4double cameraDistance = fVP.GetCameraDistance(sceneRadius);
445 G4Point3D cameraPosition = targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
446 G4Vector3D lightFrontToCameraDirection = cameraPosition - lightFrontCentre;
447 G4double lightFrontCentreDistance = lightFrontToCameraDirection.mag();
448 /*
449 G4cout << "cameraPosition: " << cameraPosition
450 << ", lightFrontCentre: " << lightFrontCentre
451 << ", lightFrontRadius: " << lightFrontRadius
452 << ", lightFrontCentreDistance: " << lightFrontCentreDistance
453 << ", dot: " << lightFrontToCameraDirection * fVP.GetViewpointDirection()
454 << G4endl;
455 */
456 if (lightFrontToCameraDirection * fVP.GetViewpointDirection() > 0. && lightFrontRadius < lightFrontCentreDistance) {
457 // Light front in front of camera...
458 G4double sineHorizonAngle = lightFrontRadius / lightFrontCentreDistance;
459 circleCentre = lightFrontCentre + (lightFrontRadius * sineHorizonAngle) * lightFrontToCameraDirection.unit();
460 circleRadius = lightFrontRadius * std::sqrt(1. - std::pow(sineHorizonAngle, 2));
461 /*
462 G4cout << "sineHorizonAngle: " << sineHorizonAngle
463 << ", circleCentre: " << circleCentre
464 << ", circleRadius: " << circleRadius
465 << G4endl;
466 */
467 } else {
468 circleRadius = -1.;
469 }
470 }
471 if (circleRadius > 0.) {
472 G4Circle lightFront(circleCentre);
473 lightFront.SetWorldRadius(circleRadius);
475 (fDisplayLightFrontRed,
476 fDisplayLightFrontGreen,
477 fDisplayLightFrontBlue));
478 lightFront.SetVisAttributes(visAtts);
479 AddPrimitiveForASingleFrame(lightFront);
480 }
481 }
482 }
483}
484
486{
487 // We don't want this to get into a display list or a TODL or a PODL so
488 // use the fMemoryForDisplayLists flag.
489 fG4OpenGLStoredSceneHandler.fDoNotUseDisplayList = true;
490 fG4OpenGLStoredSceneHandler.G4OpenGLStoredSceneHandler::AddPrimitive(text);
491 fG4OpenGLStoredSceneHandler.fDoNotUseDisplayList = false;
492}
493
495{
496 // We don't want this to get into a display list or a TODL or a PODL so
497 // use the fMemoryForDisplayLists flag.
498 fG4OpenGLStoredSceneHandler.fDoNotUseDisplayList = true;
499 fG4OpenGLStoredSceneHandler.G4OpenGLStoredSceneHandler::AddPrimitive(circle);
500 fG4OpenGLStoredSceneHandler.fDoNotUseDisplayList = false;
501}
#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)
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
std::vector< G4Plane3D > G4Planes
G4double GetBlue() const
Definition G4Colour.hh:172
G4double GetAlpha() const
Definition G4Colour.hh:173
G4double GetRed() const
Definition G4Colour.hh:170
G4double GetGreen() const
Definition G4Colour.hh:171
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()
friend class G4OpenGLStoredSceneHandler
G4OpenGLViewer(G4OpenGLSceneHandler &scene)
G4bool transparency_enabled
void g4GlOrtho(GLdouble left, GLdouble right, GLdouble bottom, GLdouble top, GLdouble near, GLdouble far)
void SetWorldRadius(G4double)
void SetScreenSize(G4double)
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:268
void NeedKernelVisit()
Definition G4VViewer.cc:85
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:273
G4ViewParameters fVP
Definition G4VViewer.hh:272
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition G4VViewer.cc:49
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsMarkerNotHidden() const
G4double GetGlobalLineWidthScale() const
const G4Colour & GetBackgroundColour() const
G4bool IsSection() const
G4bool IsPicking() const
G4bool IsDotsSmooth() 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
G4int GetTransparencyByDepthOption() const
G4double GetGlobalMarkerScale() const
G4bool IsCullingInvisible() const
const G4VisAttributes * GetDefaultVisAttributes() const
G4bool IsDensityCulling() const
G4double GetVisibleDensity() const
SMROption GetSpecialMeshRenderingOption() const
G4double GetDotsSize() const
G4bool IsCullingCovered() const
const G4Plane3D & GetSectionPlane() const
DrawingStyle GetDrawingStyle() const
G4double GetTransparencyByDepth() const
G4bool IsAuxEdgeVisible() const
const G4Colour & GetColour() const
static constexpr G4double fVeryLongTime
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:108
BasicVector3D< T > unit() const