Geant4 11.1.1
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)
 
virtual ~G4VtkSceneHandler ()=default
 
void AddPrimitive (const G4Polyline &)
 
void AddPrimitive (const G4Text &)
 
void AddPrimitive (const G4Circle &)
 
void AddPrimitive (const G4Square &)
 
void AddPrimitive (const G4Polyhedron &)
 
void AddPrimitiveTensorGlyph (const G4Polyhedron &)
 
void AddPrimitiveBakedTransform (const G4Polyhedron &)
 
void AddPrimitive (const G4Polymarker &polymarker)
 
void AddSolid (const G4Box &box)
 
void AddCompound (const G4Mesh &mesh)
 
void Modified ()
 
void Clear ()
 
- 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 Attributes

std::map< std::size_t, const G4VisAttributes * > polylineVisAttributesMap
 
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polylineDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polylineLineMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polylinePolyDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > polylinePolyDataMapperMap
 
std::map< std::size_t, vtkSmartPointer< vtkActor > > polylinePolyDataActorMap
 
std::map< std::size_t, const G4VisAttributes * > circleVisAttributesMap
 
std::map< std::size_t, vtkSmartPointer< vtkPoints > > circleDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > circlePolyDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > circleFilterMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > circlePolyDataMapperMap
 
std::map< std::size_t, vtkSmartPointer< vtkActor > > circlePolyDataActorMap
 
std::map< std::size_t, const G4VisAttributes * > squareVisAttributesMap
 
std::map< std::size_t, vtkSmartPointer< vtkPoints > > squareDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > squarePolyDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > squareFilterMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > squarePolyDataMapperMap
 
std::map< std::size_t, vtkSmartPointer< vtkActor > > squarePolyDataActorMap
 
std::map< std::size_t, const G4VisAttributes * > polyhedronVisAttributesMap
 
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polyhedronDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polyhedronPolyMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polyhedronPolyDataMap
 
std::map< std::size_t, std::size_t > polyhedronPolyDataCountMap
 
std::map< std::size_t, vtkSmartPointer< vtkPoints > > instancePositionMap
 
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceRotationMap
 
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceColoursMap
 
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > instancePolyDataMap
 
std::map< std::size_t, vtkSmartPointer< vtkTensorGlyphColor > > instanceTensorGlyphMap
 
std::map< std::size_t, vtkSmartPointer< vtkActor > > instanceActorMap
 
- 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 G4VtkViewer
 

Additional Inherited Members

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

Detailed Description

Definition at line 76 of file G4VtkSceneHandler.hh.

Constructor & Destructor Documentation

◆ G4VtkSceneHandler()

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

Definition at line 118 of file G4VtkSceneHandler.cc.

119 :
120 G4VSceneHandler(system, fSceneIdCount++, name)
121{}
static G4int fSceneIdCount

◆ ~G4VtkSceneHandler()

virtual G4VtkSceneHandler::~G4VtkSceneHandler ( )
virtualdefault

Member Function Documentation

◆ AddCompound()

void G4VtkSceneHandler::AddCompound ( const G4Mesh mesh)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 763 of file G4VtkSceneHandler.cc.

