Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VSceneHandler Class Referenceabstract

#include <G4VSceneHandler.hh>

+ Inheritance diagram for G4VSceneHandler:

Classes

struct  NameAndVisAtts
 
class  PseudoSceneFor3DRectMeshPositions
 
class  PseudoSceneForTetVertices
 

Public Types

enum  MarkerSizeType { world , screen }
 

Public Member Functions

 G4VSceneHandler (G4VGraphicsSystem &system, G4int id, const G4String &name="")
 
virtual ~G4VSceneHandler ()
 
virtual void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &)
 
virtual void PostAddSolid ()
 
virtual void AddSolid (const G4Box &)
 
virtual void AddSolid (const G4Cons &)
 
virtual void AddSolid (const G4Orb &)
 
virtual void AddSolid (const G4Para &)
 
virtual void AddSolid (const G4Sphere &)
 
virtual void AddSolid (const G4Torus &)
 
virtual void AddSolid (const G4Trap &)
 
virtual void AddSolid (const G4Trd &)
 
virtual void AddSolid (const G4Tubs &)
 
virtual void AddSolid (const G4Ellipsoid &)
 
virtual void AddSolid (const G4Polycone &)
 
virtual void AddSolid (const G4Polyhedra &)
 
virtual void AddSolid (const G4TessellatedSolid &)
 
virtual void AddSolid (const G4VSolid &)
 
virtual void AddCompound (const G4VTrajectory &)
 
virtual void AddCompound (const G4VHit &)
 
virtual void AddCompound (const G4VDigi &)
 
virtual void AddCompound (const G4THitsMap< G4double > &)
 
virtual void AddCompound (const G4THitsMap< G4StatDouble > &)
 
virtual void AddCompound (const G4Mesh &)
 
virtual void BeginModeling ()
 
virtual void EndModeling ()
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation=G4Transform3D())
 
virtual void EndPrimitives ()
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation=G4Transform3D())
 
virtual void EndPrimitives2D ()
 
virtual void AddPrimitive (const G4Polyline &)=0
 
virtual void AddPrimitive (const G4Text &)=0
 
virtual void AddPrimitive (const G4Circle &)=0
 
virtual void AddPrimitive (const G4Square &)=0
 
virtual void AddPrimitive (const G4Polymarker &)
 
virtual void AddPrimitive (const G4Polyhedron &)=0
 
virtual void AddPrimitive (const G4Plotter &)
 
virtual const G4VisExtentGetExtent () const
 
const G4StringGetName () const
 
G4int GetSceneHandlerId () const
 
G4int GetViewCount () const
 
G4VGraphicsSystemGetGraphicsSystem () const
 
G4SceneGetScene () const
 
const G4ViewerListGetViewerList () const
 
G4VModelGetModel () const
 
G4VViewerGetCurrentViewer () const
 
G4bool GetMarkForClearingTransientStore () const
 
G4bool IsReadyForTransients () const
 
G4bool GetTransientsDrawnThisEvent () const
 
G4bool GetTransientsDrawnThisRun () const
 
const G4Transform3DGetObjectTransformation () const
 
void SetName (const G4String &)
 
void SetCurrentViewer (G4VViewer *)
 
virtual void SetScene (G4Scene *)
 
G4ViewerListSetViewerList ()
 
void SetModel (G4VModel *)
 
void SetMarkForClearingTransientStore (G4bool)
 
void SetTransientsDrawnThisEvent (G4bool)
 
void SetTransientsDrawnThisRun (G4bool)
 
void SetObjectTransformation (const G4Transform3D &)
 
const G4ColourGetColour ()
 
const G4ColourGetColor ()
 
const G4ColourGetColour (const G4Visible &)
 
const G4ColourGetColor (const G4Visible &)
 
const G4ColourGetTextColour (const G4Text &)
 
const G4ColourGetTextColor (const G4Text &)
 
G4double GetLineWidth (const G4VisAttributes *)
 
G4ViewParameters::DrawingStyle GetDrawingStyle (const G4VisAttributes *)
 
G4int GetNumberOfCloudPoints (const G4VisAttributes *) const
 
G4bool GetAuxEdgeVisible (const G4VisAttributes *)
 
G4int GetNoOfSides (const G4VisAttributes *)
 
G4double GetMarkerSize (const G4VMarker &, MarkerSizeType &)
 
G4double GetMarkerDiameter (const G4VMarker &, MarkerSizeType &)
 
G4double GetMarkerRadius (const G4VMarker &, MarkerSizeType &)
 
G4ModelingParametersCreateModelingParameters ()
 
void DrawEvent (const G4Event *)
 
void DrawEndOfRunModels ()
 
template<class T >
void AddSolidT (const T &solid)
 
template<class T >
void AddSolidWithAuxiliaryEdges (const T &solid)
 
G4int IncrementViewCount ()
 
virtual void ClearStore ()
 
virtual void ClearTransientStore ()
 
void AddViewerToList (G4VViewer *pView)
 
void RemoveViewerFromList (G4VViewer *pView)
 
- Public Member Functions inherited from G4VGraphicsScene
 G4VGraphicsScene ()
 
virtual ~G4VGraphicsScene ()
 
virtual void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
 
virtual void PostAddSolid ()=0
 
virtual void AddSolid (const G4Box &)=0
 
virtual void AddSolid (const G4Cons &)=0
 
virtual void AddSolid (const G4Orb &)=0
 
virtual void AddSolid (const G4Para &)=0
 
virtual void AddSolid (const G4Sphere &)=0
 
virtual void AddSolid (const G4Torus &)=0
 
virtual void AddSolid (const G4Trap &)=0
 
virtual void AddSolid (const G4Trd &)=0
 
virtual void AddSolid (const G4Tubs &)=0
 
virtual void AddSolid (const G4Ellipsoid &)=0
 
virtual void AddSolid (const G4Polycone &)=0
 
virtual void AddSolid (const G4Polyhedra &)=0
 
virtual void AddSolid (const G4TessellatedSolid &)=0
 
virtual void AddSolid (const G4VSolid &)=0
 
virtual void AddCompound (const G4VTrajectory &)=0
 
virtual void AddCompound (const G4VHit &)=0
 
virtual void AddCompound (const G4VDigi &)=0
 
virtual void AddCompound (const G4THitsMap< G4double > &)=0
 
virtual void AddCompound (const G4THitsMap< G4StatDouble > &)=0
 
virtual void AddCompound (const G4Mesh &)=0
 
virtual void BeginPrimitives (const G4Transform3D &objectTransformation=G4Transform3D())=0
 
virtual void EndPrimitives ()=0
 
virtual void BeginPrimitives2D (const G4Transform3D &objectTransformation=G4Transform3D())=0
 
virtual void EndPrimitives2D ()=0
 
virtual void AddPrimitive (const G4Polyline &)=0
 
virtual void AddPrimitive (const G4Text &)=0
 
virtual void AddPrimitive (const G4Circle &)=0
 
virtual void AddPrimitive (const G4Square &)=0
 
virtual void AddPrimitive (const G4Polymarker &)=0
 
virtual void AddPrimitive (const G4Polyhedron &)=0
 
virtual void AddPrimitive (const G4Plotter &)=0
 
virtual const G4VisExtentGetExtent () const
 

Protected Member Functions

virtual void ProcessScene ()
 
virtual void RequestPrimitives (const G4VSolid &solid)
 
virtual G4DisplacedSolidCreateSectionSolid ()
 
virtual G4DisplacedSolidCreateCutawaySolid ()
 
void LoadAtts (const G4Visible &, G4AttHolder *)
 
void StandardSpecialMeshRendering (const G4Mesh &)
 
void Draw3DRectMeshAsDots (const G4Mesh &)
 
void Draw3DRectMeshAsSurfaces (const G4Mesh &)
 
void DrawTetMeshAsDots (const G4Mesh &)
 
void DrawTetMeshAsSurfaces (const G4Mesh &)
 
G4ThreeVector GetPointInBox (const G4ThreeVector &pos, G4double halfX, G4double halfY, G4double halfZ) const
 
G4ThreeVector GetPointInTet (const std::vector< G4ThreeVector > &vertices) const
 

Protected Attributes

G4VGraphicsSystemfSystem
 
const G4int fSceneHandlerId
 
G4String fName
 
G4int fViewCount
 
G4ViewerList fViewerList
 
G4VViewerfpViewer
 
G4ScenefpScene
 
G4bool fMarkForClearingTransientStore
 
G4bool fReadyForTransients
 
G4bool fTransientsDrawnThisEvent
 
G4bool fTransientsDrawnThisRun
 
G4bool fProcessingSolid
 
G4bool fProcessing2D
 
G4VModelfpModel
 
G4Transform3D fObjectTransformation
 
G4int fNestingDepth
 
const G4VisAttributesfpVisAttribs
 
const G4Transform3D fIdentityTransformation
 

Friends

class G4VViewer
 
std::ostream & operator<< (std::ostream &os, const G4VSceneHandler &s)
 

Detailed Description

Definition at line 52 of file G4VSceneHandler.hh.

Member Enumeration Documentation

◆ MarkerSizeType

Enumerator
world 
screen 

Definition at line 59 of file G4VSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4VSceneHandler()

G4VSceneHandler::G4VSceneHandler ( G4VGraphicsSystem system,
G4int  id,
const G4String name = "" 
)

Definition at line 99 of file G4VSceneHandler.cc.

99 :
100 fSystem (system),
101 fSceneHandlerId (id),
102 fViewCount (0),
103 fpViewer (0),
104 fpScene (0),
105 fMarkForClearingTransientStore (true), // Ready for first
106 // ClearTransientStoreIfMarked(),
107 // e.g., at end of run (see
108 // G4VisManager.cc).
109 fReadyForTransients (true), // Only false while processing scene.
110 fProcessingSolid (false),
111 fProcessing2D (false),
112 fpModel (0),
113 fNestingDepth (0),
114 fpVisAttribs (0)
115{
117 fpScene = pVMan -> GetCurrentScene ();
118 if (name == "") {
119 std::ostringstream ost;
120 ost << fSystem.GetName () << '-' << fSceneHandlerId;
121 fName = ost.str();
122 }
123 else {
124 fName = name;
125 }
128}
const G4String & GetName() const
G4bool fTransientsDrawnThisEvent
const G4int fSceneHandlerId
G4bool fTransientsDrawnThisRun
G4VViewer * fpViewer
G4bool fMarkForClearingTransientStore
const G4VisAttributes * fpVisAttribs
G4VGraphicsSystem & fSystem
G4bool GetTransientsDrawnThisEvent() const
G4bool GetTransientsDrawnThisRun() const
static G4VisManager * GetInstance()
const char * name(G4int ptype)

◆ ~G4VSceneHandler()

G4VSceneHandler::~G4VSceneHandler ( )
virtual

Definition at line 130 of file G4VSceneHandler.cc.

130 {
131 G4VViewer* last;
132 while( ! fViewerList.empty() ) {
133 last = fViewerList.back();
134 fViewerList.pop_back();
135 delete last;
136 }
137}
G4ViewerList fViewerList

Member Function Documentation

◆ AddCompound() [1/6]

void G4VSceneHandler::AddCompound ( const G4Mesh mesh)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4VtkSceneHandler.

Definition at line 431 of file G4VSceneHandler.cc.

432{
433 G4warn <<
434 "There has been an attempt to draw a mesh with option \""
436 << "\":\n" << mesh
437 << "but it is not of a recognised type or is not implemented"
438 "\nby the current graphics driver. Instead we draw its"
439 "\ncontainer \"" << mesh.GetContainerVolume()->GetName() << "\"."
440 << G4endl;
441 const auto& pv = mesh.GetContainerVolume();
442 const auto& lv = pv->GetLogicalVolume();
443 const auto& solid = lv->GetSolid();
444 const auto& transform = mesh.GetTransform();
445 // Make sure container is visible
446 G4VisAttributes tmpVisAtts; // Visible, white, not forced.
447 const auto& saveVisAtts = lv->GetVisAttributes();
448 if (saveVisAtts) {
449 tmpVisAtts = *saveVisAtts;
450 tmpVisAtts.SetVisibility(true);
451 auto colour = saveVisAtts->GetColour();
452 colour.SetAlpha(1.);
453 tmpVisAtts.SetColour(colour);
454 }
455 // Draw container
456 PreAddSolid(transform,tmpVisAtts);
457 solid->DescribeYourselfTo(*this);
458 PostAddSolid();
459 // Restore vis attributes
460 lv->SetVisAttributes(saveVisAtts);
461}
#define G4warn
Definition: G4Scene.cc:41
#define G4endl
Definition: G4ios.hh:57
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetContainerVolume() const
Definition: G4Mesh.hh:73
const G4Transform3D & GetTransform() const
Definition: G4Mesh.hh:77
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
virtual void PostAddSolid()
const G4ViewParameters & GetViewParameters() const
SMROption GetSpecialMeshRenderingOption() const
void SetColour(const G4Colour &)
void SetVisibility(G4bool=true)

