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

#include <G4ToolsSGSceneHandler.hh>

Inheritance diagram for G4ToolsSGSceneHandler:

Classes

class  Messenger

Public Member Functions

virtual void AddPrimitive (const G4Polyline &)
virtual void AddPrimitive (const G4Text &)
virtual void AddPrimitive (const G4Circle &)
virtual void AddPrimitive (const G4Square &)
virtual void AddPrimitive (const G4Polymarker &)
virtual void AddPrimitive (const G4Polyhedron &)
virtual void AddPrimitive (const G4Plotter &)
virtual void AddCompound (const G4Mesh &)
virtual void ClearStore ()
virtual void ClearTransientStore ()
 G4ToolsSGSceneHandler (G4VGraphicsSystem &system, const G4String &name)
virtual ~G4ToolsSGSceneHandler ()
tools::sg::separator & GetTransient2DObjects ()
tools::sg::separator & GetPersistent2DObjects ()
tools::sg::separator & GetTransient3DObjects ()
tools::sg::separator & GetPersistent3DObjects ()
void TouchPlotters (tools::sg::node &)
void SetMarkerScale (double)
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 > &)
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 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 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 ()
void AddViewerToList (G4VViewer *pView)
void RemoveViewerFromList (G4VViewer *pView)
Public Member Functions inherited from G4VGraphicsScene
 G4VGraphicsScene ()
virtual ~G4VGraphicsScene ()
G4int GetMaxGeometryDepth () const
void SetMaxGeometryDepth (G4int maxDepth)

Protected Types

using Region_h1 = std::pair<unsigned int,int>
using Region_h2 = std::pair<unsigned int,int>

Protected Member Functions

 G4ToolsSGSceneHandler (const G4ToolsSGSceneHandler &)
G4ToolsSGSceneHandleroperator= (const G4ToolsSGSceneHandler &)
void CreateSG ()
void EstablishBaseNodes ()
tools::sg::separator * GetOrCreateNode ()
void SetPlotterHistograms (tools::sg::plots &)
Protected Member Functions inherited from G4VSceneHandler
virtual void ProcessScene ()
virtual void ProcessTransients ()
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

tools::sg::separator fpTransient2DObjects
tools::sg::separator fpPersistent2DObjects
tools::sg::separator fpTransient3DObjects
tools::sg::separator fpPersistent3DObjects
std::vector< G4ToolsSGNode * > fpPhysicalVolumeObjects
tools::sg::base_freetype * fFreetypeNode
double fMarkerScale
std::vector< Region_h1fRegionH1s
std::vector< Region_h2fRegionH2s
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
std::map< G4VPhysicalVolume *, G4StringfProblematicVolumes
Protected Attributes inherited from G4VGraphicsScene
G4int fMaxGeometryDepth = 0

Static Protected Attributes

static G4int fSceneIdCount = 0

Additional Inherited Members

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

Detailed Description

Definition at line 44 of file G4ToolsSGSceneHandler.hh.

Member Typedef Documentation

◆ Region_h1

using G4ToolsSGSceneHandler::Region_h1 = std::pair<unsigned int,int>
protected

Definition at line 97 of file G4ToolsSGSceneHandler.hh.

◆ Region_h2

using G4ToolsSGSceneHandler::Region_h2 = std::pair<unsigned int,int>
protected

Definition at line 98 of file G4ToolsSGSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4ToolsSGSceneHandler() [1/2]

G4ToolsSGSceneHandler::G4ToolsSGSceneHandler ( G4VGraphicsSystem & system,
const G4String & name )

Definition at line 76 of file G4ToolsSGSceneHandler.cc.

78:parent(system, fSceneIdCount++, name)
81{
82 //::printf("debug : G4ToolsSGSceneHandler : %lu, %s\n",this,name.c_str());
84#if defined(TOOLS_USE_FREETYPE)
85 fFreetypeNode = new toolx::sg::text_freetype();
86 fFreetypeNode->add_embedded_font(tools::sg::font_lato_regular_ttf(),tools::font::lato_regular_ttf);
87 fFreetypeNode->add_embedded_font(tools::sg::font_roboto_bold_ttf(),tools::font::roboto_bold_ttf);
88#else
89 fFreetypeNode = new tools::sg::dummy_freetype();
90#endif
92}
tools::sg::base_freetype * fFreetypeNode

Referenced by G4ToolsSGSceneHandler(), operator=(), and G4ToolsSGSceneHandler::Messenger::SetNewValue().

◆ ~G4ToolsSGSceneHandler()

G4ToolsSGSceneHandler::~G4ToolsSGSceneHandler ( )
virtual

Definition at line 94 of file G4ToolsSGSceneHandler.cc.

95{
96 //::printf("debug : ~G4ToolsSGSceneHandler : %lu\n",this);
97 //WARNING : nodes may refer graphics managers (as tools/sg/[GL_manager,gl2ps_manager,zb_manager]
98 // used by viewers) to handle gstos (for GPU) or textures, then we have to delete them first.
99 // It is assumed that we pass here BEFORE the attached/managed viewers are deleted.
100 fpTransient2DObjects.clear();
101 fpPersistent2DObjects.clear();
102 fpTransient3DObjects.clear();
103 fpPersistent3DObjects.clear();
104 delete fFreetypeNode;
105}
tools::sg::separator fpPersistent2DObjects
tools::sg::separator fpTransient3DObjects
tools::sg::separator fpTransient2DObjects
tools::sg::separator fpPersistent3DObjects

