]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveTrackCounter.cxx
94b59ede3a40889789f852c2d5564b4e961ac8fa
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveTrackCounter.cxx
1 // @(#)root/eve:$Id$
2 // Author: Matevz Tadel 2007
3
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  **************************************************************************/
9
10 #include "AliEveTrackCounter.h"
11
12 #include "TEveManager.h"
13 #include "AliEveEventManager.h"
14 #include "AliEveTrack.h"
15 #include "AliEveTracklet.h"
16
17 #include "AliESDEvent.h"
18 #include "AliESDtrack.h"
19 #include "AliMultiplicity.h"
20
21 #include <TEveGedEditor.h>
22
23 #include <TFile.h>
24
25 //==============================================================================
26 // AliEveTrackCounter
27 //==============================================================================
28
29 //______________________________________________________________________________
30 //
31 // Provides event-based method for tagging of good / bad (or primary /
32 // secondary) tracks. A report can be written into a text file.
33 //
34 // AliEveTrack status is toggled by using secondary-selection / ctrl-click
35 // functionality of the GL viewer.
36 //
37 // Some of the functionality is implemented in AliEveTrackCounterEditor
38 // class.
39
40 ClassImp(AliEveTrackCounter)
41
42 //______________________________________________________________________________
43 AliEveTrackCounter* AliEveTrackCounter::fgInstance = 0;
44
45 //______________________________________________________________________________
46 AliEveTrackCounter::AliEveTrackCounter(const Text_t* name, const Text_t* title) :
47   TEveElement(),
48   TNamed(name, title),
49
50   fBadLineStyle (6),
51   fClickAction  (kCA_ToggleTrack),
52   fEventId      (-1),
53   fAllTracks    (0), fGoodTracks   (0),
54   fAllTracklets (0), fGoodTracklets(0),
55   fTrackLists   (),  fTrackletLists()
56 {
57   // Constructor.
58   // Connects to global signal "AliEveTrack", "SecSelected(AliEveTrack*)".
59
60   if (fgInstance == 0) fgInstance = this;
61
62   TQObject::Connect("AliEveTrack", "SecSelectedTrack(AliEveTrack*)",
63                     "AliEveTrackCounter", this, "DoTrackAction(AliEveTrack*)");
64   TQObject::Connect("AliEveTracklet", "SecSelectedTracklet(AliEveTracklet*)",
65                     "AliEveTrackCounter", this, "DoTrackletAction(AliEveTracklet*)");
66 }
67
68 //______________________________________________________________________________
69 AliEveTrackCounter::~AliEveTrackCounter()
70 {
71   // Destructor.
72   // Disconnect from the global track signals.
73
74   TQObject::Disconnect("AliEveTrack", "DoTrackAction(AliEveTrack*)");
75   TQObject::Disconnect("AliEveTracklet", "DoTrackletAction(AliEveTracklet*)");
76   if (fgInstance == this) fgInstance = 0;
77 }
78
79 /******************************************************************************/
80
81 //______________________________________________________________________________
82 void AliEveTrackCounter::Reset()
83 {
84   // Reset internal track-counters and track-list.
85
86   fAllTracks     = fGoodTracks    = 0;
87   fAllTracklets  = fGoodTracklets = 0;
88   {
89     TIter next(&fTrackLists);
90     TEveTrackList* tlist;
91     while ((tlist = dynamic_cast<TEveTrackList*>(next())))
92       tlist->DecDenyDestroy();
93     fTrackLists.Clear("nodelete");
94   }
95   {
96     TIter next(&fTrackletLists);
97     TEveTrackList* tlist;
98     while ((tlist = dynamic_cast<TEveTrackList*>(next())))
99       tlist->DecDenyDestroy();
100     fTrackletLists.Clear("nodelete");
101   }
102 }
103
104 //______________________________________________________________________________
105 void AliEveTrackCounter::RegisterTracks(TEveTrackList* tlist, Bool_t goodTracks)
106 {
107   // Register tracks from tlist and tlist itself.
108   // If goodTracks is true, they are considered as primary/good
109   // tracks.
110
111   tlist->IncDenyDestroy();
112   fTrackLists.Add(tlist);
113
114   List_i i = tlist->BeginChildren();
115   while (i != tlist->EndChildren())
116   {
117     AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
118     if (t != 0)
119     {
120       if (goodTracks)
121       {
122         ++fGoodTracks;
123         t->GetESDTrack()->SetLabel(3);
124       } else {
125         t->SetLineStyle(fBadLineStyle);
126         t->GetESDTrack()->SetLabel(0);
127       }
128       ++fAllTracks;
129     }
130     ++i;
131   }
132 }
133
134 //______________________________________________________________________________
135 void AliEveTrackCounter::RegisterTracklets(TEveTrackList* tlist, Bool_t goodTracks)
136 {
137   // Register tracklets from tlist and tlist itself.
138   // If goodTracks is true, they are considered as primary/good
139   // tracks.
140
141   AliESDEvent     *esd = AliEveEventManager::AssertESD();
142   AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
143
144   tlist->IncDenyDestroy();
145   fTrackletLists.Add(tlist);
146
147   List_i i = tlist->BeginChildren();
148   while (i != tlist->EndChildren())
149   {
150     AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
151     if (t != 0)
152     {
153       if (goodTracks && mul->GetLabel(t->GetIndex(), 0) == 3)
154       {
155         ++fGoodTracklets;
156       } else {
157         t->SetLineStyle(fBadLineStyle);
158       }
159       ++fAllTracklets;
160     }
161     ++i;
162   }
163 }
164
165 //______________________________________________________________________________
166 void AliEveTrackCounter::DoTrackAction(AliEveTrack* track)
167 {
168    // Slot called when track is secondary selected.
169    //
170    // No check is done if track actually belongs to one of the
171    // registered track-lists.
172    //
173    // Probably it would be safer to copy good/bad tracks into special
174    // sub-containers.
175    // In this case one should also override RemoveElementLocal.
176
177    static const TEveException eh("AliEveTrackCounter::DoTrackAction ");
178
179    switch (fClickAction)
180    {
181
182       case kCA_PrintTrackInfo:
183       {
184          printf("AliEveTrack '%s'\n", track->GetObject(eh)->GetName());
185          const TEveVector &v = track->GetVertex();
186          const TEveVector &p = track->GetMomentum();;
187          printf("  Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n",
188                 v.fX, v.fY, v.fZ, p.Perp(), p.fZ, TMath::RadToDeg()*p.Phi());
189          printf("  <other information should be printed ... full AliESDtrack>\n");
190          break;
191       }
192
193       case kCA_ToggleTrack:
194       {
195          AliESDtrack *esdt = track->GetESDTrack();
196          if (track->GetLineStyle() == 1)
197          {
198             track->SetLineStyle(fBadLineStyle);
199             esdt->SetLabel(esdt->GetLabel() & ~1);
200             --fGoodTracks;
201          } else {
202             track->SetLineStyle(1);
203             esdt->SetLabel(esdt->GetLabel() | 1);
204             ++fGoodTracks;
205          }
206          track->ElementChanged();
207          gEve->Redraw3D();
208
209          printf("AliEveTrackCounter::DoTrackAction All=%d, Good=%d, Bad=%d\n",
210                 fAllTracks, fGoodTracks, fAllTracks-fGoodTracks);
211
212          if (gEve->GetEditor()->GetModel() == GetObject(eh))
213             gEve->EditElement(this);
214
215          break;
216       }
217
218    } // end switch fClickAction
219 }
220
221 //______________________________________________________________________________
222 void AliEveTrackCounter::DoTrackletAction(AliEveTracklet* track)
223 {
224    // Slot called when tracklet is secondary selected.
225    //
226    // No check is done if track actually belongs to one of the
227    // registered track-lists.
228    //
229    // Probably it would be safer to copy good/bad tracks into special
230    // sub-containers.
231    // In this case one should also override RemoveElementLocal.
232
233    static const TEveException eh("AliEveTrackCounter::DoTrackletAction ");
234
235    switch (fClickAction)
236    {
237
238       case kCA_PrintTrackInfo:
239       {
240          printf("AliEveTracklet '%s'\n", track->GetObject(eh)->GetName());
241          const TEveVector &v = track->GetVertex();
242          const TEveVector &p = track->GetMomentum();;
243          printf("  Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n",
244                 v.fX, v.fY, v.fZ, p.Perp(), p.fZ, TMath::RadToDeg()*p.Phi());
245          printf("  <other information should be printed ... full AliESDtrack>\n");
246          break;
247       }
248
249       case kCA_ToggleTrack:
250       {
251          AliESDEvent     *esd = AliEveEventManager::AssertESD();
252          AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
253
254          if (track->GetLineStyle() == 1)
255          {
256             track->SetLineStyle(fBadLineStyle);
257             mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) & ~1);
258             --fGoodTracklets;
259          } else {
260             track->SetLineStyle(1);
261             mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) | 1);
262             ++fGoodTracklets;
263          }
264          track->ElementChanged();
265          gEve->Redraw3D();
266
267          printf("AliEveTrackCounter::DoTrackletAction All=%d, Good=%d, Bad=%d\n",
268                 fAllTracklets, fGoodTracklets, fAllTracklets-fGoodTracklets);
269
270          if (gEve->GetEditor()->GetModel() == GetObject(eh))
271             gEve->EditElement(this);
272
273          break;
274       }
275
276    } // end switch fClickAction
277 }
278
279 /******************************************************************************/
280
281 //______________________________________________________________________________
282 void AliEveTrackCounter::OutputEventTracks()
283 {
284   // Print good-track summary into a plain-text file by iteration
285   // through all registered track-lists.
286   // State of each track is determined by its line-style, it is
287   // considered a good track if it's line style is solid.
288
289   {
290     TFile *f = TFile::Open("scan_results.root", "UPDATE");
291
292     AliESDEvent     *esd = AliEveEventManager::AssertESD();
293     TClonesArray    *trk = static_cast<TClonesArray*>   (esd->GetList()->FindObject("Tracks"));
294     AliMultiplicity *mul = const_cast <AliMultiplicity*>(esd->GetMultiplicity());
295
296     trk->Write(TString::Format("Tracks_%04d",    fEventId), kWriteDelete | kSingleKey);
297     mul->Write(TString::Format("Tracklets_%04d", fEventId), kWriteDelete);
298
299     f->Close();
300   }
301 }
302
303 //______________________________________________________________________________
304 void AliEveTrackCounter::PrintEventTracks()
305 {
306   // Print good-track summary to stdout 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.
310
311   FILE* out = stdout;
312
313   fprintf(out, "AliEveTrackCounter::PrintEventTracks()\n");
314
315   fprintf(out, "Event=%d\n", fEventId);
316   fprintf(out, "GoodTracks=%d  AllTracks=%d\n", fGoodTracks, fAllTracks);
317
318   {
319     TIter tlists(&fTrackLists);
320     TEveTrackList* tlist;
321     Int_t cnt = 0;
322     while ((tlist = (TEveTrackList*) tlists()) != 0)
323     {
324       List_i i = tlist->BeginChildren();
325       while (i != tlist->EndChildren())
326       {
327         AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
328         if (t != 0 && t->GetLineStyle() == 1)
329         {
330           ++cnt;
331           fprintf(out, " %2d: chg=%+2d  pt=%8.5f  eta=%+8.5f  phi=%+8.5f\n",
332                   cnt, t->GetCharge(), t->GetMomentum().Perp(),
333                   t->GetMomentum().Eta(), t->GetMomentum().Phi());
334         }
335         ++i;
336       }
337     }
338   }
339
340   fprintf(out, "GoodTracklets=%d  AllTracklets=%d\n", fGoodTracklets, fAllTracklets);
341   {
342     TIter tlists(&fTrackletLists);
343     TEveTrackList* tlist;
344     Int_t cnt = 0;
345     while ((tlist = (TEveTrackList*) tlists()) != 0)
346     {
347       List_i i = tlist->BeginChildren();
348       while (i != tlist->EndChildren())
349       {
350         AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
351         if (t != 0 && t->GetLineStyle() == 1)
352         {
353           ++cnt;
354           fprintf(out, " %2d: theta=%+8.5f  eta=%+8.5f  phi=%+8.5f\n",
355                   cnt, t->GetMomentum().Theta(), t->GetMomentum().Eta(), t->GetMomentum().Phi());
356         }
357         ++i;
358       }
359     }
360   }
361 }