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

#include <G4VtkSceneHandler.hh>

Inheritance diagram for G4VtkSceneHandler:

Public Member Functions

 G4VtkSceneHandler (G4VGraphicsSystem &system, const G4String &name)
 ~G4VtkSceneHandler () override=default
void AddPrimitive (const G4Polyline &) override
void AddPrimitive (const G4Text &) override
void AddPrimitive (const G4Circle &) override
void AddPrimitive (const G4Square &) override
void AddPrimitive (const G4Polyhedron &) override
void AddPrimitive (const G4Polymarker &polymarker) override
void AddSolid (const G4Box &box) override
void AddCompound (const G4Mesh &mesh) override
void Modified ()
void ClearStore () override
void ClearTransientStore () override
G4VtkVisContext MakeDefaultVisContext ()
G4VtkStoreGetStore ()
G4VtkStoreGetTransientStore ()
virtual void Print ()
void SetPolyhedronPipeline (const G4String &str)
virtual void AddPrimitive (const G4Plotter &)
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 > &)
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 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 Attributes

G4VtkStore store = G4VtkStore("perm")
G4VtkStore transientStore = G4VtkStore("trans")
G4String polyhedronPipelineType
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

Friends

class G4VtkViewer

Additional Inherited Members

Public Types inherited from G4VSceneHandler
enum  MarkerSizeType { world , screen }
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

Detailed Description

Definition at line 82 of file G4VtkSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4VtkSceneHandler()

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

Definition at line 64 of file G4VtkSceneHandler.cc.

65 : G4VSceneHandler(system, fSceneIdCount++, name), polyhedronPipelineType(G4String("tensor"))
66{}
G4VSceneHandler(G4VGraphicsSystem &system, G4int id, const G4String &name="")

Referenced by G4VtkQtSceneHandler::G4VtkQtSceneHandler().

◆ ~G4VtkSceneHandler()

G4VtkSceneHandler::~G4VtkSceneHandler ( )
overridedefault

Member Function Documentation

◆ AddCompound() [1/6]

void G4VtkSceneHandler::AddCompound ( const G4Mesh & mesh)
overridevirtual

Reimplemented from G4VSceneHandler.

Definition at line 259 of file G4VtkSceneHandler.cc.

260{
261#ifdef G4VTKDEBUG
262 G4cout << "G4VtkSceneHandler::AddCompound> mesh type " << mesh.GetMeshType() << " "
263 << fpViewer->GetViewParameters().GetSpecialMeshRenderingOption() << G4endl;
264#endif
265
266 if(fpViewer->GetViewParameters().GetSpecialMeshRenderingOption() == G4ViewParameters::meshAsDefault)
267 {
268 auto vc = MakeDefaultVisContext();
269
271 transientStore.AddCompound(mesh, vc);
272 else
273 store.AddCompound(mesh, vc);
274 }
275 else {
277 }
278}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
MeshType GetMeshType() const
Definition G4Mesh.hh:75
void StandardSpecialMeshRendering(const G4Mesh &)
G4VtkVisContext MakeDefaultVisContext()

◆ 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
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 G4VtkSceneHandler::AddPrimitive ( const G4Circle & circle)
overridevirtual

Implements G4VSceneHandler.

Definition at line 95 of file G4VtkSceneHandler.cc.

96{
97#ifdef G4VTKDEBUG
98 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle)" << G4endl;
99#endif
100
101 auto vc = MakeDefaultVisContext();
103 vc.fSize = GetMarkerSize(circle, sizeType);
104
106 transientStore.AddPrimitive(circle, vc);
107 else
108 store.AddPrimitive(circle, vc);
109}
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)

◆ AddPrimitive() [2/7]

void G4VSceneHandler::AddPrimitive ( const G4Plotter & )
virtual

Reimplemented from G4VSceneHandler.

Definition at line 195 of file G4VSceneHandler.cc.

508 {
509 G4warn << "WARNING: Plotter not implemented for " << fSystem.GetName() << G4endl;
510 G4warn << " Open a plotter-aware graphics system or remove plotter with" << G4endl;
511 G4warn << " /vis/scene/removeModel Plotter" << G4endl;
512}
#define G4warn
Definition G4Scene.cc:41
G4VGraphicsSystem & fSystem

◆ AddPrimitive() [3/7]

void G4VtkSceneHandler::AddPrimitive ( const G4Polyhedron & polyhedron)
overridevirtual