◆ G4ToolsSGSceneHandler() [2/2]

G4ToolsSGSceneHandler::G4ToolsSGSceneHandler ( const G4ToolsSGSceneHandler & )
protected

Member Function Documentation

◆ AddCompound() [1/6]

void G4ToolsSGSceneHandler::AddCompound ( const G4Mesh & mesh)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 593 of file G4ToolsSGSceneHandler.cc.

594{
596}
void StandardSpecialMeshRendering(const G4Mesh &)

◆ AddCompound() [2/6]

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

Reimplemented from G4VSceneHandler.

Definition at line 133 of file G4VSceneHandler.cc.

345 {
346 using MeshScoreMap = G4VScoringMesh::MeshScoreMap;
347 //G4cout << "AddCompound: hits: " << &hits << G4endl;
348 G4bool scoreMapHits = false;
349 G4ScoringManager* scoringManager = G4ScoringManager::GetScoringManagerIfExist();
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(const auto& i : scoreMap) {
358 const G4String& scoreMapName = i.first;
359 if (scoreMapName == mapNam) {
360 G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
361 scoreMapHits = true;
362 mesh->DrawMesh(scoreMapName, &colorMap);
363 }
364 }
365 }
366 }
367 }
368 if (scoreMapHits) {
369 static G4bool first = true;
370 if (first) {
371 first = false;
372 G4cout <<
373 "Scoring map drawn with default parameters."
374 "\n To get gMocren file for gMocren browser:"
375 "\n /vis/open gMocrenFile"
376 "\n /vis/viewer/flush"
377 "\n Many other options available with /score/draw... commands."
378 "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
379 << G4endl;
380 }
381 } else { // Not score map hits. Just call DrawAllHits.
382 // Cast away const because DrawAllHits is non-const!!!!
383 const_cast<G4THitsMap<G4double>&>(hits).DrawAllHits();
384 }
385}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VScoringMesh * GetMesh(G4int i) const
std::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

Reimplemented from G4VSceneHandler.

Definition at line 134 of file G4VSceneHandler.cc.

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

◆ AddCompound() [4/6]

void G4VSceneHandler::AddCompound ( const G4VDigi & digi)
virtual

Reimplemented from G4VSceneHandler.

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

Reimplemented from G4VSceneHandler.

Definition at line 131 of file G4VSceneHandler.cc.

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

◆ AddCompound() [6/6]

void G4VSceneHandler::AddCompound ( const G4VTrajectory & traj)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 130 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)
virtual void DrawTrajectory() const

◆ AddPrimitive() [1/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Circle & a_circle)
virtual

Implements G4VSceneHandler.

Definition at line 421 of file G4ToolsSGSceneHandler.cc.

422{
423 G4Polymarker oneCircle(a_circle);
424 oneCircle.push_back(a_circle.GetPosition());
425 oneCircle.SetMarkerType(G4Polymarker::circles);
426 // Call this AddPrimitive to avoid re-doing sub-class code.
428}
virtual void AddPrimitive(const G4Polyline &)
G4Point3D GetPosition() const

◆ AddPrimitive() [2/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Plotter & a_plotter)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 820 of file G4ToolsSGSceneHandler.cc.

821{
822//G4cout << "debug : G4ToolsSGSceneHandler::AddPrimitive : 004" << std::endl;
823 if(!fpViewer) return;
824
825 auto currentNode = GetOrCreateNode();
826 if (!currentNode) return; // Node not available
827
828 currentNode->add(new tools::sg::light_off());
829
830 tools::sg::plots* _plots = new tools::sg::plots(*fFreetypeNode);
831 currentNode->add(_plots);
832
833 _plots->view_border = false;
834 _plots->set_regions(a_plotter.GetColumns(),a_plotter.GetRows());
835
836 {tools::sg::event_dispatcher* dpt = new tools::sg::event_dispatcher;
837 dpt->add_callback(new plots_cbk(*_plots));
838 currentNode->add(dpt);}
839
840 SetPlotterStyles(*_plots,a_plotter.GetStyles(),a_plotter.GetRegionStyles());
841
842 tools::sg::cmaps_t _cmaps;
843 _cmaps["default"] = tools::sg::style_default_colormap();
844 _cmaps["ROOT"] = tools::sg::style_ROOT_colormap();
845
846 SetPlotterParameters(_cmaps,*_plots,a_plotter.GetRegionParameters());
847
848 fRegionH1s = a_plotter.GetRegionH1s();
849 fRegionH2s = a_plotter.GetRegionH2s();
850
851 SetPlotterHistograms(*_plots);
852}
void SetPlotterParameters(tools::sg::cmaps_t &a_cmaps, tools::sg::plots &a_plots, const std::vector< G4Plotter::RegionParameter > &a_region_parameters)
void SetPlotterStyles(tools::sg::plots &a_plots, const std::vector< G4String > &a_plotter_styles, const std::vector< G4Plotter::RegionStyle > &a_region_styles)
unsigned int GetRows() const
Definition G4Plotter.hh:68
const std::vector< RegionParameter > & GetRegionParameters() const
Definition G4Plotter.hh:71
const std::vector< G4String > & GetStyles() const
Definition G4Plotter.hh:69
const std::vector< Region_h2 > & GetRegionH2s() const
Definition G4Plotter.hh:76
const std::vector< RegionStyle > & GetRegionStyles() const
Definition G4Plotter.hh:70
const std::vector< Region_h1 > & GetRegionH1s() const
Definition G4Plotter.hh:75
unsigned int GetColumns() const
Definition G4Plotter.hh:67
void SetPlotterHistograms(tools::sg::plots &)
std::vector< Region_h1 > fRegionH1s
tools::sg::separator * GetOrCreateNode()
std::vector< Region_h2 > fRegionH2s

