Line data Source code
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 "TMath.h"
32 : #include "TH1F.h"
33 : #include "TTree.h"
34 : #include "TRandom.h"
35 : #include "TObjArray.h"
36 :
37 : #include "AliLoader.h"
38 : #include "AliLog.h"
39 : #include "AliDigitizationInput.h"
40 : #include "AliRunLoader.h"
41 : #include "AliRun.h"
42 :
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"
52 : #include "AliTOF.h"
53 :
54 : extern TRandom *gRandom;
55 :
56 : extern AliRun *gAlice;
57 :
58 :
59 26 : ClassImp(AliTOFDigitizer)
60 :
61 : //___________________________________________
62 : AliTOFDigitizer::AliTOFDigitizer() :
63 0 : AliDigitizer(),
64 0 : fDigits(new TClonesArray("AliTOFdigit",4000)),
65 0 : fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
66 0 : fhitMap(0x0),
67 0 : fCalib(new AliTOFcalib())
68 0 : {
69 : // Default ctor - don't use it
70 0 : InitDecalibration();
71 0 : }
72 :
73 : //___________________________________________
74 : AliTOFDigitizer::AliTOFDigitizer(AliDigitizationInput* digInput):
75 1 : AliDigitizer(digInput),
76 3 : fDigits(new TClonesArray("AliTOFdigit",4000)),
77 3 : fSDigitsArray(new TClonesArray("AliTOFSDigit",1000)),
78 1 : fhitMap(0x0),
79 3 : fCalib(new AliTOFcalib())
80 5 : {
81 : //ctor with RunDigitizer
82 1 : InitDecalibration();
83 2 : }
84 :
85 : //------------------------------------------------------------------------
86 : AliTOFDigitizer::AliTOFDigitizer(const AliTOFDigitizer &source):
87 0 : AliDigitizer(source),
88 0 : fDigits(source.fDigits),
89 0 : fSDigitsArray(source.fSDigitsArray),
90 0 : fhitMap(source.fhitMap),
91 0 : fCalib(source.fCalib)
92 0 : {
93 : // copy constructor
94 0 : }
95 :
96 : //------------------------------------------------------------------------
97 : AliTOFDigitizer& AliTOFDigitizer::operator=(const AliTOFDigitizer &source)
98 : {
99 : // ass. op.
100 :
101 0 : if (this == &source)
102 0 : return *this;
103 :
104 0 : AliDigitizer::operator=(source);
105 0 : fDigits=source.fDigits;
106 0 : fSDigitsArray=source.fSDigitsArray;
107 0 : fhitMap=source.fhitMap;
108 0 : fCalib=source.fCalib;
109 0 : return *this;
110 :
111 0 : }
112 :
113 : //------------------------------------------------------------------------
114 : AliTOFDigitizer::~AliTOFDigitizer()
115 6 : {
116 : // Destructor
117 2 : delete fCalib;
118 1 : if (fDigits){
119 1 : fDigits->Delete();
120 2 : delete fDigits;
121 1 : fDigits=0x0;
122 1 : }
123 1 : if (fSDigitsArray){
124 1 : fSDigitsArray->Delete();
125 2 : delete fSDigitsArray;
126 1 : fSDigitsArray=0x0;
127 1 : }
128 3 : }
129 :
130 : //---------------------------------------------------------------------
131 :
132 : void AliTOFDigitizer::Digitize(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 16 : AliDebug(1, "");
145 :
146 :
147 : // get the ptr to TOF detector
148 4 : AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
149 :
150 : //Make branches
151 :
152 : const Int_t kSize = 20;
153 4 : char branchname[kSize];
154 4 : snprintf(branchname,kSize,"%s", tof->GetName ());
155 :
156 4 : AliRunLoader* outrl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
157 4 : if (outrl == 0x0)
158 : {
159 0 : AliError("Can not find Run Loader in output folder.");
160 0 : return;
161 : }
162 :
163 4 : AliLoader* outgime = outrl->GetLoader("TOFLoader");
164 4 : if (outgime == 0x0)
165 : {
166 0 : AliError("Can not get TOF Loader from Output Run Loader.");
167 0 : return;
168 : }
169 :
170 4 : TTree* treeD = outgime->TreeD();
171 4 : if (treeD == 0x0)
172 : {
173 4 : outgime->MakeTree("D");
174 4 : treeD = outgime->TreeD();
175 4 : }
176 : //Make branch for digits (to be created in Init())
177 4 : tof->MakeBranchInTree(treeD,branchname,&fDigits,4000);
178 :
179 : // container for all summed sdigits (to be created in Init())
180 : //fSDigitsArray=new TClonesArray("AliTOFSDigit",1000);
181 :
182 : // create hit map (to be created in Init())
183 8 : fhitMap = new AliTOFHitMap(fSDigitsArray);
184 :
185 : // Loop over files to digitize
186 :
187 16 : for (Int_t inputFile=0; inputFile<fDigInput->GetNinputs();
188 4 : inputFile++) {
189 4 : ReadSDigit(inputFile);
190 : }
191 :
192 : // create digits
193 4 : CreateDigits();
194 :
195 : // free used memory for Hit Map in current event
196 8 : delete fhitMap;
197 4 : fSDigitsArray->Clear();
198 :
199 4 : treeD->Fill();
200 :
201 12 : AliDebug(2,"----------------------------------------");
202 12 : AliDebug(1,Form("%d digits have been created", fDigits->GetEntriesFast()));
203 12 : AliDebug(2,"----------------------------------------");
204 :
205 4 : outgime->WriteDigits("OVERWRITE");
206 4 : outgime->UnloadDigits();
207 4 : fDigits->Clear();
208 :
209 8 : }
210 :
211 : //---------------------------------------------------------------------
212 :
213 : void AliTOFDigitizer::CreateDigits()
214 : {
215 : // loop on sdigits container to fill the AliTOFdigit TClonesArray
216 : // start digitizing all the collected sdigits
217 :
218 : Int_t ndump=0; // dump the first ndump created digits for each event
219 :
220 : // get the total number of collected sdigits
221 8 : Int_t ndig = fSDigitsArray->GetEntriesFast();
222 :
223 4 : Int_t vol[5]={-1,-1,-1,-1,-1}; // location for a digit
224 4 : Int_t digit[4] = {0,0,0,0}; // TOF digit variables
225 4 : Int_t tracknum[AliTOFSDigit::kMAXDIGITS]; // contributing tracks for the current slot
226 32 : for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
227 :
228 208 : for (Int_t k = 0; k < ndig; k++) {
229 :
230 1200 : for (Int_t i=0; i<5; i++) vol[i] = -1;
231 :
232 : // Get the information for this digit
233 100 : AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigitsArray->UncheckedAt(k);
234 :
235 100 : Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
236 : // for current sdigit
237 :
238 : // TOF sdigit volumes (always the same for all slots)
239 100 : Int_t sector = tofsdigit->GetSector(); // range [0-17]
240 100 : Int_t plate = tofsdigit->GetPlate(); // range [0- 4]
241 100 : Int_t strip = tofsdigit->GetStrip(); // range [0-14/18/19]
242 100 : Int_t padz = tofsdigit->GetPadz(); // range [0- 1]
243 100 : Int_t padx = tofsdigit->GetPadx(); // range [0-47]
244 :
245 100 : vol[0] = sector;
246 100 : vol[1] = plate;
247 100 : vol[2] = strip;
248 100 : vol[3] = padx;
249 100 : vol[4] = padz;
250 :
251 : //--------------------- QA section ----------------------
252 : // in the while, I perform QA
253 300 : Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
254 :
255 100 : if (isSDigitBad)
256 0 : AliFatal(Form("strange sdigit found %2d %1d %2d %1d %2d", sector, plate, strip, padz, padx));
257 : //-------------------------------------------------------
258 :
259 : //------------------- Dump section ----------------------
260 100 : if (k<ndump) {
261 0 : AliInfo(Form("%2d-th digit: Sector %2d | Plate %1d | Strip %2d | PadZ %1d | PadX %2d ", k, sector, plate, strip, padz, padx));
262 0 : AliInfo("----------------------------------------------------");
263 0 : }
264 : // ------------------------------------------------------
265 :
266 : // start loop on number of slots for current sdigit
267 400 : for (Int_t islot = 0; islot < nslot; islot++) {
268 1000 : for (Int_t aa=0; aa<4; aa++) digit[aa] = 0; // TOF digit variables
269 800 : for (Int_t aa=0; aa<AliTOFSDigit::kMAXDIGITS; aa++) tracknum[aa] = -1;
270 :
271 100 : Int_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
272 100 : Int_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
273 :
274 : //if (tdc>=8192) continue;//AdC
275 :
276 100 : tracknum[0]=tofsdigit->GetTrack(islot,0);
277 100 : tracknum[1]=tofsdigit->GetTrack(islot,1);
278 100 : tracknum[2]=tofsdigit->GetTrack(islot,2);
279 :
280 : // new with placement must be used
281 : // adding a TOF digit for each slot
282 100 : TClonesArray &aDigits = *fDigits;
283 100 : Int_t last=fDigits->GetEntriesFast();
284 100 : new (aDigits[last]) AliTOFdigit(tracknum, vol, digit);
285 :
286 : }
287 :
288 : } // end loop on sdigits - end digitizing all collected sdigits
289 :
290 : //Insert Decalibration
291 12 : AliDebug(2,"in digitizer, create digits");
292 4 : DecalibrateTOFSignal();
293 4 : }
294 :
295 : //---------------------------------------------------------------------
296 :
297 : void AliTOFDigitizer::ReadSDigit(Int_t inputFile )
298 : {
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.
304 :
305 : // creating the TClonesArray to store the digits
306 11 : static TClonesArray sdigitsClonesArray("AliTOFSDigit", 1000);
307 4 : sdigitsClonesArray.Clear();
308 :
309 : // get the treeS from digInput
310 4 : AliRunLoader* rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
311 4 : if (rl == 0x0)
312 : {
313 0 : AliError(Form("Can not find Run Loader in input %d folder.",inputFile));
314 0 : return;
315 : }
316 :
317 4 : AliLoader* gime = rl->GetLoader("TOFLoader");
318 4 : if (gime == 0x0)
319 : {
320 0 : AliError(Form("Can not get TOF Loader from Input %d Run Loader.",inputFile));
321 0 : return;
322 : }
323 :
324 4 : TTree* currentTreeS=gime->TreeS();
325 4 : if (currentTreeS == 0x0)
326 : {
327 1 : Int_t retval = gime->LoadSDigits();
328 1 : if (retval)
329 : {
330 0 : AliError(Form("Error occured while loading S. Digits for Input %d",inputFile));
331 0 : return;
332 : }
333 1 : currentTreeS=gime->TreeS();
334 1 : if (currentTreeS == 0x0)
335 : {
336 0 : AliError(Form("Can not get S. Digits Tree for Input %d",inputFile));
337 0 : return;
338 : }
339 1 : }
340 : // get the branch TOF inside the treeS
341 4 : TClonesArray * sdigitsDummyContainer=&sdigitsClonesArray;
342 : // check if the branch exist
343 4 : TBranch* tofBranch=currentTreeS->GetBranch("TOF");
344 :
345 4 : if(!tofBranch){
346 0 : AliFatal(Form("TOF branch not found for input %d",inputFile));
347 0 : return;
348 : }
349 :
350 4 : tofBranch->SetAddress(&sdigitsDummyContainer);
351 :
352 4 : Int_t nEntries = (Int_t)tofBranch->GetEntries();
353 :
354 : // Loop through all entries in the tree
355 : Int_t nbytes = 0;
356 :
357 4 : Int_t vol[5]; // location for a sdigit
358 48 : for (Int_t i=0; i<5; i++) vol[i] = -1;
359 :
360 16 : for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
361 :
362 : // Import the tree
363 4 : nbytes += tofBranch->GetEvent(iEntry);
364 :
365 : // Get the number of sdigits
366 4 : Int_t ndig = sdigitsDummyContainer->GetEntriesFast();
367 :
368 208 : for (Int_t k=0; k<ndig; k++) {
369 100 : AliTOFSDigit *tofSdigit= (AliTOFSDigit*) sdigitsDummyContainer->UncheckedAt(k);
370 :
371 1200 : for (Int_t i=0; i<5; i++) vol[i] = -1;
372 :
373 : // check the sdigit volume
374 100 : vol[0] = tofSdigit->GetSector();
375 100 : vol[1] = tofSdigit->GetPlate();
376 100 : vol[2] = tofSdigit->GetStrip();
377 100 : vol[3] = tofSdigit->GetPadx();
378 100 : vol[4] = tofSdigit->GetPadz();
379 :
380 100 : if (fhitMap->TestHit(vol) != kEmpty) {
381 0 : AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(fhitMap->GetHit(vol));
382 0 : sdig->Update(tofSdigit);
383 :
384 0 : } else {
385 :
386 100 : CollectSDigit(tofSdigit); // collect the current sdigit
387 100 : fhitMap->SetHit(vol); // update the hitmap for location vol
388 :
389 : } // if (hitMap->TestHit(vol) != kEmpty)
390 :
391 : } // for (Int_t k=0; k<ndig; k++)
392 :
393 : } // end loop on entries
394 :
395 12 : }
396 :
397 :
398 : //_____________________________________________________________________________
399 : void AliTOFDigitizer::CollectSDigit(const AliTOFSDigit * const sdigit)
400 : {
401 : //
402 : // Add a TOF sdigit in container
403 : // new with placement must be used
404 200 : TClonesArray &aSDigitsArray = *fSDigitsArray;
405 100 : Int_t last=fSDigitsArray->GetEntriesFast();
406 : // make a copy of the current sdigit and
407 : // put it into tmp array
408 100 : new (aSDigitsArray[last]) AliTOFSDigit(*sdigit);
409 100 : }
410 :
411 : //_____________________________________________________________________________
412 : void AliTOFDigitizer::InitDecalibration() const {
413 : //
414 : // Initialize TOF digits decalibration
415 : //
416 :
417 2 : fCalib->Init();
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 1 : }
424 : //---------------------------------------------------------------------
425 : void AliTOFDigitizer::DecalibrateTOFSignal() {
426 : //
427 : // Decalibrate TOF signals according to OCDB parameters
428 : //
429 :
430 : Double_t time=0., tot=0., corr=0.;
431 : Int_t deltaBC=0, l0l1=0, tdcBin=0;
432 : Int_t index = -1;
433 8 : Int_t detId[5] ={-1,-1,-1,-1,-1};
434 : UInt_t timestamp=0;
435 :
436 4 : Int_t ndigits = fDigits->GetEntriesFast();
437 : // Loop on TOF Digits
438 208 : for (Int_t i=0;i<ndigits;i++){
439 100 : AliTOFdigit * dig = (AliTOFdigit*)fDigits->At(i);
440 100 : detId[0] = dig->GetSector();
441 100 : detId[1] = dig->GetPlate();
442 100 : detId[2] = dig->GetStrip();
443 100 : detId[3] = dig->GetPadz();
444 100 : detId[4] = dig->GetPadx();
445 100 : dig->SetTdcND(dig->GetTdc()); // save the non decalibrated time
446 :
447 100 : index = AliTOFGeometry::GetIndex(detId); // The channel index
448 :
449 : // Read Calibration parameters from the CDB
450 : // get digit info
451 100 : time = dig->GetTdc() * AliTOFGeometry::TdcBinWidth(); /* ps */
452 100 : tot = dig->GetToT() * AliTOFGeometry::ToTBinWidth() * 1.e-3; /* ns */
453 : deltaBC = 0;//dig->GetDeltaBC();
454 : l0l1 = 0;//dig->GetL0L1Latency();
455 :
456 : // get correction
457 100 : corr = fCalib->GetTimeCorrection(index, tot, deltaBC, l0l1, timestamp); /* ps */
458 300 : 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));
460 :
461 : // apply time correction
462 100 : time += corr;
463 :
464 : // convert in TDC bins and set digit
465 : //tdcBin = (Int_t)(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
466 100 : tdcBin = TMath::Nint(time / AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
467 100 : dig->SetTdc(tdcBin);
468 :
469 : }
470 :
471 12 : AliDebug(1,"Simulating miscalibrated digits");
472 :
473 : return;
474 4 : }
475 :
476 : //---------------------------------------------------------------------
477 : /*
478 : void AliTOFDigitizer::DecalibrateTOFSignal(){ // Old implementation
479 :
480 : // Read Calibration parameters from the CDB
481 :
482 : TObjArray * calOffline= fCalib->GetTOFCalArrayOffline();
483 :
484 : AliDebug(2,Form("Size of array for Offline Calibration = %i",calOffline->GetEntries()));
485 :
486 : // Initialize Quantities to Simulate ToT Spectra
487 :
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
495 :
496 : for (Int_t ii=nbins; ii>0; ii--){
497 : if (hToT->GetBinContent(ii)!= 0) {
498 : max = maxch - (nbins-ii-1)*delta;
499 : maxbin = ii;
500 : break;}
501 : }
502 : for (Int_t j=1; j<nbins; j++){
503 : if (hToT->GetBinContent(j)!= 0) {
504 : min = minch + (j-1)*delta;
505 : minbin = j;
506 : break;}
507 : }
508 :
509 : Float_t maxToT=max;
510 : Float_t minToT=min;
511 :
512 : Float_t maxToTDistr=hToT->GetMaximum();
513 :
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));
517 :
518 : // Loop on TOF Digits
519 :
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);
526 : Int_t detId[5];
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
533 : if(isToTSimulated){
534 :
535 : //A realistic ToT Spectrum was found in input,
536 : //decalibrated TOF Digits likely to be simulated....
537 :
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
541 :
542 : //check whether we actually ask for miscalibration
543 :
544 : for (Int_t j = 0; j<6; j++){
545 : par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
546 : if(par[j]!=0)misCalibPars=kTRUE;
547 : }
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]));
549 :
550 : // Now generate Realistic ToT distribution from TestBeam Data.
551 : // Tot is in ns, assuming a Matching Window of 10 ns.
552 :
553 : Float_t simToT = 0;
554 : Float_t trix = 0;
555 : Float_t triy = 0;
556 : Double_t timeCorr;
557 : Double_t tToT;
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);
565 : }
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
575 : //convert to ps
576 : timeCorr*=1E3;
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));
580 :
581 : // now update the digit info
582 :
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)
590 : dig->SetToT(tot);
591 : AliDebug(2,Form(" Final Time and ToT (counts): %d: , %d:",dig->GetTdc(),dig->GetToT()));
592 : if(tdcCorr<0){
593 : AliWarning (Form(" The bad Slewed Time(TDC counts)= %d ", tdcCorr));
594 : AliWarning(Form(" The bad ToT (TDC counts)= %d ", tot));
595 : }
596 : }
597 : else{
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)
600 : }
601 : }
602 :
603 : if(!isToTSimulated)
604 : AliDebug(1,"Standard Production, no miscalibrated digits");
605 : else
606 : if(!misCalibPars)
607 : AliDebug(1,"Standard Production, no miscalibrated digits");
608 : else
609 : AliDebug(1,"Simulating miscalibrated digits");
610 :
611 : return;
612 : }
613 : */
|