764{
765#ifdef G4VTKDEBUG
766 G4cout << "G4VtkSceneHandler::AddCompound" << G4endl;
767#endif
768
769 if(mesh.GetMeshType() != G4Mesh::rectangle || mesh.GetMeshDepth() != 3)
770 {
772 }
773
774 auto container = mesh.GetContainerVolume();
775 auto rep1 = container->GetLogicalVolume()->GetDaughter(0);
776 auto rep2 = rep1->GetLogicalVolume()->GetDaughter(0);
777 auto rep3 = rep2->GetLogicalVolume()->GetDaughter(0);
778
779 EAxis rep1_axis, rep2_axis, rep3_axis;
780 G4int rep1_nReplicas, rep2_nReplicas, rep3_nReplicas;
781 G4double rep1_width, rep2_width, rep3_width;
782 G4double rep1_offset, rep2_offset, rep3_offset;
783 G4bool rep1_consuming, rep2_consuming, rep3_consuming;
784
785 rep1->GetReplicationData(rep1_axis,rep1_nReplicas, rep1_width, rep1_offset, rep1_consuming);
786 rep2->GetReplicationData(rep2_axis,rep2_nReplicas, rep2_width, rep2_offset, rep2_consuming);
787 rep3->GetReplicationData(rep3_axis,rep3_nReplicas, rep3_width, rep3_offset, rep3_consuming);
788
789 // Instantiate a temporary G4PhysicalVolumeModel
791 tmpMP.SetCulling(false); // This avoids drawing transparent...
792 tmpMP.SetCullingInvisible(false); // ... or invisble volumes.
793 const G4bool useFullExtent = false; // To avoid calculating the extent
794
796 G4Transform3D(), &tmpMP, useFullExtent);
797
798 // Instantiate a pseudo scene so that we can make a "private" descent and fill vtkImageData
799 vtkSmartPointer<vtkImageData> imagedata = vtkSmartPointer<vtkImageData>::New();
800 imagedata->SetDimensions(rep1_nReplicas+1,rep2_nReplicas+1,rep3_nReplicas+1);
801 imagedata->SetSpacing(rep1_width,rep2_width,rep2_width);
802 imagedata->SetOrigin(0,0,0);
803 imagedata->AllocateScalars(VTK_DOUBLE, 1);
804
805 G4double halfX = 0., halfY = 0., halfZ = 0.;
806 struct PseudoScene : public G4PseudoScene
807 {
808 PseudoScene(
809 G4PhysicalVolumeModel* pvModel, // input...the following are output
810 vtkImageData *vtkId,
811 G4int nx, G4int ny, G4int nz,
812 G4double& halfX, G4double& halfY, G4double& halfZ) : fpPVModel(pvModel)
813 , fVtkId(vtkId)
814 , fNx(nx)
815 , fNy(ny)
816 , fNz(nz)
817 , fHalfX(halfX)
818 , fHalfY(halfY)
819 , fHalfZ(halfZ)
820 {}
821 using G4PseudoScene::AddSolid; // except for...
822 void AddSolid(const G4Box& box)
823 {
824 const G4Colour& colour = fpPVModel->GetCurrentLV()->GetVisAttributes()->GetColour();
825 // G4double density = fpPVModel->GetCurrentLV()->GetMaterial()->GetDensity();
826 const G4ThreeVector& position = fpCurrentObjectTransformation->getTranslation();
827 fHalfX = box.GetXHalfLength();
828 fHalfY = box.GetYHalfLength();
829 fHalfZ = box.GetZHalfLength();
830 fVtkId->SetScalarComponentFromDouble(int(position.x()/(2*fHalfX))+fNx,
831 int(position.y()/(2*fHalfY))+fNy,
832 int(position.z()/(2*fHalfZ))+fNz,
833 0,
834 (colour.GetRed() + colour.GetBlue() + colour.GetGreen())/3.0);
835 }
836 G4PhysicalVolumeModel* fpPVModel;
837 vtkImageData *fVtkId;
838 G4int fNx, fNy, fNz;
839 G4double &fHalfX, &fHalfY, &fHalfZ;
840 }
841 // construct the pseudoScene
842 pseudoScene(&tmpPVModel, imagedata, rep1_nReplicas, rep2_nReplicas, rep3_nReplicas, halfX, halfY, halfZ);
843
844 // Make private descent into the nested parameterisation
845 tmpPVModel.DescribeYourselfTo(pseudoScene);
846
847 vtkSmartPointer<vtkOpenGLGPUVolumeRayCastMapper> volumeMapper = vtkSmartPointer<vtkOpenGLGPUVolumeRayCastMapper>::New();
848 volumeMapper->SetInputData(imagedata);
849
850 vtkNew<vtkVolume> volume;
851 volume->SetMapper(volumeMapper);
852
853 vtkSmartPointer<vtkMatrix4x4> vtkMatrix = vtkSmartPointer<vtkMatrix4x4>::New();
854
855 vtkMatrix->SetElement(0,0,fObjectTransformation.xx());
856 vtkMatrix->SetElement(0,1,fObjectTransformation.xy());
857 vtkMatrix->SetElement(0,2,fObjectTransformation.xz());
858
859 vtkMatrix->SetElement(1,0,fObjectTransformation.yx());
860 vtkMatrix->SetElement(1,1,fObjectTransformation.yy());
861 vtkMatrix->SetElement(1,2,fObjectTransformation.yz());
862
863 vtkMatrix->SetElement(2,0,fObjectTransformation.zx());
864 vtkMatrix->SetElement(2,1,fObjectTransformation.zy());
865 vtkMatrix->SetElement(2,2,fObjectTransformation.zz());
866
867 vtkMatrix->SetElement(3,0, fObjectTransformation.dx());
868 vtkMatrix->SetElement(3,1, fObjectTransformation.dy());
869 vtkMatrix->SetElement(3,2, fObjectTransformation.dz());
870 vtkMatrix->SetElement(3,3, 1);
871
872 volume->SetUserMatrix(vtkMatrix);
873
874 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
875 pVtkViewer->renderer->AddVolume(volume);
876}
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition: G4Colour.cc:155
G4double GetBlue() const
Definition: G4Colour.hh:154
G4double GetRed() const
Definition: G4Colour.hh:152
G4double GetGreen() const
Definition: G4Colour.hh:153
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
MeshType GetMeshType() const
Definition: G4Mesh.hh:75
G4VPhysicalVolume * GetContainerVolume() const
Definition: G4Mesh.hh:73
@ rectangle
Definition: G4Mesh.hh:53
G4int GetMeshDepth() const
Definition: G4Mesh.hh:76
void SetCulling(G4bool)
void SetCullingInvisible(G4bool)
void AddSolid(const G4Box &solid)
G4LogicalVolume * GetLogicalVolume() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
G4Transform3D fObjectTransformation
G4VViewer * fpViewer
virtual void AddCompound(const G4VTrajectory &)
vtkNew< vtkRenderer > renderer
Definition: G4VtkViewer.hh:186
double dy() const
Definition: Transform3D.h:287
double zz() const
Definition: Transform3D.h:281
double yz() const
Definition: Transform3D.h:272
double dz() const
Definition: Transform3D.h:290
double dx() const
Definition: Transform3D.h:284
double xy() const
Definition: Transform3D.h:260
double zx() const
Definition: Transform3D.h:275
double yx() const
Definition: Transform3D.h:266
double zy() const
Definition: Transform3D.h:278
double xx() const
Definition: Transform3D.h:257
double yy() const
Definition: Transform3D.h:269
double xz() const
Definition: Transform3D.h:263
EAxis
Definition: geomdefs.hh:54

◆ AddPrimitive() [1/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Circle circle)
virtual

Implements G4VSceneHandler.

Definition at line 285 of file G4VtkSceneHandler.cc.

