Geant4 10.7.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
45// Not yet used: G4VisAttributes::LineStyle G4VVisCommand::fCurrentLineStyle = G4VisAttributes::unbroken;
46// Not yet used: G4VMarker::FillStyle G4VVisCommand::fCurrentFillStyle = G4VMarker::filled;
47// Not yet used: G4VMarker::SizeType G4VVisCommand::fCurrentSizeType = G4VMarker::screen;
50std::vector<G4PhysicalVolumesSearchScene::Findings> G4VVisCommand::fCurrrentPVFindingsForField;
51
53
55
57
59(G4double x, G4double y, const char * unitName)
60{
61 G4double uv = G4UIcommand::ValueOf(unitName);
62
63 std::ostringstream oss;
64 oss << x/uv << " " << y/uv << " " << unitName;
65 return oss.str();
66}
67
69 G4double& xval,
70 G4double& yval)
71{
72 G4double x, y;
73 G4String unit;
74
75 std::istringstream is(paramString);
76 is >> x >> y >> unit;
77
79 xval = x*G4UIcommand::ValueOf(unit);
80 yval = y*G4UIcommand::ValueOf(unit);
81 } else {
83 if (verbosity >= G4VisManager::errors) {
84 G4cout << "ERROR: Unrecognised unit" << G4endl;
85 }
86 return false;
87 }
88
89 return true;
90}
91
93{
94 static G4String guidance
95 ("Accepts (a) RGB triplet. e.g., \".3 .4 .5\", or"
96 "\n (b) string such as \"white\", \"black\", \"grey\", \"red\"...or"
97 "\n (c) an additional number for opacity, e.g., \".3 .4 .5 .6\""
98 "\n or \"grey ! ! .6\" (note \"!\"'s for unused parameters).");
99 return guidance;
100}
101
103(G4Colour& colour,
104 const G4String& redOrString, G4double green, G4double blue, G4double opacity)
105{
106 // Note: colour is supplied by the caller and some or all of its components
107 // may act as default.
108 //
109 // Note: redOrString is either a number or string. If a string it must be
110 // one of the recognised colours.
111 //
112 // Thus the arguments can be, for example:
113 // (colour,"red",...,...,0.5): will give the colour red with opacity 0.5 (the
114 // third and fourth arguments are ignored), or
115 // (1.,0.,0.,0.5): this also will be red with opacity 0.5.
116
118
119 const size_t iPos0 = 0;
120 if (std::isalpha(redOrString[iPos0])) {
121
122 // redOrString is probably alphabetic characters defining the colour
123 if (!G4Colour::GetColour(redOrString, colour)) {
124 // Not a recognised string
125 if (verbosity >= G4VisManager::warnings) {
126 G4cout << "WARNING: Colour \"" << redOrString
127 << "\" not found. Defaulting to " << colour
128 << G4endl;
129 }
130 return;
131 } else {
132 // It was a recognised string. Now add opacity.
133 colour.SetAlpha(opacity);
134 return;
135 }
136
137 } else {
138
139 // redOrString is probably numeric defining the red component
140 std::istringstream iss(redOrString);
141 G4double red;
142 iss >> red;
143 if (iss.fail()) {
144 if (verbosity >= G4VisManager::warnings) {
145 G4cout << "WARNING: String \"" << redOrString
146 << "\" cannot be parsed. Defaulting to " << colour
147 << G4endl;
148 }
149 return;
150 } else {
151 colour = G4Colour(red,green,blue,opacity);
152 return;
153 }
154
155 }
156}
157
159(const G4String& where,
160 const G4String& unit,
161 const G4String& category,
162 G4double& value)
163{
164 // Return false if there's a problem
165
167
168 G4bool success = true;
170 if (verbosity >= G4VisManager::warnings) {
171 G4cerr << where
172 << "\n Unit \"" << unit << "\" not defined"
173 << G4endl;
174 }
175 success = false;
176 } else if (G4UnitDefinition::GetCategory(unit) != category) {
177 if (verbosity >= G4VisManager::warnings) {
178 G4cerr << where
179 << "\n Unit \"" << unit << "\" not a unit of " << category;
180 if (category == "Volumic Mass") G4cerr << " (density)";
181 G4cerr << G4endl;
182 }
183 success = false;
184 } else {
185 value = G4UnitDefinition::GetValueOf(unit);
186 }
187 return success;
188}
189
191{
193
194 if (!pScene) {
195 if (verbosity >= G4VisManager::warnings) {
196 G4cout << "WARNING: Scene pointer is null."
197 << G4endl;
198 }
199 return;
200 }
201
202 G4VSceneHandler* pSceneHandler = fpVisManager -> GetCurrentSceneHandler();
203 if (!pSceneHandler) {
204 if (verbosity >= G4VisManager::warnings) {
205 G4cout << "WARNING: Scene handler not found." << G4endl;
206 }
207 return;
208 }
209
210 // Scene has changed. If it is the scene of the currrent scene handler
211 // refresh viewers of all scene handlers using this scene. If not, it may be
212 // a scene that the user is building up before attaching to a scene handler,
213 // so do nothing.
214 if (pScene == pSceneHandler->GetScene()) {
215 G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/scene/notifyHandlers");
216 }
217
218}
219
221{
223 G4VViewer* viewer = fpVisManager -> GetCurrentViewer ();
224
225 if (!viewer) {
226 if (verbosity >= G4VisManager::errors) {
227 G4cerr <<
228 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
229 << G4endl;
230 }
231 return false;
232 }
233
234 return true;
235}
236
238(G4VisManager::Verbosity verbosity) {
239 // Some frequently used error printing...
240 if (verbosity >= G4VisManager::warnings) {
241 G4cout <<
242 "WARNING: For some reason, possibly mentioned above, it has not been"
243 "\n possible to add to the scene."
244 << G4endl;
245 }
246}
247
249(G4VViewer* viewer, const G4ViewParameters& viewParams) {
250 viewer->SetViewParameters(viewParams);
251 RefreshIfRequired(viewer);
252}
253
256 G4VSceneHandler* sceneHandler = viewer->GetSceneHandler();
257 const G4ViewParameters& viewParams = viewer->GetViewParameters();
258 if (sceneHandler && sceneHandler->GetScene()) {
259 if (viewParams.IsAutoRefresh()) {
260 G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/refresh");
261 }
262 else {
263 if (verbosity >= G4VisManager::warnings) {
264 G4cout << "Issue /vis/viewer/refresh or flush to see effect."
265 << G4endl;
266 }
267 }
268 }
269}
270
272(G4VViewer* currentViewer,
273 std::vector<G4ViewParameters> viewVector,
274 const G4int nInterpolationPoints,
275 const G4int waitTimePerPointmilliseconds,
276 const G4String exportString)
277{
278 const G4int safety = viewVector.size()*nInterpolationPoints;
279 G4int safetyCount = 0;
280 do {
281 G4ViewParameters* vp =
282 G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
283 if (!vp) break; // Finished.
284 currentViewer->SetViewParameters(*vp);
285 currentViewer->RefreshView();
286 if (exportString == "export" &&
287 currentViewer->GetName().contains("OpenGL")) {
288 G4UImanager::GetUIpointer()->ApplyCommand("/vis/ogl/export");
289 }
290 // File-writing viewers need to close the file
291 currentViewer->ShowView();
292 if (waitTimePerPointmilliseconds > 0)
293 std::this_thread::sleep_for(std::chrono::milliseconds(waitTimePerPointmilliseconds));
294 } while (safetyCount++ < safety); // Loop checking, 16.02.2016, J.Allison
295}
296
298(G4VViewer* currentViewer,
299 const G4ViewParameters& oldVP,
300 const G4ViewParameters& newVP,
301 const G4int nInterpolationPoints,
302 const G4int waitTimePerPointmilliseconds,
303 const G4String exportString)
304{
305 std::vector<G4ViewParameters> viewVector;
306 viewVector.push_back(oldVP);
307 viewVector.push_back(oldVP);
308 viewVector.push_back(newVP);
309 viewVector.push_back(newVP);
310
312 (currentViewer,
313 viewVector,
314 nInterpolationPoints,
315 waitTimePerPointmilliseconds,
316 exportString);
317}
318
320(const G4UIcommand* fromCmd, G4UIcommand* toCmd, G4int startLine)
321{
322 if (fromCmd && toCmd) {
323 const G4int nGuideEntries = fromCmd->GetGuidanceEntries();
324 for (G4int i = startLine; i < nGuideEntries; ++i) {
325 const G4String& guidance = fromCmd->GetGuidanceLine(i);
326 toCmd->SetGuidance(guidance);
327 }
328 }
329}
330
332(const G4UIcommand* fromCmd, G4UIcommand* toCmd)
333{
334 if (fromCmd && toCmd) {
335 const G4int nParEntries = fromCmd->GetParameterEntries();
336 for (G4int i = 0; i < nParEntries; ++i) {
337 G4UIparameter* parameter = new G4UIparameter(*(fromCmd->GetParameter(i)));
338 toCmd->SetParameter(parameter);
339 }
340 }
341}
342
344{
345 if (fpVisManager) {
346 const G4double halfX = (extent.GetXmax() - extent.GetXmin()) / 2.;
347 const G4double halfY = (extent.GetYmax() - extent.GetYmin()) / 2.;
348 const G4double halfZ = (extent.GetZmax() - extent.GetZmin()) / 2.;
349 if (halfX > 0. && halfY > 0. && halfZ > 0.) {
350 const G4Box box("vis_extent",halfX,halfY,halfZ);
351 const G4VisAttributes visAtts(G4Color::Red());
352 const G4Point3D& centre = extent.GetExtentCenter();
353 fpVisManager->Draw(box,visAtts,G4Translate3D(centre));
354 }
355 }
356}
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Box.hh:56
static G4Colour White()
Definition: G4Colour.hh:154
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition: G4Colour.cc:162
void SetAlpha(G4double)
Definition: G4Colour.cc:70
static G4Colour Red()
Definition: G4Colour.hh:159
static G4Colour Blue()
Definition: G4Colour.hh:161
G4bool contains(const std::string &) const
Layout
Definition: G4Text.hh:76
@ left
Definition: G4Text.hh:76
std::size_t GetParameterEntries() const
Definition: G4UIcommand.hh:138
const G4String & GetGuidanceLine(G4int i) const
Definition: G4UIcommand.hh:132
G4UIparameter * GetParameter(G4int i) const
Definition: G4UIcommand.hh:139
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:348
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
std::size_t GetGuidanceEntries() const
Definition: G4UIcommand.hh:128
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
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:119
void RefreshView()
virtual void ShowView()
Definition: G4VViewer.cc:102
G4VSceneHandler * GetSceneHandler() const
static G4double fCurrentTextSize
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
static G4Colour fCurrentTextColour
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 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
virtual ~G4VVisCommand()
static G4double fCurrentLineWidth
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="")
static G4ViewParameters * CatmullRomCubicSplineInterpolation(const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)
G4bool IsAutoRefresh() 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()