◆ AddPrimitive() [3/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Polyhedron & a_polyhedron)
virtual

Implements G4VSceneHandler.

Definition at line 439 of file G4ToolsSGSceneHandler.cc.

440{
441 if (a_polyhedron.GetNoFacets() == 0) return;
442
443 //::printf("debug : G4ToolsSGSceneHandler::AddPrimitive(const G4Polyhedron&) : %d\n",a_polyhedron.GetNoFacets());
444
445 fpVisAttribs = fpViewer->GetApplicableVisAttributes(a_polyhedron.GetVisAttributes());
446
447 // Roll out vertices and normals for the faces. Note that this means vertices
448 // are duplicated. For example a box has 8 vertices, but to define 6 faces
449 // you need 12 triangles and 36 vertices. If it was just a matter of vertices
450 // we could restrict the number to 8 and use the indices to define the
451 // triangles, but we also have to consider the normals. A vertex can be have
452 // more than one normal, depending on which face it is being used to define.
453 // So we roll out all the vertices and normals for each triangle.
454 std::vector<G4Point3D> vertices;
455 std::vector<G4Normal3D> normals;
456
457 // Also roll out edges (as lines) for wireframe. Avoid duplicate lines,
458 // including those that differ only in the order of vertices.
459 typedef std::pair<G4Point3D,G4Point3D> Line;
460 std::vector<Line> lines;
461 auto insertIfNew = [&lines](const Line& newLine) {
462// for (const auto& line: lines) {
463// if ((newLine.first==line.first && newLine.second==line.second) ||
464// (newLine.first==line.second && newLine.second==line.first))
465// return;
466// }
467 lines.push_back(newLine);
468 };
469
470 G4bool isAuxilaryEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible();
471 G4bool notLastFace;
472 do {
473 G4int nEdges;
474 G4Point3D vertex [4];
475 G4int edgeFlag[4];
476 G4Normal3D normal [4];
477 notLastFace = a_polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normal);
478 vertices.push_back(vertex[0]);
479 vertices.push_back(vertex[1]);
480 vertices.push_back(vertex[2]);
481 normals.push_back(normal[0]);
482 normals.push_back(normal[1]);
483 normals.push_back(normal[2]);
484 if(isAuxilaryEdgeVisible||edgeFlag[0]>0)insertIfNew(Line(vertex[0],vertex[1]));
485 if(isAuxilaryEdgeVisible||edgeFlag[1]>0)insertIfNew(Line(vertex[1],vertex[2]));
486 if (nEdges == 3) {
487 // Face is a triangle
488 // One more line for wireframe, triangles for surfaces are complete
489 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[0]));
490 } else if (nEdges == 4) {
491 // Face is a quadrilateral
492 // Create another triangle for surfaces, add two more lines for wireframe
493 vertices.push_back(vertex[2]);
494 vertices.push_back(vertex[3]);
495 vertices.push_back(vertex[0]);
496 normals.push_back(normal[2]);
497 normals.push_back(normal[3]);
498 normals.push_back(normal[0]);
499 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[3]));
500 if(isAuxilaryEdgeVisible||edgeFlag[3]>0)insertIfNew(Line(vertex[3],vertex[0]));
501 } else {
502 G4cerr
503 << "G4ToolsSGSceneHandler::AddPrimitive(G4Polyhedron): WARNING:"
504 << "\n G4Polyhedron facet with unexpected number of edges (" << nEdges << ")."
505 << G4endl;
506 return;
507 }
508 } while (notLastFace);
509
511 switch (drawing_style) {
513 //vertices.clear();
514 break;
516 break;
518 //lines.clear();
519 break;
521 break;
523 // Shouldn't happen in this function (it's a polyhedron!) - ignore
524 return;
525 }
526
527 auto currentNode = GetOrCreateNode();
528 if (!currentNode) return; // Node not available
529
530 tools::sg::separator* sep = new tools::sg::separator;
531 currentNode->add(sep);
532
533 // Transformation
534 {tools::sg::matrix* mtx = new tools::sg::matrix;
536 mtx->mtx.value().set_matrix(elem(0,0),elem(0,1),elem(0,2),elem(0,3),
537 elem(1,0),elem(1,1),elem(1,2),elem(1,3),
538 elem(2,0),elem(2,1),elem(2,2),elem(2,3),
539 0, 0, 0, 1);
540 sep->add(mtx);}
541
542 {const auto& colour = GetColour(a_polyhedron);
543 tools::sg::rgba* mat = new tools::sg::rgba();
544 mat->color =
545 tools::colorf(float(colour.GetRed()),
546 float(colour.GetGreen()),
547 float(colour.GetBlue()),
548 float(colour.GetAlpha()));
549 sep->add(mat);}
550
551 if (drawing_style == G4ViewParameters::hlr ||
552 drawing_style == G4ViewParameters::hsr ||
553 drawing_style == G4ViewParameters::hlhsr) {
554
555 {tools::sg::draw_style* ds = new tools::sg::draw_style;
556 ds->style = tools::sg::draw_filled;
557 //ds->cull_face = true;
558 sep->add(ds);}
559
560 tools::sg::atb_vertices* vtxs = new tools::sg::atb_vertices;
561 vtxs->mode = tools::gl::triangles();
562 sep->add(vtxs);
563
564 const auto nVerts = vertices.size();
565 for (size_t i = 0; i < nVerts; i++) {
566 vtxs->add(float(vertices[i].x()),float(vertices[i].y()),float(vertices[i].z()));
567 vtxs->add_normal(float(normals[i].x()),float(normals[i].y()),float(normals[i].z()));
568 }
569 }
570
571 if (drawing_style == G4ViewParameters::wireframe ||
572 drawing_style == G4ViewParameters::hlr ||
573 drawing_style == G4ViewParameters::hlhsr) {
574
575 {tools::sg::draw_style* ds = new tools::sg::draw_style;
576 ds->style = tools::sg::draw_lines;
577 auto visAtts = GetCurrentViewer()->GetApplicableVisAttributes(a_polyhedron.GetVisAttributes());
578 ds->line_width = GetLineWidth(visAtts); // Multiplies by GetGlobalLineWidthScale
579 sep->add(ds);}
580
581 tools::sg::vertices* vtxs = new tools::sg::vertices;
582 vtxs->mode = tools::gl::lines(); //segments
583 sep->add(vtxs);
584
585 for (const auto& line: lines) {
586 vtxs->add(float(line.first.x()),float(line.first.y()),float(line.first.z()));
587 vtxs->add(float(line.second.x()),float(line.second.y()),float(line.second.z()));
588 }
589
590 }
591}
HepGeom::Normal3D< G4double > G4Normal3D
Definition G4Normal3D.hh:34
#define elem(i, j)
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
HepGeom::Transform3D G4Transform3D
G4GLOB_DLL std::ostream G4cerr
const G4Colour & GetColour()
G4Transform3D fObjectTransformation
G4VViewer * GetCurrentViewer() const
const G4VisAttributes * fpVisAttribs
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4double GetLineWidth(const G4VisAttributes *)
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4VisAttributes * GetVisAttributes() const
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const

