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

#include <G4Qt3DSceneHandler.hh>

+ Inheritance diagram for G4Qt3DSceneHandler:

Public Member Functions

 G4Qt3DSceneHandler (G4VGraphicsSystem &system, const G4String &name)
 
virtual ~G4Qt3DSceneHandler ()
 
void PreAddSolid (const G4Transform3D &objectTransformation, const G4VisAttributes &)
 
void PostAddSolid ()
 
void BeginPrimitives2D (const G4Transform3D &objectTransformation)
 
void EndPrimitives2D ()
 
void BeginPrimitives (const G4Transform3D &objectTransformation)
 
void EndPrimitives ()
 
void AddPrimitive (const G4Polyline &)
 
void AddPrimitive (const G4Polymarker &)
 
void AddPrimitive (const G4Text &)
 
void AddPrimitive (const G4Circle &)
 
void AddPrimitive (const G4Square &)
 
void AddPrimitive (const G4Polyhedron &)
 
void AddCompound (const G4Mesh &)
 
void ClearStore ()
 
void ClearTransientStore ()
 
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 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 &)
 
- Public Member Functions inherited from G4VSceneHandler
 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

void EstablishG4Qt3DQEntities ()
 
G4Qt3DQEntityCreateNewNode ()
 
- Protected Member Functions inherited from G4VSceneHandler
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

Qt3DCore::QEntity * fpQt3DScene
 
Qt3DCore::QEntity * fpTransientObjects
 
Qt3DCore::QEntity * fpPersistentObjects
 
std::vector< G4Qt3DQEntity * > fpPhysicalVolumeObjects
 
- Protected Attributes inherited from G4VSceneHandler
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
 

Static Protected Attributes

static G4int fSceneIdCount = 0
 

Friends

class G4Qt3DViewer
 

Additional Inherited Members

- Public Types inherited from G4VSceneHandler
enum  MarkerSizeType { world , screen }
 

Detailed Description

Definition at line 40 of file G4Qt3DSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4Qt3DSceneHandler()

G4Qt3DSceneHandler::G4Qt3DSceneHandler ( G4VGraphicsSystem system,
const G4String name 
)

Definition at line 70 of file G4Qt3DSceneHandler.cc.

72: G4VSceneHandler(system, fSceneIdCount++, name)
73{
74#ifdef G4QT3DDEBUG
75 G4cout << "G4Qt3DSceneHandler::G4Qt3DSceneHandler called" << G4endl;
76#endif
77 fpQt3DScene = new Qt3DCore::QEntity;
78 fpQt3DScene->setObjectName("G4Qt3DSceneRoot");
80}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Qt3DCore::QEntity * fpQt3DScene
static G4int fSceneIdCount

◆ ~G4Qt3DSceneHandler()

G4Qt3DSceneHandler::~G4Qt3DSceneHandler ( )
virtual

Definition at line 82 of file G4Qt3DSceneHandler.cc.

83{
84 // Doesn't like this - it gives BAD_ACCESS in delete_entity_recursively.
85 // Curiously the delete traceback shows three calls to this recursively:
86 /*#1 0x0000000100411906 in (anonymous namespace)::delete_entity_recursively(Qt3DCore::QNode*) at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:60
87 #2 0x0000000100411840 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:169
88 #3 0x0000000100411fc5 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:168
89 #4 0x0000000100411fe9 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:168
90 #5 0x0000000101032510 in G4VisManager::~G4VisManager() at /Users/johna/Geant4/geant4-dev/source/visualization/management/src/G4VisManager.cc:214
91 #6 0x0000000100013885 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
92 #7 0x00000001000119a5 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
93 #8 0x00000001000119c9 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
94 #9 0x00000001000117dd in main at /Users/johna/Geant4/geant4-dev/examples/basic/B1/exampleB1.cc:108
95 */
96 //if (fpQt3DScene) delete_entity_recursively(fpQt3DScene);
97}

Member Function Documentation

◆ AddCompound() [1/7]

void G4VSceneHandler::AddCompound ( const G4Mesh mesh)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 133 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
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetContainerVolume() const
Definition: G4Mesh.hh:73
const G4Transform3D & GetTransform() const
Definition: G4Mesh.hh:77
void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4VViewer * fpViewer
const G4ViewParameters & GetViewParameters() const
SMROption GetSpecialMeshRenderingOption() const
void SetColour(const G4Colour &)
void SetVisibility(G4bool=true)

◆ AddCompound() [2/7]

void G4Qt3DSceneHandler::AddCompound ( const G4Mesh mesh)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 1056 of file G4Qt3DSceneHandler.cc.

1057{
1059}
void StandardSpecialMeshRendering(const G4Mesh &)

◆ AddCompound() [3/7]

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

Reimplemented from G4VSceneHandler.

Definition at line 131 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
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() [4/7]

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

Reimplemented from G4VSceneHandler.

Definition at line 132 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() [5/7]

void G4VSceneHandler::AddCompound ( const G4VDigi digi)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 130 of file G4VSceneHandler.cc.

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

◆ AddCompound() [6/7]

void G4VSceneHandler::AddCompound ( const G4VHit hit)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 129 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() [7/7]

void G4VSceneHandler::AddCompound ( const G4VTrajectory traj)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 128 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

◆ AddPrimitive() [1/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Circle circle)
virtual

Implements G4VSceneHandler.

Definition at line 560 of file G4Qt3DSceneHandler.cc.

