1 ////////////////////////////////////////////////////////////////////////
3 // name: AliHits2SDigits
5 // last update: 4.4.2002
6 // author: Jiri Chudoba
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
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
24 // 04.04.02 - first version
26 ////////////////////////////////////////////////////////////////////////
28 #if !defined(__CINT__) || defined(__MAKECINT__)
32 #include "TDirectory.h"
35 #include "TParticle.h"
36 #include "TPC/AliTPCDigitsArray.h"
37 #include "AliHeader.h"
38 #include "TGeometry.h"
39 #include "TObjArray.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"
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);
58 Bool_t gSameFiles = kFALSE;
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,
72 fileSDigits = Init(fileNameSDigits, fileNameHits);
73 if (!fileSDigits) return 1;
75 AliCopy(gFileHits,fileSDigits);
82 ITS = (AliITS*) gAlice->GetModule("ITS");
85 cerr<<"AliITS object not found on file." << endl;
86 } else if (!ITS->GetITSgeom()) {
87 cerr<<"AliITSgeom not found." << endl;
95 TPC = (AliTPC*)gAlice->GetDetector("TPC");
98 cerr<<"AliTPC object not found"<<endl;
103 AliTRDdigitizer *sdTRD;
105 sdTRD = InitTRDdigitizer();
110 AliPHOSSDigitizer *sdPHOS;
112 sdPHOS = new AliPHOSSDigitizer(fileNameHits.Data());
122 for (Int_t iEvent = firstEvent;iEvent<firstEvent+nEvents;iEvent++){
123 gAlice->GetEvent(iEvent);
124 gAlice->MakeTree("S",fileSDigits);
128 if (gDEBUG) {cout<<" Create ITS sdigits: ";}
129 ITS->MakeBranch("S");
130 ITS->SetTreeAddress();
132 if (gDEBUG) {cout<<"done"<<endl;}
137 if (gDEBUG) {cout<<" Create TPC sdigits: ";}
138 TPC->SetActiveSectors(1);
139 TPC->Hits2SDigits2(iEvent);
140 if (gDEBUG) {cout<<"done"<<endl;}
145 if (gDEBUG) {cout<<" Create TRD sdigits: ";}
146 sdTRD->InitOutput(fileSDigits, iEvent);
148 sdTRD->WriteDigits();
149 if (gDEBUG) {cout<<"done"<<endl;}
152 } // end of loop over events
154 // PHOS processes always all events
156 sdPHOS->ExecuteTask("deb all");
167 sdTRD->GetParameter()->Write();
171 fileSDigits->Close();
181 ////////////////////////////////////////////////////////////////////////
182 TFile* Init(TString fileNameSDigits, TString fileNameHits) {
183 // open input file, read in gAlice, prepare output file
184 if (gAlice) delete gAlice;
187 Bool_t gSameFiles = kFALSE;
188 if (fileNameSDigits == fileNameHits || fileNameSDigits == "") gSameFiles = kTRUE;
190 TString fileMode = "read";
191 if (gSameFiles) fileMode = "update";
193 gFileHits = TFile::Open(fileNameHits.Data(),fileMode.Data());
194 if (!gFileHits->IsOpen()) {
195 cerr<<"Can't open "<<fileNameHits.Data()<<" !\n";
198 if (!ImportgAlice(gFileHits)) return 0;
199 if (!gSameFiles) return gAlice->InitTreeFile("S",fileNameSDigits.Data());
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";
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;
222 ////////////////////////////////////////////////////////////////////////
223 AliTRDdigitizer *InitTRDdigitizer() {
224 // initialization of TRD digitizer
225 AliTRDdigitizer *sdTRD = new AliTRDdigitizer("TRDdigitizer"
226 ,"TRD digitizer class");
228 sdTRD->SetSDigits(kTRUE);
229 AliTRDparameter *TRDparam = new AliTRDparameter("TRDparameter"
230 ,"TRD parameter class");
232 sdTRD->SetParameter(TRDparam);
233 sdTRD->InitDetector();
236 ////////////////////////////////////////////////////////////////////////
237 void AliCopy(TFile *inputFile, TFile *outputFile) {
241 if (gDEBUG) cout<<"Copy gAlice: ";
244 if (gDEBUG) cout<<"done"<<endl;
246 TTree *treeE = gAlice->TreeE();
248 cerr<<"No TreeE found "<<endl;
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();
259 if (gDEBUG) cout<<"done"<<endl;
262 if (gDEBUG) cout<<"Copy AliceGeom: ";
263 TGeometry *AliceGeom = static_cast<TGeometry*>(inputFile->Get("AliceGeom"));
265 cerr<<"AliceGeom was not found in the input file "<<endl;
269 if (gDEBUG) cout<<"done"<<endl;