285 {
286
287 MarkerSizeType sizeType;
288 G4double size= GetMarkerSize(circle, sizeType);
289 if(fProcessing2D) {sizeType = screen;}
290 else {sizeType = world;}
291
292 // Get vis attributes - pick up defaults if none.
294 G4Color colour = pVA->GetColour();
295 G4double opacity = colour.GetAlpha();
296 // G4bool isVisible = pVA->IsVisible();
297
298#ifdef G4VTKDEBUG
299 G4cout << "=================================" << G4endl;
300 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> " << " radius:" << size << " sizeType:" << sizeType << G4endl;
301 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> colour: " << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
302 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> alpha: " << colour.GetAlpha() << G4endl;
303#endif
304
305 if (sizeType == world) {
306 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
307 if (circleVisAttributesMap.find(hash) == circleVisAttributesMap.end()) {
308 circleVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
309
310 vtkSmartPointer<vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
311 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
312 vtkSmartPointer<vtkVertexGlyphFilter> filter = vtkSmartPointer<vtkVertexGlyphFilter>::New();
313 vtkSmartPointer<vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
314 vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
315
316 polyData->SetPoints(data);
317 filter->SetInputData(polyData);
318 mapper->SetInputConnection(filter->GetOutputPort());
319 actor->SetMapper(mapper);
320
321 // Setup actor and mapper
322 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
323 actor->GetProperty()->SetOpacity(opacity);
324 actor->SetVisibility(1);
325 actor->GetProperty()->SetRenderPointsAsSpheres(true);
326 actor->GetProperty()->SetPointSize(size*5);
327
328 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
329 pVtkViewer->renderer->AddActor(actor);
330
331 circleDataMap.insert(std::pair<std::size_t,vtkSmartPointer<vtkPoints>>(hash, data));
332 circlePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyData>>(hash, polyData));
333 circleFilterMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>>(hash, filter));
334 circlePolyDataMapperMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyDataMapper>>(hash, mapper));
335 circlePolyDataActorMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkActor>>(hash, actor));
336 }
337
338 // Data data point
340 G4Point3D posPrime = rot*circle.GetPosition();
341 circleDataMap[hash]->InsertNextPoint(fObjectTransformation.dx() + posPrime.x(),
342 fObjectTransformation.dy() + posPrime.y(),
343 fObjectTransformation.dz() + posPrime.z());
344 }
345 else if (sizeType == screen) {
346
347 }
348}
G4double GetAlpha() const
Definition: G4Colour.hh:155
G4Point3D GetPosition() const
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4Colour & GetColour() const
const G4VisAttributes * GetVisAttributes() const
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > circleFilterMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > circlePolyDataMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > circleDataMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > circlePolyDataMapperMap
std::map< std::size_t, const G4VisAttributes * > circleVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > circlePolyDataActorMap
CLHEP::HepRotation getRotation() const

◆ AddPrimitive() [2/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Polyhedron polyhedron)
virtual

Implements G4VSceneHandler.

Definition at line 417 of file G4VtkSceneHandler.cc.

417 {
418 AddPrimitiveTensorGlyph(polyhedron);
419}
void AddPrimitiveTensorGlyph(const G4Polyhedron &)

◆ AddPrimitive() [3/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Polyline polyline)
virtual

Implements G4VSceneHandler.

Definition at line 144 of file G4VtkSceneHandler.cc.

144 {
145
147 // GetMarkerSize(text,sizeType);
148 if(fProcessing2D) {sizeType = screen;}
149 else {sizeType = world;}
150
151 // Get vis attributes - pick up defaults if none.
152 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(polyline.GetVisAttributes());
153 G4Color colour = pVA->GetColour();
154 G4double opacity = colour.GetAlpha();
155 G4double lineWidth = pVA->GetLineWidth();
156
157#ifdef G4VTKDEBUG
158 G4cout << "=================================" << G4endl;
159 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> vis hash: " << sizeType << " " << std::hash<G4VisAttributes>{}(*pVA) << G4endl;
160 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> sizeType: " << sizeType << G4endl;
161 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> isVisible: " << pVA->IsVisible() << G4endl;
162 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> isDaughtersInvisible: " << pVA->IsDaughtersInvisible() << G4endl;
163 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> colour: " << colour.GetRed() << " " << colour.GetGreen() << " " << colour.GetBlue() << G4endl;
164 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> alpha: " << colour.GetAlpha() << G4endl;
165 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> lineWidth: " << lineWidth << G4endl;
166#endif
167
168 if(sizeType == world) {
169 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
170 if (polylineVisAttributesMap.find(hash) == polylineVisAttributesMap.end()) {
171 polylineVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
172
173 vtkSmartPointer <vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
174 vtkSmartPointer <vtkCellArray> lines = vtkSmartPointer<vtkCellArray>::New();
175 vtkSmartPointer <vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
176 vtkSmartPointer <vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
177 vtkSmartPointer <vtkActor> actor = vtkSmartPointer<vtkActor>::New();
178
179 polyData->SetPoints(data);
180 polyData->SetLines(lines);
181 mapper->SetInputData(polyData);
182 actor->SetMapper(mapper);
183
184 // Setup actor and mapper
185 actor->GetProperty()->SetLineWidth(lineWidth);
186 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(),
187 colour.GetBlue());
188 actor->GetProperty()->SetOpacity(opacity);
189 actor->SetVisibility(1);
190 // actor->GetProperty()->BackfaceCullingOn();
191 // actor->GetProperty()->FrontfaceCullingOn();
192
193 auto pVtkViewer = dynamic_cast<G4VtkViewer*>(fpViewer);
194 pVtkViewer->renderer->AddActor(actor);
195
196 polylineDataMap.insert(
197 std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, data));
198 polylineLineMap.insert(
199 std::pair<std::size_t, vtkSmartPointer<vtkCellArray>>(hash, lines));
200 polylinePolyDataMap.insert(
201 std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash, polyData));
203 std::pair<std::size_t, vtkSmartPointer<vtkPolyDataMapper>>(hash,
204 mapper));
206 std::pair<std::size_t, vtkSmartPointer<vtkActor>>(hash, actor));
207 }
208
209 // Data data
210 const size_t nLines = polyline.size();
211
212 for (size_t i = 0; i < nLines; ++i) {
213 auto id = polylineDataMap[hash]->InsertNextPoint(polyline[i].x(),
214 polyline[i].y(),
215 polyline[i].z());
216 if (i < nLines - 1) {
217 vtkSmartPointer <vtkLine> line = vtkSmartPointer<vtkLine>::New();
218 line->GetPointIds()->SetId(0, id);
219 line->GetPointIds()->SetId(1, id + 1);
220 polylineLineMap[hash]->InsertNextCell(line);
221 }
222 }
223 }
224 else if (sizeType == screen ) {
225
226 }
227}
G4double GetLineWidth() const
G4bool IsDaughtersInvisible() const
G4bool IsVisible() const
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > polylinePolyDataMapperMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polylineLineMap
std::map< std::size_t, const G4VisAttributes * > polylineVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polylinePolyDataMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > polylinePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polylineDataMap

