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

#include <G4VtkImagePipeline.hh>

+ Inheritance diagram for G4VtkImagePipeline:

Public Member Functions

 G4VtkImagePipeline (const G4String &name, const G4VtkVisContext &vc)
 
 ~G4VtkImagePipeline () override=default
 
void AddFilter (vtkSmartPointer< vtkImageAlgorithm > f)
 
vtkSmartPointer< vtkImageAlgorithm > GetFilter (G4int iFilter)
 
vtkSmartPointer< vtkImageAlgorithm > GetFinalFilter ()
 
void AddChildPipeline (G4VtkImagePipeline *child)
 
G4VtkImagePipelineGetChildPipeline (G4int iPipeline)
 
virtual void SetImage (const G4String &fileName)
 
virtual void SetTransformation (const G4Transform3D &transformation)
 
void Print () override
 
void Modified () override
 
void Clear () override
 
virtual vtkSmartPointer< vtkImageActor > GetActor ()
 
void Disable () override
 
void Enable () override
 

Detailed Description

Definition at line 46 of file G4VtkImagePipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkImagePipeline()

G4VtkImagePipeline::G4VtkImagePipeline ( const G4String & name,
const G4VtkVisContext & vc )

Definition at line 44 of file G4VtkImagePipeline.cc.

45 : G4VVtkPipeline(nameIn, "G4VtkImagePipeline", vcIn, false, vcIn.fViewer->renderer)
46{
48
49 actor->GetProperty()->SetOpacity(vc.alpha);
50
51 auto transform = vtkSmartPointer<vtkMatrix4x4>::New();
52 double transformArray[16] = {vc.fTransform.xx(),
64 0.,
65 0.,
66 0.,
67 1.};
68 transform->DeepCopy(transformArray);
69 actor->SetUserMatrix(transform);
70 actor->GetProperty()->SetOpacity(vc.alpha);
71
72 vc.fViewer->renderer->AddActor(GetActor());
73}
G4VtkVisContext vc
virtual vtkSmartPointer< vtkImageActor > GetActor()
vtkNew< vtkRenderer > renderer
const G4Transform3D & fTransform
const G4VtkViewer * fViewer
double dy() const
double zz() const
double yz() const
double dz() const
double dx() const
double xy() const
double zx() const
double yx() const
double zy() const
double xx() const
double yy() const
double xz() const

◆ ~G4VtkImagePipeline()

G4VtkImagePipeline::~G4VtkImagePipeline ( )
overridedefault

Member Function Documentation

◆ AddChildPipeline()

void G4VtkImagePipeline::AddChildPipeline ( G4VtkImagePipeline * child)
inline

Definition at line 55 of file G4VtkImagePipeline.hh.

56 {
57 childPipelines.push_back(child);
58 if (child->GetDisableParent()) {
59 Disable();
60 }
61 }

◆ AddFilter()

void G4VtkImagePipeline::AddFilter ( vtkSmartPointer< vtkImageAlgorithm > f)
inline

Definition at line 51 of file G4VtkImagePipeline.hh.

51{ filters.push_back(f); }

Referenced by SetImage().

◆ Clear()

void G4VtkImagePipeline::Clear ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 130 of file G4VtkImagePipeline.cc.

131{
132 if (actor != nullptr) renderer->RemoveActor(actor);
133
134 for (auto c : childPipelines)
135 c->Clear();
136}
vtkSmartPointer< vtkRenderer > renderer

◆ Disable()

void G4VtkImagePipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 138 of file G4VtkImagePipeline.cc.

139{
140 GetActor()->SetVisibility(0);
141}

Referenced by AddChildPipeline().

◆ Enable()

void G4VtkImagePipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 143 of file G4VtkImagePipeline.cc.

144{
145 GetActor()->SetVisibility(1);
146}

◆ GetActor()

virtual vtkSmartPointer< vtkImageActor > G4VtkImagePipeline::GetActor ( )
inlinevirtual

Definition at line 71 of file G4VtkImagePipeline.hh.

71{ return actor; }

Referenced by Disable(), Enable(), and G4VtkImagePipeline().

◆ GetChildPipeline()

G4VtkImagePipeline * G4VtkImagePipeline::GetChildPipeline ( G4int iPipeline)
inline

Definition at line 62 of file G4VtkImagePipeline.hh.

62{ return childPipelines[iPipeline]; }

◆ GetFilter()

vtkSmartPointer< vtkImageAlgorithm > G4VtkImagePipeline::GetFilter ( G4int iFilter)
inline

Definition at line 52 of file G4VtkImagePipeline.hh.

52{ return filters[iFilter]; }

◆ GetFinalFilter()

vtkSmartPointer< vtkImageAlgorithm > G4VtkImagePipeline::GetFinalFilter ( )
inline

Definition at line 53 of file G4VtkImagePipeline.hh.

53{ return *filters.end(); }

◆ Modified()

void G4VtkImagePipeline::Modified ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 122 of file G4VtkImagePipeline.cc.

123{
124 actor->Update();
125
126 for (auto c : childPipelines)
127 c->Modified();
128}

◆ Print()

void G4VtkImagePipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 116 of file G4VtkImagePipeline.cc.

117{
118 for (auto c : childPipelines)
119 c->Print();
120}

◆ SetImage()

void G4VtkImagePipeline::SetImage ( const G4String & fileName)
virtual

Definition at line 75 of file G4VtkImagePipeline.cc.

76{
77 vtkNew<vtkImageReader2Factory> readerFactory;
79 imageReader.TakeReference(readerFactory->CreateImageReader2(fileName.c_str()));
80 imageReader->SetFileName(fileName.c_str());
81 imageReader->Update();
82
83 vtkNew<vtkImageFlip> flip;
84 flip->SetInputConnection(imageReader->GetOutputPort());
85 flip->SetFilteredAxis(1);
86
87 actor->GetMapper()->SetInputConnection(flip->GetOutputPort());
88 AddFilter(imageReader);
89}
void AddFilter(vtkSmartPointer< vtkImageAlgorithm > f)

◆ SetTransformation()

void G4VtkImagePipeline::SetTransformation ( const G4Transform3D & transformation)
virtual

Definition at line 91 of file G4VtkImagePipeline.cc.

92{
94 t->SetElement(0, 0, transformation.xx());
95 t->SetElement(0, 1, transformation.xy());
96 t->SetElement(0, 2, transformation.xz());
97 t->SetElement(0, 3, transformation.dx());
98
99 t->SetElement(1, 0, transformation.yx());
100 t->SetElement(1, 1, transformation.yy());
101 t->SetElement(1, 2, transformation.yz());
102 t->SetElement(1, 3, transformation.dy());
103
104 t->SetElement(2, 0, transformation.zx());
105 t->SetElement(2, 1, transformation.zy());
106 t->SetElement(2, 2, transformation.zz());
107 t->SetElement(2, 3, transformation.dz());
108
109 t->SetElement(3, 0, 0);
110 t->SetElement(3, 1, 0);
111 t->SetElement(3, 2, 0);
112 t->SetElement(3, 3, 1);
113 actor->SetUserMatrix(t);
114}

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