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

#include <G4VtkViewer.hh>

+ Inheritance diagram for G4VtkViewer:

Public Member Functions

 G4VtkViewer (G4VSceneHandler &, const G4String &name)
 
void Initialise ()
 
virtual ~G4VtkViewer ()
 
void SetView ()
 
void ClearView ()
 
void DrawView ()
 
void ShowView ()
 
void FinishView ()
 
void DrawViewHUD ()
 
void DrawShadows ()
 
void ExportScreenShot (G4String, G4String)
 
void ExportOBJScene (G4String)
 
void ExportVRMLScene (G4String)
 
void ExportVTPScene (G4String)
 
void ExportView ()
 
void SetGeant4View ()
 
- Public Member Functions inherited from G4VViewer
 G4VViewer (G4VSceneHandler &, G4int id, const G4String &name="")
 
virtual ~G4VViewer ()
 
virtual void Initialise ()
 
virtual void ResetView ()
 
virtual void SetView ()=0
 
virtual void ClearView ()=0
 
virtual void DrawView ()=0
 
void RefreshView ()
 
virtual void ShowView ()
 
virtual void FinishView ()
 
std::vector< G4ThreeVectorComputeFlyThrough (G4Vector3D *)
 
const G4StringGetName () const
 
const G4StringGetShortName () const
 
void SetName (const G4String &)
 
G4int GetViewId () const
 
G4VSceneHandlerGetSceneHandler () const
 
const G4ViewParametersGetViewParameters () const
 
const G4ViewParametersGetDefaultViewParameters () const
 
G4double GetKernelVisitElapsedTimeSeconds () const
 
virtual const std::vector< G4ModelingParameters::VisAttributesModifier > * GetPrivateVisAttributesModifiers () const
 
void SetViewParameters (const G4ViewParameters &vp)
 
void SetDefaultViewParameters (const G4ViewParameters &vp)
 
const G4VisAttributesGetApplicableVisAttributes (const G4VisAttributes *) const
 
void SetNeedKernelVisit (G4bool need)
 
void NeedKernelVisit ()
 
void ProcessView ()
 

Public Attributes

vtkNew< vtkTextActor > infoTextActor
 
vtkNew< vtkInfoCallbackinfoCallback
 
vtkNew< vtkGeant4Callbackgeant4Callback
 
vtkSmartPointer< vtkLight > light
 
vtkNew< vtkCamera > camera
 
vtkNew< vtkRenderer > renderer
 
vtkRenderWindow * _renderWindow
 
vtkRenderWindowInteractor * renderWindowInteractor
 

Additional Inherited Members

- Protected Member Functions inherited from G4VViewer
void SetTouchable (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath)
 
void TouchableSetVisibility (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, G4bool visibility)
 
void TouchableSetColour (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, const G4Colour &)
 
- Protected Attributes inherited from G4VViewer
G4VSceneHandlerfSceneHandler
 
G4int fViewId
 
G4String fName
 
G4String fShortName
 
G4ViewParameters fVP
 
G4ViewParameters fDefaultVP
 
G4double fKernelVisitElapsedTimeSeconds = 999.
 
G4bool fNeedKernelVisit
 

Detailed Description

Definition at line 158 of file G4VtkViewer.hh.

Constructor & Destructor Documentation

◆ G4VtkViewer()

G4VtkViewer::G4VtkViewer ( G4VSceneHandler sceneHandler,
const G4String name 
)

Definition at line 55 of file G4VtkViewer.cc.

56 : G4VViewer(sceneHandler, sceneHandler.IncrementViewCount(), name)
57{
58 vtkObject::GlobalWarningDisplayOff();
59
60 // Set default and current view parameters
61 fVP.SetAutoRefresh(true);
63}
G4int IncrementViewCount()
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:221
G4ViewParameters fVP
Definition: G4VViewer.hh:220
void SetAutoRefresh(G4bool)

◆ ~G4VtkViewer()

