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

#include <G4VtkMessenger.hh>

+ Inheritance diagram for G4VtkMessenger:

Public Member Functions

 ~G4VtkMessenger () override
 
G4String GetCurrentValue (G4UIcommand *command) override
 
void SetNewValue (G4UIcommand *command, G4String newValue) override
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
G4bool CommandsShouldBeInMaster () const
 

Static Public Member Functions

static G4VtkMessengerGetInstance ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (const G4String &s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T>
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 39 of file G4VtkMessenger.hh.

Constructor & Destructor Documentation

◆ ~G4VtkMessenger()

G4VtkMessenger::~G4VtkMessenger ( )
override

Definition at line 274 of file G4VtkMessenger.cc.

275{
276 delete fpDirectory;
277 delete fpDirectorySet;
278 delete fpDirectoryAdd;
279 delete fpCommandExport;
280 delete fpCommandExportCutter;
281 delete fpCommandWarnings;
282 delete fpCommandDebugPrint;
283 delete fpCommandPolyhedronPipeline;
284 delete fpCommandImageOverlay;
285 delete fpCommandGeometryOverlay;
286}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VtkMessenger::GetCurrentValue ( G4UIcommand * command)
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 288 of file G4VtkMessenger.cc.

289{
290 return G4String();
291}

◆ GetInstance()

G4VtkMessenger * G4VtkMessenger::GetInstance ( )
static

Definition at line 42 of file G4VtkMessenger.cc.

43{
44 if (fpInstance == nullptr) fpInstance = new G4VtkMessenger;
45 return fpInstance;
46}

Referenced by G4Vtk::G4Vtk(), and G4VtkOffscreen::G4VtkOffscreen().

◆ SetNewValue()

void G4VtkMessenger::SetNewValue ( G4UIcommand * command,
G4String newValue )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 293 of file G4VtkMessenger.cc.

294{
295 G4VisManager* pVisManager = G4VisManager::GetInstance();
296
297 G4VViewer* pViewer = pVisManager->GetCurrentViewer();
298 if (pViewer == nullptr) {
299 G4cout << "G4VtkMessenger::SetNewValue: No current viewer.\n"
300 << "\"/vis/open\", or similar, to get one." << G4endl;
301 return;
302 }
303
304 auto* pVtkViewer = dynamic_cast<G4VtkViewer*>(pViewer);
305 if (pVtkViewer == nullptr) {
306 G4cout << "G4VtkMessenger::SetNewValue: Current viewer is not of type VTK. \n"
307 << "(It is \"" << pViewer->GetName() << "\".)\n"
308 << "Use \"/vis/viewer/select\" or \"/vis/open\"." << G4endl;
309 return;
310 }
311
312 if (command == fpCommandClearNonG4) {
313 auto sceneHandler = dynamic_cast<G4VtkSceneHandler*>(pViewer->GetSceneHandler());
314 auto transientStore = sceneHandler->GetTransientStore();
315
316 transientStore.ClearNonG4();
317 }
318 else if (command == fpCommandExport) {
319 G4String format, name;
320
321 std::istringstream iss(newValue);
322 iss >> format >> name;
323
324 if (format == "jpg" || format == "tiff" || format == "png" || format == "bmp" || format == "pnm"
325 || format == "ps")
326 pVtkViewer->ExportScreenShot(name, format);
327 else if (format == "obj")
328 pVtkViewer->ExportOBJScene(name);
329 else if (format == "vrml")
330 pVtkViewer->ExportVRMLScene(name);
331 else if (format == "vtp")
332 pVtkViewer->ExportVTPScene(name);
333 else if (format == "gltf")
334 pVtkViewer->ExportGLTFScene(name);
335 else if (format == "x3d")
336 pVtkViewer->ExportX3DScene(name);
337 else
338 G4cout << "Unknown /vis/vtk/export file format" << G4endl;
339 }
340 else if (command == fpCommandExportCutter) {
341 std::istringstream iss(newValue);
342
343 G4String fileName;
344 iss >> fileName;
345 pVtkViewer->ExportVTPCutter(fileName);
346 }
347 else if (command == fpCommandWarnings) {
348 if (G4UIcommand::ConvertToBool(newValue)) {
349 vtkObject::GlobalWarningDisplayOn();
350 }
351 else {
352 vtkObject::GlobalWarningDisplayOff();
353 }
354 }
355 else if (command == fpCommandHUD) {
356 if (G4UIcommand::ConvertToBool(newValue)) {
357 pVtkViewer->EnableHUD();
358 }
359 else {
360 pVtkViewer->DisableHUD();
361 }
362 }
363 else if (command == fpCameraOrientation) {
364 G4cout << newValue << G4endl;
365 if (G4UIcommand::ConvertToBool(newValue)) {
366 pVtkViewer->EnableCameraOrientationWidget();
367 }
368 else {
369 pVtkViewer->DisableCameraOrientationWidget();
370 }
371 }
372 else if (command == fpCommandDebugPrint) {
373 pVtkViewer->Print();
374 }
375 else if (command == fpCommandPolyhedronPipeline) {
376 G4String temp;
377
378 std::istringstream iss(newValue);
379
380 G4String pipelineType;
381 iss >> pipelineType;
382
383 pVtkViewer->SetPolyhedronPipeline(pipelineType);
384 }
385 else if (command == fpCommandImageOverlay) {
386 G4String temp;
387
388 G4String fileName;
389 G4double imageBottomLeft[2] = {0, 0};
390 G4double imageTopRight[2] = {1, 1};
391 G4double worldBottomLeft[2] = {0, 0};
392 G4double worldTopRight[2] = {1, 1};
393 G4double rotation[3] = {0, 0, 0};
394 G4double translation[3] = {0, 0, 0};
395 G4double alpha = 0.5;
396
397 std::istringstream iss(newValue);
398
399 iss >> fileName >> imageBottomLeft[0] >> imageBottomLeft[1] >> imageTopRight[0]
400 >> imageTopRight[1] >> worldBottomLeft[0] >> worldBottomLeft[1] >> worldTopRight[0]
401 >> worldTopRight[1] >> rotation[0] >> rotation[1] >> rotation[2] >> translation[0]
402 >> translation[1] >> translation[2] >> alpha;
403
404 pVtkViewer->AddImageOverlay(fileName, alpha, imageBottomLeft, worldBottomLeft, imageTopRight,
405 worldTopRight, rotation, translation);
406 }
407 else if (command == fpCommandGeometryOverlay) {
408 G4String temp;
409
410 G4String fileName;
411 G4double scale[3] = {1, 1, 1};
412 G4double rotation[3] = {0,0,0};
413 G4double translation[3] = {0,0,0};
414 G4double colour[3] = {0.5, 0.5, 0.5};
415 G4double alpha = 1.0;
416 G4String representation = "w";
417
418 std::istringstream iss(newValue);
419
420 iss >> fileName
421 >> scale[0] >> scale[1] >> scale[2]
422 >> rotation[0] >> rotation[1] >> rotation[2]
423 >> translation[0] >> translation[1] >> translation[2]
424 >> alpha >> representation;
425
426 G4cout << fileName << G4endl;
427 pVtkViewer->AddGeometryOverlay(fileName, colour, alpha, representation,
428 scale, rotation, translation);
429 }
430
431
432 else if (command == fpCommandClipper) {
433 pVtkViewer->EnableClipper(G4Plane3D(), true);
434 }
435 else if (command == fpCommandCutter) {
436 pVtkViewer->EnableCutter(G4Plane3D(), true);
437 }
438 else if (command == fpCommandShadow) {
439 pVtkViewer->EnableShadows();
440 }
441 else if (command == fpCommandInteractorStart) {
442 pVtkViewer->StartInteractor();
443 }
444}
HepGeom::Plane3D< G4double > G4Plane3D
Definition G4Plane3D.hh:34
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4bool ConvertToBool(const char *st)
const G4String & GetName() const
G4VSceneHandler * GetSceneHandler() const
G4VViewer * GetCurrentViewer() const
static G4VisManager * GetInstance()
const char * name(G4int ptype)

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