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