◆ AddPrimitive() [4/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Polyline & a_polyline)
virtual

Implements G4VSceneHandler.

Definition at line 215 of file G4ToolsSGSceneHandler.cc.

216{
217 //G4cout << "debug : G4ToolsSGSceneHandler::AddPrimitive(const G4Polyline&) : \n" << a_polyline << G4endl;
218 if (a_polyline.size() == 0) return;
219
220 tools::sg::separator* parentNode = 0;
221 if(fProcessing2D) {
222 parentNode = new tools::sg::separator;
224 fpTransient2DObjects.add(parentNode);
225 } else {
226 fpPersistent2DObjects.add(parentNode);
227 }
228
229 } else {
230 parentNode = GetOrCreateNode();
231 if(!parentNode) return;
232
233 tools::sg::matrix* mtx = new tools::sg::matrix;
235 mtx->mtx.value().set_matrix(elem(0,0),elem(0,1),elem(0,2),elem(0,3),
236 elem(1,0),elem(1,1),elem(1,2),elem(1,3),
237 elem(2,0),elem(2,1),elem(2,2),elem(2,3),
238 0, 0, 0, 1);
239 parentNode->add(mtx);
240 }
241
242 {const auto& colour = GetColour(a_polyline);
243 tools::sg::rgba* mat = new tools::sg::rgba();
244 mat->color =
245 tools::colorf(float(colour.GetRed()),
246 float(colour.GetGreen()),
247 float(colour.GetBlue()),
248 float(colour.GetAlpha()));
249 parentNode->add(mat);}
250
251 {tools::sg::draw_style* ds = new tools::sg::draw_style;
252 ds->style = tools::sg::draw_lines;
253 auto visAtts = GetCurrentViewer()->GetApplicableVisAttributes(a_polyline.GetVisAttributes());
254 ds->line_width = GetLineWidth(visAtts); // Multiplies by GetGlobalLineWidthScale
255 parentNode->add(ds);}
256
257 tools::sg::vertices* vtxs = new tools::sg::vertices;
258 vtxs->mode = tools::gl::line_strip(); //polyline
259 parentNode->add(vtxs);
260
261 {for (const auto& i : a_polyline) {
262 vtxs->add(float(i.x()),float(i.y()),float(i.z()));
263 }}
264
265}

