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

#include <G4VtkCutterPipeline.hh>

+ Inheritance diagram for G4VtkCutterPipeline:

Public Member Functions

 G4VtkCutterPipeline (G4String name, const G4VtkVisContext &vc, vtkSmartPointer< vtkPolyDataAlgorithm > filter, G4bool useVcColour=false)
 
 ~G4VtkCutterPipeline () override=default
 
void SetPlane (const G4Plane3D &plane)
 
void SetPlane (G4double x, G4double y, G4double z, G4double nx, G4double ny, G4double nz)
 
void TransformPlane (G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22)
 
void Enable () override
 
void Disable () override
 
void Print () override
 
void Modified () override
 
void Clear () override
 
vtkSmartPointer< vtkActor > GetActor ()
 
- Public Member Functions inherited from G4VVtkPipeline
 G4VVtkPipeline ()
 
 G4VVtkPipeline (G4String nameIn, G4String typeIn, const G4VtkVisContext &vcIn, G4bool disableParentIn=false, vtkSmartPointer< vtkRenderer > rendererIn=nullptr)
 
virtual ~G4VVtkPipeline ()
 
void SetDisableParent (G4bool disableParentIn)
 
bool GetDisableParent ()
 
void SetName (G4String nameIn)
 
const G4String GetName ()
 
void SetTypeName (G4String typeNameIn)
 
G4String GetTypeName ()
 
G4VtkVisContextGetVtkVisContext ()
 
void AddChildPipeline (G4VVtkPipeline *child)
 
G4VVtkPipelineGetChildPipeline (G4int iPipeline)
 
G4int GetNumberOfChildPipelines ()
 
std::vector< G4VVtkPipeline * > GetChildPipelines ()
 
void ClearChildPipeline ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VVtkPipeline
G4String name
 
G4String type
 
G4bool disableParent
 
std::vector< G4VVtkPipeline * > childPipelines
 
vtkSmartPointer< vtkRenderer > renderer
 
G4VtkVisContext vc
 

Detailed Description

Definition at line 39 of file G4VtkCutterPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkCutterPipeline()

G4VtkCutterPipeline::G4VtkCutterPipeline ( G4String name,
const G4VtkVisContext & vc,
vtkSmartPointer< vtkPolyDataAlgorithm > filter,
G4bool useVcColour = false )

Definition at line 40 of file G4VtkCutterPipeline.cc.

43 : G4VVtkPipeline(nameIn, G4String("G4VtkCutterPipeline"), vcIn, false, vcIn.fViewer->renderer)
44{
45 // cutter plane
47 plane->SetOrigin(0, 0, 0);
48 plane->SetNormal(0, 1, 0);
49
50 // cutter
52 cutter->SetCutFunction(plane);
53 cutter->SetInputConnection(filter->GetOutputPort());
54
55 // mapper
57 mapper->SetInputConnection(cutter->GetOutputPort());
58 mapper->SetColorModeToDirectScalars();
59 mapper->ScalarVisibilityOn();
60
61 // add to actor
63 actor->SetMapper(mapper);
64 actor->SetVisibility(1);
65
66 // want the actor not to have shading
67 actor->GetProperty()->SetAmbient(1);
68 actor->GetProperty()->SetDiffuse(1);
69 actor->GetProperty()->SetSpecular(0.0);
70 actor->GetProperty()->SetSpecularPower(0);
71
72 actor->GetProperty()->SetLineWidth(10);
73
74 // colour parameters
75 if (useVcColour) {
76 actor->GetProperty()->SetOpacity(vc.alpha);
77 actor->GetProperty()->SetColor(vc.red, vc.green, vc.blue);
78 }
79
80 // add to renderer
81 vc.fViewer->renderer->AddActor(actor);
82}
G4VtkVisContext vc
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkCutterPipeline()

G4VtkCutterPipeline::~G4VtkCutterPipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkCutterPipeline::Clear ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 57 of file G4VtkCutterPipeline.hh.

58 {
59 if (renderer != nullptr) {
60 renderer->RemoveActor(actor);
61 }
63 }
virtual void Clear()
vtkSmartPointer< vtkRenderer > renderer

◆ Disable()

void G4VtkCutterPipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 115 of file G4VtkCutterPipeline.cc.

116{
117 actor->SetVisibility(0);
118}

◆ Enable()

void G4VtkCutterPipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 111 of file G4VtkCutterPipeline.cc.

112{
113 actor->SetVisibility(1);
114}

◆ GetActor()

vtkSmartPointer< vtkActor > G4VtkCutterPipeline::GetActor ( )
inline

Definition at line 65 of file G4VtkCutterPipeline.hh.

65{ return actor; }

◆ Modified()

void G4VtkCutterPipeline::Modified ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 56 of file G4VtkCutterPipeline.hh.

virtual void Modified()

Referenced by SetPlane().

◆ Print()

void G4VtkCutterPipeline::Print ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 55 of file G4VtkCutterPipeline.hh.

virtual void Print()

◆ SetPlane() [1/2]

void G4VtkCutterPipeline::SetPlane ( const G4Plane3D & plane)

Definition at line 84 of file G4VtkCutterPipeline.cc.

85{
86 auto normal = planeIn.normal();
87 auto point = planeIn.point();
88 this->SetPlane(point.x(), point.y(), point.z(), normal.x(), normal.y(), normal.z());
89}
void SetPlane(const G4Plane3D &plane)

Referenced by SetPlane(), TransformPlane(), and G4VtkStore::UpdateCutter().

◆ SetPlane() [2/2]

void G4VtkCutterPipeline::SetPlane ( G4double x,
G4double y,
G4double z,
G4double nx,
G4double ny,
G4double nz )

Definition at line 91 of file G4VtkCutterPipeline.cc.

93{
94 plane->SetOrigin(x, y, z);
95 plane->SetNormal(nx, ny, nz);
96 Modified();
97}

◆ TransformPlane()

void G4VtkCutterPipeline::TransformPlane ( G4double dx,
G4double dy,
G4double dz,
G4double r00,
G4double r01,
G4double r02,
G4double r10,
G4double r11,
G4double r12,
G4double r20,
G4double r21,
G4double r22 )

Definition at line 99 of file G4VtkCutterPipeline.cc.

102{
103 auto o = plane->GetOrigin();
104 auto n = plane->GetNormal();
105
106 SetPlane(r00 * o[0] + r01 * o[1] + r02 * o[2] + dx, r10 * o[0] + r11 * o[1] + r12 * o[2] + dy,
107 r20 * o[0] + r21 * o[1] + r22 * o[2] + dz, r00 * n[0] + r01 * n[1] + r02 * n[2],
108 r10 * n[0] + r11 * n[1] + r12 * n[2], r20 * n[0] + r21 * n[1] + r22 * n[2]);
109}

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