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;
218 if (fSDigitsManagerList) {
219 fSDigitsManagerList->Delete();
220 delete fSDigitsManagerList;
221 fSDigitsManagerList = 0;
231 //_____________________________________________________________________________
232 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
235 // Assignment operator
238 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
243 //_____________________________________________________________________________
244 void AliTRDdigitizer::Copy(TObject &d)
250 ((AliTRDdigitizer &) d).fRunLoader = 0;
251 ((AliTRDdigitizer &) d).fDigitsManager = 0;
252 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
253 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
254 ((AliTRDdigitizer &) d).fTRD = 0;
255 ((AliTRDdigitizer &) d).fGeo = 0;
256 ((AliTRDdigitizer &) d).fPar = 0;
257 ((AliTRDdigitizer &) d).fEvent = 0;
258 ((AliTRDdigitizer &) d).fMasks = 0;
259 ((AliTRDdigitizer &) d).fCompress = fCompress;
260 ((AliTRDdigitizer &) d).fDebug = fDebug ;
261 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
262 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
263 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
264 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
265 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
269 //_____________________________________________________________________________
270 void AliTRDdigitizer::Exec(Option_t* option)
273 // Executes the merging
278 AliTRDdigitsManager *sdigitsManager;
280 TString optionString = option;
281 if (optionString.Contains("deb")) {
283 if (optionString.Contains("2")) {
286 printf("<AliTRDdigitizer::Exec> ");
287 printf("Called with debug option %d\n",fDebug);
290 // The AliRoot file is already connected by the manager
296 printf("<AliTRDdigitizer::Exec> ");
297 printf("AliRun object found on file.\n");
301 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
303 gAlice = inrl->GetAliRun();
306 printf("<AliTRDdigitizer::Exec> ");
307 printf("Could not find AliRun object.\n");
312 Int_t nInput = fManager->GetNinputs();
313 fMasks = new Int_t[nInput];
314 for (iInput = 0; iInput < nInput; iInput++) {
315 fMasks[iInput] = fManager->GetMask(iInput);
320 AliRunLoader* orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
321 if (InitDetector()) {
322 AliLoader* ogime = orl->GetLoader("TRDLoader");
327 //if we produce SDigits
328 tree = ogime->TreeS();
331 ogime->MakeTree("S");
332 tree = ogime->TreeS();
336 {//if we produce Digits
337 tree = ogime->TreeD();
340 ogime->MakeTree("D");
341 tree = ogime->TreeD();
347 for (iInput = 0; iInput < nInput; iInput++) {
350 printf("<AliTRDdigitizer::Exec> ");
351 printf("Add input stream %d\n",iInput);
354 // check if the input tree exists
355 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
356 AliLoader* gime = inrl->GetLoader("TRDLoader");
358 TTree * treees = gime->TreeS();
361 if (gime->LoadSDigits())
363 Error("Exec","Error Occured while loading S. Digits for input %d.",iInput);
366 treees = gime->TreeS();
370 printf("<AliTRDdigitizer::Exec> ");
371 printf("Input stream %d does not exist\n",iInput);
375 // Read the s-digits via digits manager
376 sdigitsManager = new AliTRDdigitsManager();
377 sdigitsManager->SetDebug(fDebug);
378 sdigitsManager->SetSDigits(kTRUE);
380 AliRunLoader* rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
381 AliLoader* gimme = rl->GetLoader("TRDLoader");
382 if (!gimme->TreeS()) gimme->LoadSDigits();
383 sdigitsManager->ReadDigits(gimme->TreeS());
385 // Add the s-digits to the input list
386 AddSDigitsManager(sdigitsManager);
390 // Convert the s-digits to normal digits
392 printf("<AliTRDdigitizer::Exec> ");
393 printf("Do the conversion\n");
399 printf("<AliTRDdigitizer::Exec> ");
400 printf("Write the digits\n");
403 fDigitsManager->WriteDigits();
407 if (!gFile->Get("TRDparameter")) GetParameter()->Write();
410 printf("<AliTRDdigitizer::Exec> ");
414 DeleteSDigitsManager();
418 //_____________________________________________________________________________
419 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
422 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
425 // Connect the AliRoot file containing Geometry, Kine, and Hits
428 TString evfoldname = AliConfig::GetDefaultEventFolderName();
429 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
431 fRunLoader = AliRunLoader::Open(file,AliConfig::GetDefaultEventFolderName(),
436 Error("Open","Can not open session for file %s.",file);
440 if (!fRunLoader->GetAliRun()) fRunLoader->LoadgAlice();
441 gAlice = fRunLoader->GetAliRun();
445 printf("<AliTRDdigitizer::Open> ");
446 printf("AliRun object found on file.\n");
450 printf("<AliTRDdigitizer::Open> ");
451 printf("Could not find AliRun object.\n");
457 // Import the Trees for the event nEvent in the file
458 fRunLoader->GetEvent(fEvent);
460 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
463 Error("Open","Can not get TRD loader from Run Loader");
467 if (InitDetector()) {
471 //if we produce SDigits
472 tree = loader->TreeS();
475 loader->MakeTree("S");
476 tree = loader->TreeS();
480 {//if we produce Digits
483 loader->MakeTree("D");
484 tree = loader->TreeD();
487 return MakeBranch(tree);
495 //_____________________________________________________________________________
496 Bool_t AliTRDdigitizer::InitDetector()
499 // Sets the pointer to the TRD detector and the geometry
502 // Get the pointer to the detector class and check for version 1
503 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
505 printf("<AliTRDdigitizer::InitDetector> ");
506 printf("No TRD module found\n");
509 if (fTRD->IsVersion() != 1) {
510 printf("<AliTRDdigitizer::InitDetector> ");
511 printf("TRD must be version 1 (slow simulator).\n");
516 fGeo = fTRD->GetGeometry();
518 printf("<AliTRDdigitizer::InitDetector> ");
519 printf("Geometry version %d\n",fGeo->IsVersion());
522 // Create a digits manager
523 delete fDigitsManager;
524 fDigitsManager = new AliTRDdigitsManager();
525 fDigitsManager->SetSDigits(fSDigits);
526 fDigitsManager->CreateArrays();
527 fDigitsManager->SetEvent(fEvent);
528 fDigitsManager->SetDebug(fDebug);
530 // The list for the input s-digits manager to be merged
531 if (fSDigitsManagerList) {
532 fSDigitsManagerList->Delete();
534 fSDigitsManagerList = new TList();
541 //_____________________________________________________________________________
542 Bool_t AliTRDdigitizer::MakeBranch(TTree* tree) const
545 // Create the branches for the digits array
548 return fDigitsManager->MakeBranch(tree);
552 //_____________________________________________________________________________
553 Bool_t AliTRDdigitizer::MakeDigits()
559 ///////////////////////////////////////////////////////////////
561 ///////////////////////////////////////////////////////////////
563 // Converts number of electrons to fC
564 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
566 ///////////////////////////////////////////////////////////////
568 // Number of pads included in the pad response
569 const Int_t kNpad = 3;
571 // Number of track dictionary arrays
572 const Int_t kNDict = AliTRDdigitsManager::kNDict;
574 // Half the width of the amplification region
575 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
577 Int_t iRow, iCol, iTime, iPad;
581 Int_t totalSizeDigits = 0;
582 Int_t totalSizeDict0 = 0;
583 Int_t totalSizeDict1 = 0;
584 Int_t totalSizeDict2 = 0;
586 Int_t timeTRDbeg = 0;
587 Int_t timeTRDend = 1;
592 Float_t padSignal[kNpad];
593 Float_t signalOld[kNpad];
595 AliTRDdataArrayF *signals = 0;
596 AliTRDdataArrayI *digits = 0;
597 AliTRDdataArrayI *dictionary[kNDict];
599 // Create a default parameter class if none is defined
601 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
603 printf("<AliTRDdigitizer::MakeDigits> ");
604 printf("Create the default parameter object\n");
608 // Create a container for the amplitudes
609 AliTRDsegmentArray *signalsArray
610 = new AliTRDsegmentArray("AliTRDdataArrayF"
611 ,AliTRDgeometry::Ndet());
614 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
615 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
617 printf("<AliTRDdigitizer::MakeDigits> ");
618 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
622 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
625 printf("<AliTRDdigitizer::MakeDigits> ");
626 printf("No geometry defined\n");
631 printf("<AliTRDdigitizer::MakeDigits> ");
632 printf("Start creating digits.\n");
635 AliLoader* gimme = fRunLoader->GetLoader("TRDLoader");
636 if (!gimme->TreeH()) gimme->LoadHits();
637 TTree* hitTree = gimme->TreeH();
640 Error("MakeDigits","Can not get TreeH");
643 fTRD->SetTreeAddress();
645 // Get the number of entries in the hit tree
646 // (Number of primary particles creating a hit somewhere)
649 nTrack = (Int_t) hitTree->GetEntries();
651 printf("<AliTRDdigitizer::MakeDigits> ");
652 printf("Found %d primary particles\n",nTrack);
656 Int_t detectorOld = -1;
660 printf("<AliTRDdigitizer::MakeDigits> ");
661 printf("Driftvelocity = %.2f, Sampling = %.2fns\n",fPar->GetDriftVelocity(),
662 fPar->GetTimeBinSize() / fPar->GetDriftVelocity() * 1000.);
663 printf("<AliTRDdigitizer::MakeDigits> ");
664 printf("Gain = %d, Noise = %d\n",(Int_t)fPar->GetGasGain(),(Int_t)fPar->GetNoise());
665 if (fPar->TimeStructOn()) {
666 printf("<AliTRDdigitizer::MakeDigits> ");
667 printf("Time Structure of drift cells implemented.\n");
668 if (fPar->StaggeringOn()) {
669 printf("<AliTRDdigitizer::MakeDigits> ");
670 printf("Staggered geometry.\n");
672 printf("<AliTRDdigitizer::MakeDigits> ");
673 printf("Non-Staggered geometry.\n");
676 printf("<AliTRDdigitizer::MakeDigits> ");
677 printf("Constant drift velocity in drift cells.\n");
681 // Loop through all entries in the tree
682 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
686 nBytes += hitTree->GetEvent(iTrack);
689 // Loop through the TRD hits
691 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
700 Float_t q = hit->GetCharge();
701 Int_t track = hit->Track();
702 Int_t detector = hit->GetDetector();
703 Int_t plane = fGeo->GetPlane(detector);
704 Int_t sector = fGeo->GetSector(detector);
705 Int_t chamber = fGeo->GetChamber(detector);
706 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
707 Int_t nColMax = fPar->GetColMax(plane);
708 Int_t nTimeMax = fPar->GetTimeMax();
709 Int_t nTimeBefore = fPar->GetTimeBefore();
710 Int_t nTimeAfter = fPar->GetTimeAfter();
711 Int_t nTimeTotal = fPar->GetTimeTotal();
712 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
713 Float_t col0 = fPar->GetCol0(plane);
714 Float_t time0 = fPar->GetTime0(plane);
715 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
716 Float_t colPadSize = fPar->GetColPadSize(plane);
717 Float_t timeBinSize = fPar->GetTimeBinSize();
718 Float_t divideRow = 1.0 / rowPadSize;
719 Float_t divideCol = 1.0 / colPadSize;
720 Float_t divideTime = 1.0 / timeBinSize;
723 printf("Analyze hit no. %d ",iHit);
724 printf("-----------------------------------------------------------\n");
726 printf("plane = %d, sector = %d, chamber = %d\n"
727 ,plane,sector,chamber);
728 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
729 ,nRowMax,nColMax,nTimeMax);
730 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
731 ,nTimeBefore,nTimeAfter,nTimeTotal);
732 printf("row0 = %f, col0 = %f, time0 = %f\n"
734 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
735 ,rowPadSize,colPadSize,timeBinSize);
738 // Don't analyze test hits and switched off detectors
739 if ((CheckDetector(plane,chamber,sector)) &&
740 (((Int_t) q) != 0)) {
742 if (detector != detectorOld) {
745 printf("<AliTRDdigitizer::MakeDigits> ");
746 printf("Get new container. New det = %d, Old det = %d\n"
747 ,detector,detectorOld);
749 // Compress the old one if enabled
750 if ((fCompress) && (detectorOld > -1)) {
752 printf("<AliTRDdigitizer::MakeDigits> ");
753 printf("Compress the old container ...");
755 signals->Compress(1,0);
756 for (iDict = 0; iDict < kNDict; iDict++) {
757 dictionary[iDict]->Compress(1,0);
759 if (fDebug > 1) printf("done\n");
761 // Get the new container
762 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
763 if (signals->GetNtime() == 0) {
764 // Allocate a new one if not yet existing
766 printf("<AliTRDdigitizer::MakeDigits> ");
767 printf("Allocate a new container ... ");
769 signals->Allocate(nRowMax,nColMax,nTimeTotal);
771 else if (fSimpleSim) {
772 // Clear an old one for the simple simulation
774 printf("<AliTRDdigitizer::MakeDigits> ");
775 printf("Clear a old container ... ");
780 // Expand an existing one
783 printf("<AliTRDdigitizer::MakeDigits> ");
784 printf("Expand an existing container ... ");
789 // The same for the dictionary
791 for (iDict = 0; iDict < kNDict; iDict++) {
792 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
793 if (dictionary[iDict]->GetNtime() == 0) {
794 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
797 if (fCompress) dictionary[iDict]->Expand();
801 if (fDebug > 1) printf("done\n");
802 detectorOld = detector;
805 // Rotate the sectors on top of each other
812 fGeo->Rotate(detector,pos,rot);
815 // The driftlength. It is negative if the hit is in the
816 // amplification region.
817 Float_t driftlength = time0 - rot[0];
819 // Take also the drift in the amplification region into account
820 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
821 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
823 // Loop over all electrons of this hit
824 // TR photons produce hits with negative charge
825 Int_t nEl = ((Int_t) TMath::Abs(q));
826 for (Int_t iEl = 0; iEl < nEl; iEl++) {
832 // Stupid patch to take care of TR photons that are absorbed
833 // outside the chamber volume. A real fix would actually need
834 // a more clever implementation of the TR hit generation
835 Float_t zz = xyz[2] - row0;
836 if ((zz < 0.0) || (zz > rowPadSize*nRowMax)) {
838 printf("<AliTRDdigitizer::MakeDigits> ");
839 printf("Hit outside of sensitive volume, row (Q = %d)\n",((Int_t) q));
843 Int_t tt = ((Int_t) (10*(driftlength + 2.0*kAmWidth)));
846 printf("<AliTRDdigitizer::MakeDigits> ");
847 printf("Hit outside of sensitive volume, time (Q = %d)\n",((Int_t) q));
852 // Electron attachment
853 if (fPar->ElAttachOn()) {
854 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
858 // Apply the diffusion smearing
859 if (fPar->DiffusionOn()) {
860 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
863 // Apply E x B effects (depends on drift direction)
865 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
868 // Apply the drift time correction
869 if (fPar->TimeStructOn()) {
870 // Get z-position with respect to anode wire:
871 Float_t Z = xyz[2] - row0 + fPar->GetAnodeWireOffset();
872 while (Z>0.5) Z -= 0.5;
873 if (Z>0.25) Z = 0.5-Z;
874 if (!(fPar->TimeStruct(driftlength+2*kAmWidth,Z,xyz))) continue;
877 // The electron position after diffusion and ExB in pad coordinates
878 // The pad row (z-direction)
879 Float_t rowDist = xyz[2] - row0;
880 Int_t rowE = ((Int_t) (rowDist * divideRow));
881 if ((rowE < 0) || (rowE >= nRowMax)) continue;
882 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
884 // The pad column (rphi-direction)
885 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
886 Float_t colDist = xyz[1] - col0tilt;
887 Int_t colE = ((Int_t) (colDist * divideCol));
888 if ((colE < 0) || (colE >= nColMax)) continue;
889 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
891 // The time bin (negative for hits in the amplification region)
892 // In the amplification region the electrons drift from both sides
893 // to the middle (anode wire plane)
894 Float_t timeDist = time0 - xyz[0];
895 Float_t timeOffset = 0;
899 timeE = ((Int_t) (timeDist * divideTime));
900 // The distance of the position to the middle of the timebin
901 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
904 // Difference between half of the amplification gap width and
905 // the distance to the anode wire
906 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
908 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
909 // The distance of the position to the middle of the timebin
910 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
913 // Apply the gas gain including fluctuations
914 Float_t ggRndm = 0.0;
916 ggRndm = gRandom->Rndm();
917 } while (ggRndm <= 0);
918 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
920 // Apply the pad response
922 // The distance of the electron to the center of the pad
923 // in units of pad width
924 Float_t dist = - colOffset * divideCol;
925 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
929 padSignal[1] = signal;
933 // Sample the time response inside the drift region
934 // + additional time bins before and after.
935 // The sampling is done always in the middle of the time bin
936 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
937 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
940 // Apply the time response
941 Float_t timeResponse = 1.0;
942 Float_t crossTalk = 0.0;
943 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
945 timeResponse = fPar->TimeResponse(time);
948 crossTalk = fPar->CrossTalk(time);
955 for (iPad = 0; iPad < kNpad; iPad++) {
957 Int_t colPos = colE + iPad - 1;
958 if (colPos < 0) continue;
959 if (colPos >= nColMax) break;
962 // Note: The time bin number is shifted by nTimeBefore to avoid negative
963 // time bins. This has to be subtracted later.
964 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
965 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
966 if( colPos != colE ) {
967 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
970 signalOld[iPad] += padSignal[iPad] * timeResponse;
972 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
974 // Store the track index in the dictionary
975 // Note: We store index+1 in order to allow the array to be compressed
976 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
977 for (iDict = 0; iDict < kNDict; iDict++) {
978 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
981 if (oldTrack == track+1) break;
983 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
993 } // Loop: electrons of a single hit
995 } // If: detector and test hit
997 hit = (AliTRDhit *) fTRD->NextHit();
999 } // Loop: hits of one primary track
1001 } // Loop: primary tracks
1004 printf("<AliTRDdigitizer::MakeDigits> ");
1005 printf("Finished analyzing %d hits\n",countHits);
1008 // The coupling factor
1009 Float_t coupling = fPar->GetPadCoupling()
1010 * fPar->GetTimeCoupling();
1012 // The conversion factor
1013 Float_t convert = kEl2fC
1014 * fPar->GetChipGain();
1016 // Loop through all chambers to finalize the digits
1018 Int_t iDetEnd = AliTRDgeometry::Ndet();
1020 iDetBeg = fSimpleDet;
1021 iDetEnd = iDetBeg + 1;
1023 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
1025 Int_t plane = fGeo->GetPlane(iDet);
1026 Int_t sector = fGeo->GetSector(iDet);
1027 Int_t chamber = fGeo->GetChamber(iDet);
1028 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1029 Int_t nColMax = fPar->GetColMax(plane);
1030 Int_t nTimeMax = fPar->GetTimeMax();
1031 Int_t nTimeTotal = fPar->GetTimeTotal();
1033 Double_t *inADC = new Double_t[nTimeTotal];
1034 Double_t *outADC = new Double_t[nTimeTotal];
1037 printf("<AliTRDdigitizer::MakeDigits> ");
1038 printf("Digitization for chamber %d\n",iDet);
1041 // Add a container for the digits of this detector
1042 digits = fDigitsManager->GetDigits(iDet);
1043 // Allocate memory space for the digits buffer
1044 if (digits->GetNtime() == 0) {
1045 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1047 else if (fSimpleSim) {
1051 // Get the signal container
1052 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1053 if (signals->GetNtime() == 0) {
1054 // Create missing containers
1055 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1058 // Expand the container if neccessary
1059 if (fCompress) signals->Expand();
1061 // Create the missing dictionary containers
1063 for (iDict = 0; iDict < kNDict; iDict++) {
1064 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1065 if (dictionary[iDict]->GetNtime() == 0) {
1066 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1073 // Don't create noise in detectors that are switched off
1074 if (CheckDetector(plane,chamber,sector)) {
1076 // Create the digits for this chamber
1077 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1078 for (iCol = 0; iCol < nColMax; iCol++ ) {
1080 // Create summable digits
1083 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1084 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1085 signalAmp *= fSDigitsScale;
1086 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1087 Int_t adc = (Int_t) signalAmp;
1088 if (adc > 0) nDigits++;
1089 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1093 // Create normal digits
1096 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1097 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1098 // Pad and time coupling
1099 signalAmp *= coupling;
1100 // Add the noise, starting from minus ADC baseline in electrons
1101 Double_t baselineEl = fPar->GetADCbaseline() * (fPar->GetADCinRange()
1102 / fPar->GetADCoutRange())
1104 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise())
1107 signalAmp *= convert;
1108 // Add ADC baseline in mV
1109 signalAmp += fPar->GetADCbaseline() * (fPar->GetADCinRange()
1110 / fPar->GetADCoutRange());
1111 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1112 // signal is larger than fADCinRange
1114 if (signalAmp >= fPar->GetADCinRange()) {
1115 adc = ((Int_t) fPar->GetADCoutRange());
1118 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1119 / fPar->GetADCinRange())));
1122 outADC[iTime] = adc;
1125 // Apply the tail cancelation via the digital filter
1127 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1130 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1131 // Store the amplitude of the digit if above threshold
1132 if (outADC[iTime] > fPar->GetADCthreshold()) {
1134 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1135 ,iRow,iCol,iTime,outADC[iTime]);
1138 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1149 // Compress the arrays
1151 digits->Compress(1,0);
1152 for (iDict = 0; iDict < kNDict; iDict++) {
1153 dictionary[iDict]->Compress(1,0);
1156 totalSizeDigits += digits->GetSize();
1157 totalSizeDict0 += dictionary[0]->GetSize();
1158 totalSizeDict1 += dictionary[1]->GetSize();
1159 totalSizeDict2 += dictionary[2]->GetSize();
1161 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1163 printf("<AliTRDdigitizer::MakeDigits> ");
1164 printf("Found %d digits in detector %d (%3.0f).\n"
1166 ,100.0 * ((Float_t) nDigits) / nPixel);
1169 if (fCompress) signals->Compress(1,0);
1179 delete signalsArray;
1184 printf("<AliTRDdigitizer::MakeDigits> ");
1185 printf("Total number of analyzed hits = %d\n",countHits);
1187 printf("<AliTRDdigitizer::MakeDigits> ");
1188 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1199 //_____________________________________________________________________________
1200 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1203 // Add a digits manager for s-digits to the input list.
1206 fSDigitsManagerList->Add(man);
1210 //_____________________________________________________________________________
1211 void AliTRDdigitizer::DeleteSDigitsManager()
1214 // Removes digits manager from the input list.
1217 fSDigitsManagerList->Delete();
1221 //_____________________________________________________________________________
1222 Bool_t AliTRDdigitizer::ConvertSDigits()
1225 // Converts s-digits to normal digits
1228 // Number of track dictionary arrays
1229 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1231 // Converts number of electrons to fC
1232 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1240 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1242 printf("<AliTRDdigitizer::ConvertSDigits> ");
1243 printf("Create the default parameter object\n");
1247 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1248 Double_t noise = fPar->GetNoise();
1249 Double_t padCoupling = fPar->GetPadCoupling();
1250 Double_t timeCoupling = fPar->GetTimeCoupling();
1251 Double_t chipGain = fPar->GetChipGain();
1252 Double_t coupling = padCoupling * timeCoupling;
1253 Double_t convert = kEl2fC * chipGain;
1254 Double_t adcInRange = fPar->GetADCinRange();
1255 Double_t adcOutRange = fPar->GetADCoutRange();
1256 Int_t adcThreshold = fPar->GetADCthreshold();
1257 Int_t adcBaseline = fPar->GetADCbaseline();
1259 AliTRDdataArrayI *digitsIn;
1260 AliTRDdataArrayI *digitsOut;
1261 AliTRDdataArrayI *dictionaryIn[kNDict];
1262 AliTRDdataArrayI *dictionaryOut[kNDict];
1264 // Loop through the detectors
1265 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1268 printf("<AliTRDdigitizer::ConvertSDigits> ");
1269 printf("Convert detector %d to digits.\n",iDet);
1272 Int_t plane = fGeo->GetPlane(iDet);
1273 Int_t sector = fGeo->GetSector(iDet);
1274 Int_t chamber = fGeo->GetChamber(iDet);
1275 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1276 Int_t nColMax = fPar->GetColMax(plane);
1277 Int_t nTimeTotal = fPar->GetTimeTotal();
1279 Double_t *inADC = new Double_t[nTimeTotal];
1280 Double_t *outADC = new Double_t[nTimeTotal];
1282 digitsIn = fSDigitsManager->GetDigits(iDet);
1284 digitsOut = fDigitsManager->GetDigits(iDet);
1285 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1286 for (iDict = 0; iDict < kNDict; iDict++) {
1287 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1288 dictionaryIn[iDict]->Expand();
1289 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1290 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1293 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1294 for (iCol = 0; iCol < nColMax; iCol++ ) {
1296 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1297 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1298 signal *= sDigitsScale;
1299 // Pad and time coupling
1301 // Add the noise, starting from minus ADC baseline in electrons
1302 Double_t baselineEl = adcBaseline * (adcInRange / adcOutRange) / convert;
1303 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),-baselineEl);
1306 // add ADC baseline in mV
1307 signal += adcBaseline * (adcInRange / adcOutRange);
1308 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1309 // signal is larger than adcInRange
1311 if (signal >= adcInRange) {
1312 adc = ((Int_t) adcOutRange);
1315 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1318 outADC[iTime] = adc;
1321 // Apply the tail cancelation via the digital filter
1323 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1326 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1327 // Store the amplitude of the digit if above threshold
1328 if (outADC[iTime] > adcThreshold) {
1329 digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1330 // Copy the dictionary
1331 for (iDict = 0; iDict < kNDict; iDict++) {
1332 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1333 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1342 digitsIn->Compress(1,0);
1343 digitsOut->Compress(1,0);
1344 for (iDict = 0; iDict < kNDict; iDict++) {
1345 dictionaryIn[iDict]->Compress(1,0);
1346 dictionaryOut[iDict]->Compress(1,0);
1359 //_____________________________________________________________________________
1360 Bool_t AliTRDdigitizer::MergeSDigits()
1363 // Merges the input s-digits:
1364 // - The amplitude of the different inputs are summed up.
1365 // - Of the track IDs from the input dictionaries only one is
1366 // kept for each input. This works for maximal 3 different merged inputs.
1369 // Number of track dictionary arrays
1370 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1373 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1375 printf("<AliTRDdigitizer::MergeSDigits> ");
1376 printf("Create the default parameter object\n");
1383 AliTRDdataArrayI *digitsA;
1384 AliTRDdataArrayI *digitsB;
1385 AliTRDdataArrayI *dictionaryA[kNDict];
1386 AliTRDdataArrayI *dictionaryB[kNDict];
1388 // Get the first s-digits
1389 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1390 if (!fSDigitsManager) return kFALSE;
1392 // Loop through the other sets of s-digits
1393 AliTRDdigitsManager *mergeSDigitsManager;
1394 mergeSDigitsManager = (AliTRDdigitsManager *)
1395 fSDigitsManagerList->After(fSDigitsManager);
1398 if (mergeSDigitsManager) {
1399 printf("<AliTRDdigitizer::MergeSDigits> ");
1400 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1403 printf("<AliTRDdigitizer::MergeSDigits> ");
1404 printf("Only one input file.\n");
1409 while (mergeSDigitsManager) {
1413 // Loop through the detectors
1414 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1416 Int_t plane = fGeo->GetPlane(iDet);
1417 Int_t sector = fGeo->GetSector(iDet);
1418 Int_t chamber = fGeo->GetChamber(iDet);
1419 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1420 Int_t nColMax = fPar->GetColMax(plane);
1421 Int_t nTimeTotal = fPar->GetTimeTotal();
1423 // Loop through the pixels of one detector and add the signals
1424 digitsA = fSDigitsManager->GetDigits(iDet);
1425 digitsB = mergeSDigitsManager->GetDigits(iDet);
1428 for (iDict = 0; iDict < kNDict; iDict++) {
1429 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1430 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1431 dictionaryA[iDict]->Expand();
1432 dictionaryB[iDict]->Expand();
1435 // Merge only detectors that contain a signal
1436 Bool_t doMerge = kTRUE;
1437 if (fMergeSignalOnly) {
1438 if (digitsA->GetOverThreshold(0) == 0) {
1446 printf("<AliTRDdigitizer::MergeSDigits> ");
1447 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1450 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1451 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1452 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1454 // Add the amplitudes of the summable digits
1455 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1456 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1458 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1460 // Add the mask to the track id if defined.
1461 for (iDict = 0; iDict < kNDict; iDict++) {
1462 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1463 if ((fMasks) && (trackB > 0)) {
1464 for (jDict = 0; jDict < kNDict; jDict++) {
1465 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1467 trackA = trackB + fMasks[iMerge];
1468 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1481 digitsA->Compress(1,0);
1482 digitsB->Compress(1,0);
1483 for (iDict = 0; iDict < kNDict; iDict++) {
1484 dictionaryA[iDict]->Compress(1,0);
1485 dictionaryB[iDict]->Compress(1,0);
1491 // The next set of s-digits
1492 mergeSDigitsManager = (AliTRDdigitsManager *)
1493 fSDigitsManagerList->After(mergeSDigitsManager);
1501 //_____________________________________________________________________________
1502 Bool_t AliTRDdigitizer::SDigits2Digits()
1505 // Merges the input s-digits and converts them to normal digits
1508 if (!MergeSDigits()) return kFALSE;
1510 return ConvertSDigits();
1514 //_____________________________________________________________________________
1515 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1518 // Checks whether a detector is enabled
1521 if (fSimpleSim) return kTRUE;
1523 if ((fTRD->GetSensChamber() >= 0) &&
1524 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1525 if ((fTRD->GetSensPlane() >= 0) &&
1526 (fTRD->GetSensPlane() != plane)) return kFALSE;
1527 if ( fTRD->GetSensSector() >= 0) {
1528 Int_t sens1 = fTRD->GetSensSector();
1529 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1530 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1531 * AliTRDgeometry::Nsect();
1532 if (sens1 < sens2) {
1533 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1536 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1544 //_____________________________________________________________________________
1545 Bool_t AliTRDdigitizer::WriteDigits() const
1548 // Writes out the TRD-digits and the dictionaries
1552 fRunLoader->CdGAFile();
1553 if (!gFile->Get("TRDparameter")) GetParameter()->Write();
1555 // Store the digits and the dictionary in the tree
1556 return fDigitsManager->WriteDigits();
1560 //_____________________________________________________________________________
1561 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1562 , Int_t n, Int_t nexp)
1565 // Does the deconvolution by the digital filter.
1567 // Author: Marcus Gutfleisch, KIP Heidelberg
1568 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1569 // Pad-ground capacitance = 25 pF
1570 // Pad-pad cross talk capacitance = 6 pF
1571 // For 10 MHz digitization, corresponding to 20 time bins
1572 // in the drift region
1576 Double_t coefficients[2];
1578 /* initialize (coefficient = alpha, rates = lambda) */
1581 rates[0] = 0.466998;
1583 coefficients[0] = 1.0;
1586 rates[0] = 0.8988162;
1587 coefficients[0] = 0.11392069;
1588 rates[1] = 0.3745688;
1589 coefficients[1] = 0.8860793;
1591 Float_t sumc = coefficients[0]+coefficients[1];
1592 coefficients[0] /= sumc;
1593 coefficients[1] /= sumc;
1597 Double_t reminder[2];
1598 Double_t correction, result;
1600 /* attention: computation order is important */
1602 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1604 for ( i = 0; i < n; i++ ) {
1605 result = ( source[i] - correction ); /* no rescaling */
1608 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1609 * ( reminder[k] + coefficients[k] * result);
1612 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1617 //_____________________________________________________________________________
1618 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1621 // Initializes the output branches
1628 Error("InitOutput","Run Loader is NULL");
1631 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
1634 Error("Open","Can not get TRD loader from Run Loader");
1642 //if we produce SDigits
1643 tree = loader->TreeS();
1646 loader->MakeTree("S");
1647 tree = loader->TreeS();
1651 {//if we produce Digits
1652 tree = loader->TreeD();
1655 loader->MakeTree("D");
1656 tree = loader->TreeD();
1659 fDigitsManager->SetEvent(iEvent);
1660 fDigitsManager->MakeBranch(tree);