Referenced by AddPrimitive(), and AddPrimitive().

◆ AddPrimitive() [5/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Polymarker & a_polymarker)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 267 of file G4ToolsSGSceneHandler.cc.

268{
269 //::printf("debug G4ToolsSGSceneHandler::AddPrimitive(const G4Polymarker&) : %lu, type %d\n",
270 // a_polymarker.size(),a_polymarker.GetMarkerType());
271 if (a_polymarker.size() == 0) return;
272 auto currentNode = GetOrCreateNode();
273 if (!currentNode) return; // Node not available
274
275 // Transformation
276 {tools::sg::matrix* mtx = new tools::sg::matrix;
278 mtx->mtx.value().set_matrix(elem(0,0),elem(0,1),elem(0,2),elem(0,3),
279 elem(1,0),elem(1,1),elem(1,2),elem(1,3),
280 elem(2,0),elem(2,1),elem(2,2),elem(2,3),
281 0, 0, 0, 1);
282 currentNode->add(mtx);}
283
284 {const auto& colour = GetColour(a_polymarker);
285 tools::sg::rgba* mat = new tools::sg::rgba();
286 mat->color =
287 tools::colorf(float(colour.GetRed()),
288 float(colour.GetGreen()),
289 float(colour.GetBlue()),
290 float(colour.GetAlpha()));
291 currentNode->add(mat);}
292
293 MarkerSizeType markerSizeType;
294 G4double markerSize = GetMarkerSize(a_polymarker, markerSizeType);
295
296 switch (a_polymarker.GetMarkerType()) {
297 default:
298 case G4Polymarker::dots:{
299 tools::sg::draw_style* ds = new tools::sg::draw_style;
300 ds->style = tools::sg::draw_points;
301 ds->point_size = markerSize;
302 ds->point_smooth = fpViewer->GetViewParameters().IsDotsSmooth()?true:false;
303 currentNode->add(ds);
304
305 tools::sg::vertices* vtxs = new tools::sg::vertices;
306 vtxs->mode = tools::gl::points();
307 {for (const auto& i : a_polymarker) {
308 vtxs->add(float(i.x()),float(i.y()),float(i.z()));
309 }}
310 currentNode->add(vtxs);
311 }break;
313 {tools::sg::markers* markers = new tools::sg::markers;
314 G4double diameter = markerSize; // OK for "screen-size" (the usual case)
315 if (markerSizeType == G4VSceneHandler::world ) {
316 const G4double scale = 200.; // Roughly pixels per scene
317 diameter *= fpScene->GetExtent().GetExtentRadius()/scale;
318 }
319 markers->size = diameter;
320 markers->style = tools::sg::marker_circle_filled;
321 for (const auto& i : a_polymarker) {
322 markers->add(float(i.x()),float(i.y()),float(i.z()));
323 }
324 currentNode->add(markers);}
325 }break;
327 {tools::sg::markers* markers = new tools::sg::markers;
328 G4double side = markerSize; // OK for "screen-size" (the usual case)
329 if (markerSizeType == G4VSceneHandler::world ) {
330 const G4double scale = 200.; // Roughly pixels per scene
331 side *= fpScene->GetExtent().GetExtentRadius()/scale;
332 }
333 markers->size = side;
334 markers->style = tools::sg::marker_square_filled;
335 for (const auto& i : a_polymarker) {
336 markers->add(float(i.x()),float(i.y()),float(i.z()));
337 }
338 currentNode->add(markers);}
339 }break;
340 }
341}
double G4double
Definition G4Types.hh:83
MarkerType GetMarkerType() const
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)

◆ AddPrimitive() [6/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Square & a_square)
virtual

Implements G4VSceneHandler.

Definition at line 430 of file G4ToolsSGSceneHandler.cc.

431{
432 G4Polymarker oneSquare(a_square);
433 oneSquare.push_back(a_square.GetPosition());
434 oneSquare.SetMarkerType(G4Polymarker::squares);
435 // Call this AddPrimitive to avoid re-doing sub-class code.
437}

◆ AddPrimitive() [7/7]

void G4ToolsSGSceneHandler::AddPrimitive ( const G4Text & a_text)
virtual

Implements G4VSceneHandler.

Definition at line 346 of file G4ToolsSGSceneHandler.cc.

