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 // Creates and handles digits from TRD hits //
21 // Author: C. Blume (C.Blume@gsi.de) //
23 // The following effects are included: //
26 // - Gas gain including fluctuations //
27 // - Pad-response (simple Gaussian approximation) //
29 // - Electronics noise //
30 // - Electronics gain //
33 // The corresponding parameter can be adjusted via the various //
34 // Set-functions. If these parameters are not explicitly set, default //
35 // values are used (see Init-function). //
36 // As an example on how to use this class to produce digits from hits //
37 // have a look at the macro hits2digits.C //
38 // The production of summable digits is demonstrated in hits2sdigits.C //
39 // and the subsequent conversion of the s-digits into normal digits is //
40 // explained in sdigits2digits.C. //
42 ///////////////////////////////////////////////////////////////////////////////
57 #include "AliRunLoader.h"
58 #include "AliLoader.h"
59 #include "AliConfig.h"
61 #include "AliRunDigitizer.h"
62 #include "AliRunLoader.h"
63 #include "AliLoader.h"
66 #include "AliTRDhit.h"
67 #include "AliTRDdigitizer.h"
68 #include "AliTRDdataArrayI.h"
69 #include "AliTRDdataArrayF.h"
70 #include "AliTRDsegmentArray.h"
71 #include "AliTRDdigitsManager.h"
72 #include "AliTRDgeometry.h"
73 #include "AliTRDparameter.h"
75 ClassImp(AliTRDdigitizer)
77 //_____________________________________________________________________________
78 AliTRDdigitizer::AliTRDdigitizer()
81 // AliTRDdigitizer default constructor
85 fSDigitsManagerList = 0;
96 fMergeSignalOnly = kFALSE;
102 //_____________________________________________________________________________
103 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
104 :AliDigitizer(name,title)
107 // AliTRDdigitizer constructor
110 fDigitsManager = NULL;
111 fSDigitsManager = NULL;
112 fSDigitsManagerList = NULL;
116 //NewIO: These data members probably are not needed anymore
118 fSDigitsManagerList = 0;
129 fMergeSignalOnly = kFALSE;
133 // For the summable digits
134 fSDigitsScale = 100.;
138 //_____________________________________________________________________________
139 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
140 , const Text_t *name, const Text_t *title)
141 :AliDigitizer(manager,name,title)
144 // AliTRDdigitizer constructor
148 fSDigitsManagerList = 0;
158 fMergeSignalOnly = kFALSE;
162 // For the summable digits
163 fSDigitsScale = 100.;
167 //_____________________________________________________________________________
168 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
169 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
172 // AliTRDdigitizer constructor
178 fSDigitsManagerList = 0;
189 fMergeSignalOnly = kFALSE;
193 // For the summable digits
194 fSDigitsScale = 100.;
198 //_____________________________________________________________________________
199 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
202 // AliTRDdigitizer copy constructor
205 ((AliTRDdigitizer &) d).Copy(*this);
209 //_____________________________________________________________________________
210 AliTRDdigitizer::~AliTRDdigitizer()
213 // AliTRDdigitizer destructor
217 if (fDigitsManager) {
218 delete fDigitsManager;
222 if (fSDigitsManager) {
223 delete fSDigitsManager;
227 if (fSDigitsManagerList) {
228 delete fSDigitsManagerList;
229 fSDigitsManagerList = 0;
239 //_____________________________________________________________________________
240 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
243 // Assignment operator
246 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
251 //_____________________________________________________________________________
252 void AliTRDdigitizer::Copy(TObject &d)
258 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
259 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
260 ((AliTRDdigitizer &) d).fDigitsManager = 0;
261 ((AliTRDdigitizer &) d).fTRD = 0;
262 ((AliTRDdigitizer &) d).fGeo = 0;
263 ((AliTRDdigitizer &) d).fMasks = 0;
264 ((AliTRDdigitizer &) d).fEvent = 0;
265 ((AliTRDdigitizer &) d).fPar = 0;
266 ((AliTRDdigitizer &) d).fCompress = fCompress;
267 ((AliTRDdigitizer &) d).fDebug = fDebug ;
268 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
269 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
270 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
271 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
272 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
276 //_____________________________________________________________________________
277 void AliTRDdigitizer::Exec(Option_t* option)
280 // Executes the merging
285 AliTRDdigitsManager *sdigitsManager;
287 TString optionString = option;
288 if (optionString.Contains("deb")) {
290 if (optionString.Contains("2")) {
293 printf("<AliTRDdigitizer::Exec> ");
294 printf("Called with debug option %d\n",fDebug);
297 // The AliRoot file is already connected by the manager
303 printf("<AliTRDdigitizer::Exec> ");
304 printf("AliRun object found on file.\n");
308 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
310 gAlice = inrl->GetAliRun();
313 printf("<AliTRDdigitizer::Exec> ");
314 printf("Could not find AliRun object.\n");
319 Int_t nInput = fManager->GetNinputs();
320 fMasks = new Int_t[nInput];
321 for (iInput = 0; iInput < nInput; iInput++) {
322 fMasks[iInput] = fManager->GetMask(iInput);
328 for (iInput = 0; iInput < nInput; iInput++) {
331 printf("<AliTRDdigitizer::Exec> ");
332 printf("Add input stream %d\n",iInput);
335 // check if the input tree exists
336 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
337 AliLoader* gime = inrl->GetLoader("TRDLoader");
339 TTree * treees = gime->TreeS();
342 if (gime->LoadSDigits())
344 Error("Exec","Error Occured while loading S. Digits for input %d.",iInput);
347 treees = gime->TreeS();
351 printf("<AliTRDdigitizer::Exec> ");
352 printf("Input stream %d does not exist\n",iInput);
356 // Read the s-digits via digits manager
357 sdigitsManager = new AliTRDdigitsManager();
358 sdigitsManager->SetDebug(fDebug);
359 sdigitsManager->SetSDigits(kTRUE);
361 AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
362 AliLoader* gimme = rl->GetLoader("TRDLoader");
363 if (!gimme->TreeS()) gimme->LoadSDigits();
364 sdigitsManager->ReadDigits(gimme->TreeS());
366 // Add the s-digits to the input list
367 AddSDigitsManager(sdigitsManager);
371 // Convert the s-digits to normal digits
373 printf("<AliTRDdigitizer::Exec> ");
374 printf("Do the conversion\n");
380 printf("<AliTRDdigitizer::Exec> ");
381 printf("Write the digits\n");
384 AliRunLoader* orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
385 AliLoader* ogime = orl->GetLoader("TRDLoader");
387 fDigitsManager->MakeBranch(ogime->TreeD());
389 fDigitsManager->WriteDigits();
392 printf("<AliTRDdigitizer::Exec> ");
396 DeleteSDigitsManager();
400 //_____________________________________________________________________________
401 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
404 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
407 // Connect the AliRoot file containing Geometry, Kine, and Hits
409 fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,"UPDATE");
413 Error("Open","Can not open session for file %s.",file);
417 fRunLoader->LoadgAlice();
418 gAlice = fRunLoader->GetAliRun();
422 printf("<AliTRDdigitizer::Open> ");
423 printf("AliRun object found on file.\n");
427 printf("<AliTRDdigitizer::Open> ");
428 printf("Could not find AliRun object.\n");
434 // Import the Trees for the event nEvent in the file
435 fRunLoader->GetEvent(fEvent);
437 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
440 Error("Open","Can not get TRD loader from Run Loader");
444 if (InitDetector()) {
448 //if we produce SDigits
449 tree = loader->TreeS();
452 loader->MakeTree("S");
453 tree = loader->TreeS();
457 {//if we produce Digits
458 tree = loader->TreeD();
461 loader->MakeTree("D");
462 tree = loader->TreeD();
465 return MakeBranch(tree);
473 //_____________________________________________________________________________
474 Bool_t AliTRDdigitizer::InitDetector()
477 // Sets the pointer to the TRD detector and the geometry
480 // Get the pointer to the detector class and check for version 1
481 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
483 printf("<AliTRDdigitizer::InitDetector> ");
484 printf("No TRD module found\n");
487 if (fTRD->IsVersion() != 1) {
488 printf("<AliTRDdigitizer::InitDetector> ");
489 printf("TRD must be version 1 (slow simulator).\n");
494 fGeo = fTRD->GetGeometry();
496 printf("<AliTRDdigitizer::InitDetector> ");
497 printf("Geometry version %d\n",fGeo->IsVersion());
500 // Create a digits manager
501 fDigitsManager = new AliTRDdigitsManager();
502 fDigitsManager->SetSDigits(fSDigits);
503 fDigitsManager->CreateArrays();
504 fDigitsManager->SetEvent(fEvent);
505 fDigitsManager->SetDebug(fDebug);
507 // The list for the input s-digits manager to be merged
508 fSDigitsManagerList = new TList();
514 //_____________________________________________________________________________
515 Bool_t AliTRDdigitizer::MakeBranch(TTree* tree) const
518 // Create the branches for the digits array
521 return fDigitsManager->MakeBranch(tree);
525 //_____________________________________________________________________________
526 Bool_t AliTRDdigitizer::MakeDigits()
532 ///////////////////////////////////////////////////////////////
534 ///////////////////////////////////////////////////////////////
536 // Converts number of electrons to fC
537 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
539 ///////////////////////////////////////////////////////////////
541 // Number of pads included in the pad response
542 const Int_t kNpad = 3;
544 // Number of track dictionary arrays
545 const Int_t kNDict = AliTRDdigitsManager::kNDict;
547 // Half the width of the amplification region
548 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
550 Int_t iRow, iCol, iTime, iPad;
554 Int_t totalSizeDigits = 0;
555 Int_t totalSizeDict0 = 0;
556 Int_t totalSizeDict1 = 0;
557 Int_t totalSizeDict2 = 0;
559 Int_t timeTRDbeg = 0;
560 Int_t timeTRDend = 1;
565 Float_t padSignal[kNpad];
566 Float_t signalOld[kNpad];
568 AliTRDdataArrayF *signals = 0;
569 AliTRDdataArrayI *digits = 0;
570 AliTRDdataArrayI *dictionary[kNDict];
572 // Create a default parameter class if none is defined
574 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
576 printf("<AliTRDdigitizer::MakeDigits> ");
577 printf("Create the default parameter object\n");
581 // Create a container for the amplitudes
582 AliTRDsegmentArray *signalsArray
583 = new AliTRDsegmentArray("AliTRDdataArrayF"
584 ,AliTRDgeometry::Ndet());
587 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
588 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
590 printf("<AliTRDdigitizer::MakeDigits> ");
591 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
595 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
598 printf("<AliTRDdigitizer::MakeDigits> ");
599 printf("No geometry defined\n");
604 printf("<AliTRDdigitizer::MakeDigits> ");
605 printf("Start creating digits.\n");
608 AliLoader* gimme = fRunLoader->GetLoader("TRDLoader");
609 if (!gimme->TreeH()) gimme->LoadHits();
610 TTree* hitTree = gimme->TreeH();
613 Error("MakeDigits","Can not get TreeH");
616 fTRD->SetTreeAddress();
618 // Get the number of entries in the hit tree
619 // (Number of primary particles creating a hit somewhere)
622 nTrack = (Int_t) hitTree->GetEntries();
624 printf("<AliTRDdigitizer::MakeDigits> ");
625 printf("Found %d primary particles\n",nTrack);
629 Int_t detectorOld = -1;
632 // Loop through all entries in the tree
633 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
637 nBytes += hitTree->GetEvent(iTrack);
640 // Loop through the TRD hits
642 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
651 Float_t q = hit->GetCharge();
652 Int_t track = hit->Track();
653 Int_t detector = hit->GetDetector();
654 Int_t plane = fGeo->GetPlane(detector);
655 Int_t sector = fGeo->GetSector(detector);
656 Int_t chamber = fGeo->GetChamber(detector);
657 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
658 Int_t nColMax = fPar->GetColMax(plane);
659 Int_t nTimeMax = fPar->GetTimeMax();
660 Int_t nTimeBefore = fPar->GetTimeBefore();
661 Int_t nTimeAfter = fPar->GetTimeAfter();
662 Int_t nTimeTotal = fPar->GetTimeTotal();
663 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
664 Float_t col0 = fPar->GetCol0(plane);
665 Float_t time0 = fPar->GetTime0(plane);
666 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
667 Float_t colPadSize = fPar->GetColPadSize(plane);
668 Float_t timeBinSize = fPar->GetTimeBinSize();
669 Float_t divideRow = 1.0 / rowPadSize;
670 Float_t divideCol = 1.0 / colPadSize;
671 Float_t divideTime = 1.0 / timeBinSize;
674 printf("Analyze hit no. %d ",iHit);
675 printf("-----------------------------------------------------------\n");
677 printf("plane = %d, sector = %d, chamber = %d\n"
678 ,plane,sector,chamber);
679 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
680 ,nRowMax,nColMax,nTimeMax);
681 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
682 ,nTimeBefore,nTimeAfter,nTimeTotal);
683 printf("row0 = %f, col0 = %f, time0 = %f\n"
685 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
686 ,rowPadSize,colPadSize,timeBinSize);
689 // Don't analyze test hits and switched off detectors
690 if ((CheckDetector(plane,chamber,sector)) &&
691 (((Int_t) q) != 0)) {
693 if (detector != detectorOld) {
696 printf("<AliTRDdigitizer::MakeDigits> ");
697 printf("Get new container. New det = %d, Old det = %d\n"
698 ,detector,detectorOld);
700 // Compress the old one if enabled
701 if ((fCompress) && (detectorOld > -1)) {
703 printf("<AliTRDdigitizer::MakeDigits> ");
704 printf("Compress the old container ...");
706 signals->Compress(1,0);
707 for (iDict = 0; iDict < kNDict; iDict++) {
708 dictionary[iDict]->Compress(1,0);
710 if (fDebug > 1) printf("done\n");
712 // Get the new container
713 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
714 if (signals->GetNtime() == 0) {
715 // Allocate a new one if not yet existing
717 printf("<AliTRDdigitizer::MakeDigits> ");
718 printf("Allocate a new container ... ");
720 signals->Allocate(nRowMax,nColMax,nTimeTotal);
722 else if (fSimpleSim) {
723 // Clear an old one for the simple simulation
725 printf("<AliTRDdigitizer::MakeDigits> ");
726 printf("Clear a old container ... ");
731 // Expand an existing one
734 printf("<AliTRDdigitizer::MakeDigits> ");
735 printf("Expand an existing container ... ");
740 // The same for the dictionary
742 for (iDict = 0; iDict < kNDict; iDict++) {
743 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
744 if (dictionary[iDict]->GetNtime() == 0) {
745 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
748 if (fCompress) dictionary[iDict]->Expand();
752 if (fDebug > 1) printf("done\n");
753 detectorOld = detector;
756 // Rotate the sectors on top of each other
763 fGeo->Rotate(detector,pos,rot);
766 // The driftlength. It is negative if the hit is in the
767 // amplification region.
768 Float_t driftlength = time0 - rot[0];
770 // Take also the drift in the amplification region into account
771 // The drift length is at the moment still the same, regardless of
772 // the position relativ to the wire. This non-isochronity needs still
773 // to be implemented.
774 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
775 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
777 // Loop over all electrons of this hit
778 // TR photons produce hits with negative charge
779 Int_t nEl = ((Int_t) TMath::Abs(q));
780 for (Int_t iEl = 0; iEl < nEl; iEl++) {
786 // Electron attachment
787 if (fPar->ElAttachOn()) {
788 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
792 // Apply the diffusion smearing
793 if (fPar->DiffusionOn()) {
794 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
797 // Apply E x B effects (depends on drift direction)
799 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
802 // The electron position after diffusion and ExB in pad coordinates
803 // The pad row (z-direction)
804 Float_t rowDist = xyz[2] - row0;
805 Int_t rowE = ((Int_t) (rowDist * divideRow));
806 if ((rowE < 0) || (rowE >= nRowMax)) continue;
807 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
809 // The pad column (rphi-direction)
810 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
811 Float_t colDist = xyz[1] - col0tilt;
812 Int_t colE = ((Int_t) (colDist * divideCol));
813 if ((colE < 0) || (colE >= nColMax)) continue;
814 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
816 // The time bin (negative for hits in the amplification region)
817 // In the amplification region the electrons drift from both sides
818 // to the middle (anode wire plane)
819 Float_t timeDist = time0 - xyz[0];
820 Float_t timeOffset = 0;
824 timeE = ((Int_t) (timeDist * divideTime));
825 // The distance of the position to the middle of the timebin
826 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
829 // Difference between half of the amplification gap width and
830 // the distance to the anode wire
831 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
833 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
834 // The distance of the position to the middle of the timebin
835 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
838 // Apply the gas gain including fluctuations
839 Float_t ggRndm = 0.0;
841 ggRndm = gRandom->Rndm();
842 } while (ggRndm <= 0);
843 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
845 // Apply the pad response
847 // The distance of the electron to the center of the pad
848 // in units of pad width
849 Float_t dist = - colOffset * divideCol;
850 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
854 padSignal[1] = signal;
858 // Sample the time response inside the drift region
859 // + additional time bins before and after.
860 // The sampling is done always in the middle of the time bin
861 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
862 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
865 // Apply the time response
866 Float_t timeResponse = 1.0;
867 Float_t crossTalk = 0.0;
868 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
870 timeResponse = fPar->TimeResponse(time);
873 crossTalk = fPar->CrossTalk(time);
880 for (iPad = 0; iPad < kNpad; iPad++) {
882 Int_t colPos = colE + iPad - 1;
883 if (colPos < 0) continue;
884 if (colPos >= nColMax) break;
887 // Note: The time bin number is shifted by nTimeBefore to avoid negative
888 // time bins. This has to be subtracted later.
889 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
890 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
891 if( colPos != colE ) {
892 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
895 signalOld[iPad] += padSignal[iPad] * timeResponse;
897 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
899 // Store the track index in the dictionary
900 // Note: We store index+1 in order to allow the array to be compressed
901 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
902 for (iDict = 0; iDict < kNDict; iDict++) {
903 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
906 if (oldTrack == track+1) break;
908 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
918 } // Loop: electrons of a single hit
920 } // If: detector and test hit
922 hit = (AliTRDhit *) fTRD->NextHit();
924 } // Loop: hits of one primary track
926 } // Loop: primary tracks
929 printf("<AliTRDdigitizer::MakeDigits> ");
930 printf("Finished analyzing %d hits\n",countHits);
933 // The coupling factor
934 Float_t coupling = fPar->GetPadCoupling()
935 * fPar->GetTimeCoupling();
937 // The conversion factor
938 Float_t convert = kEl2fC
939 * fPar->GetChipGain();
941 // Loop through all chambers to finalize the digits
943 Int_t iDetEnd = AliTRDgeometry::Ndet();
945 iDetBeg = fSimpleDet;
946 iDetEnd = iDetBeg + 1;
948 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
950 Int_t plane = fGeo->GetPlane(iDet);
951 Int_t sector = fGeo->GetSector(iDet);
952 Int_t chamber = fGeo->GetChamber(iDet);
953 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
954 Int_t nColMax = fPar->GetColMax(plane);
955 Int_t nTimeMax = fPar->GetTimeMax();
956 Int_t nTimeTotal = fPar->GetTimeTotal();
958 Double_t *inADC = new Double_t[nTimeTotal];
959 Double_t *outADC = new Double_t[nTimeTotal];
962 printf("<AliTRDdigitizer::MakeDigits> ");
963 printf("Digitization for chamber %d\n",iDet);
966 // Add a container for the digits of this detector
967 digits = fDigitsManager->GetDigits(iDet);
968 // Allocate memory space for the digits buffer
969 if (digits->GetNtime() == 0) {
970 digits->Allocate(nRowMax,nColMax,nTimeTotal);
972 else if (fSimpleSim) {
976 // Get the signal container
977 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
978 if (signals->GetNtime() == 0) {
979 // Create missing containers
980 signals->Allocate(nRowMax,nColMax,nTimeTotal);
983 // Expand the container if neccessary
984 if (fCompress) signals->Expand();
986 // Create the missing dictionary containers
988 for (iDict = 0; iDict < kNDict; iDict++) {
989 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
990 if (dictionary[iDict]->GetNtime() == 0) {
991 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
998 // Don't create noise in detectors that are switched off
999 if (CheckDetector(plane,chamber,sector)) {
1001 // Create the digits for this chamber
1002 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1003 for (iCol = 0; iCol < nColMax; iCol++ ) {
1005 // Create summable digits
1008 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1009 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1010 signalAmp *= fSDigitsScale;
1011 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1012 Int_t adc = (Int_t) signalAmp;
1013 if (adc > 0) nDigits++;
1014 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1018 // Create normal digits
1021 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1022 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1023 // Pad and time coupling
1024 signalAmp *= coupling;
1025 // Add the noise, starting from minus ADC baseline in electrons
1026 Double_t baselineEl = fPar->GetADCbaseline() * (fPar->GetADCinRange()
1027 / fPar->GetADCoutRange())
1029 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise())
1032 signalAmp *= convert;
1033 // Add ADC baseline in mV
1034 signalAmp += fPar->GetADCbaseline() * (fPar->GetADCinRange()
1035 / fPar->GetADCoutRange());
1036 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1037 // signal is larger than fADCinRange
1039 if (signalAmp >= fPar->GetADCinRange()) {
1040 adc = ((Int_t) fPar->GetADCoutRange());
1043 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1044 / fPar->GetADCinRange())));
1047 outADC[iTime] = adc;
1050 // Apply the tail cancelation via the digital filter
1052 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1055 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1056 // Store the amplitude of the digit if above threshold
1057 if (outADC[iTime] > fPar->GetADCthreshold()) {
1059 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1060 ,iRow,iCol,iTime,outADC[iTime]);
1063 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1074 // Compress the arrays
1076 digits->Compress(1,0);
1077 for (iDict = 0; iDict < kNDict; iDict++) {
1078 dictionary[iDict]->Compress(1,0);
1081 totalSizeDigits += digits->GetSize();
1082 totalSizeDict0 += dictionary[0]->GetSize();
1083 totalSizeDict1 += dictionary[1]->GetSize();
1084 totalSizeDict2 += dictionary[2]->GetSize();
1086 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1088 printf("<AliTRDdigitizer::MakeDigits> ");
1089 printf("Found %d digits in detector %d (%3.0f).\n"
1091 ,100.0 * ((Float_t) nDigits) / nPixel);
1094 if (fCompress) signals->Compress(1,0);
1104 delete signalsArray;
1109 printf("<AliTRDdigitizer::MakeDigits> ");
1110 printf("Total number of analyzed hits = %d\n",countHits);
1112 printf("<AliTRDdigitizer::MakeDigits> ");
1113 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1124 //_____________________________________________________________________________
1125 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1128 // Add a digits manager for s-digits to the input list.
1131 fSDigitsManagerList->Add(man);
1135 //_____________________________________________________________________________
1136 void AliTRDdigitizer::DeleteSDigitsManager()
1139 // Removes digits manager from the input list.
1142 fSDigitsManagerList->Delete();
1146 //_____________________________________________________________________________
1147 Bool_t AliTRDdigitizer::ConvertSDigits()
1150 // Converts s-digits to normal digits
1153 // Number of track dictionary arrays
1154 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1156 // Converts number of electrons to fC
1157 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1165 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1167 printf("<AliTRDdigitizer::ConvertSDigits> ");
1168 printf("Create the default parameter object\n");
1172 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1173 Double_t noise = fPar->GetNoise();
1174 Double_t padCoupling = fPar->GetPadCoupling();
1175 Double_t timeCoupling = fPar->GetTimeCoupling();
1176 Double_t chipGain = fPar->GetChipGain();
1177 Double_t coupling = padCoupling * timeCoupling;
1178 Double_t convert = kEl2fC * chipGain;
1179 Double_t adcInRange = fPar->GetADCinRange();
1180 Double_t adcOutRange = fPar->GetADCoutRange();
1181 Int_t adcThreshold = fPar->GetADCthreshold();
1182 Int_t adcBaseline = fPar->GetADCbaseline();
1184 AliTRDdataArrayI *digitsIn;
1185 AliTRDdataArrayI *digitsOut;
1186 AliTRDdataArrayI *dictionaryIn[kNDict];
1187 AliTRDdataArrayI *dictionaryOut[kNDict];
1189 // Loop through the detectors
1190 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1193 printf("<AliTRDdigitizer::ConvertSDigits> ");
1194 printf("Convert detector %d to digits.\n",iDet);
1197 Int_t plane = fGeo->GetPlane(iDet);
1198 Int_t sector = fGeo->GetSector(iDet);
1199 Int_t chamber = fGeo->GetChamber(iDet);
1200 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1201 Int_t nColMax = fPar->GetColMax(plane);
1202 Int_t nTimeTotal = fPar->GetTimeTotal();
1204 Double_t *inADC = new Double_t[nTimeTotal];
1205 Double_t *outADC = new Double_t[nTimeTotal];
1207 digitsIn = fSDigitsManager->GetDigits(iDet);
1209 digitsOut = fDigitsManager->GetDigits(iDet);
1210 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1211 for (iDict = 0; iDict < kNDict; iDict++) {
1212 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1213 dictionaryIn[iDict]->Expand();
1214 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1215 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1218 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1219 for (iCol = 0; iCol < nColMax; iCol++ ) {
1221 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1222 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1223 signal *= sDigitsScale;
1224 // Pad and time coupling
1226 // Add the noise, starting from minus ADC baseline in electrons
1227 Double_t baselineEl = adcBaseline * (adcInRange / adcOutRange) / convert;
1228 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),-baselineEl);
1231 // add ADC baseline in mV
1232 signal += adcBaseline * (adcInRange / adcOutRange);
1233 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1234 // signal is larger than adcInRange
1236 if (signal >= adcInRange) {
1237 adc = ((Int_t) adcOutRange);
1240 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1243 outADC[iTime] = adc;
1246 // Apply the tail cancelation via the digital filter
1248 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1251 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1252 // Store the amplitude of the digit if above threshold
1253 if (outADC[iTime] > adcThreshold) {
1254 digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1255 // Copy the dictionary
1256 for (iDict = 0; iDict < kNDict; iDict++) {
1257 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1258 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1267 digitsIn->Compress(1,0);
1268 digitsOut->Compress(1,0);
1269 for (iDict = 0; iDict < kNDict; iDict++) {
1270 dictionaryIn[iDict]->Compress(1,0);
1271 dictionaryOut[iDict]->Compress(1,0);
1284 //_____________________________________________________________________________
1285 Bool_t AliTRDdigitizer::MergeSDigits()
1288 // Merges the input s-digits:
1289 // - The amplitude of the different inputs are summed up.
1290 // - Of the track IDs from the input dictionaries only one is
1291 // kept for each input. This works for maximal 3 different merged inputs.
1294 // Number of track dictionary arrays
1295 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1298 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1300 printf("<AliTRDdigitizer::MergeSDigits> ");
1301 printf("Create the default parameter object\n");
1308 AliTRDdataArrayI *digitsA;
1309 AliTRDdataArrayI *digitsB;
1310 AliTRDdataArrayI *dictionaryA[kNDict];
1311 AliTRDdataArrayI *dictionaryB[kNDict];
1313 // Get the first s-digits
1314 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1315 if (!fSDigitsManager) return kFALSE;
1317 // Loop through the other sets of s-digits
1318 AliTRDdigitsManager *mergeSDigitsManager;
1319 mergeSDigitsManager = (AliTRDdigitsManager *)
1320 fSDigitsManagerList->After(fSDigitsManager);
1323 if (mergeSDigitsManager) {
1324 printf("<AliTRDdigitizer::MergeSDigits> ");
1325 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1328 printf("<AliTRDdigitizer::MergeSDigits> ");
1329 printf("Only one input file.\n");
1334 while (mergeSDigitsManager) {
1338 // Loop through the detectors
1339 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1341 Int_t plane = fGeo->GetPlane(iDet);
1342 Int_t sector = fGeo->GetSector(iDet);
1343 Int_t chamber = fGeo->GetChamber(iDet);
1344 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1345 Int_t nColMax = fPar->GetColMax(plane);
1346 Int_t nTimeTotal = fPar->GetTimeTotal();
1348 // Loop through the pixels of one detector and add the signals
1349 digitsA = fSDigitsManager->GetDigits(iDet);
1350 digitsB = mergeSDigitsManager->GetDigits(iDet);
1353 for (iDict = 0; iDict < kNDict; iDict++) {
1354 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1355 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1356 dictionaryA[iDict]->Expand();
1357 dictionaryB[iDict]->Expand();
1360 // Merge only detectors that contain a signal
1361 Bool_t doMerge = kTRUE;
1362 if (fMergeSignalOnly) {
1363 if (digitsA->GetOverThreshold(0) == 0) {
1371 printf("<AliTRDdigitizer::MergeSDigits> ");
1372 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1375 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1376 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1377 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1379 // Add the amplitudes of the summable digits
1380 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1381 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1383 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1385 // Add the mask to the track id if defined.
1386 for (iDict = 0; iDict < kNDict; iDict++) {
1387 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1388 if ((fMasks) && (trackB > 0)) {
1389 for (jDict = 0; jDict < kNDict; jDict++) {
1390 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1392 trackA = trackB + fMasks[iMerge];
1393 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1406 digitsA->Compress(1,0);
1407 digitsB->Compress(1,0);
1408 for (iDict = 0; iDict < kNDict; iDict++) {
1409 dictionaryA[iDict]->Compress(1,0);
1410 dictionaryB[iDict]->Compress(1,0);
1416 // The next set of s-digits
1417 mergeSDigitsManager = (AliTRDdigitsManager *)
1418 fSDigitsManagerList->After(mergeSDigitsManager);
1426 //_____________________________________________________________________________
1427 Bool_t AliTRDdigitizer::SDigits2Digits()
1430 // Merges the input s-digits and converts them to normal digits
1433 if (!MergeSDigits()) return kFALSE;
1435 return ConvertSDigits();
1439 //_____________________________________________________________________________
1440 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1443 // Checks whether a detector is enabled
1446 if (fSimpleSim) return kTRUE;
1448 if ((fTRD->GetSensChamber() >= 0) &&
1449 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1450 if ((fTRD->GetSensPlane() >= 0) &&
1451 (fTRD->GetSensPlane() != plane)) return kFALSE;
1452 if ( fTRD->GetSensSector() >= 0) {
1453 Int_t sens1 = fTRD->GetSensSector();
1454 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1455 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1456 * AliTRDgeometry::Nsect();
1457 if (sens1 < sens2) {
1458 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1461 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1469 //_____________________________________________________________________________
1470 Bool_t AliTRDdigitizer::WriteDigits() const
1473 // Writes out the TRD-digits and the dictionaries
1476 // Store the digits and the dictionary in the tree
1477 return fDigitsManager->WriteDigits();
1481 //_____________________________________________________________________________
1482 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1483 , Int_t n, Int_t nexp)
1486 // Does the deconvolution by the digital filter.
1488 // Author: Marcus Gutfleisch, KIP Heidelberg
1489 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1490 // Pad-ground capacitance = 25 pF
1491 // Pad-pad cross talk capacitance = 6 pF
1492 // For 10 MHz digitization, corresponding to 20 time bins
1493 // in the drift region
1497 Double_t coefficients[2];
1499 /* initialize (coefficient = alpha, rates = lambda) */
1502 rates[0] = 0.466998;
1504 coefficients[0] = 1.0;
1507 rates[0] = 0.8988162;
1508 coefficients[0] = 0.11392069;
1509 rates[1] = 0.3745688;
1510 coefficients[1] = 0.8860793;
1512 Float_t sumc = coefficients[0]+coefficients[1];
1513 coefficients[0] /= sumc;
1514 coefficients[1] /= sumc;
1518 Double_t reminder[2];
1519 Double_t correction, result;
1521 /* attention: computation order is important */
1523 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1525 for ( i = 0; i < n; i++ ) {
1526 result = ( source[i] - correction ); /* no rescaling */
1529 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1530 * ( reminder[k] + coefficients[k] * result);
1533 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1538 //_____________________________________________________________________________
1539 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1542 // Initializes the output branches
1549 Error("InitOutput","Run Loader is NULL");
1552 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
1555 Error("Open","Can not get TRD loader from Run Loader");
1563 //if we produce SDigits
1564 tree = loader->TreeS();
1567 loader->MakeTree("S");
1568 tree = loader->TreeS();
1572 {//if we produce Digits
1573 tree = loader->TreeD();
1576 loader->MakeTree("D");
1577 tree = loader->TreeD();
1580 fDigitsManager->SetEvent(iEvent);
1581 fDigitsManager->MakeBranch(tree);