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