]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveKineTools.cxx
Merge changes from branches/dev/EVE. This branch was following development in ROOT...
[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. Preliminary version.
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       }
72       if (recurse)
73         SetDaughterPathMarks(track, stack, recurse);
74     }
75     ++iter;
76   }
77 }
78
79 /**************************************************************************/
80
81 void AliEveKineTools::SetTrackReferences(TEveElement* cont, TTree* treeTR, Bool_t recurse)
82 {
83   // Set decay and track reference path-marks.
84
85   static const TEveException kEH("AliEveKineTools::ImportPathMarks");
86
87   TrackMap_t map;
88   MapTracks(map, cont, recurse);
89
90   Int_t nPrimaries = (Int_t) treeTR->GetEntries();
91   TIter next(treeTR->GetListOfBranches());
92   TBranchElement* el;
93   Bool_t isRef = kTRUE;
94
95   while ((el = (TBranchElement*) next()))
96   {
97     if (strcmp("AliRun",el->GetName()) == 0)
98       isRef = kFALSE;
99
100     TClonesArray* arr = 0;
101     el->SetAddress(&arr);
102     for (Int_t iPrimPart = 0; iPrimPart<nPrimaries; iPrimPart++)
103     {
104       el->GetEntry(iPrimPart);
105
106       Int_t last_label = -1;
107       TrackMap_t::iterator iter = map.end();
108       Int_t Nent =  arr->GetEntriesFast();
109       for (Int_t iTrackRef = 0; iTrackRef < Nent; iTrackRef++)
110       {
111         AliTrackReference* atr = (AliTrackReference*)arr->UncheckedAt(iTrackRef);
112
113         Int_t label = atr->GetTrack();
114         if (label < 0)
115           throw(kEH + Form("negative label for entry %d in branch %s.",
116                           iTrackRef, el->GetName()));
117
118         if (label != last_label)
119         {
120           iter       = map.find(label);
121           last_label = label;
122         }
123
124         if (iter != map.end())
125         {
126           iter->second->AddPathMark
127             (TEvePathMark(isRef ? TEvePathMark::kReference : TEvePathMark::kDecay,
128                           TEveVector(atr->X(),  atr->Y(),  atr->Z()),
129                           TEveVector(atr->Px(), atr->Py(), atr->Pz()),
130                           atr->GetTime()));
131         }
132       } // loop track refs
133     } // loop primaries in clones arrays
134     delete arr;
135   } // end loop through top branches
136 }
137
138 /**************************************************************************/
139
140 void AliEveKineTools::SortPathMarks(TEveElement* el, Bool_t recurse)
141 {
142   // Sort path-marks for all tracks by time.
143
144   // !!!!! MT ... this is one level deep only!
145
146   TEveTrack* track = dynamic_cast<TEveTrack*>(el);
147   if (track) track->SortPathMarksByTime();
148
149   TEveElement::List_i i = el->BeginChildren();
150   while (i != el->EndChildren() && recurse)
151   {
152     track = dynamic_cast<TEveTrack*>(el);
153     if (track) track->SortPathMarksByTime();
154     ++i;
155   }
156 }