2 // Author: Matevz Tadel 2007
4 /**************************************************************************
5 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
7 * full copyright notice. *
8 **************************************************************************/
10 #include "AliEveTrackCounter.h"
12 #include "TEveManager.h"
13 #include "AliEveEventManager.h"
14 #include "AliEveTrack.h"
15 #include "AliEveTracklet.h"
17 #include "AliESDEvent.h"
18 #include "AliESDtrack.h"
19 #include "AliMultiplicity.h"
21 #include <TEveGedEditor.h>
25 //==============================================================================
27 //==============================================================================
29 //______________________________________________________________________________
31 // Provides event-based method for tagging of good / bad (or primary /
32 // secondary) tracks. A report can be written into a text file.
34 // AliEveTrack status is toggled by using secondary-selection / ctrl-click
35 // functionality of the GL viewer.
37 // Some of the functionality is implemented in AliEveTrackCounterEditor
40 ClassImp(AliEveTrackCounter)
42 //______________________________________________________________________________
43 AliEveTrackCounter* AliEveTrackCounter::fgInstance = 0;
45 Bool_t AliEveTrackCounter::IsActive()
47 // Check if instance exists and is active.
49 return fgInstance && fgInstance->fActive;
52 //______________________________________________________________________________
53 AliEveTrackCounter::AliEveTrackCounter(const Text_t* name, const Text_t* title) :
58 fClickAction (kCA_ToggleTrack),
60 fAllTracks (0), fGoodTracks (0),
61 fAllTracklets (0), fGoodTracklets(0),
62 fTrackLists (), fTrackletLists(),
66 // Connects to global signal "AliEveTrack", "SecSelected(AliEveTrack*)".
68 if (fgInstance == 0) fgInstance = this;
70 TQObject::Connect("AliEveTrack", "SecSelectedTrack(AliEveTrack*)",
71 "AliEveTrackCounter", this, "DoTrackAction(AliEveTrack*)");
72 TQObject::Connect("AliEveTracklet", "SecSelectedTracklet(AliEveTracklet*)",
73 "AliEveTrackCounter", this, "DoTrackletAction(AliEveTracklet*)");
75 AliEveEventManager::GetMaster()->Connect("NewEventDataLoaded()", "AliEveTrackCounter", this, "Reset()");
78 //______________________________________________________________________________
79 AliEveTrackCounter::~AliEveTrackCounter()
82 // Disconnect from the global track signals.
84 AliEveEventManager::GetMaster()->Disconnect("NewEventDataLoaded()", this);
86 TQObject::Disconnect("AliEveTrack", "DoTrackAction(AliEveTrack*)");
87 TQObject::Disconnect("AliEveTracklet", "DoTrackletAction(AliEveTracklet*)");
88 if (fgInstance == this) fgInstance = 0;
91 /******************************************************************************/
93 //______________________________________________________________________________
94 void AliEveTrackCounter::Reset()
96 // Reset internal track-counters and track-list.
98 fAllTracks = fGoodTracks = 0;
99 fAllTracklets = fGoodTracklets = 0;
101 TIter next(&fTrackLists);
102 TEveTrackList* tlist;
103 while ((tlist = dynamic_cast<TEveTrackList*>(next())))
104 tlist->DecDenyDestroy();
105 //fTrackLists.Clear("nodelete");
109 TIter next(&fTrackletLists);
110 TEveTrackList* tlist;
111 while ((tlist = dynamic_cast<TEveTrackList*>(next())))
112 tlist->DecDenyDestroy();
113 //fTrackletLists.Clear("nodelete");
114 fTrackletLists.Clear();
117 fEventId = AliEveEventManager::GetMaster()->GetEventId();
120 //______________________________________________________________________________
121 void AliEveTrackCounter::RegisterTracks(TEveTrackList* tlist, Bool_t goodTracks)
123 // Register tracks from tlist and tlist itself.
124 // If goodTracks is true, they are considered as primary/good
127 tlist->IncDenyDestroy();
128 fTrackLists.Add(tlist);
130 List_i i = tlist->BeginChildren();
131 while (i != tlist->EndChildren())
133 AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
139 t->GetESDTrack()->SetLabel(3);
141 t->SetLineStyle(fBadLineStyle);
142 t->GetESDTrack()->SetLabel(0);
150 //______________________________________________________________________________
151 void AliEveTrackCounter::RegisterTracklets(TEveTrackList* tlist, Bool_t goodTracks)
153 // Register tracklets from tlist and tlist itself.
154 // If goodTracks is true, they are considered as primary/good
157 AliESDEvent *esd = AliEveEventManager::AssertESD();
158 AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
160 tlist->IncDenyDestroy();
161 fTrackletLists.Add(tlist);
163 List_i i = tlist->BeginChildren();
164 while (i != tlist->EndChildren())
166 AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
171 mul->SetLabel(t->GetIndex(), 0, 3);
174 mul->SetLabel(t->GetIndex(), 0, 0);
175 t->SetLineStyle(fBadLineStyle);
183 //______________________________________________________________________________
184 void AliEveTrackCounter::DoTrackAction(AliEveTrack* track)
186 // Slot called when track is secondary selected.
188 // No check is done if track actually belongs to one of the
189 // registered track-lists.
191 // Probably it would be safer to copy good/bad tracks into special
193 // In this case one should also override RemoveElementLocal.
195 static const TEveException eh("AliEveTrackCounter::DoTrackAction ");
200 switch (fClickAction)
203 case kCA_PrintTrackInfo:
205 printf("AliEveTrack '%s'\n", track->GetObject(eh)->GetName());
206 const TEveVector &v = track->GetVertex();
207 const TEveVector &p = track->GetMomentum();;
208 printf(" Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n",
209 v.fX, v.fY, v.fZ, p.Perp(), p.fZ, TMath::RadToDeg()*p.Phi());
213 case kCA_ToggleTrack:
215 AliESDtrack *esdt = track->GetESDTrack();
216 if (track->GetLineStyle() == 1)
218 track->SetLineStyle(fBadLineStyle);
219 esdt->SetLabel(esdt->GetLabel() & ~1);
222 track->SetLineStyle(1);
223 esdt->SetLabel(esdt->GetLabel() | 1);
226 track->ElementChanged();
229 //printf("AliEveTrackCounter::DoTrackAction All=%d, Good=%d, Bad=%d\n",
230 // fAllTracks, fGoodTracks, fAllTracks-fGoodTracks);
232 if (gEve->GetEditor()->GetModel() == GetObject(eh))
233 gEve->EditElement(this);
238 } // end switch fClickAction
241 //______________________________________________________________________________
242 void AliEveTrackCounter::DoTrackletAction(AliEveTracklet* track)
244 // Slot called when tracklet is secondary selected.
246 // No check is done if track actually belongs to one of the
247 // registered track-lists.
249 // Probably it would be safer to copy good/bad tracks into special
251 // In this case one should also override RemoveElementLocal.
253 static const TEveException eh("AliEveTrackCounter::DoTrackletAction ");
258 switch (fClickAction)
261 case kCA_PrintTrackInfo:
263 printf("AliEveTracklet '%s'\n", track->GetObject(eh)->GetName());
264 const TEveVector &v = track->GetVertex();
265 const TEveVector &p = track->GetMomentum();;
266 printf(" Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n",
267 v.fX, v.fY, v.fZ, p.Perp(), p.fZ, TMath::RadToDeg()*p.Phi());
271 case kCA_ToggleTrack:
273 AliESDEvent *esd = AliEveEventManager::AssertESD();
274 AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
276 if (track->GetLineStyle() == 1)
278 track->SetLineStyle(fBadLineStyle);
279 mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) & ~1);
282 track->SetLineStyle(1);
283 mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) | 1);
286 track->ElementChanged();
289 // printf("AliEveTrackCounter::DoTrackletAction All=%d, Good=%d, Bad=%d\n",
290 // fAllTracklets, fGoodTracklets, fAllTracklets-fGoodTracklets);
292 if (gEve->GetEditor()->GetModel() == GetObject(eh))
293 gEve->EditElement(this);
298 } // end switch fClickAction
301 /******************************************************************************/
303 //______________________________________________________________________________
304 void AliEveTrackCounter::OutputEventTracks()
306 // Print good-track summary into a plain-text file by iteration
307 // through all registered track-lists.
308 // State of each track is determined by its line-style, it is
309 // considered a good track if it's line style is solid.
312 TFile *f = TFile::Open("scan_results.root", "UPDATE");
314 AliESDEvent *esd = AliEveEventManager::AssertESD();
315 TClonesArray *trk = static_cast<TClonesArray*> (esd->GetList()->FindObject("Tracks"));
316 AliMultiplicity *mul = const_cast <AliMultiplicity*>(esd->GetMultiplicity());
318 trk->Write(TString::Format("Tracks_%04d", fEventId), kWriteDelete | kSingleKey);
319 mul->Write(TString::Format("Tracklets_%04d", fEventId), kWriteDelete);
321 esd->GetPrimaryVertexTracks()->Write(TString::Format("PrimVertTracks_%04d", fEventId), kWriteDelete);
322 esd->GetPrimaryVertexTPC() ->Write(TString::Format("PrimVertTPC_%04d", fEventId), kWriteDelete);
323 esd->GetPrimaryVertexSPD() ->Write(TString::Format("PrimVertSPD_%04d", fEventId), kWriteDelete);
329 //______________________________________________________________________________
330 void AliEveTrackCounter::PrintEventTracks()
332 // Print good-track summary to stdout by iteration
333 // through all registered track-lists.
334 // State of each track is determined by its line-style, it is
335 // considered a good track if it's line style is solid.
339 fprintf(out, "AliEveTrackCounter::PrintEventTracks()\n");
341 fprintf(out, "Event=%d\n", fEventId);
342 fprintf(out, "GoodTracks=%d AllTracks=%d\n", fGoodTracks, fAllTracks);
345 TIter tlists(&fTrackLists);
346 TEveTrackList* tlist;
348 while ((tlist = (TEveTrackList*) tlists()) != 0)
350 List_i i = tlist->BeginChildren();
351 while (i != tlist->EndChildren())
353 AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
354 if (t != 0 && t->GetLineStyle() == 1)
357 fprintf(out, " %2d: chg=%+2d pt=%8.5f eta=%+8.5f phi=%+8.5f\n",
358 cnt, t->GetCharge(), t->GetMomentum().Perp(),
359 t->GetMomentum().Eta(), t->GetMomentum().Phi());
366 fprintf(out, "GoodTracklets=%d AllTracklets=%d\n", fGoodTracklets, fAllTracklets);
368 TIter tlists(&fTrackletLists);
369 TEveTrackList* tlist;
371 while ((tlist = (TEveTrackList*) tlists()) != 0)
373 List_i i = tlist->BeginChildren();
374 while (i != tlist->EndChildren())
376 AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
377 if (t != 0 && t->GetLineStyle() == 1)
380 fprintf(out, " %2d: theta=%+8.5f eta=%+8.5f phi=%+8.5f\n",
381 cnt, t->GetMomentum().Theta(), t->GetMomentum().Eta(), t->GetMomentum().Phi());