Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VtkClipOpenPipeline.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25
27
28#include "G4VtkViewer.hh"
29#include "G4VtkVisContext.hh"
30
31#include <vtkActor.h>
32#include <vtkClipPolyData.h>
33#include <vtkPlane.h>
34#include <vtkPolyDataAlgorithm.h>
35#include <vtkPolyDataMapper.h>
36#include <vtkPolyDataNormals.h>
37#include <vtkProperty.h>
38#include <vtkSmartPointer.h>
39#include <vtkClipClosedSurface.h>
40
43 G4bool useVcColour)
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}
100
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}
107
109 G4double nz)
110{
111 plane->SetOrigin(x, y, z);
112 plane->SetNormal(nx, ny, nz);
113 Modified();
114}
115
117 G4double r01, G4double r02, G4double r10, G4double r11,
118 G4double r12, G4double r20, G4double r21, G4double r22)
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}
127
129{
130 actor->SetVisibility(1);
131}
133{
134 actor->SetVisibility(0);
135}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4VtkVisContext vc
G4VtkClipOpenPipeline(G4String name, const G4VtkVisContext &vc, vtkSmartPointer< vtkPolyDataAlgorithm > filter, G4bool useVcColour=false)
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 SetPlane(const G4Plane3D &plane)
vtkNew< vtkRenderer > renderer
G4ViewParameters::DrawingStyle fDrawingStyle
const G4VtkViewer * fViewer
Point3D< T > point(const Point3D< T > &p) const
Definition Plane3D.h:115
Normal3D< T > normal() const
Definition Plane3D.h:97