◆ AddCompound() [2/6]

void G4VSceneHandler::AddCompound ( const G4THitsMap< G4double > &  hits)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4GMocrenFileSceneHandler.

Definition at line 345 of file G4VSceneHandler.cc.

345 {
346 using MeshScoreMap = G4VScoringMesh::MeshScoreMap;
347 //G4cout << "AddCompound: hits: " << &hits << G4endl;
348 G4bool scoreMapHits = false;
350 if (scoringManager) {
351 std::size_t nMeshes = scoringManager->GetNumberOfMesh();
352 for (std::size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
353 G4VScoringMesh* mesh = scoringManager->GetMesh((G4int)iMesh);
354 if (mesh && mesh->IsActive()) {
355 MeshScoreMap scoreMap = mesh->GetScoreMap();
356 const G4String& mapNam = const_cast<G4THitsMap<G4double>&>(hits).GetName();
357 for(MeshScoreMap::const_iterator i = scoreMap.cbegin();
358 i != scoreMap.cend(); ++i) {
359 const G4String& scoreMapName = i->first;
360 if (scoreMapName == mapNam) {
361 G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
362 scoreMapHits = true;
363 mesh->DrawMesh(scoreMapName, &colorMap);
364 }
365 }
366 }
367 }
368 }
369 if (scoreMapHits) {
370 static G4bool first = true;
371 if (first) {
372 first = false;
373 G4cout <<
374 "Scoring map drawn with default parameters."
375 "\n To get gMocren file for gMocren browser:"
376 "\n /vis/open gMocrenFile"
377 "\n /vis/viewer/flush"
378 "\n Many other options available with /score/draw... commands."
379 "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
380 << G4endl;
381 }
382 } else { // Not score map hits. Just call DrawAllHits.
383 // Cast away const because DrawAllHits is non-const!!!!
384 const_cast<G4THitsMap<G4double>&>(hits).DrawAllHits();
385 }
386}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4VScoringMesh * GetMesh(G4int i) const
size_t GetNumberOfMesh() const
static G4ScoringManager * GetScoringManagerIfExist()
const G4String & GetName() const
G4bool IsActive() const
std::map< G4String, RunScore * > MeshScoreMap
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
MeshScoreMap GetScoreMap() const

◆ AddCompound() [3/6]

void G4VSceneHandler::AddCompound ( const G4THitsMap< G4StatDouble > &  hits)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4GMocrenFileSceneHandler.

Definition at line 388 of file G4VSceneHandler.cc.

388 {
389 using MeshScoreMap = G4VScoringMesh::MeshScoreMap;
390 //G4cout << "AddCompound: hits: " << &hits << G4endl;
391 G4bool scoreMapHits = false;
393 if (scoringManager) {
394 std::size_t nMeshes = scoringManager->GetNumberOfMesh();
395 for (std::size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
396 G4VScoringMesh* mesh = scoringManager->GetMesh((G4int)iMesh);
397 if (mesh && mesh->IsActive()) {
398 MeshScoreMap scoreMap = mesh->GetScoreMap();
399 for(MeshScoreMap::const_iterator i = scoreMap.cbegin();
400 i != scoreMap.cend(); ++i) {
401 const G4String& scoreMapName = i->first;
402 const G4THitsMap<G4StatDouble>* foundHits = i->second;
403 if (foundHits == &hits) {
404 G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
405 scoreMapHits = true;
406 mesh->DrawMesh(scoreMapName, &colorMap);
407 }
408 }
409 }
410 }
411 }
412 if (scoreMapHits) {
413 static G4bool first = true;
414 if (first) {
415 first = false;
416 G4cout <<
417 "Scoring map drawn with default parameters."
418 "\n To get gMocren file for gMocren browser:"
419 "\n /vis/open gMocrenFile"
420 "\n /vis/viewer/flush"
421 "\n Many other options available with /score/draw... commands."
422 "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
423 << G4endl;
424 }
425 } else { // Not score map hits. Just call DrawAllHits.
426 // Cast away const because DrawAllHits is non-const!!!!
427 const_cast<G4THitsMap<G4StatDouble>&>(hits).DrawAllHits();
428 }
429}

◆ AddCompound() [4/6]

void G4VSceneHandler::AddCompound ( const G4VDigi digi)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4GMocrenFileSceneHandler.

Definition at line 340 of file G4VSceneHandler.cc.

340 {
341 // Cast away const because Draw is non-const!!!!
342 const_cast<G4VDigi&>(digi).Draw();
343}

◆ AddCompound() [5/6]

void G4VSceneHandler::AddCompound ( const G4VHit hit)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4GMocrenFileSceneHandler.

Definition at line 335 of file G4VSceneHandler.cc.

335 {
336 // Cast away const because Draw is non-const!!!!
337 const_cast<G4VHit&>(hit).Draw();
338}
Definition: G4VHit.hh:48

◆ AddCompound() [6/6]

void G4VSceneHandler::AddCompound ( const G4VTrajectory traj)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4ToolsSGSceneHandler, G4VRML2FileSceneHandler, and G4GMocrenFileSceneHandler.

Definition at line 323 of file G4VSceneHandler.cc.

323 {
324 G4TrajectoriesModel* trajectoriesModel =
325 dynamic_cast<G4TrajectoriesModel*>(fpModel);
326 if (trajectoriesModel)
327 traj.DrawTrajectory();
328 else {
330 ("G4VSceneHandler::AddCompound(const G4VTrajectory&)",
331 "visman0105", FatalException, "Not a G4TrajectoriesModel.");
332 }
333}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
virtual void DrawTrajectory() const

Referenced by G4VtkSceneHandler::AddCompound(), G4HepRepFileSceneHandler::AddCompound(), G4OpenGLSceneHandler::AddCompound(), G4GMocrenFileSceneHandler::AddCompound(), and StandardSpecialMeshRendering().

◆ AddPrimitive() [1/7]

◆ AddPrimitive() [2/7]

void G4VSceneHandler::AddPrimitive ( const G4Plotter )
virtual

Implements G4VGraphicsScene.

Reimplemented in G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLSceneHandler, G4OpenGLStoredSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4RayTracerSceneHandler, G4ToolsSGSceneHandler, G4VTreeSceneHandler, and G4VRML2FileSceneHandler.

Definition at line 510 of file G4VSceneHandler.cc.

510 {
511 G4warn << "WARNING: Plotter not implemented for " << fSystem.GetName() << G4endl;
512 G4warn << " Open a plotter-aware graphics system or remove plotter with" << G4endl;
513 G4warn << " /vis/scene/removeModel Plotter" << G4endl;
514}

◆ AddPrimitive() [3/7]

◆ AddPrimitive() [4/7]

◆ AddPrimitive() [5/7]

void G4VSceneHandler::AddPrimitive ( const G4Polymarker polymarker)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLSceneHandler, G4OpenGLSceneHandler, G4OpenGLStoredSceneHandler, G4OpenGLStoredSceneHandler, G4OpenInventorSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, G4Qt3DSceneHandler, G4RayTracerSceneHandler, G4RayTracerSceneHandler, G4ToolsSGSceneHandler, G4VTreeSceneHandler, G4VTreeSceneHandler, G4VRML2FileSceneHandler, and G4VtkSceneHandler.

Definition at line 467 of file G4VSceneHandler.cc.

467 {
468 switch (polymarker.GetMarkerType()) {
469 default:
471 {
472 G4Circle dot (polymarker);
473 dot.SetWorldSize (0.);
474 dot.SetScreenSize (0.1); // Very small circle.
475 for (std::size_t iPoint = 0; iPoint < polymarker.size (); ++iPoint) {
476 dot.SetPosition (polymarker[iPoint]);
477 AddPrimitive (dot);
478 }
479 }
480 break;
482 {
483 G4Circle circle (polymarker); // Default circle
484 for (std::size_t iPoint = 0; iPoint < polymarker.size (); ++iPoint) {
485 circle.SetPosition (polymarker[iPoint]);
486 AddPrimitive (circle);
487 }
488 }
489 break;
491 {
492 G4Square square (polymarker); // Default square
493 for (std::size_t iPoint = 0; iPoint < polymarker.size (); ++iPoint) {
494 square.SetPosition (polymarker[iPoint]);
495 AddPrimitive (square);
496 }
497 }
498 break;
499 }
500}
MarkerType GetMarkerType() const
virtual void AddPrimitive(const G4Polyline &)=0

◆ AddPrimitive() [6/7]

◆ AddPrimitive() [7/7]

◆ AddSolid() [1/14]

void G4VSceneHandler::AddSolid ( const G4Box box)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, G4VRML2FileSceneHandler, G4VRML2FileSceneHandler, G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, and G4VtkSceneHandler.

Definition at line 252 of file G4VSceneHandler.cc.

252 {
253 AddSolidT (box);
254 // If your graphics system is sophisticated enough to handle a
255 // particular solid shape as a primitive, in your derived class write a
256 // function to override this.
257 // Your function might look like this...
258 // void G4MySceneHandler::AddSolid (const G4Box& box) {
259 // Get and check applicable vis attributes.
260 // fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
261 // Do not draw if not visible.
262 // if (fpVisAttribs->IsVisible()) {
263 // Get parameters of appropriate object, e.g.:
264 // G4double dx = box.GetXHalfLength ();
265 // G4double dy = box.GetYHalfLength ();
266 // G4double dz = box.GetZHalfLength ();
267 // ...
268 // and Draw or Store in your display List.
269}
void AddSolidT(const T &solid)

Referenced by G4VtkSceneHandler::AddSolid(), G4HepRepFileSceneHandler::AddSolid(), and G4GMocrenFileSceneHandler::AddSolid().

◆ AddSolid() [2/14]

void G4VSceneHandler::AddSolid ( const G4Cons cons)
virtual

◆ AddSolid() [3/14]

void G4VSceneHandler::AddSolid ( const G4Ellipsoid ellipsoid)
virtual

Implements G4VGraphicsScene.

Reimplemented in G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, G4HepRepFileSceneHandler, G4OpenGLSceneHandler, and G4VRML2FileSceneHandler.

Definition at line 303 of file G4VSceneHandler.cc.

303 {
304 AddSolidWithAuxiliaryEdges (ellipsoid);
305}
void AddSolidWithAuxiliaryEdges(const T &solid)

◆ AddSolid() [4/14]

void G4VSceneHandler::AddSolid ( const G4Orb orb)
virtual

◆ AddSolid() [5/14]

void G4VSceneHandler::AddSolid ( const G4Para para)
virtual

◆ AddSolid() [6/14]

void G4VSceneHandler::AddSolid ( const G4Polycone polycone)
virtual

◆ AddSolid() [7/14]

void G4VSceneHandler::AddSolid ( const G4Polyhedra polyhedra)
virtual

◆ AddSolid() [8/14]

◆ AddSolid() [9/14]

void G4VSceneHandler::AddSolid ( const G4TessellatedSolid tess)
virtual

◆ AddSolid() [10/14]

◆ AddSolid() [11/14]

void G4VSceneHandler::AddSolid ( const G4Trap trap)
virtual

◆ AddSolid() [12/14]

void G4VSceneHandler::AddSolid ( const G4Trd trd)
virtual

◆ AddSolid() [13/14]

void G4VSceneHandler::AddSolid ( const G4Tubs tubs)
virtual

◆ AddSolid() [14/14]

void G4VSceneHandler::AddSolid ( const G4VSolid solid)
virtual

◆ AddSolidT()

template<class T >
void G4VSceneHandler::AddSolidT ( const T &  solid)

Definition at line 225 of file G4VSceneHandler.cc.

227{
228 // Get and check applicable vis attributes.
230 RequestPrimitives (solid);
231}
virtual void RequestPrimitives(const G4VSolid &solid)
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const

Referenced by AddSolid().

◆ AddSolidWithAuxiliaryEdges()

template<class T >
void G4VSceneHandler::AddSolidWithAuxiliaryEdges ( const T &  solid)

Definition at line 233 of file G4VSceneHandler.cc.

