Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VViewer.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//
29// John Allison 27th March 1996
30// Abstract interface class for graphics views.
31
32#include "G4VViewer.hh"
33
34#include "G4Timer.hh"
35
36#include "G4ios.hh"
37#include <sstream>
38
39#include "G4VisManager.hh"
40#include "G4VGraphicsSystem.hh"
41#include "G4VSceneHandler.hh"
42#include "G4Scene.hh"
44#include "G4VPhysicalVolume.hh"
45#include "G4Transform3D.hh"
46#include "G4UImanager.hh"
47
49 G4int id, const G4String& name):
50fSceneHandler (sceneHandler),
51fViewId (id),
52//fModified (true),
53fNeedKernelVisit (true)
54{
55 if (name == "") {
56 std::ostringstream ost;
57 ost << fSceneHandler.GetName () << '-' << fViewId;
58 fName = ost.str();
59 }
60 else {
61 fName = name;
62 }
63 fShortName = fName.substr(0, fName.find (' '));
65
68}
69
72}
73
74void G4VViewer::SetName (const G4String& name) {
75 fName = name;
76 fShortName = fName.substr(0, fName.find (' '));
78}
79
81
82 fNeedKernelVisit = true;
83
84 // At one time I thought we'd better notify all viewers. But I guess
85 // each viewer can take care of itself, so the following code is
86 // redundant (but keep it commented out for now). (John Allison)
87 // Notify all viewers that a kernel visit is required.
88 // const G4ViewerList& viewerList = fSceneHandler.GetViewerList ();
89 // G4ViewerListConstIterator i;
90 // for (i = viewerList.begin(); i != viewerList.end(); i++) {
91 // (*i) -> SetNeedKernelVisit ();
92 // }
93 // ??...but, there's a problem in OpenGL Stored which seems to
94 // require *all* viewers to revisit the kernel, so...
95 // const G4ViewerList& viewerList = fSceneHandler.GetViewerList ();
96 // G4ViewerListConstIterator i;
97 // for (i = viewerList.begin(); i != viewerList.end(); i++) {
98 // (*i) -> SetNeedKernelVisit (true);
99 // }
100 // Feb 2005 - commented out. Let's fix OpenGL if necessary.
101}
102
104
106
108{
109 // If the scene has changed, or if the concrete viewer has decided
110 // that it necessary to visit the kernel, perhaps because the view
111 // parameters have changed significantly (this should be done in the
112 // concrete viewer's DrawView)...
113 if (fNeedKernelVisit) {
114 // Reset flag. This must be done before ProcessScene to prevent
115 // recursive calls when recomputing transients...
116 G4Timer timer;
117 timer.Start();
118 fNeedKernelVisit = false;
121 timer.Stop();
123 }
124}
125
127 fVP = vp;
128}
129
131(const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& fullPath)
132{
133 // Set the touchable for /vis/touchable/set/... commands.
134 std::ostringstream oss;
135 const auto& pvStore = G4PhysicalVolumeStore::GetInstance();
136 for (const auto& pvNodeId: fullPath) {
137 const auto& pv = pvNodeId.GetPhysicalVolume();
138 auto iterator = find(pvStore->cbegin(),pvStore->cend(),pv);
139 if (iterator == pvStore->cend()) {
141 ed << "Volume no longer in physical volume store.";
142 G4Exception("G4VViewer::SetTouchable", "visman0501", JustWarning, ed);
143 } else {
144 oss
145 << ' ' << pvNodeId.GetPhysicalVolume()->GetName()
146 << ' ' << pvNodeId.GetCopyNo();
147 }
148 }
149 G4UImanager::GetUIpointer()->ApplyCommand("/vis/set/touchable" + oss.str());
150}
151
153(const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& fullPath,
154 G4bool visibiity)
155{
156 // Changes the Vis Attribute Modifiers WITHOUT triggering a rebuild.
157
158 // The following is equivalent to
159 // G4UImanager::GetUIpointer()->ApplyCommand("/vis/touchable/set/visibility ...");
160 // (assuming the touchable has already been set), but avoids view rebuild.
161
162 // Instantiate a working copy of a G4VisAttributes object...
163 G4VisAttributes workingVisAtts;
164 // and set the visibility.
165 workingVisAtts.SetVisibility(visibiity);
166
169 (workingVisAtts,
172 // G4ModelingParameters::VASVisibility (VAS = Vis Attribute Signifier)
173 // signifies that it is the visibility that should be picked out
174 // and merged with the touchable's normal vis attributes.
175}
176
178(const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& fullPath,
179 const G4Colour& colour)
180{
181 // Changes the Vis Attribute Modifiers WITHOUT triggering a rebuild.
182
183 // The following is equivalent to
184 // G4UImanager::GetUIpointer()->ApplyCommand("/vis/touchable/set/colour ...");
185 // (assuming the touchable has already been set), but avoids view rebuild.
186
187 // Instantiate a working copy of a G4VisAttributes object...
188 G4VisAttributes workingVisAtts;
189 // and set the colour.
190 workingVisAtts.SetColour(colour);
191
194 (workingVisAtts,
197 // G4ModelingParameters::VASColour (VAS = Vis Attribute Signifier)
198 // signifies that it is the colour that should be picked out
199 // and merged with the touchable's normal vis attributes.
200}
201
202std::vector <G4ThreeVector> G4VViewer::ComputeFlyThrough(G4Vector3D* /*aVect*/)
203{
204 enum CurveType {
205 Bezier,
206 G4SplineTest};
207
208 // Choose a curve type (for testing)
209// int myCurveType = Bezier;
210
211 // number if step points
212 G4int stepPoints = 500;
213
214
215 G4Spline spline;
216
217
218 // At the moment we don't use the aVect parameters, but build it here :
219 // Good step points for exampleB5
220 spline.AddSplinePoint(G4Vector3D(0,1000,-14000));
221 spline.AddSplinePoint(G4Vector3D(0,1000,0));
222 spline.AddSplinePoint(G4Vector3D(-4000,1000,4000));
223
224
225 std::vector <G4ThreeVector> viewVect;
226
227// if(myCurveType == Bezier) {
228
229
230 // Draw the spline
231
232 for (G4int i = 0; i < stepPoints; ++i) {
233 G4float t = (G4float)i / (G4float)stepPoints;
234 G4Vector3D cameraPosition = spline.GetInterpolatedSplinePoint(t);
235 // G4Vector3D targetPoint = spline.GetInterpolatedSplinePoint(t);
236
237 // viewParam->SetViewAndLights(G4ThreeVector (cameraPosition.x(), cameraPosition.y(), cameraPosition.z()));
238 // viewParam->SetCurrentTargetPoint(targetPoint);
239 G4cout << "FLY CR("<< i << "):" << cameraPosition << G4endl;
240 viewVect.push_back(G4ThreeVector (cameraPosition.x(), cameraPosition.y(), cameraPosition.z()));
241 }
242
243// } else if (myCurveType == G4SplineTest) {
244 /*
245 This method is a inspire from a Bezier curve. The problem of the Bezier curve is that the path does not go straight between two waypoints.
246 This method add "stay straight" parameter which could be between 0 and 1 where the pass will follow exactly the line between the waypoints
247 Ex : stay straight = 50%
248 m1 = 3*(P1+P0)/2
249
250 Ex : stay straight = 0%
251 m1 = (P1+P0)/2
252
253 P1
254 / \
255 / \
256 a--x--b
257 / ° ° \
258 / ° ° \
259 m1 m2
260 / \
261 / \
262 / \
263 / \
264 P0 P2
265
266 */
267// G4Vector3D a;
268// G4Vector3D b;
269// G4Vector3D m1;
270// G4Vector3D m2;
271// G4Vector3D P0;
272// G4Vector3D P1;
273// G4Vector3D P2;
274// G4double stayStraight = 0;
275// G4double bezierSpeed = 0.4; // Spend 40% time in bezier curve (time between m1-m2 is 40% of time between P0-P1)
276//
277// G4Vector3D firstPoint;
278// G4Vector3D lastPoint;
279//
280// float nbBezierSteps = (stepPoints * bezierSpeed*(1-stayStraight)) * (2./spline.GetNumPoints());
281// float nbFirstSteps = ((stepPoints/2-nbBezierSteps/2) /(1+stayStraight)) * (2./spline.GetNumPoints());
282//
283// // First points
284// firstPoint = spline.GetPoint(0);
285// lastPoint = (firstPoint + spline.GetPoint(1))/2;
286//
287// for( float j=0; j<1; j+= 1/nbFirstSteps) {
288// G4ThreeVector pt = firstPoint + (lastPoint - firstPoint) * j;
289// viewVect.push_back(pt);
290// G4cout << "FLY Bezier A1("<< viewVect.size()<< "):" << pt << G4endl;
291// }
292//
293// for (int i = 0; i < spline.GetNumPoints()-2; i++) {
294// P0 = spline.GetPoint(i);
295// P1 = spline.GetPoint(i+1);
296// P2 = spline.GetPoint(i+2);
297//
298// m1 = P1 - (P1-P0)*(1-stayStraight)/2;
299// m2 = P1 + (P2-P1)*(1-stayStraight)/2;
300//
301// // We have to get straight path from (middile of P0-P1) to (middile of P0-P1 + (dist P0-P1) * stayStraight/2)
302// if (stayStraight >0) {
303//
304// firstPoint = (P0 + P1)/2;
305// lastPoint = (P0 + P1)/2 + (P1-P0)*stayStraight/2;
306//
307// for( float j=0; j<1; j+= 1/(nbFirstSteps*stayStraight)) {
308// G4ThreeVector pt = firstPoint + (lastPoint - firstPoint)* j;
309// viewVect.push_back(pt);
310// G4cout << "FLY Bezier A2("<< viewVect.size()<< "):" << pt << G4endl;
311// }
312// }
313// // Compute Bezier curve
314// for( float delta = 0 ; delta < 1 ; delta += 1/nbBezierSteps)
315// {
316// // The Green Line
317// a = m1 + ( (P1 - m1) * delta );
318// b = P1 + ( (m2 - P1) * delta );
319//
320// // Final point
321// G4ThreeVector pt = a + ((b-a) * delta );
322// viewVect.push_back(pt);
323// G4cout << "FLY Bezier("<< viewVect.size()<< "):" << pt << G4endl;
324// }
325//
326// // We have to get straight path
327// if (stayStraight >0) {
328// firstPoint = (P1 + P2)/2 - (P2-P1)*stayStraight/2;
329// lastPoint = (P1 + P2)/2;
330//
331// for( float j=0; j<1; j+= 1/(nbFirstSteps*stayStraight)) {
332// G4ThreeVector pt = firstPoint + (lastPoint - firstPoint)* j;
333// viewVect.push_back(pt);
334// G4cout << "FLY Bezier B1("<< viewVect.size()<< "):" << pt << G4endl;
335// }
336// }
337// }
338//
339// // last points
340// firstPoint = spline.GetPoint(spline.GetNumPoints()-2);
341// lastPoint = spline.GetPoint(spline.GetNumPoints()-1);
342// for( float j=1; j>0; j-= 1/nbFirstSteps) {
343// G4ThreeVector pt = lastPoint - ((lastPoint-firstPoint)*((1-stayStraight)/2) * j );
344// viewVect.push_back(pt);
345// G4cout << "FLY Bezier B2("<< viewVect.size()<< "):" << pt << G4endl;
346// }
347// }
348 return viewVect;
349}
350
351
352#ifdef G4MULTITHREADED
353
354void G4VViewer::DoneWithMasterThread () {
355 // G4cout << "G4VViewer::DoneWithMasterThread" << G4endl;
356}
357
358void G4VViewer::MovingToMasterThread () {
359 // G4cout << "G4VViewer::MovingToMasterThread" << G4endl;
360}
361
362void G4VViewer::SwitchToVisSubThread () {
363 // G4cout << "G4VViewer::SwitchToVisSubThread" << G4endl;
364}
365
366void G4VViewer::DoneWithVisSubThread () {
367 // G4cout << "G4VViewer::DoneWithVisSubThread" << G4endl;
368}
369
370void G4VViewer::MovingToVisSubThread () {
371 // G4cout << "G4VViewer::MovingToVisSubThread" << G4endl;
372}
373
374void G4VViewer::SwitchToMasterThread () {
375 // G4cout << "G4VViewer::SwitchToMasterThread" << G4endl;
376}
377
378#endif
379
380std::ostream& operator << (std::ostream& os, const G4VViewer& v) {
381 os << "View " << v.fName << ":\n";
382 os << v.fVP;
383 return os;
384}
385
386
387// ===== G4Spline class =====
388
390: vp(), delta_t(0)
391{
392}
393
394
396{}
397
398// Solve the Catmull-Rom parametric equation for a given time(t) and vector quadruple (p1,p2,p3,p4)
400{
401 G4float t2 = t * t;
402 G4float t3 = t2 * t;
403
404 G4float b1 = .5 * ( -t3 + 2*t2 - t);
405 G4float b2 = .5 * ( 3*t3 - 5*t2 + 2);
406 G4float b3 = .5 * (-3*t3 + 4*t2 + t);
407 G4float b4 = .5 * ( t3 - t2 );
408
409 return (p1*b1 + p2*b2 + p3*b3 + p4*b4);
410}
411
413{
414 vp.push_back(v);
415 delta_t = (G4float)1 / (G4float)vp.size();
416}
417
418
420{
421 return vp[a];
422}
423
425{
426 return (G4int)vp.size();
427}
428
430{
431 // Find out in which interval we are on the spline
432 G4int p = (G4int)(t / delta_t);
433 // Compute local control point indices
434#define BOUNDS(pp) { if (pp < 0) pp = 0; else if (pp >= (G4int)vp.size()-1) pp = (G4int)vp.size() - 1; }
435 G4int p0 = p - 1; BOUNDS(p0);
436 G4int p1 = p; BOUNDS(p1);
437 G4int p2 = p + 1; BOUNDS(p2);
438 G4int p3 = p + 2; BOUNDS(p3);
439 // Relative (local) time
440 G4float lt = (t - delta_t*p) / delta_t;
441 // Interpolate
442 return CatmullRom_Eq(lt, vp[p0], vp[p1], vp[p2], vp[p3]);
443}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
float G4float
Definition: G4Types.hh:84
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::ostream & operator<<(std::ostream &os, const G4VViewer &v)
Definition: G4VViewer.cc:380
#define BOUNDS(pp)
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath(const std::vector< G4PhysicalVolumeNodeID > &)
static G4PhysicalVolumeStore * GetInstance()
void Stop()
void Start()
G4double GetRealElapsed() const
Definition: G4Timer.cc:113
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:495
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
virtual void ProcessScene()
const G4String & GetName() const
void RemoveViewerFromList(G4VViewer *pView)
virtual void ClearStore()
G4Vector3D CatmullRom_Eq(G4float t, const G4Vector3D &p1, const G4Vector3D &p2, const G4Vector3D &p3, const G4Vector3D &p4)
Definition: G4VViewer.cc:399
G4Vector3D GetPoint(int)
Definition: G4VViewer.cc:419
void AddSplinePoint(const G4Vector3D &v)
Definition: G4VViewer.cc:412
G4Vector3D GetInterpolatedSplinePoint(G4float t)
Definition: G4VViewer.cc:429
void SetTouchable(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath)
Definition: G4VViewer.cc:131
G4bool fNeedKernelVisit
Definition: G4VViewer.hh:227
void SetName(const G4String &)
Definition: G4VViewer.cc:74
void ProcessView()
Definition: G4VViewer.cc:107
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:216
G4String fShortName
Definition: G4VViewer.hh:219
G4double fKernelVisitElapsedTimeSeconds
Definition: G4VViewer.hh:222
virtual ~G4VViewer()
Definition: G4VViewer.cc:70
G4String fName
Definition: G4VViewer.hh:218
void NeedKernelVisit()
Definition: G4VViewer.cc:80
std::vector< G4ThreeVector > ComputeFlyThrough(G4Vector3D *)
Definition: G4VViewer.cc:202
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:221
G4int fViewId
Definition: G4VViewer.hh:217
void TouchableSetVisibility(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, G4bool visibility)
Definition: G4VViewer.cc:153
G4ViewParameters fVP
Definition: G4VViewer.hh:220
virtual void FinishView()
Definition: G4VViewer.cc:103
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition: G4VViewer.cc:48
void SetViewParameters(const G4ViewParameters &vp)
Definition: G4VViewer.cc:126
void TouchableSetColour(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, const G4Colour &)
Definition: G4VViewer.cc:178
virtual void ShowView()
Definition: G4VViewer.cc:105
void AddVisAttributesModifier(const G4ModelingParameters::VisAttributesModifier &)
void SetColour(const G4Colour &)
void SetVisibility(G4bool=true)
const G4ViewParameters & GetDefaultViewParameters() const
static G4VisManager * GetInstance()
void strip(G4String &str, char ch=' ')
Remove leading and trailing characters from string.