◆ AddPrimitive() [4/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Polymarker polymarker)
inlinevirtual

Reimplemented from G4VSceneHandler.

Definition at line 97 of file G4VtkSceneHandler.hh.

98 {
99#ifdef G4VTKDEBUG
100 G4cout << "=================================" << G4endl;
101 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polymarker& polymarker) called." << G4endl;
102#endif
104 }
virtual void AddPrimitive(const G4Polyline &)=0

◆ AddPrimitive() [5/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Square square)
virtual

Implements G4VSceneHandler.

Definition at line 350 of file G4VtkSceneHandler.cc.

350 {
351
352 MarkerSizeType sizeType;
353 G4double size = GetMarkerSize (square, sizeType);
354
355 // Get vis attributes - pick up defaults if none.
356 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(square.GetVisAttributes());
357 G4Color colour = pVA->GetColour();
358 G4double opacity = colour.GetAlpha();
359 // G4bool isVisible = pVA->IsVisible();
360
361 // Draw in world coordinates.
362 vtkSmartPointer<vtkRegularPolygonSource> polygonSource = vtkSmartPointer<vtkRegularPolygonSource>::New();
363 polygonSource->SetNumberOfSides(4);
364 polygonSource->SetRadius(size);
365
366
367#ifdef G4VTKDEBUG
368 G4cout << "=================================" << G4endl;
369 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Square& square) called" << G4endl;
370 G4cout << square.GetPosition().x() << " " << square.GetPosition().y() << " " << square.GetPosition().z() << G4endl;
371 //PrintThings();
372#endif
373
374 if (sizeType == world) {
375 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
376 if (squareVisAttributesMap.find(hash) == squareVisAttributesMap.end()) {
377 squareVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
378
379 vtkSmartPointer<vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
380 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
381 vtkSmartPointer<vtkVertexGlyphFilter> filter = vtkSmartPointer<vtkVertexGlyphFilter>::New();
382 vtkSmartPointer<vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
383 vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
384
385 polyData->SetPoints(data);
386 filter->SetInputData(polyData);
387 mapper->SetInputConnection(filter->GetOutputPort());
388 actor->SetMapper(mapper);
389
390 // Setup actor and mapper
391 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
392 actor->GetProperty()->SetOpacity(opacity);
393 actor->SetVisibility(1);
394 actor->GetProperty()->SetPointSize(size*5);
395
396 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
397 pVtkViewer->renderer->AddActor(actor);
398
399 squareDataMap.insert(std::pair<std::size_t,vtkSmartPointer<vtkPoints>>(hash, data));
400 squarePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyData>>(hash, polyData));
401 squareFilterMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>>(hash, filter));
402 squarePolyDataMapperMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyDataMapper>>(hash, mapper));
403 squarePolyDataActorMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkActor>>(hash, actor));
404 }
405
406 // Data data point
408 G4Point3D posPrime = rot*square.GetPosition();
409 squareDataMap[hash]->InsertNextPoint(fObjectTransformation.dx() + posPrime.x(),
410 fObjectTransformation.dy() + posPrime.y(),
411 fObjectTransformation.dz() + posPrime.z());
412 }
413 else if (sizeType == screen) {
414
415 }
416}
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > squareFilterMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > squareDataMap
std::map< std::size_t, const G4VisAttributes * > squareVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > squarePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > squarePolyDataMapperMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > squarePolyDataMap

◆ AddPrimitive() [6/6]

void G4VtkSceneHandler::AddPrimitive ( const G4Text text)
virtual

Implements G4VSceneHandler.

Definition at line 229 of file G4VtkSceneHandler.cc.

229 {
230
232 if(fProcessing2D) {sizeType = screen;}
233 else {sizeType = world;}
234
235 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(text.GetVisAttributes());
236 G4Color colour = pVA->GetColour();
237 G4double opacity = colour.GetAlpha();
238 // G4Text::Layout layout = text.GetLayout();
239 // G4double xOffset = text.GetXOffset();
240 // G4double yOffset = text.GetYOffset();
241
242 double x = text.GetPosition().x();
243 double y = text.GetPosition().y();
244 double z = text.GetPosition().z();
245
246#ifdef G4VTKDEBUG
247 G4cout << "=================================" << G4endl;
248 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> text: " << text.GetText() << " sizeType:" << sizeType << " " << fProcessing2D << G4endl;
249 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> colour: " << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
250 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> alpha: " << colour.GetAlpha() << G4endl;
251 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> position: " << x << " " << y << " " << z << G4endl;
252#endif
253
254 switch (sizeType) {
255 default:
256 case (screen): {
257 vtkSmartPointer <vtkTextActor> actor = vtkSmartPointer<vtkTextActor>::New();
258 actor->SetInput(text.GetText().c_str());
259 actor->GetPositionCoordinate()->SetCoordinateSystemToNormalizedViewport();
260 // actor->SetTextScaleModeToViewport();
261 actor->SetPosition((x+1.)/2.0, (y+1.)/2.);
262 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
263 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetBlue(), colour.GetGreen());
264 actor->GetTextProperty()->SetOpacity(opacity);
265
266 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
267 pVtkViewer->renderer->AddActor(actor);
268 break;
269 }
270 case world: {
271 vtkSmartPointer <vtkBillboardTextActor3D> actor = vtkSmartPointer<vtkBillboardTextActor3D>::New();
272 actor->SetInput(text.GetText().c_str());
273 actor->SetPosition(x, y, z);
274 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
275 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetBlue(), colour.GetGreen());
276 actor->GetTextProperty()->SetOpacity(opacity);
277
278 auto *pVtkViewer = dynamic_cast<G4VtkViewer*>(fpViewer);
279 pVtkViewer->renderer->AddActor(actor);
280 break;
281 }
282 }
283}
G4String GetText() const
G4double GetScreenSize() const

