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

#include <G4VtkClipClosedSurfacePipeline.hh>

+ Inheritance diagram for G4VtkClipClosedSurfacePipeline:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4VtkClipClosedSurfacePipeline()

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

Definition at line 42 of file G4VtkClipClosedSurfacePipeline.cc.

45 : G4VVtkPipeline(nameIn, G4String("G4VtkClipClosedSurfacePipeline"), vcIn, true,
46 vcIn.fViewer->renderer)
47{
48 // create implicit function for clipping
50 plane->SetOrigin(0, 0, 0);
51 plane->SetNormal(0, 1, 0);
52
53 vtkNew<vtkPlaneCollection> planes;
54 planes->AddItem(plane);
55
56 // clipper
58 clipper->SetClippingPlanes(planes);
59 clipper->SetInputConnection(filter->GetOutputPort());
60 // clipper->SetScalarModeToColors();
61 clipper->PassPointDataOn();
62
63 // calculate normals with a feature angle of 45 degrees
64 auto filterNormals = vtkSmartPointer<vtkPolyDataNormals>::New();
65 filterNormals->SetFeatureAngle(45);
66 filterNormals->SetInputConnection(clipper->GetOutputPort());
67
68 // mapper
70 mapper->SetInputConnection(filterNormals->GetOutputPort());
71 mapper->SetColorModeToDirectScalars();
72 mapper->ScalarVisibilityOn();
73 mapper->SetResolveCoincidentTopologyToPolygonOffset();
74 mapper->SetRelativeCoincidentTopologyPolygonOffsetParameters(1, -vc.fDepth);
75
76 // add to actor
78 actor->SetMapper(mapper);
79 actor->SetVisibility(1);
80
81 // colour parameters
82 if (useVcColour) {
83 actor->GetProperty()->SetOpacity(vc.alpha);
84 actor->GetProperty()->SetColor(vc.red, vc.green, vc.blue);
85 }
86
87 // set actor properties from vis context
89 actor->GetProperty()->SetRepresentationToSurface();
90 }
92 actor->GetProperty()->SetRepresentationToSurface();
93 }
95 actor->GetProperty()->SetRepresentationToWireframe();
96 }
97
98 // shading parameters
99 actor->GetProperty()->SetAmbient(0.2);
100 actor->GetProperty()->SetDiffuse(0.7);
101 actor->GetProperty()->SetSpecular(0.1);
102 actor->GetProperty()->SetSpecularPower(1);
103
104 // add to renderer
105 vc.fViewer->renderer->AddActor(actor);
106}
G4VtkVisContext vc
vtkNew< vtkRenderer > renderer
G4ViewParameters::DrawingStyle fDrawingStyle
const G4VtkViewer * fViewer

◆ ~G4VtkClipClosedSurfacePipeline()

G4VtkClipClosedSurfacePipeline::~G4VtkClipClosedSurfacePipeline ( )
overridedefault

Member Function Documentation

◆ Clear()

void G4VtkClipClosedSurfacePipeline::Clear ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 60 of file G4VtkClipClosedSurfacePipeline.hh.

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

◆ Disable()

void G4VtkClipClosedSurfacePipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 140 of file G4VtkClipClosedSurfacePipeline.cc.

141{
142 actor->SetVisibility(0);
143}

◆ Enable()

void G4VtkClipClosedSurfacePipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 136 of file G4VtkClipClosedSurfacePipeline.cc.

137{
138 actor->SetVisibility(1);
139}

◆ GetActor()

vtkSmartPointer< vtkActor > G4VtkClipClosedSurfacePipeline::GetActor ( )
inline

Definition at line 68 of file G4VtkClipClosedSurfacePipeline.hh.

68{ return actor; }

◆ Modified()

void G4VtkClipClosedSurfacePipeline::Modified ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 59 of file G4VtkClipClosedSurfacePipeline.hh.

virtual void Modified()

Referenced by SetPlane().

◆ Print()

void G4VtkClipClosedSurfacePipeline::Print ( )
inlineoverridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 58 of file G4VtkClipClosedSurfacePipeline.hh.

virtual void Print()

◆ SetPlane() [1/2]

void G4VtkClipClosedSurfacePipeline::SetPlane ( const G4Plane3D & plane)

Definition at line 108 of file G4VtkClipClosedSurfacePipeline.cc.

109{
110 auto normal = planeIn.normal();
111 auto point = planeIn.point();
112 this->SetPlane(point.x(), point.y(), point.z(), normal.x(), normal.y(), normal.z());
113}

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

◆ SetPlane() [2/2]

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

Definition at line 115 of file G4VtkClipClosedSurfacePipeline.cc.

117{
118 plane->SetOrigin(x, y, z);
119 plane->SetNormal(nx, ny, nz);
120 this->Modified();
121}

◆ TransformPlane()

void G4VtkClipClosedSurfacePipeline::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 123 of file G4VtkClipClosedSurfacePipeline.cc.

127{
128 auto o = plane->GetOrigin();
129 auto n = plane->GetNormal();
130
131 SetPlane(r00 * o[0] + r01 * o[1] + r02 * o[2] + dx, r10 * o[0] + r11 * o[1] + r12 * o[2] + dy,
132 r20 * o[0] + r21 * o[1] + r22 * o[2] + dz, r00 * n[0] + r01 * n[1] + r02 * n[2],
133 r10 * n[0] + r11 * n[1] + r12 * n[2], r20 * n[0] + r21 * n[1] + r22 * n[2]);
134}

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