]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliTrackMapper.cxx
Default constructor added.
[u/mrichter/AliRoot.git] / STEER / AliTrackMapper.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 */
19
20 ////////////////////////////////////////////////////////////////////////
21 //
22 // AliTrackMapper.cxx
23 // description: 
24 //        create an AliTrackMap - a relation between track label and 
25 //        it's index in TreeH. Check all branches with hits.
26 //        Output is in the root file.
27 //
28 ////////////////////////////////////////////////////////////////////////
29
30 #include <iostream.h>
31
32 #include "TTree.h"
33 #include "TROOT.h"
34 #include "TFile.h"
35 #include "TBenchmark.h"
36 #include "TStopwatch.h"
37
38 #include "AliDetector.h"
39 #include "AliTrackMapper.h"
40 #include "AliTrackMap.h"
41 #include "AliRun.h"
42 #include "AliHit.h"
43
44 ClassImp(AliTrackMapper)
45
46 ////////////////////////////////////////////////////////////////////////
47 void AliTrackMapper::CreateMap(Int_t nEvents, Int_t firstEventNr,
48                             const char* fnMap, const char* fnHits)
49 {
50 //
51 // method to create a track map for a given number of events starting 
52 // from event firstEventNr. This method just opens and closes files and
53 // loops over events, the creation of map is delegated to 
54 // CreateMap(Int_t eventNr, TFile* fileMap) method
55 //
56   TStopwatch timer;
57   timer.Start();
58   
59   TFile *fileMap=TFile::Open(fnMap,"new");
60   if (!fileMap->IsOpen()) {cerr<<"Can't open output file "<<fnMap<<"!\n"; return;}
61
62   TFile *fileHits=TFile::Open(fnHits);
63   if (!fileHits->IsOpen()) {cerr<<"Can't open input file "<<fnHits<<"!\n"; return;}
64   if (!(gAlice=(AliRun*)fileHits->Get("gAlice"))) {
65     cerr<<"gAlice have not been found on galice.root !\n";
66     return;
67   }
68
69   for (Int_t eventNr = firstEventNr; eventNr < firstEventNr+nEvents;
70        eventNr++) {
71     CreateMap(eventNr,fileMap);
72   } // end loop over events
73
74   delete gAlice;
75   gAlice = 0;
76   fileHits->Close();
77   delete fileHits;
78   fileMap->Close();
79   delete fileMap;
80   timer.Stop();
81   if (fDEBUG > 0) timer.Print();
82 }
83
84 ////////////////////////////////////////////////////////////////////////
85 Int_t  AliTrackMapper::CreateMap(Int_t eventNr, TFile* fileMap) {
86 //
87 // create an AliTrackMap for a given event
88 // correct gAlice must be already present in memory
89 //
90   Int_t nGenPrimPlusSecParticles = gAlice->GetEvent(eventNr);
91   if (fDEBUG > 1) cout<<"nGenPrimPlusSecParticles = "<<nGenPrimPlusSecParticles<<endl;
92   if (nGenPrimPlusSecParticles < 1) {
93     cerr<<"No primary particles found in event "<<eventNr<<endl;
94     return -1;
95   }
96
97   TTree *treeK = gAlice->TreeK();
98   if (!treeK) {
99     cerr<<"Error: Event "<<eventNr<<", treeK not found."<<endl;
100     return -1;
101   }
102   Int_t nAllParticles = static_cast<Int_t>(treeK->GetEntries());
103   Int_t *trackMap = new Int_t[nAllParticles];
104   for (Int_t i = 0; i<nAllParticles; i++) {trackMap[i] = kNoEntry;}
105
106   TTree *treeH = gAlice->TreeH();
107   if (!treeH) {
108     cerr<<"Error: Event "<<eventNr<<", treeH not found."<<endl;
109     return -1;
110   }
111   Int_t treeHEntries = static_cast<Int_t>(treeH->GetEntries());
112   if (fDEBUG > 1) cout<<"treeHEntries "<<treeHEntries<<endl;
113
114
115   TObjArray *modules = gAlice->Detectors();
116   if (!modules) {
117     cerr<<"TObjArray with modules not found."<<endl;
118     return -1;
119   }
120   Int_t nModules = static_cast<Int_t>(modules->GetEntries());
121   for (Int_t iModule = 0; iModule < nModules; iModule++) {
122 //  for (Int_t iModule = nModules-1; iModule >= 0; iModule--) {
123     AliDetector * detector = dynamic_cast<AliDetector*>(modules->At(iModule));
124     if (!detector) continue;
125 // process only detectors with shunt = 0
126     if (detector->GetIshunt()) continue; 
127     if (fDEBUG > 0) cerr<<"Processing detector "<<detector->GetName()<<endl;
128
129     AliHit* hit;
130     for (Int_t treeHIndex = 0; treeHIndex < treeHEntries; treeHIndex++) {
131       detector->ResetHits();
132       treeH->GetEvent(treeHIndex);
133       hit=(AliHit*)detector->FirstHit(-1);
134       Int_t lastLabel=-1, label;
135       for( ; hit; hit=(AliHit*)detector->NextHit() ) {
136         label=hit->Track();     
137         if (lastLabel != label) {
138           if (label < 0 || label >= nAllParticles) {
139             cerr<<"Error: label out of range. ";
140             cerr<<"Event "<<eventNr<<" treeHIndex "<<treeHIndex<<" label = "<<label<<endl;
141             return -2;
142           }
143           if (trackMap[label] >=0 && trackMap[label] != treeHIndex) {
144             cerr<<"Error: different treeHIndex for label "<<label
145                 <<" indeces: "<<trackMap[label]<<" != "<<treeHIndex;
146             cerr<<" event "<<eventNr<<" detector "<<detector->GetName()<<endl;
147             return -3;
148           }
149           trackMap[label] = treeHIndex;
150           if (fDEBUG > 2) cout<<"treeHIndex, label = "<<treeHIndex<<" "<<label<<endl;
151           lastLabel = label;
152         }
153       }
154     }
155   }
156
157   if (fDEBUG > 0) {
158     for (Int_t i = 0; i < nAllParticles; i++) {
159       cout<<eventNr<<"\t"<<i<<"\t"<<trackMap[i]<<endl;
160     }
161   }
162   fileMap->cd();
163   AliTrackMap* trackMapObject = new AliTrackMap(nAllParticles, trackMap);
164   trackMapObject->SetEventNr(eventNr);
165   trackMapObject->Write();
166   delete trackMapObject;
167
168   delete [] trackMap;
169   return 0;
170 }
171   
172 ////////////////////////////////////////////////////////////////////////
173 AliTrackMap* AliTrackMapper::LoadTrackMap(Int_t eventNr, const char* fnMap, TFile* &fileMap) {
174 //
175 // read an AliTrackMap object for the given event eventNr from 
176 // the file fileMap
177 //
178   fileMap=TFile::Open(fnMap);
179   if (!fileMap->IsOpen()) {cerr<<"Can't open file "<<fnMap<<" with map!\n"; return 0;}
180   char mapName[20];
181   sprintf(mapName,"AliTrackMap_%5.5d",eventNr);
182   AliTrackMap* trackMapObject = (AliTrackMap*)fileMap->Get(mapName);
183   if (!trackMapObject) {
184     cerr<<"Error: map named "<<mapName<<" not found."<<endl;
185     return 0;
186   }
187   return trackMapObject;
188 }
189 ////////////////////////////////////////////////////////////////////////
190 void AliTrackMapper::CheckTrackMap(Int_t eventNr, const char* fnMap) {
191 //
192 // 
193 //
194   TFile *fileMap;
195   AliTrackMap* trackMapObject = LoadTrackMap(eventNr, fnMap, fileMap);
196   if (!trackMapObject) return;
197
198   trackMapObject->PrintValues();
199
200   delete trackMapObject;
201   fileMap->Close();
202   delete fileMap;
203 }
204 ////////////////////////////////////////////////////////////////////////
205