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