Geant4 11.1.1
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 ()
 
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 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 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
 

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 fpCommandVolumeForField = new G4UIcmdWithABool("/vis/touchable/volumeForField",this);
124 fpCommandVolumeForField->SetGuidance("Set volume for field.");
125 fpCommandVolumeForField->SetGuidance("If parameter == true, also draw.");
126 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
127 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandVolumeForField,1);
128 fpCommandVolumeForField->SetParameterName("draw", omitable = true);
129 fpCommandVolumeForField->SetDefaultValue(false);
130}
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:157
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)

◆ ~G4VisCommandsTouchable()

G4VisCommandsTouchable::~G4VisCommandsTouchable ( )
virtual

Definition at line 132 of file G4VisCommandsTouchable.cc.

132 {
133 delete fpCommandVolumeForField;
134 delete fpCommandShowExtent;
135 delete fpCommandLocalAxes;
136 delete fpCommandFindPath;
137 delete fpCommandExtentForField;
138 delete fpCommandDump;
139 delete fpCommandDraw;
140 delete fpCommandCentreAndZoomInOn;
141 delete fpCommandCentreOn;
142}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandsTouchable::GetCurrentValue ( G4UIcommand command)
virtual

Reimplemented from G4UImessenger.

Definition at line 144 of file G4VisCommandsTouchable.cc.

144 {
145 return "";
146}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 148 of file G4VisCommandsTouchable.cc.

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