◆ AddPrimitiveBakedTransform()

void G4VtkSceneHandler::AddPrimitiveBakedTransform ( const G4Polyhedron polyhedron)

Definition at line 589 of file G4VtkSceneHandler.cc.

590{
591 // only have a single polydata object for each LV but bake in the PV
592 // transformation so clippers and cutters can be implemented
593 AddPrimitiveTensorGlyph(polyhedron);
594}

◆ AddPrimitiveTensorGlyph()

void G4VtkSceneHandler::AddPrimitiveTensorGlyph ( const G4Polyhedron polyhedron)

Definition at line 421 of file G4VtkSceneHandler.cc.

421 {
422 // only have a single polydata object for each LV but bake in the PV
423 // transformation so clippers and cutters can be implemented
424
425 //Get colour, etc..
426 if (polyhedron.GetNoFacets() == 0) {return;}
427
428 // Get vis attributes - pick up defaults if none.
430 G4Color colour = pVA->GetColour();
431 // G4bool isVisible = pVA->IsVisible();
432 // G4double lineWidth = pVA->GetLineWidth();
433 // G4VisAttributes::LineStyle lineStyle = pVA->GetLineStyle();
434
435 // G4double lineWidthScale = drawing_style.GetGlobalLineWidthScale();
436
437 // Get view parameters that the user can force through the vis attributes, thereby over-riding the current view parameter.
439 //G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
440
441#ifdef G4VTKDEBUG
442 G4cout << "=================================" << G4endl;
443 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> " << G4endl;
444 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> colour:" << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
445 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> alpha:" << colour.GetAlpha() << G4endl;
446 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> lineWidth:" << lineWidth << G4endl;
447 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> lineStyle:" << lineStyle << G4endl;
448#endif
449
450 //std::size_t vhash = std::hash<G4VisAttributes>{}(*pVA);
451
452 std::size_t vhash = 0;
453 std::hash_combine(vhash,static_cast<int>(drawing_style));
454
455 std::size_t phash = std::hash<G4Polyhedron>{}(polyhedron);
456
457 std::size_t hash = 0;
458 std::hash_combine(hash, phash);
459 std::hash_combine(hash, vhash);
460
461 if (polyhedronPolyDataMap.find(hash) == polyhedronPolyDataMap.end()) {
462
463 vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
464 vtkSmartPointer<vtkCellArray> polys = vtkSmartPointer<vtkCellArray>::New();
465 vtkSmartPointer<vtkPolyData> polydata = vtkSmartPointer<vtkPolyData>::New();
466
467 G4bool notLastFace;
468 int iVert = 0;
469 do {
470 G4Point3D vertex[4];
471 G4int edgeFlag[4];
472 G4Normal3D normals[4];
473 G4int nEdges;
474 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals);
475
476 vtkSmartPointer<vtkIdList> poly = vtkSmartPointer<vtkIdList>::New();
477 // loop over vertices
478 for (int i = 0; i < nEdges; i++) {
479 points->InsertNextPoint(vertex[i].x(), vertex[i].y(), vertex[i].z());
480 poly->InsertNextId(iVert);
481 iVert++;
482 }
483 polys->InsertNextCell(poly);
484
485 } while (notLastFace);
486
487 polydata->SetPoints(points);
488 polydata->SetPolys(polys);
489
490 polyhedronDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, points));
491 polyhedronPolyMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkCellArray>>(hash, polys));
492 polyhedronPolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash, polydata));
493 polyhedronPolyDataCountMap.insert(std::pair<std::size_t, std::size_t>(hash, 0));
494
495 vtkSmartPointer<vtkPoints> instancePosition = vtkSmartPointer<vtkPoints>::New();
496 vtkSmartPointer<vtkDoubleArray> instanceRotation = vtkSmartPointer<vtkDoubleArray>::New();
497 vtkSmartPointer<vtkDoubleArray> instanceColors = vtkSmartPointer<vtkDoubleArray>::New();
498 vtkSmartPointer<vtkPolyDataMapper> instanceMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
499 vtkSmartPointer<vtkActor> instanceActor = vtkSmartPointer<vtkActor>::New();
500 instanceColors->SetName("colors");
501
502 instanceColors->SetNumberOfComponents(4);
503 instanceRotation->SetNumberOfComponents(9);
504
505 vtkSmartPointer<vtkPolyData> instancePolyData = vtkSmartPointer<vtkPolyData>::New();
506 instancePolyData->SetPoints(instancePosition);
507 instancePolyData->GetPointData()->SetTensors(instanceRotation);
508 instancePolyData->GetPointData()->SetVectors(instanceColors);
509 instancePolyData->GetPointData()->SetScalars(instanceColors);
510
511 vtkSmartPointer<vtkCleanPolyData> filterClean = vtkSmartPointer<vtkCleanPolyData>::New();
512 filterClean->PointMergingOn();
513 filterClean->AddInputData(polydata);
514
515 vtkSmartPointer<vtkTriangleFilter> filterTriangle = vtkSmartPointer<vtkTriangleFilter>::New();
516 filterTriangle->SetInputConnection(filterClean->GetOutputPort());
517
518 vtkSmartPointer<vtkPolyDataNormals> filterNormals = vtkSmartPointer<vtkPolyDataNormals>::New();
519 filterNormals->SetFeatureAngle(45);
520 filterNormals->SetInputConnection(filterTriangle->GetOutputPort());
521
522 vtkSmartPointer<vtkFeatureEdges> filterEdge = vtkSmartPointer<vtkFeatureEdges>::New();
523 filterEdge->SetFeatureEdges(1);
524 filterEdge->SetManifoldEdges(0);
525 filterEdge->SetBoundaryEdges(0);
526 filterEdge->SetFeatureAngle(45); // TODO need to have a function and command to set this
527 filterEdge->SetInputConnection(filterTriangle->GetOutputPort());
528
529 vtkSmartPointer<vtkTensorGlyphColor> tensorGlyph = vtkSmartPointer<vtkTensorGlyphColor>::New();
530 tensorGlyph->SetInputData(instancePolyData);
531 tensorGlyph->SetSourceConnection(filterNormals->GetOutputPort());
532 tensorGlyph->ColorGlyphsOn();
533 tensorGlyph->ScalingOff();
534 tensorGlyph->ThreeGlyphsOff();
535 tensorGlyph->ExtractEigenvaluesOff();
536 tensorGlyph->SetColorModeToScalars();
537 tensorGlyph->Update();
538
539 instanceMapper->SetInputData(tensorGlyph->GetOutput());
540 instanceMapper->SetColorModeToDirectScalars();
541
542 instanceActor->SetMapper(instanceMapper);
543 // instanceActor->GetProperty()->SetLineWidth(10);
544 instanceActor->SetVisibility(1);
545
546 if(drawing_style == G4ViewParameters::hsr) {
547
548 }
549
550 if(drawing_style == G4ViewParameters::hlr) {
551 }
552
553 if(drawing_style == G4ViewParameters::wireframe) {
554 instanceActor->GetProperty()->SetRepresentationToWireframe();
555 }
556
557 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
558 pVtkViewer->renderer->AddActor(instanceActor);
559
560 instancePositionMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, instancePosition));
561 instanceRotationMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkDoubleArray>>(hash, instanceRotation));
562 instanceColoursMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkDoubleArray>>(hash, instanceColors));
563 instancePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash,instancePolyData));
564 instanceActorMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkActor>>(hash, instanceActor));
565 instanceTensorGlyphMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkTensorGlyphColor>>(hash, tensorGlyph));
566 }
567
569
570 double red = colour.GetRed();
571 double green = colour.GetGreen();
572 double blue = colour.GetBlue();
573 double alpha = colour.GetAlpha();
574
575 instanceColoursMap[hash]->InsertNextTuple4(red, green, blue, alpha);
576
577 instancePositionMap[hash]->InsertNextPoint(fObjectTransformation.dx(),
580
582
583 instanceRotationMap[hash]->InsertNextTuple9(fInvObjTrans.xx(), fInvObjTrans.xy(),fInvObjTrans.xz(),
584 fInvObjTrans.yx(), fInvObjTrans.yy(),fInvObjTrans.yz(),
585 fInvObjTrans.zx(), fInvObjTrans.zy(),fInvObjTrans.zz());
586
587}
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > instancePolyDataMap
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceColoursMap
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceRotationMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > instancePositionMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polyhedronPolyDataMap
std::map< std::size_t, std::size_t > polyhedronPolyDataCountMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > instanceActorMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polyhedronPolyMap
std::map< std::size_t, vtkSmartPointer< vtkTensorGlyphColor > > instanceTensorGlyphMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polyhedronDataMap
Transform3D inverse() const
Definition: Transform3D.cc:141
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const
void hash_combine(std::size_t)

