]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - TOF/AliTOFDigitizer.cxx
- Correct initialisation of sRandom.
[u/mrichter/AliRoot.git] / TOF / AliTOFDigitizer.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16//_________________________________________________________________________
17// This is a TTask that makes TOF-Digits out of TOF-SDigits.
18// The simulation of the detector is performed at sdigits level:
19// during digitization the unique task is the sum of all sdigits in the
20// same pad.
21// Digits are written to TreeD in branch "TOF".
22//
23// -- Author : F. Pierella (Bologna University) pierella@bo.infn.it
24//////////////////////////////////////////////////////////////////////////////
25
26#include <TTree.h>
27#include <TVector.h>
28#include <TObjArray.h>
29#include <TFile.h>
30#include <TDirectory.h>
31#include <TRandom.h>
32
33
34#include "AliTOFDigitizer.h"
35#include "AliTOF.h"
36#include "AliTOFSDigitizer.h"
37#include "AliTOFhit.h"
38#include "AliTOFdigit.h"
39#include "AliTOFSDigit.h"
40#include "AliTOFHitMap.h"
41#include "AliDigitizer.h"
42#include "AliRunDigitizer.h"
43
44#include "AliRun.h"
45#include "AliPDG.h"
46
47#include <stdlib.h>
48#include <Riostream.h>
49#include <Riostream.h>
50
51ClassImp(AliTOFDigitizer)
52
53//___________________________________________
54 AliTOFDigitizer::AliTOFDigitizer() :AliDigitizer()
55{
56 // Default ctor - don't use it
57 fDigits=0;
58 fSDigitsArray=0;
59 fhitMap=0;
60}
61
62//___________________________________________
63AliTOFDigitizer::AliTOFDigitizer(AliRunDigitizer* manager)
64 :AliDigitizer(manager)
65{
66 fDigits=0;
67 fSDigitsArray=0;
68 fhitMap=0;
69}
70
71//------------------------------------------------------------------------
72AliTOFDigitizer::~AliTOFDigitizer()
73{
74 // Destructor
75}
76
77//---------------------------------------------------------------------
78
79void AliTOFDigitizer::Exec(Option_t* option)
80{
81 //
82 // Perform digitization and merging.
83 // The algorithm is the following:
84 // - a hitmap is created to check if a pad is already activated;
85 // - an sdigits container is created to collect all sdigits from
86 // different files;
87 // - sdigits are summed using the hitmap;
88 // - the sdigits container is used to create the array of AliTOFdigit.
89 //
90
91 if(strstr(option,"deb")) cout<<"AliTOFDigitizer::Exec\n";
92
93
94 // get the ptr to TOF detector
95 AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
96
97 //Make branches
98 char branchname[20];
99 sprintf (branchname, "%s", tof->GetName ());
100
101 fDigits=new TClonesArray("AliTOFdigit",4000);
102
103 TTree* treeD = fManager->GetTreeD();
104 //Make branch for digits (to be created in Init())
105 tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
106
107 // container for all summed sdigits (to be created in Init())
108 fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
109
110 // create hit map (to be created in Init())
111 fhitMap = new AliTOFHitMap(fSDigitsArray);
112
113 // Loop over files to digitize
114
115 for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
116 inputFile++) {
117 ReadSDigit(inputFile);
118 }
119
120 // create digits
121 CreateDigits();
122
123 // free used memory for Hit Map in current event
124 delete fhitMap;
125 fSDigitsArray->Delete();
126 treeD->Fill();
127
128 fManager->GetTreeD()->AutoSave(); // to fit with the framework
129 fDigits->Delete();
130 delete fDigits;
131
132}
133
134//---------------------------------------------------------------------
135
136void AliTOFDigitizer::CreateDigits()
137{
138 // loop on sdigits container to fill the AliTOFdigit TClonesArray
139 // start digitizing all the collected sdigits
140
141 Int_t ndump=15; // dump the first ndump created digits for each event
142
143 // get the total number of collected sdigits
144 Int_t ndig = fSDigitsArray->GetEntriesFast();
145
146 for (Int_t k = 0; k < ndig; k++) {
147
148 Int_t vol[5]; // location for a digit
149
150 // Get the information for this digit
151 AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
152
153 Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
154 // for current sdigit
155
156 // TOF sdigit volumes (always the same for all slots)
157 Int_t sector = tofsdigit->GetSector(); // range [1-18]
158 Int_t plate = tofsdigit->GetPlate(); // range [1- 5]
159 Int_t strip = tofsdigit->GetStrip(); // range [1-20]
160 Int_t padz = tofsdigit->GetPadz(); // range [1- 2]
161 Int_t padx = tofsdigit->GetPadx(); // range [1-48]
162
163 vol[0] = sector;
164 vol[1] = plate;
165 vol[2] = strip;
166 vol[3] = padx;
167 vol[4] = padz;
168
169 //--------------------- QA section ----------------------
170 // in the while, I perform QA
171 Bool_t isSDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
172
173 if (isSDigitBad) {
174 cout << "<AliTOFSDigits2Digits> strange sdigit found" << endl;
175 abort();
176 }
177 //-------------------------------------------------------
178
179 //------------------- Dump section ----------------------
180 if(k<ndump){
181 cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
182 cout << k << "-th sdigit" << endl;
183 cout << "----------------------------------------------------"<< endl;
184 }
185 // ------------------------------------------------------
186
187 // start loop on number of slots for current sdigit
188 for (Int_t islot = 0; islot < nslot; islot++) {
189 Float_t digit[2]; // TOF digit variables
190 Int_t tracknum[kMAXDIGITS]; // contributing tracks for the current slot
191
192 Float_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
193 Float_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
194
195 tracknum[0]=tofsdigit->GetTrack(islot,0);
196 tracknum[1]=tofsdigit->GetTrack(islot,1);
197 tracknum[2]=tofsdigit->GetTrack(islot,2);
198
199 // new with placement must be used
200 // adding a TOF digit for each slot
201 TClonesArray &aDigits = *fDigits;
202 Int_t last=fDigits->GetEntriesFast();
203 new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
204
205 }
206
207 } // end loop on sdigits - end digitizing all collected sdigits
208
209}
210
211//---------------------------------------------------------------------
212
213void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
214{
215 // Read sdigits for current event and inputFile;
216 // store them into the sdigits container
217 // and update the hit map
218 // SDigits from different files are assumed to
219 // be created with the same simulation parameters.
220
221 // get the treeS from manager
222 TTree* currentTreeS=fManager->GetInputTreeS(inputFile);
223
224 // get the branch TOF inside the treeS
225 TClonesArray * sdigitsDummyContainer= new TClonesArray("AliTOFSDigit", 1000);
226
227 // check if the branch exist
228 TBranch* tofBranch=currentTreeS->GetBranch("TOF");
229
230 if(!tofBranch){
231 Fatal("ReadSDigit","TOF branch not found for input %d",inputFile);
232 }
233
234 tofBranch->SetAddress(&sdigitsDummyContainer);
235
236 Int_t nEntries = (Int_t)tofBranch->GetEntries();
237
238 // Loop through all entries in the tree
239 Int_t nbytes;
240
241 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
242
243 // Import the tree
244 nbytes += tofBranch->GetEvent(iEntry);
245
246 // Get the number of sdigits
247 Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
248
249 for (Int_t k=0; k<ndig; k++) {
250 AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
251
252 Int_t vol[5]; // location for a sdigit
253
254 // check the sdigit volume
255 vol[0] = tofSdigit->GetSector();
256 vol[1] = tofSdigit->GetPlate();
257 vol[2] = tofSdigit->GetStrip();
258 vol[3] = tofSdigit->GetPadx();
259 vol[4] = tofSdigit->GetPadz();
260
261 if (fhitMap->TestHit(vol) != kEmpty) {
262 AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
263 sdig->Update(tofSdigit);
264
265 } else {
266
267 CollectSDigit(tofSdigit); // collect the current sdigit
268 fhitMap->SetHit(vol); // update the hitmap for location vol
269
270 } // if (hitMap->TestHit(vol) != kEmpty)
271
272 } // for (Int_t k=0; k<ndig; k++)
273
274 } // end loop on entries
275
276 sdigitsDummyContainer->Delete();
277 sdigitsDummyContainer=0;
278
279}
280
281
282//_____________________________________________________________________________
283void AliTOFDigitizer::CollectSDigit(AliTOFSDigit * sdigit)
284{
285 //
286 // Add a TOF sdigit in container
287 // new with placement must be used
288 TClonesArray &aSDigitsArray = *fSDigitsArray;
289 Int_t last=fSDigitsArray->GetEntriesFast();
290 // make a copy of the current sdigit and
291 // put it into tmp array
292 new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
293}