]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/EveBase/AliEveKineTools.cxx
From Antonin: switch to the right formula of pseudo-rapidity in detailed view display.
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveKineTools.cxx
CommitLineData
d810d0de 1// $Id$
2// Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
c63fcc1d 3
d810d0de 4/**************************************************************************
5 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
51346b82 7 * full copyright notice. *
d810d0de 8 **************************************************************************/
9
10#include "AliEveKineTools.h"
c63fcc1d 11
12#include <TObject.h>
13#include <TTree.h>
14#include <TBranchElement.h>
7704f096 15#include <TClonesArray.h>
16
c63fcc1d 17#include <AliStack.h>
18#include <AliTrackReference.h>
19
84aff7a4 20#include <TEveTrack.h>
21#include <TEveElement.h>
c63fcc1d 22
7704f096 23#include <map>
c63fcc1d 24
57ffa5fb 25//______________________________________________________________________________
d810d0de 26// AliEveKineTools
c63fcc1d 27//
d4369ead 28// Tools for import of kinematics.
40790e5b 29//
d810d0de 30
d810d0de 31ClassImp(AliEveKineTools)
c63fcc1d 32
40790e5b 33namespace {
34
2142305c 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;
40790e5b 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);
2142305c 49 map.insert(std::make_pair(track->GetLabel(), track));
40790e5b 50 if (recurse)
51 MapTracks(map, track, recurse);
52 ++i;
53 }
54 }
55}
56
57/**************************************************************************/
76461dad 58
d810d0de 59void AliEveKineTools::SetDaughterPathMarks(TEveElement* cont, AliStack* stack, Bool_t recurse)
c63fcc1d 60{
76461dad 61 // Import daughters birth points.
62
84aff7a4 63 TEveElement::List_i iter = cont->BeginChildren();
c63fcc1d 64 while(iter != cont->EndChildren())
65 {
51346b82 66 TEveTrack* track = dynamic_cast<TEveTrack*>(*iter);
c63fcc1d 67 TParticle* p = stack->Particle(track->GetLabel());
a15e6d7d 68 if (p->GetNDaughters())
69 {
c63fcc1d 70 Int_t d0 = p->GetDaughter(0), d1 = p->GetDaughter(1);
a15e6d7d 71 for(int d = d0; d > 0 && d <= d1; ++d)
51346b82 72 {
c63fcc1d 73 TParticle* dp = stack->Particle(d);
a15e6d7d 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()));
d4369ead 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());
c63fcc1d 80 }
804bb570 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
27aa96b0 92 if (recurse)
93 SetDaughterPathMarks(track, stack, recurse);
c63fcc1d 94 }
95 ++iter;
96 }
97}
98
40790e5b 99/**************************************************************************/
974767e6 100
d810d0de 101void AliEveKineTools::SetTrackReferences(TEveElement* cont, TTree* treeTR, Bool_t recurse)
974767e6 102{
40790e5b 103 // Set decay and track reference path-marks.
974767e6 104
a15e6d7d 105 static const TEveException kEH("AliEveKineTools::ImportPathMarks");
974767e6 106
40790e5b 107 TrackMap_t map;
108 MapTracks(map, cont, recurse);
51346b82 109
d4369ead 110 TClonesArray* arr = 0;
111 treeTR->SetBranchAddress("TrackReferences", &arr);
c63fcc1d 112
d4369ead 113 Int_t nTreeEntries = (Int_t) treeTR->GetEntries();
114
115 for (Int_t te = 0; te < nTreeEntries; ++te)
c63fcc1d 116 {
d4369ead 117 treeTR->GetEntry(te);
118
119 Int_t last_label = -1;
2142305c 120 std::pair<TrackMap_ci, TrackMap_ci> range;
d4369ead 121 Int_t nArrEntries = arr->GetEntriesFast();
122
123 // printf("tree-entry %d, n-arr-entries %d\n", te, nArrEntries);
c63fcc1d 124
d4369ead 125 for (Int_t ae = 0; ae < nArrEntries; ++ae)
c63fcc1d 126 {
d4369ead 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 {
2142305c 141 range = map.equal_range(label);
d4369ead 142 last_label = label;
143 }
1d74e2b0 144
2142305c 145 for (TrackMap_ci i = range.first; i != range.second; ++i)
c63fcc1d 146 {
2142305c 147 i->second->AddPathMark
d4369ead 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;
27aa96b0 156}
157
40790e5b 158/**************************************************************************/
159
160void AliEveKineTools::SortPathMarks(TEveElement* el, Bool_t recurse)
27aa96b0 161{
d4369ead 162 // Sort path-marks for track by time.
163 // If recurse is true, descends down all track children.
a15e6d7d 164
40790e5b 165 TEveTrack* track = dynamic_cast<TEveTrack*>(el);
d4369ead 166 if (track)
167 track->SortPathMarksByTime();
c63fcc1d 168
d4369ead 169 if (recurse)
a15e6d7d 170 {
d4369ead 171 for (TEveElement::List_i i = el->BeginChildren(); i != el->EndChildren(); ++i)
172 SortPathMarks(*i, kTRUE);
c63fcc1d 173 }
174}