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
87 fSDigitsManagerList = 0;
97 fMergeSignalOnly = kFALSE;
103 //_____________________________________________________________________________
104 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
105 :AliDigitizer(name,title)
108 // AliTRDdigitizer constructor
113 //NewIO: These data members probably are not needed anymore
116 fSDigitsManagerList = 0;
126 fSDigitsScale = 100.; // For the summable digits
127 fMergeSignalOnly = kFALSE;
134 //_____________________________________________________________________________
135 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
136 , const Text_t *name, const Text_t *title)
137 :AliDigitizer(manager,name,title)
140 // AliTRDdigitizer constructor
146 fSDigitsManagerList = 0;
155 fSDigitsScale = 100.; // For the summable digits
156 fMergeSignalOnly = kFALSE;
163 //_____________________________________________________________________________
164 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
165 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
168 // AliTRDdigitizer constructor
175 fSDigitsManagerList = 0;
184 fSDigitsScale = 100.; // For the summable digits
185 fMergeSignalOnly = kFALSE;
192 //_____________________________________________________________________________
193 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d):AliDigitizer(d)
196 // AliTRDdigitizer copy constructor
199 ((AliTRDdigitizer &) d).Copy(*this);
203 //_____________________________________________________________________________
204 AliTRDdigitizer::~AliTRDdigitizer()
207 // AliTRDdigitizer destructor
211 if (fDigitsManager) {
212 delete fDigitsManager;
216 if (fSDigitsManager) {
217 delete fSDigitsManager;
221 if (fSDigitsManagerList) {
222 delete fSDigitsManagerList;
223 fSDigitsManagerList = 0;
233 //_____________________________________________________________________________
234 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
237 // Assignment operator
240 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
245 //_____________________________________________________________________________
246 void AliTRDdigitizer::Copy(TObject &d)
252 ((AliTRDdigitizer &) d).fRunLoader = 0;
253 ((AliTRDdigitizer &) d).fDigitsManager = 0;
254 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
255 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
256 ((AliTRDdigitizer &) d).fTRD = 0;
257 ((AliTRDdigitizer &) d).fGeo = 0;
258 ((AliTRDdigitizer &) d).fPar = 0;
259 ((AliTRDdigitizer &) d).fEvent = 0;
260 ((AliTRDdigitizer &) d).fMasks = 0;
261 ((AliTRDdigitizer &) d).fCompress = fCompress;
262 ((AliTRDdigitizer &) d).fDebug = fDebug ;
263 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
264 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
265 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
266 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
267 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
271 //_____________________________________________________________________________
272 void AliTRDdigitizer::Exec(Option_t* option)
275 // Executes the merging
280 AliTRDdigitsManager *sdigitsManager;
282 TString optionString = option;
283 if (optionString.Contains("deb")) {
285 if (optionString.Contains("2")) {
288 printf("<AliTRDdigitizer::Exec> ");
289 printf("Called with debug option %d\n",fDebug);
292 // The AliRoot file is already connected by the manager
298 printf("<AliTRDdigitizer::Exec> ");
299 printf("AliRun object found on file.\n");
303 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
305 gAlice = inrl->GetAliRun();
308 printf("<AliTRDdigitizer::Exec> ");
309 printf("Could not find AliRun object.\n");
314 Int_t nInput = fManager->GetNinputs();
315 fMasks = new Int_t[nInput];
316 for (iInput = 0; iInput < nInput; iInput++) {
317 fMasks[iInput] = fManager->GetMask(iInput);
322 AliRunLoader* orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
323 if (InitDetector()) {
324 AliLoader* ogime = orl->GetLoader("TRDLoader");
329 //if we produce SDigits
330 tree = ogime->TreeS();
333 ogime->MakeTree("S");
334 tree = ogime->TreeS();
338 {//if we produce Digits
339 tree = ogime->TreeD();
342 ogime->MakeTree("D");
343 tree = ogime->TreeD();
349 for (iInput = 0; iInput < nInput; iInput++) {
352 printf("<AliTRDdigitizer::Exec> ");
353 printf("Add input stream %d\n",iInput);
356 // check if the input tree exists
357 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
358 AliLoader* gime = inrl->GetLoader("TRDLoader");
360 TTree * treees = gime->TreeS();
363 if (gime->LoadSDigits())
365 Error("Exec","Error Occured while loading S. Digits for input %d.",iInput);
368 treees = gime->TreeS();
372 printf("<AliTRDdigitizer::Exec> ");
373 printf("Input stream %d does not exist\n",iInput);
377 // Read the s-digits via digits manager
378 sdigitsManager = new AliTRDdigitsManager();
379 sdigitsManager->SetDebug(fDebug);
380 sdigitsManager->SetSDigits(kTRUE);
382 AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
383 AliLoader* gimme = rl->GetLoader("TRDLoader");
384 if (!gimme->TreeS()) gimme->LoadSDigits();
385 sdigitsManager->ReadDigits(gimme->TreeS());
387 // Add the s-digits to the input list
388 AddSDigitsManager(sdigitsManager);
392 // Convert the s-digits to normal digits
394 printf("<AliTRDdigitizer::Exec> ");
395 printf("Do the conversion\n");
401 printf("<AliTRDdigitizer::Exec> ");
402 printf("Write the digits\n");
405 fDigitsManager->WriteDigits();
409 GetParameter()->Write();
412 printf("<AliTRDdigitizer::Exec> ");
416 DeleteSDigitsManager();
420 //_____________________________________________________________________________
421 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
424 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
427 // Connect the AliRoot file containing Geometry, Kine, and Hits
429 fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
434 Error("Open","Can not open session for file %s.",file);
438 fRunLoader->LoadgAlice();
439 gAlice = fRunLoader->GetAliRun();
443 printf("<AliTRDdigitizer::Open> ");
444 printf("AliRun object found on file.\n");
448 printf("<AliTRDdigitizer::Open> ");
449 printf("Could not find AliRun object.\n");
455 // Import the Trees for the event nEvent in the file
456 fRunLoader->GetEvent(fEvent);
458 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
461 Error("Open","Can not get TRD loader from Run Loader");
465 if (InitDetector()) {
469 //if we produce SDigits
470 tree = loader->TreeS();
473 loader->MakeTree("S");
474 tree = loader->TreeS();
478 {//if we produce Digits
479 tree = loader->TreeD();
482 loader->MakeTree("D");
483 tree = loader->TreeD();
486 return MakeBranch(tree);
494 //_____________________________________________________________________________
495 Bool_t AliTRDdigitizer::InitDetector()
498 // Sets the pointer to the TRD detector and the geometry
501 // Get the pointer to the detector class and check for version 1
502 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
504 printf("<AliTRDdigitizer::InitDetector> ");
505 printf("No TRD module found\n");
508 if (fTRD->IsVersion() != 1) {
509 printf("<AliTRDdigitizer::InitDetector> ");
510 printf("TRD must be version 1 (slow simulator).\n");
515 fGeo = fTRD->GetGeometry();
517 printf("<AliTRDdigitizer::InitDetector> ");
518 printf("Geometry version %d\n",fGeo->IsVersion());
521 // Create a digits manager
522 fDigitsManager = new AliTRDdigitsManager();
523 fDigitsManager->SetSDigits(fSDigits);
524 fDigitsManager->CreateArrays();
525 fDigitsManager->SetEvent(fEvent);
526 fDigitsManager->SetDebug(fDebug);
528 // The list for the input s-digits manager to be merged
529 fSDigitsManagerList = new TList();
535 //_____________________________________________________________________________
536 Bool_t AliTRDdigitizer::MakeBranch(TTree* tree) const
539 // Create the branches for the digits array
542 return fDigitsManager->MakeBranch(tree);
546 //_____________________________________________________________________________
547 Bool_t AliTRDdigitizer::MakeDigits()
553 ///////////////////////////////////////////////////////////////
555 ///////////////////////////////////////////////////////////////
557 // Converts number of electrons to fC
558 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
560 ///////////////////////////////////////////////////////////////
562 // Number of pads included in the pad response
563 const Int_t kNpad = 3;
565 // Number of track dictionary arrays
566 const Int_t kNDict = AliTRDdigitsManager::kNDict;
568 // Half the width of the amplification region
569 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
571 Int_t iRow, iCol, iTime, iPad;
575 Int_t totalSizeDigits = 0;
576 Int_t totalSizeDict0 = 0;
577 Int_t totalSizeDict1 = 0;
578 Int_t totalSizeDict2 = 0;
580 Int_t timeTRDbeg = 0;
581 Int_t timeTRDend = 1;
586 Float_t padSignal[kNpad];
587 Float_t signalOld[kNpad];
589 AliTRDdataArrayF *signals = 0;
590 AliTRDdataArrayI *digits = 0;
591 AliTRDdataArrayI *dictionary[kNDict];
593 // Create a default parameter class if none is defined
595 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
597 printf("<AliTRDdigitizer::MakeDigits> ");
598 printf("Create the default parameter object\n");
602 // Create a container for the amplitudes
603 AliTRDsegmentArray *signalsArray
604 = new AliTRDsegmentArray("AliTRDdataArrayF"
605 ,AliTRDgeometry::Ndet());
608 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
609 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
611 printf("<AliTRDdigitizer::MakeDigits> ");
612 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
616 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
619 printf("<AliTRDdigitizer::MakeDigits> ");
620 printf("No geometry defined\n");
625 printf("<AliTRDdigitizer::MakeDigits> ");
626 printf("Start creating digits.\n");
629 AliLoader* gimme = fRunLoader->GetLoader("TRDLoader");
630 if (!gimme->TreeH()) gimme->LoadHits();
631 TTree* hitTree = gimme->TreeH();
634 Error("MakeDigits","Can not get TreeH");
637 fTRD->SetTreeAddress();
639 // Get the number of entries in the hit tree
640 // (Number of primary particles creating a hit somewhere)
643 nTrack = (Int_t) hitTree->GetEntries();
645 printf("<AliTRDdigitizer::MakeDigits> ");
646 printf("Found %d primary particles\n",nTrack);
650 Int_t detectorOld = -1;
653 // Loop through all entries in the tree
654 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
658 nBytes += hitTree->GetEvent(iTrack);
661 // Loop through the TRD hits
663 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
672 Float_t q = hit->GetCharge();
673 Int_t track = hit->Track();
674 Int_t detector = hit->GetDetector();
675 Int_t plane = fGeo->GetPlane(detector);
676 Int_t sector = fGeo->GetSector(detector);
677 Int_t chamber = fGeo->GetChamber(detector);
678 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
679 Int_t nColMax = fPar->GetColMax(plane);
680 Int_t nTimeMax = fPar->GetTimeMax();
681 Int_t nTimeBefore = fPar->GetTimeBefore();
682 Int_t nTimeAfter = fPar->GetTimeAfter();
683 Int_t nTimeTotal = fPar->GetTimeTotal();
684 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
685 Float_t col0 = fPar->GetCol0(plane);
686 Float_t time0 = fPar->GetTime0(plane);
687 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
688 Float_t colPadSize = fPar->GetColPadSize(plane);
689 Float_t timeBinSize = fPar->GetTimeBinSize();
690 Float_t divideRow = 1.0 / rowPadSize;
691 Float_t divideCol = 1.0 / colPadSize;
692 Float_t divideTime = 1.0 / timeBinSize;
695 printf("Analyze hit no. %d ",iHit);
696 printf("-----------------------------------------------------------\n");
698 printf("plane = %d, sector = %d, chamber = %d\n"
699 ,plane,sector,chamber);
700 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
701 ,nRowMax,nColMax,nTimeMax);
702 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
703 ,nTimeBefore,nTimeAfter,nTimeTotal);
704 printf("row0 = %f, col0 = %f, time0 = %f\n"
706 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
707 ,rowPadSize,colPadSize,timeBinSize);
710 // Don't analyze test hits and switched off detectors
711 if ((CheckDetector(plane,chamber,sector)) &&
712 (((Int_t) q) != 0)) {
714 if (detector != detectorOld) {
717 printf("<AliTRDdigitizer::MakeDigits> ");
718 printf("Get new container. New det = %d, Old det = %d\n"
719 ,detector,detectorOld);
721 // Compress the old one if enabled
722 if ((fCompress) && (detectorOld > -1)) {
724 printf("<AliTRDdigitizer::MakeDigits> ");
725 printf("Compress the old container ...");
727 signals->Compress(1,0);
728 for (iDict = 0; iDict < kNDict; iDict++) {
729 dictionary[iDict]->Compress(1,0);
731 if (fDebug > 1) printf("done\n");
733 // Get the new container
734 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
735 if (signals->GetNtime() == 0) {
736 // Allocate a new one if not yet existing
738 printf("<AliTRDdigitizer::MakeDigits> ");
739 printf("Allocate a new container ... ");
741 signals->Allocate(nRowMax,nColMax,nTimeTotal);
743 else if (fSimpleSim) {
744 // Clear an old one for the simple simulation
746 printf("<AliTRDdigitizer::MakeDigits> ");
747 printf("Clear a old container ... ");
752 // Expand an existing one
755 printf("<AliTRDdigitizer::MakeDigits> ");
756 printf("Expand an existing container ... ");
761 // The same for the dictionary
763 for (iDict = 0; iDict < kNDict; iDict++) {
764 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
765 if (dictionary[iDict]->GetNtime() == 0) {
766 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
769 if (fCompress) dictionary[iDict]->Expand();
773 if (fDebug > 1) printf("done\n");
774 detectorOld = detector;
777 // Rotate the sectors on top of each other
784 fGeo->Rotate(detector,pos,rot);
787 // The driftlength. It is negative if the hit is in the
788 // amplification region.
789 Float_t driftlength = time0 - rot[0];
791 // Take also the drift in the amplification region into account
792 // The drift length is at the moment still the same, regardless of
793 // the position relativ to the wire. This non-isochronity needs still
794 // to be implemented.
795 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
796 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
798 // Loop over all electrons of this hit
799 // TR photons produce hits with negative charge
800 Int_t nEl = ((Int_t) TMath::Abs(q));
801 for (Int_t iEl = 0; iEl < nEl; iEl++) {
807 // Electron attachment
808 if (fPar->ElAttachOn()) {
809 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
813 // Apply the diffusion smearing
814 if (fPar->DiffusionOn()) {
815 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
818 // Apply E x B effects (depends on drift direction)
820 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
823 // The electron position after diffusion and ExB in pad coordinates
824 // The pad row (z-direction)
825 Float_t rowDist = xyz[2] - row0;
826 Int_t rowE = ((Int_t) (rowDist * divideRow));
827 if ((rowE < 0) || (rowE >= nRowMax)) continue;
828 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
830 // The pad column (rphi-direction)
831 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
832 Float_t colDist = xyz[1] - col0tilt;
833 Int_t colE = ((Int_t) (colDist * divideCol));
834 if ((colE < 0) || (colE >= nColMax)) continue;
835 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
837 // The time bin (negative for hits in the amplification region)
838 // In the amplification region the electrons drift from both sides
839 // to the middle (anode wire plane)
840 Float_t timeDist = time0 - xyz[0];
841 Float_t timeOffset = 0;
845 timeE = ((Int_t) (timeDist * divideTime));
846 // The distance of the position to the middle of the timebin
847 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
850 // Difference between half of the amplification gap width and
851 // the distance to the anode wire
852 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
854 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
855 // The distance of the position to the middle of the timebin
856 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
859 // Apply the gas gain including fluctuations
860 Float_t ggRndm = 0.0;
862 ggRndm = gRandom->Rndm();
863 } while (ggRndm <= 0);
864 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
866 // Apply the pad response
868 // The distance of the electron to the center of the pad
869 // in units of pad width
870 Float_t dist = - colOffset * divideCol;
871 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
875 padSignal[1] = signal;
879 // Sample the time response inside the drift region
880 // + additional time bins before and after.
881 // The sampling is done always in the middle of the time bin
882 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
883 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
886 // Apply the time response
887 Float_t timeResponse = 1.0;
888 Float_t crossTalk = 0.0;
889 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
891 timeResponse = fPar->TimeResponse(time);
894 crossTalk = fPar->CrossTalk(time);
901 for (iPad = 0; iPad < kNpad; iPad++) {
903 Int_t colPos = colE + iPad - 1;
904 if (colPos < 0) continue;
905 if (colPos >= nColMax) break;
908 // Note: The time bin number is shifted by nTimeBefore to avoid negative
909 // time bins. This has to be subtracted later.
910 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
911 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
912 if( colPos != colE ) {
913 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
916 signalOld[iPad] += padSignal[iPad] * timeResponse;
918 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
920 // Store the track index in the dictionary
921 // Note: We store index+1 in order to allow the array to be compressed
922 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
923 for (iDict = 0; iDict < kNDict; iDict++) {
924 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
927 if (oldTrack == track+1) break;
929 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
939 } // Loop: electrons of a single hit
941 } // If: detector and test hit
943 hit = (AliTRDhit *) fTRD->NextHit();
945 } // Loop: hits of one primary track
947 } // Loop: primary tracks
950 printf("<AliTRDdigitizer::MakeDigits> ");
951 printf("Finished analyzing %d hits\n",countHits);
954 // The coupling factor
955 Float_t coupling = fPar->GetPadCoupling()
956 * fPar->GetTimeCoupling();
958 // The conversion factor
959 Float_t convert = kEl2fC
960 * fPar->GetChipGain();
962 // Loop through all chambers to finalize the digits
964 Int_t iDetEnd = AliTRDgeometry::Ndet();
966 iDetBeg = fSimpleDet;
967 iDetEnd = iDetBeg + 1;
969 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
971 Int_t plane = fGeo->GetPlane(iDet);
972 Int_t sector = fGeo->GetSector(iDet);
973 Int_t chamber = fGeo->GetChamber(iDet);
974 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
975 Int_t nColMax = fPar->GetColMax(plane);
976 Int_t nTimeMax = fPar->GetTimeMax();
977 Int_t nTimeTotal = fPar->GetTimeTotal();
979 Double_t *inADC = new Double_t[nTimeTotal];
980 Double_t *outADC = new Double_t[nTimeTotal];
983 printf("<AliTRDdigitizer::MakeDigits> ");
984 printf("Digitization for chamber %d\n",iDet);
987 // Add a container for the digits of this detector
988 digits = fDigitsManager->GetDigits(iDet);
989 // Allocate memory space for the digits buffer
990 if (digits->GetNtime() == 0) {
991 digits->Allocate(nRowMax,nColMax,nTimeTotal);
993 else if (fSimpleSim) {
997 // Get the signal container
998 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
999 if (signals->GetNtime() == 0) {
1000 // Create missing containers
1001 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1004 // Expand the container if neccessary
1005 if (fCompress) signals->Expand();
1007 // Create the missing dictionary containers
1009 for (iDict = 0; iDict < kNDict; iDict++) {
1010 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1011 if (dictionary[iDict]->GetNtime() == 0) {
1012 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1019 // Don't create noise in detectors that are switched off
1020 if (CheckDetector(plane,chamber,sector)) {
1022 // Create the digits for this chamber
1023 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1024 for (iCol = 0; iCol < nColMax; iCol++ ) {
1026 // Create summable digits
1029 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1030 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1031 signalAmp *= fSDigitsScale;
1032 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1033 Int_t adc = (Int_t) signalAmp;
1034 if (adc > 0) nDigits++;
1035 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1039 // Create normal digits
1042 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1043 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1044 // Pad and time coupling
1045 signalAmp *= coupling;
1046 // Add the noise, starting from minus ADC baseline in electrons
1047 Double_t baselineEl = fPar->GetADCbaseline() * (fPar->GetADCinRange()
1048 / fPar->GetADCoutRange())
1050 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise())
1053 signalAmp *= convert;
1054 // Add ADC baseline in mV
1055 signalAmp += fPar->GetADCbaseline() * (fPar->GetADCinRange()
1056 / fPar->GetADCoutRange());
1057 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1058 // signal is larger than fADCinRange
1060 if (signalAmp >= fPar->GetADCinRange()) {
1061 adc = ((Int_t) fPar->GetADCoutRange());
1064 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1065 / fPar->GetADCinRange())));
1068 outADC[iTime] = adc;
1071 // Apply the tail cancelation via the digital filter
1073 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1076 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1077 // Store the amplitude of the digit if above threshold
1078 if (outADC[iTime] > fPar->GetADCthreshold()) {
1080 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1081 ,iRow,iCol,iTime,outADC[iTime]);
1084 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1095 // Compress the arrays
1097 digits->Compress(1,0);
1098 for (iDict = 0; iDict < kNDict; iDict++) {
1099 dictionary[iDict]->Compress(1,0);
1102 totalSizeDigits += digits->GetSize();
1103 totalSizeDict0 += dictionary[0]->GetSize();
1104 totalSizeDict1 += dictionary[1]->GetSize();
1105 totalSizeDict2 += dictionary[2]->GetSize();
1107 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1109 printf("<AliTRDdigitizer::MakeDigits> ");
1110 printf("Found %d digits in detector %d (%3.0f).\n"
1112 ,100.0 * ((Float_t) nDigits) / nPixel);
1115 if (fCompress) signals->Compress(1,0);
1125 delete signalsArray;
1130 printf("<AliTRDdigitizer::MakeDigits> ");
1131 printf("Total number of analyzed hits = %d\n",countHits);
1133 printf("<AliTRDdigitizer::MakeDigits> ");
1134 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1145 //_____________________________________________________________________________
1146 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1149 // Add a digits manager for s-digits to the input list.
1152 fSDigitsManagerList->Add(man);
1156 //_____________________________________________________________________________
1157 void AliTRDdigitizer::DeleteSDigitsManager()
1160 // Removes digits manager from the input list.
1163 fSDigitsManagerList->Delete();
1167 //_____________________________________________________________________________
1168 Bool_t AliTRDdigitizer::ConvertSDigits()
1171 // Converts s-digits to normal digits
1174 // Number of track dictionary arrays
1175 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1177 // Converts number of electrons to fC
1178 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1186 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1188 printf("<AliTRDdigitizer::ConvertSDigits> ");
1189 printf("Create the default parameter object\n");
1193 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1194 Double_t noise = fPar->GetNoise();
1195 Double_t padCoupling = fPar->GetPadCoupling();
1196 Double_t timeCoupling = fPar->GetTimeCoupling();
1197 Double_t chipGain = fPar->GetChipGain();
1198 Double_t coupling = padCoupling * timeCoupling;
1199 Double_t convert = kEl2fC * chipGain;
1200 Double_t adcInRange = fPar->GetADCinRange();
1201 Double_t adcOutRange = fPar->GetADCoutRange();
1202 Int_t adcThreshold = fPar->GetADCthreshold();
1203 Int_t adcBaseline = fPar->GetADCbaseline();
1205 AliTRDdataArrayI *digitsIn;
1206 AliTRDdataArrayI *digitsOut;
1207 AliTRDdataArrayI *dictionaryIn[kNDict];
1208 AliTRDdataArrayI *dictionaryOut[kNDict];
1210 // Loop through the detectors
1211 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1214 printf("<AliTRDdigitizer::ConvertSDigits> ");
1215 printf("Convert detector %d to digits.\n",iDet);
1218 Int_t plane = fGeo->GetPlane(iDet);
1219 Int_t sector = fGeo->GetSector(iDet);
1220 Int_t chamber = fGeo->GetChamber(iDet);
1221 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1222 Int_t nColMax = fPar->GetColMax(plane);
1223 Int_t nTimeTotal = fPar->GetTimeTotal();
1225 Double_t *inADC = new Double_t[nTimeTotal];
1226 Double_t *outADC = new Double_t[nTimeTotal];
1228 digitsIn = fSDigitsManager->GetDigits(iDet);
1230 digitsOut = fDigitsManager->GetDigits(iDet);
1231 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1232 for (iDict = 0; iDict < kNDict; iDict++) {
1233 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1234 dictionaryIn[iDict]->Expand();
1235 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1236 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1239 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1240 for (iCol = 0; iCol < nColMax; iCol++ ) {
1242 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1243 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1244 signal *= sDigitsScale;
1245 // Pad and time coupling
1247 // Add the noise, starting from minus ADC baseline in electrons
1248 Double_t baselineEl = adcBaseline * (adcInRange / adcOutRange) / convert;
1249 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),-baselineEl);
1252 // add ADC baseline in mV
1253 signal += adcBaseline * (adcInRange / adcOutRange);
1254 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1255 // signal is larger than adcInRange
1257 if (signal >= adcInRange) {
1258 adc = ((Int_t) adcOutRange);
1261 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1264 outADC[iTime] = adc;
1267 // Apply the tail cancelation via the digital filter
1269 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1272 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1273 // Store the amplitude of the digit if above threshold
1274 if (outADC[iTime] > adcThreshold) {
1275 digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1276 // Copy the dictionary
1277 for (iDict = 0; iDict < kNDict; iDict++) {
1278 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1279 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1288 digitsIn->Compress(1,0);
1289 digitsOut->Compress(1,0);
1290 for (iDict = 0; iDict < kNDict; iDict++) {
1291 dictionaryIn[iDict]->Compress(1,0);
1292 dictionaryOut[iDict]->Compress(1,0);
1305 //_____________________________________________________________________________
1306 Bool_t AliTRDdigitizer::MergeSDigits()
1309 // Merges the input s-digits:
1310 // - The amplitude of the different inputs are summed up.
1311 // - Of the track IDs from the input dictionaries only one is
1312 // kept for each input. This works for maximal 3 different merged inputs.
1315 // Number of track dictionary arrays
1316 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1319 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1321 printf("<AliTRDdigitizer::MergeSDigits> ");
1322 printf("Create the default parameter object\n");
1329 AliTRDdataArrayI *digitsA;
1330 AliTRDdataArrayI *digitsB;
1331 AliTRDdataArrayI *dictionaryA[kNDict];
1332 AliTRDdataArrayI *dictionaryB[kNDict];
1334 // Get the first s-digits
1335 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1336 if (!fSDigitsManager) return kFALSE;
1338 // Loop through the other sets of s-digits
1339 AliTRDdigitsManager *mergeSDigitsManager;
1340 mergeSDigitsManager = (AliTRDdigitsManager *)
1341 fSDigitsManagerList->After(fSDigitsManager);
1344 if (mergeSDigitsManager) {
1345 printf("<AliTRDdigitizer::MergeSDigits> ");
1346 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1349 printf("<AliTRDdigitizer::MergeSDigits> ");
1350 printf("Only one input file.\n");
1355 while (mergeSDigitsManager) {
1359 // Loop through the detectors
1360 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1362 Int_t plane = fGeo->GetPlane(iDet);
1363 Int_t sector = fGeo->GetSector(iDet);
1364 Int_t chamber = fGeo->GetChamber(iDet);
1365 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1366 Int_t nColMax = fPar->GetColMax(plane);
1367 Int_t nTimeTotal = fPar->GetTimeTotal();
1369 // Loop through the pixels of one detector and add the signals
1370 digitsA = fSDigitsManager->GetDigits(iDet);
1371 digitsB = mergeSDigitsManager->GetDigits(iDet);
1374 for (iDict = 0; iDict < kNDict; iDict++) {
1375 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1376 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1377 dictionaryA[iDict]->Expand();
1378 dictionaryB[iDict]->Expand();
1381 // Merge only detectors that contain a signal
1382 Bool_t doMerge = kTRUE;
1383 if (fMergeSignalOnly) {
1384 if (digitsA->GetOverThreshold(0) == 0) {
1392 printf("<AliTRDdigitizer::MergeSDigits> ");
1393 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1396 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1397 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1398 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1400 // Add the amplitudes of the summable digits
1401 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1402 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1404 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1406 // Add the mask to the track id if defined.
1407 for (iDict = 0; iDict < kNDict; iDict++) {
1408 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1409 if ((fMasks) && (trackB > 0)) {
1410 for (jDict = 0; jDict < kNDict; jDict++) {
1411 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1413 trackA = trackB + fMasks[iMerge];
1414 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1427 digitsA->Compress(1,0);
1428 digitsB->Compress(1,0);
1429 for (iDict = 0; iDict < kNDict; iDict++) {
1430 dictionaryA[iDict]->Compress(1,0);
1431 dictionaryB[iDict]->Compress(1,0);
1437 // The next set of s-digits
1438 mergeSDigitsManager = (AliTRDdigitsManager *)
1439 fSDigitsManagerList->After(mergeSDigitsManager);
1447 //_____________________________________________________________________________
1448 Bool_t AliTRDdigitizer::SDigits2Digits()
1451 // Merges the input s-digits and converts them to normal digits
1454 if (!MergeSDigits()) return kFALSE;
1456 return ConvertSDigits();
1460 //_____________________________________________________________________________
1461 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1464 // Checks whether a detector is enabled
1467 if (fSimpleSim) return kTRUE;
1469 if ((fTRD->GetSensChamber() >= 0) &&
1470 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1471 if ((fTRD->GetSensPlane() >= 0) &&
1472 (fTRD->GetSensPlane() != plane)) return kFALSE;
1473 if ( fTRD->GetSensSector() >= 0) {
1474 Int_t sens1 = fTRD->GetSensSector();
1475 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1476 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1477 * AliTRDgeometry::Nsect();
1478 if (sens1 < sens2) {
1479 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1482 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1490 //_____________________________________________________________________________
1491 Bool_t AliTRDdigitizer::WriteDigits() const
1494 // Writes out the TRD-digits and the dictionaries
1497 // Store the digits and the dictionary in the tree
1498 return fDigitsManager->WriteDigits();
1502 //_____________________________________________________________________________
1503 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1504 , Int_t n, Int_t nexp)
1507 // Does the deconvolution by the digital filter.
1509 // Author: Marcus Gutfleisch, KIP Heidelberg
1510 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1511 // Pad-ground capacitance = 25 pF
1512 // Pad-pad cross talk capacitance = 6 pF
1513 // For 10 MHz digitization, corresponding to 20 time bins
1514 // in the drift region
1518 Double_t coefficients[2];
1520 /* initialize (coefficient = alpha, rates = lambda) */
1523 rates[0] = 0.466998;
1525 coefficients[0] = 1.0;
1528 rates[0] = 0.8988162;
1529 coefficients[0] = 0.11392069;
1530 rates[1] = 0.3745688;
1531 coefficients[1] = 0.8860793;
1533 Float_t sumc = coefficients[0]+coefficients[1];
1534 coefficients[0] /= sumc;
1535 coefficients[1] /= sumc;
1539 Double_t reminder[2];
1540 Double_t correction, result;
1542 /* attention: computation order is important */
1544 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1546 for ( i = 0; i < n; i++ ) {
1547 result = ( source[i] - correction ); /* no rescaling */
1550 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1551 * ( reminder[k] + coefficients[k] * result);
1554 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1559 //_____________________________________________________________________________
1560 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1563 // Initializes the output branches
1570 Error("InitOutput","Run Loader is NULL");
1573 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
1576 Error("Open","Can not get TRD loader from Run Loader");
1584 //if we produce SDigits
1585 tree = loader->TreeS();
1588 loader->MakeTree("S");
1589 tree = loader->TreeS();
1593 {//if we produce Digits
1594 tree = loader->TreeD();
1597 loader->MakeTree("D");
1598 tree = loader->TreeD();
1601 fDigitsManager->SetEvent(iEvent);
1602 fDigitsManager->MakeBranch(tree);