G4VtkViewer::~G4VtkViewer ( )
virtual

Definition at line 123 of file G4VtkViewer.cc.

123{}

Member Function Documentation

◆ ClearView()

void G4VtkViewer::ClearView ( void  )
virtual

Implements G4VViewer.

Definition at line 219 of file G4VtkViewer.cc.

219 {
220 vtkActorCollection *actors = renderer->GetActors();
221 vtkActor *actor = actors->GetLastActor();
222
223 while(actor) {
224#ifdef G4VTKDEBUG
225 G4cout << "G4VtkViewer::ClearView() remove actor " << actor << G4endl;
226#endif
227 renderer->RemoveActor(actor);
228 actor = actors->GetLastActor();
229 }
230
231 vtkPropCollection *props = renderer->GetViewProps();
232 vtkProp *prop = props->GetLastProp();
233
234 while(prop) {
235#ifdef G4VTKDEBUG
236 G4cout << "G4VtkViewer::ClearView() remove prop " << prop << G4endl;
237#endif
238 renderer->RemoveViewProp(prop);
239 prop = props->GetLastProp();
240 }
241}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
vtkNew< vtkRenderer > renderer
Definition: G4VtkViewer.hh:186

◆ DrawShadows()

void G4VtkViewer::DrawShadows ( )

Definition at line 276 of file G4VtkViewer.cc.

277{
278 _renderWindow->SetMultiSamples(0);
279
280 vtkNew<vtkShadowMapPass> shadows;
281 vtkNew<vtkSequencePass> seq;
282
283 vtkNew<vtkRenderPassCollection> passes;
284 passes->AddItem(shadows->GetShadowMapBakerPass());
285 passes->AddItem(shadows);
286 seq->SetPasses(passes);
287
288 vtkNew<vtkCameraPass> cameraP;
289 cameraP->SetDelegatePass(seq);
290
291 // tell the renderer to use our render pass pipeline
292 vtkOpenGLRenderer* glrenderer = dynamic_cast<vtkOpenGLRenderer*>(renderer.GetPointer());
293 glrenderer->SetPass(cameraP);
294}
vtkRenderWindow * _renderWindow
Definition: G4VtkViewer.hh:187

◆ DrawView()

void G4VtkViewer::DrawView ( )
virtual

Implements G4VViewer.

Definition at line 243 of file G4VtkViewer.cc.

243 {
244 // First, a view should decide when to re-visit the G4 kernel.
245 // Sometimes it might not be necessary, e.g., if the scene is stored
246 // in a graphical database (e.g., OpenGL's display lists) and only
247 // the viewing angle has changed. But graphics systems without a
248 // graphical database will always need to visit the G4 kernel.
249
250 NeedKernelVisit(); // Default is - always visit G4 kernel.
251 // Note: this routine sets the fNeedKernelVisit flag of *all* the
252 // views of the scene.
253
254 ProcessView(); // The basic logic is here.
255
256 // Add HUD
257 DrawViewHUD();
258
259 // ...before finally...
260 FinishView(); // Flush streams and/or swap buffers.
261}
void ProcessView()
Definition: G4VViewer.cc:107
void NeedKernelVisit()
Definition: G4VViewer.cc:80
void FinishView()
Definition: G4VtkViewer.cc:321
void DrawViewHUD()
Definition: G4VtkViewer.cc:263

◆ DrawViewHUD()

void G4VtkViewer::DrawViewHUD ( )

Definition at line 263 of file G4VtkViewer.cc.

