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