]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFSDigits2Digits.C
Fixed a bug related to xmin and xmax values for probability histos (F.Pierella)
[u/mrichter/AliRoot.git] / TOF / AliTOFSDigits2Digits.C
CommitLineData
ca4508a6 1Int_t AliTOFSDigits2Digits(TString infileNameSDigits, Int_t firstEvent=0,Int_t nEvents=1, Int_t ndump=15)
2{
3 // Create TOF digits out of TOF sdigits (no merging implemented here).
4 // Input: infileNameSDigits --> TOF sdigits filename
5 // firstEvent --> first event to digitize
6 // nEvents --> number of events to digitize
7 // including the first and the last ones
8 // if nEvents==0 we sdigitize all events in infileNameSDigits
9 // ndump --> number of dumped digits
10
11 // Author: F. Pierella (Bologna University)
12 // report problem to pierella@bo.infn.it
13
14 // Use case: (start root)
15 //root [0] .L AliTOFSDigits2Digits.C
16 //root [1] AliTOFSDigits2Digits("fileWithTOFSdigits.root")
17
18
19 // number 3 is a legacy from AliDigit object
20 const Int_t kMAXDIGITS = 3;
21
22
23 // Dynamically link some shared libs
24 if (gClassTable->GetID("AliRun") < 0) {
25 gROOT->LoadMacro("loadlibs.C");
26 loadlibs();
27 } else {
28 delete gAlice;
29 gAlice = 0;
30 }
31
32 Int_t rc=0;
33
34
35 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(infileNameSDigits.Data());
36 if(file){
37 cout<<"headerFile already open \n";
38 }
39 else {
40 // file is open in "update" mode
41 // in order to have a writable digits tree
42 if(!file)file=TFile::Open(infileNameSDigits.Data(),"update");
43 }
44
45 // Get AliRun object from file
46 if (!gAlice) {
47 gAlice = (AliRun*)file->Get("gAlice");
48 if (gAlice) printf("AliRun object found on file\n");
49 }
50
51 if (nEvents == 0) nEvents = (Int_t) gAlice->TreeE()->GetEntries();
52
53 AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
54
55 if (!tof) {
56 cout << "<AliTOFSDigits2Digits> No TOF detector found" << endl;
57 rc = 2;
58 return rc;
59 }
60
61 Int_t totndig=0; // total number of digits
62 Int_t tottracks=0; // total number of tracks contributing to totndig
63 Int_t upperBoundEvNum=firstEvent+nEvents;
64 for (Int_t ievent = firstEvent; ievent < upperBoundEvNum; ievent++) {
65
66 gAlice->GetEvent(ievent);
67 if (gAlice->TreeD () == 0)
68 gAlice->MakeTree ("D");
69
70 //Make branches
71 char branchname[20];
72 sprintf (branchname, "%s", tof->GetName ());
73 //Make branch for digits
74 tof->MakeBranch ("D");
75
76
77 // get the TOF branch in TreeS for current event
78 char tname[100]; sprintf(tname,"TreeS%d",ievent);
79
80 TTree *sdigitstree=(TTree*)file->Get(tname);
81
82 TClonesArray * fSDigits= new TClonesArray("AliTOFSDigit", 1000);
83 sdigitstree->GetBranch("TOF")->SetAddress(&fSDigits);
84
85 Int_t nEntries = sdigitstree->GetEntries();
86 //cout << nEntries << endl;
87
88 // Loop through all entries in the tree
89 Int_t nbytes;
90
91 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
92
93 // Import the tree
94 nbytes += sdigitstree->GetEvent(iEntry);
95
96 // Get the number of sdigits
97 Int_t ndig = fSDigits->GetEntriesFast();
98 cout << "------------------<AliTOFSDigits2Digits>------------------" << endl;
99 cout << "Found " << ndig << " TOF SDigits for event " << ievent << endl;
100 cout << "------------------------------------------------------------" << endl;
101
102 // start loop on sdigits
103 for (Int_t k = 0; k < ndig; k++) {
104
105 Int_t vol[5]; // location for a digit
106
107 // Get the information for this digit
108 AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigits->UncheckedAt(k);
109
110 Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
111 // for current sdigit
112
113 // TOF sdigit volumes (always the same for all slots)
114 Int_t sector = tofsdigit->GetSector(); // range [1-18]
115 Int_t plate = tofsdigit->GetPlate(); // range [1- 5]
116 Int_t strip = tofsdigit->GetStrip(); // range [1-20]
117 Int_t padz = tofsdigit->GetPadz(); // range [1- 2]
118 Int_t padx = tofsdigit->GetPadx(); // range [1-48]
119
120 vol[0] = sector;
121 vol[1] = plate;
122 vol[2] = strip;
123 vol[3] = padx;
124 vol[4] = padz;
125
126 //--------------------- QA section ----------------------
127 // in the while, I perform QA
128 Bool_t isSDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
129
130 if (isSDigitBad) {
131 cout << "<AliTOFSDigits2Digits> strange sdigit found" << endl;
132 rc = 3;
133 return rc;
134 }
135 //-------------------------------------------------------
136
137 //------------------- Dump section ----------------------
138 if(k<ndump){
139 cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
140 cout << k << "-th sdigit" << endl;
141 cout << "----------------------------------------------------"<< endl;
142 }
143 // ------------------------------------------------------
144
145 // start loop on number of slots for current sdigit
146 for (Int_t islot = 0; islot < nslot; islot++) {
147 Float_t digit[2]; // TOF digit variables
148 Int_t tracknum[kMAXDIGITS]; // contributing tracks for the current slot
149
150 Float_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
151 Float_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
152
153 tracknum[0]=tofsdigit->GetTrack(islot,0);
154 tracknum[1]=tofsdigit->GetTrack(islot,1);
155 tracknum[2]=tofsdigit->GetTrack(islot,2);
156
157 for (Int_t i = 0; i < kMAXDIGITS; i++) {
158 tottracks++;
159 // search for the first empty location
160 if(tracknum[i]==-1){
161 tottracks--;
162 break;
163 }
164 }
165
166 // adding a TOF digit for each slot
167 tof->AddDigit(tracknum, vol, digit);
168 totndig++;
169 }
170
171
172 } // end loop on sdigits
173
174 } // end loop on entries
175
176 // free used memory
177 fSDigits->Clear();
178 fSDigits=0;
179
180 gAlice->TreeD()->Reset();
181 gAlice->TreeD()->Fill();
182 //gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
183 gAlice->TreeD()->AutoSave();
184
185 } // end loop on events
186
187 cout << "----------------------------------------------------------" << endl;
188 cout << "<AliTOFSDigits2Digits> Summary" << endl;
189 cout << "contributing tracks to " << totndig << " digits: " << tottracks << endl;
190 cout << "----------------------------------------------------------" << endl;
191
192 return rc;
193
194
195}