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