Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VisCommandViewerCentreOn Class Reference

#include <G4VisCommandsViewer.hh>

Inheritance diagram for G4VisCommandViewerCentreOn:

Public Member Functions

 G4VisCommandViewerCentreOn ()
virtual ~G4VisCommandViewerCentreOn ()
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 56 of file G4VisCommandsViewer.hh.

Constructor & Destructor Documentation

◆ G4VisCommandViewerCentreOn()

G4VisCommandViewerCentreOn::G4VisCommandViewerCentreOn ( )

Definition at line 142 of file G4VisCommandsViewer.cc.

142 {
143 G4bool omitable;
144 fpCommandCentreAndZoomInOn = new G4UIcommand ("/vis/viewer/centreAndZoomInOn", this);
145 fpCommandCentreAndZoomInOn->SetGuidance
146 ("Centre and zoom in on the given physical volume.");
147 fpCommandCentreAndZoomInOn->SetGuidance
148 ("The names of all volumes in all worlds are matched against pv-name. If"
149 "\ncopy-no is supplied, it matches the copy number too. If pv-name is of the"
150 "\nform \"/regexp/\", where regexp is a regular expression (see C++ regex),"
151 "\nthe match uses the usual rules of regular expression matching."
152 "\nOtherwise an exact match is required."
153 "\nFor example, \"/Shap/\" matches \"Shape1\" and \"Shape2\".");
154 fpCommandCentreAndZoomInOn->SetGuidance
155 ("It may help to see a textual representation of the geometry hierarchy of"
156 "\nthe worlds. Try \"/vis/drawTree [worlds]\"");
157 fpCommandCentreAndZoomInOn->SetGuidance
158 ("If there are more than one matching physical volumes they will all be"
159 "\nincluded. If this is not what you want, and what you want is to centre on a"
160 "\nparticular touchable, then select the touchable (\"/vis/set/touchable\") and"
161 "\nuse \"/vis/touchable/centreOn\". (You may need \"/vis/touchable/findPath\".)");
162 G4UIparameter* parameter;
163 parameter = new G4UIparameter("pv-name",'s',omitable = false);
164 parameter->SetGuidance ("Physical volume name.");
165 fpCommandCentreAndZoomInOn->SetParameter(parameter);
166 parameter = new G4UIparameter("copy-no",'i',omitable = true);
167 parameter->SetDefaultValue (-1);
168 parameter->SetGuidance ("Copy number. -1 means any or all copy numbers");
169 fpCommandCentreAndZoomInOn->SetParameter(parameter);
170
171 fpCommandCentreOn = new G4UIcommand ("/vis/viewer/centreOn", this);
172 fpCommandCentreOn->SetGuidance ("Centre the view on the given physical volume.");
173 // Pick up additional guidance from /vis/viewer/centreAndZoomInOn
174 CopyGuidanceFrom(fpCommandCentreAndZoomInOn,fpCommandCentreOn,1);
175 // Pick up parameters from /vis/viewer/centreAndZoomInOn
176 CopyParametersFrom(fpCommandCentreAndZoomInOn,fpCommandCentreOn);
177}
bool G4bool
Definition G4Types.hh:86
void SetDefaultValue(const char *theDefaultValue)
void SetGuidance(const char *theGuidance)
void CopyParametersFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd)
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)

◆ ~G4VisCommandViewerCentreOn()

G4VisCommandViewerCentreOn::~G4VisCommandViewerCentreOn ( )
virtual

Definition at line 179 of file G4VisCommandsViewer.cc.

179 {
180 delete fpCommandCentreAndZoomInOn;
181 delete fpCommandCentreOn;
182}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandViewerCentreOn::GetCurrentValue ( G4UIcommand * command)
virtual

Reimplemented from G4UImessenger.

Definition at line 184 of file G4VisCommandsViewer.cc.

184 {
185 return "";
186}

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 188 of file G4VisCommandsViewer.cc.