Implements G4VSceneHandler.

Definition at line 127 of file G4VtkSceneHandler.cc.

128{
129#ifdef G4VTKDEBUG
130 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron)" << G4endl;
131#endif
132
133 auto vc = MakeDefaultVisContext();
134 auto visAtt = vc.fViewer->GetApplicableVisAttributes(polyhedron.GetVisAttributes());
135 auto colour = visAtt->GetColour();
136
137 vc.fDrawingStyle = GetDrawingStyle(visAtt);
138 vc.alpha = colour.GetAlpha();
139 vc.red = colour.GetRed();
140 vc.green = colour.GetGreen();
141 vc.blue = colour.GetBlue();
142
143 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
144 if (pPVModel != nullptr) {
145 vc.fDepth = pPVModel->GetCurrentDepth();
146 vc.fDescription = pPVModel->GetCurrentDescription();
147 }
148
150 if (polyhedronPipelineType == "tensor")
151 transientStore.AddPrimitiveTensorGlyph(polyhedron, vc);
152 else if (polyhedronPipelineType == "append")
153 transientStore.AddPrimitiveAppend(polyhedron, vc);
154 else if (polyhedronPipelineType == "bake")
155 transientStore.AddPrimitiveTransformBake(polyhedron, vc);
156 else if (polyhedronPipelineType == "separate")
157 transientStore.AddPrimitiveSeparate(polyhedron, vc);
158 }
159 else {
160 if (polyhedronPipelineType == "tensor")
161 store.AddPrimitiveTensorGlyph(polyhedron, vc);
162 else if (polyhedronPipelineType == "append")
163 store.AddPrimitiveAppend(polyhedron, vc);
164 else if (polyhedronPipelineType == "bake")
165 store.AddPrimitiveTransformBake(polyhedron, vc);
166 else if (polyhedronPipelineType == "separate")
167 store.AddPrimitiveSeparate(polyhedron, vc);
168 }
169}
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
const G4VisAttributes * GetVisAttributes() const

◆ AddPrimitive() [4/7]

void G4VtkSceneHandler::AddPrimitive ( const G4Polyline & polyline)
overridevirtual

Implements G4VSceneHandler.

Definition at line 68 of file G4VtkSceneHandler.cc.

69{
70#ifdef G4VTKDEBUG
71 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline&)" << G4endl;
72#endif
73 auto vc = MakeDefaultVisContext();
74
76 transientStore.AddPrimitive(polyline, vc);
77 else
78 store.AddPrimitive(polyline, vc);
79}

Referenced by G4VtkQtSceneHandler::AddPrimitive(), G4VtkQtSceneHandler::AddPrimitive(), G4VtkQtSceneHandler::AddPrimitive(), G4VtkQtSceneHandler::AddPrimitive(), G4VtkQtSceneHandler::AddPrimitive(), and G4VtkQtSceneHandler::~G4VtkQtSceneHandler().

◆ AddPrimitive() [5/7]

void G4VtkSceneHandler::AddPrimitive ( const G4Polymarker & polymarker)
inlineoverridevirtual

Reimplemented from G4VSceneHandler.

Definition at line 101 of file G4VtkSceneHandler.hh.

102 {
104 }
virtual void AddPrimitive(const G4Polyline &)=0

◆ AddPrimitive() [6/7]

void G4VtkSceneHandler::AddPrimitive ( const G4Square & square)
overridevirtual

Implements G4VSceneHandler.

Definition at line 111 of file G4VtkSceneHandler.cc.

112{
113#ifdef G4VTKDEBUG
114 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Square& square)" << G4endl;
115#endif
116
117 auto vc = MakeDefaultVisContext();
119 vc.fSize = GetMarkerSize(square, sizeType);
120
122 transientStore.AddPrimitive(square, vc);
123 else
124 store.AddPrimitive(square, vc);
125}

◆ AddPrimitive() [7/7]

void G4VtkSceneHandler::AddPrimitive ( const G4Text & text)
overridevirtual

Implements G4VSceneHandler.

Definition at line 81 of file G4VtkSceneHandler.cc.

82{
83#ifdef G4VTKDEBUG
84 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Text& text)" << G4endl;
85#endif
86
87 auto vc = MakeDefaultVisContext();
88
90 transientStore.AddPrimitive(text, vc);
91 else
92 store.AddPrimitive(text, vc);
93}

◆ AddSolid() [1/14]

