1 /**************************************************************************
2 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //_________________________________________________________________________//
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 //
22 // Digits are written to TreeD in branch "TOF". //
24 // -- Author : F. Pierella (Bologna University) pierella@bo.infn.it //
26 //_________________________________________________________________________//
28 //#include "Riostream.h"
35 #include "TObjArray.h"
37 #include "AliLoader.h"
39 #include "AliRunDigitizer.h"
40 #include "AliRunLoader.h"
43 #include "AliTOFcalib.h"
44 //#include "AliTOFChannelOnlineArray.h"
45 //#include "AliTOFChannelOnlineStatusArray.h"
46 //#include "AliTOFChannelOffline.h"
47 #include "AliTOFDigitizer.h"
48 #include "AliTOFdigit.h"
49 #include "AliTOFHitMap.h"
50 #include "AliTOFGeometry.h"
51 #include "AliTOFSDigit.h"
54 extern TDirectory *gDirectory;
55 //extern TFile *gFile;
56 extern TRandom *gRandom;
58 extern AliRun *gAlice;
61 ClassImp(AliTOFDigitizer)
63 //___________________________________________
64 AliTOFDigitizer::AliTOFDigitizer() :
66 fDigits(new TClonesArray("AliTOFdigit",4000)),
67 fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
69 fCalib(new AliTOFcalib())
71 // Default ctor - don't use it
75 //___________________________________________
76 AliTOFDigitizer::AliTOFDigitizer(AliRunDigitizer* manager):
77 AliDigitizer(manager),
78 fDigits(new TClonesArray("AliTOFdigit",4000)),
79 fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
81 fCalib(new AliTOFcalib())
83 //ctor with RunDigitizer
87 //------------------------------------------------------------------------
88 AliTOFDigitizer::AliTOFDigitizer(const AliTOFDigitizer &source):
90 fDigits(source.fDigits),
91 fSDigitsArray(source.fSDigitsArray),
92 fhitMap(source.fhitMap),
98 //------------------------------------------------------------------------
99 AliTOFDigitizer& AliTOFDigitizer::operator=(const AliTOFDigitizer &source)
106 AliDigitizer::operator=(source);
107 fDigits=source.fDigits;
108 fSDigitsArray=source.fSDigitsArray;
109 fhitMap=source.fhitMap;
110 fCalib=source.fCalib;
115 //------------------------------------------------------------------------
116 AliTOFDigitizer::~AliTOFDigitizer()
126 fSDigitsArray->Delete();
127 delete fSDigitsArray;
132 //---------------------------------------------------------------------
134 void AliTOFDigitizer::Exec(Option_t* /*option*/)
137 // Perform digitization and merging.
138 // The algorithm is the following:
139 // - a hitmap is created to check if a pad is already activated;
140 // - an sdigits container is created to collect all sdigits from
142 // - sdigits are summed using the hitmap;
143 // - the sdigits container is used to create the array of AliTOFdigit.
149 // get the ptr to TOF detector
150 AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
154 const Int_t kSize = 20;
155 char branchname[kSize];
156 snprintf(branchname,kSize,"%s", tof->GetName ());
158 AliRunLoader* outrl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
161 AliError("Can not find Run Loader in output folder.");
165 AliLoader* outgime = outrl->GetLoader("TOFLoader");
168 AliError("Can not get TOF Loader from Output Run Loader.");
172 TTree* treeD = outgime->TreeD();
175 outgime->MakeTree("D");
176 treeD = outgime->TreeD();
178 //Make branch for digits (to be created in Init())
179 tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
181 // container for all summed sdigits (to be created in Init())
182 //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
184 // create hit map (to be created in Init())
185 fhitMap = new AliTOFHitMap(fSDigitsArray);
187 // Loop over files to digitize
189 for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
191 ReadSDigit(inputFile);
197 // free used memory for Hit Map in current event
199 fSDigitsArray->Clear();
203 AliDebug(2,"----------------------------------------");
204 AliDebug(1,Form("%d digits have been created", fDigits->GetEntriesFast()));
205 AliDebug(2,"----------------------------------------");
207 outgime->WriteDigits("OVERWRITE");
208 outgime->UnloadDigits();
213 //---------------------------------------------------------------------
215 void AliTOFDigitizer::CreateDigits()
217 // loop on sdigits container to fill the AliTOFdigit TClonesArray
218 // start digitizing all the collected sdigits
220 Int_t ndump=0; // dump the first ndump created digits for each event
222 // get the total number of collected sdigits
223 Int_t ndig = fSDigitsArray->GetEntriesFast();
225 Int_t vol[5]={-1,-1,-1,-1,-1}; // location for a digit
226 Int_t digit[4] = {0,0,0,0}; // TOF digit variables
227 Int_t tracknum[AliTOFSDigit::kMAXDIGITS]; // contributing tracks for the current slot
228 for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
230 for (Int_t k = 0; k < ndig; k++) {
232 for (Int_t i=0; i<5; i++) vol[i] = -1;
234 // Get the information for this digit
235 AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
237 Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
238 // for current sdigit
240 // TOF sdigit volumes (always the same for all slots)
241 Int_t sector = tofsdigit->GetSector(); // range [0-17]
242 Int_t plate = tofsdigit->GetPlate(); // range [0- 4]
243 Int_t strip = tofsdigit->GetStrip(); // range [0-14/18/19]
244 Int_t padz = tofsdigit->GetPadz(); // range [0- 1]
245 Int_t padx = tofsdigit->GetPadx(); // range [0-47]
253 //--------------------- QA section ----------------------
254 // in the while, I perform QA
255 Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
258 AliFatal(Form("strange sdigit found %2d %1d %2d %1d %2d", sector, plate, strip, padz, padx));
259 //-------------------------------------------------------
261 //------------------- Dump section ----------------------
263 AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
264 AliInfo("----------------------------------------------------");
266 // ------------------------------------------------------
268 // start loop on number of slots for current sdigit
269 for (Int_t islot = 0; islot < nslot; islot++) {
270 for (Int_t aa=0; aa<4; aa++) digit[aa] = 0; // TOF digit variables
271 for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
273 Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
274 Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
276 //if (tdc>=8192) continue;//AdC
278 tracknum[0]=tofsdigit->GetTrack(islot,0);
279 tracknum[1]=tofsdigit->GetTrack(islot,1);
280 tracknum[2]=tofsdigit->GetTrack(islot,2);
282 // new with placement must be used
283 // adding a TOF digit for each slot
284 TClonesArray &aDigits = *fDigits;
285 Int_t last=fDigits->GetEntriesFast();
286 new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
290 } // end loop on sdigits - end digitizing all collected sdigits
292 //Insert Decalibration
293 AliDebug(2,"in digitizer, create digits");
294 DecalibrateTOFSignal();
297 //---------------------------------------------------------------------
299 void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
301 // Read sdigits for current event and inputFile;
302 // store them into the sdigits container
303 // and update the hit map
304 // SDigits from different files are assumed to
305 // be created with the same simulation parameters.
307 // creating the TClonesArray to store the digits
308 static TClonesArray sdigitsClonesArray("AliTOFSDigit", 1000);
309 sdigitsClonesArray.Clear();
311 // get the treeS from manager
312 AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
315 AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
319 AliLoader* gime = rl->GetLoader("TOFLoader");
322 AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
326 TTree* currentTreeS=gime->TreeS();
327 if (currentTreeS == 0x0)
329 Int_t retval = gime->LoadSDigits();
332 AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
335 currentTreeS=gime->TreeS();
336 if (currentTreeS == 0x0)
338 AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
342 // get the branch TOF inside the treeS
343 TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
344 // check if the branch exist
345 TBranch* tofBranch=currentTreeS->GetBranch("TOF");
348 AliFatal(Form("TOF branch not found for input %d",inputFile));
352 tofBranch->SetAddress(&sdigitsDummyContainer);
354 Int_t nEntries = (Int_t)tofBranch->GetEntries();
356 // Loop through all entries in the tree
359 Int_t vol[5]; // location for a sdigit
360 for (Int_t i=0; i<5; i++) vol[i] = -1;
362 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
365 nbytes += tofBranch->GetEvent(iEntry);
367 // Get the number of sdigits
368 Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
370 for (Int_t k=0; k<ndig; k++) {
371 AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
373 for (Int_t i=0; i<5; i++) vol[i] = -1;
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();
382 if (fhitMap->TestHit(vol) != kEmpty) {
383 AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
384 sdig->Update(tofSdigit);
388 CollectSDigit(tofSdigit); // collect the current sdigit
389 fhitMap->SetHit(vol); // update the hitmap for location vol
391 } // if (hitMap->TestHit(vol) != kEmpty)
393 } // for (Int_t k=0; k<ndig; k++)
395 } // end loop on entries
400 //_____________________________________________________________________________
401 void AliTOFDigitizer::CollectSDigit(const AliTOFSDigit * const sdigit)
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);
413 //_____________________________________________________________________________
414 void AliTOFDigitizer::InitDecalibration() const {
416 // Initialize TOF digits decalibration
421 fCalib->CreateCalArrays();
422 fCalib->ReadSimHistoFromCDB("TOF/Calib", -1); // use AliCDBManager's number
423 fCalib->ReadParOfflineFromCDB("TOF/Calib", -1); // use AliCDBManager's number
426 //---------------------------------------------------------------------
427 void AliTOFDigitizer::DecalibrateTOFSignal() {
429 // Decalibrate TOF signals according to OCDB parameters
432 Double_t time=0., tot=0., corr=0.;
433 Int_t deltaBC=0, l0l1=0, tdcBin=0;
435 Int_t detId[5] ={-1,-1,-1,-1,-1};
438 Int_t ndigits = fDigits->GetEntriesFast();
439 // Loop on TOF Digits
440 for (Int_t i=0;i<ndigits;i++){
441 AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
442 detId[0] = dig->GetSector();
443 detId[1] = dig->GetPlate();
444 detId[2] = dig->GetStrip();
445 detId[3] = dig->GetPadz();
446 detId[4] = dig->GetPadx();
447 dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
449 index = AliTOFGeometry::GetIndex(detId); // The channel index
451 // Read Calibration parameters from the CDB
453 time = dig->GetTdc() * AliTOFGeometry::TdcBinWidth(); /* ps */
454 tot = dig->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
455 deltaBC = 0;//dig->GetDeltaBC();
456 l0l1 = 0;//dig->GetL0L1Latency();
459 corr = fCalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
460 AliDebug(2, Form("calibrate index %d: time=%f (ps) tot=%f (ns) deltaBC=%d l0l1=%d timestamp=%d corr=%f (ps)",
461 index, time, tot, deltaBC, l0l1, timestamp, corr));
463 // apply time correction
466 // convert in TDC bins and set digit
467 //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
468 tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
473 AliDebug(1,"Simulating miscalibrated digits");
478 //---------------------------------------------------------------------
480 void AliTOFDigitizer::DecalibrateTOFSignal(){ // Old implementation
482 // Read Calibration parameters from the CDB
484 TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
486 AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
488 // Initialize Quantities to Simulate ToT Spectra
490 TH1F * hToT= fCalib->GetTOFSimToT();
491 Int_t nbins = hToT->GetNbinsX();
492 Float_t delta = hToT->GetBinWidth(1);
493 Float_t maxch = hToT->GetBinLowEdge(nbins)+delta;
494 Float_t minch = hToT->GetBinLowEdge(1);
495 Float_t max=0,min=0; //maximum and minimum value of the distribution
496 Int_t maxbin=0,minbin=0; //maximum and minimum bin of the distribution
498 for (Int_t ii=nbins; ii>0; ii--){
499 if (hToT->GetBinContent(ii)!= 0) {
500 max = maxch - (nbins-ii-1)*delta;
504 for (Int_t j=1; j<nbins; j++){
505 if (hToT->GetBinContent(j)!= 0) {
506 min = minch + (j-1)*delta;
514 Float_t maxToTDistr=hToT->GetMaximum();
516 AliDebug (1, Form(" The minimum ToT = %f", minToT));
517 AliDebug (1, Form(" The maximum ToT = %f", maxToT));
518 AliDebug (1, Form(" The maximum peak in ToT = %f", maxToTDistr));
520 // Loop on TOF Digits
522 Bool_t isToTSimulated=kFALSE;
523 Bool_t misCalibPars=kFALSE;
524 if(hToT->GetEntries()>0)isToTSimulated=kTRUE;
525 Int_t ndigits = fDigits->GetEntriesFast();
526 for (Int_t i=0;i<ndigits;i++){
527 AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
529 detId[0] = dig->GetSector();
530 detId[1] = dig->GetPlate();
531 detId[2] = dig->GetStrip();
532 detId[3] = dig->GetPadz();
533 detId[4] = dig->GetPadx();
534 dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
537 //A realistic ToT Spectrum was found in input,
538 //decalibrated TOF Digits likely to be simulated....
540 Int_t index = AliTOFGeometry::GetIndex(detId); // The channel index
541 AliTOFChannelOffline *calChannelOffline = (AliTOFChannelOffline *)calOffline->At(index); //retrieve the info for time slewing
542 Double_t par[6]; // time slewing parameters
544 //check whether we actually ask for miscalibration
546 for (Int_t j = 0; j<6; j++){
547 par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
548 if(par[j]!=0)misCalibPars=kTRUE;
550 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]));
552 // Now generate Realistic ToT distribution from TestBeam Data.
553 // Tot is in ns, assuming a Matching Window of 10 ns.
560 while (simToT <= triy){
561 trix = gRandom->Rndm(i);
562 triy = gRandom->Rndm(i);
563 trix = (maxToT-minToT)*trix + minToT;
564 triy = maxToTDistr*triy;
565 Int_t binx=hToT->FindBin(trix);
566 simToT=hToT->GetBinContent(binx);
568 // the generated ToT (ns)
569 tToT= (Double_t) trix; // to apply slewing we start from ns..
570 // transform TOF signal in ns
571 AliDebug(2,Form(" The Initial Time (counts): %i: ",dig->GetTdc()));
572 AliDebug(2,Form(" Time before miscalibration (ps) %e: ",dig->GetTdc()*(Double_t)AliTOFGeometry::TdcBinWidth()));
573 // add slewing effect
574 timeCorr=par[0] + tToT*(par[1] +tToT*(par[2] +tToT*(par[3] +tToT*(par[4] +tToT*par[5]))));
575 AliDebug(2,Form(" The Time slewing + delay (ns): %f: ",timeCorr));
576 // add global time shift
579 Double_t timeMis = (Double_t)(dig->GetTdc())*(Double_t)AliTOFGeometry::TdcBinWidth();
580 timeMis = timeMis+timeCorr;
581 AliDebug(2,Form(" The Miscalibrated time (ps): %e: ",timeMis));
583 // now update the digit info
585 Int_t tdcCorr= (Int_t)(timeMis/AliTOFGeometry::TdcBinWidth());
586 AliDebug(2,Form(" Final Time (counts): %i: ",tdcCorr));
587 // Setting Decalibrated Time signal (TDC counts)
588 dig->SetTdc(tdcCorr);
589 // Setting realistic ToT signal (TDC counts)
590 tToT*=1E3; //back to ps
591 Int_t tot=(Int_t)(tToT/AliTOFGeometry::ToTBinWidth());//(factor 1E3 as input ToT is in ns)
593 AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
595 AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr));
596 AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot));
600 // For Data with no Miscalibration, set ToT signal == Adc
601 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)
606 AliDebug(1,"Standard Production, no miscalibrated digits");
609 AliDebug(1,"Standard Production, no miscalibrated digits");
611 AliDebug(1,"Simulating miscalibrated digits");