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