1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.2.2.1 2002/10/14 09:45:57 hristov
19 Updating VirtualMC to v3-09-02
21 Revision 1.2 2002/10/01 08:47:04 jchudoba
22 Change loop order to run faster.
24 Revision 1.1 2002/09/17 08:37:12 jchudoba
25 Classes to create and store tracks maps - correcpondence between track label and entry number in the TreeH
29 ////////////////////////////////////////////////////////////////////////
33 // create an AliTrackMap - a relation between track label and
34 // it's index in TreeH. Check all branches with hits.
35 // Output is in the root file.
37 ////////////////////////////////////////////////////////////////////////
44 #include "TBenchmark.h"
45 #include "TStopwatch.h"
47 #include "AliDetector.h"
48 #include "AliTrackMapper.h"
49 #include "AliTrackMap.h"
53 ClassImp(AliTrackMapper)
55 ////////////////////////////////////////////////////////////////////////
56 void AliTrackMapper::CreateMap(Int_t nEvents, Int_t firstEventNr,
57 const char* fnMap, const char* fnHits)
60 // method to create a track map for a given number of events starting
61 // from event firstEventNr. This method just opens and closes files and
62 // loops over events, the creation of map is delegated to
63 // CreateMap(Int_t eventNr, TFile* fileMap) method
68 TFile *fileMap=TFile::Open(fnMap,"new");
69 if (!fileMap->IsOpen()) {cerr<<"Can't open output file "<<fnMap<<"!\n"; return;}
71 TFile *fileHits=TFile::Open(fnHits);
72 if (!fileHits->IsOpen()) {cerr<<"Can't open input file "<<fnHits<<"!\n"; return;}
73 if (!(gAlice=(AliRun*)fileHits->Get("gAlice"))) {
74 cerr<<"gAlice have not been found on galice.root !\n";
78 for (Int_t eventNr = firstEventNr; eventNr < firstEventNr+nEvents;
80 CreateMap(eventNr,fileMap);
81 } // end loop over events
90 if (fDEBUG > 0) timer.Print();
93 ////////////////////////////////////////////////////////////////////////
94 Int_t AliTrackMapper::CreateMap(Int_t eventNr, TFile* fileMap) {
96 // create an AliTrackMap for a given event
97 // correct gAlice must be already present in memory
100 Int_t nGenPrimPlusSecParticles = gAlice->GetEvent(eventNr);
101 if (fDEBUG > 1) cout<<"nGenPrimPlusSecParticles = "<<nGenPrimPlusSecParticles<<endl;
102 if (nGenPrimPlusSecParticles < 1) {
103 cerr<<"No primary particles found in event "<<eventNr<<endl;
107 TTree *treeK = gAlice->TreeK();
109 cerr<<"Error: Event "<<eventNr<<", treeK not found."<<endl;
112 Int_t nAllParticles = static_cast<Int_t>(treeK->GetEntries());
113 Int_t *trackMap = new Int_t[nAllParticles];
114 for (Int_t i = 0; i<nAllParticles; i++) {trackMap[i] = kNoEntry;}
116 TTree *treeH = gAlice->TreeH();
118 cerr<<"Error: Event "<<eventNr<<", treeH not found."<<endl;
121 Int_t treeHEntries = static_cast<Int_t>(treeH->GetEntries());
122 if (fDEBUG > 1) cout<<"treeHEntries "<<treeHEntries<<endl;
125 TObjArray *modules = gAlice->Detectors();
127 cerr<<"TObjArray with modules not found."<<endl;
130 Int_t nModules = static_cast<Int_t>(modules->GetEntries());
132 for (Int_t treeHIndex = 0; treeHIndex < treeHEntries; treeHIndex++) {
134 treeH->GetEvent(treeHIndex);
135 for (Int_t iModule = 0; iModule < nModules; iModule++) {
136 AliDetector * detector = dynamic_cast<AliDetector*>(modules->At(iModule));
137 if (!detector) continue;
138 // process only detectors with shunt = 0
139 if (detector->GetIshunt()) continue;
141 hit=(AliHit*)detector->FirstHit(-1);
142 Int_t lastLabel=-1, label;
143 for( ; hit; hit=(AliHit*)detector->NextHit() ) {
145 if (lastLabel != label) {
146 if (label < 0 || label >= nAllParticles) {
147 cerr<<"Error: label out of range. ";
148 cerr<<"Event "<<eventNr<<" treeHIndex "<<treeHIndex<<" label = "<<label<<endl;
151 if (trackMap[label] >=0 && trackMap[label] != treeHIndex) {
152 cerr<<"Error: different treeHIndex for label "<<label
153 <<" indeces: "<<trackMap[label]<<" != "<<treeHIndex;
154 cerr<<" event "<<eventNr<<" detector "<<detector->GetName()<<endl;
157 trackMap[label] = treeHIndex;
158 if (fDEBUG > 2) cout<<"treeHIndex, label = "<<treeHIndex<<" "<<label<<endl;
166 for (Int_t i = 0; i < nAllParticles; i++) {
167 cout<<eventNr<<"\t"<<i<<"\t"<<trackMap[i]<<endl;
171 AliTrackMap* trackMapObject = new AliTrackMap(nAllParticles, trackMap);
172 trackMapObject->SetEventNr(eventNr);
173 trackMapObject->Write();
174 delete trackMapObject;
180 ////////////////////////////////////////////////////////////////////////
181 AliTrackMap* AliTrackMapper::LoadTrackMap(Int_t eventNr, const char* fnMap, TFile* &fileMap) {
183 // read an AliTrackMap object for the given event eventNr from
186 fileMap=TFile::Open(fnMap);
187 if (!fileMap->IsOpen()) {cerr<<"Can't open file "<<fnMap<<" with map!\n"; return 0;}
189 sprintf(mapName,"AliTrackMap_%5.5d",eventNr);
190 AliTrackMap* trackMapObject = (AliTrackMap*)fileMap->Get(mapName);
191 if (!trackMapObject) {
192 cerr<<"Error: map named "<<mapName<<" not found."<<endl;
195 return trackMapObject;
197 ////////////////////////////////////////////////////////////////////////
198 void AliTrackMapper::CheckTrackMap(Int_t eventNr, const char* fnMap) {
203 AliTrackMap* trackMapObject = LoadTrackMap(eventNr, fnMap, fileMap);
204 if (!trackMapObject) return;
206 trackMapObject->PrintValues();
208 delete trackMapObject;
212 ////////////////////////////////////////////////////////////////////////