347{
348 //::printf("debug : G4ToolsSGSceneHandler::AddPrimitive(const G4Text&) : 000 : \"%s\"\n",a_text.GetText().c_str());
349 //::printf("debug : G4ToolsSGSceneHandler::AddPrimitive(const G4Text&) : 2D ? %d\n",fProcessing2D);
350 auto pos = a_text.GetPosition();
351 //::printf("debug : Add Text : pos %g %g %g\n",pos.x(),pos.y(),pos.z());
352
353 tools::sg::separator* parentNode = 0;
354 if(fProcessing2D) {
355 parentNode = new tools::sg::separator;
357 fpTransient2DObjects.add(parentNode);
358 } else {
359 fpPersistent2DObjects.add(parentNode);
360 }
361
362 tools::sg::matrix* mtx = new tools::sg::matrix;
363 mtx->set_translate(pos.x(),pos.y(),pos.z());
364 parentNode->add(mtx);
365
366 } else {
367 parentNode = GetOrCreateNode();
368 if (!parentNode) return;
369
370 tools::sg::matrix* mtx = new tools::sg::matrix;
372 mtx->mtx.value().set_matrix(elem(0,0),elem(0,1),elem(0,2),elem(0,3),
373 elem(1,0),elem(1,1),elem(1,2),elem(1,3),
374 elem(2,0),elem(2,1),elem(2,2),elem(2,3),
375 0, 0, 0, 1);
376 parentNode->add(mtx);
377 }
378
379 MarkerSizeType sizeType;
380 G4double size = GetMarkerSize(a_text, sizeType);
381 size *= fMarkerScale;
382
383 {const auto& colour = GetTextColour(a_text);
384 tools::sg::rgba* mat = new tools::sg::rgba();
385 mat->color =
386 tools::colorf(float(colour.GetRed()),
387 float(colour.GetGreen()),
388 float(colour.GetBlue()),
389 float(colour.GetAlpha()));
390 parentNode->add(mat);}
391
392#ifdef TOOLS_USE_FREETYPE
393 toolx::sg::text_freetype_marker* text = new toolx::sg::text_freetype_marker;
394 text->add_embedded_font(tools::sg::font_lato_regular_ttf(),tools::font::lato_regular_ttf);
395 text->font = tools::sg::font_lato_regular_ttf();
396 text->front_face = tools::sg::winding_cw;
397//text->modeling = tools::sg::font_pixmap; //problem with Qt/GL. It slows rendering!
398#else
399 tools::sg::text_hershey_marker* text = new tools::sg::text_hershey_marker;
400//text->encoding.value(a_encoding);
401#endif
402 text->height = float(size); //pixels
403 text->strings.add(a_text.GetText());
404 {switch (a_text.GetLayout()) {
405 default:
406 case G4Text::left:
407 text->hjust = tools::sg::left;
408 break;
409 case G4Text::centre:
410 text->hjust = tools::sg::center;
411 break;
412 case G4Text::right:
413 text->hjust = tools::sg::right;
414 break;
415 }}
416//text->vjust.value(a_vjust);
417 parentNode->add(text);
418
419}
HepGeom::Translate3D G4Translate3D
Layout GetLayout() const
G4String GetText() const
@ centre
Definition G4Text.hh:76
@ right
Definition G4Text.hh:76
@ left
Definition G4Text.hh:76
const G4Colour & GetTextColour(const G4Text &)

◆ ClearStore()

void G4ToolsSGSceneHandler::ClearStore ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 194 of file G4ToolsSGSceneHandler.cc.

195{
196 fpTransient2DObjects.clear();
197 fpPersistent2DObjects.clear();
198 fpTransient3DObjects.clear();
199 fpPersistent3DObjects.clear();
201}

◆ ClearTransientStore()

void G4ToolsSGSceneHandler::ClearTransientStore ( )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 205 of file G4ToolsSGSceneHandler.cc.

206{
207 fpTransient2DObjects.clear();
208 fpTransient3DObjects.clear();
209 if(fpViewer) {
210 G4ToolsSG_viewer* tsg_viewer = dynamic_cast<G4ToolsSG_viewer*>(fpViewer);
211 if(tsg_viewer) tsg_viewer->EmitWindowRender();
212 }
213}
virtual void EmitWindowRender()

◆ CreateSG()

void G4ToolsSGSceneHandler::CreateSG ( )
protected

◆ EstablishBaseNodes()

void G4ToolsSGSceneHandler::EstablishBaseNodes ( )
protected

Definition at line 107 of file G4ToolsSGSceneHandler.cc.

108{
109 // Physical volume objects for each world hang from POs
110 G4TransportationManager* transportationManager = G4TransportationManager::GetTransportationManager ();
111 size_t nWorlds = transportationManager->GetNoWorlds();
112 std::vector<G4VPhysicalVolume*>::iterator iterWorld = transportationManager->GetWorldsIterator();
113 fpPhysicalVolumeObjects.resize(nWorlds);
114 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
115 G4VPhysicalVolume* _world = (*iterWorld);
116 auto entity = new G4ToolsSGNode;
117 fpPersistent3DObjects.add(entity);
118 entity->SetPVNodeID(G4PhysicalVolumeModel::G4PhysicalVolumeNodeID(_world));
119 fpPhysicalVolumeObjects[i] = entity;
120 }
121}
std::vector< G4ToolsSGNode * > fpPhysicalVolumeObjects
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const

Referenced by ClearStore(), and G4ToolsSGSceneHandler().

◆ GetOrCreateNode()

tools::sg::separator * G4ToolsSGSceneHandler::GetOrCreateNode ( )
protected

Definition at line 123 of file G4ToolsSGSceneHandler.cc.