235{
236 // Get and check applicable vis attributes.
238 // Draw with auxiliary edges unless otherwise specified.
240 // Create a vis atts object for the modified vis atts.
241 // It is static so that we may return a reliable pointer to it.
242 static G4VisAttributes visAttsWithAuxEdges;
243 // Initialise it with the current vis atts and reset the pointer.
244 visAttsWithAuxEdges = *fpVisAttribs;
245 // Force auxiliary edges visible.
246 visAttsWithAuxEdges.SetForceAuxEdgeVisible();
247 fpVisAttribs = &visAttsWithAuxEdges;
248 }
249 RequestPrimitives (solid);
250}
void SetForceAuxEdgeVisible(G4bool=true)
G4bool IsForceAuxEdgeVisible() const

Referenced by AddSolid().

◆ AddViewerToList()

void G4VSceneHandler::AddViewerToList ( G4VViewer pView)

Definition at line 463 of file G4VSceneHandler.cc.

463 {
464 fViewerList.push_back (pViewer);
465}

◆ BeginModeling()

◆ BeginPrimitives()

void G4VSceneHandler::BeginPrimitives ( const G4Transform3D objectTransformation = G4Transform3D())
virtual

Implements G4VGraphicsScene.

Reimplemented in G4DAWNFILESceneHandler, G4GMocrenFileSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLSceneHandler, G4OpenGLStoredSceneHandler, G4OpenInventorSceneHandler, G4Qt3DSceneHandler, and G4VRML2FileSceneHandler.

Definition at line 165 of file G4VSceneHandler.cc.

