]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliHits2SDigits.C
Introducing Header instead of Log
[u/mrichter/AliRoot.git] / STEER / AliHits2SDigits.C
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // name: AliHits2SDigits
4 // date: 4.4.2002
5 // last update: 4.4.2002
6 // author: Jiri Chudoba
7 // version: 1.0
8 //
9 // description: 
10 //       creates sdigits for several detectors
11 //       stores sdigits in separate file (or in the source file
12 //       with hits). Stores gAlice object and copies TE to the
13 //       file with sdigits
14 //
15 // input:
16 //       TString fileNameSDigits ... output file with sdigits
17 //       TString fileNameHits ... input file with hits
18 //       Int_t nEvents  ... how many events to proceed
19 //       Int_t firstEvent  ... first event number
20 //       Int_t ITS, TPC, ...   many flags for diff. detectors
21 //
22 // History:
23 //
24 // 04.04.02 - first version
25 // 
26 ////////////////////////////////////////////////////////////////////////
27
28 #if !defined(__CINT__) || defined(__MAKECINT__)
29 #include "iostream.h"
30 #include "TTree.h"
31 #include "TBranch.h"
32 #include "TDirectory.h"
33 #include "TFile.h"
34 #include "AliRun.h"
35 #include "TParticle.h"
36 #include "TPC/AliTPCDigitsArray.h"
37 #include "AliHeader.h"
38 #include "TGeometry.h"
39 #include "TObjArray.h"
40 #include "TString.h"
41 #include "ITS/AliITS.h"
42 #include "TPC/AliTPC.h"
43 #include "PHOS/AliPHOSSDigitizer.h"
44 #include "TRD/AliTRDdigitizer.h"
45 #include "TStopwatch.h"
46 #include "TRD/AliTRDparameter.h"
47 #include "TOF/AliTOFSDigitizer.h"
48 #endif
49
50 TFile* Init(TString fileNameSDigits, TString fileNameHits);
51 TFile* OpenFile(TString fileName);
52 Bool_t ImportgAlice(TFile *file);
53 AliTRDdigitizer *InitTRDdigitizer();
54 void AliCopy(TFile *inputFile, TFile *outputFile);
55
56 // global variables
57
58 TFile *gFileHits = 0;
59 Bool_t gSameFiles = kFALSE;
60 Int_t gDEBUG = 1;
61
62
63 Int_t AliHits2SDigits(TString fileNameSDigits="sdigits.root", 
64                       TString fileNameHits="rfio:galice.root", 
65                       Int_t nEvents = 1, Int_t firstEvent = 0, Int_t iITS = 0,
66                       Int_t iTPC = 0, Int_t iTRD = 0,Int_t iPHOS = 0, 
67                       Int_t iTOF = 0, Int_t iCopy = 1)
68 {
69 //
70 // Initialization
71 //
72   TFile *fileSDigits;
73   fileSDigits = Init(fileNameSDigits, fileNameHits);
74   if (!fileSDigits) return 1;
75   if (iCopy) {
76     AliCopy(gFileHits,fileSDigits);
77     gFileHits->cd();
78   }  
79
80 // ITS
81   AliITS *ITS;
82   if (iITS) {
83     ITS  = (AliITS*) gAlice->GetModule("ITS");
84     if (!ITS) {
85       iITS = 0;
86       cerr<<"AliITS object not found on file." << endl;
87     } else if (!ITS->GetITSgeom()) {
88       iITS = 0;
89       cerr<<"AliITSgeom not found." << endl;
90     }
91   }
92
93 // TPC
94   AliTPC *TPC;
95   if (iTPC) {
96     TPC = (AliTPC*)gAlice->GetDetector("TPC");
97     if (!TPC) {
98       iTPC = 0;
99       cerr<<"AliTPC object not found"<<endl;
100     }
101   }
102
103 // TRD
104   AliTRDdigitizer *sdTRD;
105   if (iTRD) {
106     sdTRD = InitTRDdigitizer();
107   }
108
109
110 // PHOS
111   AliPHOSSDigitizer *sdPHOS;
112   if (iPHOS) {
113     sdPHOS = new AliPHOSSDigitizer(fileNameHits.Data());
114   }
115
116 // TOF
117   AliTOFSDigitizer *sdTOF;
118   if (iTOF) {
119     sdTOF = new AliTOFSDigitizer(fileNameHits.Data(),firstEvent,nEvents);
120   }
121
122 //
123 // loop over events
124 //
125   TStopwatch timer;
126   timer.Start();
127   for (Int_t iEvent = firstEvent;iEvent<firstEvent+nEvents;iEvent++){
128     gAlice->GetEvent(iEvent);
129     if (!gAlice->TreeS()) gAlice->MakeTree("S",fileSDigits);
130     
131 // ITS
132     if (iITS) {
133       if (gDEBUG) {cout<<"  Create ITS sdigits: ";}
134       ITS->MakeBranch("S");
135       ITS->SetTreeAddress();
136       ITS->Hits2SDigits();
137       if (gDEBUG) {cout<<"done"<<endl;}
138     }
139
140 // TPC
141     if (iTPC) {
142       if (gDEBUG) {cout<<"  Create TPC sdigits: ";}
143       if (iTPC == 1) {
144 // process all sectors
145         TPC->SetActiveSectors(1);
146         if (gDEBUG) {cout<<"All TPC sectors set active."<<endl;}
147       } else if (iTPC == 2) {
148 // process only sectors with hits
149         TPC->SetActiveSectors(0);
150         if (gDEBUG) {
151           printf("\nActive sectors\n");
152           Int_t iActive = 0;
153           for (Int_t i=0;i<72;i++) {
154             if (TPC->IsSectorActive(i)) {
155               printf("%2d ",i);
156               iActive++;
157               if (iActive%10 == 0) printf("\n");
158             }
159           }
160           printf("\n");
161         }
162       }    
163       TPC->Hits2SDigits2(iEvent);
164       if (gDEBUG) {cout<<"done"<<endl;}
165     }
166
167 // TRD
168     if (iTRD) {
169       if (gDEBUG) {cout<<"  Create TRD sdigits: ";}
170       sdTRD->InitOutput(fileSDigits, iEvent);
171       sdTRD->MakeDigits();
172       sdTRD->WriteDigits();
173       if (gDEBUG) {cout<<"done"<<endl;}
174     }
175     
176   } // end of loop over events
177
178 // PHOS processes always all events
179   if (iPHOS) {
180     sdPHOS->ExecuteTask("deb all");
181   }
182
183 // TOF does its own loop
184   if (iTOF) {
185     sdTOF->Exec("");
186   }
187
188
189 //
190 // finish 
191 //
192   timer.Stop(); 
193   timer.Print();
194
195   if (iTRD) { 
196     fileSDigits->cd();
197     sdTRD->GetParameter()->Write();
198   }
199
200   gFileHits->cd();
201   delete gAlice;
202   gAlice = 0;
203   if (!gSameFiles) {
204     gFileHits->Close();
205     delete gFileHits;
206   }
207
208 }
209  
210
211 ////////////////////////////////////////////////////////////////////////
212 TFile* Init(TString fileNameSDigits, TString fileNameHits) {
213 // open input file, read in gAlice, prepare output file
214   if (gAlice) delete gAlice;
215   gAlice = 0;
216
217   Bool_t gSameFiles = kFALSE;
218   if (fileNameSDigits == fileNameHits || fileNameSDigits == "") gSameFiles = kTRUE;
219
220   TString fileMode = "read";
221   if (gSameFiles) fileMode = "update";
222
223   gFileHits =  TFile::Open(fileNameHits.Data(),fileMode.Data());
224   if (!gFileHits->IsOpen()) {
225     cerr<<"Can't open "<<fileNameHits.Data()<<" !\n";
226     return 0;
227   }
228   if (!ImportgAlice(gFileHits)) return 0;
229   if (!gSameFiles) return gAlice->InitTreeFile("S",fileNameSDigits.Data());
230   return gFileHits;
231
232 }
233
234 ////////////////////////////////////////////////////////////////////////
235 TFile* OpenFile(TString fileName) {
236 // open file fileName
237   TFile *file = TFile::Open(fileName.Data());
238   if (!file->IsOpen()) {
239     cerr<<"Can't open "<<fileName.Data()<<" !\n";
240     return 0;
241   }
242   return file;
243 }
244
245 ////////////////////////////////////////////////////////////////////////
246 Bool_t ImportgAlice(TFile *file) {
247 // read in gAlice object from the file
248   gAlice = (AliRun*)file->Get("gAlice");
249   if (!gAlice)  return kFALSE;
250   return kTRUE;
251 }
252 ////////////////////////////////////////////////////////////////////////
253 AliTRDdigitizer *InitTRDdigitizer() {
254 // initialization of TRD digitizer
255   AliTRDdigitizer *sdTRD = new AliTRDdigitizer("TRDdigitizer"
256                                                ,"TRD digitizer class");
257   sdTRD->SetDebug(0);
258   sdTRD->SetSDigits(kTRUE);
259   AliTRDparameter *TRDparam = new AliTRDparameter("TRDparameter"
260                                                   ,"TRD parameter class");
261
262   sdTRD->SetParameter(TRDparam);
263   sdTRD->InitDetector();
264   if (!sdTRD->MakeBranch()) {
265     cerr<<"Problems with TRD digitizer initialization."<<endl;
266   }
267   return sdTRD;
268 }
269 ////////////////////////////////////////////////////////////////////////
270 void AliCopy(TFile *inputFile, TFile *outputFile) {
271 // copy some objects
272
273 // copy gAlice
274   if (gDEBUG) cout<<"Copy gAlice: ";
275   outputFile->cd();
276   gAlice->Write();
277   if (gDEBUG) cout<<"done"<<endl;
278
279   TTree *treeE  = gAlice->TreeE();
280   if (!treeE) {
281     cerr<<"No TreeE found "<<endl;
282     return;
283   }      
284
285 // copy TreeE
286   if (gDEBUG) cout<<"Copy TreeE: ";
287   AliHeader *header = new AliHeader();
288   treeE->SetBranchAddress("Header", &header);
289   treeE->SetBranchStatus("*",1);
290   TTree *treeENew =  treeE->CloneTree();
291   treeENew->Write();
292   if (gDEBUG) cout<<"done"<<endl;
293
294 // copy AliceGeom
295   if (gDEBUG) cout<<"Copy AliceGeom: ";
296   TGeometry *AliceGeom = static_cast<TGeometry*>(inputFile->Get("AliceGeom"));
297   if (!AliceGeom) {
298     cerr<<"AliceGeom was not found in the input file "<<endl;
299     return;
300   }
301   AliceGeom->Write();
302   if (gDEBUG) cout<<"done"<<endl;
303
304 }