561{
562#ifdef G4QT3DDEBUG
563 G4cout <<
564 "G4Qt3DSceneHandler::AddPrimitive(const G4Circle& circle) called.\n"
565 << circle
566 << G4endl;
567#endif
568
569#ifdef G4QT3DDEBUG
570 MarkerSizeType sizeType;
571 G4double size = GetMarkerSize (circle, sizeType);
572 switch (sizeType) {
573 default:
574 case screen:
575 // Draw in screen coordinates.
576 G4cout << "screen";
577 break;
578 case world:
579 // Draw in world coordinates.
580 G4cout << "world";
581 break;
582 }
583 G4cout << " size: " << size << G4endl;
584#endif
585
586 auto currentNode = CreateNewNode();
587 if (!currentNode) {
588 static G4bool first = true;
589 if (first) {
590 first = false;
591 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Circle&)",
592 "qt3d-0003", JustWarning,
593 "No available node!");
594 }
595 return;
596 }
597
599
602
603 const auto& colour = fpVisAttribs->GetColour();
604 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
605 material->setObjectName("materialForCircle");
606 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
607 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
608
609 auto sphereMesh = new Qt3DExtras::QSphereMesh;
610 sphereMesh->setObjectName("sphereMesh");
611 G4double radius;
612 if (circle.GetSizeType() == G4VMarker::world ) {
613 radius =circle.GetWorldRadius();
614 } else { // Screen-size or none
615 // Not figured out how to do screen-size, so use scene extent
616 const G4double scale = 200.; // Roughly pixles per scene
617 radius = circle.GetScreenRadius()*fpScene->GetExtent().GetExtentRadius()/scale;
618 }
619 sphereMesh->setRadius(radius);
620
621 auto currentEntity = new Qt3DCore::QEntity(currentNode);
622 currentEntity->addComponent(material);
623 currentEntity->addComponent(transform);
624 currentEntity->addComponent(sphereMesh);
625}
@ JustWarning
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:83
G4Qt3DQEntity * CreateNewNode()
const G4VisExtent & GetExtent() const
SizeType GetSizeType() const
Definition: G4VMarker.cc:79
G4double GetScreenRadius() const
G4Point3D GetPosition() const
G4double GetWorldRadius() const
G4Transform3D fObjectTransformation
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
const G4VisAttributes * fpVisAttribs
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4Colour & GetColour() const
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
const G4VisAttributes * GetVisAttributes() const
Qt3DCore::QTransform * CreateQTransformFrom(const G4Transform3D &)
Definition: G4Qt3DUtils.cc:33
QColor ConvertToQColor(const G4Colour &c)
Definition: G4Qt3DUtils.cc:46

◆ AddPrimitive() [2/13]

virtual void G4VSceneHandler::AddPrimitive ( const G4Circle )
virtual

Implements G4VSceneHandler.

◆ AddPrimitive() [3/13]

void G4VSceneHandler::AddPrimitive ( const G4Plotter )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 193 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}
const G4String & GetName() const
G4VGraphicsSystem & fSystem

◆ AddPrimitive() [4/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Polyhedron polyhedron)
virtual

Implements G4VSceneHandler.

Definition at line 696 of file G4Qt3DSceneHandler.cc.

