Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkClipClosedSurfacePipeline.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 <vtkPlaneCollection.h>
35#include <vtkPolyDataAlgorithm.h>
36#include <vtkPolyDataMapper.h>
37#include <vtkPolyDataNormals.h>
38#include <vtkProperty.h>
39#include <vtkSmartPointer.h>
40#include <vtkClipClosedSurface.h>
41
44 G4bool useVcColour)
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}
107
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}
114
116 G4double ny, G4double nz)
117{
118 plane->SetOrigin(x, y, z);
119 plane->SetNormal(nx, ny, nz);
120 this->Modified();
121}
122
124 G4double r00, G4double r01, G4double r02,
125 G4double r10, G4double r11, G4double r12,
126 G4double r20, G4double r21, G4double r22)
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}
135
137{
138 actor->SetVisibility(1);
139}
141{
142 actor->SetVisibility(0);
143}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4VtkVisContext vc
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)
G4VtkClipClosedSurfacePipeline(G4String name, const G4VtkVisContext &vc, vtkSmartPointer< vtkPolyDataAlgorithm > filter, G4bool useVcColour=false)
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