]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveKineTools.cxx
* show_scan_results.C
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveKineTools.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 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 "AliEveKineTools.h"
11 #include "AliEveTrack.h"
12
13 #include <AliStack.h>
14 #include <AliTrackReference.h>
15
16 #include <TTree.h>
17 #include <TBranchElement.h>
18 #include <TClonesArray.h>
19
20 #include <map>
21
22 //______________________________________________________________________________
23 // AliEveKineTools
24 //
25 // Tools for import of kinematics.
26 //
27
28 ClassImp(AliEveKineTools)
29
30 namespace {
31
32   // Map to store label-to-track association.
33   //
34   // multimap is used as there are cases when initial particles (in
35   // particular resonances) are not assigned proper status-codes
36   // and can thus be found several times in the eve-track-list.
37
38   typedef std::multimap<Int_t, AliEveTrack*>                 TrackMap_t;
39   typedef std::multimap<Int_t, AliEveTrack*>::const_iterator TrackMap_ci;
40
41   void MapTracks(TrackMap_t& map, TEveElement* cont, Bool_t recurse)
42   {
43     TEveElement::List_i i = cont->BeginChildren();
44     while (i != cont->EndChildren()) {
45       AliEveTrack* track = dynamic_cast<AliEveTrack*>(*i);
46       map.insert(std::make_pair(track->GetLabel(), track));
47       if (recurse)
48         MapTracks(map, track, recurse);
49       ++i;
50     }
51   }
52 }
53
54 /**************************************************************************/
55
56 void AliEveKineTools::SetDaughterPathMarks(TEveElement* cont, AliStack* stack, Bool_t recurse)
57 {
58   // Import daughters birth points.
59
60   TEveElement::List_i  iter = cont->BeginChildren();
61   while(iter != cont->EndChildren())
62   {
63     AliEveTrack* track = dynamic_cast<AliEveTrack*>(*iter);
64     TParticle* p = stack->Particle(track->GetLabel());
65     if (p->GetNDaughters())
66     {
67       Int_t d0 = p->GetDaughter(0), d1 = p->GetDaughter(1);
68       for(int d = d0; d > 0 && d <= d1; ++d)
69       {
70         TParticle* dp = stack->Particle(d);
71         track->AddPathMark(TEvePathMark(TEvePathMark::kDaughter,
72                                         TEveVector(dp->Vx(), dp->Vy(), dp->Vz()),
73                                         TEveVector(dp->Px(), dp->Py(), dp->Pz()),
74                                         dp->T()));
75         // printf("Daughter path-mark for %d, %d, t=%e, r=%f,%f,%f\n",
76         //        track->GetLabel(), d, dp->T(), dp->Vx(), dp->Vy(), dp->Vz());
77       }
78
79       // Check last process, set decay if needed.
80       Int_t lp = stack->Particle(d1)->GetUniqueID();
81       if (lp != kPBrem && lp != kPDeltaRay && lp < kPCerenkov)
82       {
83         TParticle* dp = stack->Particle(d1);
84         track->AddPathMark(TEvePathMark(TEvePathMark::kDecay,
85                                         TEveVector(dp->Vx(), dp->Vy(), dp->Vz()),
86                                         TEveVector(0, 0,0),  dp->T()));
87       }
88
89       if (recurse)
90         SetDaughterPathMarks(track, stack, recurse);
91     }
92     ++iter;
93   }
94 }
95
96 /**************************************************************************/
97
98 void AliEveKineTools::SetTrackReferences(TEveElement* cont, TTree* treeTR, Bool_t recurse)
99 {
100   // Set decay and track reference path-marks.
101
102   static const TEveException kEH("AliEveKineTools::ImportPathMarks");
103
104   TrackMap_t map;
105   MapTracks(map, cont, recurse);
106
107   TClonesArray* arr = 0;
108   treeTR->SetBranchAddress("TrackReferences", &arr);
109
110   Int_t nTreeEntries = (Int_t) treeTR->GetEntries();
111
112   for (Int_t te = 0; te < nTreeEntries; ++te)
113   {
114     treeTR->GetEntry(te);
115
116     Int_t last_label = -1;
117     std::pair<TrackMap_ci, TrackMap_ci> range;
118     Int_t nArrEntries = arr->GetEntriesFast();
119
120     // printf("tree-entry %d, n-arr-entries %d\n", te, nArrEntries);
121
122     for (Int_t ae = 0; ae < nArrEntries; ++ae)
123     {
124       AliTrackReference* atr = (AliTrackReference*)arr->UncheckedAt(ae);
125       Bool_t isRef = (atr->DetectorId() != -1);
126       Int_t  label = atr->GetTrack();
127
128       // printf("    arr-entry %d, label %d, detid %d, len=%f, t=%e r=%f,%f,%f\n",
129       //        ae, label, atr->DetectorId(), 
130       //        atr->GetLength(), atr->GetTime(), atr->X(), atr->Y(), atr->Z());
131
132       if (label < 0)
133         throw(kEH + Form("negative label for array-entry %d in tree-entry %d.",
134                          ae, te));
135
136       if (label != last_label)
137       {
138         range      = map.equal_range(label);
139         last_label = label;
140       }
141
142       for (TrackMap_ci i = range.first; i != range.second; ++i)
143       {
144         i->second->AddPathMark
145           (TEvePathMark(isRef ? TEvePathMark::kReference : TEvePathMark::kDecay,
146                         TEveVector(atr->X(),  atr->Y(),  atr->Z()),
147                         TEveVector(atr->Px(), atr->Py(), atr->Pz()),
148                         atr->GetTime()));
149       }
150     }
151   }
152   delete arr;
153 }
154
155 /**************************************************************************/
156
157 void AliEveKineTools::SortPathMarks(TEveElement* el, Bool_t recurse)
158 {
159   // Sort path-marks for track by time.
160   // If recurse is true, descends down all track children.
161
162   AliEveTrack* track = dynamic_cast<AliEveTrack*>(el);
163   if (track)
164     track->SortPathMarksByTime();
165
166   if (recurse)
167   {
168     for (TEveElement::List_i i = el->BeginChildren(); i != el->EndChildren(); ++i)
169       SortPathMarks(*i, kTRUE);
170   }
171 }