]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSDDLRawData.C
Added option TGeo in the constructor in order to initialize some parameters directly...
[u/mrichter/AliRoot.git] / ITS / AliITSDDLRawData.C
1 #if !defined(__CINT__)
2 #include <Riostream.h>
3 #include "AliITSDDLRawData.h"
4 #include "AliRunLoader.h"
5 #include "AliLoader.h"
6 #include "AliITS.h"
7 #endif
8
9 /*
10 Before running this macro it is necessary to comment the following line of the method
11 AddDigit in the class AliITSsimulationSDD
12 //if( fResponse->Do10to8() ) signal = Convert8to10( signal ); 
13 In this way the amplitude value for signal coming from SDD takes only 8 bits and not 10.
14 */
15 //DigitsFile is the input file that contains digits
16
17 void AliITSDDLRawData(Int_t eventNumber=0){
18
19   Int_t spdLDCs=2;
20   Int_t sddLDCs=4;
21   Int_t ssdLDCs=2;
22   Int_t eventn=0;
23   const char * inFile_new = "galice.root";
24   AliRunLoader *rl = AliRunLoader::Open(inFile_new,"Event","read");
25   rl->LoadgAlice();
26   gAlice=rl->GetAliRun();
27   Int_t nevents=rl->GetNumberOfEvents();
28   cout<<"Number of Events:"<<nevents<<endl;
29   while (eventNumber<=0 || eventNumber>nevents){
30     cout<<"Insert the event number:";
31     cin>>eventNumber;
32     cout<<endl;
33   }
34   rl->GetEvent(eventNumber-1);
35   AliLoader *itsloader=rl->GetLoader("ITSLoader");
36   itsloader->LoadDigits();
37   TTree *TD=itsloader->TreeD();
38   gAlice=rl->GetAliRun(); 
39   if(!gAlice){
40     cout<<"gAlice is null"<<endl;
41     return;
42   } 
43   AliITS *ITS  = (AliITS*)gAlice->GetDetector("ITS");
44
45   Int_t nmodules;
46   ITS->InitModules(-1,nmodules);
47   ITS->GetDetTypeSim()->SetTreeAddressD(TD,"ITS");
48
49     AliITSDDLRawData *util=new AliITSDDLRawData();
50     //Verbose level
51     // 0: Silent
52     // 1: cout messages
53     // 2: txt files with digits 
54     //BE CAREFUL, verbose level 2 MUST be used only for debugging and
55     //it is highly suggested to use this mode only for debugging digits files
56     //reasonably small, because otherwise the size of the txt files can reach
57     //quickly several MB wasting time and disk space.
58     util->SetVerbose(0);
59     TStopwatch timer;
60     
61     //SILICON PIXEL DETECTOR
62     cout<<"Formatting data for SPD"<<endl;
63     timer.Start();
64     util->RawDataSPD(ITS,TD,spdLDCs,eventNumber);
65     timer.Stop();
66     timer.Print();
67     //ONLY FOR DEBUGGING 
68     //    util->TestFormat(eventNumber);
69     
70     //SILICON DRIFT DETECTOR
71     cout<<"Formatting data for SDD"<<endl;
72     timer.Start();
73     util->RawDataSDD(ITS,TD,sddLDCs,eventNumber);
74     timer.Stop();
75     timer.Print();
76     
77     //SILICON STRIP DETECTOR
78     cout<<"Formatting data for SSD"<<endl;
79     timer.Start();
80     util->RawDataSSD(ITS,TD,ssdLDCs,eventNumber);
81     timer.Stop();
82     timer.Print();
83     
84     delete util;
85
86     return;
87 }