697{
698 auto currentNode = CreateNewNode();
699 if (!currentNode) {
700 static G4bool first = true;
701 if (first) {
702 first = false;
703 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polyhedron&)",
704 "qt3d-0003", JustWarning,
705 "No available node!");
706 }
707 return;
708 }
709
710 if (polyhedron.GetNoFacets() == 0) return;
711
713
714 // Roll out vertices and normals for the faces. Note that this means vertices
715 // are duplicated. For example a box has 8 vertices, but to define 6 faces
716 // you need 12 triangles and 36 vertices. If it was just a matter of vertices
717 // we could restrict the number to 8 and use the indices to define the
718 // triangles, but we also have to consider the normals. A vertex can be have
719 // more than one normal, depending on which face it is being used to define.
720 // So we roll out all the vertices and normals for each triangle.
721 std::vector<G4Point3D> vertices;
722 std::vector<G4Normal3D> normals;
723
724 // Also roll out edges (as lines) for wireframe. Avoid duplicate lines,
725 // including those that differ only in the order of vertices.
726 typedef std::pair<G4Point3D,G4Point3D> Line;
727 std::vector<Line> lines;
728 auto insertIfNew = [&lines](const Line& newLine) {
729 // For a large polyhedron, eliminating lines like this is prohibitively
730 // expensive. Comment out for now, and maybe unwind altogether in future.
731 // Allow the graphics-reps utilities to optimise things like this.
732// for (const auto& line: lines) {
733// if ((newLine.first==line.first && newLine.second==line.second) ||
734// (newLine.first==line.second && newLine.second==line.first))
735// return;
736// }
737 lines.push_back(newLine);
738 };
739
740 G4bool isAuxilaryEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible();
741 G4bool notLastFace;
742 do {
743 G4int nEdges;
744 G4Point3D vertex [4];
745 G4int edgeFlag[4];
746 G4Normal3D normal [4];
747 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normal);
748 vertices.push_back(vertex[0]);
749 vertices.push_back(vertex[1]);
750 vertices.push_back(vertex[2]);
751 normals.push_back(normal[0]);
752 normals.push_back(normal[1]);
753 normals.push_back(normal[2]);
754 if(isAuxilaryEdgeVisible||edgeFlag[0]>0)insertIfNew(Line(vertex[0],vertex[1]));
755 if(isAuxilaryEdgeVisible||edgeFlag[1]>0)insertIfNew(Line(vertex[1],vertex[2]));
756 if (nEdges == 3) {
757 // Face is a triangle
758 // One more line for wireframe, triangles for surfaces are complete
759 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[0]));
760 } else if (nEdges == 4) {
761 // Face is a quadrilateral
762 // Create another triangle for surfaces, add two more lines for wireframe
763 vertices.push_back(vertex[2]);
764 vertices.push_back(vertex[3]);
765 vertices.push_back(vertex[0]);
766 normals.push_back(normal[2]);
767 normals.push_back(normal[3]);
768 normals.push_back(normal[0]);
769 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[3]));
770 if(isAuxilaryEdgeVisible||edgeFlag[3]>0)insertIfNew(Line(vertex[3],vertex[0]));
771 } else {
772 G4warn
773 << "ERROR: polyhedron face with unexpected number of edges (" << nEdges << ')'
774 << "\n Tag: " << fpModel->GetCurrentTag()
775 << G4endl;
776 return;
777 }
778 } while (notLastFace);
779 const auto nVerts = vertices.size();
780 const auto nLines = lines.size();
781
782 // Now put stuff into Qt objects
783
785 transform->setObjectName("transform");
786
787 Qt3DCore::QEntity* wireframeEntity = nullptr;
788 Qt3DCore::QEntity* surfaceEntity = nullptr;
789 static G4int errorCount = 0;
791 switch (drawing_style) {
793 wireframeEntity = new Qt3DCore::QEntity(currentNode);
794 wireframeEntity->addComponent(transform);
795 break;
797 wireframeEntity = new Qt3DCore::QEntity(currentNode);
798 wireframeEntity->addComponent(transform);
799 surfaceEntity = new Qt3DCore::QEntity(currentNode);
800 surfaceEntity->addComponent(transform);
801 break;
803 surfaceEntity = new Qt3DCore::QEntity(currentNode);
804 surfaceEntity->addComponent(transform);
805 break;
807 wireframeEntity = new Qt3DCore::QEntity(currentNode);
808 wireframeEntity->addComponent(transform);
809 surfaceEntity = new Qt3DCore::QEntity(currentNode);
810 surfaceEntity->addComponent(transform);
811 break;
813 // Shouldn't happen in this function (it's a polyhedron!)
814 if (errorCount == 0) {
815 ++errorCount;
816 G4warn << "WARNING: Qt3D: cloud drawing not implemented" << G4endl;
817 }
818 return;
819 break;
820 }
821
822 const auto vertexByteSize = 3*sizeof(PRECISION);
823
824 Qt3DRender::QGeometry* vertexGeometry = nullptr;
825 Qt3DRender::QGeometry* lineGeometry = nullptr;
826
827 Qt3DRender::QAttribute* positionAtt = nullptr;
828 Qt3DRender::QAttribute* normalAtt = nullptr;
829 Qt3DRender::QAttribute* lineAtt = nullptr;
830
831 Qt3DRender::QBuffer* vertexBuffer = nullptr;
832 if (drawing_style == G4ViewParameters::hlr ||
833 drawing_style == G4ViewParameters::hsr ||
834 drawing_style == G4ViewParameters::hlhsr) {
835
836 // Put vertices, normals into QByteArray
837 // Accomodates both vertices and normals - hence 2*
838 QByteArray vertexByteArray;
839 const auto vertexBufferByteSize = 2*nVerts*vertexByteSize;
840 vertexByteArray.resize((G4int)vertexBufferByteSize);
841 auto vertexBufferArray = reinterpret_cast<PRECISION*>(vertexByteArray.data());
842 G4int i1 = 0;
843 for (std::size_t i = 0; i < nVerts; ++i) {
844 vertexBufferArray[i1++] = vertices[i].x();
845 vertexBufferArray[i1++] = vertices[i].y();
846 vertexBufferArray[i1++] = vertices[i].z();
847 vertexBufferArray[i1++] = normals[i].x();
848 vertexBufferArray[i1++] = normals[i].y();
849 vertexBufferArray[i1++] = normals[i].z();
850 }
851 // Vertex buffer (vertices and normals)
852 vertexGeometry = new Qt3DRender::QGeometry();
853 vertexGeometry->setObjectName("vertexGeometry");
854 vertexBuffer = new Qt3DRender::QBuffer(vertexGeometry);
855 vertexBuffer->setObjectName("Vertex buffer");
856 vertexBuffer->setData(vertexByteArray);
857
858 // Position attribute
859 positionAtt = new Qt3DRender::QAttribute;
860 positionAtt->setObjectName("Position attribute");
861 positionAtt->setName(Qt3DRender::QAttribute::defaultPositionAttributeName());
862 positionAtt->setBuffer(vertexBuffer);
863 positionAtt->setAttributeType(Qt3DRender::QAttribute::VertexAttribute);
864 positionAtt->setVertexBaseType(BASETYPE);
865 positionAtt->setVertexSize(3);
866 positionAtt->setCount((G4int)nVerts);
867 positionAtt->setByteOffset(0);
868 positionAtt->setByteStride(2*vertexByteSize);
869
870 // Normal attribute
871 normalAtt = new Qt3DRender::QAttribute;
872 normalAtt->setObjectName("Normal attribute");
873 normalAtt->setName(Qt3DRender::QAttribute::defaultNormalAttributeName());
874 normalAtt->setBuffer(vertexBuffer);
875 normalAtt->setAttributeType(Qt3DRender::QAttribute::VertexAttribute);
876 normalAtt->setVertexBaseType(BASETYPE);
877 normalAtt->setVertexSize(3);
878 normalAtt->setCount((G4int)nVerts);
879 normalAtt->setByteOffset(vertexByteSize);
880 normalAtt->setByteStride(2*vertexByteSize);
881 }
882
883 Qt3DRender::QBuffer* lineBuffer = nullptr;
884 if (drawing_style == G4ViewParameters::wireframe ||
885 drawing_style == G4ViewParameters::hlr ||
886 drawing_style == G4ViewParameters::hlhsr) {
887
888 // Put lines into a QByteArray
889 QByteArray lineByteArray;
890 const auto lineBufferByteSize = 2*nLines*vertexByteSize;
891 lineByteArray.resize((G4int)lineBufferByteSize);
892 auto lineBufferArray = reinterpret_cast<PRECISION*>(lineByteArray.data());
893 G4int i2 = 0;
894 for (const auto& line: lines) {
895 lineBufferArray[i2++] = line.first.x();
896 lineBufferArray[i2++] = line.first.y();
897 lineBufferArray[i2++] = line.first.z();
898 lineBufferArray[i2++] = line.second.x();
899 lineBufferArray[i2++] = line.second.y();
900 lineBufferArray[i2++] = line.second.z();
901 }
902 // Line loop buffer
903 lineGeometry = new Qt3DRender::QGeometry();
904 lineGeometry->setObjectName("lineGeometry");
905 lineBuffer = new Qt3DRender::QBuffer(lineGeometry);
906 lineBuffer->setObjectName("Line buffer");
907 lineBuffer->setData(lineByteArray);
908
909 // Line attribute
910 lineAtt = new Qt3DRender::QAttribute;
911 lineAtt->setObjectName("Position attribute");
912 lineAtt->setName(Qt3DRender::QAttribute::defaultPositionAttributeName());
913 lineAtt->setBuffer(lineBuffer);
914 lineAtt->setAttributeType(Qt3DRender::QAttribute::VertexAttribute);
915 lineAtt->setVertexBaseType(BASETYPE);
916 lineAtt->setVertexSize(3);
917 lineAtt->setCount((G4int)nLines);
918 lineAtt->setByteOffset(0);
919 lineAtt->setByteStride(vertexByteSize);
920 }
921
922 // Create material and renderer(s)...
923
924 const auto& colour = fpVisAttribs->GetColour();
925 Qt3DExtras::QDiffuseSpecularMaterial* material;
926 Qt3DRender::QGeometryRenderer* renderer;
927 switch (drawing_style) {
928
930
931 lineGeometry->addAttribute(lineAtt);
932
933 material = new Qt3DExtras::QDiffuseSpecularMaterial();
934 material->setObjectName("materialForWireframe");
935 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
936 material->setShininess(0.);
937 material->setSpecular(0.);
938 wireframeEntity->addComponent(material);
939
940 renderer = new Qt3DRender::QGeometryRenderer;
941 renderer->setObjectName("polyhedronWireframeRenderer");
942 renderer->setGeometry(lineGeometry);
943 renderer->setVertexCount(2*(G4int)nLines);
944 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
945 wireframeEntity->addComponent(renderer);
946
947 break;
948
950
951 // Surfaces with background colour to hide the edges
952
953 vertexGeometry->addAttribute(positionAtt);
954 vertexGeometry->addAttribute(normalAtt);
955
956 material = new Qt3DExtras::QDiffuseSpecularMaterial();
957 material->setObjectName("materialForHiddenLines");
958 material->setAmbient(Qt::white); // White for now (should be from fVP)
959 material->setShininess(0.);
960 material->setSpecular(0.);
961 surfaceEntity->addComponent(material);
962
963 renderer = new Qt3DRender::QGeometryRenderer;
964 renderer->setObjectName("polyhedronSurfaceRenderer");
965 renderer->setGeometry(vertexGeometry);
966 renderer->setVertexCount((G4int)nVerts);
967 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
968 surfaceEntity->addComponent(renderer);
969
970 // Edges
971
972 lineGeometry->addAttribute(lineAtt);
973
974 material = new Qt3DExtras::QDiffuseSpecularMaterial();
975 material->setObjectName("materialForWireFrame");
976 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
977 material->setShininess(0.);
978 material->setSpecular(0.);
979 wireframeEntity->addComponent(material);
980
981 renderer = new Qt3DRender::QGeometryRenderer;
982 renderer->setObjectName("polyhedronWireframeRenderer");
983 renderer->setGeometry(lineGeometry);
984 renderer->setVertexCount(2*(G4int)nLines);
985 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
986 wireframeEntity->addComponent(renderer);
987
988 break;
989
991
992 vertexGeometry->addAttribute(positionAtt);
993 vertexGeometry->addAttribute(normalAtt);
994
995 material = new Qt3DExtras::QDiffuseSpecularMaterial();
996 material->setObjectName("materialForSurface");
997 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
998 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
999 surfaceEntity->addComponent(material);
1000
1001 renderer = new Qt3DRender::QGeometryRenderer;
1002 renderer->setObjectName("polyhedronSurfaceRenderer");
1003 renderer->setGeometry(vertexGeometry);
1004 renderer->setVertexCount((G4int)nVerts);
1005 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
1006 surfaceEntity->addComponent(renderer);
1007
1008 break;
1009
1011
1012 // Surfaces
1013
1014 vertexGeometry->addAttribute(positionAtt);
1015 vertexGeometry->addAttribute(normalAtt);
1016
1017 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1018 material->setObjectName("materialForSurface");
1019 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1020 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
1021 surfaceEntity->addComponent(material);
1022
1023 renderer = new Qt3DRender::QGeometryRenderer;
1024 renderer->setObjectName("polyhedronSurfaceRenderer");
1025 renderer->setGeometry(vertexGeometry);
1026 renderer->setVertexCount((G4int)nVerts);
1027 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
1028 surfaceEntity->addComponent(renderer);
1029
1030 // Edges
1031
1032 lineGeometry->addAttribute(lineAtt);
1033
1034 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1035 material->setObjectName("materialForWireframe");
1036 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1037 material->setShininess(0.);
1038 material->setSpecular(0.);
1039 wireframeEntity->addComponent(material);
1040
1041 renderer = new Qt3DRender::QGeometryRenderer;
1042 renderer->setObjectName("polyhedronSurfaceRenderer");
1043 renderer->setGeometry(lineGeometry);
1044 renderer->setVertexCount(2*(G4int)nLines);
1045 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
1046 wireframeEntity->addComponent(renderer);
1047
1048 break;
1049
1051 // Case trapped at start of function, so no need to implement
1052 break;
1053 }
1054}
#define PRECISION
#define BASETYPE
virtual G4String GetCurrentTag() const
Definition: G4VModel.cc:46
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4bool IsAuxEdgeVisible() const
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const