Referenced by AddPrimitive(), and AddPrimitiveBakedTransform().

◆ AddSolid()

void G4VtkSceneHandler::AddSolid ( const G4Box box)
virtual

Reimplemented from G4VSceneHandler.

Definition at line 724 of file G4VtkSceneHandler.cc.

724 {
725
727
728#if 0
729 return;
730
731 const G4VModel* pv_model = GetModel();
732 if (!pv_model) { return ; }
733
734 G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
735 if (!pPVModel) { return ; }
736
737 //-- debug information
738 if(1) {
739 G4VPhysicalVolume *pv = pPVModel->GetCurrentPV();
741 G4cout << "name=" << box.GetName() << " volumeType=" << pv->VolumeType() << " pvName=" << pv->GetName() << " lvName=" << lv->GetName() << " multiplicity=" << pv->GetMultiplicity() << " isparametrised=" << pv->IsParameterised() << " isreplicated=" << pv->IsReplicated() << " parametrisation=" << pv->GetParameterisation() << G4endl;
742 }
743
744 if(0) {
745 G4Material *mat = pPVModel->GetCurrentMaterial();
746 G4String name = mat->GetName();
747 G4double dens = mat->GetDensity()/(g/cm3);
748 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
749 G4int depth = pPVModel->GetCurrentDepth();
750 G4cout << " name : " << box.GetName() << G4endl;
751 G4cout << " copy no.: " << copyNo << G4endl;
752 G4cout << " depth : " << depth << G4endl;
753 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
754 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
755 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
756 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
757 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
758 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
759 }
760#endif
761}
const G4String & GetName() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4String & GetName() const
Definition: G4Material.hh:172
G4VPhysicalVolume * GetCurrentPV() const
G4VPhysicalVolume * GetTopPhysicalVolume() const
G4Material * GetCurrentMaterial() const
virtual G4bool IsReplicated() const =0
virtual EVolume VolumeType() const =0
virtual G4int GetMultiplicity() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
G4VModel * GetModel() const
virtual void AddSolid(const G4Box &)
G4String GetName() const
const char * name(G4int ptype)

