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