◆ AddPrimitive() [5/13]

virtual void G4VSceneHandler::AddPrimitive ( const G4Polyhedron )
virtual

Implements G4VSceneHandler.

◆ AddPrimitive() [6/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Polyline polyline)
virtual

Implements G4VSceneHandler.

Definition at line 258 of file G4Qt3DSceneHandler.cc.

259{
260#ifdef G4QT3DDEBUG
261 G4cout <<
262 "G4Qt3DSceneHandler::AddPrimitive(const G4Polyline& polyline) called.\n"
263 << polyline
264 << G4endl;
265#endif
266
267 if (polyline.size() == 0) return;
268
269 auto currentNode = CreateNewNode();
270 if (!currentNode) {
271 static G4bool first = true;
272 if (first) {
273 first = false;
274 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polyline&)",
275 "qt3d-0003", JustWarning,
276 "No available node!");
277 }
278 return;
279 }
280
282
284 transform->setObjectName("transform");
285
286 auto polylineEntity = new Qt3DCore::QEntity(currentNode);
287 polylineEntity->addComponent(transform);
288
289 const auto vertexByteSize = 3*sizeof(PRECISION);
290
291 const std::size_t nLines = polyline.size() - 1;
292 QByteArray polylineByteArray;
293 const auto polylineBufferByteSize = 2*nLines*vertexByteSize;
294 polylineByteArray.resize((G4int)polylineBufferByteSize);
295 auto polylineBufferArray = reinterpret_cast<PRECISION*>(polylineByteArray.data());
296 G4int iLine = 0;
297 for (std::size_t i = 0; i < nLines; ++i) {
298 polylineBufferArray[iLine++] = polyline[i].x();
299 polylineBufferArray[iLine++] = polyline[i].y();
300 polylineBufferArray[iLine++] = polyline[i].z();
301 polylineBufferArray[iLine++] = polyline[i+1].x();
302 polylineBufferArray[iLine++] = polyline[i+1].y();
303 polylineBufferArray[iLine++] = polyline[i+1].z();
304 }
305 auto polylineGeometry = new Qt3DRender::QGeometry();
306 polylineGeometry->setObjectName("polylineGeometry");
307 auto polylineBuffer = new Qt3DRender::QBuffer(polylineGeometry);
308 polylineBuffer->setObjectName("Polyline buffer");
309 polylineBuffer->setData(polylineByteArray);
310
311 auto polylineAtt = new Qt3DRender::QAttribute;
312 polylineAtt->setObjectName("Position attribute");
313 polylineAtt->setName(Qt3DRender::QAttribute::defaultPositionAttributeName());
314 polylineAtt->setBuffer(polylineBuffer);
315 polylineAtt->setAttributeType(Qt3DRender::QAttribute::VertexAttribute);
316 polylineAtt->setVertexBaseType(BASETYPE);
317 polylineAtt->setVertexSize(3);
318 polylineAtt->setCount((G4int)nLines);
319 polylineAtt->setByteOffset(0);
320 polylineAtt->setByteStride(vertexByteSize);
321
322 const auto& colour = fpVisAttribs->GetColour();
323
324 polylineGeometry->addAttribute(polylineAtt);
325
326 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
327 material->setObjectName("materialForPolyline");
328 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
329 material->setShininess(0.);
330 material->setSpecular(0.);
331 polylineEntity->addComponent(material);
332
333 auto renderer = new Qt3DRender::QGeometryRenderer;
334 renderer->setObjectName("polylineWireframeRenderer");
335 renderer->setGeometry(polylineGeometry);
336 renderer->setVertexCount(2*(G4int)nLines);
337 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
338 polylineEntity->addComponent(renderer);
339}