124{ // Retrieve or create a G4ToolsSGNode node suitable for next solid or primitive
125
126 if (fReadyForTransients) { // All transients hang from this node
127 tools::sg::separator* sep = new tools::sg::separator;
128 fpTransient3DObjects.add(sep);
129 return sep;
130 }
131
132 auto* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
133
134 if (!pPVModel) { // Persistent objects (e.g., axes)
135 tools::sg::separator* sep = new tools::sg::separator;
136 fpPersistent3DObjects.add(sep);
137 return sep;
138 }
139
140 // So this is a G4PhysicalVolumeModel
141 typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
142 typedef std::vector<PVNodeID> PVPath;
143 //const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
144 const PVPath& fullPVPath = pPVModel->GetFullPVPath();
145 //G4int currentDepth = pPVModel->GetCurrentDepth();
146 //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
147 //G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
148 //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
149 // Note: pCurrentMaterial may be zero (parallel world).
150
151 // Find appropriate root
152 const size_t nWorlds = fpPhysicalVolumeObjects.size();
153 size_t iWorld = 0;
154 for (; iWorld < nWorlds; ++iWorld) {
155 if (fullPVPath[0].GetPhysicalVolume() ==
156 fpPhysicalVolumeObjects[iWorld]->GetPVNodeID().GetPhysicalVolume()) break;
157 }
158 if (iWorld == nWorlds) {
159 G4Exception("G4ToolsSGSceneHandler::GetOrCreateNode", "ToolsSG-0000", FatalException,
160 "World mis-match - not possible(!?)");
161 }
162
163 // (Re-)establish pv path of root entity
164 G4ToolsSGNode* _world = fpPhysicalVolumeObjects[iWorld];
165 _world->SetPVNodeID(fullPVPath[0]);
166
167 // Provide nodes as required - may be a new node or a pre-existing node
168 G4ToolsSGNode* node = _world; // Working variable - default to world
169 const size_t depth = fullPVPath.size();
170 size_t iDepth = 1;
171 while (iDepth < depth) {
172 const auto& children = node->children();
173 const G4int nChildren = (G4int)children.size();
174 G4int iChild = 0;
175 G4ToolsSGNode* child = nullptr;
176 for (; iChild < nChildren; ++iChild) {
177 child = static_cast<G4ToolsSGNode*>(children[iChild]);
178 if (child->GetPVNodeID() == fullPVPath[iDepth]) break;
179 }
180 if (iChild != nChildren) { // Existing node found
181 node = child; // Must be the ancestor of new node (subsequent iteration)
182 } else {
183 // Add a new node as child of node
184 G4ToolsSGNode* newNode = new G4ToolsSGNode;
185 node->add(newNode);
186 newNode->SetPVNodeID(fullPVPath[iDepth]);
187 node = newNode;
188 }
189 ++iDepth;
190 }
191 return node;
192}
const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID & GetPVNodeID() const
void SetPVNodeID(const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID &id)

Referenced by AddPrimitive(), AddPrimitive(), AddPrimitive(), AddPrimitive(), and AddPrimitive().

◆ GetPersistent2DObjects()

tools::sg::separator & G4ToolsSGSceneHandler::GetPersistent2DObjects ( )
inline

Definition at line 65 of file G4ToolsSGSceneHandler.hh.

◆ GetPersistent3DObjects()

tools::sg::separator & G4ToolsSGSceneHandler::GetPersistent3DObjects ( )
inline

Definition at line 67 of file G4ToolsSGSceneHandler.hh.

◆ GetTransient2DObjects()

tools::sg::separator & G4ToolsSGSceneHandler::GetTransient2DObjects ( )
inline

Definition at line 64 of file G4ToolsSGSceneHandler.hh.

◆ GetTransient3DObjects()

tools::sg::separator & G4ToolsSGSceneHandler::GetTransient3DObjects ( )
inline

Definition at line 66 of file G4ToolsSGSceneHandler.hh.

◆ operator=()

G4ToolsSGSceneHandler & G4ToolsSGSceneHandler::operator= ( const G4ToolsSGSceneHandler & )
inlineprotected

Definition at line 75 of file G4ToolsSGSceneHandler.hh.

75{return *this;}

◆ SetMarkerScale()

void G4ToolsSGSceneHandler::SetMarkerScale ( double a_scale)

Definition at line 344 of file G4ToolsSGSceneHandler.cc.

344{fMarkerScale = a_scale;}

◆ SetPlotterHistograms()

void G4ToolsSGSceneHandler::SetPlotterHistograms ( tools::sg::plots & a_plots)
protected

Definition at line 712 of file G4ToolsSGSceneHandler.cc.

