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