void G4VtkSceneHandler::AddSolid ( const G4Box & box)
overridevirtual

Reimplemented from G4VSceneHandler.

Definition at line 215 of file G4VtkSceneHandler.cc.

216{
218
219 return;
220
221 const G4VModel* pv_model = GetModel();
222 if (pv_model == nullptr) {
223 return;
224 }
225
226 auto pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
227 if (pPVModel == nullptr) {
228 return;
229 }
230
231 //-- debug information
232#ifdef G4VTKDEBUG
233 G4VPhysicalVolume* pv = pPVModel->GetCurrentPV();
234 G4LogicalVolume* lv = pv->GetLogicalVolume();
235 G4cout << "name=" << box.GetName() << " volumeType=" << pv->VolumeType()
236 << " pvName=" << pv->GetName() << " lvName=" << lv->GetName()
237 << " multiplicity=" << pv->GetMultiplicity() << " isparametrised=" << pv->IsParameterised()
238 << " isreplicated=" << pv->IsReplicated()
239 << " parametrisation=" << pv->GetParameterisation()
240 << G4endl
241
242 G4Material* mat = pPVModel->GetCurrentMaterial();
243 G4String name = mat->GetName();
244 G4double dens = mat->GetDensity() / (g / cm3);
245 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
246 G4int depth = pPVModel->GetCurrentDepth();
247 G4cout << " name : " << box.GetName() << G4endl;
248 G4cout << " copy no.: " << copyNo << G4endl;
249 G4cout << " depth : " << depth << G4endl;
250 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
251 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
252 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
253 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
254 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
255 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
256#endif
257}
double G4double
Definition G4Types.hh:83
const G4String & GetName() const
virtual G4bool IsReplicated() const =0
virtual EVolume VolumeType() const =0
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4bool IsParameterised() const =0
G4VModel * GetModel() const
virtual void AddSolid(const G4Box &)
G4String GetName() const
const char * name(G4int ptype)

◆ AddSolid() [2/14]

void G4VSceneHandler::AddSolid ( const G4Cons & cons)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 102 of file G4VSceneHandler.cc.

271 {
272 AddSolidT (cons);
273}
void AddSolidT(const T &solid)

◆ AddSolid() [3/14]

void G4VSceneHandler::AddSolid ( const G4Ellipsoid & ellipsoid)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 112 of file G4VSceneHandler.cc.

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

◆ AddSolid() [4/14]

void G4VSceneHandler::AddSolid ( const G4Orb & orb)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 103 of file G4VSceneHandler.cc.

275 {
277}

◆ AddSolid() [5/14]

void G4VSceneHandler::AddSolid ( const G4Para & para)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 104 of file G4VSceneHandler.cc.

279 {
280 AddSolidT (para);
281}

◆ AddSolid() [6/14]

void G4VSceneHandler::AddSolid ( const G4Polycone & polycone)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 113 of file G4VSceneHandler.cc.

307 {
308 AddSolidT (polycone);
309}

◆ AddSolid() [7/14]

void G4VSceneHandler::AddSolid ( const G4Polyhedra & polyhedra)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 114 of file G4VSceneHandler.cc.

311 {
312 AddSolidT (polyhedra);
313}

◆ AddSolid() [8/14]

void G4VSceneHandler::AddSolid ( const G4Sphere & sphere)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 105 of file G4VSceneHandler.cc.

283 {
285}

◆ AddSolid() [9/14]

void G4VSceneHandler::AddSolid ( const G4TessellatedSolid & tess)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 115 of file G4VSceneHandler.cc.

315 {
316 AddSolidT (tess);
317}

◆ AddSolid() [10/14]

void G4VSceneHandler::AddSolid ( const G4Torus & torus)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 106 of file G4VSceneHandler.cc.

287 {
289}

◆ AddSolid() [11/14]

void G4VSceneHandler::AddSolid ( const G4Trap & trap)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 107 of file G4VSceneHandler.cc.

291 {
292 AddSolidT (trap);
293}

◆ AddSolid() [12/14]

void G4VSceneHandler::AddSolid ( const G4Trd & trd)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 108 of file G4VSceneHandler.cc.

295 {
296 AddSolidT (trd);
297}

◆ AddSolid() [13/14]

void G4VSceneHandler::AddSolid ( const G4Tubs & tubs)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 109 of file G4VSceneHandler.cc.

299 {
300 AddSolidT (tubs);
301}

