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

#include <G4VtkClipOpenPipeline.hh>

+ Inheritance diagram for G4VtkClipOpenPipeline:

Public Member Functions

 G4VtkClipOpenPipeline (G4String name, const G4VtkVisContext &vc, vtkSmartPointer< vtkPolyDataAlgorithm > filter, G4bool useVcColour=false)
 
 ~G4VtkClipOpenPipeline () 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 41 of file G4VtkClipOpenPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkClipOpenPipeline()

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

Definition at line 41 of file G4VtkClipOpenPipeline.cc.

44 : G4VVtkPipeline(nameIn, G4String("G4VtkClipOpenPipeline"), vcIn, true, vcIn.fViewer->renderer)
45{
46 // create implicit function for clipping
48 plane->SetOrigin(0, 0, 0);
49 plane->SetNormal(-1, 0, 0);
50
51 // clipper
53 clipper->SetClipFunction(plane);
54 clipper->SetInputConnection(filter->GetOutputPort());
55 // clipper->SetScalarModeToColors();
56 // clipper->PassPointDataOn();
57
58 // calculate normals with a feature angle of 45 degrees
59 auto filterNormals = vtkSmartPointer<vtkPolyDataNormals>::New();
60 filterNormals->SetFeatureAngle(45);
61 filterNormals->SetInputConnection(clipper->GetOutputPort());
62
63 // mapper
65 mapper->SetInputConnection(filterNormals->GetOutputPort());
66 mapper->SetColorModeToDirectScalars();
67 mapper->ScalarVisibilityOn();
68
69 // add to actor
71 actor->SetMapper(mapper);
72 actor->SetVisibility(1);
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 // set actor properties from vis context
82 actor->GetProperty()->SetRepresentationToSurface();
83 }
85 actor->GetProperty()->SetRepresentationToSurface();
86 }
88 actor->GetProperty()->SetRepresentationToWireframe();
89 }
90
91 // shading parameters
92 actor->GetProperty()->SetAmbient(0.2);
93 actor->GetProperty()->SetDiffuse(0.7);
94 actor->GetProperty()->SetSpecular(0.1);
95 actor->GetProperty()->SetSpecularPower(1);
96
97 // add to renderer
98 vc.fViewer->renderer->AddActor(actor);
99}
G4VtkVisContext vc
vtkNew< vtkRenderer > renderer
G4ViewParameters::DrawingStyle fDrawingStyle
const G4VtkViewer * fViewer

◆ ~G4VtkClipOpenPipeline()

G4VtkClipOpenPipeline::~G4VtkClipOpenPipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkClipOpenPipeline::Clear ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 59 of file G4VtkClipOpenPipeline.hh.

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

◆ Disable()

void G4VtkClipOpenPipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 132 of file G4VtkClipOpenPipeline.cc.

133{
134 actor->SetVisibility(0);
135}

◆ Enable()

void G4VtkClipOpenPipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 128 of file G4VtkClipOpenPipeline.cc.

129{
130 actor->SetVisibility(1);
131}

◆ GetActor()

vtkSmartPointer< vtkActor > G4VtkClipOpenPipeline::GetActor ( )
inline

Definition at line 67 of file G4VtkClipOpenPipeline.hh.

67{ return actor; }

◆ Modified()

void G4VtkClipOpenPipeline::Modified ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 58 of file G4VtkClipOpenPipeline.hh.

virtual void Modified()

Referenced by SetPlane().

◆ Print()

void G4VtkClipOpenPipeline::Print ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 57 of file G4VtkClipOpenPipeline.hh.

virtual void Print()

◆ SetPlane() [1/2]

void G4VtkClipOpenPipeline::SetPlane ( const G4Plane3D & plane)

Definition at line 101 of file G4VtkClipOpenPipeline.cc.

102{
103 auto normal = planeIn.normal();
104 auto point = planeIn.point();
105 this->SetPlane(point.x(), point.y(), point.z(), normal.x(), normal.y(), normal.z());
106}
void SetPlane(const G4Plane3D &plane)

Referenced by SetPlane(), and TransformPlane().

◆ SetPlane() [2/2]

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

Definition at line 108 of file G4VtkClipOpenPipeline.cc.

110{
111 plane->SetOrigin(x, y, z);
112 plane->SetNormal(nx, ny, nz);
113 Modified();
114}

◆ TransformPlane()

void G4VtkClipOpenPipeline::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 116 of file G4VtkClipOpenPipeline.cc.

119{
120 auto o = plane->GetOrigin();
121 auto n = plane->GetNormal();
122
123 SetPlane(r00 * o[0] + r01 * o[1] + r02 * o[2] + dx, r10 * o[0] + r11 * o[1] + r12 * o[2] + dy,
124 r20 * o[0] + r21 * o[1] + r22 * o[2] + dz, r00 * n[0] + r01 * n[1] + r02 * n[2],
125 r10 * n[0] + r11 * n[1] + r12 * n[2], r20 * n[0] + r21 * n[1] + r22 * n[2]);
126}

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