◆ AddPrimitive() [7/13]

virtual void G4VSceneHandler::AddPrimitive ( const G4Polyline )
virtual

Implements G4VSceneHandler.

◆ AddPrimitive() [8/13]

void G4VSceneHandler::AddPrimitive ( const G4Polymarker polymarker)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 191 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
void AddPrimitive(const G4Polyline &)

◆ AddPrimitive() [9/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Polymarker polymarker)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 341 of file G4Qt3DSceneHandler.cc.

342{
343 if (polymarker.size() == 0) return;
344
345 auto currentNode = CreateNewNode();
346 if (!currentNode) {
347 static G4bool first = true;
348 if (first) {
349 first = false;
350 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polymarker&)",
351 "qt3d-0003", JustWarning,
352 "No available node!");
353 }
354 return;
355 }
356
358
359 MarkerSizeType markerSizeType;
360 G4double markerSize = GetMarkerSize(polymarker, markerSizeType);
361
362 switch (polymarker.GetMarkerType()) {
363 default:
365 {
366 const std::size_t nDots = polymarker.size();
367
369 transform->setObjectName("transform");
370
371 auto polymarkerEntity = new Qt3DCore::QEntity(currentNode);
372 polymarkerEntity->addComponent(transform);
373
374 const auto vertexByteSize = 3*sizeof(PRECISION);
375
376 QByteArray polymarkerByteArray;
377 const auto polymarkerBufferByteSize = nDots*vertexByteSize;
378 polymarkerByteArray.resize((G4int)polymarkerBufferByteSize);
379 auto polymarkerBufferArray = reinterpret_cast<PRECISION*>(polymarkerByteArray.data());
380 G4int iMarker = 0;
381 for (std::size_t i = 0; i < polymarker.size(); ++i) {
382 polymarkerBufferArray[iMarker++] = polymarker[i].x();
383 polymarkerBufferArray[iMarker++] = polymarker[i].y();
384 polymarkerBufferArray[iMarker++] = polymarker[i].z();
385 }
386 auto polymarkerGeometry = new Qt3DRender::QGeometry();
387 polymarkerGeometry->setObjectName("polymarkerGeometry");
388 auto polymarkerBuffer = new Qt3DRender::QBuffer(polymarkerGeometry);
389 polymarkerBuffer->setObjectName("Polymarker buffer");
390 polymarkerBuffer->setData(polymarkerByteArray);
391
392 auto polymarkerAtt = new Qt3DRender::QAttribute;
393 polymarkerAtt->setObjectName("Position attribute");
394 polymarkerAtt->setName(Qt3DRender::QAttribute::defaultPositionAttributeName());
395 polymarkerAtt->setBuffer(polymarkerBuffer);
396 polymarkerAtt->setAttributeType(Qt3DRender::QAttribute::VertexAttribute);
397 polymarkerAtt->setVertexBaseType(BASETYPE);
398 polymarkerAtt->setVertexSize(3);
399 polymarkerAtt->setCount((G4int)nDots);
400 polymarkerAtt->setByteOffset(0);
401 polymarkerAtt->setByteStride(vertexByteSize);
402
403 const auto& colour = fpVisAttribs->GetColour();
404
405 polymarkerGeometry->addAttribute(polymarkerAtt);
406
407 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
408 material->setObjectName("materialForPolymarker");
409 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
410 material->setShininess(0.);
411 material->setSpecular(0.);
412 polymarkerEntity->addComponent(material);
413
414 auto renderer = new Qt3DRender::QGeometryRenderer;
415 renderer->setObjectName("polymarkerWireframeRenderer");
416 renderer->setGeometry(polymarkerGeometry);
417 renderer->setVertexCount((G4int)nDots);
418 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Points);
419 polymarkerEntity->addComponent(renderer);
420 }
421 break;
423 {
424 G4Circle circle (polymarker); // Default circle
425
426 const auto& colour = fpVisAttribs->GetColour();
427 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
428 material->setObjectName("materialForCircle");
429 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
430 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
431
432 auto sphereMesh = new Qt3DExtras::QSphereMesh;
433 sphereMesh->setObjectName("sphereMesh");
434 G4double radius = markerSize/2.;
435 if (markerSizeType == G4VSceneHandler::screen ) {
436 // Not figured out how to do screen-size, so use scene extent
437 const G4double scale = 200.; // Roughly pixels per scene
438 radius *= fpScene->GetExtent().GetExtentRadius()/scale;
439 }
440 sphereMesh->setRadius(radius);
441// sphereMesh->setInstanceCount(polymarker.size()); // Not undertood instancing yet
442
443// auto currentEntity = new Qt3DCore::QEntity(currentNode); // Not undertood instancing yet
444 for (std::size_t iPoint = 0; iPoint < polymarker.size(); ++iPoint) {
445 auto position = fObjectTransformation*G4Translate3D(polymarker[iPoint]);
447 auto currentEntity = new Qt3DCore::QEntity(currentNode); // Not undertood instancing yet
448 currentEntity->addComponent(material);
449 currentEntity->addComponent(transform);
450 currentEntity->addComponent(sphereMesh);
451 }
452 }
453 break;
455 {
456 G4Square square (polymarker); // Default square
457
458 const auto& colour = fpVisAttribs->GetColour();
459 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
460 material->setObjectName("materialForSquare");
461 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
462 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
463
464 auto boxMesh = new Qt3DExtras::QCuboidMesh();
465 boxMesh->setObjectName("boxMesh");
466 G4double side = markerSize;
467 if (markerSizeType == G4VSceneHandler::screen ) {
468 // Not figured out how to do screen-size, so use scene extent
469 const G4double scale = 200.; // Roughly pixles per scene
470 side *= fpScene->GetExtent().GetExtentRadius()/scale;
471 }
472 boxMesh->setXExtent(side);
473 boxMesh->setYExtent(side);
474 boxMesh->setZExtent(side);
475
476 for (std::size_t iPoint = 0; iPoint < polymarker.size(); ++iPoint) {
477 auto position = fObjectTransformation*G4Translate3D(polymarker[iPoint]);
479 auto currentEntity = new Qt3DCore::QEntity(currentNode);
480 currentEntity->addComponent(material);
481 currentEntity->addComponent(transform);
482 currentEntity->addComponent(boxMesh);
483 }
484 }
485 break;
486 }
487}

