]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliHits2SDigits.C
8a74e6693f14587cc5c3a9edff7b630b9abf9f15
[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 #endif
48
49 TFile* Init(TString fileNameSDigits, TString fileNameHits);
50 TFile* OpenFile(TString fileName);
51 Bool_t ImportgAlice(TFile *file);
52 AliTRDdigitizer *InitTRDdigitizer();
53 void AliCopy(TFile *inputFile, TFile *outputFile);
54
55 // global variables
56
57 TFile *gFileHits = 0;
58 Bool_t gSameFiles = kFALSE;
59 Int_t gDEBUG = 1;
60
61
62 Int_t AliHits2SDigits(TString fileNameSDigits="sdigits.root", 
63                       TString fileNameHits="rfio:galice.root", 
64                       Int_t nEvents = 1, Int_t firstEvent = 0, Int_t iITS = 0,
65                       Int_t iTPC = 0, Int_t iTRD = 0,Int_t iPHOS = 0, 
66                       Int_t iCopy = 1)
67 {
68 //
69 // Initialization
70 //
71   TFile *fileSDigits;
72   fileSDigits = Init(fileNameSDigits, fileNameHits);
73   if (!fileSDigits) return 1;
74   if (iCopy) {
75     AliCopy(gFileHits,fileSDigits);
76     gFileHits->cd();
77   }  
78
79 // ITS
80   AliITS *ITS;
81   if (iITS) {
82     ITS  = (AliITS*) gAlice->GetModule("ITS");
83     if (!ITS) {
84       iITS = 0;
85       cerr<<"AliITS object not found on file." << endl;
86     } else if (!ITS->GetITSgeom()) {
87       cerr<<"AliITSgeom not found." << endl;
88       iITS = 0;
89     }
90   }
91
92 // TPC
93   AliTPC *TPC;
94   if (iTPC) {
95     TPC = (AliTPC*)gAlice->GetDetector("TPC");
96     if (!TPC) {
97       iTPC = 0;
98       cerr<<"AliTPC object not found"<<endl;
99     }
100   }
101
102 // TRD
103   AliTRDdigitizer *sdTRD;
104   if (iTRD) {
105     sdTRD = InitTRDdigitizer();
106   }
107
108
109 // PHOS
110   AliPHOSSDigitizer *sdPHOS;
111   if (iPHOS) {
112     sdPHOS = new AliPHOSSDigitizer(fileNameHits.Data());
113   }
114
115
116
117 //
118 // loop over events
119 //
120   TStopwatch timer;
121   timer.Start();
122   for (Int_t iEvent = firstEvent;iEvent<firstEvent+nEvents;iEvent++){
123     gAlice->GetEvent(iEvent);
124     gAlice->MakeTree("S",fileSDigits);
125     
126 // ITS
127     if (iITS) {
128       if (gDEBUG) {cout<<"  Create ITS sdigits: ";}
129       ITS->MakeBranch("S");
130       ITS->SetTreeAddress();
131       ITS->Hits2SDigits();
132       if (gDEBUG) {cout<<"done"<<endl;}
133     }
134
135 // TPC
136     if (iTPC) {
137       if (gDEBUG) {cout<<"  Create TPC sdigits: ";}
138       TPC->SetActiveSectors(1);
139       TPC->Hits2SDigits2(iEvent);
140       if (gDEBUG) {cout<<"done"<<endl;}
141     }
142
143 // TRD
144     if (iTRD) {
145       if (gDEBUG) {cout<<"  Create TRD sdigits: ";}
146       sdTRD->InitOutput(fileSDigits, iEvent);
147       sdTRD->MakeDigits();
148       sdTRD->WriteDigits();
149       if (gDEBUG) {cout<<"done"<<endl;}
150     }
151     
152   } // end of loop over events
153
154 // PHOS processes always all events
155   if (iPHOS) {
156     sdPHOS->ExecuteTask("deb all");
157   }
158
159 //
160 // finish 
161 //
162   timer.Stop(); 
163   timer.Print();
164
165   if (iTRD) { 
166     fileSDigits->cd();
167     sdTRD->GetParameter()->Write();
168     gFileHits->cd();
169   }
170
171   fileSDigits->Close();
172   delete fileSDigits;
173   if (!gSameFiles) {
174     gFileHits->Close();
175     delete gFileHits;
176   }
177
178 }
179  
180
181 ////////////////////////////////////////////////////////////////////////
182 TFile* Init(TString fileNameSDigits, TString fileNameHits) {
183 // open input file, read in gAlice, prepare output file
184   if (gAlice) delete gAlice;
185   gAlice = 0;
186
187   Bool_t gSameFiles = kFALSE;
188   if (fileNameSDigits == fileNameHits || fileNameSDigits == "") gSameFiles = kTRUE;
189
190   TString fileMode = "read";
191   if (gSameFiles) fileMode = "update";
192
193   gFileHits =  TFile::Open(fileNameHits.Data(),fileMode.Data());
194   if (!gFileHits->IsOpen()) {
195     cerr<<"Can't open "<<fileNameHits.Data()<<" !\n";
196     return 0;
197   }
198   if (!ImportgAlice(gFileHits)) return 0;
199   if (!gSameFiles) return gAlice->InitTreeFile("S",fileNameSDigits.Data());
200   return gFileHits;
201
202 }
203
204 ////////////////////////////////////////////////////////////////////////
205 TFile* OpenFile(TString fileName) {
206 // open file fileName
207   TFile *file = TFile::Open(fileName.Data());
208   if (!file->IsOpen()) {
209     cerr<<"Can't open "<<fileName.Data()<<" !\n";
210     return 0;
211   }
212   return file;
213 }
214
215 ////////////////////////////////////////////////////////////////////////
216 Bool_t ImportgAlice(TFile *file) {
217 // read in gAlice object from the file
218   gAlice = (AliRun*)file->Get("gAlice");
219   if (!gAlice)  return kFALSE;
220   return kTRUE;
221 }
222 ////////////////////////////////////////////////////////////////////////
223 AliTRDdigitizer *InitTRDdigitizer() {
224 // initialization of TRD digitizer
225   AliTRDdigitizer *sdTRD = new AliTRDdigitizer("TRDdigitizer"
226                                                ,"TRD digitizer class");
227   sdTRD->SetDebug(0);
228   sdTRD->SetSDigits(kTRUE);
229   AliTRDparameter *TRDparam = new AliTRDparameter("TRDparameter"
230                                                   ,"TRD parameter class");
231
232   sdTRD->SetParameter(TRDparam);
233   sdTRD->InitDetector();
234   return sdTRD;
235 }
236 ////////////////////////////////////////////////////////////////////////
237 void AliCopy(TFile *inputFile, TFile *outputFile) {
238 // copy some objects
239
240 // copy gAlice
241   if (gDEBUG) cout<<"Copy gAlice: ";
242   outputFile->cd();
243   gAlice->Write();
244   if (gDEBUG) cout<<"done"<<endl;
245
246   TTree *treeE  = gAlice->TreeE();
247   if (!treeE) {
248     cerr<<"No TreeE found "<<endl;
249     return;
250   }      
251
252 // copy TreeE
253   if (gDEBUG) cout<<"Copy TreeE: ";
254   AliHeader *header = new AliHeader();
255   treeE->SetBranchAddress("Header", &header);
256   treeE->SetBranchStatus("*",1);
257   TTree *treeENew =  treeE->CloneTree();
258   treeENew->Write();
259   if (gDEBUG) cout<<"done"<<endl;
260
261 // copy AliceGeom
262   if (gDEBUG) cout<<"Copy AliceGeom: ";
263   TGeometry *AliceGeom = static_cast<TGeometry*>(inputFile->Get("AliceGeom"));
264   if (!AliceGeom) {
265     cerr<<"AliceGeom was not found in the input file "<<endl;
266     return;
267   }
268   AliceGeom->Write();
269   if (gDEBUG) cout<<"done"<<endl;
270
271 }