Geant4 10.7.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 ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VVisCommand
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, 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="")
 
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 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 DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (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 = 0
 
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
 

Detailed Description

Definition at line 39 of file G4VisCommandsTouchable.hh.

Constructor & Destructor Documentation

◆ G4VisCommandsTouchable()

G4VisCommandsTouchable::G4VisCommandsTouchable ( )

Definition at line 43 of file G4VisCommandsTouchable.cc.

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

◆ ~G4VisCommandsTouchable()

G4VisCommandsTouchable::~G4VisCommandsTouchable ( )
virtual

Definition at line 120 of file G4VisCommandsTouchable.cc.

120 {
121 delete fpCommandVolumeForField;
122 delete fpCommandShowExtent;
123 delete fpCommandFindPath;
124 delete fpCommandExtentForField;
125 delete fpCommandDump;
126 delete fpCommandDraw;
127 delete fpCommandCentreAndZoomInOn;
128 delete fpCommandCentreOn;
129}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandsTouchable::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 131 of file G4VisCommandsTouchable.cc.

131 {
132 return "";
133}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 135 of file G4VisCommandsTouchable.cc.

137{
139 G4bool warn = verbosity >= G4VisManager::warnings;
140
142
143 G4TransportationManager* transportationManager =
145
146 size_t nWorlds = transportationManager->GetNoWorlds();
147
148 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
149 if (!world) {
150 if (verbosity >= G4VisManager::errors) {
151 G4cerr <<
152 "ERROR: G4VisCommandsTouchable::SetNewValue:"
153 "\n No world. Maybe the geometry has not yet been defined."
154 "\n Try \"/run/initialize\""
155 << G4endl;
156 }
157 return;
158 }
159
160 G4VViewer* currentViewer = fpVisManager -> GetCurrentViewer ();
161 if (!currentViewer) {
162 if (verbosity >= G4VisManager::errors) {
163 G4cerr <<
164 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
165 << G4endl;
166 }
167 return;
168 }
169
170 G4Scene* currentScene = fpVisManager->GetCurrentScene();
171 if (!currentScene) {
172 if (verbosity >= G4VisManager::errors) {
173 G4cerr <<
174 "ERROR: No current scene - \"/vis/scene/list\" to see possibilities."
175 << G4endl;
176 }
177 return;
178 }
179
180 if (command == fpCommandCentreOn || command == fpCommandCentreAndZoomInOn) {
181
184 if (properties.fpTouchablePV) {
185 // To handle parameterisations, set copy number
186 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
187 G4PhysicalVolumeModel tempPVModel
188 (properties.fpTouchablePV,
190 properties.fTouchableGlobalTransform,
191 nullptr, // Modelling parameters (not used)
192 true, // use full extent (prevents calculating own extent, which crashes)
193 properties.fTouchableBaseFullPVPath);
194 // Use a temporary scene in order to find vis extent
195 G4Scene tempScene("Centre Scene");
196 G4bool successful = tempScene.AddRunDurationModel(&tempPVModel,warn);
197 if (successful) {
198 if (verbosity >= G4VisManager::confirmations) {
199 G4cout
201 << ",\n has been added to temporary scene \"" << tempScene.GetName() << "\"."
202 << G4endl;
203 }
204 }
205 const G4VisExtent& newExtent = tempScene.GetExtent();
206 const G4ThreeVector& newTargetPoint = newExtent.GetExtentCentre();
207 G4ViewParameters saveVP = currentViewer->GetViewParameters();
208 G4ViewParameters newVP = saveVP;
209 if (command == fpCommandCentreAndZoomInOn) {
210 // Calculate the new zoom factor
211 const G4double zoomFactor
212 = currentScene->GetExtent().GetExtentRadius()/newExtent.GetExtentRadius();
213 newVP.SetZoomFactor(zoomFactor);
214 }
215 // Change the target point
216 const G4Point3D& standardTargetPoint = currentScene->GetStandardTargetPoint();
217 newVP.SetCurrentTargetPoint(newTargetPoint - standardTargetPoint);
218 // Interpolate
219 InterpolateToNewView(currentViewer, saveVP, newVP);
220 if (verbosity >= G4VisManager::confirmations) {
221 G4cout
222 << "Viewer \"" << currentViewer->GetName()
223 << "\" centred ";
224 if (fpCommandCentreAndZoomInOn) {
225 G4cout << "and zoomed in";
226 }
228 << G4endl;
229 }
230 SetViewParameters(currentViewer, newVP);
231 } else {
232 G4cout << "Touchable not found." << G4endl;
233 }
234 return;
235
236 } else if (command == fpCommandDraw) {
237
240 if (properties.fpTouchablePV) {
241 // To handle paramaterisations we have to set the copy number
242 properties.fpTouchablePV->SetCopyNo(properties.fCopyNo);
244 (properties.fpTouchablePV,
246 properties.fTouchableGlobalTransform,
247 nullptr, // Modelling parameters (not used)
248 true, // use full extent (prevents calculating own extent, which crashes)
249 properties.fTouchableBaseFullPVPath);
250
251 G4int keepVerbose = UImanager->GetVerboseLevel();
252 G4int newVerbose(0);
253 if (keepVerbose >= 2 || verbosity >= G4VisManager::confirmations)
254 newVerbose = 2;
255 UImanager->SetVerboseLevel(newVerbose);
256 UImanager->ApplyCommand("/vis/scene/create");
257 currentScene = fpVisManager->GetCurrentScene(); // New current scene
258 G4bool successful = currentScene->AddRunDurationModel(pvModel,warn);
259 UImanager->ApplyCommand("/vis/sceneHandler/attach");
260 UImanager->SetVerboseLevel(keepVerbose);
261
262 if (successful) {
263 if (verbosity >= G4VisManager::confirmations) {
264 G4cout << "\"" << properties.fpTouchablePV->GetName()
265 << "\", copy no. " << properties.fCopyNo << " drawn"
266 << G4endl;
267 }
268 } else {
270 }
271 } else {
272 G4cout << "Touchable not found." << G4endl;
273 }
274 return;
275
276 } else if (command == fpCommandDump) {
277
280 if (properties.fpTouchablePV) {
281 // To handle paramaterisations we have to set the 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 const std::map<G4String,G4AttDef>* attDefs = tempPVModel.GetAttDefs();
291 std::vector<G4AttValue>* attValues = tempPVModel.CreateCurrentAttValues();
292 G4cout << G4AttCheck(attValues,attDefs);
293 delete attValues;
294 G4Polyhedron* polyhedron =
296 G4cout << "\nLocal polyhedron coordinates:\n" << *polyhedron;
297 G4Transform3D* transform = tempPVModel.GetCurrentTransform();
298 polyhedron->Transform(*transform);
299 G4cout << "\nGlobal polyhedron coordinates:\n" << *polyhedron;
300 } else {
301 G4cout << "Touchable not found." << G4endl;
302 }
303 return;
304
305 } else if (command == fpCommandExtentForField) {
306
309 if (properties.fpTouchablePV) {
310 G4VisExtent extent
311 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
312 extent.Transform(properties.fTouchableGlobalTransform);
313 fCurrentExtentForField = extent;
315 if (verbosity >= G4VisManager::confirmations) {
316 G4cout << "Extent for field set to " << extent
317 << "\nVolume for field has been cleared."
318 << G4endl;
319 }
320 if (fpCommandExtentForField->GetNewBoolValue(newValue)) {
321 DrawExtent(extent);
322 }
323 } else {
324 G4cout << "Touchable not found." << G4endl;
325 }
326 return;
327
328 } else if (command == fpCommandFindPath) {
329
330 G4String pvName;
331 G4int copyNo;
332 std::istringstream iss(newValue);
333 iss >> pvName >> copyNo;
334 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
335 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
336 transportationManager->GetWorldsIterator();
337 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
338 G4PhysicalVolumeModel searchModel (*iterWorld); // Unlimited depth.
339 G4ModelingParameters mp; // Default - no culling.
340 searchModel.SetModelingParameters (&mp);
341 G4PhysicalVolumesSearchScene searchScene (&searchModel, pvName, copyNo);
342 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
343 for (const auto& findings: searchScene.GetFindings()) {
344 findingsVector.push_back(findings);
345 }
346 }
347 for (const auto& findings: findingsVector) {
348 G4cout
349 << findings.fFoundBasePVPath
350 << ' ' << findings.fpFoundPV->GetName()
351 << ' ' << findings.fFoundPVCopyNo
352 << " (mother logical volume: "
353 << findings.fpFoundPV->GetMotherLogical()->GetName()
354 << ')'
355 << G4endl;
356 }
357 if (findingsVector.size()) {
358 G4cout
359 << "Use this to set a particular touchable with \"/vis/set/touchable <path>\""
360 << "\nor to see overlaps: \"/vis/drawLogicalVolume <mother-logical-volume-name>\""
361 << G4endl;
362 } else {
363 G4cout << pvName;
364 if (copyNo >= 0) G4cout << ':' << copyNo;
365 G4cout << " not found" << G4endl;
366 }
367
368 } else if (command == fpCommandShowExtent) {
369
372 if (properties.fpTouchablePV) {
373 G4VisExtent extent
374 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
375 extent.Transform(properties.fTouchableGlobalTransform);
376 G4cout << extent << G4endl;
377 if (fpCommandShowExtent->GetNewBoolValue(newValue)) DrawExtent(extent);
378 } else {
379 G4cout << "Touchable not found." << G4endl;
380 }
381 return;
382
383 } else if (command == fpCommandVolumeForField) {
384
387 if (properties.fpTouchablePV) {
388 G4VisExtent extent
389 = properties.fpTouchablePV->GetLogicalVolume()->GetSolid()->GetExtent();
390 extent.Transform(properties.fTouchableGlobalTransform);
391 fCurrentExtentForField = extent;
395 if (verbosity >= G4VisManager::confirmations) {
396 G4cout
397 << "Volume for field set to " << properties.fpTouchablePV->GetName()
398 << ':' << properties.fCopyNo
399 << " at " << properties.fTouchableBaseFullPVPath
400 << G4endl;
401 }
402 if (fpCommandVolumeForField->GetNewBoolValue(newValue)) {
403 DrawExtent(extent);
404 }
405 } else {
406 G4cout << "Touchable not found." << G4endl;
407 }
408 return;
409
410 } else {
411
412 if (verbosity >= G4VisManager::errors) {
413 G4cerr <<
414 "ERROR: G4VisCommandsTouchable::SetNewValue: unrecognised command."
415 << G4endl;
416 }
417 return;
418 }
419}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VSolid * GetSolid() const
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition: G4Scene.cc:50
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const
static G4bool GetNewBoolValue(const char *paramString)
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
G4int GetVerboseLevel() const
Definition: G4UImanager.hh:193
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
void SetVerboseLevel(G4int val)
Definition: G4UImanager.hh:192
virtual void SetCopyNo(G4int CopyNo)=0
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4VisExtent GetExtent() const
Definition: G4VSolid.cc:670
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:693
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 SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
void SetZoomFactor(G4double zoomFactor)
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
static Verbosity GetVerbosity()
G4PhysicalVolumeModel::TouchableProperties FindTouchableProperties(G4ModelingParameters::PVNameCopyNoPath path)
G4ModelingParameters::PVNameCopyNoPath fTouchablePath
std::vector< G4PhysicalVolumeNodeID > fTouchableBaseFullPVPath

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