166 {
167 //static G4int count = 0;
168 //G4cout << "G4VSceneHandler::BeginPrimitives: " << count++ << G4endl;
170 if (fNestingDepth > 1)
172 ("G4VSceneHandler::BeginPrimitives",
173 "visman0101", FatalException,
174 "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
175 fObjectTransformation = objectTransformation;
176}
G4Transform3D fObjectTransformation

Referenced by G4GMocrenFileSceneHandler::BeginPrimitives(), G4OpenGLSceneHandler::BeginPrimitives(), G4OpenInventorSceneHandler::BeginPrimitives(), G4Qt3DSceneHandler::BeginPrimitives(), Draw3DRectMeshAsDots(), Draw3DRectMeshAsSurfaces(), DrawTetMeshAsDots(), DrawTetMeshAsSurfaces(), RequestPrimitives(), and StandardSpecialMeshRendering().

◆ BeginPrimitives2D()

void G4VSceneHandler::BeginPrimitives2D ( const G4Transform3D objectTransformation = G4Transform3D())
virtual

Implements G4VGraphicsScene.

Reimplemented in G4HepRepFileSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLSceneHandler, G4OpenGLStoredSceneHandler, and G4Qt3DSceneHandler.

Definition at line 189 of file G4VSceneHandler.cc.

190 {
192 if (fNestingDepth > 1)
194 ("G4VSceneHandler::BeginPrimitives2D",
195 "visman0103", FatalException,
196 "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
197 fObjectTransformation = objectTransformation;
198 fProcessing2D = true;
199}

Referenced by G4HepRepFileSceneHandler::BeginPrimitives2D(), G4OpenGLSceneHandler::BeginPrimitives2D(), and G4Qt3DSceneHandler::BeginPrimitives2D().

◆ ClearStore()

◆ ClearTransientStore()

◆ CreateCutawaySolid()

G4DisplacedSolid * G4VSceneHandler::CreateCutawaySolid ( )
protectedvirtual

Reimplemented in G4OpenGLSceneHandler.

Definition at line 938 of file G4VSceneHandler.cc.

939{
941 if (vp.IsCutaway()) {
942
943 std::vector<G4DisplacedSolid*> cutaway_solids;
944
946 G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
947 G4VSolid* cutawayBox =
948 new G4Box("_cutaway_box", safe, safe, safe); // world box...
949
950 for (int plane_no = 0; plane_no < int(vp.GetCutawayPlanes().size()); plane_no++){
951
952 const G4Normal3D originalNormal(0,0,1); // ...so this is original normal.
953
954 const G4Plane3D& sp = vp.GetCutawayPlanes()[plane_no]; //];
955 const G4double& a = sp.a();
956 const G4double& b = sp.b();
957 const G4double& c = sp.c();
958 const G4double& d = sp.d();
959 const G4Normal3D newNormal(-a,-b,-c); // Convention: keep a*x+b*y+c*z+d>=0
960 // Not easy to see why the above gives the right convention, but it has been
961 // arrived at by trial and error to agree with the OpenGL implementation
962 // of clipping planes.
963
964 G4Transform3D requiredTransform; // Null transform
965 // Calculate the rotation
966 // If newNormal is (0,0,1), no need to do anything
967 // Treat (0,0,-1) as a special case, since cannot define axis in this case
968 if (newNormal == G4Normal3D(0,0,-1)) {
969 requiredTransform = G4Rotate3D(pi,G4Vector3D(1,0,0));
970 } else if (newNormal != originalNormal) {
971 const G4double& angle = std::acos(newNormal.dot(originalNormal));
972 const G4Vector3D& axis = originalNormal.cross(newNormal);
973 requiredTransform = G4Rotate3D(angle, axis);
974 }
975 // Translation
976 requiredTransform = requiredTransform * G4TranslateZ3D(d + safe);
977 cutaway_solids.push_back
978 (new G4DisplacedSolid("_displaced_cutaway_box", cutawayBox, requiredTransform));
979 }
980
981 if (cutaway_solids.size() == 1){
982 return (G4DisplacedSolid*) cutaway_solids[0];
984 G4UnionSolid* union2 =
985 new G4UnionSolid("_union_2", cutaway_solids[0], cutaway_solids[1]);
986 if (cutaway_solids.size() == 2)
987 return (G4DisplacedSolid*)union2;
988 else
989 return (G4DisplacedSolid*)
990 new G4UnionSolid("_union_3", union2, cutaway_solids[2]);
992 G4IntersectionSolid* intersection2 =
993 new G4IntersectionSolid("_intersection_2", cutaway_solids[0], cutaway_solids[1]);
994 if (cutaway_solids.size() == 2)
995 return (G4DisplacedSolid*)intersection2;
996 else
997 return (G4DisplacedSolid*)
998 new G4IntersectionSolid("_intersection_3", intersection2, cutaway_solids[2]);
999 }
1000 }
1001
1002 return 0;
1003}
HepGeom::Normal3D< G4double > G4Normal3D
Definition: G4Normal3D.hh:34
HepGeom::TranslateZ3D G4TranslateZ3D
HepGeom::Rotate3D G4Rotate3D
double G4double
Definition: G4Types.hh:83
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
Definition: G4Box.hh:56
const G4VisExtent & GetExtent() const
CutawayMode GetCutawayMode() const
G4bool IsCutaway() const
const G4Planes & GetCutawayPlanes() const
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
const G4Point3D & GetExtentCentre() const
Definition: G4VisExtent.cc:65
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

Referenced by CreateModelingParameters().

◆ CreateModelingParameters()

G4ModelingParameters * G4VSceneHandler::CreateModelingParameters ( )

Definition at line 832 of file G4VSceneHandler.cc.

833{
834 // Create modeling parameters from View Parameters...
835 if (!fpViewer) return NULL;
836
837 const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
838
839 // Convert drawing styles...
840 G4ModelingParameters::DrawingStyle modelDrawingStyle =
842 switch (vp.GetDrawingStyle ()) {
843 default:
845 modelDrawingStyle = G4ModelingParameters::wf;
846 break;
848 modelDrawingStyle = G4ModelingParameters::hlr;
849 break;
851 modelDrawingStyle = G4ModelingParameters::hsr;
852 break;
854 modelDrawingStyle = G4ModelingParameters::hlhsr;
855 break;
857 modelDrawingStyle = G4ModelingParameters::cloud;
858 break;
859 }
860
861 // Decide if covered daughters are really to be culled...
862 G4bool reallyCullCovered =
863 vp.IsCullingCovered() // Culling daughters depends also on...
864 && !vp.IsSection () // Sections (DCUT) not requested.
865 && !vp.IsCutaway () // Cutaways not requested.
866 ;
867
868 G4ModelingParameters* pModelingParams = new G4ModelingParameters
870 modelDrawingStyle,
871 vp.IsCulling (),
872 vp.IsCullingInvisible (),
873 vp.IsDensityCulling (),
874 vp.GetVisibleDensity (),
875 reallyCullCovered,
876 vp.GetNoOfSides ()
877 );
878
879 pModelingParams->SetNumberOfCloudPoints(vp.GetNumberOfCloudPoints());
880 pModelingParams->SetWarning
882
883 pModelingParams->SetCBDAlgorithmNumber(vp.GetCBDAlgorithmNumber());
884 pModelingParams->SetCBDParameters(vp.GetCBDParameters());
885
886 pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
887 pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
888
889 pModelingParams->SetSectionSolid(CreateSectionSolid());
890 pModelingParams->SetCutawaySolid(CreateCutawaySolid());
891 // The polyhedron objects are deleted in the modeling parameters destructor.
892
894
895 pModelingParams->SetSpecialMeshRendering(vp.IsSpecialMeshRendering());
896 pModelingParams->SetSpecialMeshVolumes(vp.GetSpecialMeshVolumes());
897
898 return pModelingParams;
899}
void SetCBDParameters(const std::vector< G4double > &)
void SetWarning(G4bool)
void SetNumberOfCloudPoints(G4int)
void SetCBDAlgorithmNumber(G4int)
void SetExplodeFactor(G4double explodeFactor)
void SetVisAttributesModifiers(const std::vector< VisAttributesModifier > &)
void SetExplodeCentre(const G4Point3D &explodeCentre)
void SetCutawaySolid(G4DisplacedSolid *pCutawaySolid)
void SetSectionSolid(G4DisplacedSolid *pSectionSolid)
void SetSpecialMeshVolumes(const std::vector< PVNameCopyNo > &)
void SetSpecialMeshRendering(G4bool)
virtual G4DisplacedSolid * CreateSectionSolid()
virtual G4DisplacedSolid * CreateCutawaySolid()
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsSection() const
G4bool IsCulling() const
const std::vector< G4double > & GetCBDParameters() const
G4int GetCBDAlgorithmNumber() const
const std::vector< G4ModelingParameters::PVNameCopyNo > & GetSpecialMeshVolumes() const
G4bool IsCullingInvisible() const
const G4VisAttributes * GetDefaultVisAttributes() const
G4bool IsDensityCulling() const
G4double GetVisibleDensity() const
const G4Point3D & GetExplodeCentre() const
G4bool IsCullingCovered() const
DrawingStyle GetDrawingStyle() const
static Verbosity GetVerbosity()

Referenced by DrawEndOfRunModels(), DrawEvent(), G4VisManager::DrawGeometry(), and ProcessScene().

◆ CreateSectionSolid()

G4DisplacedSolid * G4VSceneHandler::CreateSectionSolid ( )
protectedvirtual

Reimplemented in G4OpenGLSceneHandler.

Definition at line 901 of file G4VSceneHandler.cc.

902{
903 G4DisplacedSolid* sectioner = 0;
904
906 if (vp.IsSection () ) {
907
909 G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
910 G4VSolid* sectionBox =
911 new G4Box("_sectioner", safe, safe, 1.e-5 * radius); // Thin in z-plane...
912 const G4Normal3D originalNormal(0,0,1); // ...so this is original normal.
913
914 const G4Plane3D& sp = vp.GetSectionPlane ();
915 const G4double& a = sp.a();
916 const G4double& b = sp.b();
917 const G4double& c = sp.c();
918 const G4double& d = sp.d();
919 const G4Normal3D newNormal(a,b,c);
920
921 G4Transform3D requiredTransform;
922 // Rotate
923 if (newNormal != originalNormal) {
924 const G4double& angle = std::acos(newNormal.dot(originalNormal));
925 const G4Vector3D& axis = originalNormal.cross(newNormal);
926 requiredTransform = G4Rotate3D(angle, axis);
927 }
928 // Translate
929 requiredTransform = requiredTransform * G4TranslateZ3D(-d);
930
931 sectioner = new G4DisplacedSolid
932 ("_displaced_sectioning_box", sectionBox, requiredTransform);
933 }
934
935 return sectioner;
936}
const G4Plane3D & GetSectionPlane() const

Referenced by CreateModelingParameters(), and G4OpenGLSceneHandler::CreateSectionSolid().

◆ Draw3DRectMeshAsDots()

void G4VSceneHandler::Draw3DRectMeshAsDots ( const G4Mesh mesh)
protected

Definition at line 1331 of file G4VSceneHandler.cc.

1334{
1335 // Check
1336 if (mesh.GetMeshType() != G4Mesh::rectangle &&
1339 ed << "Called with a mesh that is not rectangular:" << mesh;
1340 G4Exception("G4VSceneHandler::Draw3DRectMeshAsDots","visman0108",JustWarning,ed);
1341 return;
1342 }
1343
1344 static G4bool firstPrint = true;
1345 const auto& verbosity = G4VisManager::GetVerbosity();
1346 G4bool print = firstPrint && verbosity >= G4VisManager::errors;
1347 if (print) {
1348 G4cout
1349 << "Special case drawing of 3D rectangular G4VNestedParameterisation as dots:"
1350 << '\n' << mesh
1351 << G4endl;
1352 }
1353
1354 const auto& container = mesh.GetContainerVolume();
1355
1356 // This map is static so that once filled it stays filled.
1357 static std::map<G4String,std::map<const G4Material*,G4Polymarker>> dotsByMaterialAndMesh;
1358 auto& dotsByMaterial = dotsByMaterialAndMesh[mesh.GetContainerVolume()->GetName()];
1359
1360 // Fill map if not already filled
1361 if (dotsByMaterial.empty()) {
1362
1363 // Get positions and material one cell at a time (using PseudoSceneFor3DRectMeshPositions).
1364 // The pseudo scene allows a "private" descent into the parameterisation.
1365 // Instantiate a temporary G4PhysicalVolumeModel
1367 tmpMP.SetCulling(true); // This avoids drawing transparent...
1368 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
1369 const G4bool useFullExtent = true; // To avoid calculating the extent
1370 G4PhysicalVolumeModel tmpPVModel
1371 (container,
1373 G4Transform3D(), // so that positions are in local coordinates
1374 &tmpMP,
1375 useFullExtent);
1376 // Accumulate information in temporary maps by material
1377 std::multimap<const G4Material*,const G4ThreeVector> positionByMaterial;
1378 std::map<const G4Material*,G4VSceneHandler::NameAndVisAtts> nameAndVisAttsByMaterial;
1379 // Instantiate the pseudo scene
1380 PseudoSceneFor3DRectMeshPositions pseudoScene
1381 (&tmpPVModel,mesh.GetMeshDepth(),positionByMaterial,nameAndVisAttsByMaterial);
1382 // Make private descent into the parameterisation
1383 tmpPVModel.DescribeYourselfTo(pseudoScene);
1384 // Now we have a map of positions by material.
1385 // Also a map of name and colour by material.
1386
1387 const auto& prms = mesh.GetThreeDRectParameters();
1388 const auto& halfX = prms.fHalfX;
1389 const auto& halfY = prms.fHalfY;
1390 const auto& halfZ = prms.fHalfZ;
1391
1392 // Fill the permanent (static) map of dots by material
1393 G4int nDotsTotal = 0;
1394 for (const auto& entry: nameAndVisAttsByMaterial) {
1395 G4int nDots = 0;
1396 const auto& material = entry.first;
1397 const auto& nameAndVisAtts = nameAndVisAttsByMaterial[material];
1398 const auto& name = nameAndVisAtts.fName;
1399 const auto& visAtts = nameAndVisAtts.fVisAtts;
1400 G4Polymarker dots;
1401 dots.SetInfo(name);
1402 dots.SetVisAttributes(visAtts);
1404 dots.SetSize(G4VMarker::screen,1.);
1405 // Enter empty polymarker into the map
1406 dotsByMaterial[material] = dots;
1407 // Now fill it in situ
1408 auto& dotsInMap = dotsByMaterial[material];
1409 const auto& range = positionByMaterial.equal_range(material);
1410 for (auto posByMat = range.first; posByMat != range.second; ++posByMat) {
1411 dotsInMap.push_back(GetPointInBox(posByMat->second, halfX, halfY, halfZ));
1412 ++nDots;
1413 }
1414
1415 if (print) {
1416 G4cout
1417 << std::setw(30) << std::left << name.substr(0,30) << std::right
1418 << ": " << std::setw(7) << nDots << " dots"
1419 << ": colour " << std::fixed << std::setprecision(2)
1420 << visAtts.GetColour() << std::defaultfloat
1421 << G4endl;
1422 }
1423
1424 nDotsTotal += nDots;
1425 }
1426
1427 if (print) {
1428 G4cout << "Total number of dots: " << nDotsTotal << G4endl;
1429 }
1430 }
1431
1432 // Some subsequent expressions apply only to G4PhysicalVolumeModel
1433 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1434
1435 G4String parameterisationName;
1436 if (pPVModel) {
1437 parameterisationName = pPVModel->GetFullPVPath().back().GetPhysicalVolume()->GetName();
1438 }
1439
1440 // Draw the dots by material
1441 // Ensure they are "hidden", i.e., use the z-buffer as non-marker primitives do
1442 auto keepVP = fpViewer->GetViewParameters();
1443 auto vp = fpViewer->GetViewParameters();
1444 vp.SetMarkerHidden();
1446 // Now we transform to world coordinates
1448 for (const auto& entry: dotsByMaterial) {
1449 const auto& dots = entry.second;
1450 // The current "leaf" node in the PVPath is the parameterisation. Here it has
1451 // been converted into polymarkers by material. So...temporarily...change
1452 // its name to that of the material (whose name has been stored in Info)
1453 // so that its appearance in the scene tree of, e.g., G4OpenGLQtViewer, has
1454 // an appropriate name and its visibility and colour may be changed.
1455 if (pPVModel) {
1456 const auto& fullPVPath = pPVModel->GetFullPVPath();
1457 auto leafPV = fullPVPath.back().GetPhysicalVolume();
1458 leafPV->SetName(dots.GetInfo());
1459 }
1460 // Add dots to the scene
1461 AddPrimitive(dots);
1462 }
1463 EndPrimitives ();
1464 // Restore view parameters
1465 fpViewer->SetViewParameters(keepVP);
1466 // Restore parameterisation name
1467 if (pPVModel) {
1468 pPVModel->GetFullPVPath().back().GetPhysicalVolume()->SetName(parameterisationName);
1469 }
1470
1471 firstPrint = false;
1472 return;
1473}
@ JustWarning
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
HepGeom::Transform3D G4Transform3D
void print(G4double elem)
MeshType GetMeshType() const
Definition: G4Mesh.hh:75
@ rectangle
Definition: G4Mesh.hh:53
@ nested3DRectangular
Definition: G4Mesh.hh:54
const ThreeDRectangleParameters & GetThreeDRectParameters() const
Definition: G4Mesh.hh:78
G4int GetMeshDepth() const
Definition: G4Mesh.hh:76
void SetCulling(G4bool)
void SetCullingInvisible(G4bool)
void SetMarkerType(MarkerType)
void SetSize(SizeType, G4double)
Definition: G4VMarker.cc:86
virtual void EndPrimitives()
G4ThreeVector GetPointInBox(const G4ThreeVector &pos, G4double halfX, G4double halfY, G4double halfZ) const
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
void SetViewParameters(const G4ViewParameters &vp)
Definition: G4VViewer.cc:126
void SetMarkerHidden()
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:98
virtual void SetInfo(const G4String &info)
virtual const G4String & GetInfo() const

Referenced by StandardSpecialMeshRendering().

◆ Draw3DRectMeshAsSurfaces()

void G4VSceneHandler::Draw3DRectMeshAsSurfaces ( const G4Mesh mesh)
protected

Definition at line 1475 of file G4VSceneHandler.cc.

1478{
1479 // Check
1480 if (mesh.GetMeshType() != G4Mesh::rectangle &&
1483 ed << "Called with a mesh that is not rectangular:" << mesh;
1484 G4Exception("G4VSceneHandler::Draw3DRectMeshAsSurfaces","visman0108",JustWarning,ed);
1485 return;
1486 }
1487
1488 static G4bool firstPrint = true;
1489 const auto& verbosity = G4VisManager::GetVerbosity();
1490 G4bool print = firstPrint && verbosity >= G4VisManager::errors;
1491 if (print) {
1492 G4cout
1493 << "Special case drawing of 3D rectangular G4VNestedParameterisation as surfaces:"
1494 << '\n' << mesh
1495 << G4endl;
1496 }
1497
1498 const auto& container = mesh.GetContainerVolume();
1499
1500 // This map is static so that once filled it stays filled.
1501 static std::map<G4String,std::map<const G4Material*,G4Polyhedron>> boxesByMaterialAndMesh;
1502 auto& boxesByMaterial = boxesByMaterialAndMesh[mesh.GetContainerVolume()->GetName()];
1503
1504 // Fill map if not already filled
1505 if (boxesByMaterial.empty()) {
1506
1507 // Get positions and material one cell at a time (using PseudoSceneFor3DRectMeshPositions).
1508 // The pseudo scene allows a "private" descent into the parameterisation.
1509 // Instantiate a temporary G4PhysicalVolumeModel
1511 tmpMP.SetCulling(true); // This avoids drawing transparent...
1512 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
1513 const G4bool useFullExtent = true; // To avoid calculating the extent
1514 G4PhysicalVolumeModel tmpPVModel
1515 (container,
1517 G4Transform3D(), // so that positions are in local coordinates
1518 &tmpMP,
1519 useFullExtent);
1520 // Accumulate information in temporary maps by material
1521 std::multimap<const G4Material*,const G4ThreeVector> positionByMaterial;
1522 std::map<const G4Material*,G4VSceneHandler::NameAndVisAtts> nameAndVisAttsByMaterial;
1523 // Instantiate the pseudo scene
1524 PseudoSceneFor3DRectMeshPositions pseudoScene
1525 (&tmpPVModel,mesh.GetMeshDepth(),positionByMaterial,nameAndVisAttsByMaterial);
1526 // Make private descent into the parameterisation
1527 tmpPVModel.DescribeYourselfTo(pseudoScene);
1528 // Now we have a map of positions by material.
1529 // Also a map of name and colour by material.
1530
1531 const auto& prms = mesh.GetThreeDRectParameters();
1532 const auto& sizeX = 2.*prms.fHalfX;
1533 const auto& sizeY = 2.*prms.fHalfY;
1534 const auto& sizeZ = 2.*prms.fHalfZ;
1535
1536 // Fill the permanent (static) map of boxes by material
1537 G4int nBoxesTotal = 0, nFacetsTotal = 0;
1538 for (const auto& entry: nameAndVisAttsByMaterial) {
1539 G4int nBoxes = 0;
1540 const auto& material = entry.first;
1541 const auto& nameAndVisAtts = nameAndVisAttsByMaterial[material];
1542 const auto& name = nameAndVisAtts.fName;
1543 const auto& visAtts = nameAndVisAtts.fVisAtts;
1544 // Transfer positions into a vector ready for creating polyhedral surface
1545 std::vector<G4ThreeVector> positionsForPolyhedron;
1546 const auto& range = positionByMaterial.equal_range(material);
1547 for (auto posByMat = range.first; posByMat != range.second; ++posByMat) {
1548 const auto& position = posByMat->second;
1549 positionsForPolyhedron.push_back(position);
1550 ++nBoxes;
1551 }
1552 // The polyhedron will be in local coordinates
1553 // Add an empty place-holder to the map and get a reference to it
1554 auto& polyhedron = boxesByMaterial[material];
1555 // Replace with the desired polyhedron (uses efficient "move assignment")
1556 polyhedron = G4PolyhedronBoxMesh(sizeX,sizeY,sizeZ,positionsForPolyhedron);
1557 polyhedron.SetVisAttributes(visAtts);
1558 polyhedron.SetInfo(name);
1559
1560 if (print) {
1561 G4cout
1562 << std::setw(30) << std::left << name.substr(0,30) << std::right
1563 << ": " << std::setw(7) << nBoxes << " boxes"
1564 << " (" << std::setw(7) << 6*nBoxes << " faces)"
1565 << ": reduced to " << std::setw(7) << polyhedron.GetNoFacets() << " facets ("
1566 << std::setw(2) << std::fixed << std::setprecision(2) << 100*polyhedron.GetNoFacets()/(6*nBoxes)
1567 << "%): colour " << std::fixed << std::setprecision(2)
1568 << visAtts.GetColour() << std::defaultfloat
1569 << G4endl;
1570 }
1571
1572 nBoxesTotal += nBoxes;
1573 nFacetsTotal += polyhedron.GetNoFacets();
1574 }
1575
1576 if (print) {
1577 G4cout << "Total number of boxes: " << nBoxesTotal << " (" << 6*nBoxesTotal << " faces)"
1578 << ": reduced to " << nFacetsTotal << " facets ("
1579 << std::setw(2) << std::fixed << std::setprecision(2) << 100*nFacetsTotal/(6*nBoxesTotal) << "%)"
1580 << G4endl;
1581 }
1582 }
1583
1584 // Some subsequent expressions apply only to G4PhysicalVolumeModel
1585 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1586
1587 G4String parameterisationName;
1588 if (pPVModel) {
1589 parameterisationName = pPVModel->GetFullPVPath().back().GetPhysicalVolume()->GetName();
1590 }
1591
1592 // Draw the boxes by material
1593 // Now we transform to world coordinates
1595 for (const auto& entry: boxesByMaterial) {
1596 const auto& poly = entry.second;
1597 // The current "leaf" node in the PVPath is the parameterisation. Here it has
1598 // been converted into polyhedra by material. So...temporarily...change
1599 // its name to that of the material (whose name has been stored in Info)
1600 // so that its appearance in the scene tree of, e.g., G4OpenGLQtViewer, has
1601 // an appropriate name and its visibility and colour may be changed.
1602 if (pPVModel) {
1603 const auto& fullPVPath = pPVModel->GetFullPVPath();
1604 auto leafPV = fullPVPath.back().GetPhysicalVolume();
1605 leafPV->SetName(poly.GetInfo());
1606 }
1607 AddPrimitive(poly);
1608 }
1609 EndPrimitives ();
1610 // Restore parameterisation name
1611 if (pPVModel) {
1612 pPVModel->GetFullPVPath().back().GetPhysicalVolume()->SetName(parameterisationName);
1613 }
1614
1615 firstPrint = false;
1616 return;
1617}

Referenced by StandardSpecialMeshRendering().

◆ DrawEndOfRunModels()

void G4VSceneHandler::DrawEndOfRunModels ( )

Definition at line 811 of file G4VSceneHandler.cc.

812{
813 const std::vector<G4Scene::Model>& EORModelList =
814 fpScene -> GetEndOfRunModelList ();
815 std::size_t nModels = EORModelList.size();
816 if (nModels) {
818 pMP->SetEvent(0);
819 for (std::size_t i = 0; i < nModels; ++i) {
820 if (EORModelList[i].fActive) {
821 fpModel = EORModelList[i].fpModel;
822 fpModel -> SetModelingParameters(pMP);
823 fpModel -> DescribeYourselfTo (*this);
824 fpModel -> SetModelingParameters(0);
825 }
826 }
827 fpModel = 0;
828 delete pMP;
829 }
830}
void SetEvent(const G4Event *pEvent)
G4ModelingParameters * CreateModelingParameters()

Referenced by ProcessScene().

◆ DrawEvent()

void G4VSceneHandler::DrawEvent ( const G4Event event)

Definition at line 790 of file G4VSceneHandler.cc.

791{
792 const std::vector<G4Scene::Model>& EOEModelList =
793 fpScene -> GetEndOfEventModelList ();
794 std::size_t nModels = EOEModelList.size();
795 if (nModels) {
797 pMP->SetEvent(event);
798 for (std::size_t i = 0; i < nModels; ++i) {
799 if (EOEModelList[i].fActive) {
800 fpModel = EOEModelList[i].fpModel;
801 fpModel -> SetModelingParameters(pMP);
802 fpModel -> DescribeYourselfTo (*this);
803 fpModel -> SetModelingParameters(0);
804 }
805 }
806 fpModel = 0;
807 delete pMP;
808 }
809}

Referenced by ProcessScene().

◆ DrawTetMeshAsDots()

void G4VSceneHandler::DrawTetMeshAsDots ( const G4Mesh mesh)
protected

Definition at line 1619 of file G4VSceneHandler.cc.

1622{
1623 // Check
1624 if (mesh.GetMeshType() != G4Mesh::tetrahedron) {
1626 ed << "Called with mesh that is not a tetrahedron mesh:" << mesh;
1627 G4Exception("G4VSceneHandler::DrawTetMeshAsDots","visman0108",JustWarning,ed);
1628 return;
1629 }
1630
1631 static G4bool firstPrint = true;
1632 const auto& verbosity = G4VisManager::GetVerbosity();
1633 G4bool print = firstPrint && verbosity >= G4VisManager::errors;
1634
1635 if (print) {
1636 G4cout
1637 << "Special case drawing of tetrahedron mesh as dots"
1638 << '\n' << mesh
1639 << G4endl;
1640 }
1641
1642 const auto& container = mesh.GetContainerVolume();
1643
1644 // This map is static so that once filled it stays filled.
1645 static std::map<G4String,std::map<const G4Material*,G4Polymarker>> dotsByMaterialAndMesh;
1646 auto& dotsByMaterial = dotsByMaterialAndMesh[mesh.GetContainerVolume()->GetName()];
1647
1648 // Fill map if not already filled
1649 if (dotsByMaterial.empty()) {
1650
1651 // Get vertices and colour one cell at a time (using PseudoSceneForTetVertices).
1652 // The pseudo scene allows a "private" descent into the parameterisation.
1653 // Instantiate a temporary G4PhysicalVolumeModel
1655 tmpMP.SetCulling(true); // This avoids drawing transparent...
1656 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
1657 const G4bool useFullExtent = true; // To avoid calculating the extent
1658 G4PhysicalVolumeModel tmpPVModel
1659 (container,
1661 G4Transform3D(), // so that positions are in local coordinates
1662 &tmpMP,
1663 useFullExtent);
1664 // Accumulate information in temporary maps by material
1665 std::multimap<const G4Material*,std::vector<G4ThreeVector>> verticesByMaterial;
1666 std::map<const G4Material*,G4VSceneHandler::NameAndVisAtts> nameAndVisAttsByMaterial;
1667 // Instantiate a pseudo scene
1668 PseudoSceneForTetVertices pseudoScene
1669 (&tmpPVModel,mesh.GetMeshDepth(),verticesByMaterial,nameAndVisAttsByMaterial);
1670 // Make private descent into the parameterisation
1671 tmpPVModel.DescribeYourselfTo(pseudoScene);
1672 // Now we have a map of vertices by material.
1673 // Also a map of name and colour by material.
1674
1675 // Fill the permanent (static) map of dots by material
1676 G4int nDotsTotal = 0;
1677 for (const auto& entry: nameAndVisAttsByMaterial) {
1678 G4int nDots = 0;
1679 const auto& material = entry.first;
1680 const auto& nameAndVisAtts = nameAndVisAttsByMaterial[material];
1681 const auto& name = nameAndVisAtts.fName;
1682 const auto& visAtts = nameAndVisAtts.fVisAtts;
1683 G4Polymarker dots;
1684 dots.SetVisAttributes(visAtts);
1686 dots.SetSize(G4VMarker::screen,1.);
1687 dots.SetInfo(name);
1688 // Enter empty polymarker into the map
1689 dotsByMaterial[material] = dots;
1690 // Now fill it in situ
1691 auto& dotsInMap = dotsByMaterial[material];
1692 const auto& range = verticesByMaterial.equal_range(material);
1693 for (auto vByMat = range.first; vByMat != range.second; ++vByMat) {
1694 dotsInMap.push_back(GetPointInTet(vByMat->second));
1695 ++nDots;
1696 }
1697
1698 if (print) {
1699 G4cout
1700 << std::setw(30) << std::left << name.substr(0,30) << std::right
1701 << ": " << std::setw(7) << nDots << " dots"
1702 << ": colour " << std::fixed << std::setprecision(2)
1703 << visAtts.GetColour() << std::defaultfloat
1704 << G4endl;
1705 }
1706
1707 nDotsTotal += nDots;
1708 }
1709
1710 if (print) {
1711 G4cout << "Total number of dots: " << nDotsTotal << G4endl;
1712 }
1713 }
1714
1715 // Some subsequent expressions apply only to G4PhysicalVolumeModel
1716 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1717
1718 G4String parameterisationName;
1719 if (pPVModel) {
1720 parameterisationName = pPVModel->GetFullPVPath().back().GetPhysicalVolume()->GetName();
1721 }
1722
1723 // Draw the dots by material
1724 // Ensure they are "hidden", i.e., use the z-buffer as non-marker primitives do
1725 auto keepVP = fpViewer->GetViewParameters();
1726 auto vp = fpViewer->GetViewParameters();
1727 vp.SetMarkerHidden();
1729
1730 // Now we transform to world coordinates
1732 for (const auto& entry: dotsByMaterial) {
1733 const auto& dots = entry.second;
1734 // The current "leaf" node in the PVPath is the parameterisation. Here it has
1735 // been converted into polymarkers by material. So...temporarily...change
1736 // its name to that of the material (whose name has been stored in Info)
1737 // so that its appearance in the scene tree of, e.g., G4OpenGLQtViewer, has
1738 // an appropriate name and its visibility and colour may be changed.
1739 if (pPVModel) {
1740 const auto& fullPVPath = pPVModel->GetFullPVPath();
1741 auto leafPV = fullPVPath.back().GetPhysicalVolume();
1742 leafPV->SetName(dots.GetInfo());
1743 }
1744 AddPrimitive(dots);
1745 }
1746 EndPrimitives ();
1747
1748 // Restore view parameters
1749 fpViewer->SetViewParameters(keepVP);
1750 // Restore parameterisation name
1751 if (pPVModel) {
1752 pPVModel->GetFullPVPath().back().GetPhysicalVolume()->SetName(parameterisationName);
1753 }
1754
1755 firstPrint = false;
1756 return;
1757}
@ tetrahedron
Definition: G4Mesh.hh:57
G4ThreeVector GetPointInTet(const std::vector< G4ThreeVector > &vertices) const

Referenced by StandardSpecialMeshRendering().

◆ DrawTetMeshAsSurfaces()

void G4VSceneHandler::DrawTetMeshAsSurfaces ( const G4Mesh mesh)
protected

Definition at line 1759 of file G4VSceneHandler.cc.

1762{
1763 // Check
1764 if (mesh.GetMeshType() != G4Mesh::tetrahedron) {
1766 ed << "Called with mesh that is not a tetrahedron mesh:" << mesh;
1767 G4Exception("G4VSceneHandler::DrawTetMeshAsSurfaces","visman0108",JustWarning,ed);
1768 return;
1769 }
1770
1771 static G4bool firstPrint = true;
1772 const auto& verbosity = G4VisManager::GetVerbosity();
1773 G4bool print = firstPrint && verbosity >= G4VisManager::errors;
1774
1775 if (print) {
1776 G4cout
1777 << "Special case drawing of tetrahedron mesh as surfaces"
1778 << '\n' << mesh
1779 << G4endl;
1780 }
1781
1782 // This map is static so that once filled it stays filled.
1783 static std::map<G4String,std::map<const G4Material*,G4Polyhedron>> surfacesByMaterialAndMesh;
1784 auto& surfacesByMaterial = surfacesByMaterialAndMesh[mesh.GetContainerVolume()->GetName()];
1785
1786 // Fill map if not already filled
1787 if (surfacesByMaterial.empty()) {
1788
1789 // Get vertices and colour one cell at a time (using PseudoSceneForTetVertices).
1790 // The pseudo scene allows a "private" descent into the parameterisation.
1791 // Instantiate a temporary G4PhysicalVolumeModel
1793 tmpMP.SetCulling(true); // This avoids drawing transparent...
1794 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
1795 const G4bool useFullExtent = true; // To avoid calculating the extent
1796 G4PhysicalVolumeModel tmpPVModel
1797 (mesh.GetContainerVolume(),
1799 G4Transform3D(), // so that positions are in local coordinates
1800 &tmpMP,
1801 useFullExtent);
1802 // Accumulate information in temporary maps by material
1803 std::multimap<const G4Material*,std::vector<G4ThreeVector>> verticesByMaterial;
1804 std::map<const G4Material*,G4VSceneHandler::NameAndVisAtts> nameAndVisAttsByMaterial;
1805 // Instantiate a pseudo scene
1806 PseudoSceneForTetVertices pseudoScene
1807 (&tmpPVModel,mesh.GetMeshDepth(),verticesByMaterial,nameAndVisAttsByMaterial);
1808 // Make private descent into the parameterisation
1809 tmpPVModel.DescribeYourselfTo(pseudoScene);
1810 // Now we have a map of vertices by material.
1811 // Also a map of name and colour by material.
1812
1813 // Fill the permanent (static) map of surfaces by material
1814 G4int nTetsTotal = 0, nFacetsTotal = 0;
1815 for (const auto& entry: nameAndVisAttsByMaterial) {
1816 G4int nTets = 0;
1817 const auto& material = entry.first;
1818 const auto& nameAndVisAtts = nameAndVisAttsByMaterial[material];
1819 const auto& name = nameAndVisAtts.fName;
1820 const auto& visAtts = nameAndVisAtts.fVisAtts;
1821 // Transfer vertices into a vector ready for creating polyhedral surface
1822 std::vector<G4ThreeVector> verticesForPolyhedron;
1823 const auto& range = verticesByMaterial.equal_range(material);
1824 for (auto vByMat = range.first; vByMat != range.second; ++vByMat) {
1825 const std::vector<G4ThreeVector>& vertices = vByMat->second;
1826 for (const auto& vertex: vertices)
1827 verticesForPolyhedron.push_back(vertex);
1828 ++nTets;
1829 }
1830 // The polyhedron will be in local coordinates
1831 // Add an empty place-holder to the map and get a reference to it
1832 auto& polyhedron = surfacesByMaterial[material];
1833 // Replace with the desired polyhedron (uses efficient "move assignment")
1834 polyhedron = G4PolyhedronTetMesh(verticesForPolyhedron);
1835 polyhedron.SetVisAttributes(visAtts);
1836 polyhedron.SetInfo(name);
1837
1838 if (print) {
1839 G4cout
1840 << std::setw(30) << std::left << name.substr(0,30) << std::right
1841 << ": " << std::setw(7) << nTets << " tetrahedra"
1842 << " (" << std::setw(7) << 4*nTets << " faces)"
1843 << ": reduced to " << std::setw(7) << polyhedron.GetNoFacets() << " facets ("
1844 << std::setw(2) << std::fixed << std::setprecision(2) << 100*polyhedron.GetNoFacets()/(4*nTets)
1845 << "%): colour " << std::fixed << std::setprecision(2)
1846 << visAtts.GetColour() << std::defaultfloat
1847 << G4endl;
1848 }
1849
1850 nTetsTotal += nTets;
1851 nFacetsTotal += polyhedron.GetNoFacets();
1852 }
1853
1854 if (print) {
1855 G4cout << "Total number of tetrahedra: " << nTetsTotal << " (" << 4*nTetsTotal << " faces)"
1856 << ": reduced to " << nFacetsTotal << " facets ("
1857 << std::setw(2) << std::fixed << std::setprecision(2) << 100*nFacetsTotal/(4*nTetsTotal) << "%)"
1858 << G4endl;
1859 }
1860 }
1861
1862 // Some subsequent expressions apply only to G4PhysicalVolumeModel
1863 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1864
1865 G4String parameterisationName;
1866 if (pPVModel) {
1867 parameterisationName = pPVModel->GetFullPVPath().back().GetPhysicalVolume()->GetName();
1868 }
1869
1870 // Draw the surfaces by material
1871 // Now we transform to world coordinates
1873 for (const auto& entry: surfacesByMaterial) {
1874 const auto& poly = entry.second;
1875 // The current "leaf" node in the PVPath is the parameterisation. Here it has
1876 // been converted into polyhedra by material. So...temporarily...change
1877 // its name to that of the material (whose name has been stored in Info)
1878 // so that its appearance in the scene tree of, e.g., G4OpenGLQtViewer, has
1879 // an appropriate name and its visibility and colour may be changed.
1880 if (pPVModel) {
1881 const auto& fullPVPath = pPVModel->GetFullPVPath();
1882 auto leafPV = fullPVPath.back().GetPhysicalVolume();
1883 leafPV->SetName(poly.GetInfo());
1884 }
1885 AddPrimitive(poly);
1886 }
1887 EndPrimitives ();
1888
1889 // Restore parameterisation name
1890 if (pPVModel) {
1891 pPVModel->GetFullPVPath().back().GetPhysicalVolume()->SetName(parameterisationName);
1892 }
1893
1894 firstPrint = false;
1895 return;
1896}

Referenced by StandardSpecialMeshRendering().

◆ EndModeling()

◆ EndPrimitives()

◆ EndPrimitives2D()

void G4VSceneHandler::EndPrimitives2D ( )
virtual

Implements G4VGraphicsScene.

Reimplemented in G4HepRepFileSceneHandler, G4OpenGLImmediateSceneHandler, G4OpenGLSceneHandler, G4OpenGLStoredSceneHandler, and G4Qt3DSceneHandler.

Definition at line 201 of file G4VSceneHandler.cc.

201 {
202 if (fNestingDepth <= 0)
203 G4Exception("G4VSceneHandler::EndPrimitives2D",
204 "visman0104", FatalException, "Nesting error.");
209 }
210 fProcessing2D = false;
211}