264{
265 // make sure text is always visible
267 infoTextActor->GetTextProperty()->SetColor(std::fmod(colour.GetRed() + 0.5, 1.0),
268 std::fmod(colour.GetGreen() + 0.5, 1.0),
269 std::fmod(colour.GetBlue() + 0.5, 1.0));
270 infoTextActor->GetTextProperty()->SetFontSize(20);
271 infoCallback->SetTextActor(infoTextActor);
272 renderer->AddObserver(vtkCommand::EndEvent, infoCallback);
273 renderer->AddActor(infoTextActor);
274}
G4double GetBlue() const
Definition: G4Colour.hh:154
G4double GetRed() const
Definition: G4Colour.hh:152
G4double GetGreen() const
Definition: G4Colour.hh:153
const G4Colour & GetBackgroundColour() const
vtkNew< vtkTextActor > infoTextActor
Definition: G4VtkViewer.hh:181
vtkNew< vtkInfoCallback > infoCallback
Definition: G4VtkViewer.hh:182

Referenced by DrawView().

◆ ExportOBJScene()

void G4VtkViewer::ExportOBJScene ( G4String  path)

Definition at line 377 of file G4VtkViewer.cc.

378{
379 vtkSmartPointer<vtkRenderWindow> _rw1 = vtkSmartPointer<vtkRenderWindow>::New();
380 _rw1->AddRenderer(_renderWindow->GetRenderers()->GetFirstRenderer());
381 vtkSmartPointer<vtkOBJExporter> exporter = vtkSmartPointer<vtkOBJExporter>::New();
382 exporter->SetRenderWindow(_rw1);
383 exporter->SetFilePrefix(path.c_str());
384 exporter->Write();
385}

◆ ExportScreenShot()

void G4VtkViewer::ExportScreenShot ( G4String  path,
G4String  format 
)

Definition at line 331 of file G4VtkViewer.cc.

332{
333
334 vtkImageWriter *imWriter = nullptr;
335
336 if(format == "bmp") {
337 imWriter = vtkBMPWriter::New();
338 }
339 else if (format == "jpg") {
340 imWriter = vtkJPEGWriter::New();
341 }
342 else if (format == "pnm") {
343 imWriter = vtkPNMWriter::New();
344 }
345 else if (format == "png") {
346 imWriter = vtkPNGWriter::New();
347 }
348 else if (format == "tiff") {
349 imWriter = vtkTIFFWriter::New();
350 }
351 else if (format == "ps") {
352 imWriter = vtkPostScriptWriter::New();
353 }
354 else {
355 imWriter = vtkPNGWriter::New();
356 }
357
358 _renderWindow->Render();
359
360 vtkSmartPointer<vtkWindowToImageFilter> winToImage = vtkSmartPointer<vtkWindowToImageFilter>::New();
361 winToImage->SetInput(_renderWindow);
362 winToImage->SetScale(1);
363 if(format == "ps")
364 {
365 winToImage->SetInputBufferTypeToRGB();
366 winToImage->ReadFrontBufferOff();
367 winToImage->Update();
368 }
369 else
370 {winToImage->SetInputBufferTypeToRGBA();}
371
372 imWriter->SetFileName((path+"."+format).c_str());
373 imWriter->SetInputConnection(winToImage->GetOutputPort());
374 imWriter->Write();
375}

◆ ExportView()

void G4VtkViewer::ExportView ( )
inline

Definition at line 178 of file G4VtkViewer.hh.

178{};

◆ ExportVRMLScene()

void G4VtkViewer::ExportVRMLScene ( G4String  path)

Definition at line 387 of file G4VtkViewer.cc.

388{
389 vtkSmartPointer<vtkRenderWindow> _rw1 = vtkSmartPointer<vtkRenderWindow>::New();
390 _rw1->AddRenderer(_renderWindow->GetRenderers()->GetFirstRenderer());
391 vtkSmartPointer<vtkVRMLExporter> exporter = vtkSmartPointer<vtkVRMLExporter>::New();
392 exporter->SetRenderWindow(_rw1);
393 exporter->SetFileName((path+".vrml").c_str());
394 exporter->Write();
395}

◆ ExportVTPScene()

void G4VtkViewer::ExportVTPScene ( G4String  path)

Definition at line 397 of file G4VtkViewer.cc.