◆ Clear()

void G4VtkSceneHandler::Clear ( )

Definition at line 679 of file G4VtkSceneHandler.cc.

679 {
680
682 polylineDataMap.clear();
683 polylineLineMap.clear();
684 polylinePolyDataMap.clear();
687
688 for (auto& v : polylineDataMap)
689 {v.second->Reset();}
690 for (auto& v : polylineLineMap)
691 {v.second->Reset();}
692
694 circleDataMap.clear();
695 circlePolyDataMap.clear();
696 circleFilterMap.clear();
699
701 squareDataMap.clear();
702 squarePolyDataMap.clear();
703 squareFilterMap.clear();
706
707 for (auto& v : squareDataMap)
708 {v.second->Reset();}
709
711 polyhedronDataMap.clear();
712 polyhedronPolyMap.clear();
713 polyhedronPolyDataMap.clear();
715
716 instancePositionMap.clear();
717 instanceRotationMap.clear();
718 instanceColoursMap.clear();
719 instancePolyDataMap.clear();
721 instanceActorMap.clear();
722}
std::map< std::size_t, const G4VisAttributes * > polyhedronVisAttributesMap

◆ Modified()

void G4VtkSceneHandler::Modified ( )

Definition at line 595 of file G4VtkSceneHandler.cc.

595 {
596 for (auto it = polylineDataMap.begin(); it != polylineDataMap.end(); it++)
597 {
598 it->second->Modified();
599 }
600
601 for (auto it = polylineLineMap.begin(); it != polylineLineMap.end(); it++)
602 {
603 it->second->Modified();
604 }
605 for (auto it = circleDataMap.begin(); it != circleDataMap.end(); it++)
606 {
607 it->second->Modified();
608 }
609
610 for (auto it = squareDataMap.begin(); it != squareDataMap.end(); it++)
611 {
612 it->second->Modified();
613 }
614
615 for(auto it = instancePositionMap.begin(); it != instancePositionMap.end(); it++)
616 {
617 it->second->Modified();
618 }
619
620 for(auto it = instanceRotationMap.begin(); it != instanceRotationMap.end(); it++)
621 {
622 it->second->Modified();
623 }
624
625 for(auto it = instanceColoursMap.begin(); it != instanceColoursMap.end(); it++)
626 {
627 it->second->Modified();
628 }
629
630 for(auto it = instancePolyDataMap.begin(); it != instancePolyDataMap.end(); it++)
631 {
632 it->second->Modified();
633 }
634
635 for(auto it = instanceActorMap.begin(); it != instanceActorMap.end(); it++)
636 {
637 it->second->Modified();
638 }
639
640 for(auto it = instanceTensorGlyphMap.begin(); it != instanceTensorGlyphMap.end(); it++)
641 {
642 it->second->Update();
643 }
644
645#ifdef G4VTKDEBUG
646 G4cout << "G4VtkSceneHandler::Modified() polyline styles: " << polylineVisAttributesMap.size() << G4endl;
647 for (auto it = polylineLineMap.begin(); it != polylineLineMap.end(); it++)
648 {
649 G4cout << "G4VtkSceneHandler::Modified() polyline segments: "
650 << it->second->GetNumberOfCells() << G4endl;
651 }
652 G4cout << "G4VtkSceneHandler::Modified() circle styles: " << circleVisAttributesMap.size() << G4endl;
653 for (auto it = circleDataMap.begin(); it != circleDataMap.end(); it++)
654 {
655 G4cout << "G4VtkSceneHandler::Modified() circles: "
656 << it->second->GetNumberOfPoints() << G4endl;
657 }
658
659 G4cout << "G4VtkSceneHandler::Modified() square styles: " << squareVisAttributesMap.size() << G4endl;
660 for (auto it = squareDataMap.begin(); it != squareDataMap.end(); it++)
661 {
662 G4cout << "G4VtkSceneHanler::Modified() squares: "
663 << it->second->GetNumberOfPoints() << G4endl;
664 }
665
666 G4cout << "G4VtkSceneHandler::Modified() unique polyhedra: " << polyhedronDataMap.size() << G4endl;
667
668 int nPlacements = 0;
669 int nCells = 0;
670 for (auto it = polyhedronPolyDataMap.begin(); it != polyhedronPolyDataMap.end(); it++) {
671 G4cout << "G4VtkSceneHandler::Modified() polyhedronPolyData: " << it->second->GetPoints()->GetNumberOfPoints() << " " << it->second->GetPolys()->GetNumberOfCells() << " " << polyhedronPolyDataCountMap[it->first] <<G4endl;
672 nCells += it->second->GetPolys()->GetNumberOfCells()*polyhedronPolyDataCountMap[it->first];
673 nPlacements += polyhedronPolyDataCountMap[it->first];
674 }
675 G4cout << "G4VtkSceneHandler::Modified() polyhedronPolyData: " << nPlacements << " " << nCells << G4endl;
676#endif
677}

Referenced by G4VtkViewer::FinishView(), and G4VtkViewer::ShowView().

Friends And Related Function Documentation

◆ G4VtkViewer

friend class G4VtkViewer
friend

Definition at line 77 of file G4VtkSceneHandler.hh.

Member Data Documentation

◆ circleDataMap

std::map<std::size_t, vtkSmartPointer<vtkPoints> > G4VtkSceneHandler::circleDataMap
protected

Definition at line 134 of file G4VtkSceneHandler.hh.

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

