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 "AliDigitizationInput.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 TRandom *gRandom;
56 extern AliRun *gAlice;
59 ClassImp(AliTOFDigitizer)
61 //___________________________________________
62 AliTOFDigitizer::AliTOFDigitizer() :
64 fDigits(new TClonesArray("AliTOFdigit",4000)),
65 fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
67 fCalib(new AliTOFcalib())
69 // Default ctor - don't use it
73 //___________________________________________
74 AliTOFDigitizer::AliTOFDigitizer(AliDigitizationInput* digInput):
75 AliDigitizer(digInput),
76 fDigits(new TClonesArray("AliTOFdigit",4000)),
77 fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
79 fCalib(new AliTOFcalib())
81 //ctor with RunDigitizer
85 //------------------------------------------------------------------------
86 AliTOFDigitizer::AliTOFDigitizer(const AliTOFDigitizer &source):
88 fDigits(source.fDigits),
89 fSDigitsArray(source.fSDigitsArray),
90 fhitMap(source.fhitMap),
96 //------------------------------------------------------------------------
97 AliTOFDigitizer& AliTOFDigitizer::operator=(const AliTOFDigitizer &source)
104 AliDigitizer::operator=(source);
105 fDigits=source.fDigits;
106 fSDigitsArray=source.fSDigitsArray;
107 fhitMap=source.fhitMap;
108 fCalib=source.fCalib;
113 //------------------------------------------------------------------------
114 AliTOFDigitizer::~AliTOFDigitizer()
124 fSDigitsArray->Delete();
125 delete fSDigitsArray;
130 //---------------------------------------------------------------------
132 void AliTOFDigitizer::Digitize(Option_t* /*option*/)
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
140 // - sdigits are summed using the hitmap;
141 // - the sdigits container is used to create the array of AliTOFdigit.
147 // get the ptr to TOF detector
148 AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
152 const Int_t kSize = 20;
153 char branchname[kSize];
154 snprintf(branchname,kSize,"%s", tof->GetName ());
156 AliRunLoader* outrl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
159 AliError("Can not find Run Loader in output folder.");
163 AliLoader* outgime = outrl->GetLoader("TOFLoader");
166 AliError("Can not get TOF Loader from Output Run Loader.");
170 TTree* treeD = outgime->TreeD();
173 outgime->MakeTree("D");
174 treeD = outgime->TreeD();
176 //Make branch for digits (to be created in Init())
177 tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
179 // container for all summed sdigits (to be created in Init())
180 //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
182 // create hit map (to be created in Init())
183 fhitMap = new AliTOFHitMap(fSDigitsArray);
185 // Loop over files to digitize
187 for (Int_t inputFile=0; inputFile<fDigInput->GetNinputs();
189 ReadSDigit(inputFile);
195 // free used memory for Hit Map in current event
197 fSDigitsArray->Clear();
201 AliDebug(2,"----------------------------------------");
202 AliDebug(1,Form("%d digits have been created", fDigits->GetEntriesFast()));
203 AliDebug(2,"----------------------------------------");
205 outgime->WriteDigits("OVERWRITE");
206 outgime->UnloadDigits();
211 //---------------------------------------------------------------------
213 void AliTOFDigitizer::CreateDigits()
215 // loop on sdigits container to fill the AliTOFdigit TClonesArray
216 // start digitizing all the collected sdigits
218 Int_t ndump=0; // dump the first ndump created digits for each event
220 // get the total number of collected sdigits
221 Int_t ndig = fSDigitsArray->GetEntriesFast();
223 Int_t vol[5]={-1,-1,-1,-1,-1}; // location for a digit
224 Int_t digit[4] = {0,0,0,0}; // TOF digit variables
225 Int_t tracknum[AliTOFSDigit::kMAXDIGITS]; // contributing tracks for the current slot
226 for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
228 for (Int_t k = 0; k < ndig; k++) {
230 for (Int_t i=0; i<5; i++) vol[i] = -1;
232 // Get the information for this digit
233 AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
235 Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
236 // for current sdigit
238 // TOF sdigit volumes (always the same for all slots)
239 Int_t sector = tofsdigit->GetSector(); // range [0-17]
240 Int_t plate = tofsdigit->GetPlate(); // range [0- 4]
241 Int_t strip = tofsdigit->GetStrip(); // range [0-14/18/19]
242 Int_t padz = tofsdigit->GetPadz(); // range [0- 1]
243 Int_t padx = tofsdigit->GetPadx(); // range [0-47]
251 //--------------------- QA section ----------------------
252 // in the while, I perform QA
253 Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
256 AliFatal(Form("strange sdigit found %2d %1d %2d %1d %2d", sector, plate, strip, padz, padx));
257 //-------------------------------------------------------
259 //------------------- Dump section ----------------------
261 AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
262 AliInfo("----------------------------------------------------");
264 // ------------------------------------------------------
266 // start loop on number of slots for current sdigit
267 for (Int_t islot = 0; islot < nslot; islot++) {
268 for (Int_t aa=0; aa<4; aa++) digit[aa] = 0; // TOF digit variables
269 for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
271 Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
272 Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
274 //if (tdc>=8192) continue;//AdC
276 tracknum[0]=tofsdigit->GetTrack(islot,0);
277 tracknum[1]=tofsdigit->GetTrack(islot,1);
278 tracknum[2]=tofsdigit->GetTrack(islot,2);
280 // new with placement must be used
281 // adding a TOF digit for each slot
282 TClonesArray &aDigits = *fDigits;
283 Int_t last=fDigits->GetEntriesFast();
284 new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
288 } // end loop on sdigits - end digitizing all collected sdigits
290 //Insert Decalibration
291 AliDebug(2,"in digitizer, create digits");
292 DecalibrateTOFSignal();
295 //---------------------------------------------------------------------
297 void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
299 // Read sdigits for current event and inputFile;
300 // store them into the sdigits container
301 // and update the hit map
302 // SDigits from different files are assumed to
303 // be created with the same simulation parameters.
305 // creating the TClonesArray to store the digits
306 static TClonesArray sdigitsClonesArray("AliTOFSDigit", 1000);
307 sdigitsClonesArray.Clear();
309 // get the treeS from digInput
310 AliRunLoader* rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
313 AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
317 AliLoader* gime = rl->GetLoader("TOFLoader");
320 AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
324 TTree* currentTreeS=gime->TreeS();
325 if (currentTreeS == 0x0)
327 Int_t retval = gime->LoadSDigits();
330 AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
333 currentTreeS=gime->TreeS();
334 if (currentTreeS == 0x0)
336 AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
340 // get the branch TOF inside the treeS
341 TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
342 // check if the branch exist
343 TBranch* tofBranch=currentTreeS->GetBranch("TOF");
346 AliFatal(Form("TOF branch not found for input %d",inputFile));
350 tofBranch->SetAddress(&sdigitsDummyContainer);
352 Int_t nEntries = (Int_t)tofBranch->GetEntries();
354 // Loop through all entries in the tree
357 Int_t vol[5]; // location for a sdigit
358 for (Int_t i=0; i<5; i++) vol[i] = -1;
360 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
363 nbytes += tofBranch->GetEvent(iEntry);
365 // Get the number of sdigits
366 Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
368 for (Int_t k=0; k<ndig; k++) {
369 AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
371 for (Int_t i=0; i<5; i++) vol[i] = -1;
373 // check the sdigit volume
374 vol[0] = tofSdigit->GetSector();
375 vol[1] = tofSdigit->GetPlate();
376 vol[2] = tofSdigit->GetStrip();
377 vol[3] = tofSdigit->GetPadx();
378 vol[4] = tofSdigit->GetPadz();
380 if (fhitMap->TestHit(vol) != kEmpty) {
381 AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
382 sdig->Update(tofSdigit);
386 CollectSDigit(tofSdigit); // collect the current sdigit
387 fhitMap->SetHit(vol); // update the hitmap for location vol
389 } // if (hitMap->TestHit(vol) != kEmpty)
391 } // for (Int_t k=0; k<ndig; k++)
393 } // end loop on entries
398 //_____________________________________________________________________________
399 void AliTOFDigitizer::CollectSDigit(const AliTOFSDigit * const sdigit)
402 // Add a TOF sdigit in container
403 // new with placement must be used
404 TClonesArray &aSDigitsArray = *fSDigitsArray;
405 Int_t last=fSDigitsArray->GetEntriesFast();
406 // make a copy of the current sdigit and
407 // put it into tmp array
408 new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
411 //_____________________________________________________________________________
412 void AliTOFDigitizer::InitDecalibration() const {
414 // Initialize TOF digits decalibration
419 fCalib->CreateCalArrays();
420 fCalib->ReadSimHistoFromCDB("TOF/Calib", -1); // use AliCDBManager's number
421 fCalib->ReadParOfflineFromCDB("TOF/Calib", -1); // use AliCDBManager's number
424 //---------------------------------------------------------------------
425 void AliTOFDigitizer::DecalibrateTOFSignal() {
427 // Decalibrate TOF signals according to OCDB parameters
430 Double_t time=0., tot=0., corr=0.;
431 Int_t deltaBC=0, l0l1=0, tdcBin=0;
433 Int_t detId[5] ={-1,-1,-1,-1,-1};
436 Int_t ndigits = fDigits->GetEntriesFast();
437 // Loop on TOF Digits
438 for (Int_t i=0;i<ndigits;i++){
439 AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
440 detId[0] = dig->GetSector();
441 detId[1] = dig->GetPlate();
442 detId[2] = dig->GetStrip();
443 detId[3] = dig->GetPadz();
444 detId[4] = dig->GetPadx();
445 dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
447 index = AliTOFGeometry::GetIndex(detId); // The channel index
449 // Read Calibration parameters from the CDB
451 time = dig->GetTdc() * AliTOFGeometry::TdcBinWidth(); /* ps */
452 tot = dig->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
453 deltaBC = 0;//dig->GetDeltaBC();
454 l0l1 = 0;//dig->GetL0L1Latency();
457 corr = fCalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
458 AliDebug(2, Form("calibrate index %d: time=%f (ps) tot=%f (ns) deltaBC=%d l0l1=%d timestamp=%d corr=%f (ps)",
459 index, time, tot, deltaBC, l0l1, timestamp, corr));
461 // apply time correction
464 // convert in TDC bins and set digit
465 //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
466 tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
471 AliDebug(1,"Simulating miscalibrated digits");
476 //---------------------------------------------------------------------
478 void AliTOFDigitizer::DecalibrateTOFSignal(){ // Old implementation
480 // Read Calibration parameters from the CDB
482 TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
484 AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
486 // Initialize Quantities to Simulate ToT Spectra
488 TH1F * hToT= fCalib->GetTOFSimToT();
489 Int_t nbins = hToT->GetNbinsX();
490 Float_t delta = hToT->GetBinWidth(1);
491 Float_t maxch = hToT->GetBinLowEdge(nbins)+delta;
492 Float_t minch = hToT->GetBinLowEdge(1);
493 Float_t max=0,min=0; //maximum and minimum value of the distribution
494 Int_t maxbin=0,minbin=0; //maximum and minimum bin of the distribution
496 for (Int_t ii=nbins; ii>0; ii--){
497 if (hToT->GetBinContent(ii)!= 0) {
498 max = maxch - (nbins-ii-1)*delta;
502 for (Int_t j=1; j<nbins; j++){
503 if (hToT->GetBinContent(j)!= 0) {
504 min = minch + (j-1)*delta;
512 Float_t maxToTDistr=hToT->GetMaximum();
514 AliDebug (1, Form(" The minimum ToT = %f", minToT));
515 AliDebug (1, Form(" The maximum ToT = %f", maxToT));
516 AliDebug (1, Form(" The maximum peak in ToT = %f", maxToTDistr));
518 // Loop on TOF Digits
520 Bool_t isToTSimulated=kFALSE;
521 Bool_t misCalibPars=kFALSE;
522 if(hToT->GetEntries()>0)isToTSimulated=kTRUE;
523 Int_t ndigits = fDigits->GetEntriesFast();
524 for (Int_t i=0;i<ndigits;i++){
525 AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
527 detId[0] = dig->GetSector();
528 detId[1] = dig->GetPlate();
529 detId[2] = dig->GetStrip();
530 detId[3] = dig->GetPadz();
531 detId[4] = dig->GetPadx();
532 dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
535 //A realistic ToT Spectrum was found in input,
536 //decalibrated TOF Digits likely to be simulated....
538 Int_t index = AliTOFGeometry::GetIndex(detId); // The channel index
539 AliTOFChannelOffline *calChannelOffline = (AliTOFChannelOffline *)calOffline->At(index); //retrieve the info for time slewing
540 Double_t par[6]; // time slewing parameters
542 //check whether we actually ask for miscalibration
544 for (Int_t j = 0; j<6; j++){
545 par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
546 if(par[j]!=0)misCalibPars=kTRUE;
548 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]));
550 // Now generate Realistic ToT distribution from TestBeam Data.
551 // Tot is in ns, assuming a Matching Window of 10 ns.
558 while (simToT <= triy){
559 trix = gRandom->Rndm(i);
560 triy = gRandom->Rndm(i);
561 trix = (maxToT-minToT)*trix + minToT;
562 triy = maxToTDistr*triy;
563 Int_t binx=hToT->FindBin(trix);
564 simToT=hToT->GetBinContent(binx);
566 // the generated ToT (ns)
567 tToT= (Double_t) trix; // to apply slewing we start from ns..
568 // transform TOF signal in ns
569 AliDebug(2,Form(" The Initial Time (counts): %i: ",dig->GetTdc()));
570 AliDebug(2,Form(" Time before miscalibration (ps) %e: ",dig->GetTdc()*(Double_t)AliTOFGeometry::TdcBinWidth()));
571 // add slewing effect
572 timeCorr=par[0] + tToT*(par[1] +tToT*(par[2] +tToT*(par[3] +tToT*(par[4] +tToT*par[5]))));
573 AliDebug(2,Form(" The Time slewing + delay (ns): %f: ",timeCorr));
574 // add global time shift
577 Double_t timeMis = (Double_t)(dig->GetTdc())*(Double_t)AliTOFGeometry::TdcBinWidth();
578 timeMis = timeMis+timeCorr;
579 AliDebug(2,Form(" The Miscalibrated time (ps): %e: ",timeMis));
581 // now update the digit info
583 Int_t tdcCorr= (Int_t)(timeMis/AliTOFGeometry::TdcBinWidth());
584 AliDebug(2,Form(" Final Time (counts): %i: ",tdcCorr));
585 // Setting Decalibrated Time signal (TDC counts)
586 dig->SetTdc(tdcCorr);
587 // Setting realistic ToT signal (TDC counts)
588 tToT*=1E3; //back to ps
589 Int_t tot=(Int_t)(tToT/AliTOFGeometry::ToTBinWidth());//(factor 1E3 as input ToT is in ns)
591 AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
593 AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr));
594 AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot));
598 // For Data with no Miscalibration, set ToT signal == Adc
599 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)
604 AliDebug(1,"Standard Production, no miscalibrated digits");
607 AliDebug(1,"Standard Production, no miscalibrated digits");
609 AliDebug(1,"Simulating miscalibrated digits");