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

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