◆ circleFilterMap

std::map<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter> > G4VtkSceneHandler::circleFilterMap
protected

Definition at line 136 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ circlePolyDataActorMap

std::map<std::size_t, vtkSmartPointer<vtkActor> > G4VtkSceneHandler::circlePolyDataActorMap
protected

Definition at line 138 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ circlePolyDataMap

std::map<std::size_t, vtkSmartPointer<vtkPolyData> > G4VtkSceneHandler::circlePolyDataMap
protected

Definition at line 135 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ circlePolyDataMapperMap

std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper> > G4VtkSceneHandler::circlePolyDataMapperMap
protected

Definition at line 137 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ circleVisAttributesMap

std::map<std::size_t, const G4VisAttributes*> G4VtkSceneHandler::circleVisAttributesMap
protected

Definition at line 133 of file G4VtkSceneHandler.hh.

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

◆ fSceneIdCount

G4int G4VtkSceneHandler::fSceneIdCount = 0
staticprotected

Definition at line 122 of file G4VtkSceneHandler.hh.

◆ instanceActorMap

std::map<std::size_t, vtkSmartPointer<vtkActor> > G4VtkSceneHandler::instanceActorMap
protected

Definition at line 161 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ instanceColoursMap

std::map<std::size_t, vtkSmartPointer<vtkDoubleArray> > G4VtkSceneHandler::instanceColoursMap
protected

Definition at line 158 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ instancePolyDataMap

std::map<std::size_t, vtkSmartPointer<vtkPolyData> > G4VtkSceneHandler::instancePolyDataMap
protected

Definition at line 159 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ instancePositionMap

std::map<std::size_t, vtkSmartPointer<vtkPoints> > G4VtkSceneHandler::instancePositionMap
protected

Definition at line 156 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ instanceRotationMap

std::map<std::size_t, vtkSmartPointer<vtkDoubleArray> > G4VtkSceneHandler::instanceRotationMap
protected

Definition at line 157 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ instanceTensorGlyphMap

std::map<std::size_t, vtkSmartPointer<vtkTensorGlyphColor> > G4VtkSceneHandler::instanceTensorGlyphMap
protected

Definition at line 160 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ polyhedronDataMap

std::map<std::size_t, vtkSmartPointer<vtkPoints> > G4VtkSceneHandler::polyhedronDataMap
protected

Definition at line 150 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ polyhedronPolyDataCountMap

std::map<std::size_t, std::size_t> G4VtkSceneHandler::polyhedronPolyDataCountMap
protected

Definition at line 153 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ polyhedronPolyDataMap

std::map<std::size_t, vtkSmartPointer<vtkPolyData> > G4VtkSceneHandler::polyhedronPolyDataMap
protected

Definition at line 152 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), Clear(), and Modified().

◆ polyhedronPolyMap

std::map<std::size_t, vtkSmartPointer<vtkCellArray> > G4VtkSceneHandler::polyhedronPolyMap
protected

Definition at line 151 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitiveTensorGlyph(), and Clear().

◆ polyhedronVisAttributesMap

std::map<std::size_t, const G4VisAttributes*> G4VtkSceneHandler::polyhedronVisAttributesMap
protected

Definition at line 149 of file G4VtkSceneHandler.hh.

Referenced by Clear().

◆ polylineDataMap

std::map<std::size_t, vtkSmartPointer<vtkPoints> > G4VtkSceneHandler::polylineDataMap
protected

Definition at line 126 of file G4VtkSceneHandler.hh.

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

◆ polylineLineMap

std::map<std::size_t, vtkSmartPointer<vtkCellArray> > G4VtkSceneHandler::polylineLineMap
protected

Definition at line 127 of file G4VtkSceneHandler.hh.

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

◆ polylinePolyDataActorMap

std::map<std::size_t, vtkSmartPointer<vtkActor> > G4VtkSceneHandler::polylinePolyDataActorMap
protected

Definition at line 130 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ polylinePolyDataMap

std::map<std::size_t, vtkSmartPointer<vtkPolyData> > G4VtkSceneHandler::polylinePolyDataMap
protected

Definition at line 128 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ polylinePolyDataMapperMap

std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper> > G4VtkSceneHandler::polylinePolyDataMapperMap
protected

Definition at line 129 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ polylineVisAttributesMap

std::map<std::size_t, const G4VisAttributes*> G4VtkSceneHandler::polylineVisAttributesMap
protected

Definition at line 125 of file G4VtkSceneHandler.hh.

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

◆ squareDataMap

std::map<std::size_t, vtkSmartPointer<vtkPoints> > G4VtkSceneHandler::squareDataMap
protected

Definition at line 142 of file G4VtkSceneHandler.hh.

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

◆ squareFilterMap

std::map<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter> > G4VtkSceneHandler::squareFilterMap
protected

Definition at line 144 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ squarePolyDataActorMap

std::map<std::size_t, vtkSmartPointer<vtkActor> > G4VtkSceneHandler::squarePolyDataActorMap
protected

Definition at line 146 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ squarePolyDataMap

std::map<std::size_t, vtkSmartPointer<vtkPolyData> > G4VtkSceneHandler::squarePolyDataMap
protected

Definition at line 143 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ squarePolyDataMapperMap

std::map<std::size_t, vtkSmartPointer<vtkPolyDataMapper> > G4VtkSceneHandler::squarePolyDataMapperMap
protected

Definition at line 145 of file G4VtkSceneHandler.hh.

Referenced by AddPrimitive(), and Clear().

◆ squareVisAttributesMap

std::map<std::size_t, const G4VisAttributes*> G4VtkSceneHandler::squareVisAttributesMap
protected

Definition at line 141 of file G4VtkSceneHandler.hh.

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


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