Geant4-11
G4VisCommandsTouchable.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//
26//
27
28// /vis/touchable commands - John Allison 14th May 2014
29
31
32#include "G4UImanager.hh"
33#include "G4UIcmdWithAString.hh"
35#include "G4UIcmdWithABool.hh"
37#include "G4TouchableUtils.hh"
39#include "G4AttDef.hh"
40#include "G4AttValue.hh"
41#include "G4AttCheck.hh"
42#include "G4AxesModel.hh"
43
45{
46 G4bool omitable;
47
48 fpCommandCentreAndZoomInOn = new G4UIcmdWithoutParameter("/vis/touchable/centreAndZoomInOn",this);
49 fpCommandCentreAndZoomInOn->SetGuidance ("Centre and zoom in on the current touchable.");
51 ("Use \"/vis/set/touchable\" to set current touchable.");
53 ("You may also need \"/vis/touchable/findPath\".");
55 ("Use \"/vis/touchable/set\" to set attributes.");
56
57 fpCommandCentreOn = new G4UIcmdWithoutParameter("/vis/touchable/centreOn",this);
58 fpCommandCentreOn->SetGuidance ("Centre the view on the current touchable.");
59 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
61
62 fpCommandDraw = new G4UIcmdWithoutParameter("/vis/touchable/draw",this);
63 fpCommandDraw->SetGuidance("Draw touchable.");
64 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
66
67 fpCommandDump = new G4UIcmdWithoutParameter("/vis/touchable/dump",this);
68 fpCommandDump->SetGuidance("Dump touchable attributes.");
69 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
71
72 fpCommandExtentForField = new G4UIcmdWithABool("/vis/touchable/extentForField",this);
73 fpCommandExtentForField->SetGuidance("Set extent for field.");
74 fpCommandExtentForField->SetGuidance("If parameter == true, also draw.");
75 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
77 fpCommandExtentForField->SetParameterName("draw", omitable = true);
79
80 fpCommandFindPath = new G4UIcommand("/vis/touchable/findPath",this);
82 ("Prints the path to touchable and its logical volume mother"
83 "\ngiven a physical volume name and copy no.");
84 fpCommandFindPath -> SetGuidance
85 ("A search of all worlds is made and all physical volume names are"
86 "\nmatched against the argument of this command. If this is of the"
87 "\nform \"/regexp/\", where regexp is a regular expression (see C++ regex),"
88 "\nthe physical volume name is matched against regexp by the usual rules"
89 "\nof regular expression matching. Otherwise an exact match is required."
90 "\nFor example, \"/Shap/\" matches \"Shape1\" and \"Shape2\".");
91 fpCommandFindPath -> SetGuidance
92 ("It may help to see a textual representation of the geometry hierarchy of"
93 "\nthe worlds. Try \"/vis/drawTree [worlds]\" or one of the driver/browser"
94 "\ncombinations that have the required functionality, e.g., HepRep.");
95 G4UIparameter* parameter;
96 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
97 parameter -> SetDefaultValue ("world");
98 fpCommandFindPath -> SetParameter (parameter);
99 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
100 parameter -> SetGuidance ("If negative, matches any copy no.");
101 parameter -> SetDefaultValue (-1);
102 fpCommandFindPath -> SetParameter (parameter);
103
104 fpCommandLocalAxes = new G4UIcmdWithoutParameter("/vis/touchable/localAxes",this);
105 fpCommandLocalAxes->SetGuidance("Draw local axes.");
106 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
108
109 fpCommandShowExtent = new G4UIcmdWithABool("/vis/touchable/showExtent",this);
110 fpCommandShowExtent->SetGuidance("Print extent of touchable.");
111 fpCommandShowExtent->SetGuidance("If parameter == true, also draw.");
112 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
114 fpCommandShowExtent->SetParameterName("draw", omitable = true);
116
117 fpCommandVolumeForField = new G4UIcmdWithABool("/vis/touchable/volumeForField",this);
118 fpCommandVolumeForField->SetGuidance("Set volume for field.");
119 fpCommandVolumeForField->SetGuidance("If parameter == true, also draw.");
120 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
122 fpCommandVolumeForField->SetParameterName("draw", omitable = true);
124}
125
128 delete fpCommandShowExtent;
129 delete fpCommandLocalAxes;
130 delete fpCommandFindPath;
132 delete fpCommandDump;
133 delete fpCommandDraw;
135 delete fpCommandCentreOn;
136}
137
139 return "";
140}
141
143(G4UIcommand* command, G4String newValue)
144{
146 G4bool warn = verbosity >= G4VisManager::warnings;
147
149
150 G4TransportationManager* transportationManager =
152
153 size_t nWorlds = transportationManager->GetNoWorlds();
154
155 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
156 if (!world) {
157 if (verbosity >= G4VisManager::errors) {
158 G4cerr <<
159 "ERROR: G4VisCommandsTouchable::SetNewValue:"
160 "\n No world. Maybe the geometry has not yet been defined."
161 "\n Try \"/run/initialize\""
162 << G4endl;
163 }
164 return;
165 }
166
167 G4VViewer* currentViewer = fpVisManager -> GetCurrentViewer ();
168 if (!currentViewer) {
169 if (verbosity >= G4VisManager::errors) {
170 G4cerr <<
171 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
172 << G4endl;
173 }
174 return;
175 }
176
177 G4Scene* currentScene = fpVisManager->GetCurrentScene();
178 if (!currentScene) {
179 if (verbosity >= G4VisManager::errors) {
180 G4cerr <<
181 "ERROR: No current scene - \"/vis/scene/list\" to see possibilities."
182 << G4endl;
183 }
184 return;
185 }
186
187 if (command == fpCommandCentreOn || command == fpCommandCentreAndZoomInOn) {
188
191 if (properties.fpTouchablePV) {
192 // To handle parameterisations, set copy number
193 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
194 G4PhysicalVolumeModel tempPVModel
195 (properties.fpTouchablePV,
197 properties.fTouchableGlobalTransform,
198 nullptr, // Modelling parameters (not used)
199 true, // use full extent (prevents calculating own extent, which crashes)
200 properties.fTouchableBaseFullPVPath);
201 // Use a temporary scene in order to find vis extent
202 G4Scene tempScene("Centre Scene");
203 G4bool successful = tempScene.AddRunDurationModel(&tempPVModel,warn);
204 if (successful) {
205 if (verbosity >= G4VisManager::confirmations) {
206 G4cout
208 << ",\n has been added to temporary scene \"" << tempScene.GetName() << "\"."
209 << G4endl;
210 }
211 }
212 const G4VisExtent& newExtent = tempScene.GetExtent();
213 const G4ThreeVector& newTargetPoint = newExtent.GetExtentCentre();
214 G4ViewParameters saveVP = currentViewer->GetViewParameters();
215 G4ViewParameters newVP = saveVP;
216 if (command == fpCommandCentreAndZoomInOn) {
217 // Calculate the new zoom factor
218 const G4double zoomFactor
219 = currentScene->GetExtent().GetExtentRadius()/newExtent.GetExtentRadius();
220 newVP.SetZoomFactor(zoomFactor);
221 }
222 // Change the target point
223 const G4Point3D& standardTargetPoint = currentScene->GetStandardTargetPoint();
224 newVP.SetCurrentTargetPoint(newTargetPoint - standardTargetPoint);
225 // Interpolate
226 InterpolateToNewView(currentViewer, saveVP, newVP);
227 if (verbosity >= G4VisManager::confirmations) {
228 G4cout
229 << "Viewer \"" << currentViewer->GetName()
230 << "\" centred ";
232 G4cout << "and zoomed in";
233 }
235 << G4endl;
236 }
237 SetViewParameters(currentViewer, newVP);
238 } else {
239 G4cout << "Touchable not found." << G4endl;
240 }
241 return;
242
243 } else if (command == fpCommandDraw) {
244
247 if (properties.fpTouchablePV) {
248 // To handle paramaterisations we have to set the copy number
249 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
251 (properties.fpTouchablePV,
253 properties.fTouchableGlobalTransform,
254 nullptr, // Modelling parameters (not used)
255 true, // use full extent (prevents calculating own extent, which crashes)
256 properties.fTouchableBaseFullPVPath);
257
258 G4int keepVerbose = UImanager->GetVerboseLevel();
259 G4int newVerbose(0);
260 if (keepVerbose >= 2 || verbosity >= G4VisManager::confirmations)
261 newVerbose = 2;
262 UImanager->SetVerboseLevel(newVerbose);
263 UImanager->ApplyCommand("/vis/scene/create");
264 currentScene = fpVisManager->GetCurrentScene(); // New current scene
265 G4bool successful = currentScene->AddRunDurationModel(pvModel,warn);
266 UImanager->ApplyCommand("/vis/sceneHandler/attach");
267 UImanager->SetVerboseLevel(keepVerbose);
268
269 if (successful) {
270 if (verbosity >= G4VisManager::confirmations) {
271 G4cout << "\"" << properties.fpTouchablePV->GetName()
272 << "\", copy no. " << properties.fCopyNo << " drawn"
273 << G4endl;
274 }
275 } else {
277 }
278 } else {
279 G4cout << "Touchable not found." << G4endl;
280 }
281 return;
282
283 } else if (command == fpCommandDump) {
284
287 if (properties.fpTouchablePV) {
288 // To handle paramaterisations we have to set the copy number
289 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
290 G4PhysicalVolumeModel tempPVModel
291 (properties.fpTouchablePV,
293 properties.fTouchableGlobalTransform,
294 nullptr, // Modelling parameters (not used)
295 true, // use full extent (prevents calculating own extent, which crashes)
296 properties.fTouchableBaseFullPVPath);
297 const std::map<G4String,G4AttDef>* attDefs = tempPVModel.GetAttDefs();
298 std::vector<G4AttValue>* attValues = tempPVModel.CreateCurrentAttValues();
299 G4cout << G4AttCheck(attValues,attDefs);
300 delete attValues;
301 G4Polyhedron* polyhedron =
303 G4cout << "\nLocal polyhedron coordinates:\n" << *polyhedron;
304 const G4Transform3D& transform = tempPVModel.GetCurrentTransform();
305 polyhedron->Transform(transform);
306 G4cout << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
307 } else {
308 G4cout << "Touchable not found." << G4endl;
309 }
310 return;
311
312 } else if (command == fpCommandExtentForField) {
313
316 if (properties.fpTouchablePV) {
317 G4VisExtent extent
318 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
319 extent.Transform(properties.fTouchableGlobalTransform);
320 fCurrentExtentForField = extent;
322 if (verbosity >= G4VisManager::confirmations) {
323 G4cout << "Extent for field set to " << extent
324 << "\nVolume for field has been cleared."
325 << G4endl;
326 }
328 DrawExtent(extent);
329 }
330 } else {
331 G4cout << "Touchable not found." << G4endl;
332 }
333 return;
334
335 } else if (command == fpCommandFindPath) {
336
337 G4String pvName;
338 G4int copyNo;
339 std::istringstream iss(newValue);
340 iss >> pvName >> copyNo;
341 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
342 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
343 transportationManager->GetWorldsIterator();
344 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
345 G4PhysicalVolumeModel searchModel (*iterWorld); // Unlimited depth.
346 G4ModelingParameters mp; // Default - no culling.
347 searchModel.SetModelingParameters (&mp);
348 // Find all instances at any position in the tree
349 G4PhysicalVolumesSearchScene searchScene (&searchModel, pvName, copyNo);
350 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
351 for (const auto& findings: searchScene.GetFindings()) {
352 findingsVector.push_back(findings);
353 }
354 }
355 for (const auto& findings: findingsVector) {
356 G4cout
357 << findings.fFoundBasePVPath
358 << ' ' << findings.fpFoundPV->GetName()
359 << ' ' << findings.fFoundPVCopyNo
360 << " (mother logical volume: "
361 << findings.fpFoundPV->GetMotherLogical()->GetName()
362 << ')'
363 << G4endl;
364 }
365 if (findingsVector.size()) {
366 G4cout
367 << "Use this to set a particular touchable with \"/vis/set/touchable <path>\""
368 << "\nor to see overlaps: \"/vis/drawLogicalVolume <mother-logical-volume-name>\""
369 << G4endl;
370 } else {
371 G4cout << pvName;
372 if (copyNo >= 0) G4cout << ':' << copyNo;
373 G4cout << " not found" << G4endl;
374 }
375
376 } else if (command == fpCommandLocalAxes) {
377
380 const G4double lengthMax = extent.GetExtentRadius()/2.;
381 const G4double intLog10LengthMax = std::floor(std::log10(lengthMax));
382 G4double length = std::pow(10,intLog10LengthMax);
383 if (5.*length < lengthMax) length *= 5.;
384 else if (2.*length < lengthMax) length *= 2.;
385 G4AxesModel axesModel(0.,0.,0.,length,transform);
386 axesModel.SetGlobalTag("LocalAxesModel");
388
389 } else if (command == fpCommandShowExtent) {
390
393 if (properties.fpTouchablePV) {
394 G4VisExtent extent
395 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
396 extent.Transform(properties.fTouchableGlobalTransform);
397 G4cout << extent << G4endl;
398 if (fpCommandShowExtent->GetNewBoolValue(newValue)) DrawExtent(extent);
399 } else {
400 G4cout << "Touchable not found." << G4endl;
401 }
402 return;
403
404 } else if (command == fpCommandVolumeForField) {
405
408 if (properties.fpTouchablePV) {
409 G4VisExtent extent
410 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
411 extent.Transform(properties.fTouchableGlobalTransform);
412 fCurrentExtentForField = extent;
416 if (verbosity >= G4VisManager::confirmations) {
417 G4cout
418 << "Volume for field set to " << properties.fpTouchablePV->GetName()
419 << ':' << properties.fCopyNo
420 << " at " << properties.fTouchableBaseFullPVPath
421 << G4endl;
422 }
424 DrawExtent(extent);
425 }
426 } else {
427 G4cout << "Touchable not found." << G4endl;
428 }
429 return;
430
431 } else {
432
433 if (verbosity >= G4VisManager::errors) {
434 G4cerr <<
435 "ERROR: G4VisCommandsTouchable::SetNewValue: unrecognised command."
436 << G4endl;
437 }
438 return;
439 }
440}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void DescribeYourselfTo(G4VGraphicsScene &) override
Definition: G4AxesModel.cc:208
G4VSolid * GetSolid() const
const G4Transform3D & GetCurrentTransform() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
void DescribeYourselfTo(G4VGraphicsScene &)
const std::map< G4String, G4AttDef > * GetAttDefs() const
const std::vector< Findings > & GetFindings() const
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition: G4Scene.cc:158
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
const G4String & GetName() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const
static G4bool GetNewBoolValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
G4int GetVerboseLevel() const
Definition: G4UImanager.hh:200
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
void SetVerboseLevel(G4int val)
Definition: G4UImanager.hh:199
void SetModelingParameters(const G4ModelingParameters *)
void SetGlobalTag(const G4String &)
virtual void SetCopyNo(G4int CopyNo)=0
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:682
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:705
const G4String & GetName() const
const G4ViewParameters & GetViewParameters() const
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
static std::vector< G4PhysicalVolumesSearchScene::Findings > fCurrrentPVFindingsForField
static G4VisManager * fpVisManager
static G4VisExtent fCurrentExtentForField
void DrawExtent(const G4VisExtent &)
void InterpolateToNewView(G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String exportString="")
void SetViewParameters(G4VViewer *viewer, const G4ViewParameters &viewParams)
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
void SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
void SetZoomFactor(G4double zoomFactor)
G4UIcmdWithoutParameter * fpCommandDraw
G4UIcmdWithABool * fpCommandExtentForField
G4String GetCurrentValue(G4UIcommand *command)
G4UIcmdWithoutParameter * fpCommandDump
G4UIcmdWithoutParameter * fpCommandLocalAxes
G4UIcmdWithABool * fpCommandVolumeForField
void SetNewValue(G4UIcommand *command, G4String newValue)
G4UIcmdWithoutParameter * fpCommandCentreAndZoomInOn
G4UIcmdWithoutParameter * fpCommandCentreOn
G4UIcmdWithABool * fpCommandShowExtent
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
G4VisExtent & Transform(const G4Transform3D &)
Definition: G4VisExtent.cc:102
const G4Point3D & GetExtentCentre() const
Definition: G4VisExtent.cc:65
G4Scene * GetCurrentScene() const
G4VSceneHandler * GetCurrentSceneHandler() const
static Verbosity GetVerbosity()
HepPolyhedron & Transform(const G4Transform3D &t)
G4PhysicalVolumeModel::TouchableProperties FindTouchableProperties(G4ModelingParameters::PVNameCopyNoPath path)
G4bool transform(G4String &input, const G4String &type)
G4ModelingParameters::PVNameCopyNoPath fTouchablePath
std::vector< G4PhysicalVolumeNodeID > fTouchableBaseFullPVPath