KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVEventViewer.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Fri Apr 12 15:48:03 2013
2 //Author: John Frankland,,,
3 
4 #include "KVEventViewer.h"
5 #include "TGeoMatrix.h"
6 #include "TGLViewer.h"
7 #include "TGLLightSet.h"
8 #include "TVirtualPad.h"
9 #include "TSystem.h"
10 #include "TRandom.h"
11 #include "TGLPerspectiveCamera.h"
12 #include <TBranch.h>
13 #include "TTree.h"
14 
16 
17 
18 
26 KVEventViewer::KVEventViewer(Int_t protoncolor, Int_t neutroncolor, Int_t highlightprotoncolor, Int_t highlightneutroncolor, Double_t freenucleonradius, Double_t nucleonradius, Double_t nuclearradiusparameter)
27  : fproton_color(protoncolor), fneutron_color(neutroncolor),
28  fProton_color(highlightprotoncolor), fNeutron_color(highlightneutroncolor),
29  free_nucleon_radius(freenucleonradius),
30  nucleon_radius(nucleonradius),
31  nuclear_radius_parameter(nuclearradiusparameter)
32 {
33  // protoncolor = color of protons
34  // neutroncolor = color of neutrons
35  // highlight...color = color of protons/neutrons in highlighted nuclei
36  // freenucleonradius = radius of free nucleons
37  // nucleonradius = radius of nucleons in nuclei
38  // nuclearradiusparameter = r0 in expression r=r0*A**1/3
39 
40  ivol = 0;
41 
42  geom = 0;
43 
44  fMaxVelocity = -1.;
45  fFixSeed = kFALSE;
46  fRefresh = 0;
47  fSeed = 863167;
48  fXYMode = kFALSE;
49  fMomentumSpace = kFALSE;
50  fScaleFactor = 1.;
51 
52  fSavePicture = kFALSE;
53  textInput = kFALSE;
54  theEvent = 0;
55 
56  fHighlightMode = kNoHighlight;
57  fAxesMode = kFALSE;
58 
59  fFixPerspective = kFALSE;
60 }
61 
62 
63 
66 
68 {
69  // Destructor
70 }
71 
72 
73 
76 
77 void KVEventViewer::DrawNucleus(KVNucleus* nucleus, const Char_t* frame)
78 {
79  // Draw nucleus
80 
81  Int_t N = nucleus->GetN();
82  Int_t Z = nucleus->GetZ();
83 
84  TVector3 V = (fMomentumSpace ? nucleus->GetFrame(frame, kFALSE)->GetMomentum() : nucleus->GetFrame(frame, kFALSE)->GetV());
85 
86  Bool_t Highlight = SetHighlight(nucleus);
87 
88  if (nucleus->GetA() == 1) {
89  TGeoVolume* ball = 0;
90  if (Z == 0) ball = geom->MakeSphere("n", Nuc, 0., free_nucleon_radius);
91  else if (Z == 1) ball = geom->MakeSphere("p", Nuc, 0., free_nucleon_radius);
92  int color = 0;
93  if (Z == 0) color = (Highlight ? fNeutron_color : fneutron_color);
94  else if (Z == 1) color = (Highlight ? fProton_color : fproton_color);
95  ball->SetLineColor(color);
96  top->AddNode(ball, ivol++, new TGeoTranslation(fScaleFactor * V.X(), fScaleFactor * V.Y(), fScaleFactor * V.Z()));
97  return;
98  }
99 
100  double sph_rad = pow(nuclear_radius_parameter, 3.) * (N + Z);
101 
102  for (int i = 0; i < N; i++) {
103  TGeoVolume* ball = geom->MakeSphere("n", Nuc, 0., nucleon_radius);
104  ball->SetLineColor(Highlight ? fNeutron_color : fneutron_color);
105  double rvec = pow(gRandom->Uniform(0., sph_rad), .33333);
106  double x, y, z;
107  gRandom->Sphere(x, y, z, rvec);
108  top->AddNode(ball, ivol++, new TGeoTranslation(fScaleFactor * V.X() + x, fScaleFactor * V.Y() + y, fScaleFactor * V.Z() + z));
109  }
110  for (int i = 0; i < Z; i++) {
111  TGeoVolume* ball = geom->MakeSphere("p", Nuc, 0., nucleon_radius);
112  ball->SetLineColor(Highlight ? fProton_color : fproton_color);
113  double rvec = pow(gRandom->Uniform(0., sph_rad), .33333);
114  double x, y, z;
115  gRandom->Sphere(x, y, z, rvec);
116  top->AddNode(ball, ivol++, new TGeoTranslation(fScaleFactor * V.X() + x, fScaleFactor * V.Y() + y, fScaleFactor * V.Z() + z));
117  }
118 }
119 
120 
121 
124 
126 {
127  // Draw all particles in event which are "ok"
128 
129  if (geom) delete geom;
130  geom = new TGeoManager(Form("Event%d", event->GetNumber()), "Event view");
131  geom->SetNsegments(500);
132  matEmptySpace = new TGeoMaterial("EmptySpace", 0, 0, 0);
133  matNuc = new TGeoMaterial("Nuc", .938, 1., 10000.);
134  matBox = new TGeoMaterial("Box", .938, 1., 10000.);
135 
136  EmptySpace = new TGeoMedium("Empty", 1, matEmptySpace);
137  Nuc = new TGeoMedium("Nuc", 1, matNuc);
138  Box = new TGeoMedium("Box", 1, matBox);
139 
140  top = geom->MakeBox("WORLD", EmptySpace, 30, 30, 30);
142 
143  // Find biggest velocity/momentum & biggest fragment
144  maxV = 0;
145  maxZ = 0;
146  KVNucleus* nuc;
147  if (fMomentumSpace) fScaleFactor = 1.e-2;
148  else fScaleFactor = 1.;
149  while ((nuc = event->GetNextParticle("ok"))) {
151  nuc->GetFrame(frame, kFALSE)->GetV().Mag());
152  if (v > maxV) maxV = v;
153  if (nuc->GetZ() > maxZ) maxZ = nuc->GetZ();
154  }
155  // scale down
156  maxV *= 0.5;
157  if (fFixPerspective) {
158  // fixed perspective mode
159  // calculate maxV once (first event), store it in fMaxVelocity, then use
160  // for all other events
161  if (fMaxVelocity > 0) maxV = fMaxVelocity;
162  else fMaxVelocity = maxV;
163  }
164  else if (fMaxVelocity > 0) maxV = fMaxVelocity;
165 
166  TGeoVolume* box = geom->MakeSphere("box", Box, 0.99 * maxV, maxV);
167  box->SetLineColor(kBlack);
168  box->SetTransparency(100);
169  top->AddNode(box, 1000);
170 
171  // set volume (nucleon) counter to 0
172  ivol = 0;
173 
174  if (fFixSeed) {
175  if (!(event->GetNumber() % fRefresh)) fSeed += 7;
177  }
178 
179  while ((nuc = event->GetNextParticle("ok"))) DrawNucleus(nuc, frame);
180 
181  geom->CloseGeometry();
182 
184  top->Draw("ogl");
185  TGLViewer* view = (TGLViewer*)gPad->GetViewer3D();
186  if (fAxesMode) view->SetGuideState(2, 1, 0, 0);
187  view->SetClearColor(kWhite);
188  view->SetSmoothLines(1);
189  view->SetSmoothPoints(1);
197  ((TGLOrthoCamera&)view->CurrentCamera()).SetEnableRotate(kTRUE);
198  if (fXYMode) view->SetOrthoCamera(TGLViewer::kCameraOrthoZOY, 0.9, 0, 0, 0, TMath::Pi() / 2.); //TMath::Pi()/8.,TMath::Pi()/8.);
199  else view->SetOrthoCamera(TGLViewer::kCameraOrthoZOY, 0.9, 0, 0, TMath::Pi() / 8., TMath::Pi() / 8.);
200 // if(fSavePicture) view->SavePicture(Form("Event%d.gif",event->GetNumber()));
201  if (fSavePicture) view->SavePicture(Form("Event-%d.png", event->GetNumber()));
202 }
203 
204 
205 
207 
208 void KVEventViewer::DrawEvent(KVEvent* event, const Char_t* frame) const
209 {
210  const_cast<KVEventViewer*>(this)->DrawEvent(event, frame);
211 }
212 
213 
214 
217 
219 {
220  // Read events from branch of a TTree
221 
223  eventbranch->SetAddress(&theEvent);
224  textInput = kFALSE;
225  theTree = eventbranch->GetTree();
227  eventCounter = 0;
228 }
229 
230 
231 
234 
235 void KVEventViewer::SetInput(const char* filename)
236 {
237  // Read events from file
238 
239  eventFile.open(filename);
240  if (!eventFile.good()) {
241  Error("SetInput", "Problem opening %s", filename);
242  return;
243  }
244  eventCounter = 0;
245  textInput = kTRUE;
246 }
247 
248 
249 
252 
254 {
255  // Read an event from input source
256 
257  if (textInput) {
258  if (!theEvent) theEvent = new KVEvent;
259  ReadTextEvent();
260  }
261  else {
262  ReadTreeEvent();
263  }
264 
265 }
266 
267 
268 
275 
277 {
278  // Read event from text input
279  // We assume a file containing, for each event,
280  // the multiplicity, N,
281  // followed by N lines containing Z Vx Vy Vz
282  // the atomic number and velocity components (in cm/ns)
283 
284  theEvent->Clear();
285 
286  Double_t V[3];
287  Int_t Z;
288 
289  Int_t nballs;
290  eventFile >> nballs;
291  eventCounter++;
292 
294 
295  for (int i = 0; i < nballs; i++) {
296 
297  eventFile >> Z >> V[0] >> V[1] >> V[2];
298 
299  KVNucleus* nuc = theEvent->AddParticle();
300  nuc->SetZ(Z);
301  TVector3 v(V);
302  nuc->SetVelocity(v);
303 
304  }
305 }
306 
307 
308 
310 
312 {
313  if (eventCounter == (Int_t)NtreeEntries) {
314  Warning("ReadTreeEvent", "All events have been read");
315  return;
316  }
318 }
319 
320 
321 
324 
326 {
327  // Draw next event read from input source
328 
329  ReadEvent();
331 }
332 
333 
334 
337 
339 {
340  // Decide whether or not to highlight this nucleus in the event display
341 
342  if (fHighlightMode == kHighlightZmax && n->GetZ() == maxZ) return kTRUE;
343  else if (fHighlightMode == kHighlightParameter && n->GetParameters()->IsValue("EventViewer.Highlight", 1)) return kTRUE;
344  return kFALSE;
345 }
346 
347 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define SafeDelete(p)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
kBlack
kWhite
#define N
double pow(double, double)
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
#define gPad
virtual void SetNumber(UInt_t num)
Definition: KVBase.h:209
Draw events in 3D using OpenGL.
Definition: KVEventViewer.h:52
Long64_t NtreeEntries
Definition: KVEventViewer.h:90
virtual Bool_t SetHighlight(KVNucleus *n)
Decide whether or not to highlight this nucleus in the event display.
Double_t nucleon_radius
radius of nucleons in nuclei
Definition: KVEventViewer.h:59
void ReadEvent()
Read an event from input source.
void DrawNextEvent()
Draw next event read from input source.
Double_t fMaxVelocity
Definition: KVEventViewer.h:62
void DrawNucleus(KVNucleus *, const Char_t *frame="")
Draw nucleus.
Int_t fproton_color
proton colour
Definition: KVEventViewer.h:54
Double_t fScaleFactor
Definition: KVEventViewer.h:68
Int_t fProton_color
proton colour for highlighted nuclei
Definition: KVEventViewer.h:56
Int_t eventCounter
Definition: KVEventViewer.h:86
TGeoVolume * top
Definition: KVEventViewer.h:80
virtual ~KVEventViewer()
Destructor.
Int_t fHighlightMode
Definition: KVEventViewer.h:95
TGeoMaterial * matNuc
Definition: KVEventViewer.h:74
Double_t nuclear_radius_parameter
to calculate radii of nuclei
Definition: KVEventViewer.h:60
void DrawEvent(KVEvent *, const Char_t *frame="")
Draw all particles in event which are "ok".
TGeoMedium * Nuc
Definition: KVEventViewer.h:78
void SetInput(TBranch *eventbranch)
Read events from branch of a TTree.
TGeoManager * geom
Definition: KVEventViewer.h:81
Int_t fNeutron_color
neutron colour for highlighted nuclei
Definition: KVEventViewer.h:57
ifstream eventFile
Definition: KVEventViewer.h:85
TTree * theTree
Definition: KVEventViewer.h:89
Double_t maxV
largest velocity of event
Definition: KVEventViewer.h:92
KVEvent * theEvent
Definition: KVEventViewer.h:83
TGeoMedium * Box
Definition: KVEventViewer.h:79
Double_t free_nucleon_radius
radius of free nucleons
Definition: KVEventViewer.h:58
TGeoMedium * EmptySpace
Definition: KVEventViewer.h:77
Int_t ivol
geovolume counter
Definition: KVEventViewer.h:70
Bool_t fSavePicture
kTRUE to save in GIF file
Definition: KVEventViewer.h:71
Int_t fneutron_color
neutron colour
Definition: KVEventViewer.h:55
Int_t maxZ
largest Z of event
Definition: KVEventViewer.h:93
Bool_t fAxesMode
Definition: KVEventViewer.h:96
Bool_t fMomentumSpace
kTRUE to show event in momentum space
Definition: KVEventViewer.h:67
Bool_t fFixPerspective
Definition: KVEventViewer.h:98
Bool_t textInput
Definition: KVEventViewer.h:87
TGeoMaterial * matBox
Definition: KVEventViewer.h:75
TGeoMaterial * matEmptySpace
Definition: KVEventViewer.h:73
Base class container for multi-particle events.
Definition: KVEvent.h:176
virtual void Clear(Option_t *opt="")
Definition: KVEvent.cpp:200
KVNucleus * AddParticle()
Definition: KVEvent.cpp:166
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Int_t GetA() const
Definition: KVNucleus.cpp:799
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:704
Int_t GetN() const
Return the number of neutron.
Definition: KVNucleus.cpp:781
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:770
TVector3 GetMomentum() const
Definition: KVParticle.h:565
void SetVelocity(const TVector3 &)
Set velocity of particle (in cm/ns units)
KVParticle const * GetFrame(const Char_t *frame, Bool_t warn_and_return_null_if_unknown=kTRUE) const
Definition: KVParticle.cpp:952
TVector3 GetV() const
Definition: KVParticle.h:632
TTree * GetTree() const
virtual void SetAddress(void *add)
void SetLight(ELight light, Bool_t on)
void SetGuideState(Int_t axesType, Bool_t axesDepthTest, Bool_t referenceOn, const Double_t *referencePos)
void SetCurrentCamera(ECameraType camera)
Bool_t SavePicture()
void SetSmoothPoints(Bool_t s)
void SetOrthoCamera(ECameraType camera, Double_t zoom, Double_t dolly, Double_t center[3], Double_t hRotate, Double_t vRotate)
void SetSmoothLines(Bool_t s)
TGLCamera & CurrentCamera() const
TGLLightSet * GetLightSet() const
void SetClearColor(Color_t col)
void CloseGeometry(Option_t *option="d")
TGeoVolume * MakeBox(const char *name, TGeoMedium *medium, Double_t dx, Double_t dy, Double_t dz)
TGeoVolume * MakeSphere(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t themin=0, Double_t themax=180, Double_t phimin=0, Double_t phimax=360)
void SetTopVolume(TGeoVolume *vol)
void SetNsegments(Int_t nseg)
virtual void Draw(Option_t *option="")
virtual void SetVisContainers(Bool_t flag=kTRUE)
virtual void SetLineColor(Color_t lcolor)
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void SetSeed(ULong_t seed=0)
virtual Double_t Uniform(Double_t x1, Double_t x2)
virtual void Sphere(Double_t &x, Double_t &y, Double_t &z, Double_t r)
virtual Int_t GetEntry(Long64_t entry, Int_t getall=0)
virtual Long64_t GetEntries() const
Double_t Z() const
Double_t Y() const
Double_t X() const
Double_t Mag() const
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Double_t y[n]
Double_t x[n]
const Int_t n
constexpr Double_t Pi()
v