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

#include <G4VisCommandsTouchable.hh>

Inheritance diagram for G4VisCommandsTouchable:

Public Member Functions

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

Additional Inherited Members

Static Public Member Functions inherited from G4VVisCommand
static G4VisManagerGetVisManager ()
static void SetVisManager (G4VisManager *pVisManager)
static const G4ColourGetCurrentTextColour ()
Protected Member Functions inherited from G4VVisCommand
void SetViewParameters (G4VViewer *viewer, const G4ViewParameters &viewParams)
void RefreshIfRequired (G4VViewer *viewer)
void InterpolateViews (G4VViewer *currentViewer, const std::vector< G4ViewParameters > &viewVector, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
void InterpolateToNewView (G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
void Twinkle (G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
const G4StringConvertToColourGuidance ()
void ConvertToColour (G4Colour &colour, const G4String &redOrString, G4double green, G4double blue, G4double opacity)
G4bool ProvideValueOfUnit (const G4String &where, const G4String &unit, const G4String &category, G4double &value)
void CopyCameraParameters (G4ViewParameters &target, const G4ViewParameters &from)
void CheckSceneAndNotifyHandlers (G4Scene *=nullptr)
G4bool CheckView ()
void G4VisCommandsSceneAddUnsuccessful (G4VisManager::Verbosity verbosity)
void CopyGuidanceFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
void CopyParametersFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd)
void DrawExtent (const G4VisExtent &)
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)
Static Protected Member Functions inherited from G4VVisCommand
static G4String ConvertToString (G4double x, G4double y, const char *unitName)
static G4bool ConvertToDoublePair (const G4String &paramString, G4double &xval, G4double &yval)
Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
G4String baseDirName = ""
G4bool commandsShouldBeInMaster = false
Static Protected Attributes inherited from G4VVisCommand
static G4VisManagerfpVisManager = nullptr
static G4int fCurrentArrow3DLineSegmentsPerCircle = 6
static G4Colour fCurrentColour = G4Colour::White()
static G4double fCurrentLineWidth = 1.
static G4Colour fCurrentTextColour = G4Colour::Blue()
static G4Text::Layout fCurrentTextLayout = G4Text::left
static G4double fCurrentTextSize = 12.
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
static G4VisExtent fCurrentExtentForField
static std::vector< G4PhysicalVolumesSearchScene::FindingsfCurrrentPVFindingsForField
static G4bool fThereWasAViewer = false
static G4ViewParameters fExistingVP
static G4SceneTreeItem fExistingSceneTree

Detailed Description

Definition at line 39 of file G4VisCommandsTouchable.hh.

Constructor & Destructor Documentation

◆ G4VisCommandsTouchable()

G4VisCommandsTouchable::G4VisCommandsTouchable ( )

Definition at line 46 of file G4VisCommandsTouchable.cc.

47{
48 G4bool omitable;
49
50 fpCommandCentreAndZoomInOn = new G4UIcmdWithoutParameter("/vis/touchable/centreAndZoomInOn",this);
51 fpCommandCentreAndZoomInOn->SetGuidance ("Centre and zoom in on the current touchable.");
52 fpCommandCentreAndZoomInOn->SetGuidance
53 ("Use \"/vis/set/touchable\" to set current touchable.");
54 fpCommandCentreAndZoomInOn->SetGuidance
55 ("You may also need \"/vis/touchable/findPath\".");
56 fpCommandCentreAndZoomInOn->SetGuidance
57 ("Use \"/vis/touchable/set\" to set attributes.");
58
59 fpCommandCentreOn = new G4UIcmdWithoutParameter("/vis/touchable/centreOn",this);
60 fpCommandCentreOn->SetGuidance ("Centre the view on the current touchable.");
61 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
62 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandCentreOn,1);
63
64 fpCommandDraw = new G4UIcmdWithABool("/vis/touchable/draw",this);
65 fpCommandDraw->SetGuidance("Draw touchable.");
66 fpCommandDraw->SetGuidance
67 ("If parameter == true, also draw extent as a white wireframe box.");
68 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
69 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandDraw,1);
70 fpCommandDraw->SetParameterName("extent", omitable = true);
71 fpCommandDraw->SetDefaultValue(false);
72
73 fpCommandDump = new G4UIcmdWithoutParameter("/vis/touchable/dump",this);
74 fpCommandDump->SetGuidance("Dump touchable attributes.");
75 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
76 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandDump,1);
77
78 fpCommandExtentForField = new G4UIcmdWithABool("/vis/touchable/extentForField",this);
79 fpCommandExtentForField->SetGuidance("Set extent for field.");
80 fpCommandExtentForField->SetGuidance("If parameter == true, also draw.");
81 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
82 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandExtentForField,1);
83 fpCommandExtentForField->SetParameterName("draw", omitable = true);
84 fpCommandExtentForField->SetDefaultValue(false);
85
86 fpCommandFindPath = new G4UIcommand("/vis/touchable/findPath",this);
87 fpCommandFindPath->SetGuidance
88 ("Prints the path to touchable and its logical volume mother"
89 "\ngiven a physical volume name and copy no.");
90 fpCommandFindPath -> SetGuidance
91 ("A search of all worlds is made and all physical volume names are"
92 "\nmatched against the argument of this command. If this is of the"
93 "\nform \"/regexp/\", where regexp is a regular expression (see C++ regex),"
94 "\nthe physical volume name is matched against regexp by the usual rules"
95 "\nof regular expression matching. Otherwise an exact match is required."
96 "\nFor example, \"/Shap/\" matches \"Shape1\" and \"Shape2\".");
97 fpCommandFindPath -> SetGuidance
98 ("It may help to see a textual representation of the geometry hierarchy of"
99 "\nthe worlds. Try \"/vis/drawTree [worlds]\".");
100 G4UIparameter* parameter;
101 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
102 parameter -> SetDefaultValue ("world");
103 fpCommandFindPath -> SetParameter (parameter);
104 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
105 parameter -> SetGuidance ("If negative, matches any copy no.");
106 parameter -> SetDefaultValue (-1);
107 fpCommandFindPath -> SetParameter (parameter);
108
109 fpCommandLocalAxes = new G4UIcmdWithoutParameter("/vis/touchable/localAxes",this);
110 fpCommandLocalAxes->SetGuidance("Draw local axes.");
111 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
112 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandLocalAxes,1);
113
114 fpCommandShowExtent = new G4UIcmdWithABool("/vis/touchable/showExtent",this);
115 fpCommandShowExtent->SetGuidance("Print extent of touchable.");
116 fpCommandShowExtent->SetGuidance("If parameter == true, also draw.");
117 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
118 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandShowExtent,1);
119 fpCommandShowExtent->SetParameterName("draw", omitable = true);
120 fpCommandShowExtent->SetDefaultValue(false);
121
122 fpCommandTwinkle = new G4UIcmdWithoutParameter("/vis/touchable/twinkle",this);
123 fpCommandTwinkle->SetGuidance("Cause touchable to twinkle.");
124 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
125 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandTwinkle,1);
126
127 fpCommandVolumeForField = new G4UIcmdWithABool("/vis/touchable/volumeForField",this);
128 fpCommandVolumeForField->SetGuidance("Set volume for field.");
129 fpCommandVolumeForField->SetGuidance("If parameter == true, also draw.");
130 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
131 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandVolumeForField,1);
132 fpCommandVolumeForField->SetParameterName("draw", omitable = true);
133 fpCommandVolumeForField->SetDefaultValue(false);
134}
bool G4bool
Definition G4Types.hh:86
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)

