]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveTrackCounter.cxx
Improved quitting Event Display with close button.
[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         fTrackLists.Clear();
107   }
108   {
109     TIter next(&fTrackletLists);
110     TEveTrackList* tlist;
111     while ((tlist = dynamic_cast<TEveTrackList*>(next())))
112       tlist->DecDenyDestroy();
113     //fTrackletLists.Clear("nodelete");
114         fTrackletLists.Clear();
115   }
116
117   fEventId = AliEveEventManager::GetMaster()->GetEventId();
118 }
119
120 //______________________________________________________________________________
121 void AliEveTrackCounter::RegisterTracks(TEveTrackList* tlist, Bool_t goodTracks)
122 {
123   // Register tracks from tlist and tlist itself.
124   // If goodTracks is true, they are considered as primary/good
125   // tracks.
126
127   tlist->IncDenyDestroy();
128   fTrackLists.Add(tlist);
129
130   List_i i = tlist->BeginChildren();
131   while (i != tlist->EndChildren())
132   {
133     AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
134     if (t != 0)
135     {
136       if (goodTracks)
137       {
138         ++fGoodTracks;
139         t->GetESDTrack()->SetLabel(3);
140       } else {
141         t->SetLineStyle(fBadLineStyle);
142         t->GetESDTrack()->SetLabel(0);
143       }
144       ++fAllTracks;
145     }
146     ++i;
147   }
148 }
149
150 //______________________________________________________________________________
151 void AliEveTrackCounter::RegisterTracklets(TEveTrackList* tlist, Bool_t goodTracks)
152 {
153   // Register tracklets from tlist and tlist itself.
154   // If goodTracks is true, they are considered as primary/good
155   // tracks.
156
157   AliESDEvent     *esd = AliEveEventManager::AssertESD();
158   AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
159
160   tlist->IncDenyDestroy();
161   fTrackletLists.Add(tlist);
162
163   List_i i = tlist->BeginChildren();
164   while (i != tlist->EndChildren())
165   {
166     AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
167     if (t != 0)
168     {
169       if (goodTracks)
170       {
171         mul->SetLabel(t->GetIndex(), 0, 3);
172         ++fGoodTracklets;
173       } else {
174         mul->SetLabel(t->GetIndex(), 0, 0);
175         t->SetLineStyle(fBadLineStyle);
176       }
177       ++fAllTracklets;
178     }
179     ++i;
180   }
181 }
182
183 //______________________________________________________________________________
184 void AliEveTrackCounter::DoTrackAction(AliEveTrack* track)
185 {
186    // Slot called when track is secondary selected.
187    //
188    // No check is done if track actually belongs to one of the
189    // registered track-lists.
190    //
191    // Probably it would be safer to copy good/bad tracks into special
192    // sub-containers.
193    // In this case one should also override RemoveElementLocal.
194
195    static const TEveException eh("AliEveTrackCounter::DoTrackAction ");
196
197    if (!fActive)
198      return;
199
200    switch (fClickAction)
201    {
202
203       case kCA_PrintTrackInfo:
204       {
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());
210          break;
211       }
212
213       case kCA_ToggleTrack:
214       {
215          AliESDtrack *esdt = track->GetESDTrack();
216          if (track->GetLineStyle() == 1)
217          {
218             track->SetLineStyle(fBadLineStyle);
219             esdt->SetLabel(esdt->GetLabel() & ~1);
220             --fGoodTracks;
221          } else {
222             track->SetLineStyle(1);
223             esdt->SetLabel(esdt->GetLabel() | 1);
224             ++fGoodTracks;
225          }
226          track->ElementChanged();
227          gEve->Redraw3D();
228
229          //printf("AliEveTrackCounter::DoTrackAction All=%d, Good=%d, Bad=%d\n",
230          //       fAllTracks, fGoodTracks, fAllTracks-fGoodTracks);
231
232          if (gEve->GetEditor()->GetModel() == GetObject(eh))
233             gEve->EditElement(this);
234
235          break;
236       }
237
238    } // end switch fClickAction
239 }
240
241 //______________________________________________________________________________
242 void AliEveTrackCounter::DoTrackletAction(AliEveTracklet* track)
243 {
244    // Slot called when tracklet is secondary selected.
245    //
246    // No check is done if track actually belongs to one of the
247    // registered track-lists.
248    //
249    // Probably it would be safer to copy good/bad tracks into special
250    // sub-containers.
251    // In this case one should also override RemoveElementLocal.
252
253    static const TEveException eh("AliEveTrackCounter::DoTrackletAction ");
254
255    if (!fActive)
256      return;
257
258    switch (fClickAction)
259    {
260
261       case kCA_PrintTrackInfo:
262       {
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());
268          break;
269       }
270
271       case kCA_ToggleTrack:
272       {
273          AliESDEvent     *esd = AliEveEventManager::AssertESD();
274          AliMultiplicity *mul = const_cast<AliMultiplicity*>(esd->GetMultiplicity());
275
276          if (track->GetLineStyle() == 1)
277          {
278             track->SetLineStyle(fBadLineStyle);
279             mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) & ~1);
280             --fGoodTracklets;
281          } else {
282             track->SetLineStyle(1);
283             mul->SetLabel(track->GetIndex(), 0, mul->GetLabel(track->GetIndex(), 0) | 1);
284             ++fGoodTracklets;
285          }
286          track->ElementChanged();
287          gEve->Redraw3D();
288
289          // printf("AliEveTrackCounter::DoTrackletAction All=%d, Good=%d, Bad=%d\n",
290          //        fAllTracklets, fGoodTracklets, fAllTracklets-fGoodTracklets);
291
292          if (gEve->GetEditor()->GetModel() == GetObject(eh))
293             gEve->EditElement(this);
294
295          break;
296       }
297
298    } // end switch fClickAction
299 }
300
301 /******************************************************************************/
302
303 //______________________________________________________________________________
304 void AliEveTrackCounter::OutputEventTracks()
305 {
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.
310
311   {
312     TFile *f = TFile::Open("scan_results.root", "UPDATE");
313
314     AliESDEvent     *esd = AliEveEventManager::AssertESD();
315     TClonesArray    *trk = static_cast<TClonesArray*>   (esd->GetList()->FindObject("Tracks"));
316     AliMultiplicity *mul = const_cast <AliMultiplicity*>(esd->GetMultiplicity());
317
318     trk->Write(TString::Format("Tracks_%04d",    fEventId), kWriteDelete | kSingleKey);
319     mul->Write(TString::Format("Tracklets_%04d", fEventId), kWriteDelete);
320
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);
324
325     f->Close();
326   }
327 }
328
329 //______________________________________________________________________________
330 void AliEveTrackCounter::PrintEventTracks()
331 {
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.
336
337   FILE* out = stdout;
338
339   fprintf(out, "AliEveTrackCounter::PrintEventTracks()\n");
340
341   fprintf(out, "Event=%d\n", fEventId);
342   fprintf(out, "GoodTracks=%d  AllTracks=%d\n", fGoodTracks, fAllTracks);
343
344   {
345     TIter tlists(&fTrackLists);
346     TEveTrackList* tlist;
347     Int_t cnt = 0;
348     while ((tlist = (TEveTrackList*) tlists()) != 0)
349     {
350       List_i i = tlist->BeginChildren();
351       while (i != tlist->EndChildren())
352       {
353         AliEveTrack* t = dynamic_cast<AliEveTrack*>(*i);
354         if (t != 0 && t->GetLineStyle() == 1)
355         {
356           ++cnt;
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());
360         }
361         ++i;
362       }
363     }
364   }
365
366   fprintf(out, "GoodTracklets=%d  AllTracklets=%d\n", fGoodTracklets, fAllTracklets);
367   {
368     TIter tlists(&fTrackletLists);
369     TEveTrackList* tlist;
370     Int_t cnt = 0;
371     while ((tlist = (TEveTrackList*) tlists()) != 0)
372     {
373       List_i i = tlist->BeginChildren();
374       while (i != tlist->EndChildren())
375       {
376         AliEveTracklet* t = dynamic_cast<AliEveTracklet*>(*i);
377         if (t != 0 && t->GetLineStyle() == 1)
378         {
379           ++cnt;
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());
382         }
383         ++i;
384       }
385     }
386   }
387 }