◆ AddSolid() [14/14]

void G4VSceneHandler::AddSolid ( const G4VSolid & solid)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 118 of file G4VSceneHandler.cc.

319 {
320 AddSolidT (solid);
321}

◆ ClearStore()

void G4VtkSceneHandler::ClearStore ( )
overridevirtual

Reimplemented from G4VSceneHandler.

Definition at line 181 of file G4VtkSceneHandler.cc.

182{
183#ifdef G4VTKDEBUG
184 G4cout << "G4VtkSceneHandler::ClearStore()" << G4endl;
185#endif
186 store.Clear();
187 transientStore.Clear();
188}

◆ ClearTransientStore()

void G4VtkSceneHandler::ClearTransientStore ( )
overridevirtual

Reimplemented from G4VSceneHandler.

Definition at line 190 of file G4VtkSceneHandler.cc.

191{
192#ifdef G4VTKDEBUG
193 G4cout << "G4VtkSceneHandler::ClearTransientStore()" << G4endl;
194#endif
195 transientStore.Clear();
196}

◆ GetStore()

G4VtkStore & G4VtkSceneHandler::GetStore ( )
inline

Definition at line 116 of file G4VtkSceneHandler.hh.

116{ return store; }

◆ GetTransientStore()

G4VtkStore & G4VtkSceneHandler::GetTransientStore ( )
inline

Definition at line 117 of file G4VtkSceneHandler.hh.

117{ return transientStore; }

Referenced by G4VtkMessenger::SetNewValue().

◆ MakeDefaultVisContext()

G4VtkVisContext G4VtkSceneHandler::MakeDefaultVisContext ( )

Definition at line 198 of file G4VtkSceneHandler.cc.

199{
200 auto vc = G4VtkVisContext(dynamic_cast<G4VtkViewer*>(fpViewer), fpVisAttribs, fProcessing2D,
202
203 if (fpVisAttribs != nullptr) {
204 G4Colour c = fpVisAttribs->GetColour();
205 vc.red = c.GetRed();
206 vc.green = c.GetGreen();
207 vc.blue = c.GetBlue();
208 vc.alpha = c.GetAlpha();
209 vc.fDrawingStyle = fpViewer->GetViewParameters().GetDrawingStyle();
210 }
211
212 return vc;
213}
G4double GetBlue() const
Definition G4Colour.hh:172
G4double GetAlpha() const
Definition G4Colour.hh:173
G4double GetRed() const
Definition G4Colour.hh:170
G4double GetGreen() const
Definition G4Colour.hh:171
G4Transform3D fObjectTransformation
const G4VisAttributes * fpVisAttribs

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

◆ Modified()

void G4VtkSceneHandler::Modified ( )

Definition at line 171 of file G4VtkSceneHandler.cc.

172{
173#ifdef G4VTKDEBUG
174 G4cout << "G4VtkSceneHandler::Modified()" << G4endl;
175#endif
176
177 store.Modified();
178 transientStore.Modified();
179}

Referenced by G4VtkOffscreenViewer::FinishView(), G4VtkQtViewer::FinishView(), and G4VtkViewer::FinishView().

◆ Print()

void G4VtkSceneHandler::Print ( )
virtual

Definition at line 280 of file G4VtkSceneHandler.cc.

280{}

◆ SetPolyhedronPipeline()

void G4VtkSceneHandler::SetPolyhedronPipeline ( const G4String & str)

Definition at line 282 of file G4VtkSceneHandler.cc.

283{
285}

Referenced by G4VtkViewer::SetPolyhedronPipeline().

◆ G4VtkViewer

friend class G4VtkViewer
friend

Definition at line 132 of file G4VtkSceneHandler.hh.

Referenced by G4VtkViewer, and MakeDefaultVisContext().

Member Data Documentation

◆ fSceneIdCount

G4int G4VtkSceneHandler::fSceneIdCount = 0
staticprotected

Definition at line 124 of file G4VtkSceneHandler.hh.

Referenced by G4VtkSceneHandler().

◆ polyhedronPipelineType

G4String G4VtkSceneHandler::polyhedronPipelineType
protected

Definition at line 129 of file G4VtkSceneHandler.hh.

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

◆ store

G4VtkStore G4VtkSceneHandler::store = G4VtkStore("perm")
protected

◆ transientStore

G4VtkStore G4VtkSceneHandler::transientStore = G4VtkStore("trans")
protected

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