]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFDigitizer.cxx
AliESDtrack flags update during the TOF matching procedure. Now, (1) when a TPC/TRD...
[u/mrichter/AliRoot.git] / TOF / AliTOFDigitizer.cxx
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 //                                                                         //
18 // This is a TTask that makes TOF-Digits out of TOF-SDigits.               //
19 // The simulation of the detector is performed at sdigits level:           //
20 // during digitization the unique task is the sum of all sdigits in the    //
21 // same pad.                                                               //
22 // Digits are written to TreeD in branch "TOF".                            //
23 //                                                                         //
24 // -- Author :  F. Pierella (Bologna University) pierella@bo.infn.it       //
25 //                                                                         //
26 //_________________________________________________________________________//
27
28 //#include "Riostream.h"
29
30 //#include "TFile.h"
31 #include "TH1F.h"
32 #include "TTree.h"
33 #include "TRandom.h"
34 #include "TObjArray.h"
35
36 #include "AliLoader.h"
37 #include "AliLog.h"
38 #include "AliRunDigitizer.h"
39 #include "AliRunLoader.h"
40 #include "AliRun.h"
41
42 #include "AliTOFcalib.h"
43 //#include "AliTOFChannelOnline.h"
44 #include "AliTOFChannelOffline.h"
45 #include "AliTOFDigitizer.h"
46 #include "AliTOFdigit.h"
47 #include "AliTOFHitMap.h"
48 #include "AliTOFGeometry.h"
49 #include "AliTOFSDigit.h"
50 #include "AliTOF.h"
51
52 extern TDirectory *gDirectory;
53 //extern TFile *gFile;
54 extern TRandom *gRandom;
55
56 extern AliRun *gAlice;
57
58
59 ClassImp(AliTOFDigitizer)
60
61 //___________________________________________
62   AliTOFDigitizer::AliTOFDigitizer()  :
63     AliDigitizer(),
64     fDigits(new TClonesArray("AliTOFdigit",4000)),
65     fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
66   fhitMap(0x0),
67   fCalib(new AliTOFcalib())
68 {
69   // Default ctor - don't use it
70   InitDecalibration();
71 }
72
73 //___________________________________________
74 AliTOFDigitizer::AliTOFDigitizer(AliRunDigitizer* manager): 
75   AliDigitizer(manager), 
76   fDigits(new TClonesArray("AliTOFdigit",4000)),
77   fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
78   fhitMap(0x0),
79   fCalib(new AliTOFcalib())
80 {
81   //ctor with RunDigitizer
82   InitDecalibration();
83 }
84
85 //------------------------------------------------------------------------
86 AliTOFDigitizer::AliTOFDigitizer(const AliTOFDigitizer &source):
87   AliDigitizer(source),
88   fDigits(source.fDigits),
89   fSDigitsArray(source.fSDigitsArray),
90   fhitMap(source.fhitMap),
91   fCalib(source.fCalib)
92 {
93   // copy constructor
94 }
95
96 //------------------------------------------------------------------------
97   AliTOFDigitizer& AliTOFDigitizer::operator=(const AliTOFDigitizer &source)
98 {
99   // ass. op.
100   
101   if (this == &source)
102     return *this;
103
104   AliDigitizer::operator=(source);
105   fDigits=source.fDigits;
106   fSDigitsArray=source.fSDigitsArray;
107   fhitMap=source.fhitMap;
108   fCalib=source.fCalib;
109   return *this;
110
111 }
112
113 //------------------------------------------------------------------------
114 AliTOFDigitizer::~AliTOFDigitizer()
115 {
116   // Destructor
117   delete fCalib;
118   if (fDigits){
119     fDigits->Delete();
120     delete fDigits;
121     fDigits=0x0;
122   }
123   if (fSDigitsArray){
124     fSDigitsArray->Delete();
125     delete fSDigitsArray;
126     fSDigitsArray=0x0;
127   }
128 }
129
130 //---------------------------------------------------------------------
131
132 void AliTOFDigitizer::Exec(Option_t* /*option*/)
133 {
134   //
135   // Perform digitization and merging.
136   // The algorithm is the following:
137   // - a hitmap is created to check if a pad is already activated;
138   // - an sdigits container is created to collect all sdigits from
139   //   different files;
140   // - sdigits are summed using the hitmap;
141   // - the sdigits container is used to create the array of AliTOFdigit.
142   //
143
144   AliDebug(1, "");
145
146
147   // get the ptr to TOF detector
148   AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
149
150   //Make branches
151   char branchname[20];
152   sprintf (branchname, "%s", tof->GetName ());
153  
154   AliRunLoader* outrl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
155   if (outrl == 0x0)
156    {
157      AliError("Can not find Run Loader in output folder.");
158      return;
159    }
160    
161   AliLoader* outgime = outrl->GetLoader("TOFLoader");
162   if (outgime == 0x0)
163    {
164      AliError("Can not get TOF Loader from Output Run Loader.");
165      return;
166    }
167   
168   TTree* treeD = outgime->TreeD();
169   if (treeD == 0x0)
170    {
171      outgime->MakeTree("D");
172      treeD = outgime->TreeD();
173    }
174   //Make branch for digits (to be created in Init())
175   tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
176
177   // container for all summed sdigits (to be created in Init())
178   //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
179   
180   // create hit map (to be created in Init())
181   fhitMap = new AliTOFHitMap(fSDigitsArray);
182   
183   // Loop over files to digitize
184
185   for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
186        inputFile++) {
187     ReadSDigit(inputFile);
188    }
189
190   // create digits
191   CreateDigits();
192
193   // free used memory for Hit Map in current event
194   delete fhitMap;
195   fSDigitsArray->Clear();
196
197   treeD->Fill();
198
199   AliDebug(2,"----------------------------------------");
200   AliInfo(Form("%d digits have been created", fDigits->GetEntriesFast()));
201   AliDebug(2,"----------------------------------------");
202
203   outgime->WriteDigits("OVERWRITE");
204   outgime->UnloadDigits();
205   fDigits->Clear();
206
207 }
208
209 //---------------------------------------------------------------------
210
211 void AliTOFDigitizer::CreateDigits()
212 {
213   // loop on sdigits container to fill the AliTOFdigit TClonesArray
214   // start digitizing all the collected sdigits 
215
216   Int_t ndump=0; // dump the first ndump created digits for each event
217
218   // get the total number of collected sdigits
219   Int_t ndig = fSDigitsArray->GetEntriesFast();
220
221   for (Int_t k = 0; k < ndig; k++) {
222     
223     Int_t  vol[5];  // location for a digit
224     for (Int_t i=0; i<5; i++) vol[i] = -1;
225     
226     // Get the information for this digit
227     AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
228     
229     Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
230     // for current sdigit
231     
232     // TOF sdigit volumes (always the same for all slots)
233     Int_t sector    = tofsdigit->GetSector(); // range [0-17]
234     Int_t plate     = tofsdigit->GetPlate();  // range [0- 4]
235     Int_t strip     = tofsdigit->GetStrip();  // range [0-14/18/19]
236     Int_t padz      = tofsdigit->GetPadz();   // range [0- 1]
237     Int_t padx      = tofsdigit->GetPadx();   // range [0-47]
238     
239     vol[0] = sector;
240     vol[1] = plate;
241     vol[2] = strip;
242     vol[3] = padx;
243     vol[4] = padz;
244     
245     //--------------------- QA section ----------------------
246     // in the while, I perform QA
247     Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
248     
249     if (isSDigitBad)
250       AliFatal(Form("strange sdigit found   %2d  %1d  %2d  %1d %2d", sector, plate, strip, padz, padx));
251     //-------------------------------------------------------
252     
253     //------------------- Dump section ----------------------
254     if (k<ndump) {
255       AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
256       AliInfo("----------------------------------------------------");
257     }
258     // ------------------------------------------------------
259     
260     // start loop on number of slots for current sdigit
261     for (Int_t islot = 0; islot < nslot; islot++) {
262       Int_t  digit[4] = {-1,-1,-1,-1};     // TOF digit variables
263       Int_t tracknum[AliTOFSDigit::kMAXDIGITS];     // contributing tracks for the current slot
264       
265       Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
266       Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
267
268       //if (tdc>=8192) continue;//AdC
269
270       tracknum[0]=tofsdigit->GetTrack(islot,0);
271       tracknum[1]=tofsdigit->GetTrack(islot,1);
272       tracknum[2]=tofsdigit->GetTrack(islot,2);
273       
274       // new with placement must be used
275       // adding a TOF digit for each slot
276       TClonesArray &aDigits = *fDigits;
277       Int_t last=fDigits->GetEntriesFast();
278       new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
279
280     }
281     
282   } // end loop on sdigits - end digitizing all collected sdigits
283
284   //Insert Decalibration 
285   AliDebug(2,"in digitizer, create digits");
286   DecalibrateTOFSignal();
287 }
288
289 //---------------------------------------------------------------------
290
291 void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
292 {
293   // Read sdigits for current event and inputFile; 
294   // store them into the sdigits container
295   // and update the hit map
296   // SDigits from different files are assumed to
297   // be created with the same simulation parameters.
298   
299   // creating the TClonesArray to store the digits
300   static TClonesArray sdigitsClonesArray("AliTOFSDigit",  1000); 
301   sdigitsClonesArray.Clear();
302
303   // get the treeS from manager
304   AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
305   if (rl == 0x0)
306    {
307      AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
308      return;
309    }
310
311   AliLoader* gime = rl->GetLoader("TOFLoader");
312   if (gime == 0x0)
313    {
314      AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
315      return;
316    }
317
318   TTree* currentTreeS=gime->TreeS();
319   if (currentTreeS == 0x0)
320    {
321      Int_t retval = gime->LoadSDigits();
322      if (retval) 
323       {
324          AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
325          return;
326       }
327      currentTreeS=gime->TreeS();
328      if (currentTreeS == 0x0)
329       {
330          AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
331          return;
332       }
333    } 
334   // get the branch TOF inside the treeS
335   TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
336   // check if the branch exist
337   TBranch* tofBranch=currentTreeS->GetBranch("TOF");
338
339   if(!tofBranch){
340     AliFatal(Form("TOF branch not found for input %d",inputFile));
341   }
342   
343   tofBranch->SetAddress(&sdigitsDummyContainer);           
344   
345   Int_t nEntries = (Int_t)tofBranch->GetEntries();                                
346
347   // Loop through all entries in the tree
348   Int_t nbytes = 0;
349   
350   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
351     
352     // Import the tree
353     nbytes += tofBranch->GetEvent(iEntry);
354     
355     // Get the number of sdigits
356     Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
357     
358     for (Int_t k=0; k<ndig; k++) {
359       AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
360       
361       Int_t  vol[5]; // location for a sdigit
362       for (Int_t i=0; i<5; i++) vol[i] = -1;
363
364       // check the sdigit volume
365       vol[0] = tofSdigit->GetSector();
366       vol[1] = tofSdigit->GetPlate();
367       vol[2] = tofSdigit->GetStrip();
368       vol[3] = tofSdigit->GetPadx();
369       vol[4] = tofSdigit->GetPadz();
370       
371       if (fhitMap->TestHit(vol) != kEmpty) {
372         AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
373         sdig->Update(tofSdigit);
374
375       } else {
376
377         CollectSDigit(tofSdigit); // collect the current sdigit
378         fhitMap->SetHit(vol);     // update the hitmap for location vol
379
380       } // if (hitMap->TestHit(vol) != kEmpty)
381       
382     } // for (Int_t k=0; k<ndig; k++)
383
384   } // end loop on entries
385
386 }
387
388
389 //_____________________________________________________________________________
390 void AliTOFDigitizer::CollectSDigit(AliTOFSDigit * sdigit)
391 {
392   //
393   // Add a TOF sdigit in container
394   // new with placement must be used
395   TClonesArray &aSDigitsArray = *fSDigitsArray;
396   Int_t last=fSDigitsArray->GetEntriesFast();
397   // make a copy of the current sdigit and
398   // put it into tmp array
399   new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
400 }
401
402 //_____________________________________________________________________________
403 void AliTOFDigitizer::InitDecalibration() const {
404   //
405   //
406   //
407
408   fCalib->CreateCalArrays();
409   fCalib->ReadSimHistoFromCDB("TOF/Calib", -1); // use AliCDBManager's number
410   fCalib->ReadParOfflineFromCDB("TOF/Calib", -1); // use AliCDBManager's number
411 }
412 //---------------------------------------------------------------------
413 void AliTOFDigitizer::DecalibrateTOFSignal(){
414
415   // Read Calibration parameters from the CDB
416
417   TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
418
419   AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
420
421   // Initialize Quantities to Simulate ToT Spectra
422
423   TH1F * hToT= fCalib->GetTOFSimToT();
424   Int_t nbins = hToT->GetNbinsX();
425   Float_t delta = hToT->GetBinWidth(1);
426   Float_t maxch = hToT->GetBinLowEdge(nbins)+delta;
427   Float_t minch = hToT->GetBinLowEdge(1);
428   Float_t max=0,min=0; //maximum and minimum value of the distribution
429   Int_t maxbin=0,minbin=0; //maximum and minimum bin of the distribution
430
431   for (Int_t ii=nbins; ii>0; ii--){
432     if (hToT->GetBinContent(ii)!= 0) {
433       max = maxch - (nbins-ii-1)*delta;
434       maxbin = ii; 
435       break;}
436   }
437   for (Int_t j=1; j<nbins; j++){
438     if (hToT->GetBinContent(j)!= 0) {
439       min = minch + (j-1)*delta;
440       minbin = j; 
441       break;}
442   }
443
444   Float_t maxToT=max;
445   Float_t minToT=min;
446  
447   Float_t maxToTDistr=hToT->GetMaximum();
448
449   AliDebug (1, Form(" The minimum ToT = %f", minToT)); 
450   AliDebug (1, Form(" The maximum ToT = %f", maxToT)); 
451   AliDebug (1, Form(" The maximum peak in ToT = %f", maxToTDistr)); 
452   
453   // Loop on TOF Digits
454
455   Bool_t isToTSimulated=kFALSE;
456   Bool_t misCalibPars=kFALSE;
457   if(hToT->GetEntries()>0)isToTSimulated=kTRUE;  
458   Int_t ndigits = fDigits->GetEntriesFast();
459   for (Int_t i=0;i<ndigits;i++){
460     AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
461     Int_t detId[5];
462     detId[0] = dig->GetSector();
463     detId[1] = dig->GetPlate();
464     detId[2] = dig->GetStrip();
465     detId[3] = dig->GetPadz();
466     detId[4] = dig->GetPadx();
467     dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
468     if(isToTSimulated){  
469
470       //A realistic ToT Spectrum was found in input, 
471       //decalibrated TOF Digits likely to be simulated....
472  
473       Int_t index = AliTOFGeometry::GetIndex(detId); // The channel index    
474       AliTOFChannelOffline *calChannelOffline = (AliTOFChannelOffline *)calOffline->At(index); //retrieve the info for time slewing 
475       Double_t par[6];  // time slewing parameters
476   
477       //check whether we actually ask for miscalibration
478
479       for (Int_t j = 0; j<6; j++){
480         par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
481         if(par[j]!=0)misCalibPars=kTRUE;
482       }
483       AliDebug(2,Form(" Calib Pars = %f (0-th parameter for time slewing + time delay), %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
484
485       // Now generate Realistic ToT distribution from TestBeam Data. 
486       // Tot is in ns, assuming a Matching Window of 10 ns.
487
488       Float_t simToT = 0;
489       Float_t trix = 0;
490       Float_t triy = 0;
491       Double_t timeCorr;
492       Double_t tToT;
493       while (simToT <= triy){
494         trix = gRandom->Rndm(i);
495         triy = gRandom->Rndm(i);
496         trix = (maxToT-minToT)*trix + minToT; 
497         triy = maxToTDistr*triy;
498         Int_t binx=hToT->FindBin(trix);
499         simToT=hToT->GetBinContent(binx);
500       }
501       // the generated ToT (ns)
502       tToT= (Double_t) trix; // to apply slewing we start from ns..
503       // transform TOF signal in ns
504       AliDebug(2,Form(" The Initial Time (counts): %i: ",dig->GetTdc()));
505       AliDebug(2,Form(" Time before miscalibration (ps) %e: ",dig->GetTdc()*(Double_t)AliTOFGeometry::TdcBinWidth()));
506       // add slewing effect
507       timeCorr=par[0] + tToT*(par[1] +tToT*(par[2] +tToT*(par[3] +tToT*(par[4] +tToT*par[5])))); 
508       AliDebug(2,Form(" The Time slewing + delay (ns): %f: ",timeCorr));
509       // add global time shift
510       //convert to ps
511       timeCorr*=1E3;
512       Double_t timeMis = (Double_t)(dig->GetTdc())*(Double_t)AliTOFGeometry::TdcBinWidth();
513       timeMis = timeMis+timeCorr;
514       AliDebug(2,Form(" The Miscalibrated time (ps): %e: ",timeMis));
515
516       // now update the digit info
517  
518       Int_t tdcCorr= (Int_t)(timeMis/AliTOFGeometry::TdcBinWidth());
519       AliDebug(2,Form(" Final Time (counts): %i: ",tdcCorr));
520       // Setting Decalibrated Time signal (TDC counts)    
521       dig->SetTdc(tdcCorr);   
522       // Setting realistic ToT signal (TDC counts) 
523       tToT*=1E3; //back to ps  
524       Int_t tot=(Int_t)(tToT/AliTOFGeometry::ToTBinWidth());//(factor 1E3 as input ToT is in ns)
525       dig->SetToT(tot); 
526       AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
527       if(tdcCorr<0){
528         AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr)); 
529         AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot)); 
530       }
531     }
532     else{
533     // For Data with no Miscalibration, set ToT signal == Adc
534       dig->SetToT((Int_t)(dig->GetAdc()/AliTOFGeometry::ToTBinWidth())); //remove the factor 10^3 just to have a reasonable ToT range for raw data simulation even in the case of non-realistic ToT distribution (n.b. fAdc is practically an arbitrary quantity, and ToT has no impact on the TOF reco for non-miscalibrated digits)
535     }
536   }
537
538   if(!isToTSimulated)
539     AliDebug(1,"Standard Production, no miscalibrated digits");
540   else
541     if(!misCalibPars)
542       AliDebug(1,"Standard Production, no miscalibrated digits");
543     else
544       AliDebug(1,"Simulating miscalibrated digits");
545
546   return;
547 }
548