Referenced by G4HepRepFileSceneHandler::EndPrimitives2D(), G4OpenGLSceneHandler::EndPrimitives2D(), and G4Qt3DSceneHandler::EndPrimitives2D().

◆ GetAuxEdgeVisible()

G4bool G4VSceneHandler::GetAuxEdgeVisible ( const G4VisAttributes pVisAttribs)

Definition at line 1152 of file G4VSceneHandler.cc.

1152 {
1153 G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
1154 if (pVisAttribs -> IsForceAuxEdgeVisible()) {
1155 isAuxEdgeVisible = pVisAttribs->IsForcedAuxEdgeVisible();
1156 }
1157 return isAuxEdgeVisible;
1158}
G4bool IsAuxEdgeVisible() const
G4bool IsForcedAuxEdgeVisible() const

Referenced by G4OpenGLSceneHandler::AddPrimitive().

◆ GetColor() [1/2]

const G4Colour & G4VSceneHandler::GetColor ( )

◆ GetColor() [2/2]

const G4Colour & G4VSceneHandler::GetColor ( const G4Visible )

◆ GetColour() [1/2]

const G4Colour & G4VSceneHandler::GetColour ( )

◆ GetColour() [2/2]

const G4Colour & G4VSceneHandler::GetColour ( const G4Visible visible)