◆ ~G4VisCommandsTouchable()

G4VisCommandsTouchable::~G4VisCommandsTouchable ( )
virtual

Definition at line 136 of file G4VisCommandsTouchable.cc.

136 {
137 delete fpCommandVolumeForField;
138 delete fpCommandTwinkle;
139 delete fpCommandShowExtent;
140 delete fpCommandLocalAxes;
141 delete fpCommandFindPath;
142 delete fpCommandExtentForField;
143 delete fpCommandDump;
144 delete fpCommandDraw;
145 delete fpCommandCentreAndZoomInOn;
146 delete fpCommandCentreOn;
147}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandsTouchable::GetCurrentValue ( G4UIcommand * command)
virtual

Reimplemented from G4UImessenger.

Definition at line 149 of file G4VisCommandsTouchable.cc.

149 {
150 return "";
151}

◆ SetNewValue()

void G4VisCommandsTouchable::SetNewValue ( G4UIcommand * command,
G4String newValue )
virtual

Reimplemented from G4UImessenger.

Definition at line 153 of file G4VisCommandsTouchable.cc.

155{
156 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
157 G4bool warn = verbosity >= G4VisManager::warnings;
158
159 G4UImanager* UImanager = G4UImanager::GetUIpointer();
160
161 G4TransportationManager* transportationManager =
163
164 size_t nWorlds = transportationManager->GetNoWorlds();
165
166 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
167 if (!world) {
168 if (verbosity >= G4VisManager::errors) {
169 G4warn <<
170 "ERROR: G4VisCommandsTouchable::SetNewValue:"
171 "\n No world. Maybe the geometry has not yet been defined."
172 "\n Try \"/run/initialize\""
173 << G4endl;
174 }
175 return;
176 }
177
178 if (command == fpCommandDump) {
179
180 G4PhysicalVolumeModel::TouchableProperties properties =
182 if (properties.fpTouchablePV) {
183 // To handle parameterisations we have to set the copy number
184 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
185 G4PhysicalVolumeModel tempPVModel
186 (properties.fpTouchablePV,
188 properties.fTouchableGlobalTransform,
189 nullptr, // Modelling parameters (not used)
190 true, // use full extent (prevents calculating own extent, which crashes)
191 properties.fTouchableBaseFullPVPath);
192 const std::map<G4String,G4AttDef>* attDefs = tempPVModel.GetAttDefs();
193 std::vector<G4AttValue>* attValues = tempPVModel.CreateCurrentAttValues();
194 G4cout << G4AttCheck(attValues,attDefs);
195 delete attValues;
196 const auto lv = properties.fpTouchablePV->GetLogicalVolume();
197 const auto polyhedron = lv->GetSolid()->GetPolyhedron();
198 if (polyhedron) {
199 polyhedron->SetVisAttributes(lv->GetVisAttributes());
200 G4cout << "\nLocal polyhedron coordinates:\n" << *polyhedron;
201 const G4Transform3D& transform = tempPVModel.GetCurrentTransform();
202 polyhedron->Transform(transform);
203 G4cout << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
204 }
205 } else {
206 G4warn << "Touchable not found." << G4endl;
207 }
208 return;
209
210 } else if (command == fpCommandFindPath) {
211
212 G4String pvName;
213 G4int copyNo;
214 std::istringstream iss(newValue);
215 iss >> pvName >> copyNo;
216 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
217 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
218 transportationManager->GetWorldsIterator();
219 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
220 G4PhysicalVolumeModel searchModel (*iterWorld); // Unlimited depth.
221 G4ModelingParameters mp; // Default - no culling.
222 searchModel.SetModelingParameters (&mp);
223 // Find all instances at any position in the tree
224 G4PhysicalVolumesSearchScene searchScene (&searchModel, pvName, copyNo);
225 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
226 for (const auto& findings: searchScene.GetFindings()) {
227 findingsVector.push_back(findings);
228 }
229 }
230 for (const auto& findings: findingsVector) {
231 G4cout
232 << findings.fFoundBasePVPath
233 << ' ' << findings.fpFoundPV->GetName()
234 << ' ' << findings.fFoundPVCopyNo
235 << " (mother logical volume: "
236 << findings.fpFoundPV->GetMotherLogical()->GetName()
237 << ')'
238 << G4endl;
239 }
240 if (findingsVector.size()) {
241 G4cout
242 << "Use this to set a particular touchable with \"/vis/set/touchable <path>\""
243 << "\nor to see overlaps: \"/vis/drawLogicalVolume <mother-logical-volume-name>\""
244 << G4endl;
245 } else {
246 G4warn << pvName;
247 if (copyNo >= 0) G4warn << ':' << copyNo;
248 G4warn << " not found" << G4endl;
249 }
250 return;
251 }
252
253 G4VViewer* currentViewer = fpVisManager -> GetCurrentViewer ();
254 if (!currentViewer) {
255 if (verbosity >= G4VisManager::errors) {
256 G4warn <<
257 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
258 << G4endl;
259 }
260 return;
261 }
262
263 G4Scene* currentScene = fpVisManager->GetCurrentScene();
264 if (!currentScene) {
265 if (verbosity >= G4VisManager::errors) {
266 G4warn <<
267 "ERROR: No current scene - \"/vis/scene/list\" to see possibilities."
268 << G4endl;
269 }
270 return;
271 }
272
273 if (command == fpCommandCentreOn || command == fpCommandCentreAndZoomInOn) {
274
275 // For twinkling...
276 std::vector<std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>> touchables;
277
278 G4PhysicalVolumeModel::TouchableProperties properties =
280 if (properties.fpTouchablePV) {
281 // To handle parameterisations, set copy number
282 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
283 G4PhysicalVolumeModel tempPVModel
284 (properties.fpTouchablePV,
286 properties.fTouchableGlobalTransform,
287 nullptr, // Modelling parameters (not used)
288 true, // use full extent (prevents calculating own extent, which crashes)
289 properties.fTouchableBaseFullPVPath);
290 touchables.push_back(properties.fTouchableFullPVPath); // Only one in this case
291 // Use a temporary scene in order to find vis extent
292 G4Scene tempScene("Centre Scene");
293 G4bool successful = tempScene.AddRunDurationModel(&tempPVModel,warn);
294 if (!successful) return;
295 if (verbosity >= G4VisManager::parameters) {
296 G4cout
297 << "Touchable " << fCurrentTouchableProperties.fTouchablePath
298 << ",\n has been added to temporary scene \"" << tempScene.GetName() << "\"."
299 << G4endl;
300 }
301
302 const G4VisExtent& newExtent = tempScene.GetExtent();
303 const G4ThreeVector& newTargetPoint = newExtent.GetExtentCentre();
304 G4ViewParameters saveVP = currentViewer->GetViewParameters();
305 G4ViewParameters newVP = saveVP;
306 if (command == fpCommandCentreAndZoomInOn) {
307 // Calculate the new zoom factor
308 const G4double zoomFactor
309 = currentScene->GetExtent().GetExtentRadius()/newExtent.GetExtentRadius();
310 newVP.SetZoomFactor(zoomFactor);
311 }
312 // Change the target point
313 const G4Point3D& standardTargetPoint = currentScene->GetStandardTargetPoint();
314 newVP.SetCurrentTargetPoint(newTargetPoint - standardTargetPoint);
315
316 if (currentViewer->GetKernelVisitElapsedTimeSeconds() < 0.1) {
317 // Interpolate
318 auto keepVisVerbose = fpVisManager->GetVerbosity();
319 fpVisManager->SetVerboseLevel(G4VisManager::errors);
320 if (newVP != saveVP) InterpolateToNewView(currentViewer, saveVP, newVP);
321 // ...and twinkle
322 Twinkle(currentViewer,newVP,touchables);
323 fpVisManager->SetVerboseLevel(keepVisVerbose);
324 }
325
326 if (verbosity >= G4VisManager::confirmations) {
327 G4cout
328 << "Viewer \"" << currentViewer->GetName()
329 << "\" centred ";
330 if (fpCommandCentreAndZoomInOn) {
331 G4cout << "and zoomed in";
332 }
333 G4cout << " on touchable\n" << fCurrentTouchableProperties.fTouchablePath
334 << G4endl;
335 }
336 SetViewParameters(currentViewer, newVP);
337 } else {
338 G4warn << "Touchable not found." << G4endl;
339 }
340 return;
341
342 } else if (command == fpCommandDraw) {
343
344 G4PhysicalVolumeModel::TouchableProperties properties =
346 if (properties.fpTouchablePV) {
347 // To handle parameterisations we have to set the copy number
348 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
349 G4PhysicalVolumeModel* pvModel = new G4PhysicalVolumeModel
350 (properties.fpTouchablePV,
352 properties.fTouchableGlobalTransform,
353 nullptr, // Modelling parameters (not used)
354 true, // use full extent (prevents calculating own extent, which crashes)
355 properties.fTouchableBaseFullPVPath);
356
357 UImanager->ApplyCommand("/vis/scene/create");
358 currentScene = fpVisManager->GetCurrentScene(); // New current scene
359 G4bool successful = currentScene->AddRunDurationModel(pvModel,warn);
360 UImanager->ApplyCommand("/vis/sceneHandler/attach");
361
362 if (successful) {
363 if (fpCommandDraw->GetNewBoolValue(newValue)) {
364 const auto& extent = pvModel->GetExtent();
365 const G4double halfX = (extent.GetXmax()-extent.GetXmin())/2.;
366 const G4double halfY = (extent.GetYmax()-extent.GetYmin())/2.;
367 const G4double halfZ = (extent.GetZmax()-extent.GetZmin())/2.;
368 G4Box extentBox("extent",halfX,halfY,halfZ);
369 G4VisAttributes extentVA;
370 extentVA.SetForceWireframe();
371 fpVisManager->Draw(extentBox,extentVA,G4Translate3D(extent.GetExtentCentre()));
372 }
373 if (verbosity >= G4VisManager::confirmations) {
374 G4cout << "\"" << properties.fpTouchablePV->GetName()
375 << "\", copy no. " << properties.fCopyNo << " drawn";
376 if (fpCommandDraw->GetNewBoolValue(newValue)) {
377 G4cout << " with extent box";
378 }
379 G4cout << '.' << G4endl;
380 }
381 } else {
383 }
384 } else {
385 G4warn << "Touchable not found." << G4endl;
386 }
387 return;
388
389 } else if (command == fpCommandExtentForField) {
390
391 G4PhysicalVolumeModel::TouchableProperties properties =
393 if (properties.fpTouchablePV) {
394 G4VisExtent extent
395 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
396 extent.Transform(properties.fTouchableGlobalTransform);
397 fCurrentExtentForField = extent;
399 if (verbosity >= G4VisManager::confirmations) {
400 G4cout << "Extent for field set to " << extent
401 << "\nVolume for field has been cleared."
402 << G4endl;
403 }
404 if (fpCommandExtentForField->GetNewBoolValue(newValue)) {
405 DrawExtent(extent);
406 }
407 } else {
408 G4warn << "Touchable not found." << G4endl;
409 }
410 return;
411
412 } else if (command == fpCommandLocalAxes) {
413
414 G4PhysicalVolumeModel::TouchableProperties properties =
416 if (properties.fpTouchablePV) {
417 const auto& transform = fCurrentTouchableProperties.fTouchableGlobalTransform;
418 const auto& extent = fCurrentTouchableProperties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
419 const G4double lengthMax = extent.GetExtentRadius()/2.;
420 const G4double intLog10LengthMax = std::floor(std::log10(lengthMax));
421 G4double length = std::pow(10,intLog10LengthMax);
422 if (5.*length < lengthMax) length *= 5.;
423 else if (2.*length < lengthMax) length *= 2.;
424 G4AxesModel axesModel(0.,0.,0.,length,transform);
425 axesModel.SetGlobalTag("LocalAxesModel");
426 axesModel.DescribeYourselfTo(*fpVisManager->GetCurrentSceneHandler());
427 G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/refresh");
428 } else {
429 G4warn << "Touchable not found." << G4endl;
430 }
431 return;
432
433 } else if (command == fpCommandShowExtent) {
434
435 G4PhysicalVolumeModel::TouchableProperties properties =
437 if (properties.fpTouchablePV) {
438 G4VisExtent extent
439 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
440 extent.Transform(properties.fTouchableGlobalTransform);
441 G4cout << extent << G4endl;
442 if (fpCommandShowExtent->GetNewBoolValue(newValue)) DrawExtent(extent);
443 } else {
444 G4warn << "Touchable not found." << G4endl;
445 }
446 return;
447
448 } else if (command == fpCommandTwinkle) {
449
450 if (currentViewer->GetKernelVisitElapsedTimeSeconds() < 0.1) {
451 G4PhysicalVolumeModel::TouchableProperties properties =
453 if (properties.fpTouchablePV) {
454 std::vector<std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>> touchables;
455 touchables.push_back(properties.fTouchableFullPVPath);
456 auto keepVisVerbose = fpVisManager->GetVerbosity();
457 fpVisManager->SetVerboseLevel(G4VisManager::errors);
458 auto keepVP = currentViewer->GetViewParameters();
459 Twinkle(currentViewer,currentViewer->GetViewParameters(),touchables);
460 SetViewParameters(currentViewer, keepVP);
461 fpVisManager->SetVerboseLevel(keepVisVerbose);
462 } else {
463 G4warn << "Touchable not found." << G4endl;
464 }
465 } else {
466 G4warn << "Twinkling not available - image construction time too long." << G4endl;
467 }
468 return;
469
470 } else if (command == fpCommandVolumeForField) {
471
472 G4PhysicalVolumeModel::TouchableProperties properties =
474 if (properties.fpTouchablePV) {
475 G4VisExtent extent
476 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
477 extent.Transform(properties.fTouchableGlobalTransform);
478 fCurrentExtentForField = extent;
481 (G4PhysicalVolumesSearchScene::Findings(properties));
482 if (verbosity >= G4VisManager::confirmations) {
483 G4cout
484 << "Volume for field set to " << properties.fpTouchablePV->GetName()
485 << ':' << properties.fCopyNo
486 << " at " << properties.fTouchableBaseFullPVPath
487 << G4endl;
488 }
489 if (fpCommandVolumeForField->GetNewBoolValue(newValue)) {
490 DrawExtent(extent);
491 }
492 } else {
493 G4warn << "Touchable not found." << G4endl;
494 }
495 return;
496
497 } else {
498
499 if (verbosity >= G4VisManager::errors) {
500 G4warn <<
501 "ERROR: G4VisCommandsTouchable::SetNewValue: unrecognised command."
502 << G4endl;
503 }
504 return;
505 }
506}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
#define G4warn
Definition G4Scene.cc:41
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VSolid * GetSolid() const
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition G4Scene.cc:160
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()
const G4VisExtent & GetExtent() const
virtual void SetCopyNo(G4int CopyNo)=0
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4VisExtent GetExtent() const
Definition G4VSolid.cc:708
virtual G4Polyhedron * GetPolyhedron() const
Definition G4VSolid.cc:731
const G4String & GetName() const
const G4ViewParameters & GetViewParameters() const
G4double GetKernelVisitElapsedTimeSeconds() const
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
void InterpolateToNewView(G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
static std::vector< G4PhysicalVolumesSearchScene::Findings > fCurrrentPVFindingsForField
static G4VisManager * fpVisManager
static G4VisExtent fCurrentExtentForField
void DrawExtent(const G4VisExtent &)
void SetViewParameters(G4VViewer *viewer, const G4ViewParameters &viewParams)
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
void Twinkle(G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
void SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
void SetZoomFactor(G4double zoomFactor)
void SetForceWireframe(G4bool=true)
G4double GetExtentRadius() const
G4VisExtent & Transform(const G4Transform3D &)
const G4Point3D & GetExtentCentre() const
G4PhysicalVolumeModel::TouchableProperties FindTouchableProperties(G4ModelingParameters::PVNameCopyNoPath path)
std::vector< G4PhysicalVolumeNodeID > fTouchableFullPVPath
std::vector< G4PhysicalVolumeNodeID > fTouchableBaseFullPVPath

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