◆ AddPrimitive() [10/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Square square)
virtual

Implements G4VSceneHandler.

Definition at line 627 of file G4Qt3DSceneHandler.cc.

628{
629#ifdef G4QT3DDEBUG
630 G4cout <<
631 "G4Qt3DSceneHandler::AddPrimitive(const G4Square& square) called.\n"
632 << square
633 << G4endl;
634#endif
635
636#ifdef G4QT3DDEBUG
637 MarkerSizeType sizeType;
638 G4double size = GetMarkerSize (square, sizeType);
639 switch (sizeType) {
640 default:
641 case screen:
642 // Draw in screen coordinates.
643 G4cout << "screen";
644 break;
645 case world:
646 // Draw in world coordinates.
647 G4cout << "world";
648 break;
649 }
650 G4cout << " size: " << size << G4endl;
651#endif
652
653 auto currentNode = CreateNewNode();
654 if (!currentNode) {
655 static G4bool first = true;
656 if (first) {
657 first = false;
658 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Square&)",
659 "qt3d-0003", JustWarning,
660 "No available node!");
661 }
662 return;
663 }
664
666
669
670 const auto& colour = fpVisAttribs->GetColour();
671 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
672 material->setObjectName("materialForSquare");
673 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
674 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
675
676 auto boxMesh = new Qt3DExtras::QCuboidMesh();
677 boxMesh->setObjectName("boxMesh");
678 G4double side;
679 if (square.GetSizeType() == G4VMarker::world ) {
680 side = square.GetWorldDiameter();
681 } else { // Screen-size or none
682 // Not figured out how to do screen-size, so use scene extent
683 const G4double scale = 200.; // Roughly pixles per scene
684 side = square.GetScreenDiameter()*fpScene->GetExtent().GetExtentRadius()/scale;
685 }
686 boxMesh->setXExtent(side);
687 boxMesh->setYExtent(side);
688 boxMesh->setZExtent(side);
689
690 auto currentEntity = new Qt3DCore::QEntity(currentNode);
691 currentEntity->addComponent(material);
692 currentEntity->addComponent(transform);
693 currentEntity->addComponent(boxMesh);
694}
G4double GetWorldDiameter() const
G4double GetScreenDiameter() const

◆ AddPrimitive() [11/13]