398{
399 vtkSmartPointer<vtkRenderWindow> _rw1 = vtkSmartPointer<vtkRenderWindow>::New();
400 _rw1->AddRenderer(_renderWindow->GetRenderers()->GetFirstRenderer());
401 vtkSmartPointer<vtkSingleVTPExporter> exporter = vtkSmartPointer<vtkSingleVTPExporter>::New();
402 exporter->SetRenderWindow(_rw1);
403 exporter->SetFileName((path+".vtp").c_str());
404 exporter->Write();
405}

◆ FinishView()

void G4VtkViewer::FinishView ( void  )
virtual

Reimplemented from G4VViewer.

Definition at line 321 of file G4VtkViewer.cc.

322{
323 G4VtkSceneHandler& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
324 fVtkSceneHandler.Modified();
325
326 _renderWindow->Render();
327 _renderWindow->GetInteractor()->Initialize();
328 _renderWindow->GetInteractor()->Start();
329}
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:216

Referenced by DrawView(), and G4VtkQtViewer::FinishView().

◆ Initialise()

void G4VtkViewer::Initialise ( )
virtual

Reimplemented from G4VViewer.

Definition at line 65 of file G4VtkViewer.cc.

66{
67 _renderWindow = vtkRenderWindow::New();
68 renderWindowInteractor = vtkRenderWindowInteractor::New();
69
70#ifdef G4VTKDEBUG
71 G4cout << "G4VtkViewer::G4VtkViewer" << G4endl;
72 G4cout << "G4VtkViewer::G4VtkViewer> " << fVP.GetWindowSizeHintX() << " "
74 G4cout << "G4VtkViewer::G4VtkViewer> " << fVP.GetWindowLocationHintX() << " "
76#endif
77
78 // Need windowSizeX/Y - obtain from _renderWindow?
79 G4int screenSizeX = _renderWindow->GetScreenSize()[0];
80 G4int screenSizeY = _renderWindow->GetScreenSize()[1];
81 G4int positionX = fVP.GetWindowLocationHintX();
83 positionX = screenSizeX + positionX - fVP.GetWindowSizeHintX();
84 }
85 G4int positionY = fVP.GetWindowLocationHintY();
87 positionY = screenSizeY + positionY - fVP.GetWindowSizeHintY();
88 }
89 _renderWindow->SetPosition(positionX, positionY);
90#ifdef __APPLE__
91 // Adjust window size for Apple to make it correspond to OpenGL.
92 // Maybe it's OpenGL that shoud be adjusted.
93 const G4double pixelFactor = 2.;
94#else
95 const G4double pixelFactor = 1.;
96#endif
97 _renderWindow->SetSize
98 (pixelFactor*fVP.GetWindowSizeHintX(),pixelFactor*fVP.GetWindowSizeHintY());
99 _renderWindow->SetWindowName("Vtk viewer");
100
101 _renderWindow->AddRenderer(renderer);
102 renderWindowInteractor->SetRenderWindow(_renderWindow);
103
104 // TODO proper camera parameter settings
105 camera->SetPosition(0, 0, 1000);
106 camera->SetFocalPoint(0, 0, 0);
107 renderer->SetActiveCamera(camera);
108
109 //renderer->SetUseHiddenLineRemoval(1); // TODO needs to be an option
110 //renderer->SetUseShadows(1); // TODO needs to be an option
111
112 // Set callback to match VTK parameters to Geant4
113 geant4Callback->SetGeant4ViewParameters(&fVP);
114 renderer->AddObserver(vtkCommand::EndEvent, geant4Callback);
115
116 vtkSmartPointer<vtkInteractorStyleTrackballCamera> style =
117 vtkSmartPointer<vtkInteractorStyleTrackballCamera>::New();
118 renderWindowInteractor->SetInteractorStyle(style);
119
120 // DrawShadows();
121}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4int GetWindowLocationHintX() const
unsigned int GetWindowSizeHintX() const
G4bool IsWindowLocationHintXNegative() const
G4bool IsWindowLocationHintYNegative() const
G4int GetWindowLocationHintY() const
unsigned int GetWindowSizeHintY() const
vtkNew< vtkCamera > camera
Definition: G4VtkViewer.hh:185
vtkNew< vtkGeant4Callback > geant4Callback
Definition: G4VtkViewer.hh:183
vtkRenderWindowInteractor * renderWindowInteractor
Definition: G4VtkViewer.hh:188

