Geant4 11.3.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]\" or one of the driver/browser"
100 "\ncombinations that have the required functionality, e.g., HepRep.");
101 G4UIparameter* parameter;
102 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
103 parameter -> SetDefaultValue ("world");
104 fpCommandFindPath -> SetParameter (parameter);
105 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
106 parameter -> SetGuidance ("If negative, matches any copy no.");
107 parameter -> SetDefaultValue (-1);
108 fpCommandFindPath -> SetParameter (parameter);
109
110 fpCommandLocalAxes = new G4UIcmdWithoutParameter("/vis/touchable/localAxes",this);
111 fpCommandLocalAxes->SetGuidance("Draw local axes.");
112 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
113 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandLocalAxes,1);
114
115 fpCommandShowExtent = new G4UIcmdWithABool("/vis/touchable/showExtent",this);
116 fpCommandShowExtent->SetGuidance("Print extent of touchable.");
117 fpCommandShowExtent->SetGuidance("If parameter == true, also draw.");
118 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
119 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandShowExtent,1);
120 fpCommandShowExtent->SetParameterName("draw", omitable = true);
121 fpCommandShowExtent->SetDefaultValue(false);
122
123 fpCommandTwinkle = new G4UIcmdWithoutParameter("/vis/touchable/twinkle",this);
124 fpCommandTwinkle->SetGuidance("Cause touchable to twinkle.");
125 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
126 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandTwinkle,1);
127
128 fpCommandVolumeForField = new G4UIcmdWithABool("/vis/touchable/volumeForField",this);
129 fpCommandVolumeForField->SetGuidance("Set volume for field.");
130 fpCommandVolumeForField->SetGuidance("If parameter == true, also draw.");
131 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
132 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandVolumeForField,1);
133 fpCommandVolumeForField->SetParameterName("draw", omitable = true);
134 fpCommandVolumeForField->SetDefaultValue(false);
135}
bool G4bool
Definition G4Types.hh:86
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)

◆ ~G4VisCommandsTouchable()

G4VisCommandsTouchable::~G4VisCommandsTouchable ( )
virtual

Definition at line 137 of file G4VisCommandsTouchable.cc.

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

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandsTouchable::GetCurrentValue ( G4UIcommand * command)
virtual

Reimplemented from G4UImessenger.

Definition at line 150 of file G4VisCommandsTouchable.cc.

150 {
151 return "";
152}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 154 of file G4VisCommandsTouchable.cc.

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