1 /**************************************************************************
2 * Copyright(c) 1998-1999, 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 **************************************************************************/
18 //__________________________________________________________//
20 // This is a TTask that constructs SDigits out of Hits //
21 // A Summable Digits is the "sum" of all hits in a pad //
22 // Detector response has been simulated via the method //
23 // SimulateDetectorResponse //
25 // -- Authors: F. Pierella, A. De Caro //
26 // Use case: see AliTOFhits2sdigits.C macro in the CVS //
27 //__________________________________________________________//
29 #include <TBenchmark.h>
30 #include <TClonesArray.h>
33 #include <TParticle.h>
38 #include "AliLoader.h"
41 #include "AliRunLoader.h"
44 #include "AliTOFcalib.h"
45 #include "AliTOFRecoParam.h"
46 #include "AliTOFGeometry.h"
47 #include "AliTOFHitMap.h"
48 #include "AliTOFhitT0.h"
49 #include "AliTOFhit.h"
50 #include "AliTOFSDigitizer.h"
51 #include "AliTOFSDigit.h"
56 ClassImp(AliTOFSDigitizer)
58 //____________________________________________________________________________
59 AliTOFSDigitizer::AliTOFSDigitizer():
60 TTask("TOFSDigitizer",""),
69 fTimeResolution(100.),
93 fLogChargeSmearing(0),
99 fCalib(new AliTOFcalib())
105 //------------------------------------------------------------------------
106 AliTOFSDigitizer::AliTOFSDigitizer(const AliTOFSDigitizer &source):
116 fTimeResolution(100.),
133 fTimeWalkBoundary(0),
136 fPulseHeightSlope(0),
140 fLogChargeSmearing(0),
142 fAverageTimeFlag(-1),
146 fCalib(new AliTOFcalib())
149 //this->fTOFGeometry=source.fTOFGeometry;
153 //____________________________________________________________________________
154 AliTOFSDigitizer& AliTOFSDigitizer::operator=(const AliTOFSDigitizer &/*source*/)
161 //____________________________________________________________________________
162 AliTOFSDigitizer::AliTOFSDigitizer(const char* HeaderFile, Int_t evNumber1, Int_t nEvents):
163 TTask("TOFSDigitizer",""),
167 fHeadersFile(HeaderFile), // input filename (with hits)
170 fSelectedSector(-1), // by default we sdigitize all sectors
171 fSelectedPlate(-1), // by default we sdigitize all plates in all sectors
172 fTimeResolution(100.),
189 fTimeWalkBoundary(0),
192 fPulseHeightSlope(0),
196 fLogChargeSmearing(0),
198 fAverageTimeFlag(-1),
202 fCalib(new AliTOFcalib())
204 //ctor, reading from input file
206 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data());
208 //File was not opened yet open file and get alirun object
210 file = TFile::Open(fHeadersFile.Data(),"update") ;
211 gAlice = (AliRun *) file->Get("gAlice") ;
214 // add Task to //root/Tasks folder
215 TString evfoldname = AliConfig::GetDefaultEventFolderName();
216 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
218 fRunLoader = AliRunLoader::Open(HeaderFile);//open session and mount on default event folder
219 if (fRunLoader == 0x0)
221 AliFatal("Event is not loaded. Exiting");
226 fRunLoader->CdGAFile();
227 TDirectory *savedir=gDirectory;
228 TFile *in=(TFile*)gFile;
231 // when fTOFGeometry was needed
233 AliWarning("Geometry file is not open default TOF geometry will be used");
234 fTOFGeometry = new AliTOFGeometry();
238 fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
243 if (fRunLoader->TreeE() == 0x0) fRunLoader->LoadHeader();
245 if (evNumber1>=0) fEvent1 = evNumber1;
248 if (nEvents==0) fEvent2 = (Int_t)(fRunLoader->GetNumberOfEvents());
249 else if (nEvents>0) fEvent2 = evNumber1+nEvents;
252 if (!(fEvent2>fEvent1)) {
253 AliError(Form("fEvent2 = %d <= fEvent1 = %d", fEvent2, fEvent1));
256 AliError(Form("Correction: fEvent2 = %d <= fEvent1 = %d", fEvent2, fEvent1));
259 // init parameters for sdigitization
262 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
263 if (fTOFLoader == 0x0)
265 AliFatal("Can not find TOF loader in event. Exiting.");
268 fTOFLoader->PostSDigitizer(this);
272 //____________________________________________________________________________
273 AliTOFSDigitizer::~AliTOFSDigitizer()
276 fTOFLoader->CleanSDigitizer();
278 if (fCalib) delete fCalib;
282 //____________________________________________________________________________
283 void AliTOFSDigitizer::InitParameters()
285 // set parameters for detector simulation
289 //fTimeResolution = 80.; //120.; OLD
290 AliTOFRecoParam *recoParams = (AliTOFRecoParam*)fCalib->ReadRecParFromCDB("TOF/Calib",fRunLoader->GetRunNumber());
291 fTimeResolution = recoParams->GetTimeResolution(); // now from OCDB
292 if (fTimeResolution==0.) {
293 AliWarning("In OCDB found 0ps for TOF time resolution. It is set to 100ps.");
294 fTimeResolution = 100.;
296 AliDebug(1,Form(" TOF time resolution read from OCDB = %f ps",fTimeResolution));
297 fpadefficiency = 0.99 ;
298 //fEdgeEffect = 2 ; // edge effects according to test beam results
299 fEdgeEffect = 1 ; // edge effects according to test beam results
300 // but with fixed time resolution, i.e. fTimeResolution
306 fEffCenter = fpadefficiency;
308 fEff2Boundary = 0.90;
309 fEff3Boundary = 0.08;
310 fAddTRes = 68. ; // \sqrt{2x20^2 + 15^2 + 2x10^2 + 30^2 + 50^2} (p-p)
311 //fAddTRes = 48. ; // \sqrt{2x20^2 + 15^2 + 2x10^2 + 30^2 + 15^2} (Pb-Pb)
312 // 30^2+20^2+40^2+50^2+50^2+50^2 = 10400 ps^2 (very old value)
313 fResCenter = 35. ; //50. ; // OLD
315 fResSlope = 37. ; //40. ; // OLD
316 fTimeWalkCenter = 0. ;
317 fTimeWalkBoundary=0. ;
318 fTimeWalkSlope = 0. ;
320 fPulseHeightSlope=2.0 ;
321 fTimeDelaySlope =0.060;
322 // was fMinimumCharge = TMath::Exp(fPulseHeightSlope*fKparameter/2.);
323 fMinimumCharge = TMath::Exp(-fPulseHeightSlope*fHparameter);
324 fChargeSmearing=0.0 ;
325 fLogChargeSmearing=0.13;
326 fTimeSmearing =0.022;
329 fAdcBin = 0.25; // 1 ADC bin = 0.25 pC (or 0.03 pC)
330 fAdcMean = 50.; // ADC distribution mpv value for Landau (in bins)
331 // it corresponds to a mean value of ~100 bins
332 fAdcRms = 25.; // ADC distribution rms value (in bins)
333 // it corresponds to distribution rms ~50 bins
336 //__________________________________________________________________
337 Double_t TimeWithTail(const Double_t * const x, const Double_t * const par)
339 // sigma - par[0], alpha - par[1], part - par[2]
340 // at x<part*sigma - gauss
341 // at x>part*sigma - TMath::Exp(-x/alpha)
344 if(xx<par[0]*par[2]) {
345 f = TMath::Exp(-xx*xx/(2*par[0]*par[0]));
347 f = TMath::Exp(-(xx-par[0]*par[2])/par[1]-0.5*par[2]*par[2]);
352 //____________________________________________________________________________
353 void AliTOFSDigitizer::Exec(Option_t *verboseOption) {
354 //execute TOF sdigitization
355 if (strstr(verboseOption,"tim") || strstr(verboseOption,"all"))
356 gBenchmark->Start("TOFSDigitizer");
358 if (fEdgeTails) ftail = new TF1("tail",TimeWithTail,-2,2,3);
360 Int_t nselectedHits=0;
361 Int_t ntotalsdigits=0;
362 Int_t ntotalupdates=0;
363 Int_t nnoisesdigits=0;
364 Int_t nsignalsdigits=0;
365 Int_t nHitsFromPrim=0;
366 Int_t nHitsFromSec=0;
367 Int_t nlargeTofDiff=0;
369 Bool_t thereIsNotASelection=(fSelectedSector==-1) && (fSelectedPlate==-1);
371 if (fRunLoader->GetAliRun() == 0x0) fRunLoader->LoadgAlice();
372 gAlice = fRunLoader->GetAliRun();
374 fRunLoader->LoadKinematics();
376 AliTOF *tof = (AliTOF *) gAlice->GetDetector("TOF");
379 AliError("TOF not found");
383 fTOFLoader->LoadHits("read");
384 fTOFLoader->LoadSDigits("recreate");
386 Int_t vol[5]={-1,-1,-1,-1,-1}; // location for a digit
387 Int_t digit[2]={0,0}; // TOF digit variables
389 Int_t nselectedHitsinEv=0;
390 Int_t ntotalsdigitsinEv=0;
391 Int_t ntotalupdatesinEv=0;
392 Int_t nnoisesdigitsinEv=0;
393 Int_t nsignalsdigitsinEv=0;
395 for (Int_t iEvent=fEvent1; iEvent<fEvent2; iEvent++) {
396 //AliInfo(Form("------------------- %s -------------", GetName()));
397 //AliInfo(Form("Sdigitizing event %i", iEvent));
399 fRunLoader->GetEvent(iEvent);
401 TTree *hitTree = fTOFLoader->TreeH();
402 if (!hitTree) return;
404 if (fTOFLoader->TreeS () == 0) fTOFLoader->MakeTree ("S");
406 //Make branch for digits
407 tof->MakeBranch("S");
409 // recreate TClonesArray fSDigits - for backward compatibility
410 if (tof->SDigits() == 0) {
411 tof->CreateSDigitsArray();
413 tof->RecreateSDigitsArray();
416 tof->SetTreeAddress();
418 Int_t version=tof->IsVersion();
424 nsignalsdigitsinEv=0;
428 TClonesArray *tofHitArray = tof->Hits();
431 //AliTOFHitMap *hitMap = new AliTOFHitMap(tof->SDigits(), fTOFGeometry);
432 AliTOFHitMap *hitMap = new AliTOFHitMap(tof->SDigits());
434 TBranch * tofHitsBranch = hitTree->GetBranch("TOF");
436 Int_t ntracks = static_cast<Int_t>(hitTree->GetEntries());
437 for (Int_t track = 0; track < ntracks; track++)
439 gAlice->GetMCApp()->ResetHits();
440 tofHitsBranch->GetEvent(track);
442 AliMC *mcApplication = (AliMC*)gAlice->GetMCApp();
444 particle = (TParticle*)mcApplication->Particle(track);
445 Int_t nhits = tofHitArray->GetEntriesFast();
446 // cleaning all hits of the same track in the same pad volume
447 // it is a rare event, however it happens
449 Int_t previousTrack =-1;
450 Int_t previousSector=-1;
451 Int_t previousPlate =-1;
452 Int_t previousStrip =-1;
453 Int_t previousPadX =-1;
454 Int_t previousPadZ =-1;
456 for (Int_t hit = 0; hit < nhits; hit++) {
457 for (Int_t aa=0; aa<5;aa++) vol[aa]=-1; // location for a digit
458 for (Int_t aa=0; aa<2;aa++) digit[aa]=0; // TOF digit variables
464 // fp: really sorry for this, it is a temporary trick to have
466 if (version<6) { //(version!=6 && version!=7)
467 AliTOFhit *tofHit = (AliTOFhit *) tofHitArray->UncheckedAt(hit);
468 tracknum = tofHit->GetTrack();
469 vol[0] = tofHit->GetSector();
470 vol[1] = tofHit->GetPlate();
471 vol[2] = tofHit->GetStrip();
472 vol[3] = tofHit->GetPadx();
473 vol[4] = tofHit->GetPadz();
474 dxPad = tofHit->GetDx();
475 dzPad = tofHit->GetDz();
476 geantTime = tofHit->GetTof(); // unit [s] // already corrected per event_time smearing
478 AliTOFhitT0 *tofHit = (AliTOFhitT0 *) tofHitArray->UncheckedAt(hit);
479 tracknum = tofHit->GetTrack();
480 vol[0] = tofHit->GetSector();
481 vol[1] = tofHit->GetPlate();
482 vol[2] = tofHit->GetStrip();
483 vol[3] = tofHit->GetPadx();
484 vol[4] = tofHit->GetPadz();
485 dxPad = tofHit->GetDx();
486 dzPad = tofHit->GetDz();
487 geantTime = tofHit->GetTof(); // unit [s] // already corrected per event_time_smearing
490 geantTime *= 1.e+09; // conversion from [s] to [ns]
491 // TOF matching window (~200ns) control
492 if (geantTime>=AliTOFGeometry::MatchingWindow()*1E-3) {
493 AliDebug(2,Form("Time measurement (%f) greater than the matching window (%f)",
494 geantTime, AliTOFGeometry::MatchingWindow()*1E-3));
498 // selection case for sdigitizing only hits in a given plate of a given sector
499 if(thereIsNotASelection || (vol[0]==fSelectedSector && vol[1]==fSelectedPlate)){
501 Bool_t dummy=((tracknum==previousTrack) && (vol[0]==previousSector) && (vol[1]==previousPlate) && (vol[2]==previousStrip));
503 Bool_t isCloneOfThePrevious=dummy && ((vol[3]==previousPadX) && (vol[4]==previousPadZ));
505 Bool_t isNeighOfThePrevious=dummy && ((((vol[3]==previousPadX-1) || (vol[3]==previousPadX+1)) && (vol[4]==previousPadZ)) || ((vol[3]==previousPadX) && ((vol[4]==previousPadZ+1) || (vol[4]==previousPadZ-1))));
507 if(!isCloneOfThePrevious && !isNeighOfThePrevious){
508 // update "previous" values
509 // in fact, we are yet in the future, so the present is past
510 previousTrack=tracknum;
511 previousSector=vol[0];
512 previousPlate=vol[1];
513 previousStrip=vol[2];
519 if (particle->GetFirstMother() < 0) nHitsFromPrim++; // counts hits due to primary particles
521 Float_t xStrip=AliTOFGeometry::XPad()*(vol[3]+0.5-0.5*AliTOFGeometry::NpadX())+dxPad;
522 Float_t zStrip=AliTOFGeometry::ZPad()*(vol[4]+0.5-0.5*AliTOFGeometry::NpadZ())+dzPad;
524 Int_t nActivatedPads = 0, nFiredPads = 0;
525 Bool_t isFired[4] = {kFALSE, kFALSE, kFALSE, kFALSE};
526 Float_t tofAfterSimul[4] = {0., 0., 0., 0.};
527 Float_t qInduced[4] = {0.,0.,0.,0.};
528 Int_t nPlace[4] = {0, 0, 0, 0};
529 Float_t averageTime = 0.;
530 SimulateDetectorResponse(zStrip,xStrip,geantTime,nActivatedPads,nFiredPads,isFired,nPlace,qInduced,tofAfterSimul,averageTime);
532 for(Int_t indexOfPad=0; indexOfPad<nActivatedPads; indexOfPad++) {
533 if(isFired[indexOfPad]){ // the pad has fired
534 Float_t timediff=geantTime-tofAfterSimul[indexOfPad];
536 // TOF matching window (~200ns) control
537 if (tofAfterSimul[indexOfPad]>=AliTOFGeometry::MatchingWindow()*1E-3) {
538 AliDebug(2,Form("Time measurement (%f) greater than the matching window (%f)",
539 tofAfterSimul[indexOfPad], AliTOFGeometry::MatchingWindow()*1E-3));
543 if(timediff>=0.2) nlargeTofDiff++; // greater than 200ps
545 digit[0] = (Int_t) ((tofAfterSimul[indexOfPad]*1.e+03)/AliTOFGeometry::TdcBinWidth()); // TDC bin number (each bin -> 24.4 ps)
547 Float_t landauFactor = gRandom->Landau(fAdcMean, fAdcRms);
548 digit[1] = (Int_t) (qInduced[indexOfPad] * landauFactor); // ADC bins (each bin -> 0.25 (or 0.03) pC)
550 // recalculate the volume only for neighbouring pads
552 (nPlace[indexOfPad]<=AliTOFGeometry::NpadX()) ? vol[4] = 0 : vol[4] = 1;
553 (nPlace[indexOfPad]<=AliTOFGeometry::NpadX()) ? vol[3] = nPlace[indexOfPad] - 1 : vol[3] = nPlace[indexOfPad] - AliTOFGeometry::NpadX() - 1;
555 // check if two sdigit are on the same pad;
556 // in that case we sum the two or more sdigits
557 if (hitMap->TestHit(vol) != kEmpty) {
558 AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(hitMap->GetHit(vol));
559 Int_t tdctime = (Int_t) digit[0];
560 Int_t adccharge = (Int_t) digit[1];
561 sdig->Update(AliTOFGeometry::TdcBinWidth(),tdctime,adccharge,tracknum);
566 tof->AddSDigit(tracknum, vol, digit);
573 nsignalsdigitsinEv++;
578 } // if (hitMap->TestHit(vol) != kEmpty)
579 } // if(isFired[indexOfPad])
580 } // end loop on nActivatedPads
581 } // if(nFiredPads) i.e. if some pads has fired
582 } // close if(!isCloneOfThePrevious)
583 } // close the selection on sector and plate
584 } // end loop on hits for the current track
585 } // end loop on ntracks
589 fTOFLoader->TreeS()->Reset();
590 fTOFLoader->TreeS()->Fill();
591 fTOFLoader->WriteSDigits("OVERWRITE");
593 if (tof->SDigits()) tof->ResetSDigits();
595 if (strstr(verboseOption,"all") || strstr(verboseOption,"partial")) {
596 AliDebug(2,"----------------------------------------");
597 AliDebug(2,Form("After sdigitizing %d hits in event %d", nselectedHitsinEv, iEvent));
598 //" (" << nHitsFromPrim << " from primaries and " << nHitsFromSec << " from secondaries) TOF hits, "
599 AliDebug(1,Form("%d sdigits have been created", ntotalsdigitsinEv));
600 AliDebug(2,Form("(%d due to signals and %d due to border effect)", nsignalsdigitsinEv, nnoisesdigitsinEv));
601 AliDebug(2,Form("%d total updates of the hit map have been performed in current event", ntotalupdatesinEv));
602 AliDebug(2,"----------------------------------------");
605 } //event loop on events
607 fTOFLoader->UnloadSDigits();
608 fTOFLoader->UnloadHits();
609 fRunLoader->UnloadKinematics();
610 //fRunLoader->UnloadgAlice();
618 nHitsFromSec=nselectedHits-nHitsFromPrim;
619 if (strstr(verboseOption,"all") || strstr(verboseOption,"partial")) {
620 AliDebug(2,"----------------------------------------");
621 AliDebug(2,Form("After sdigitizing %d hits in %d events ", nselectedHits, fEvent2-fEvent1));
622 //" (" << nHitsFromPrim << " from primaries and " << nHitsFromSec << " from secondaries) TOF hits, "
623 AliDebug(2,Form("%d sdigits have been created", ntotalsdigits));
624 AliDebug(2,Form("(%d due to signals and %d due to border effect)", nsignalsdigits, nnoisesdigits));
625 AliDebug(2,Form("%d total updates of the hit map have been performed", ntotalupdates));
626 AliDebug(2,Form("in %d cases the time of flight difference is greater than 200 ps", nlargeTofDiff));
627 AliDebug(2,"----------------------------------------");
631 if(strstr(verboseOption,"tim") || strstr(verboseOption,"all")){
632 gBenchmark->Stop("TOFSDigitizer");
633 AliInfo("AliTOFSDigitizer:");
634 AliInfo(Form(" took %f seconds in order to make sdigits "
635 "%f seconds per event", gBenchmark->GetCpuTime("TOFSDigitizer"), gBenchmark->GetCpuTime("TOFSDigitizer")/(fEvent2-fEvent1)));
636 AliInfo(" +++++++++++++++++++++++++++++++++++++++++++++++++++ ");
641 //__________________________________________________________________
642 void AliTOFSDigitizer::Print(Option_t* /*opt*/)const
644 AliInfo(Form(" ------------------- %s ------------- ", GetName()));
647 //__________________________________________________________________
648 void AliTOFSDigitizer::SelectSectorAndPlate(Int_t sector, Int_t plate)
650 //Select sector and plate
651 Bool_t isaWrongSelection=(sector < 0) || (sector >= AliTOFGeometry::NSectors()) || (plate < 0) || (plate >= AliTOFGeometry::NPlates());
652 if(isaWrongSelection){
653 AliError("You have selected an invalid value for sector or plate ");
654 AliError(Form("The correct range for sector is [0,%d]", AliTOFGeometry::NSectors()-1));
655 AliError(Form("The correct range for plate is [0,%d]", AliTOFGeometry::NPlates()-1));
656 AliError("By default we continue sdigitizing all hits in all plates of all sectors");
658 fSelectedSector=sector;
659 fSelectedPlate =plate;
660 AliInfo(Form("SDigitizing only hits in plate %d of the sector %d", fSelectedPlate, fSelectedSector));
664 //__________________________________________________________________
665 void AliTOFSDigitizer::SimulateDetectorResponse(Float_t z0, Float_t x0, Float_t geantTime, Int_t& nActivatedPads, Int_t& nFiredPads, Bool_t* isFired, Int_t* nPlace, Float_t* qInduced, Float_t* tofTime, Float_t& averageTime)
668 // Input: z0, x0 - hit position in the strip system (0,0 - center of the strip), cm
669 // geantTime - time generated by Geant, ns
670 // Output: nActivatedPads - the number of pads activated by the hit (1 || 2 || 4)
671 // nFiredPads - the number of pads fired (really activated) by the hit (nFiredPads <= nActivatedPads)
672 // qInduced[iPad]- charge induced on pad, arb. units
673 // this array is initialized at zero by the caller
674 // tofAfterSimul[iPad] - time calculated with edge effect algorithm, ns
675 // this array is initialized at zero by the caller
676 // averageTime - time given by pad hited by the Geant track taking into account the times (weighted) given by the pads fired for edge effect also.
677 // The weight is given by the qInduced[iPad]/qCenterPad
678 // this variable is initialized at zero by the caller
679 // nPlace[iPad] - the number of the pad place, iPad = 0, 1, 2, 3
680 // this variable is initialized at zero by the caller
682 // Description of used variables:
683 // eff[iPad] - efficiency of the pad
684 // res[iPad] - resolution of the pad, ns
685 // timeWalk[iPad] - time walk of the pad, ns
686 // timeDelay[iPad] - time delay for neighbouring pad to hited pad, ns
687 // PadId[iPad] - Pad Identifier
688 // E | F --> PadId[iPad] = 5 | 6
689 // A | B --> PadId[iPad] = 1 | 2
690 // C | D --> PadId[iPad] = 3 | 4
691 // nTail[iPad] - the tail number, = 1 for tailA, = 2 for tailB
692 // qCenterPad - charge extimated for each pad, arb. units
693 // weightsSum - sum of weights extimated for each pad fired, arb. units
695 const Float_t kSigmaForTail[2] = {AliTOFGeometry::SigmaForTail1(),AliTOFGeometry::SigmaForTail2()}; //for tail
696 Int_t iz = 0, ix = 0;
697 Float_t dX = 0., dZ = 0., x = 0., z = 0.;
698 Float_t h = fHparameter, h2 = fH2parameter, k = fKparameter, k2 = fK2parameter;
699 Float_t effX = 0., effZ = 0., resX = 0., resZ = 0., timeWalkX = 0., timeWalkZ = 0.;
700 Float_t logOfqInd = 0.;
701 Float_t weightsSum = 0.;
702 Int_t nTail[4] = {0,0,0,0};
703 Int_t padId[4] = {0,0,0,0};
704 Float_t eff[4] = {0.,0.,0.,0.};
705 Float_t res[4] = {0.,0.,0.,0.};
706 // Float_t qCenterPad = fMinimumCharge * fMinimumCharge;
707 Float_t qCenterPad = 1.;
708 Float_t timeWalk[4] = {0.,0.,0.,0.};
709 Float_t timeDelay[4] = {0.,0.,0.,0.};
714 (z0 <= 0) ? iz = 0 : iz = 1;
715 dZ = z0 + (0.5 * AliTOFGeometry::NpadZ() - iz - 0.5) * AliTOFGeometry::ZPad(); // hit position in the pad frame, (0,0) - center of the pad
716 z = 0.5 * AliTOFGeometry::ZPad() - TMath::Abs(dZ); // variable for eff., res. and timeWalk. functions
717 iz++; // z row: 1, ..., AliTOFGeometry::NpadZ = 2
718 ix = (Int_t)((x0 + 0.5 * AliTOFGeometry::NpadX() * AliTOFGeometry::XPad()) / AliTOFGeometry::XPad());
719 dX = x0 + (0.5 * AliTOFGeometry::NpadX() - ix - 0.5) * AliTOFGeometry::XPad(); // hit position in the pad frame, (0,0) - center of the pad
720 x = 0.5 * AliTOFGeometry::XPad() - TMath::Abs(dX); // variable for eff., res. and timeWalk. functions;
721 ix++; // x row: 1, ..., AliTOFGeometry::NpadX = 48
725 nPlace[nActivatedPads-1] = (iz - 1) * AliTOFGeometry::NpadX() + ix;
726 qInduced[nActivatedPads-1] = qCenterPad;
727 padId[nActivatedPads-1] = 1;
729 if (fEdgeEffect == 0) {
730 eff[nActivatedPads-1] = fEffCenter;
731 if (gRandom->Rndm() < eff[nActivatedPads-1]) {
733 res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + fResCenter * fResCenter); // ns
734 isFired[nActivatedPads-1] = kTRUE;
735 tofTime[nActivatedPads-1] = gRandom->Gaus(geantTime + fTimeWalkCenter, res[0]);
736 averageTime = tofTime[nActivatedPads-1];
738 } else { // if (fEdgeEffet!=0)
742 effZ = fEffBoundary + (fEff2Boundary - fEffBoundary) * z / h2;
744 effZ = fEff2Boundary + (fEffCenter - fEff2Boundary) * (z - h2) / (h - h2);
747 resZ = fTimeResolution;
748 else if (fEdgeEffect==2)
749 resZ = fResBoundary + (fResCenter - fResBoundary) * z / h;
750 timeWalkZ = fTimeWalkBoundary + (fTimeWalkCenter - fTimeWalkBoundary) * z / h;
751 nTail[nActivatedPads-1] = 1;
755 resZ = fTimeResolution;
756 else if (fEdgeEffect==2)
758 timeWalkZ = fTimeWalkCenter;
763 effX = fEffBoundary + (fEff2Boundary - fEffBoundary) * x / h2;
765 effX = fEff2Boundary + (fEffCenter - fEff2Boundary) * (x - h2) / (h - h2);
768 resX = fTimeResolution;
769 else if (fEdgeEffect==2)
770 resX = fResBoundary + (fResCenter - fResBoundary) * x / h;
771 timeWalkX = fTimeWalkBoundary + (fTimeWalkCenter - fTimeWalkBoundary) * x / h;
772 nTail[nActivatedPads-1] = 1;
776 resX = fTimeResolution;
777 else if (fEdgeEffect==2)
779 timeWalkX = fTimeWalkCenter;
782 (effZ<effX) ? eff[nActivatedPads-1] = effZ : eff[nActivatedPads-1] = effX;
783 (resZ<resX) ? res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resX * resX) : res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resZ * resZ); // ns
784 (timeWalkZ<timeWalkX) ? timeWalk[nActivatedPads-1] = 0.001 * timeWalkZ : timeWalk[nActivatedPads-1] = 0.001 * timeWalkX; // ns
789 effZ = fEffBoundary - (fEffBoundary - fEff3Boundary) * (z / k2);
791 effZ = fEff3Boundary * (k - z) / (k - k2);
794 resZ = fTimeResolution;
795 else if (fEdgeEffect==2)
796 resZ = fResBoundary + fResSlope * z / k;
797 timeWalkZ = fTimeWalkBoundary + fTimeWalkSlope * z / k;
800 if( (iz == 1 && dZ > 0) || (iz == 2 && dZ < 0) ) {
802 nPlace[nActivatedPads-1] = nPlace[0] + (3 - 2 * iz) * AliTOFGeometry::NpadX();
803 eff[nActivatedPads-1] = effZ;
804 res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resZ * resZ); // ns
805 timeWalk[nActivatedPads-1] = 0.001 * timeWalkZ; // ns
806 nTail[nActivatedPads-1] = 2;
807 if (fTimeDelayFlag) {
808 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * z / 2.);
809 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * z / 2.);
810 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * z);
811 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * z, fLogChargeSmearing);
812 timeDelay[nActivatedPads-1] = gRandom->Gaus(-fTimeDelaySlope * logOfqInd, fTimeSmearing);
814 timeDelay[nActivatedPads-1] = 0.;
816 padId[nActivatedPads-1] = 2;
821 ////// Pad C, D, E, F:
823 effX = fEffBoundary - (fEffBoundary - fEff3Boundary) * (x / k2);
825 effX = fEff3Boundary * (k - x) / (k - k2);
828 resX = fTimeResolution;
829 else if (fEdgeEffect==2)
830 resX = fResBoundary + fResSlope*x/k;
831 timeWalkX = fTimeWalkBoundary + fTimeWalkSlope*x/k;
835 if(ix > 1 && dX < 0) {
837 nPlace[nActivatedPads-1] = nPlace[0] - 1;
838 eff[nActivatedPads-1] = effX;
839 res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resX * resX); // ns
840 timeWalk[nActivatedPads-1] = 0.001 * timeWalkX; // ns
841 nTail[nActivatedPads-1] = 2;
842 if (fTimeDelayFlag) {
843 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * x / 2.);
844 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * x / 2.);
845 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * x);
846 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * x, fLogChargeSmearing);
847 timeDelay[nActivatedPads-1] = gRandom->Gaus(-fTimeDelaySlope * logOfqInd, fTimeSmearing);
849 timeDelay[nActivatedPads-1] = 0.;
851 padId[nActivatedPads-1] = 3;
855 if( (iz == 1 && dZ > 0) || (iz == 2 && dZ < 0) ) {
857 nPlace[nActivatedPads-1] = nPlace[0] + (3 - 2 * iz) * AliTOFGeometry::NpadX() - 1;
858 eff[nActivatedPads-1] = effX * effZ;
859 (resZ<resX) ? res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resX * resX) : res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resZ * resZ); // ns
860 (timeWalkZ<timeWalkX) ? timeWalk[nActivatedPads-1] = 0.001 * timeWalkZ : timeWalk[nActivatedPads-1] = 0.001 * timeWalkX; // ns
862 nTail[nActivatedPads-1] = 2;
863 if (fTimeDelayFlag) {
864 if (TMath::Abs(x) < TMath::Abs(z)) {
865 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * z / 2.);
866 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * z / 2.);
867 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * z);
868 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * z, fLogChargeSmearing);
870 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * x / 2.);
871 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * x / 2.);
872 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * x);
873 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * x, fLogChargeSmearing);
875 timeDelay[nActivatedPads-1] = gRandom->Gaus(-fTimeDelaySlope * logOfqInd, fTimeSmearing);
877 timeDelay[nActivatedPads-1] = 0.;
879 padId[nActivatedPads-1] = 4;
885 if(ix < AliTOFGeometry::NpadX() && dX > 0) {
887 nPlace[nActivatedPads-1] = nPlace[0] + 1;
888 eff[nActivatedPads-1] = effX;
889 res[nActivatedPads-1] = 0.001 * (TMath::Sqrt(fAddTRes*fAddTRes + resX * resX)); // ns
890 timeWalk[nActivatedPads-1] = 0.001 * timeWalkX; // ns
891 nTail[nActivatedPads-1] = 2;
892 if (fTimeDelayFlag) {
893 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * x / 2.);
894 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * x / 2.);
895 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * x);
896 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * x, fLogChargeSmearing);
897 timeDelay[nActivatedPads-1] = gRandom->Gaus(-fTimeDelaySlope * logOfqInd, fTimeSmearing);
899 timeDelay[nActivatedPads-1] = 0.;
901 padId[nActivatedPads-1] = 5;
906 if( (iz == 1 && dZ > 0) || (iz == 2 && dZ < 0) ) {
908 nPlace[nActivatedPads - 1] = nPlace[0] + (3 - 2 * iz) * AliTOFGeometry::NpadX() + 1;
909 eff[nActivatedPads - 1] = effX * effZ;
910 (resZ<resX) ? res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resX * resX) : res[nActivatedPads-1] = 0.001 * TMath::Sqrt(fAddTRes*fAddTRes + resZ * resZ); // ns
911 (timeWalkZ<timeWalkX) ? timeWalk[nActivatedPads-1] = 0.001 * timeWalkZ : timeWalk[nActivatedPads-1] = 0.001*timeWalkX; // ns
912 nTail[nActivatedPads-1] = 2;
913 if (fTimeDelayFlag) {
914 if (TMath::Abs(x) < TMath::Abs(z)) {
915 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * z / 2.);
916 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * z / 2.);
917 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * z);
918 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * z, fLogChargeSmearing);
920 // qInduced[0] = fMinimumCharge * TMath::Exp(fPulseHeightSlope * x / 2.);
921 // qInduced[nActivatedPads-1] = fMinimumCharge * TMath::Exp(-fPulseHeightSlope * x / 2.);
922 qInduced[nActivatedPads-1] = TMath::Exp(-fPulseHeightSlope * x);
923 logOfqInd = gRandom->Gaus(-fPulseHeightSlope * x, fLogChargeSmearing);
925 timeDelay[nActivatedPads-1] = gRandom->Gaus(-fTimeDelaySlope * logOfqInd, fTimeSmearing);
927 timeDelay[nActivatedPads-1] = 0.;
929 padId[nActivatedPads-1] = 6;
936 for (Int_t iPad = 0; iPad < nActivatedPads; iPad++) {
937 if (fEdgeEffect==2 && res[iPad] < fTimeResolution) res[iPad] = fTimeResolution;
938 if(gRandom->Rndm() < eff[iPad]) {
939 isFired[iPad] = kTRUE;
942 if(nTail[iPad] == 0) {
943 tofTime[iPad] = gRandom->Gaus(geantTime + timeWalk[iPad] + timeDelay[iPad], res[iPad]);
945 ftail->SetParameters(res[iPad], 2. * res[iPad], kSigmaForTail[nTail[iPad]-1]);
946 Double_t timeAB = ftail->GetRandom();
947 tofTime[iPad] = geantTime + timeWalk[iPad] + timeDelay[iPad] + timeAB;
950 AliDebug(1,Form(" ----------------- TOF time resolution = %f",res[iPad]));
951 tofTime[iPad] = gRandom->Gaus(geantTime + timeWalk[iPad] + timeDelay[iPad], res[iPad]);
953 if (fAverageTimeFlag) {
954 averageTime += tofTime[iPad] * qInduced[iPad];
955 weightsSum += qInduced[iPad];
957 averageTime += tofTime[iPad];
962 if (weightsSum!=0) averageTime /= weightsSum;
963 } // end else (fEdgeEffect != 0)
966 //__________________________________________________________________
967 void AliTOFSDigitizer::PrintParameters()const
970 // Print parameters used for sdigitization
972 AliInfo(Form(" ------------------- %s -------------", GetName()));
973 AliInfo(" Parameters used for TOF SDigitization ");
974 // Printing the parameters
976 AliInfo(Form(" Number of events: %i ", (fEvent2-fEvent1)));
977 AliInfo(Form(" from event %i to event %i", fEvent1, (fEvent2-1)));
978 AliInfo(Form(" Time Resolution (ps) %f Pad Efficiency: %f ", fTimeResolution, fpadefficiency));
979 AliInfo(Form(" Edge Effect option: %d", fEdgeEffect));
981 AliInfo(" Boundary Effect Simulation Parameters ");
982 AliInfo(Form(" Hparameter: %f H2parameter: %f Kparameter: %f K2parameter: %f", fHparameter, fH2parameter, fKparameter, fK2parameter));
983 AliInfo(Form(" Efficiency in the central region of the pad: %f", fEffCenter));
984 AliInfo(Form(" Efficiency at the boundary region of the pad: %f", fEffBoundary));
985 AliInfo(Form(" Efficiency value at H2parameter %f", fEff2Boundary));
986 AliInfo(Form(" Efficiency value at K2parameter %f", fEff3Boundary));
987 AliInfo(Form(" Resolution (ps) in the central region of the pad: %f", fResCenter));
988 AliInfo(Form(" Resolution (ps) at the boundary of the pad : %f", fResBoundary));
989 AliInfo(Form(" Slope (ps/K) for neighbouring pad : %f", fResSlope));
990 AliInfo(Form(" Time walk (ps) in the central region of the pad : %f", fTimeWalkCenter));
991 AliInfo(Form(" Time walk (ps) at the boundary of the pad : %f", fTimeWalkBoundary));
992 AliInfo(Form(" Slope (ps/K) for neighbouring pad : %f", fTimeWalkSlope));
993 AliInfo(" Pulse Heigth Simulation Parameters ");
994 AliInfo(Form(" Flag for delay due to the PulseHeightEffect : %d", fTimeDelayFlag));
995 AliInfo(Form(" Pulse Height Slope : %f", fPulseHeightSlope));
996 AliInfo(Form(" Time Delay Slope : %f", fTimeDelaySlope));
997 AliInfo(Form(" Minimum charge amount which could be induced : %f", fMinimumCharge));
998 AliInfo(Form(" Smearing in charge in (q1/q2) vs x plot : %f", fChargeSmearing));
999 AliInfo(Form(" Smearing in log of charge ratio : %f", fLogChargeSmearing));
1000 AliInfo(Form(" Smearing in time in time vs log(q1/q2) plot : %f", fTimeSmearing));
1001 AliInfo(Form(" Flag for average time : %d", fAverageTimeFlag));
1002 AliInfo(Form(" Edge tails option : %d", fEdgeTails));