Definition at line 1071 of file G4VSceneHandler.cc.

1071 {
1072 auto pVA = visible.GetVisAttributes();
1074 return pVA->GetColour();
1075}
const G4Colour & GetColour() const
const G4VisAttributes * GetVisAttributes() const

◆ GetCurrentViewer()

G4VViewer * G4VSceneHandler::GetCurrentViewer ( ) const

◆ GetDrawingStyle()

G4ViewParameters::DrawingStyle G4VSceneHandler::GetDrawingStyle ( const G4VisAttributes pVisAttribs)

Definition at line 1092 of file G4VSceneHandler.cc.

1093 {
1094 // Drawing style is normally determined by the view parameters, but
1095 // it can be overriddden by the ForceDrawingStyle flag in the vis
1096 // attributes.
1098 const G4ViewParameters::DrawingStyle viewerStyle = vp.GetDrawingStyle();
1099 G4ViewParameters::DrawingStyle resultantStyle = viewerStyle;
1100 if (pVisAttribs -> IsForceDrawingStyle ()) {
1102 pVisAttribs -> GetForcedDrawingStyle ();
1103 // This is complicated because if hidden line and surface removal
1104 // has been requested we wish to preserve this sometimes.
1105 switch (forcedStyle) {
1107 switch (viewerStyle) {
1108 case (G4ViewParameters::hlr):
1109 resultantStyle = G4ViewParameters::hlhsr;
1110 break;
1112 resultantStyle = G4ViewParameters::hsr;
1113 break;
1115 resultantStyle = G4ViewParameters::hsr;
1116 break;
1118 case (G4ViewParameters::hsr):
1119 break;
1120 }
1121 break;
1123 resultantStyle = G4ViewParameters::cloud;
1124 break;
1126 default:
1127 // But if forced style is wireframe, do it, because one of its
1128 // main uses is in displaying the consituent solids of a Boolean
1129 // solid and their surfaces overlap with the resulting Booean
1130 // solid, making a mess if hlr is specified.
1131 resultantStyle = G4ViewParameters::wireframe;
1132 break;
1133 }
1134 }
1135 return resultantStyle;
1136}

Referenced by G4OpenGLSceneHandler::AddPrimitive(), G4Qt3DSceneHandler::AddPrimitive(), G4ToolsSGSceneHandler::AddPrimitive(), G4VtkSceneHandler::AddPrimitiveTensorGlyph(), and RequestPrimitives().

◆ GetExtent()

const G4VisExtent & G4VSceneHandler::GetExtent ( ) const
virtual

Reimplemented from G4VGraphicsScene.

Definition at line 139 of file G4VSceneHandler.cc.

140{
141 if (fpScene) {
142 return fpScene->GetExtent();
143 } else {
144 static const G4VisExtent defaultExtent = G4VisExtent();
145 return defaultExtent;
146 }
147}

Referenced by G4VtkViewer::SetView(), G4ToolsSGViewer< SG_SESSION, SG_VIEWER >::wheel_rotate(), and G4Qt3DViewer::wheelEvent().

◆ GetGraphicsSystem()

◆ GetLineWidth()

G4double G4VSceneHandler::GetLineWidth ( const G4VisAttributes pVisAttribs)

Definition at line 1083 of file G4VSceneHandler.cc.

1084{
1085 G4double lineWidth = pVisAttribs->GetLineWidth();
1086 if (lineWidth < 1.) lineWidth = 1.;
1087 lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
1088 if (lineWidth < 1.) lineWidth = 1.;
1089 return lineWidth;
1090}
G4double GetLineWidth() const

Referenced by G4OpenGLSceneHandler::AddPrimitive().

◆ GetMarkerDiameter()

G4double G4VSceneHandler::GetMarkerDiameter ( const G4VMarker ,
MarkerSizeType  
)

◆ GetMarkerRadius()

G4double G4VSceneHandler::GetMarkerRadius ( const G4VMarker ,
MarkerSizeType  
)

◆ GetMarkerSize()

G4double G4VSceneHandler::GetMarkerSize ( const G4VMarker marker,
G4VSceneHandler::MarkerSizeType markerSizeType 
)

Definition at line 1160 of file G4VSceneHandler.cc.

1163{
1164 G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
1165 const G4VMarker& defaultMarker =
1166 fpViewer -> GetViewParameters().GetDefaultMarker();
1167 G4double size = userSpecified ?
1168 marker.GetWorldSize() : defaultMarker.GetWorldSize();
1169 if (size) {
1170 // Draw in world coordinates.
1171 markerSizeType = world;
1172 }
1173 else {
1174 size = userSpecified ?
1175 marker.GetScreenSize() : defaultMarker.GetScreenSize();
1176 // Draw in screen coordinates.
1177 markerSizeType = screen;
1178 }
1179 size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
1180 if (markerSizeType == screen && size < 1.) size = 1.;
1181 return size;
1182}
G4double GetScreenSize() const
G4double GetWorldSize() const

Referenced by G4VtkSceneHandler::AddPrimitive(), G4HepRepFileSceneHandler::AddPrimitive(), G4OpenGLSceneHandler::AddPrimitive(), G4OpenInventorSceneHandler::AddPrimitive(), G4Qt3DSceneHandler::AddPrimitive(), G4ToolsSGSceneHandler::AddPrimitive(), G4OpenGLQtViewer::DrawText(), G4OpenGLViewer::DrawText(), and G4OpenGLXViewer::DrawText().

◆ GetMarkForClearingTransientStore()

G4bool G4VSceneHandler::GetMarkForClearingTransientStore ( ) const

◆ GetModel()

G4VModel * G4VSceneHandler::GetModel ( ) const

◆ GetName()

◆ GetNoOfSides()

G4int G4VSceneHandler::GetNoOfSides ( const G4VisAttributes pVisAttribs)

Definition at line 1184 of file G4VSceneHandler.cc.

1185{
1186 // No. of sides (lines segments per circle) is normally determined
1187 // by the view parameters, but it can be overriddden by the
1188 // ForceLineSegmentsPerCircle in the vis attributes.
1189 G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
1190 if (pVisAttribs) {
1191 if (pVisAttribs->IsForceLineSegmentsPerCircle())
1192 lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
1193 if (lineSegmentsPerCircle < pVisAttribs->GetMinLineSegmentsPerCircle()) {
1194 lineSegmentsPerCircle = pVisAttribs->GetMinLineSegmentsPerCircle();
1195 G4warn <<
1196 "G4VSceneHandler::GetNoOfSides: attempt to set the"
1197 "\nnumber of line segments per circle < " << lineSegmentsPerCircle
1198 << "; forced to " << pVisAttribs->GetMinLineSegmentsPerCircle() << G4endl;
1199 }
1200 }
1201 return lineSegmentsPerCircle;
1202}
G4bool IsForceLineSegmentsPerCircle() const
G4int GetForcedLineSegmentsPerCircle() const
static G4int GetMinLineSegmentsPerCircle()