188 {
189
190 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
191 G4bool warn = verbosity >= G4VisManager::warnings;
192
193 G4VViewer* currentViewer = fpVisManager -> GetCurrentViewer ();
194 if (!currentViewer) {
195 if (verbosity >= G4VisManager::errors) {
196 G4warn <<
197 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
198 << G4endl;
199 }
200 return;
201 }
202
203 G4String pvName;
204 G4int copyNo;
205 std::istringstream is (newValue);
206 is >> pvName >> copyNo;
207
208 // Find physical volumes
209 G4TransportationManager* transportationManager =
211 std::size_t nWorlds = transportationManager->GetNoWorlds();
212 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
213 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
214 transportationManager->GetWorldsIterator();
215 for (std::size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
216 G4PhysicalVolumeModel searchModel (*iterWorld); // Unlimited depth.
217 G4ModelingParameters mp; // Default - no culling.
218 searchModel.SetModelingParameters (&mp);
219 // Find all instances at any position in the tree
220 G4PhysicalVolumesSearchScene searchScene (&searchModel, pvName, copyNo);
221 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
222 for (const auto& findings: searchScene.GetFindings()) {
223 findingsVector.push_back(findings);
224 }
225 }
226
227 if (findingsVector.empty()) {
228 if (verbosity >= G4VisManager::warnings) {
229 G4warn
230 << "WARNING: Volume \"" << pvName << "\" ";
231 if (copyNo > 0) {
232 G4warn << "copy number " << copyNo;
233 }
234 G4warn << " not found." << G4endl;
235 }
236 return;
237 }
238
239 // A vector of found paths so that we can highlight (twinkle) the found volume(s).
240 std::vector<std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>> foundPaths;
241
242 // Use a temporary scene in order to find vis extent
243 G4Scene tempScene("Centre Scene");
244 G4bool successfullyAdded = true;
245 for (const auto& findings: findingsVector) {
246 // To handle paramaterisations we have to set the copy number
247 findings.fpFoundPV->SetCopyNo(findings.fFoundPVCopyNo);
248 // Create a temporary physical volume model.
249 // They have to be created on the heap because they have
250 // to hang about long enough to be conflated.
251 G4PhysicalVolumeModel* tempPVModel = new G4PhysicalVolumeModel
252 (findings.fpFoundPV,
253 0, // Only interested in top volume
254 findings.fFoundObjectTransformation,
255 0, // No modelling parameters (these are set later by the scene handler).
256 true, // Use full extent
257 findings.fFoundBasePVPath);
258 // ...and add it to the scene.
259 auto successful = tempScene.AddRunDurationModel(tempPVModel,warn);
260 if (!successful) {
261 successfullyAdded = false;
262 continue;
263 }
264 if (verbosity >= G4VisManager::parameters) {
265 G4cout << "\"" << findings.fpFoundPV->GetName()
266 << "\", copy no. " << findings.fFoundPVCopyNo
267 << ",\n found in searched volume \""
268 << findings.fpSearchPV->GetName()
269 << "\" at depth " << findings.fFoundDepth
270 << ",\n base path: \"" << findings.fFoundBasePVPath
271 << ",\n has been added to temporary scene \"" << tempScene.GetName() << "\"."
272 << G4endl;
273 }
274 foundPaths.push_back(findings.fFoundFullPVPath);
275 }
276 // Delete temporary physical volume models
277 for (const auto& sceneModel: tempScene.GetRunDurationModelList()) {
278 delete sceneModel.fpModel;
279 }
280 if (!successfullyAdded) return;
281
282 // Relevant results
283 const G4VisExtent& newExtent = tempScene.GetExtent();
284 const G4ThreeVector& newTargetPoint = newExtent.GetExtentCentre();
285
286 G4Scene* currentScene = currentViewer->GetSceneHandler()->GetScene();
287 G4ViewParameters saveVP = currentViewer->GetViewParameters();
288 G4ViewParameters newVP = saveVP;
289 if (command == fpCommandCentreAndZoomInOn) {
290 // Calculate the new zoom factor
291 const G4double zoomFactor
292 = currentScene->GetExtent().GetExtentRadius()/newExtent.GetExtentRadius();
293 newVP.SetZoomFactor(zoomFactor);
294 }
295 // Change the target point
296 const G4Point3D& standardTargetPoint = currentScene->GetStandardTargetPoint();
297 newVP.SetCurrentTargetPoint(newTargetPoint - standardTargetPoint);
298
299 // If this particular view is simple enough
300 if (currentViewer->GetKernelVisitElapsedTimeSeconds() < 0.1) {
301 // Interpolate
302 auto keepVisVerbosity = fpVisManager->GetVerbosity();
303 fpVisManager->SetVerboseLevel(G4VisManager::errors);
304 if (newVP != saveVP) InterpolateToNewView(currentViewer, saveVP, newVP);
305 // ...and twinkle
306 Twinkle(currentViewer,newVP,foundPaths);
307 fpVisManager->SetVerboseLevel(keepVisVerbosity);
308 }
309
310 if (verbosity >= G4VisManager::confirmations) {
311 G4cout
312 << "Viewer \"" << currentViewer->GetName()
313 << "\" centred ";
314 if (fpCommandCentreAndZoomInOn) {
315 G4cout << "and zoomed in";
316 }
317 G4cout << " on physical volume(s) \"" << pvName << '\"'
318 << G4endl;
319 }
320
321 SetViewParameters(currentViewer, newVP);
322}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
#define G4warn
Definition G4Scene.cc:41
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const
G4Scene * GetScene() const
const G4String & GetName() const
const G4ViewParameters & GetViewParameters() const
G4double GetKernelVisitElapsedTimeSeconds() const
G4VSceneHandler * GetSceneHandler() const
void InterpolateToNewView(G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
static G4VisManager * fpVisManager
void SetViewParameters(G4VViewer *viewer, const G4ViewParameters &viewParams)
void Twinkle(G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
void SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
void SetZoomFactor(G4double zoomFactor)
G4double GetExtentRadius() const
const G4Point3D & GetExtentCentre() const

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