Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VVisCommand.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// Base class for visualization commands - John Allison 9th August 1998
29// It is really a messenger - we have one command per messenger.
30
31#include "G4VVisCommand.hh"
32
33#include "G4UIcommand.hh"
34#include "G4UImanager.hh"
35#include "G4UnitsTable.hh"
36#include <sstream>
37#include <cctype>
38
40#include "G4LogicalVolume.hh"
41
42#define G4warn G4cout
43
50// Not yet used: G4VisAttributes::LineStyle G4VVisCommand::fCurrentLineStyle = G4VisAttributes::unbroken;
51// Not yet used: G4VMarker::FillStyle G4VVisCommand::fCurrentFillStyle = G4VMarker::filled;
52// Not yet used: G4VMarker::SizeType G4VVisCommand::fCurrentSizeType = G4VMarker::screen;
55std::vector<G4PhysicalVolumesSearchScene::Findings> G4VVisCommand::fCurrrentPVFindingsForField;
58
60
62
64
66{
67 return fpVisManager;
68}
69
71{
72 fpVisManager = pVisManager;
73}
74
76{
77 return fCurrentTextColour;
78}
79
81(G4double x, G4double y, const char * unitName)
82{
83 G4double uv = G4UIcommand::ValueOf(unitName);
84
85 std::ostringstream oss;
86 oss << x/uv << " " << y/uv << " " << unitName;
87 return oss.str();
88}
89
91 G4double& xval,
92 G4double& yval)
93{
94 G4double x, y;
95 G4String unit;
96
97 std::istringstream is(paramString);
98 is >> x >> y >> unit;
99
101 xval = x*G4UIcommand::ValueOf(unit);
102 yval = y*G4UIcommand::ValueOf(unit);
103 } else {
105 if (verbosity >= G4VisManager::errors) {
106 G4warn << "ERROR: Unrecognised unit" << G4endl;
107 }
108 return false;
109 }
110
111 return true;
112}
113
115{
116 static G4String guidance
117 ("Accepts (a) RGB triplet. e.g., \".3 .4 .5\", or"
118 "\n (b) string such as \"white\", \"black\", \"grey\", \"red\"...or"
119 "\n (c) an additional number for opacity, e.g., \".3 .4 .5 .6\""
120 "\n or \"grey ! ! .6\" (note \"!\"'s for unused parameters).");
121 return guidance;
122}
123
125(G4Colour& colour,
126 const G4String& redOrString, G4double green, G4double blue, G4double opacity)
127{
128 // Note: colour is supplied by the caller and some or all of its components
129 // may act as default.
130 //
131 // Note: redOrString is either a number or string. If a string it must be
132 // one of the recognised colours.
133 //
134 // Thus the arguments can be, for example:
135 // (colour,"red",...,...,0.5): will give the colour red with opacity 0.5 (the
136 // third and fourth arguments are ignored), or
137 // (1.,0.,0.,0.5): this also will be red with opacity 0.5.
138
140
141 const std::size_t iPos0 = 0;
142 if (std::isalpha(redOrString[iPos0])) {
143
144 // redOrString is probably alphabetic characters defining the colour
145 if (!G4Colour::GetColour(redOrString, colour)) {
146 // Not a recognised string
147 if (verbosity >= G4VisManager::warnings) {
148 G4warn << "WARNING: Colour \"" << redOrString
149 << "\" not found. Defaulting to " << colour
150 << G4endl;
151 }
152 return;
153 } else {
154 // It was a recognised string. Now add opacity.
155 colour.SetAlpha(opacity);
156 return;
157 }
158
159 } else {
160
161 // redOrString is probably numeric defining the red component
162 std::istringstream iss(redOrString);
163 G4double red;
164 iss >> red;
165 if (iss.fail()) {
166 if (verbosity >= G4VisManager::warnings) {
167 G4warn << "WARNING: String \"" << redOrString
168 << "\" cannot be parsed. Defaulting to " << colour
169 << G4endl;
170 }
171 return;
172 } else {
173 colour = G4Colour(red,green,blue,opacity);
174 return;
175 }
176
177 }
178}
179
181(const G4String& where,
182 const G4String& unit,
183 const G4String& category,
184 G4double& value)
185{
186 // Return false if there's a problem
187
189
190 G4bool success = true;
192 if (verbosity >= G4VisManager::warnings) {
193 G4warn << where
194 << "\n Unit \"" << unit << "\" not defined"
195 << G4endl;
196 }
197 success = false;
198 } else if (G4UnitDefinition::GetCategory(unit) != category) {
199 if (verbosity >= G4VisManager::warnings) {
200 G4warn << where
201 << "\n Unit \"" << unit << "\" not a unit of " << category;
202 if (category == "Volumic Mass") G4warn << " (density)";
203 G4warn << G4endl;
204 }
205 success = false;
206 } else {
207 value = G4UnitDefinition::GetValueOf(unit);
208 }
209 return success;
210}
211
213{
215
216 if (!pScene) {
217 if (verbosity >= G4VisManager::warnings) {
218 G4warn << "WARNING: Scene pointer is null."
219 << G4endl;
220 }
221 return;
222 }
223
224 G4VSceneHandler* pSceneHandler = fpVisManager -> GetCurrentSceneHandler();
225 if (!pSceneHandler) {
226 if (verbosity >= G4VisManager::warnings) {
227 G4warn << "WARNING: Scene handler not found." << G4endl;
228 }
229 return;
230 }
231
232 // Scene has changed. If it is the scene of the currrent scene handler
233 // refresh viewers of all scene handlers using this scene. If not, it may be
234 // a scene that the user is building up before attaching to a scene handler,
235 // so do nothing.
236 if (pScene == pSceneHandler->GetScene()) {
237 G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/scene/notifyHandlers");
238 }
239
240}
241
243{
245 G4VViewer* viewer = fpVisManager -> GetCurrentViewer ();
246
247 if (!viewer) {
248 if (verbosity >= G4VisManager::errors) {
249 G4warn <<
250 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
251 << G4endl;
252 }
253 return false;
254 }
255
256 return true;
257}
258
260(G4VisManager::Verbosity verbosity) {
261 // Some frequently used error printing...
262 if (verbosity >= G4VisManager::warnings) {
263 G4warn <<
264 "WARNING: For some reason, possibly mentioned above, it has not been"
265 "\n possible to add to the scene."
266 << G4endl;
267 }
268}
269
271(G4VViewer* viewer, const G4ViewParameters& viewParams) {
272 viewer->SetViewParameters(viewParams);
273 RefreshIfRequired(viewer);
274}
275
278 G4VSceneHandler* sceneHandler = viewer->GetSceneHandler();
279 const G4ViewParameters& viewParams = viewer->GetViewParameters();
280 if (sceneHandler && sceneHandler->GetScene()) {
281 if (viewParams.IsAutoRefresh()) {
282 G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/refresh");
283 }
284 else {
285 if (verbosity >= G4VisManager::warnings) {
286 G4warn << "Issue /vis/viewer/refresh or flush to see effect."
287 << G4endl;
288 }
289 }
290 }
291}
292
294(G4VViewer* currentViewer,
295 std::vector<G4ViewParameters> viewVector,
296 const G4int nInterpolationPoints,
297 const G4int waitTimePerPointmilliseconds,
298 const G4String exportString)
299{
300 const G4int safety = (G4int)viewVector.size()*nInterpolationPoints;
301 G4int safetyCount = 0;
302 do {
303 G4ViewParameters* vp =
304 G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
305 if (!vp) break; // Finished.
306 currentViewer->SetViewParameters(*vp);
307 currentViewer->RefreshView();
308 if (exportString == "export" &&
309 G4StrUtil::contains(currentViewer->GetName(), "OpenGL")) {
310 G4UImanager::GetUIpointer()->ApplyCommand("/vis/ogl/export");
311 }
312 currentViewer->ShowView();
313 if (waitTimePerPointmilliseconds > 0)
314 std::this_thread::sleep_for(std::chrono::milliseconds(waitTimePerPointmilliseconds));
315 } while (safetyCount++ < safety); // Loop checking, 16.02.2016, J.Allison
316}
317
319(G4VViewer* currentViewer,
320 const G4ViewParameters& oldVP,
321 const G4ViewParameters& newVP,
322 const G4int nInterpolationPoints,
323 const G4int waitTimePerPointmilliseconds,
324 const G4String exportString)
325{
326 std::vector<G4ViewParameters> viewVector;
327 viewVector.push_back(oldVP);
328 viewVector.push_back(oldVP);
329 viewVector.push_back(newVP);
330 viewVector.push_back(newVP);
331
333 (currentViewer,
334 viewVector,
335 nInterpolationPoints,
336 waitTimePerPointmilliseconds,
337 exportString);
338}
339
341// Twinkles the touchables in paths
342// /vis/viewer/centreOn to see its effect
343 (G4VViewer* currentViewer,
344 const G4ViewParameters& baseVP,
345 const std::vector<std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>>& paths)
346{
347 // Copy view parameters to temporary variables ready for adding VisAttributes Modifiers (VAMs)
348 auto loVP = baseVP; // For black and solid VAMs
349 auto hiVP = baseVP; // For white and solid VAMs
350
351 // Modify them with vis attribute modifiers (VAMs)
352 for (const auto& path: paths) {
353 const auto& touchable = path.back().GetPhysicalVolume();
354 auto loVisAtts
355 = *(currentViewer->GetApplicableVisAttributes
356 (touchable->GetLogicalVolume()->GetVisAttributes()));
357 auto hiVisAtts = loVisAtts;
358 loVisAtts.SetColour(G4Colour::Black());
359 loVisAtts.SetForceSolid();
360 hiVisAtts.SetColour(G4Colour::White());
361 hiVisAtts.SetForceSolid();
362 auto pvNameCopyNoPath
364
366 (loVisAtts, G4ModelingParameters::VASColour, pvNameCopyNoPath);
367 loVP.AddVisAttributesModifier(loVAMColour);
369 (loVisAtts, G4ModelingParameters::VASForceSolid, pvNameCopyNoPath);
370 loVP.AddVisAttributesModifier(loVAMStyle);
371
373 (hiVisAtts, G4ModelingParameters::VASColour, pvNameCopyNoPath);
374 hiVP.AddVisAttributesModifier(hiVAMColour);
376 (hiVisAtts, G4ModelingParameters::VASForceSolid, pvNameCopyNoPath);
377 hiVP.AddVisAttributesModifier(hiVAMStyle);
378 }
379
380 // Twinkle
381 std::vector<G4ViewParameters> viewVector;
382 viewVector.push_back(loVP);
383 viewVector.push_back(hiVP);
384 viewVector.push_back(loVP);
385 viewVector.push_back(hiVP);
386 viewVector.push_back(loVP);
387 viewVector.push_back(hiVP);
388 viewVector.push_back(loVP);
389 viewVector.push_back(hiVP);
390 viewVector.push_back(loVP);
391 viewVector.push_back(hiVP);
392 // Just 5 interpolation points for a reasonable twinkle rate
393 InterpolateViews(currentViewer,viewVector,5);
394}
395
397(const G4UIcommand* fromCmd, G4UIcommand* toCmd, G4int startLine)
398{
399 if (fromCmd && toCmd) {
400 const G4int nGuideEntries = (G4int)fromCmd->GetGuidanceEntries();
401 for (G4int i = startLine; i < nGuideEntries; ++i) {
402 const G4String& guidance = fromCmd->GetGuidanceLine(i);
403 toCmd->SetGuidance(guidance);
404 }
405 }
406}
407
409(const G4UIcommand* fromCmd, G4UIcommand* toCmd)
410{
411 if (fromCmd && toCmd) {
412 const G4int nParEntries = (G4int)fromCmd->GetParameterEntries();
413 for (G4int i = 0; i < nParEntries; ++i) {
414 G4UIparameter* parameter = new G4UIparameter(*(fromCmd->GetParameter(i)));
415 toCmd->SetParameter(parameter);
416 }
417 }
418}
419
421{
422 if (fpVisManager) {
423 const G4double halfX = (extent.GetXmax() - extent.GetXmin()) / 2.;
424 const G4double halfY = (extent.GetYmax() - extent.GetYmin()) / 2.;
425 const G4double halfZ = (extent.GetZmax() - extent.GetZmin()) / 2.;
426 if (halfX > 0. && halfY > 0. && halfZ > 0.) {
427 const G4Box box("vis_extent",halfX,halfY,halfZ);
428 const G4VisAttributes visAtts(G4Color::Red());
429 const G4Point3D& centre = extent.GetExtentCenter();
430 fpVisManager->Draw(box,visAtts,G4Translate3D(centre));
431 }
432 }
433}
434
436(G4ViewParameters& target, const G4ViewParameters& from)
437{
438 // Copy view parameters pertaining only to camera
442 target.SetUpVector (from.GetUpVector());
443 target.SetFieldHalfAngle (from.GetFieldHalfAngle());
444 target.SetZoomFactor (from.GetZoomFactor());
445 target.SetScaleFactor (from.GetScaleFactor());
447 target.SetDolly (from.GetDolly());
448}
#define G4warn
Definition: G4Scene.cc:41
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
Definition: G4Box.hh:56
static G4Colour White()
Definition: G4Colour.hh:156
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition: G4Colour.cc:155
void SetAlpha(G4double)
Definition: G4Colour.cc:70
static G4Colour Red()
Definition: G4Colour.hh:161
static G4Colour Black()
Definition: G4Colour.hh:159
static G4Colour Blue()
Definition: G4Colour.hh:163
static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath(const std::vector< G4PhysicalVolumeNodeID > &)
Layout
Definition: G4Text.hh:76
@ left
Definition: G4Text.hh:76
std::size_t GetParameterEntries() const
Definition: G4UIcommand.hh:139
const G4String & GetGuidanceLine(G4int i) const
Definition: G4UIcommand.hh:133
G4UIparameter * GetParameter(G4int i) const
Definition: G4UIcommand.hh:140
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:362
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:147
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
std::size_t GetGuidanceEntries() const
Definition: G4UIcommand.hh:129
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:495
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
static G4bool IsUnitDefined(const G4String &)
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
G4Scene * GetScene() const
const G4String & GetName() const
const G4ViewParameters & GetViewParameters() const
void SetViewParameters(const G4ViewParameters &vp)
Definition: G4VViewer.cc:126
void RefreshView()
virtual void ShowView()
Definition: G4VViewer.cc:105
G4VSceneHandler * GetSceneHandler() const
static G4double fCurrentTextSize
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
static void SetVisManager(G4VisManager *pVisManager)
static G4ViewParameters fExistingVP
static G4Colour fCurrentTextColour
void CopyCameraParameters(G4ViewParameters &target, const G4ViewParameters &from)
static G4VisManager * GetVisManager()
void CheckSceneAndNotifyHandlers(G4Scene *=nullptr)
static std::vector< G4PhysicalVolumesSearchScene::Findings > fCurrrentPVFindingsForField
G4bool CheckView()
void ConvertToColour(G4Colour &colour, const G4String &redOrString, G4double green, G4double blue, G4double opacity)
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="")
static const G4Colour & GetCurrentTextColour()
static G4bool ConvertToDoublePair(const G4String &paramString, G4double &xval, G4double &yval)
void RefreshIfRequired(G4VViewer *viewer)
void SetViewParameters(G4VViewer *viewer, const G4ViewParameters &viewParams)
const G4String & ConvertToColourGuidance()
void CopyParametersFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd)
static G4int fCurrentArrow3DLineSegmentsPerCircle
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
G4bool ProvideValueOfUnit(const G4String &where, const G4String &unit, const G4String &category, G4double &value)
static G4Text::Layout fCurrentTextLayout
static G4bool fThereWasAViewer
virtual ~G4VVisCommand()
static G4double fCurrentLineWidth
void Twinkle(G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
static G4String ConvertToString(G4double x, G4double y, const char *unitName)
static G4Colour fCurrentColour
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
void InterpolateViews(G4VViewer *currentViewer, std::vector< G4ViewParameters > viewVector, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String exportString="")
void SetViewpointDirection(const G4Vector3D &viewpointDirection)
void SetScaleFactor(const G4Vector3D &scaleFactor)
static G4ViewParameters * CatmullRomCubicSplineInterpolation(const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)
const G4Vector3D & GetScaleFactor() const
void SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
const G4Vector3D & GetLightpointDirection() const
void SetFieldHalfAngle(G4double fieldHalfAngle)
const G4Vector3D & GetViewpointDirection() const
const G4Point3D & GetCurrentTargetPoint() const
G4double GetFieldHalfAngle() const
G4double GetZoomFactor() const
void SetDolly(G4double dolly)
const G4Vector3D & GetUpVector() const
void SetZoomFactor(G4double zoomFactor)
void SetUpVector(const G4Vector3D &upVector)
void SetLightpointDirection(const G4Vector3D &lightpointDirection)
void SetLightsMoveWithCamera(G4bool moves)
G4bool IsAutoRefresh() const
G4bool GetLightsMoveWithCamera() const
G4double GetDolly() const
G4double GetYmin() const
Definition: G4VisExtent.hh:101
G4double GetXmax() const
Definition: G4VisExtent.hh:100
const G4Point3D & GetExtentCenter() const
Definition: G4VisExtent.hh:106
G4double GetYmax() const
Definition: G4VisExtent.hh:102
G4double GetZmax() const
Definition: G4VisExtent.hh:104
G4double GetZmin() const
Definition: G4VisExtent.hh:103
G4double GetXmin() const
Definition: G4VisExtent.hh:99
void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())
static Verbosity GetVerbosity()
G4bool contains(const G4String &str, std::string_view ss)
Check if a string contains a given substring.