virtual void G4VSceneHandler::AddPrimitive ( const G4Square )
virtual

Implements G4VSceneHandler.

◆ AddPrimitive() [12/13]

void G4Qt3DSceneHandler::AddPrimitive ( const G4Text )
virtual

Implements G4VSceneHandler.

Definition at line 496 of file G4Qt3DSceneHandler.cc.

496 {
497#endif
498
499 static G4bool first = true;
500 if (first) {
501 first = false;
502 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Text&)",
503 "qt3D-0002", JustWarning,
504 "Text drawing doesn't work yet");
505 } // OK. Not working, but let it execute, which it does without error.
506
507 /* But it crashes after /vis/viewer/rebuild!!!
508 auto currentNode = CreateNewNode();
509 if (!currentNode) {
510 static G4bool first = true;
511 if (first) {
512 first = false;
513 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Text&)",
514 "qt3d-0003", JustWarning,
515 "No available node!");
516 }
517 return;
518 }
519
520 fpVisAttribs = fpViewer->GetApplicableVisAttributes(text.GetVisAttributes());
521
522 auto position = fObjectTransformation*G4Translate3D(text.GetPosition());
523 auto transform = G4Qt3DUtils::CreateQTransformFrom(position);
524// transform->setScale(10);
525 transform->setScale(0.1);
526
527// auto currentEntity = new Qt3DCore::QEntity(currentNode);
528
529 // This simply does not work
530 auto qtext = new Qt3DExtras::QText2DEntity();
531 qtext->setParent(currentNode);
532// qtext->setParent(currentEntity); // ?? Doesn't help
533 qtext->setText(text.GetText().c_str());
534// qtext->setHeight(100);
535// qtext->setWidth(1000);
536 qtext->setHeight(20);
537 qtext->setWidth(100);
538 qtext->setColor(Qt::green);
539 qtext->setFont(QFont("Courier New", 10));
540 qtext->addComponent(transform);
541
542 // This produces text in 3D facing +z - not what we want
543// const auto& colour = GetTextColour(text);
544// auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
545// material->setObjectName("materialForText");
546// material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
547// if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
548//
549// auto textMesh = new Qt3DExtras::QExtrudedTextMesh();
550// textMesh->setText(text.GetText().c_str());
551// textMesh->setFont(QFont("Courier New", 10));
552// textMesh->setDepth(.01f);
553//
554// currentNode->addComponent(material);
555// currentNode->addComponent(transform);
556// currentNode->addComponent(textMesh);
557 */
558}

◆ AddPrimitive() [13/13]

virtual void G4VSceneHandler::AddPrimitive ( const G4Text )
virtual

Implements G4VSceneHandler.

◆ BeginPrimitives()

void G4Qt3DSceneHandler::BeginPrimitives ( const G4Transform3D objectTransformation)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 247 of file G4Qt3DSceneHandler.cc.

249{
250 G4VSceneHandler::BeginPrimitives(objectTransformation);
251}
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())

◆ BeginPrimitives2D()

void G4Qt3DSceneHandler::BeginPrimitives2D ( const G4Transform3D objectTransformation)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 220 of file G4Qt3DSceneHandler.cc.

221{
222// The x,y coordinates of the primitives passed to AddPrimitive are
223// intrepreted as screen coordinates, -1 < x,y < 1. The
224// z-coordinate is ignored.
225// IMPORTANT: invoke this from your polymorphic versions, e.g.:
226// void MyXXXSceneHandler::BeginPrimitives2D
227// (const G4Transform3D& objectTransformation) {
228 static G4bool first = true;
229 if (first) {
230 first = false;
231 G4Exception("G4Qt3DSceneHandler::BeginPrimitives2D", "qt3D-0001",
233 "2D drawing not yet implemented");
234 }
235 G4VSceneHandler::BeginPrimitives2D (objectTransformation);
236// ...
237}
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation=G4Transform3D())

◆ ClearStore()

void G4Qt3DSceneHandler::ClearStore ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 1061 of file G4Qt3DSceneHandler.cc.

1062{
1065}
void delete_components_and_children_of_entity_recursively(Qt3DCore::QNode *node)
Definition: G4Qt3DUtils.cc:102

◆ ClearTransientStore()

void G4Qt3DSceneHandler::ClearTransientStore ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 1067 of file G4Qt3DSceneHandler.cc.

◆ CreateNewNode()

G4Qt3DQEntity * G4Qt3DSceneHandler::CreateNewNode ( )
protected

Definition at line 122 of file G4Qt3DSceneHandler.cc.