Referenced by G4OpenGLSceneHandler::AddPrimitive(), and RequestPrimitives().

◆ GetNumberOfCloudPoints()

G4int G4VSceneHandler::GetNumberOfCloudPoints ( const G4VisAttributes pVisAttribs) const

Definition at line 1138 of file G4VSceneHandler.cc.

1139 {
1140 // Returns no of cloud points from current view parameters, unless the user
1141 // has forced through the vis attributes, thereby over-riding the
1142 // current view parameter.
1143 G4int numberOfCloudPoints = fpViewer->GetViewParameters().GetNumberOfCloudPoints();
1144 if (pVisAttribs -> IsForceDrawingStyle() &&
1145 pVisAttribs -> GetForcedDrawingStyle() == G4VisAttributes::cloud &&
1146 pVisAttribs -> GetForcedNumberOfCloudPoints() > 0) {
1147 numberOfCloudPoints = pVisAttribs -> GetForcedNumberOfCloudPoints();
1148 }
1149 return numberOfCloudPoints;
1150}

Referenced by RequestPrimitives().

◆ GetObjectTransformation()

const G4Transform3D & G4VSceneHandler::GetObjectTransformation ( ) const

◆ GetPointInBox()

G4ThreeVector G4VSceneHandler::GetPointInBox ( const G4ThreeVector pos,
G4double  halfX,
G4double  halfY,
G4double  halfZ 
) const
protected

Definition at line 1899 of file G4VSceneHandler.cc.

1903{
1904 G4double x = pos.getX() + (2.*G4QuickRand() - 1.)*halfX;
1905 G4double y = pos.getY() + (2.*G4QuickRand() - 1.)*halfY;
1906 G4double z = pos.getZ() + (2.*G4QuickRand() - 1.)*halfZ;
1907 return G4ThreeVector(x, y, z);
1908}
G4double G4QuickRand()
Definition: G4QuickRand.hh:34
CLHEP::Hep3Vector G4ThreeVector

Referenced by Draw3DRectMeshAsDots().

◆ GetPointInTet()

G4ThreeVector G4VSceneHandler::GetPointInTet ( const std::vector< G4ThreeVector > &  vertices) const
protected

Definition at line 1911 of file G4VSceneHandler.cc.

1912{
1913 G4double p = G4QuickRand();
1914 G4double q = G4QuickRand();
1915 G4double r = G4QuickRand();
1916 if (p + q > 1.)
1917 {
1918 p = 1. - p;
1919 q = 1. - q;
1920 }
1921 if (q + r > 1.)
1922 {
1923 G4double tmp = r;
1924 r = 1. - p - q;
1925 q = 1. - tmp;
1926 }
1927 else if (p + q + r > 1.)
1928 {
1929 G4double tmp = r;
1930 r = p + q + r - 1.;
1931 p = 1. - q - tmp;
1932 }
1933 G4double a = 1. - p - q - r;
1934 return vertices[0]*a + vertices[1]*p + vertices[2]*q + vertices[3]*r;
1935}

Referenced by DrawTetMeshAsDots().

◆ GetScene()

◆ GetSceneHandlerId()

G4int G4VSceneHandler::GetSceneHandlerId ( ) const

◆ GetTextColor()

const G4Colour & G4VSceneHandler::GetTextColor ( const G4Text )

◆ GetTextColour()

const G4Colour & G4VSceneHandler::GetTextColour ( const G4Text text)

◆ GetTransientsDrawnThisEvent()

G4bool G4VSceneHandler::GetTransientsDrawnThisEvent ( ) const

◆ GetTransientsDrawnThisRun()

G4bool G4VSceneHandler::GetTransientsDrawnThisRun ( ) const

◆ GetViewCount()

G4int G4VSceneHandler::GetViewCount ( ) const

◆ GetViewerList()

const G4ViewerList & G4VSceneHandler::GetViewerList ( ) const

◆ IncrementViewCount()

G4int G4VSceneHandler::IncrementViewCount ( )

◆ IsReadyForTransients()

G4bool G4VSceneHandler::IsReadyForTransients ( ) const

◆ LoadAtts()

void G4VSceneHandler::LoadAtts ( const G4Visible visible,
G4AttHolder holder 
)
protected

Definition at line 1005 of file G4VSceneHandler.cc.

1006{
1007 // Load G4Atts from G4VisAttributes, if any...
1008 const G4VisAttributes* va = visible.GetVisAttributes();
1009 if (va) {
1010 const std::map<G4String,G4AttDef>* vaDefs =
1011 va->GetAttDefs();
1012 if (vaDefs) {
1013 holder->AddAtts(visible.GetVisAttributes()->CreateAttValues(), vaDefs);
1014 }
1015 }
1016
1017 G4PhysicalVolumeModel* pPVModel =
1018 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1019 if (pPVModel) {
1020 // Load G4Atts from G4PhysicalVolumeModel...
1021 const std::map<G4String,G4AttDef>* pvDefs = pPVModel->GetAttDefs();
1022 if (pvDefs) {
1023 holder->AddAtts(pPVModel->CreateCurrentAttValues(), pvDefs);
1024 }
1025 }
1026
1027 G4TrajectoriesModel* trajModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
1028 if (trajModel) {
1029 // Load G4Atts from trajectory model...
1030 const std::map<G4String,G4AttDef>* trajModelDefs = trajModel->GetAttDefs();
1031 if (trajModelDefs) {
1032 holder->AddAtts(trajModel->CreateCurrentAttValues(), trajModelDefs);
1033 }
1034 // Load G4Atts from trajectory...
1035 const G4VTrajectory* traj = trajModel->GetCurrentTrajectory();
1036 if (traj) {
1037 const std::map<G4String,G4AttDef>* trajDefs = traj->GetAttDefs();
1038 if (trajDefs) {
1039 holder->AddAtts(traj->CreateAttValues(), trajDefs);
1040 }
1041 G4int nPoints = traj->GetPointEntries();
1042 for (G4int i = 0; i < nPoints; ++i) {
1043 G4VTrajectoryPoint* trajPoint = traj->GetPoint(i);
1044 if (trajPoint) {
1045 const std::map<G4String,G4AttDef>* pointDefs = trajPoint->GetAttDefs();
1046 if (pointDefs) {
1047 holder->AddAtts(trajPoint->CreateAttValues(), pointDefs);
1048 }
1049 }
1050 }
1051 }
1052 }
1053
1054 G4HitsModel* hitsModel = dynamic_cast<G4HitsModel*>(fpModel);
1055 if (hitsModel) {
1056 // Load G4Atts from hit...
1057 const G4VHit* hit = hitsModel->GetCurrentHit();
1058 const std::map<G4String,G4AttDef>* hitsDefs = hit->GetAttDefs();
1059 if (hitsDefs) {
1060 holder->AddAtts(hit->CreateAttValues(), hitsDefs);
1061 }
1062 }
1063}
void AddAtts(const std::vector< G4AttValue > *values, const std::map< G4String, G4AttDef > *defs)
Definition: G4AttHolder.hh:66
const G4VHit * GetCurrentHit() const
Definition: G4HitsModel.hh:57
std::vector< G4AttValue > * CreateCurrentAttValues() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4VTrajectory * GetCurrentTrajectory() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4VHit.hh:64
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
Definition: G4VHit.hh:58
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
const std::vector< G4AttValue > * CreateAttValues() const

Referenced by G4OpenInventorSceneHandler::AddPrimitive().

◆ PostAddSolid()

void G4VSceneHandler::PostAddSolid ( )
virtual

Implements G4VGraphicsScene.

Reimplemented in G4Qt3DSceneHandler, and G4VTreeSceneHandler.

Definition at line 156 of file G4VSceneHandler.cc.

156 {
157 fpVisAttribs = 0;
158 fProcessingSolid = false;
162 }
163}

Referenced by AddCompound(), and G4Qt3DSceneHandler::PostAddSolid().

◆ PreAddSolid()

void G4VSceneHandler::PreAddSolid ( const G4Transform3D objectTransformation,
const G4VisAttributes visAttribs 
)
virtual

◆ ProcessScene()

void G4VSceneHandler::ProcessScene ( )
protectedvirtual

Reimplemented in G4OpenGLSceneHandler.

Definition at line 646 of file G4VSceneHandler.cc.

647{
648 // Assumes graphics database store has already been cleared if
649 // relevant for the particular scene handler.
650
651 if(!fpScene)
652 return;
653
655 {
656 G4Exception("G4VSceneHandler::ProcessScene", "visman0106", JustWarning,
657 "The scene has no extent.");
658 }
659
661
662 if(!visManager->GetConcreteInstance())
663 return;
664
665 G4VisManager::Verbosity verbosity = visManager->GetVerbosity();
666
667 fReadyForTransients = false;
668
669 // Reset fMarkForClearingTransientStore. (Leaving
670 // fMarkForClearingTransientStore true causes problems with
671 // recomputing transients below.) Restore it again at end...
672 G4bool tmpMarkForClearingTransientStore = fMarkForClearingTransientStore;
674
675 // Traverse geometry tree and send drawing primitives to window(s).
676
677 const std::vector<G4Scene::Model>& runDurationModelList =
679
680 if(runDurationModelList.size())
681 {
682 if(verbosity >= G4VisManager::confirmations)
683 {
684 G4cout << "Traversing scene data..." << G4endl;
685 }
686
688
689 // Create modeling parameters from view parameters...
691
692 for(std::size_t i = 0; i < runDurationModelList.size(); ++i)
693 {
694 if(runDurationModelList[i].fActive)
695 {
696 fpModel = runDurationModelList[i].fpModel;
699 // To see the extents of each model represented as wireframe boxes,
700 // uncomment the next line and DrawExtent in namespace above
701 // DrawExtent(fpModel);
703 }
704 }
705
706 fpModel = 0;
707 delete pMP;
708
709 EndModeling();
710 }
711
712 fReadyForTransients = true;
713
714 // Refresh event from end-of-event model list.
715 // Allow only in Idle or GeomClosed state...
717 G4ApplicationState state = stateManager->GetCurrentState();
718 if(state == G4State_Idle || state == G4State_GeomClosed)
719 {
720 visManager->SetEventRefreshing(true);
721
722 if(visManager->GetRequestedEvent())
723 {
724 DrawEvent(visManager->GetRequestedEvent());
725 }
726 else
727 {
729 if(runManager)
730 {
731 const G4Run* run = runManager->GetCurrentRun();
732 const std::vector<const G4Event*>* events =
733 run ? run->GetEventVector() : 0;
734 std::size_t nKeptEvents = 0;
735 if(events)
736 nKeptEvents = events->size();
737 if(nKeptEvents)
738 {
740 {
741 if(verbosity >= G4VisManager::confirmations)
742 {
743 G4cout << "Refreshing event..." << G4endl;
744 }
745 const G4Event* event = 0;
746 if(events && events->size())
747 event = events->back();
748 if(event)
749 DrawEvent(event);
750 }
751 else
752 { // Accumulating events.
753
754 if(verbosity >= G4VisManager::confirmations)
755 {
756 G4cout << "Refreshing events in run..." << G4endl;
757 }
758 for(const auto& event : *events)
759 {
760 if(event)
761 DrawEvent(event);
762 }
763
765 {
766 if(verbosity >= G4VisManager::warnings)
767 {
768 G4warn << "WARNING: Cannot refresh events accumulated over more"
769 "\n than one runs. Refreshed just the last run."
770 << G4endl;
771 }
772 }
773 }
774 }
775 }
776 }
777 visManager->SetEventRefreshing(false);
778 }
779
780 // Refresh end-of-run model list.
781 // Allow only in Idle or GeomClosed state...
782 if(state == G4State_Idle || state == G4State_GeomClosed)
783 {
785 }
786
787 fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
788}
G4ApplicationState
@ G4State_Idle
@ G4State_GeomClosed
static G4RunManager * GetMasterRunManager()
const G4Run * GetCurrentRun() const
Definition: G4Run.hh:49
const std::vector< const G4Event * > * GetEventVector() const
Definition: G4Run.hh:96
const std::vector< Model > & GetRunDurationModelList() const
G4bool GetRefreshAtEndOfEvent() const
G4bool GetRefreshAtEndOfRun() const
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
void SetModelingParameters(const G4ModelingParameters *)
virtual void DescribeYourselfTo(G4VGraphicsScene &)=0
virtual void BeginModeling()
void DrawEvent(const G4Event *)
virtual void EndModeling()
static G4VVisManager * GetConcreteInstance()
static const G4VisExtent & GetNullExtent()
Definition: G4VisExtent.cc:60
void SetEventRefreshing(G4bool)
const G4Event * GetRequestedEvent() const