◆ SetGeant4View()

void G4VtkViewer::SetGeant4View ( )
inline

Definition at line 179 of file G4VtkViewer.hh.

179{};

◆ SetView()

void G4VtkViewer::SetView ( )
virtual

Implements G4VViewer.

Definition at line 125 of file G4VtkViewer.cc.

125 {
126
127 // background colour
128 const G4Colour backgroundColour = fVP.GetBackgroundColour();
129 renderer->SetBackground(backgroundColour.GetRed(), backgroundColour.GetGreen(), backgroundColour.GetBlue());
130
131 // target and camera positions
133 if(radius <= 0.)
134 {radius = 1.;}
135 G4double cameraDistance = fVP.GetCameraDistance(radius);
136 G4Point3D viewpointDirection = fVP.GetViewpointDirection();
137 G4Point3D targetPoint = fVP.GetCurrentTargetPoint();
138 G4Point3D cameraPosition =
139 targetPoint + viewpointDirection.unit() * cameraDistance;
140 renderer->GetActiveCamera()->SetFocalPoint(targetPoint.x(),
141 targetPoint.y(),
142 targetPoint.z());
143 renderer->GetActiveCamera()->SetPosition(cameraPosition.x(),
144 cameraPosition.y(),
145 cameraPosition.z());
146 renderer->GetActiveCamera()->SetParallelScale(cameraDistance);
147
148 // need to set camera distance and parallel scale on first set view
149 if(firstSetView)
150 {
151 geant4Callback->SetVtkInitialValues(cameraDistance, cameraDistance);
152 firstSetView = false;
153 }
154
155 // projection type and view angle and zoom factor
156 G4double fieldHalfAngle = fVP.GetFieldHalfAngle();
157 G4double zoomFactor = fVP.GetZoomFactor();
158 vtkCamera* activeCamera = renderer->GetActiveCamera();
159 if(fieldHalfAngle == 0) {
160 activeCamera->SetParallelProjection(1);
161 activeCamera->SetParallelScale(activeCamera->GetParallelScale()/zoomFactor);
162 }
163 else {
164 activeCamera->SetParallelProjection(0);
165 activeCamera->SetViewAngle(2*fieldHalfAngle/M_PI*180);
166 activeCamera->SetPosition(cameraPosition.x()/zoomFactor,
167 cameraPosition.y()/zoomFactor,
168 cameraPosition.z()/zoomFactor);
169 }
170
171 // camera roll
172 // renderer->GetActiveCamera()->SetRoll(0);
173
174 // camera up direction
175 const G4Vector3D upVector = fVP.GetUpVector();
176 renderer->GetActiveCamera()->SetViewUp(upVector.x(),
177 upVector.y(),
178 upVector.z());
179
180 // Light
181 const G4Vector3D lightDirection = fVP.GetLightpointDirection();
182 G4bool lightsMoveWithCamera = fVP.GetLightsMoveWithCamera();
183 G4Vector3D lightPosition =
184 targetPoint + lightDirection.unit() * cameraDistance;
185
186 vtkLightCollection* currentLights = renderer->GetLights();
187 if (currentLights->GetNumberOfItems() != 0)
188 {
189 auto currentLight = dynamic_cast<vtkLight*>(currentLights->GetItemAsObject(0));
190 if (currentLight)
191 {
192 currentLight->SetPosition(lightPosition.x(),
193 lightPosition.y(),
194 lightPosition.z());
195 if (lightsMoveWithCamera)
196 {currentLight->SetLightTypeToCameraLight();}
197 else
198 {currentLight->SetLightTypeToSceneLight();}
199 }
200 }
201
202 // Rotation style
203#if 0
206 vtkSmartPointer<vtkInteractorStyleTrackballCamera> style =
207 vtkSmartPointer<vtkInteractorStyleTrackballCamera>::New();
208 _renderWindow->GetInteractor()->SetInteractorStyle(style);
209 }
211 // camera->SetViewUp(upVector.x(), upVector.y(), upVector.z());
212 vtkSmartPointer<vtkInteractorStyleTerrain> style =
213 vtkSmartPointer<vtkInteractorStyleTerrain>::New();
214 _renderWindow->GetInteractor()->SetInteractorStyle(style);
215 }
216#endif
217}
bool G4bool
Definition: G4Types.hh:86
#define M_PI
Definition: SbMath.h:33
virtual const G4VisExtent & GetExtent() const
G4double GetCameraDistance(G4double radius) const
const G4Vector3D & GetLightpointDirection() const
const G4Vector3D & GetViewpointDirection() const
const G4Point3D & GetCurrentTargetPoint() const
G4double GetFieldHalfAngle() const
G4double GetZoomFactor() const
const G4Vector3D & GetUpVector() const
RotationStyle GetRotationStyle() const
G4bool GetLightsMoveWithCamera() const
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
BasicVector3D< T > unit() const

