]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSHits2SDR.C
Removing obsolete macros
[u/mrichter/AliRoot.git] / ITS / AliITSHits2SDR.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <iostream.h>
3 #include <TClonesArray.h>
4 #include <TFile.h>
5 #include <TParticle.h>
6 #include <TStopwatch.h>
7 #include <TTree.h>
8 #include <TFile.h>
9
10 #include "AliRun.h"
11 #include "AliHeader.h"
12 #include "AliITS.h"
13 #include "AliITSDetType.h"
14 #include "AliITSsimulationFastPoints.h"
15 #include "AliITSresponse.h"
16 #include "AliITSresponseSDD.h"
17 #include "AliITSreconstruction.h"
18 #include "AliRunDigitizer.h"
19 #include "AliITSsegmentationSDD.h"
20 #include "AliITSDigitizer.h"
21
22 #endif
23
24 TFile* AccessInpFile(TString inFile);
25 void writeAR(TFile * fin, TFile *fou);
26 void AliITSSD2D(TString inFile = "galice.root", TString outFile = "");
27 void AliITSH2FR2files (TString inFile, TString outfile);
28 void AliITSH2SD2Files(TFile *file2);
29 void copyTK(TString inFile, TString outFile);
30
31 void AliITSHits2SDR(TString inFile = "galice.root", TString File1 = "galiceS.root",TString File2 = "galiceD.root", TString File3 = "galiceR.root"){
32   ////////////////////////////////////////////////////////////////////
33   // This macro takes hits from inFile file and 
34   // produce fast RecPoints, Sdigits, Digits and slow Rec Points 
35   // for the ITS using the standard detector simulations. 
36   // The first argument is the input file which contains the tree with hits
37   // The second argument is the output file name onto which TreeS will  
38   // be written. 
39   // The third and fourth arguments are the file names onto which TreeD  
40   // and TreeR+TreeC will be written (respectively).
41   // The AliRun object and the KINE tree will be saved onto 
42   // these files as well. 
43   // This macro will process all of the events on input the root file.
44   // This macro can be compiled and it should be launched from an aliroot 
45   // session.
46   // WARNING: Fast Rec Points are written on a branch named ITSRecPointsF
47   // this macro should be used with aliroot (it does not load the libraries) 
48   ////////////////////////////////////////////////////////////////////
49  
50  
51  
52   // Produce Fast Rec Points  (on File3)
53   cout<<"Producing Fast rec points...\n";
54   AliITSH2FR2files(inFile,File3);
55
56   // TreeK is copied to File3 (for comparison purposes)
57  
58   copyTK(inFile,File3);
59  
60   // Produce SDigits
61   TFile *file = AccessInpFile(inFile);
62   if(!file){
63     cerr<<"Problems in accessing the input file\n";
64     return;
65   }
66   TFile *file2 = gAlice->InitTreeFile("S",File1);
67   file2->cd();
68   cout<<"Producing Summable digits ....\n";
69   AliITSH2SD2Files(file2);
70   // write the AliRun object to the output file
71   writeAR(file,file2);
72   delete gAlice;
73   gAlice = 0;
74   file->Close();
75   delete file;
76   file2 = 0;
77   
78   // Produce Digits 
79   cout<<"Producing Digits ...\n";
80   AliITSSD2D(File1,File2);
81
82   // Produce Slow Rec. Points 
83   TFile *fin = new TFile(File2);
84   gAlice = (AliRun*) fin->Get("gAlice");
85   const char *nulptr=0;
86   AliITSreconstruction *itsr = new AliITSreconstruction(nulptr);
87   itsr->SetOutputFile(File3);
88   itsr->Init();
89   itsr->Exec();
90   TFile *fou = gAlice->GetTreeRFile();
91   writeAR(fin,fou);
92   delete itsr;
93   
94 }
95
96 void AliITSH2SD2Files(TFile *file2){
97
98   // Function used to produce s-digits
99  
100   AliITS *ITS  = (AliITS*) gAlice->GetModule("ITS");
101   if (!ITS) {
102         cerr<<"AliITS object not found on file" << endl;
103         return;
104   }
105   if(!(ITS->GetITSgeom())){
106         cerr << " AliITSgeom not found. Can't digitize without it." << endl;
107         return ;
108   }
109
110   // these re-instantiations are necessary if the input file was created
111   // with aliroot version V3.07 or older
112   AliITSresponseSDD *resp1 = (AliITSresponseSDD*)ITS->DetType(1)->GetResponseModel();
113   resp1 = new AliITSresponseSDD();
114   ITS->SetResponseModel(1,resp1);
115   cout << "Changed response class for SDD: \n";
116   resp1->Print();
117
118   for (Int_t nevent=0; nevent<gAlice->TreeE()->GetEntries(); nevent++) {
119     gAlice->GetEvent(nevent);
120     gAlice->MakeTree("S",file2);
121     ITS->MakeBranch("S");
122     ITS->SetTreeAddress();   
123     ITS->Hits2SDigits();
124   }
125 }
126
127
128 TFile * AccessInpFile(TString inFile){
129
130   // Function used to open the input file and fetch the AliRun object
131
132   TFile *retfil = 0;
133   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(inFile);
134   if (file) {file->Close(); delete file;}
135   file = new TFile(inFile);
136   if (!file->IsOpen()) {
137         cerr<<"Can't open "<<inFile<<" !" << endl;
138         return retfil;
139   } 
140
141   // Get AliRun object from file or return if not on file
142   if (gAlice) delete gAlice; gAlice = 0;
143   gAlice = (AliRun*)file->Get("gAlice");
144   if (!gAlice) {
145         cerr << "AliRun object not found on file"<< endl;
146         return retfil;
147   } 
148   return file;
149 }
150
151 void writeAR(TFile * fin, TFile *fou) {
152   TDirectory *current = gDirectory;
153   TTree *Te;
154   TTree *TeNew;
155   AliHeader *alhe = new AliHeader();
156   Te = (TTree*)fin->Get("TE");
157   Te->SetBranchAddress("Header",&alhe);
158   Te->SetBranchStatus("*",1);
159   fou->cd();
160   TeNew = Te->CloneTree();
161   TeNew->Write(0,TObject::kOverwrite);
162   gAlice->Write(0,TObject::kOverwrite);
163   current->cd();
164   delete alhe;
165   cout<<"AliRun object written to file\n";
166 }
167
168
169 void AliITSSD2D(TString inFile, TString outFile){
170   AliRunDigitizer * manager = new AliRunDigitizer(1,1);
171   char ftmp[50];
172   sprintf(ftmp,"%s",inFile.Data());
173   manager->SetInputStream(0,ftmp);
174   if(outFile != "")manager->SetOutputFile(outFile);
175   AliITSDigitizer *dITS  = new AliITSDigitizer(manager);
176   manager->Exec("");
177   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(inFile);
178   TFile *file2= new TFile(outFile,"update");
179   writeAR(file,file2);
180   delete dITS;
181   delete manager;
182   file = (TFile*)gROOT->GetListOfFiles()->FindObject(inFile);
183   if(file){
184     file->Close();
185     delete file;
186   }
187   file2->Close();
188   delete file2;
189 }
190
191 void AliITSH2FR2files (TString inFile, TString outfile)
192 {
193
194   TFile * file = AccessInpFile(inFile);
195   if(!file){
196     cerr<<"*** AliITSH2FR2files: Problems in accessing the input file\n";
197     return;
198   }
199  
200   // open output file and create TreeR on it
201   TFile * file2 = gAlice->InitTreeFile("R",outfile);
202   file2->cd();
203   // get ITS  
204   AliITS *ITS  = (AliITS*) gAlice->GetModule("ITS");
205   if (!ITS) {
206         cerr<<"ITSH2FR2files.C : AliITS object not found on file" << endl;
207         return;
208   }  // end if !ITS
209
210   // Set the simulation model
211  
212   for (Int_t i=0;i<3;i++) {
213     ITS->SetSimulationModel(i,new AliITSsimulationFastPoints());
214   }
215  
216
217   //
218   // Event Loop
219   //
220   Int_t nsignal=25;
221   Int_t size=-1;
222   TStopwatch timer;
223   for (Int_t nevent=0; nevent<gAlice->TreeE()->GetEntries(); nevent++) {
224     gAlice->GetEvent(nevent);
225     if(!gAlice->TreeR())gAlice->MakeTree("R",file2);
226     ITS->MakeBranch("RF"); 
227     ITS->SetTreeAddress();  
228     Int_t bgr_ev=nevent/nsignal;
229     ITS->HitsToFastRecPoints(nevent,bgr_ev,size," ","All"," ");
230   }
231
232   timer.Stop(); timer.Print();
233
234   delete gAlice;
235   gAlice = 0;
236   file->Close();
237 }
238
239 void copyTK(TString inFile, TString outFile) {
240
241   TDirectory *current = gDirectory;
242   TParticle *part = new TParticle();
243   TTree *Tk;
244   TTree *TkNew;
245   TFile *fin = AccessInpFile(inFile);
246   TFile *fou = new TFile(outFile,"update");
247   char tname[30];
248   for(Int_t event= 0; event < gAlice->TreeE()->GetEntries(); event++){
249     current->cd();
250     sprintf(tname,"TreeK%d",event);
251     Tk = (TTree*)fin->Get(tname);
252     if(!Tk){
253       cerr<<"Trying to access a non-existent TTree : "<<tname<<endl;
254     }
255     else {
256       Tk->SetBranchAddress("Particles",&part);
257       Tk->SetBranchStatus("*",1);
258       fou->cd();
259       TkNew = Tk->CloneTree();
260       TkNew->Write();
261     }
262   }
263   delete gAlice;
264   gAlice = 0;
265   fin->Close();
266   fou->Close();
267   delete fin;
268   delete fou;
269   current->cd();
270
271 }