Referenced by G4OpenGLSceneHandler::ProcessScene(), and G4VViewer::ProcessView().

◆ RemoveViewerFromList()

void G4VSceneHandler::RemoveViewerFromList ( G4VViewer pView)

Definition at line 502 of file G4VSceneHandler.cc.

502 {
503 fViewerList.remove(pViewer); // Does nothing if already removed
504 // And reset current viewer
505 auto visManager = G4VisManager::GetInstance();
506 visManager->SetCurrentViewer(nullptr);
507}
void remove(G4VViewer *)
Definition: G4ViewerList.cc:30

Referenced by G4VViewer::~G4VViewer().

◆ RequestPrimitives()

void G4VSceneHandler::RequestPrimitives ( const G4VSolid solid)
protectedvirtual

Reimplemented in G4ASCIITreeSceneHandler, and G4RayTracerSceneHandler.

Definition at line 525 of file G4VSceneHandler.cc.

526{
527 // Sometimes solids that have no substance get requested. They may
528 // be part of the geometry tree but have been "spirited away", for
529 // example by a Boolean subtraction in wich the original volume
530 // is entirely inside the subtractor.
531 // The problem is that the Boolean Processor still returns a
532 // polyhedron in these cases (IMHO it should not), so the
533 // workaround is to return before the damage is done.
534 auto pSolid = &solid;
535 auto pBooleanSolid = dynamic_cast<const G4BooleanSolid*>(pSolid);
536 if (pBooleanSolid) {
537 G4ThreeVector bmin, bmax;
538 pBooleanSolid->BoundingLimits(bmin, bmax);
539 G4bool isGood = false;
540 for (G4int i=0; i<100000; ++i) {
541 G4double x = bmin.x() + (bmax.x() - bmin.x())*G4QuickRand();
542 G4double y = bmin.y() + (bmax.y() - bmin.y())*G4QuickRand();
543 G4double z = bmin.z() + (bmax.z() - bmin.z())*G4QuickRand();
544 if (pBooleanSolid->Inside(G4ThreeVector(x,y,z)) == kInside) {
545 isGood = true;
546 break;
547 }
548 }
549 if (!isGood) return;
550 }
551
554
555 switch (style) {
556 default:
561 {
562 // Use polyhedral representation
564 G4Polyhedron* pPolyhedron = solid.GetPolyhedron ();
566 if (pPolyhedron) {
567 pPolyhedron -> SetVisAttributes (fpVisAttribs);
569 AddPrimitive (*pPolyhedron);
570 EndPrimitives ();
571 break;
572 } else { // Print warnings and drop through to cloud
574 static std::set<const G4VSolid*> problematicSolids;
575 if (verbosity >= G4VisManager::errors &&
576 problematicSolids.find(&solid) == problematicSolids.end()) {
577 problematicSolids.insert(&solid);
578 G4warn <<
579 "ERROR: G4VSceneHandler::RequestPrimitives"
580 "\n Polyhedron not available for " << solid.GetName ();
581 G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
582 if (pPVModel) {
583 G4warn << "\n Touchable path: " << pPVModel->GetFullPVPath();
584 }
585 static G4bool explanation = false;
586 if (!explanation) {
587 explanation = true;
588 G4warn <<
589 "\n This means it cannot be visualized in the usual way on most systems."
590 "\n 1) The solid may not have implemented the CreatePolyhedron method."
591 "\n 2) For Boolean solids, the BooleanProcessor, which attempts to create"
592 "\n the resultant polyhedron, may have failed."
593 "\n Try RayTracer. It uses Geant4's tracking algorithms instead.";
594 }
595 G4warn << "\n Drawing solid with cloud of points.";
596 G4warn << G4endl;
597 }
598 }
599 }
600 [[fallthrough]];
601
603 {
604 // Form solid out of cloud of dots on surface of solid
605 G4Polymarker dots;
606 // Note: OpenGL has a fast implementation of polymarker so it's better
607 // to build a polymarker rather than add a succession of circles.
608 // And anyway, in Qt, in the latter case each circle would be a scene-tree
609 // entry, something we would want to avoid.
612 dots.SetSize(G4VMarker::screen,1.);
613 G4int numberOfCloudPoints = GetNumberOfCloudPoints(fpVisAttribs);
614 if (numberOfCloudPoints <= 0) numberOfCloudPoints = vp.GetNumberOfCloudPoints();
615 for (G4int i = 0; i < numberOfCloudPoints; ++i) {
617 dots.push_back(p);
618 }
620 AddPrimitive(dots);
621 EndPrimitives ();
622 break;
623 }
624 }
625}
double z() const
double x() const
double y() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
G4int GetNumberOfCloudPoints(const G4VisAttributes *) const
G4int GetNoOfSides(const G4VisAttributes *)
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4String GetName() const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:705
static void SetNumberOfRotationSteps(G4int n)
static void ResetNumberOfRotationSteps()
@ kInside
Definition: geomdefs.hh:70

Referenced by AddSolidT(), and AddSolidWithAuxiliaryEdges().

◆ SetCurrentViewer()

void G4VSceneHandler::SetCurrentViewer ( G4VViewer )

◆ SetMarkForClearingTransientStore()

◆ SetModel()

void G4VSceneHandler::SetModel ( G4VModel )

◆ SetName()

void G4VSceneHandler::SetName ( const G4String )

◆ SetObjectTransformation()

void G4VSceneHandler::SetObjectTransformation ( const G4Transform3D )

◆ SetScene()

void G4VSceneHandler::SetScene ( G4Scene pScene)
virtual

Reimplemented in G4OpenGLStoredQtSceneHandler.

Definition at line 516 of file G4VSceneHandler.cc.

516 {
517 fpScene = pScene;
518 // Notify all viewers that a kernel visit is required.
520 for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
521 (*i) -> SetNeedKernelVisit (true);
522 }
523}
std::vector< G4VViewer * >::iterator G4ViewerListIterator
Definition: G4ViewerList.hh:42

Referenced by G4OpenGLStoredQtSceneHandler::SetScene().

◆ SetTransientsDrawnThisEvent()

void G4VSceneHandler::SetTransientsDrawnThisEvent ( G4bool  )

◆ SetTransientsDrawnThisRun()

void G4VSceneHandler::SetTransientsDrawnThisRun ( G4bool  )

◆ SetViewerList()

G4ViewerList & G4VSceneHandler::SetViewerList ( )

◆ StandardSpecialMeshRendering()

void G4VSceneHandler::StandardSpecialMeshRendering ( const G4Mesh mesh)
protected

Definition at line 1271 of file G4VSceneHandler.cc.

1275{
1276 G4bool implemented = false;
1277 switch (mesh.GetMeshType()) {
1278 case G4Mesh::rectangle: [[fallthrough]];
1282 Draw3DRectMeshAsDots(mesh); // Rectangular 3-deep mesh as dots
1283 implemented = true;
1284 break;
1286 Draw3DRectMeshAsSurfaces(mesh); // Rectangular 3-deep mesh as surfaces
1287 implemented = true;
1288 break;
1289 }
1290 break;
1294 DrawTetMeshAsDots(mesh); // Tetrahedron mesh as dots
1295 implemented = true;
1296 break;
1298 DrawTetMeshAsSurfaces(mesh); // Tetrahedron mesh as surfaces
1299 implemented = true;
1300 break;
1301 }
1302 break;
1303 case G4Mesh::cylinder: [[fallthrough]];
1304 case G4Mesh::sphere: [[fallthrough]];
1305 case G4Mesh::invalid: break;
1306 }
1307 if (implemented) {
1308 // Draw container if not marked invisible...
1309 auto container = mesh.GetContainerVolume();
1310 auto containerLogical = container->GetLogicalVolume();
1311 auto containerVisAtts = containerLogical->GetVisAttributes();
1312 if (containerVisAtts == nullptr || containerVisAtts->IsVisible()) {
1313 auto solid = containerLogical->GetSolid();
1314 auto polyhedron = solid->GetPolyhedron();
1315 // Always draw as wireframe
1316 G4VisAttributes tmpVisAtts;
1317 if (containerVisAtts != nullptr) tmpVisAtts = *containerVisAtts;
1318 tmpVisAtts.SetForceWireframe();
1319 polyhedron->SetVisAttributes(tmpVisAtts);
1321 AddPrimitive(*polyhedron);
1322 EndPrimitives();
1323 }
1324 } else {
1325 // Invoke base class function
1327 }
1328 return;
1329}
const G4VisAttributes * GetVisAttributes() const
@ invalid
Definition: G4Mesh.hh:52
@ sphere
Definition: G4Mesh.hh:56
@ cylinder
Definition: G4Mesh.hh:55
void DrawTetMeshAsSurfaces(const G4Mesh &)
void Draw3DRectMeshAsDots(const G4Mesh &)
void DrawTetMeshAsDots(const G4Mesh &)
void Draw3DRectMeshAsSurfaces(const G4Mesh &)
virtual void AddCompound(const G4VTrajectory &)
void SetForceWireframe(G4bool=true)

Referenced by G4OpenGLSceneHandler::AddCompound(), G4OpenInventorSceneHandler::AddCompound(), G4Qt3DSceneHandler::AddCompound(), and G4ToolsSGSceneHandler::AddCompound().

Friends And Related Function Documentation

◆ G4VViewer

friend class G4VViewer
friend

Definition at line 54 of file G4VSceneHandler.hh.

◆ operator<<

std::ostream & operator<< ( std::ostream &  os,
const G4VSceneHandler s 
)
friend

Definition at line 1204 of file G4VSceneHandler.cc.

1204 {
1205
1206 os << "Scene handler " << sh.fName << " has "
1207 << sh.fViewerList.size () << " viewer(s):";
1208 for (std::size_t i = 0; i < sh.fViewerList.size (); ++i) {
1209 os << "\n " << *(sh.fViewerList [i]);
1210 }
1211
1212 if (sh.fpScene) {
1213 os << "\n " << *sh.fpScene;
1214 }
1215 else {
1216 os << "\n This scene handler currently has no scene.";
1217 }
1218
1219 return os;
1220}

Member Data Documentation

◆ fIdentityTransformation

const G4Transform3D G4VSceneHandler::fIdentityTransformation
protected

Definition at line 453 of file G4VSceneHandler.hh.

◆ fMarkForClearingTransientStore

G4bool G4VSceneHandler::fMarkForClearingTransientStore
protected

Definition at line 441 of file G4VSceneHandler.hh.

Referenced by ProcessScene().

◆ fName

G4String G4VSceneHandler::fName
protected

Definition at line 436 of file G4VSceneHandler.hh.

◆ fNestingDepth

G4int G4VSceneHandler::fNestingDepth
protected

Definition at line 451 of file G4VSceneHandler.hh.

Referenced by BeginPrimitives(), EndPrimitives(), and EndPrimitives2D().

◆ fObjectTransformation

◆ fpModel

◆ fProcessing2D

◆ fProcessingSolid

G4bool G4VSceneHandler::fProcessingSolid
protected

◆ fpScene

◆ fpViewer

◆ fpVisAttribs

◆ fReadyForTransients

◆ fSceneHandlerId

const G4int G4VSceneHandler::fSceneHandlerId
protected

Definition at line 435 of file G4VSceneHandler.hh.

◆ fSystem

G4VGraphicsSystem& G4VSceneHandler::fSystem
protected

Definition at line 434 of file G4VSceneHandler.hh.

Referenced by AddPrimitive().

◆ fTransientsDrawnThisEvent

G4bool G4VSceneHandler::fTransientsDrawnThisEvent
protected

Definition at line 444 of file G4VSceneHandler.hh.

Referenced by EndPrimitives(), EndPrimitives2D(), and PostAddSolid().

◆ fTransientsDrawnThisRun

G4bool G4VSceneHandler::fTransientsDrawnThisRun
protected

Definition at line 445 of file G4VSceneHandler.hh.

Referenced by EndPrimitives(), EndPrimitives2D(), and PostAddSolid().

◆ fViewCount

G4int G4VSceneHandler::fViewCount
protected

Definition at line 437 of file G4VSceneHandler.hh.

◆ fViewerList

G4ViewerList G4VSceneHandler::fViewerList
protected

Definition at line 438 of file G4VSceneHandler.hh.

Referenced by AddViewerToList(), RemoveViewerFromList(), and SetScene().


The documentation for this class was generated from the following files: