Geant4  10.03
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 // $Id: G4VViewer.cc 101714 2016-11-22 08:53:13Z gcosmo $
28 //
29 //
30 // John Allison 27th March 1996
31 // Abstract interface class for graphics views.
32 
33 #include "G4VViewer.hh"
34 
35 #include "G4ios.hh"
36 #include <sstream>
37 
38 #include "G4VisManager.hh"
39 #include "G4VGraphicsSystem.hh"
40 #include "G4VSceneHandler.hh"
41 #include "G4Scene.hh"
42 #include "G4VPhysicalVolume.hh"
43 #include "G4Transform3D.hh"
44 #include "G4UImanager.hh"
45 
47  G4int id, const G4String& name):
48 fSceneHandler (sceneHandler),
49 fViewId (id),
50 //fModified (true),
51 fNeedKernelVisit (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 (0, fName.find (' '));
62  fShortName.strip ();
63 
65  fDefaultVP = fVP;
66 }
67 
70 }
71 
73  fName = name;
74  fShortName = fName (0, fName.find (' '));
75  fShortName.strip ();
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  for (const auto& pvNodeId: fullPath) {
130  oss
131  << ' ' << pvNodeId.GetPhysicalVolume()->GetName()
132  << ' ' << pvNodeId.GetCopyNo();
133  }
134  G4UImanager::GetUIpointer()->ApplyCommand("/vis/set/touchable" + oss.str());
135 }
136 
138 (const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& fullPath,
139  G4bool visibiity)
140 {
141  // Changes the Vis Attribute Modifiers WITHOUT triggering a rebuild.
142 
143  std::ostringstream oss;
144  oss << "/vis/touchable/set/visibility ";
145  if (visibiity) oss << "true"; else oss << "false";
146 
147  // The following is equivalent to
148  // G4UImanager::GetUIpointer()->ApplyCommand(oss.str());
149  // (assuming the touchable has already been set), but avoids view rebuild.
150 
151  // Instantiate a working copy of a G4VisAttributes object...
152  G4VisAttributes workingVisAtts;
153  // and set the visibility.
154  workingVisAtts.SetVisibility(visibiity);
155 
156  fVP.AddVisAttributesModifier
158  (workingVisAtts,
160  fullPath));
161  // G4ModelingParameters::VASVisibility (VAS = Vis Attribute Signifier)
162  // signifies that it is the visibility that should be picked out
163  // and merged with the touchable's normal vis attributes.
164 
165  // Record on G4cout (with #) for information.
166  if (G4UImanager::GetUIpointer()->GetVerboseLevel() >= 2) {
167  G4cout << "# " << oss.str() << G4endl;
168  }
169 }
170 
172 (const std::vector<G4PhysicalVolumeModel::G4PhysicalVolumeNodeID>& fullPath,
173  const G4Colour& colour)
174 {
175  // Changes the Vis Attribute Modifiers WITHOUT triggering a rebuild.
176 
177  std::ostringstream oss;
178  oss << "/vis/touchable/set/colour "
179  << colour.GetRed() << ' ' << colour.GetGreen()
180  << ' ' << colour.GetBlue() << ' ' << colour.GetAlpha();
181 
182  // The following is equivalent to
183  // G4UImanager::GetUIpointer()->ApplyCommand(oss.str());
184  // (assuming the touchable has already been set), but avoids view rebuild.
185 
186  // Instantiate a working copy of a G4VisAttributes object...
187  G4VisAttributes workingVisAtts;
188  // and set the colour.
189  workingVisAtts.SetColour(colour);
190 
191  fVP.AddVisAttributesModifier
193  (workingVisAtts,
195  fullPath));
196  // G4ModelingParameters::VASColour (VAS = Vis Attribute Signifier)
197  // signifies that it is the colour that should be picked out
198  // and merged with the touchable's normal vis attributes.
199 
200  // Record on G4cout (with #) for information.
201  if (G4UImanager::GetUIpointer()->GetVerboseLevel() >= 2) {
202  G4cout << "# " << oss.str() << G4endl;
203  }
204 }
205 
206 std::vector <G4ThreeVector> G4VViewer::ComputeFlyThrough(G4Vector3D* /*aVect*/)
207 {
208  enum CurveType {
209  Bezier,
210  G4SplineTest};
211 
212  // Choose a curve type (for testing)
213 // int myCurveType = Bezier;
214 
215  // number if step points
216  int stepPoints = 500;
217 
218 
219  G4Spline spline;
220 
221 
222  // At the moment we don't use the aVect parameters, but build it here :
223  // Good step points for exampleB5
224  spline.AddSplinePoint(G4Vector3D(0,1000,-14000));
225  spline.AddSplinePoint(G4Vector3D(0,1000,0));
226  spline.AddSplinePoint(G4Vector3D(-4000,1000,4000));
227 
228 
229  std::vector <G4ThreeVector> viewVect;
230 
231 // if(myCurveType == Bezier) {
232 
233 
234  // Draw the spline
235 
236  for (int i = 0; i < stepPoints; i++) {
237  float t = (float)i / (float)stepPoints;
238  G4Vector3D cameraPosition = spline.GetInterpolatedSplinePoint(t);
239  // G4Vector3D targetPoint = spline.GetInterpolatedSplinePoint(t);
240 
241  // viewParam->SetViewAndLights(G4ThreeVector (cameraPosition.x(), cameraPosition.y(), cameraPosition.z()));
242  // viewParam->SetCurrentTargetPoint(targetPoint);
243  G4cout << "FLY CR("<< i << "):" << cameraPosition << G4endl;
244  viewVect.push_back(G4ThreeVector (cameraPosition.x(), cameraPosition.y(), cameraPosition.z()));
245  }
246 
247 // } else if (myCurveType == G4SplineTest) {
248  /*
249  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.
250  This method add "stay straight" parameter which could be between 0 and 1 where the pass will follow exactly the line between the waypoints
251  Ex : stay straight = 50%
252  m1 = 3*(P1+P0)/2
253 
254  Ex : stay straight = 0%
255  m1 = (P1+P0)/2
256 
257  P1
258  / \
259  / \
260  a--x--b
261  / ° ° \
262  / ° ° \
263  m1 m2
264  / \
265  / \
266  / \
267  / \
268  P0 P2
269 
270  */
271 // G4Vector3D a;
272 // G4Vector3D b;
273 // G4Vector3D m1;
274 // G4Vector3D m2;
275 // G4Vector3D P0;
276 // G4Vector3D P1;
277 // G4Vector3D P2;
278 // G4double stayStraight = 0;
279 // G4double bezierSpeed = 0.4; // Spend 40% time in bezier curve (time between m1-m2 is 40% of time between P0-P1)
280 //
281 // G4Vector3D firstPoint;
282 // G4Vector3D lastPoint;
283 //
284 // float nbBezierSteps = (stepPoints * bezierSpeed*(1-stayStraight)) * (2./spline.GetNumPoints());
285 // float nbFirstSteps = ((stepPoints/2-nbBezierSteps/2) /(1+stayStraight)) * (2./spline.GetNumPoints());
286 //
287 // // First points
288 // firstPoint = spline.GetPoint(0);
289 // lastPoint = (firstPoint + spline.GetPoint(1))/2;
290 //
291 // for( float j=0; j<1; j+= 1/nbFirstSteps) {
292 // G4ThreeVector pt = firstPoint + (lastPoint - firstPoint) * j;
293 // viewVect.push_back(pt);
294 // G4cout << "FLY Bezier A1("<< viewVect.size()<< "):" << pt << G4endl;
295 // }
296 //
297 // for (int i = 0; i < spline.GetNumPoints()-2; i++) {
298 // P0 = spline.GetPoint(i);
299 // P1 = spline.GetPoint(i+1);
300 // P2 = spline.GetPoint(i+2);
301 //
302 // m1 = P1 - (P1-P0)*(1-stayStraight)/2;
303 // m2 = P1 + (P2-P1)*(1-stayStraight)/2;
304 //
305 // // We have to get straight path from (middile of P0-P1) to (middile of P0-P1 + (dist P0-P1) * stayStraight/2)
306 // if (stayStraight >0) {
307 //
308 // firstPoint = (P0 + P1)/2;
309 // lastPoint = (P0 + P1)/2 + (P1-P0)*stayStraight/2;
310 //
311 // for( float j=0; j<1; j+= 1/(nbFirstSteps*stayStraight)) {
312 // G4ThreeVector pt = firstPoint + (lastPoint - firstPoint)* j;
313 // viewVect.push_back(pt);
314 // G4cout << "FLY Bezier A2("<< viewVect.size()<< "):" << pt << G4endl;
315 // }
316 // }
317 // // Compute Bezier curve
318 // for( float delta = 0 ; delta < 1 ; delta += 1/nbBezierSteps)
319 // {
320 // // The Green Line
321 // a = m1 + ( (P1 - m1) * delta );
322 // b = P1 + ( (m2 - P1) * delta );
323 //
324 // // Final point
325 // G4ThreeVector pt = a + ((b-a) * delta );
326 // viewVect.push_back(pt);
327 // G4cout << "FLY Bezier("<< viewVect.size()<< "):" << pt << G4endl;
328 // }
329 //
330 // // We have to get straight path
331 // if (stayStraight >0) {
332 // firstPoint = (P1 + P2)/2 - (P2-P1)*stayStraight/2;
333 // lastPoint = (P1 + P2)/2;
334 //
335 // for( float j=0; j<1; j+= 1/(nbFirstSteps*stayStraight)) {
336 // G4ThreeVector pt = firstPoint + (lastPoint - firstPoint)* j;
337 // viewVect.push_back(pt);
338 // G4cout << "FLY Bezier B1("<< viewVect.size()<< "):" << pt << G4endl;
339 // }
340 // }
341 // }
342 //
343 // // last points
344 // firstPoint = spline.GetPoint(spline.GetNumPoints()-2);
345 // lastPoint = spline.GetPoint(spline.GetNumPoints()-1);
346 // for( float j=1; j>0; j-= 1/nbFirstSteps) {
347 // G4ThreeVector pt = lastPoint - ((lastPoint-firstPoint)*((1-stayStraight)/2) * j );
348 // viewVect.push_back(pt);
349 // G4cout << "FLY Bezier B2("<< viewVect.size()<< "):" << pt << G4endl;
350 // }
351 // }
352  return viewVect;
353 }
354 
355 
356 #ifdef G4MULTITHREADED
357 
358 void G4VViewer::DoneWithMasterThread () {
359  // G4cout << "G4VViewer::DoneWithMasterThread" << G4endl;
360 }
361 
362 void G4VViewer::MovingToMasterThread () {
363  // G4cout << "G4VViewer::MovingToMasterThread" << G4endl;
364 }
365 
366 void G4VViewer::SwitchToVisSubThread () {
367  // G4cout << "G4VViewer::SwitchToVisSubThread" << G4endl;
368 }
369 
370 void G4VViewer::DoneWithVisSubThread () {
371  // G4cout << "G4VViewer::DoneWithVisSubThread" << G4endl;
372 }
373 
374 void G4VViewer::MovingToVisSubThread () {
375  // G4cout << "G4VViewer::MovingToVisSubThread" << G4endl;
376 }
377 
378 void G4VViewer::SwitchToMasterThread () {
379  // G4cout << "G4VViewer::SwitchToMasterThread" << G4endl;
380 }
381 
382 #endif
383 
384 std::ostream& operator << (std::ostream& os, const G4VViewer& v) {
385  os << "View " << v.fName << ":\n";
386  os << v.fVP;
387  return os;
388 }
389 
390 
391 // ===== G4Spline class =====
392 
394 : vp(), delta_t(0)
395 {
396 }
397 
398 
400 {}
401 
402 // Solve the Catmull-Rom parametric equation for a given time(t) and vector quadruple (p1,p2,p3,p4)
403 G4Vector3D G4VViewer::G4Spline::CatmullRom_Eq(float t, const G4Vector3D& p1, const G4Vector3D& p2, const G4Vector3D& p3, const G4Vector3D& p4)
404 {
405  float t2 = t * t;
406  float t3 = t2 * t;
407 
408  float b1 = .5 * ( -t3 + 2*t2 - t);
409  float b2 = .5 * ( 3*t3 - 5*t2 + 2);
410  float b3 = .5 * (-3*t3 + 4*t2 + t);
411  float b4 = .5 * ( t3 - t2 );
412 
413  return (p1*b1 + p2*b2 + p3*b3 + p4*b4);
414 }
415 
417 {
418  vp.push_back(v);
419  delta_t = (float)1 / (float)vp.size();
420 }
421 
422 
424 {
425  return vp[a];
426 }
427 
429 {
430  return vp.size();
431 }
432 
434 {
435  // Find out in which interval we are on the spline
436  int p = (int)(t / delta_t);
437  // Compute local control point indices
438 #define BOUNDS(pp) { if (pp < 0) pp = 0; else if (pp >= (int)vp.size()-1) pp = vp.size() - 1; }
439  int p0 = p - 1; BOUNDS(p0);
440  int p1 = p; BOUNDS(p1);
441  int p2 = p + 1; BOUNDS(p2);
442  int p3 = p + 2; BOUNDS(p3);
443  // Relative (local) time
444  float lt = (t - delta_t*(float)p) / delta_t;
445  // Interpolate
446  return CatmullRom_Eq(lt, vp[p0], vp[p1], vp[p2], vp[p3]);
447 }
void SetName(const G4String &)
Definition: G4VViewer.cc:72
virtual void ClearStore()
G4int fViewId
Definition: G4VViewer.hh:202
void SetColour(const G4Colour &)
G4double GetAlpha() const
Definition: G4Colour.hh:142
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
G4String strip(G4int strip_Type=trailing, char c=' ')
virtual void ShowView()
Definition: G4VViewer.cc:103
void RemoveViewerFromList(G4VViewer *pView)
const G4String & GetName() const
void SetViewParameters(const G4ViewParameters &vp)
Definition: G4VViewer.cc:120
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
G4double GetBlue() const
Definition: G4Colour.hh:141
const char * name(G4int ptype)
int G4int
Definition: G4Types.hh:78
void TouchableSetVisibility(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, G4bool visibility)
Definition: G4VViewer.cc:138
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:206
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition: G4VViewer.cc:46
std::ostream & operator<<(std::ostream &os, const G4VViewer &v)
Definition: G4VViewer.cc:384
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
void TouchableSetColour(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, const G4Colour &)
Definition: G4VViewer.cc:172
const G4ViewParameters & GetDefaultViewParameters() const
G4GLOB_DLL std::ostream G4cout
G4double GetRed() const
Definition: G4Colour.hh:139
void SetVisibility(G4bool=true)
bool G4bool
Definition: G4Types.hh:79
G4double GetGreen() const
Definition: G4Colour.hh:140
static G4VisManager * GetInstance()
G4ViewParameters fVP
Definition: G4VViewer.hh:205
G4Vector3D GetPoint(int)
Definition: G4VViewer.cc:423
void SetTouchable(const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath)
Definition: G4VViewer.cc:125
G4bool fNeedKernelVisit
Definition: G4VViewer.hh:210
G4Vector3D GetInterpolatedSplinePoint(float t)
Definition: G4VViewer.cc:433
G4String fName
Definition: G4VViewer.hh:203
#define BOUNDS(pp)
virtual ~G4VViewer()
Definition: G4VViewer.cc:68
virtual void ProcessScene()
void NeedKernelVisit()
Definition: G4VViewer.cc:78
G4Vector3D CatmullRom_Eq(float t, const G4Vector3D &p1, const G4Vector3D &p2, const G4Vector3D &p3, const G4Vector3D &p4)
Definition: G4VViewer.cc:403
#define G4endl
Definition: G4ios.hh:61
void ProcessView()
Definition: G4VViewer.cc:105
std::vector< G4ThreeVector > ComputeFlyThrough(G4Vector3D *)
Definition: G4VViewer.cc:206
G4String fShortName
Definition: G4VViewer.hh:204
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:201
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:446
void AddSplinePoint(const G4Vector3D &v)
Definition: G4VViewer.cc:416
virtual void FinishView()
Definition: G4VViewer.cc:101