Geant4 11.2.2
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;
59
61
63
65
70
72{
73 fpVisManager = pVisManager;
74}
75
80
82(G4double x, G4double y, const char * unitName)
83{
84 G4double uv = G4UIcommand::ValueOf(unitName);
85
86 std::ostringstream oss;
87 oss << x/uv << " " << y/uv << " " << unitName;
88 return oss.str();
89}
90
92 G4double& xval,
93 G4double& yval)
94{
95 G4double x, y;
96 G4String unit;
97
98 std::istringstream is(paramString);
99 is >> x >> y >> unit;
100
102 xval = x*G4UIcommand::ValueOf(unit);
103 yval = y*G4UIcommand::ValueOf(unit);
104 } else {
106 if (verbosity >= G4VisManager::errors) {
107 G4warn << "ERROR: Unrecognised unit" << G4endl;
108 }
109 return false;
110 }
111
112 return true;
113}
114
116{
117 static G4String guidance
118 ("Accepts (a) RGB triplet. e.g., \".3 .4 .5\", or"
119 "\n (b) string such as \"white\", \"black\", \"grey\", \"red\"...or"
120 "\n (c) an additional number for opacity, e.g., \".3 .4 .5 .6\""
121 "\n or \"grey ! ! .6\" (note \"!\"'s for unused parameters).");
122 return guidance;
123}
124
126(G4Colour& colour,
127 const G4String& redOrString, G4double green, G4double blue, G4double opacity)
128{
129 // Note: colour is supplied by the caller and some or all of its components
130 // may act as default.
131 //
132 // Note: redOrString is either a number or string. If a string it must be
133 // one of the recognised colours.
134 //
135 // Thus the arguments can be, for example:
136 // (colour,"red",...,...,0.5): will give the colour red with opacity 0.5 (the
137 // third and fourth arguments are ignored), or
138 // (1.,0.,0.,0.5): this also will be red with opacity 0.5.
139
141
142 const std::size_t iPos0 = 0;
143 if (std::isalpha(redOrString[iPos0])) {
144
145 // redOrString is probably alphabetic characters defining the colour
146 if (!G4Colour::GetColour(redOrString, colour)) {
147 // Not a recognised string
148 if (verbosity >= G4VisManager::warnings) {
149 G4warn << "WARNING: Colour \"" << redOrString
150 << "\" not found. Defaulting to " << colour
151 << G4endl;
152 }
153 return;
154 } else {
155 // It was a recognised string. Now add opacity.
156 colour.SetAlpha(opacity);
157 return;
158 }
159
160 } else {
161
162 // redOrString is probably numeric defining the red component
163 std::istringstream iss(redOrString);
164 G4double red;
165 iss >> red;
166 if (iss.fail()) {
167 if (verbosity >= G4VisManager::warnings) {
168 G4warn << "WARNING: String \"" << redOrString
169 << "\" cannot be parsed. Defaulting to " << colour
170 << G4endl;
171 }
172 return;
173 } else {
174 colour = G4Colour(red,green,blue,opacity);
175 return;
176 }
177
178 }
179}
180
182(const G4String& where,
183 const G4String& unit,
184 const G4String& category,
185 G4double& value)
186{
187 // Return false if there's a problem
188
190
191 G4bool success = true;
193 if (verbosity >= G4VisManager::warnings) {
194 G4warn << where
195 << "\n Unit \"" << unit << "\" not defined"
196 << G4endl;
197 }
198 success = false;
199 } else if (G4UnitDefinition::GetCategory(unit) != category) {
200 if (verbosity >= G4VisManager::warnings) {
201 G4warn << where
202 << "\n Unit \"" << unit << "\" not a unit of " << category;
203 if (category == "Volumic Mass") G4warn << " (density)";
204 G4warn << G4endl;
205 }
206 success = false;
207 } else {
208 value = G4UnitDefinition::GetValueOf(unit);
209 }
210 return success;
211}
212
214{
216
217 if (!pScene) {
218 if (verbosity >= G4VisManager::warnings) {
219 G4warn << "WARNING: Scene pointer is null."
220 << G4endl;
221 }
222 return;
223 }
224
225 G4VSceneHandler* pSceneHandler = fpVisManager -> GetCurrentSceneHandler();
226 if (!pSceneHandler) {
227 if (verbosity >= G4VisManager::warnings) {
228 G4warn << "WARNING: Scene handler not found." << G4endl;
229 }
230 return;
231 }
232
233 // Scene has changed. If it is the scene of the currrent scene handler
234 // refresh viewers of all scene handlers using this scene. If not, it may be
235 // a scene that the user is building up before attaching to a scene handler,
236 // so do nothing.
237 if (pScene == pSceneHandler->GetScene()) {
238 G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/scene/notifyHandlers");
239 }
240
241}
242
244{
246 G4VViewer* viewer = fpVisManager -> GetCurrentViewer ();
247
248 if (!viewer) {
249 if (verbosity >= G4VisManager::errors) {
250 G4warn <<
251 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
252 << G4endl;
253 }
254 return false;
255 }
256
257 return true;
258}
259
261(G4VisManager::Verbosity verbosity) {
262 // Some frequently used error printing...
263 if (verbosity >= G4VisManager::warnings) {
264 G4warn <<
265 "WARNING: For some reason, possibly mentioned above, it has not been"
266 "\n possible to add to the scene."
267 << G4endl;
268 }
269}
270
272(G4VViewer* viewer, const G4ViewParameters& viewParams) {
273 viewer->SetViewParameters(viewParams);
274 RefreshIfRequired(viewer);
275}
276
279 G4VSceneHandler* sceneHandler = viewer->GetSceneHandler();
280 const G4ViewParameters& viewParams = viewer->GetViewParameters();
281 if (sceneHandler && sceneHandler->GetScene()) {
282 if (viewParams.IsAutoRefresh()) {
283 G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/refresh");
284 }
285 else {
286 if (verbosity >= G4VisManager::warnings) {
287 G4warn << "Issue /vis/viewer/refresh or flush to see effect."
288 << G4endl;
289 }
290 }
291 }
292}
293
295(G4VViewer* currentViewer,
296 std::vector<G4ViewParameters> viewVector,
297 const G4int nInterpolationPoints,
298 const G4int waitTimePerPointmilliseconds,
299 const G4String exportString)
300{
301 const G4int safety = (G4int)viewVector.size()*nInterpolationPoints;
302 G4int safetyCount = 0;
303 do {
304 G4ViewParameters* vp =
305 G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
306 if (!vp) break; // Finished.
307 currentViewer->SetViewParameters(*vp);
308 currentViewer->RefreshView();
309 if (exportString == "export" &&
310 G4StrUtil::contains(currentViewer->GetName(), "OpenGL")) {
311 G4UImanager::GetUIpointer()->ApplyCommand("/vis/ogl/export");
312 }
313 currentViewer->ShowView();
314 if (waitTimePerPointmilliseconds > 0)
315 std::this_thread::sleep_for(std::chrono::milliseconds(waitTimePerPointmilliseconds));
316 } while (safetyCount++ < safety); // Loop checking, 16.02.2016, J.Allison
317}
318
320(G4VViewer* currentViewer,
321 const G4ViewParameters& oldVP,
322 const G4ViewParameters& newVP,
323 const G4int nInterpolationPoints,
324 const G4int waitTimePerPointmilliseconds,
325 const G4String exportString)
326{
327 std::vector<G4ViewParameters> viewVector;
328 viewVector.push_back(oldVP);
329 viewVector.push_back(oldVP);
330 viewVector.push_back(newVP);
331 viewVector.push_back(newVP);
332
334 (currentViewer,
335 viewVector,
336 nInterpolationPoints,
337 waitTimePerPointmilliseconds,
338 exportString);
339}
340
342// Twinkles the touchables in paths
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 // Just 5 twinkles is reasonable to get a human's attention
383 for (G4int i = 0; i < 5; i++) {
384 viewVector.push_back(loVP);
385 viewVector.push_back(hiVP);
386 }
387 // Just 5 interpolation points for a reasonable twinkle rate
388 InterpolateViews(currentViewer,viewVector,5);
389}
390
392(const G4UIcommand* fromCmd, G4UIcommand* toCmd, G4int startLine)
393{
394 if (fromCmd && toCmd) {
395 const G4int nGuideEntries = (G4int)fromCmd->GetGuidanceEntries();
396 for (G4int i = startLine; i < nGuideEntries; ++i) {
397 const G4String& guidance = fromCmd->GetGuidanceLine(i);
398 toCmd->SetGuidance(guidance);
399 }
400 }
401}
402
404(const G4UIcommand* fromCmd, G4UIcommand* toCmd)
405{
406 if (fromCmd && toCmd) {
407 const G4int nParEntries = (G4int)fromCmd->GetParameterEntries();
408 for (G4int i = 0; i < nParEntries; ++i) {
409 G4UIparameter* parameter = new G4UIparameter(*(fromCmd->GetParameter(i)));
410 toCmd->SetParameter(parameter);
411 }
412 }
413}
414
416{
417 if (fpVisManager) {
418 const G4double halfX = (extent.GetXmax() - extent.GetXmin()) / 2.;
419 const G4double halfY = (extent.GetYmax() - extent.GetYmin()) / 2.;
420 const G4double halfZ = (extent.GetZmax() - extent.GetZmin()) / 2.;
421 if (halfX > 0. && halfY > 0. && halfZ > 0.) {
422 const G4Box box("vis_extent",halfX,halfY,halfZ);
423 const G4VisAttributes visAtts(G4Color::Red());
424 const G4Point3D& centre = extent.GetExtentCenter();
425 fpVisManager->Draw(box,visAtts,G4Translate3D(centre));
426 }
427 }
428}
429
431(G4ViewParameters& target, const G4ViewParameters& from)
432{
433 // Copy view parameters pertaining only to camera
437 target.SetUpVector (from.GetUpVector());
438 target.SetFieldHalfAngle (from.GetFieldHalfAngle());
439 target.SetZoomFactor (from.GetZoomFactor());
440 target.SetScaleFactor (from.GetScaleFactor());
442 target.SetDolly (from.GetDolly());
443}
#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
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 > &)
@ 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:128
void RefreshView()
virtual void ShowView()
Definition G4VViewer.cc:106
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
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)
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 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
G4double GetXmax() const
const G4Point3D & GetExtentCenter() const
G4double GetYmax() const
G4double GetZmax() const
G4double GetZmin() const
G4double GetXmin() const
void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())
static Verbosity GetVerbosity()