◆ ShowView()

void G4VtkViewer::ShowView ( void  )
virtual

Reimplemented from G4VViewer.

Definition at line 296 of file G4VtkViewer.cc.

297{
298#ifdef G4VTKDEBUG
299 G4cout << "G4VtkViewer::ShowView() called." << G4endl;
300 // static_cast<G4VtkSceneHandler&>(fSceneHandler).PrintStores();
301#endif
302
303 G4VtkSceneHandler& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
304 fVtkSceneHandler.Modified();
305
306 infoTextActor->GetTextProperty()->SetFontSize(28);
308
309 // make sure text is always visible
310 infoTextActor->GetTextProperty()->SetColor(std::fmod(colour.GetRed() + 0.5, 1.0),
311 std::fmod(colour.GetGreen() + 0.5, 1.0),
312 std::fmod(colour.GetBlue() + 0.5, 1.0));
313 infoTextActor->GetTextProperty()->SetFontSize(20);
314 infoCallback->SetTextActor(infoTextActor);
315 renderer->AddObserver(vtkCommand::EndEvent, infoCallback);
316 geant4Callback->SetGeant4ViewParameters(&fVP);
317 renderer->AddObserver(vtkCommand::EndEvent, geant4Callback);
318 renderer->AddActor(infoTextActor);
319}

Member Data Documentation

◆ _renderWindow

vtkRenderWindow* G4VtkViewer::_renderWindow

◆ camera

vtkNew<vtkCamera> G4VtkViewer::camera

Definition at line 185 of file G4VtkViewer.hh.

Referenced by Initialise().

◆ geant4Callback

vtkNew<vtkGeant4Callback> G4VtkViewer::geant4Callback

Definition at line 183 of file G4VtkViewer.hh.

Referenced by G4VtkQtViewer::Initialise(), Initialise(), SetView(), and ShowView().

◆ infoCallback

vtkNew<vtkInfoCallback> G4VtkViewer::infoCallback

Definition at line 182 of file G4VtkViewer.hh.

Referenced by DrawViewHUD(), and ShowView().

◆ infoTextActor

vtkNew<vtkTextActor> G4VtkViewer::infoTextActor

Definition at line 181 of file G4VtkViewer.hh.

Referenced by DrawViewHUD(), and ShowView().

◆ light

vtkSmartPointer<vtkLight> G4VtkViewer::light

Definition at line 184 of file G4VtkViewer.hh.

◆ renderer

◆ renderWindowInteractor

vtkRenderWindowInteractor* G4VtkViewer::renderWindowInteractor

Definition at line 188 of file G4VtkViewer.hh.

Referenced by G4VtkQtViewer::Initialise(), and Initialise().


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