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