123{
124 // Create a G4Qt3DQEntity node suitable for next solid or primitive
125
126 G4Qt3DQEntity* newNode = nullptr;
127
128 if (fReadyForTransients) { // All transients hang from this node
129 newNode = new G4Qt3DQEntity(fpTransientObjects);
131 newNode->setObjectName(name.c_str());
132 return newNode;
133 }
134
135 G4PhysicalVolumeModel* pPVModel =
136 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
137
138 if (!pPVModel) { // Persistent objects (e.g., axes)
139 newNode = new G4Qt3DQEntity(fpPersistentObjects);
140 newNode->setObjectName(fpModel->GetGlobalTag().c_str());
141 return newNode;
142 }
143
144 // So this is a G4PhysicalVolumeModel
145
147 typedef std::vector<PVNodeID> PVPath;
148// const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
149 const PVPath& fullPVPath = pPVModel->GetFullPVPath();
150 //G4int currentDepth = pPVModel->GetCurrentDepth();
151 //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
152 //G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
153 //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
154 // Note: pCurrentMaterial may be zero (parallel world).
155
156#ifdef G4QTDEBUG
157 G4cout << "A: " << fullPVPath << G4endl; // DEBUG
158#endif
159
160 // Find appropriate root
161 const std::size_t nWorlds = fpPhysicalVolumeObjects.size();
162 std::size_t iWorld = 0;
163 for (; iWorld < nWorlds; ++iWorld) {
164 if (fullPVPath[0].GetPhysicalVolume() ==
165 fpPhysicalVolumeObjects[iWorld]->GetPVNodeID().GetPhysicalVolume()) break;
166 }
167 if (iWorld == nWorlds) {
168 G4Exception("G4Qt3DSceneHandler::CreateNewNode", "qt3D-0000", FatalException,
169 "World mis-match - not possible(!?)");
170 }
171
172 // (Re-)establish pv path of root entity
174 wrld->SetPVNodeID(fullPVPath[0]);
175
176 // Create nodes as required
177 G4Qt3DQEntity* node = wrld;
178 newNode = node;
179 const std::size_t depth = fullPVPath.size();
180 std::size_t iDepth = 1;
181 while (iDepth < depth) {
182 const auto& children = node->children();
183 const G4int nChildren = children.size(); // int size() (Qt covention?)
184 G4int iChild = 0;
185 G4Qt3DQEntity* child = nullptr;
186 for (; iChild < nChildren; ++iChild) {
187 child = static_cast<G4Qt3DQEntity*>(children[iChild]);
188 if (child->GetPVNodeID() == fullPVPath[iDepth]) break;
189 }
190 if (iChild != nChildren) { // Existing node found
191 node = child; // Must be the ancestor of new node (subsequent iteration)
192 } else {
193 // Add a new node as child of node
194 newNode = new G4Qt3DQEntity(node);
195 newNode->SetPVNodeID(fullPVPath[iDepth]);
196 std::ostringstream oss;
197 oss << newNode->GetPVNodeID().GetPhysicalVolume()->GetName()
198 << ':' << newNode->GetPVNodeID().GetCopyNo();
199 newNode->setObjectName(oss.str().c_str());
200 node = newNode;
201 }
202 ++iDepth;
203 }
204
205 return node;
206}
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
void SetPVNodeID(const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID &id)
const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID & GetPVNodeID() const
std::vector< G4Qt3DQEntity * > fpPhysicalVolumeObjects
Qt3DCore::QEntity * fpPersistentObjects
const G4String & GetGlobalTag() const
const char * name(G4int ptype)

Referenced by AddPrimitive().

◆ EndPrimitives()

void G4Qt3DSceneHandler::EndPrimitives ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 253 of file G4Qt3DSceneHandler.cc.

254{
256}
virtual void EndPrimitives()

◆ EndPrimitives2D()

void G4Qt3DSceneHandler::EndPrimitives2D ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 239 of file G4Qt3DSceneHandler.cc.

240{
241// IMPORTANT: invoke this from your polymorphic versions, e.g.:
242// void MyXXXSceneHandler::EndPrimitives2D () {
243// ...
245}
virtual void EndPrimitives2D()

◆ EstablishG4Qt3DQEntities()

void G4Qt3DSceneHandler::EstablishG4Qt3DQEntities ( )
protected

Definition at line 99 of file G4Qt3DSceneHandler.cc.

100{
101 fpTransientObjects = new G4Qt3DQEntity(fpQt3DScene); // Hangs from root
102 fpTransientObjects ->setObjectName("G4Qt3DTORoot");
103 fpPersistentObjects = new G4Qt3DQEntity(fpQt3DScene); // Hangs from root
104 fpPersistentObjects ->setObjectName("G4Qt3DPORoot");
105
106 // Physical volume objects for each world hang from POs
107 G4TransportationManager* transportationManager
109 std::size_t nWorlds = transportationManager->GetNoWorlds();
110 std::vector<G4VPhysicalVolume*>::iterator iterWorld
111 = transportationManager->GetWorldsIterator();
112 fpPhysicalVolumeObjects.resize(nWorlds);
113 for (std::size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
114 G4VPhysicalVolume* wrld = (*iterWorld);
115 auto entity = new G4Qt3DQEntity(fpPersistentObjects);
116 entity->setObjectName("G4Qt3DPORoot_"+QString(wrld->GetName()));
117 entity->SetPVNodeID(G4PhysicalVolumeModel::G4PhysicalVolumeNodeID(wrld));
118 fpPhysicalVolumeObjects[i] = entity;
119 }
120}
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const

Referenced by ClearStore(), and G4Qt3DSceneHandler().

◆ PostAddSolid()

void G4Qt3DSceneHandler::PostAddSolid ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 215 of file G4Qt3DSceneHandler.cc.

216{
218}
virtual void PostAddSolid()

◆ PreAddSolid()

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

Reimplemented from G4VSceneHandler.

Definition at line 208 of file G4Qt3DSceneHandler.cc.

211{
212 G4VSceneHandler::PreAddSolid(objectTransformation, visAttribs);
213}
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)

Friends And Related Function Documentation

◆ G4Qt3DViewer

friend class G4Qt3DViewer
friend

Definition at line 42 of file G4Qt3DSceneHandler.hh.

Member Data Documentation

◆ fpPersistentObjects

Qt3DCore::QEntity* G4Qt3DSceneHandler::fpPersistentObjects
protected

Definition at line 83 of file G4Qt3DSceneHandler.hh.

Referenced by CreateNewNode(), and EstablishG4Qt3DQEntities().

◆ fpPhysicalVolumeObjects

std::vector<G4Qt3DQEntity*> G4Qt3DSceneHandler::fpPhysicalVolumeObjects
protected

Definition at line 84 of file G4Qt3DSceneHandler.hh.

Referenced by CreateNewNode(), and EstablishG4Qt3DQEntities().

◆ fpQt3DScene

Qt3DCore::QEntity* G4Qt3DSceneHandler::fpQt3DScene
protected

◆ fpTransientObjects

Qt3DCore::QEntity* G4Qt3DSceneHandler::fpTransientObjects
protected

◆ fSceneIdCount

G4int G4Qt3DSceneHandler::fSceneIdCount = 0
staticprotected

Definition at line 79 of file G4Qt3DSceneHandler.hh.


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