Detect possible particle decay by the process-id of its last daughter.
[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   typedef std::map<Int_t, TEveTrack*> TrackMap_t;
36
37   void MapTracks(TrackMap_t& map, TEveElement* cont, Bool_t recurse)
38   {
39     TEveElement::List_i i = cont->BeginChildren();
40     while (i != cont->EndChildren()) {
41       TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
42       map[track->GetLabel()] = track;
43       if (recurse)
44         MapTracks(map, track, recurse);
45       ++i;
46     }
47   }
48 }
49
50 /**************************************************************************/
51
52 void AliEveKineTools::SetDaughterPathMarks(TEveElement* cont, AliStack* stack, Bool_t recurse)
53 {
54   // Import daughters birth points.
55
56   TEveElement::List_i  iter = cont->BeginChildren();
57   while(iter != cont->EndChildren())
58   {
59     TEveTrack* track = dynamic_cast<TEveTrack*>(*iter);
60     TParticle* p = stack->Particle(track->GetLabel());
61     if (p->GetNDaughters())
62     {
63       Int_t d0 = p->GetDaughter(0), d1 = p->GetDaughter(1);
64       for(int d = d0; d > 0 && d <= d1; ++d)
65       {
66         TParticle* dp = stack->Particle(d);
67         track->AddPathMark(TEvePathMark(TEvePathMark::kDaughter,
68                                         TEveVector(dp->Vx(), dp->Vy(), dp->Vz()),
69                                         TEveVector(dp->Px(), dp->Py(), dp->Pz()),
70                                         dp->T()));
71         // printf("Daughter path-mark for %d, %d, t=%e, r=%f,%f,%f\n",
72         //        track->GetLabel(), d, dp->T(), dp->Vx(), dp->Vy(), dp->Vz());
73       }
74
75       // Check last process, set decay if needed.
76       Int_t lp = stack->Particle(d1)->GetUniqueID();
77       if (lp != kPBrem && lp != kPDeltaRay && lp < kPCerenkov)
78       {
79         TParticle* dp = stack->Particle(d1);
80         track->AddPathMark(TEvePathMark(TEvePathMark::kDecay,
81                                         TEveVector(dp->Vx(), dp->Vy(), dp->Vz()),
82                                         TEveVector(0, 0,0),  dp->T()));
83       }
84
85       if (recurse)
86         SetDaughterPathMarks(track, stack, recurse);
87     }
88     ++iter;
89   }
90 }
91
92 /**************************************************************************/
93
94 void AliEveKineTools::SetTrackReferences(TEveElement* cont, TTree* treeTR, Bool_t recurse)
95 {
96   // Set decay and track reference path-marks.
97
98   static const TEveException kEH("AliEveKineTools::ImportPathMarks");
99
100   TrackMap_t map;
101   MapTracks(map, cont, recurse);
102
103   TClonesArray* arr = 0;
104   treeTR->SetBranchAddress("TrackReferences", &arr);
105
106   Int_t nTreeEntries = (Int_t) treeTR->GetEntries();
107
108   for (Int_t te = 0; te < nTreeEntries; ++te)
109   {
110     treeTR->GetEntry(te);
111
112     Int_t last_label = -1;
113     TrackMap_t::iterator iter = map.end();
114     Int_t nArrEntries = arr->GetEntriesFast();
115
116     // printf("tree-entry %d, n-arr-entries %d\n", te, nArrEntries);
117
118     for (Int_t ae = 0; ae < nArrEntries; ++ae)
119     {
120       AliTrackReference* atr = (AliTrackReference*)arr->UncheckedAt(ae);
121       Bool_t isRef = (atr->DetectorId() != -1);
122       Int_t  label = atr->GetTrack();
123
124       // printf("    arr-entry %d, label %d, detid %d, len=%f, t=%e r=%f,%f,%f\n",
125       //        ae, label, atr->DetectorId(), 
126       //        atr->GetLength(), atr->GetTime(), atr->X(), atr->Y(), atr->Z());
127
128       if (label < 0)
129         throw(kEH + Form("negative label for array-entry %d in tree-entry %d.",
130                          ae, te));
131
132       if (label != last_label)
133       {
134         iter       = map.find(label);
135         last_label = label;
136       }
137
138       if (iter != map.end())
139       {
140         iter->second->AddPathMark
141           (TEvePathMark(isRef ? TEvePathMark::kReference : TEvePathMark::kDecay,
142                         TEveVector(atr->X(),  atr->Y(),  atr->Z()),
143                         TEveVector(atr->Px(), atr->Py(), atr->Pz()),
144                         atr->GetTime()));
145       }
146     }
147   }
148   delete arr;
149 }
150
151 /**************************************************************************/
152
153 void AliEveKineTools::SortPathMarks(TEveElement* el, Bool_t recurse)
154 {
155   // Sort path-marks for track by time.
156   // If recurse is true, descends down all track children.
157
158   TEveTrack* track = dynamic_cast<TEveTrack*>(el);
159   if (track)
160     track->SortPathMarksByTime();
161
162   if (recurse)
163   {
164     for (TEveElement::List_i i = el->BeginChildren(); i != el->EndChildren(); ++i)
165       SortPathMarks(*i, kTRUE);
166   }
167 }