712 {
713 a_plots.clear();
714 G4UImanager* UI = G4UImanager::GetUIpointer();
715 if(UI==NULL) return;
716
717 for(const auto& [region, hid] : fRegionH1s) {
718 tools::sg::plotter* _plotter = a_plots.find_plotter(region);
719 if(_plotter) {
720 std::ostringstream os;
721 os << hid;
722 std::string cmd("/analysis/h1/get ");
723 cmd += std::string(os.str());
724 auto keepControlVerbose = UI->GetVerboseLevel();
725 UI->SetVerboseLevel(0);
726 G4int status = UI->ApplyCommand(cmd.c_str());
727 UI->SetVerboseLevel(keepControlVerbose);
729 G4String hexString = UI->GetCurrentValues("/analysis/h1/get");
730 if(hexString.size()) {
731 void* ptr;
732 std::istringstream is(hexString);
733 is >> ptr;
734 tools::histo::h1d* _h = (tools::histo::h1d*)ptr;
735 tools::sg::plottable* p = new tools::sg::h1d2plot_cp(*_h);
736 _plotter->add_plottable(p); //give ownership of p to sg::plotter.
737 }
738 } else {
739 G4cerr <<
740 "G4ToolsSGSceneHandler::SetPlotterHistograms: cannot get histogram - maybe doesn't exist?"
741 "\n Maybe this app does not do analysis at all?"
742 << G4endl;
743 }
744 }
745 }
746
747 for(const auto& [region, hid] : fRegionH2s) {
748 tools::sg::plotter* _plotter = a_plots.find_plotter(region);
749 if(_plotter) {
750 std::ostringstream os;
751 os << hid;
752 std::string cmd("/analysis/h2/get ");
753 cmd += std::string(os.str());
754 auto keepControlVerbose = UI->GetVerboseLevel();
755 UI->SetVerboseLevel(0);
756 G4int status = UI->ApplyCommand(cmd.c_str());
757 UI->SetVerboseLevel(keepControlVerbose);
759 G4String hexString = UI->GetCurrentValues("/analysis/h2/get");
760 if(hexString.size()) {
761 void* ptr;
762 std::istringstream is(hexString);
763 is >> ptr;
764 tools::histo::h2d* _h = (tools::histo::h2d*)ptr;
765 tools::sg::plottable* p = new tools::sg::h2d2plot_cp(*_h);
766 _plotter->add_plottable(p); //give ownership of p to sg::plotter.
767 }
768 } else {
769 G4cerr <<
770 "G4ToolsSGSceneHandler::SetPlotterHistograms: cannot get histogram - maybe doesn't exist?"
771 "\n Maybe this app does not do analysis at all?"
772 << G4endl;
773 }
774 }
775 }
776}
@ fCommandSucceeded
G4int ApplyCommand(const char *aCommand)
G4int GetVerboseLevel() const
G4String GetCurrentValues(const char *aCommand)
static G4UImanager * GetUIpointer()
void SetVerboseLevel(G4int val)

Referenced by AddPrimitive(), and TouchPlotters().

◆ TouchPlotters()

void G4ToolsSGSceneHandler::TouchPlotters ( tools::sg::node & a_sg)

Definition at line 809 of file G4ToolsSGSceneHandler.cc.

809 {
810 tools::sg::search_action sa(G4cout);
811 const tools::sg::search_action::paths_t& paths = tools::sg::find_paths<tools::sg::plots>(sa,a_sg);
812 for(const auto& path : paths) {
813 tools::sg::plots* _plots = tools::sg::tail<tools::sg::plots>(path);
814 if(_plots) {
815 SetPlotterHistograms(*_plots);
816 }
817 }
818}

Member Data Documentation

◆ fFreetypeNode

tools::sg::base_freetype* G4ToolsSGSceneHandler::fFreetypeNode
protected

◆ fMarkerScale

double G4ToolsSGSceneHandler::fMarkerScale
protected

Definition at line 95 of file G4ToolsSGSceneHandler.hh.

Referenced by AddPrimitive(), G4ToolsSGSceneHandler(), and SetMarkerScale().

◆ fpPersistent2DObjects

tools::sg::separator G4ToolsSGSceneHandler::fpPersistent2DObjects
protected

◆ fpPersistent3DObjects

tools::sg::separator G4ToolsSGSceneHandler::fpPersistent3DObjects
protected

◆ fpPhysicalVolumeObjects

std::vector<G4ToolsSGNode*> G4ToolsSGSceneHandler::fpPhysicalVolumeObjects
protected

Definition at line 92 of file G4ToolsSGSceneHandler.hh.

Referenced by EstablishBaseNodes(), and GetOrCreateNode().

◆ fpTransient2DObjects

tools::sg::separator G4ToolsSGSceneHandler::fpTransient2DObjects
protected

◆ fpTransient3DObjects

tools::sg::separator G4ToolsSGSceneHandler::fpTransient3DObjects
protected

◆ fRegionH1s

std::vector<Region_h1> G4ToolsSGSceneHandler::fRegionH1s
protected

Definition at line 99 of file G4ToolsSGSceneHandler.hh.

Referenced by AddPrimitive(), and SetPlotterHistograms().

◆ fRegionH2s

std::vector<Region_h2> G4ToolsSGSceneHandler::fRegionH2s
protected

Definition at line 100 of file G4ToolsSGSceneHandler.hh.

Referenced by AddPrimitive(), and SetPlotterHistograms().

◆ fSceneIdCount

G4int G4ToolsSGSceneHandler::fSceneIdCount = 0
staticprotected

Definition at line 84 of file G4ToolsSGSceneHandler.hh.

Referenced by G4ToolsSGSceneHandler().


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