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 //
22 // Authors: C. Blume (blume@ikf.uni-frankfurt.de) //
26 // The following effects are included: //
29 // - Gas gain including fluctuations //
30 // - Pad-response (simple Gaussian approximation) //
32 // - Electronics noise //
33 // - Electronics gain //
35 // - Zero suppression //
37 ////////////////////////////////////////////////////////////////////////////
39 #include <TGeoManager.h>
47 #include "AliRunLoader.h"
48 #include "AliLoader.h"
49 #include "AliConfig.h"
50 #include "AliDigitizationInput.h"
51 #include "AliRunLoader.h"
52 #include "AliLoader.h"
56 #include "AliTRDhit.h"
57 #include "AliTRDdigitizer.h"
58 #include "AliTRDarrayDictionary.h"
59 #include "AliTRDarrayADC.h"
60 #include "AliTRDarraySignal.h"
61 #include "AliTRDdigitsManager.h"
62 #include "AliTRDgeometry.h"
63 #include "AliTRDpadPlane.h"
64 #include "AliTRDcalibDB.h"
65 #include "AliTRDSimParam.h"
66 #include "AliTRDCommonParam.h"
67 #include "AliTRDfeeParam.h"
68 #include "AliTRDmcmSim.h"
69 #include "AliTRDdigitsParam.h"
71 #include "Cal/AliTRDCalROC.h"
72 #include "Cal/AliTRDCalDet.h"
73 #include "Cal/AliTRDCalOnlineGainTableROC.h"
75 ClassImp(AliTRDdigitizer)
77 //_____________________________________________________________________________
78 AliTRDdigitizer::AliTRDdigitizer()
83 ,fSDigitsManagerList(0)
86 ,fMcmSim(new AliTRDmcmSim)
91 ,fMergeSignalOnly(kFALSE)
94 // AliTRDdigitizer default constructor
99 //_____________________________________________________________________________
100 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
101 :AliDigitizer(name,title)
105 ,fSDigitsManagerList(0)
108 ,fMcmSim(new AliTRDmcmSim)
113 ,fMergeSignalOnly(kFALSE)
116 // AliTRDdigitizer constructor
121 //_____________________________________________________________________________
122 AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput
123 , const Text_t *name, const Text_t *title)
124 :AliDigitizer(digInput,name,title)
128 ,fSDigitsManagerList(0)
131 ,fMcmSim(new AliTRDmcmSim)
136 ,fMergeSignalOnly(kFALSE)
139 // AliTRDdigitizer constructor
144 //_____________________________________________________________________________
145 AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput)
146 :AliDigitizer(digInput,"AliTRDdigitizer","TRD digitizer")
150 ,fSDigitsManagerList(0)
153 ,fMcmSim(new AliTRDmcmSim)
158 ,fMergeSignalOnly(kFALSE)
161 // AliTRDdigitizer constructor
166 //_____________________________________________________________________________
167 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
172 ,fSDigitsManagerList(0)
175 ,fMcmSim(new AliTRDmcmSim)
178 ,fCompress(d.fCompress)
179 ,fSDigits(d.fSDigits)
180 ,fMergeSignalOnly(d.fMergeSignalOnly)
183 // AliTRDdigitizer copy constructor
188 //_____________________________________________________________________________
189 AliTRDdigitizer::~AliTRDdigitizer()
192 // AliTRDdigitizer destructor
195 delete fDigitsManager;
198 // s-digitsmanager will be deleted via list
200 fSDigitsManagerList->Delete();
201 delete fSDigitsManagerList;
202 fSDigitsManagerList = 0;
215 //_____________________________________________________________________________
216 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
219 // Assignment operator
223 ((AliTRDdigitizer &) d).Copy(*this);
230 //_____________________________________________________________________________
231 void AliTRDdigitizer::Copy(TObject &d) const
237 ((AliTRDdigitizer &) d).fRunLoader = 0;
238 ((AliTRDdigitizer &) d).fDigitsManager = 0;
239 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
240 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
241 ((AliTRDdigitizer &) d).fTRD = 0;
242 ((AliTRDdigitizer &) d).fGeo = 0;
243 ((AliTRDdigitizer &) d).fEvent = 0;
244 ((AliTRDdigitizer &) d).fMasks = 0;
245 ((AliTRDdigitizer &) d).fCompress = fCompress;
246 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
247 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
251 //_____________________________________________________________________________
252 void AliTRDdigitizer::Digitize(const Option_t* option)
255 // Executes the merging
260 AliTRDdigitsManager *sdigitsManager;
262 TString optionString = option;
263 if (optionString.Contains("deb")) {
264 AliLog::SetClassDebugLevel("AliTRDdigitizer",1);
265 AliInfo("Called with debug option");
268 // The AliRoot file is already connected by the manager
269 AliRunLoader *inrl = 0x0;
272 AliDebug(1,"AliRun object found on file.");
275 inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0));
277 gAlice = inrl->GetAliRun();
279 AliError("Could not find AliRun object.");
284 Int_t nInput = fDigInput->GetNinputs();
285 fMasks = new Int_t[nInput];
286 for (iInput = 0; iInput < nInput; iInput++) {
287 fMasks[iInput] = fDigInput->GetMask(iInput);
294 AliRunLoader *orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
296 if (InitDetector()) {
298 AliLoader *ogime = orl->GetLoader("TRDLoader");
302 // If we produce SDigits
303 tree = ogime->TreeS();
305 ogime->MakeTree("S");
306 tree = ogime->TreeS();
310 // If we produce Digits
311 tree = ogime->TreeD();
313 ogime->MakeTree("D");
314 tree = ogime->TreeD();
322 for (iInput = 0; iInput < nInput; iInput++) {
324 AliDebug(1,Form("Add input stream %d",iInput));
326 // Check if the input tree exists
327 inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
328 AliLoader *gime = inrl->GetLoader("TRDLoader");
330 TTree *treees = gime->TreeS();
332 if (gime->LoadSDigits()) {
333 AliError(Form("Error Occured while loading S. Digits for input %d.",iInput));
336 treees = gime->TreeS();
340 AliError(Form("Input stream %d does not exist",iInput));
344 // Read the s-digits via digits manager
345 sdigitsManager = new AliTRDdigitsManager();
346 sdigitsManager->SetSDigits(kTRUE);
348 AliRunLoader *rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
349 AliLoader *gimme = rl->GetLoader("TRDLoader");
352 gimme->LoadSDigits();
355 sdigitsManager->ReadDigits(gimme->TreeS());
357 // Add the s-digits to the input list
358 AddSDigitsManager(sdigitsManager);
362 // Convert the s-digits to normal digits
363 AliDebug(1,"Do the conversion");
367 AliDebug(1,"Write the digits");
368 fDigitsManager->WriteDigits();
374 DeleteSDigitsManager();
380 //_____________________________________________________________________________
381 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
384 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
386 // Connect the AliRoot file containing Geometry, Kine, and Hits
389 TString evfoldname = AliConfig::GetDefaultEventFolderName();
391 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
393 fRunLoader = AliRunLoader::Open(file,evfoldname,"UPDATE");
396 AliError(Form("Can not open session for file %s.",file));
400 if (!fRunLoader->GetAliRun()) {
401 fRunLoader->LoadgAlice();
403 gAlice = fRunLoader->GetAliRun();
406 AliDebug(1,"AliRun object found on file.");
409 AliError("Could not find AliRun object.");
415 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
417 AliError("Can not get TRD loader from Run Loader");
421 if (InitDetector()) {
424 // If we produce SDigits
425 tree = loader->TreeS();
427 loader->MakeTree("S");
428 tree = loader->TreeS();
432 // If we produce Digits
433 tree = loader->TreeD();
435 loader->MakeTree("D");
436 tree = loader->TreeD();
439 return MakeBranch(tree);
447 //_____________________________________________________________________________
448 Bool_t AliTRDdigitizer::Open(AliRunLoader * const runLoader, Int_t nEvent)
451 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
453 // Connect the AliRoot file containing Geometry, Kine, and Hits
456 fRunLoader = runLoader;
458 AliError("RunLoader does not exist");
462 if (!fRunLoader->GetAliRun()) {
463 fRunLoader->LoadgAlice();
465 gAlice = fRunLoader->GetAliRun();
468 AliDebug(1,"AliRun object found on file.");
471 AliError("Could not find AliRun object.");
477 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
479 AliError("Can not get TRD loader from Run Loader");
483 if (InitDetector()) {
486 // If we produce SDigits
487 tree = loader->TreeS();
489 loader->MakeTree("S");
490 tree = loader->TreeS();
494 // If we produce Digits
495 tree = loader->TreeD();
497 loader->MakeTree("D");
498 tree = loader->TreeD();
501 return MakeBranch(tree);
509 //_____________________________________________________________________________
510 Bool_t AliTRDdigitizer::InitDetector()
513 // Sets the pointer to the TRD detector and the geometry
516 // Get the pointer to the detector class and check for version 1
517 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
519 AliFatal("No TRD module found");
522 if (fTRD->IsVersion() != 1) {
523 AliFatal("TRD must be version 1 (slow simulator)");
528 fGeo = new AliTRDgeometry();
530 // Create a digits manager
531 if (fDigitsManager) {
532 delete fDigitsManager;
534 fDigitsManager = new AliTRDdigitsManager();
535 fDigitsManager->SetSDigits(fSDigits);
536 fDigitsManager->CreateArrays();
537 fDigitsManager->SetEvent(fEvent);
539 // The list for the input s-digits manager to be merged
540 if (fSDigitsManagerList) {
541 fSDigitsManagerList->Delete();
544 fSDigitsManagerList = new TList();
551 //_____________________________________________________________________________
552 Bool_t AliTRDdigitizer::MakeBranch(TTree *tree) const
555 // Create the branches for the digits array
558 return fDigitsManager->MakeBranch(tree);
562 //_____________________________________________________________________________
563 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
566 // Add a digits manager for s-digits to the input list.
569 fSDigitsManagerList->Add(man);
573 //_____________________________________________________________________________
574 void AliTRDdigitizer::DeleteSDigitsManager()
577 // Removes digits manager from the input list.
580 fSDigitsManagerList->Delete();
584 //_____________________________________________________________________________
585 Bool_t AliTRDdigitizer::MakeDigits()
591 AliDebug(1,"Start creating digits");
594 AliError("No geometry defined");
598 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
600 AliFatal("Could not get calibration object");
604 const Int_t kNdet = AliTRDgeometry::Ndet();
606 Float_t **hits = new Float_t*[kNdet];
607 Int_t *nhit = new Int_t[kNdet];
608 memset(nhit,0,kNdet*sizeof(Int_t));
610 AliTRDarraySignal *signals = 0x0;
612 // Check the number of time bins from simParam against OCDB,
613 // if OCDB value is not supposed to be used.
614 // As default, the value from OCDB is taken
615 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
616 if (calibration->GetNumberOfTimeBinsDCS() != AliTRDSimParam::Instance()->GetNTimeBins()) {
617 AliWarning(Form("Number of time bins is different to OCDB value [SIM=%d, OCDB=%d]"
618 ,AliTRDSimParam::Instance()->GetNTimeBins()
619 ,calibration->GetNumberOfTimeBinsDCS()));
621 // Save the values for the raw data headers
622 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
625 // Save the values for the raw data headers
626 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(calibration->GetNumberOfTimeBinsDCS());
629 // Save the values for the raw data headers
630 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
632 // Sort all hits according to detector number
633 if (!SortHits(hits,nhit)) {
634 AliError("Sorting hits failed");
640 // Loop through all detectors
641 for (Int_t det = 0; det < kNdet; det++) {
643 // Detectors that are switched off, not installed, etc.
644 if ((!calibration->IsChamberNoData(det)) &&
645 ( fGeo->ChamberInGeometry(det)) &&
648 signals = new AliTRDarraySignal();
650 // Convert the hits of the current detector to detector signals
651 if (!ConvertHits(det,hits[det],nhit[det],signals)) {
652 AliError(Form("Conversion of hits failed for detector=%d",det));
660 // Convert the detector signals to digits or s-digits
661 if (!ConvertSignals(det,signals)) {
662 AliError(Form("Conversion of signals failed for detector=%d",det));
670 // Delete the signals array
674 } // if: detector status
681 if (AliDataLoader *trklLoader
682 = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
683 if (trklLoader->Tree())
684 trklLoader->WriteData("OVERWRITE");
695 //_____________________________________________________________________________
696 Bool_t AliTRDdigitizer::SortHits(Float_t **hits, Int_t *nhit)
699 // Read all the hits and sorts them according to detector number
700 // in the output array <hits>.
703 AliDebug(1,"Start sorting hits");
705 const Int_t kNdet = AliTRDgeometry::Ndet();
706 // Size of the hit vector
707 const Int_t kNhit = 6;
712 Int_t *lhit = new Int_t[kNdet];
713 memset(lhit,0,kNdet*sizeof(Int_t));
715 for (Int_t det = 0; det < kNdet; det++) {
719 AliLoader *gimme = fRunLoader->GetLoader("TRDLoader");
720 if (!gimme->TreeH()) {
723 TTree *hitTree = gimme->TreeH();
724 if (hitTree == 0x0) {
725 AliError("Can not get TreeH");
729 fTRD->SetTreeAddress();
731 // Get the number of entries in the hit tree
732 // (Number of primary particles creating a hit somewhere)
733 Int_t nTrk = (Int_t) hitTree->GetEntries();
734 AliDebug(1,Form("Found %d tracks",nTrk));
736 // Loop through all the tracks in the tree
737 for (Int_t iTrk = 0; iTrk < nTrk; iTrk++) {
739 gAlice->GetMCApp()->ResetHits();
740 hitTree->GetEvent(iTrk);
743 AliError(Form("No hits array for track = %d",iTrk));
747 // Number of hits for this track
748 nhitTrk = fTRD->Hits()->GetEntriesFast();
751 // Loop through the TRD hits
752 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
757 // Don't analyze test hits
758 if (((Int_t) hit->GetCharge()) != 0) {
760 Int_t trk = hit->Track();
761 Int_t det = hit->GetDetector();
762 Int_t q = hit->GetCharge();
763 Float_t x = hit->X();
764 Float_t y = hit->Y();
765 Float_t z = hit->Z();
766 Float_t time = hit->GetTime();
768 if (nhit[det] == lhit[det]) {
769 // Inititialization of new detector
770 xyz = new Float_t[kNhit*(nhitTrk+lhit[det])];
772 memcpy(xyz,hits[det],sizeof(Float_t)*kNhit*lhit[det]);
775 lhit[det] += nhitTrk;
781 xyz[nhit[det]*kNhit+0] = x;
782 xyz[nhit[det]*kNhit+1] = y;
783 xyz[nhit[det]*kNhit+2] = z;
784 xyz[nhit[det]*kNhit+3] = q;
785 xyz[nhit[det]*kNhit+4] = trk;
786 xyz[nhit[det]*kNhit+5] = time;
791 hit = (AliTRDhit *) fTRD->NextHit();
793 } // for: hits of one track
803 //_____________________________________________________________________________
804 Bool_t AliTRDdigitizer::ConvertHits(Int_t det
805 , const Float_t * const hits
807 , AliTRDarraySignal *signals)
810 // Converts the detectorwise sorted hits to detector signals
813 AliDebug(1,Form("Start converting hits for detector=%d (nhits=%d)",det,nhit));
815 // Number of pads included in the pad response
816 const Int_t kNpad = 3;
817 // Number of track dictionary arrays
818 const Int_t kNdict = AliTRDdigitsManager::kNDict;
819 // Size of the hit vector
820 const Int_t kNhit = 6;
822 // Width of the amplification region
823 const Float_t kAmWidth = AliTRDgeometry::AmThick();
824 // Width of the drift region
825 const Float_t kDrWidth = AliTRDgeometry::DrThick();
826 // Drift + amplification region
827 const Float_t kDrMin = - 0.5 * kAmWidth;
828 const Float_t kDrMax = kDrWidth + 0.5 * kAmWidth;
832 Int_t timeBinTRFend = 1;
836 Double_t padSignal[kNpad];
837 Double_t signalOld[kNpad];
839 AliTRDarrayDictionary *dictionary[kNdict];
841 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
842 AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
843 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
846 AliFatal("Could not get common parameterss");
850 AliFatal("Could not get simulation parameters");
854 AliFatal("Could not get calibration object");
858 // Get the detector wise calibration objects
859 AliTRDCalROC *calVdriftROC = 0;
860 Float_t calVdriftDetValue = 0.0;
861 const AliTRDCalDet *calVdriftDet = calibration->GetVdriftDet();
862 AliTRDCalROC *calT0ROC = 0;
863 Float_t calT0DetValue = 0.0;
864 const AliTRDCalDet *calT0Det = calibration->GetT0Det();
865 Double_t calExBDetValue = 0.0;
866 const AliTRDCalDet *calExBDet = calibration->GetExBDet();
868 if (simParam->TRFOn()) {
869 timeBinTRFend = ((Int_t) (simParam->GetTRFhi()
870 * commonParam->GetSamplingFrequency())) - 1;
873 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
874 Float_t samplingRate = commonParam->GetSamplingFrequency();
875 Float_t elAttachProp = simParam->GetElAttachProp() / 100.0;
877 AliTRDpadPlane *padPlane = fGeo->GetPadPlane(det);
878 Int_t layer = fGeo->GetLayer(det); //update
879 Float_t row0 = padPlane->GetRow0ROC();
880 Int_t nRowMax = padPlane->GetNrows();
881 Int_t nColMax = padPlane->GetNcols();
883 // Create a new array for the signals
884 signals->Allocate(nRowMax,nColMax,nTimeTotal);
886 // Create a new array for the dictionary
887 for (dict = 0; dict < kNdict; dict++) {
888 dictionary[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
889 dictionary[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
892 // Loop through the hits in this detector
893 for (Int_t hit = 0; hit < nhit; hit++) {
895 pos[0] = hits[hit*kNhit+0];
896 pos[1] = hits[hit*kNhit+1];
897 pos[2] = hits[hit*kNhit+2];
898 Float_t q = hits[hit*kNhit+3];
899 Float_t hittime = hits[hit*kNhit+5];
900 Int_t track = ((Int_t) hits[hit*kNhit+4]);
904 // Find the current volume with the geo manager
905 gGeoManager->SetCurrentPoint(pos);
906 gGeoManager->FindNode();
907 if (strstr(gGeoManager->GetPath(),"/UK")) {
911 // Get the calibration objects
912 calVdriftROC = calibration->GetVdriftROC(det);
913 calVdriftDetValue = calVdriftDet->GetValue(det);
914 calT0ROC = calibration->GetT0ROC(det);
915 calT0DetValue = calT0Det->GetValue(det);
916 calExBDetValue = calExBDet->GetValue(det);
918 // Go to the local coordinate system:
919 // loc[0] - col direction in amplification or driftvolume
920 // loc[1] - row direction in amplification or driftvolume
921 // loc[2] - time direction in amplification or driftvolume
922 gGeoManager->MasterToLocal(pos,loc);
924 // Relative to middle of amplification region
925 loc[2] = loc[2] - kDrWidth/2.0 - kAmWidth/2.0;
928 // The driftlength [cm] (w/o diffusion yet !).
929 // It is negative if the hit is between pad plane and anode wires.
930 Double_t driftlength = -1.0 * loc[2];
932 // Stupid patch to take care of TR photons that are absorbed
933 // outside the chamber volume. A real fix would actually need
934 // a more clever implementation of the TR hit generation
936 if ((loc[1] < padPlane->GetRowEndROC()) ||
937 (loc[1] > padPlane->GetRow0ROC())) {
940 if ((driftlength < kDrMin) ||
941 (driftlength > kDrMax)) {
946 // Get row and col of unsmeared electron to retrieve drift velocity
947 // The pad row (z-direction)
948 Int_t rowE = padPlane->GetPadRowNumberROC(loc[1]);
952 Double_t rowOffset = padPlane->GetPadRowOffsetROC(rowE,loc[1]);
953 // The pad column (rphi-direction)
954 Double_t offsetTilt = padPlane->GetTiltOffset(rowOffset);
955 Int_t colE = padPlane->GetPadColNumber(loc[0]+offsetTilt);
959 Double_t colOffset = 0.0;
961 // Normalized drift length
962 Float_t driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
963 Double_t absdriftlength = TMath::Abs(driftlength);
964 if (commonParam->ExBOn()) {
965 absdriftlength /= TMath::Sqrt(1.0 / (1.0 + calExBDetValue*calExBDetValue));
968 // Loop over all electrons of this hit
969 // TR photons produce hits with negative charge
970 Int_t nEl = ((Int_t) TMath::Abs(q));
971 for (Int_t iEl = 0; iEl < nEl; iEl++) {
973 // Now the real local coordinate system of the ROC
974 // column direction: locC
975 // row direction: locR
976 // time direction: locT
977 // locR and locC are identical to the coordinates of the corresponding
978 // volumina of the drift or amplification region.
979 // locT is defined relative to the wire plane (i.e. middle of amplification
980 // region), meaning locT = 0, and is negative for hits coming from the
982 Double_t locC = loc[0];
983 Double_t locR = loc[1];
984 Double_t locT = loc[2];
986 // Electron attachment
987 if (simParam->ElAttachOn()) {
988 if (gRandom->Rndm() < (absdriftlength * elAttachProp)) {
993 // Apply the diffusion smearing
994 if (simParam->DiffusionOn()) {
995 if (!(Diffusion(driftvelocity,absdriftlength,calExBDetValue,locR,locC,locT))) {
1000 // Apply E x B effects (depends on drift direction)
1001 if (commonParam->ExBOn()) {
1002 locC = locC + calExBDetValue * driftlength;
1005 // The electron position after diffusion and ExB in pad coordinates.
1006 // The pad row (z-direction)
1007 rowE = padPlane->GetPadRowNumberROC(locR);
1008 if (rowE < 0) continue;
1009 rowOffset = padPlane->GetPadRowOffsetROC(rowE,locR);
1011 // The pad column (rphi-direction)
1012 offsetTilt = padPlane->GetTiltOffset(rowOffset);
1013 colE = padPlane->GetPadColNumber(locC+offsetTilt);
1014 if (colE < 0) continue;
1015 colOffset = padPlane->GetPadColOffset(colE,locC+offsetTilt);
1017 // Also re-retrieve drift velocity because col and row may have changed
1018 driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
1019 Float_t t0 = calT0DetValue + calT0ROC->GetValue(colE,rowE);
1021 // Convert the position to drift time [mus], using either constant drift velocity or
1022 // time structure of drift cells (non-isochronity, GARFIELD calculation).
1023 // Also add absolute time of hits to take pile-up events into account properly
1025 if (simParam->TimeStructOn()) {
1026 // Get z-position with respect to anode wire
1027 Double_t zz = row0 - locR + padPlane->GetAnodeWireOffset();
1028 zz -= ((Int_t)(2 * zz)) / 2.0;
1032 // Use drift time map (GARFIELD)
1033 drifttime = commonParam->TimeStruct(driftvelocity,0.5*kAmWidth-1.0*locT,zz)
1037 // Use constant drift velocity
1038 drifttime = TMath::Abs(locT) / driftvelocity
1042 // Apply the gas gain including fluctuations
1043 Double_t ggRndm = 0.0;
1045 ggRndm = gRandom->Rndm();
1046 } while (ggRndm <= 0);
1047 Double_t signal = -(simParam->GetGasGain()) * TMath::Log(ggRndm);
1049 // Apply the pad response
1050 if (simParam->PRFOn()) {
1051 // The distance of the electron to the center of the pad
1052 // in units of pad width
1053 Double_t dist = (colOffset - 0.5*padPlane->GetColSize(colE))
1054 / padPlane->GetColSize(colE);
1055 // This is a fixed parametrization, i.e. not dependent on
1056 // calibration values !
1057 if (!(calibration->PadResponse(signal,dist,layer,padSignal))) continue;
1061 padSignal[1] = signal;
1065 // The time bin (always positive), with t0 distortion
1066 Double_t timeBinIdeal = drifttime * samplingRate + t0;
1068 if (TMath::Abs(timeBinIdeal) > 2*nTimeTotal) {
1069 timeBinIdeal = 2 * nTimeTotal;
1071 Int_t timeBinTruncated = ((Int_t) timeBinIdeal);
1072 // The distance of the position to the middle of the timebin
1073 Double_t timeOffset = ((Float_t) timeBinTruncated
1074 + 0.5 - timeBinIdeal) / samplingRate;
1076 // Sample the time response inside the drift region
1077 // + additional time bins before and after.
1078 // The sampling is done always in the middle of the time bin
1079 for (Int_t iTimeBin = TMath::Max(timeBinTruncated,0)
1080 ;iTimeBin < TMath::Min(timeBinTruncated+timeBinTRFend,nTimeTotal)
1083 // Apply the time response
1084 Double_t timeResponse = 1.0;
1085 Double_t crossTalk = 0.0;
1086 Double_t time = (iTimeBin - timeBinTruncated) / samplingRate + timeOffset;
1088 if (simParam->TRFOn()) {
1089 timeResponse = simParam->TimeResponse(time);
1091 if (simParam->CTOn()) {
1092 crossTalk = simParam->CrossTalk(time);
1099 for (iPad = 0; iPad < kNpad; iPad++) {
1101 Int_t colPos = colE + iPad - 1;
1102 if (colPos < 0) continue;
1103 if (colPos >= nColMax) break;
1106 signalOld[iPad] = signals->GetData(rowE,colPos,iTimeBin);
1108 if (colPos != colE) {
1109 // Cross talk added to non-central pads
1110 signalOld[iPad] += padSignal[iPad]
1111 * (timeResponse + crossTalk);
1114 // W/o cross talk at central pad
1115 signalOld[iPad] += padSignal[iPad]
1119 signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);
1121 // Store the track index in the dictionary
1122 // Note: We store index+1 in order to allow the array to be compressed
1123 // Note2: Taking out the +1 in track
1124 if (signalOld[iPad] > 0.0) {
1125 for (dict = 0; dict < kNdict; dict++) {
1126 Int_t oldTrack = dictionary[dict]->GetData(rowE,colPos,iTimeBin);
1127 if (oldTrack == track) break;
1128 if (oldTrack == -1 ) {
1129 dictionary[dict]->SetData(rowE,colPos,iTimeBin,track);
1137 } // Loop: time bins
1139 } // Loop: electrons of a single hit
1143 AliDebug(2,Form("Finished analyzing %d hits",nhit));
1149 //_____________________________________________________________________________
1150 Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
1153 // Convert signals to digits
1156 AliDebug(1,Form("Start converting the signals for detector %d",det));
1159 // Convert the signal array to s-digits
1160 if (!Signal2SDigits(det,signals)) {
1165 // Convert the signal array to digits
1166 if (!Signal2ADC(det,signals)) {
1169 // Run digital processing for digits
1170 RunDigitalProcessing(det);
1173 // Compress the arrays
1174 CompressOutputArrays(det);
1180 //_____________________________________________________________________________
1181 Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
1184 // Converts the sampled electron signals to ADC values for a given chamber
1187 AliDebug(1,Form("Start converting signals to ADC values for detector=%d",det));
1189 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1191 AliFatal("Could not get calibration object");
1195 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1197 AliFatal("Could not get simulation parameters");
1201 // Converts number of electrons to fC
1202 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1205 Double_t coupling = simParam->GetPadCoupling()
1206 * simParam->GetTimeCoupling();
1207 // Electronics conversion factor
1208 Double_t convert = kEl2fC
1209 * simParam->GetChipGain();
1210 // ADC conversion factor
1211 Double_t adcConvert = simParam->GetADCoutRange()
1212 / simParam->GetADCinRange();
1213 // The electronics baseline in mV
1214 Double_t baseline = simParam->GetADCbaseline()
1216 // The electronics baseline in electrons
1217 Double_t baselineEl = baseline
1224 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1225 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1226 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1227 if (fSDigitsManager->GetDigitsParam()->GetNTimeBins(det)) {
1228 nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1231 AliFatal("Could not get number of time bins");
1235 // The gain factor calibration objects
1236 const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1237 AliTRDCalROC *calGainFactorROC = 0x0;
1238 Float_t calGainFactorDetValue = 0.0;
1240 AliTRDarrayADC *digits = 0x0;
1243 AliError(Form("Signals array for detector %d does not exist\n",det));
1246 if (signals->HasData()) {
1247 // Expand the container if neccessary
1251 // Create missing containers
1252 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1255 // Get the container for the digits of this detector
1256 if (fDigitsManager->HasSDigits()) {
1257 AliError("Digits manager has s-digits");
1261 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1262 // Allocate memory space for the digits buffer
1263 if (!digits->HasData()) {
1264 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1267 // Get the calibration objects
1268 calGainFactorROC = calibration->GetGainFactorROC(det);
1269 calGainFactorDetValue = calGainFactorDet->GetValue(det);
1271 // Get the online gain factors
1272 //AliTRDCalOnlineGainTableROC *onlineGainFactorROC
1273 // = calibration->GetOnlineGainTableROC(det);
1275 // Create the digits for this chamber
1276 for (row = 0; row < nRowMax; row++ ) {
1277 for (col = 0; col < nColMax; col++ ) {
1279 // Check whether pad is masked
1280 // Bridged pads are not considered yet!!!
1281 if (calibration->IsPadMasked(det,col,row) ||
1282 calibration->IsPadNotConnected(det,col,row)) {
1287 Float_t padgain = calGainFactorDetValue
1288 * calGainFactorROC->GetValue(col,row);
1290 AliError(Form("Not a valid gain %f, %d %d %d",padgain,det,col,row));
1293 for (time = 0; time < nTimeTotal; time++) {
1295 // Get the signal amplitude
1296 Float_t signalAmp = signals->GetData(row,col,time);
1297 // Pad and time coupling
1298 signalAmp *= coupling;
1300 signalAmp *= padgain;
1302 // Add the noise, starting from minus ADC baseline in electrons
1303 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,simParam->GetNoise())
1307 signalAmp *= convert;
1308 // Add ADC baseline in mV
1309 signalAmp += baseline;
1311 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1312 // signal is larger than fADCinRange
1314 if (signalAmp >= simParam->GetADCinRange()) {
1315 adc = ((Short_t) simParam->GetADCoutRange());
1318 adc = TMath::Nint(signalAmp * adcConvert);
1321 // Saving all digits
1322 digits->SetData(row,col,time,adc);
1333 //_____________________________________________________________________________
1334 Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
1337 // Converts the sampled electron signals to s-digits
1340 AliDebug(1,Form("Start converting signals to s-digits for detector=%d",det));
1342 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1344 AliFatal("Could not get calibration object");
1352 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1353 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1354 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1356 // Get the container for the digits of this detector
1357 if (!fDigitsManager->HasSDigits()) {
1358 AliError("Digits manager has no s-digits");
1362 AliTRDarraySignal *digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1363 // Allocate memory space for the digits buffer
1364 if (!digits->HasData()) {
1365 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1368 // Create the sdigits for this chamber
1369 for (row = 0; row < nRowMax; row++ ) {
1370 for (col = 0; col < nColMax; col++ ) {
1371 for (time = 0; time < nTimeTotal; time++) {
1372 digits->SetData(row,col,time,signals->GetData(row,col,time));
1381 //_____________________________________________________________________________
1382 Bool_t AliTRDdigitizer::Digits2SDigits(AliTRDdigitsManager * const manDig
1383 , AliTRDdigitsManager * const manSDig)
1386 // Converts digits into s-digits. Needed for embedding into real data.
1389 AliDebug(1,"Start converting digits to s-digits");
1392 fGeo = new AliTRDgeometry();
1395 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1397 AliFatal("Could not get calibration object");
1401 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1403 AliFatal("Could not get simulation parameters");
1407 // Converts number of electrons to fC
1408 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1411 Double_t coupling = simParam->GetPadCoupling()
1412 * simParam->GetTimeCoupling();
1413 // Electronics conversion factor
1414 Double_t convert = kEl2fC
1415 * simParam->GetChipGain();
1416 // ADC conversion factor
1417 Double_t adcConvert = simParam->GetADCoutRange()
1418 / simParam->GetADCinRange();
1419 // The electronics baseline in mV
1420 Double_t baseline = simParam->GetADCbaseline()
1422 // The electronics baseline in electrons
1423 //Double_t baselineEl = baseline
1426 // The gainfactor calibration objects
1427 //const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1428 //AliTRDCalROC *calGainFactorROC = 0;
1429 //Float_t calGainFactorDetValue = 0.0;
1435 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1437 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1438 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1439 Int_t nTimeTotal = manDig->GetDigitsParam()->GetNTimeBins(det);
1441 // Get the calibration objects
1442 //calGainFactorROC = calibration->GetGainFactorROC(det);
1443 //calGainFactorDetValue = calGainFactorDet->GetValue(det);
1446 AliTRDarrayADC *digits = (AliTRDarrayADC *) manDig->GetDigits(det);
1448 if (!manSDig->HasSDigits()) {
1449 AliError("SDigits manager has no s-digits");
1453 AliTRDarraySignal *sdigits = (AliTRDarraySignal *) manSDig->GetSDigits(det);
1454 AliTRDarrayDictionary *tracks0 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,0);
1455 AliTRDarrayDictionary *tracks1 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,1);
1456 AliTRDarrayDictionary *tracks2 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,2);
1457 // Allocate memory space for the digits buffer
1458 sdigits->Allocate(nRowMax,nColMax,nTimeTotal);
1459 tracks0->Allocate(nRowMax,nColMax,nTimeTotal);
1460 tracks1->Allocate(nRowMax,nColMax,nTimeTotal);
1461 tracks2->Allocate(nRowMax,nColMax,nTimeTotal);
1463 // Keep the digits param
1464 manSDig->GetDigitsParam()->SetNTimeBinsAll(manDig->GetDigitsParam()->GetNTimeBins(0));
1465 manSDig->GetDigitsParam()->SetADCbaselineAll(manDig->GetDigitsParam()->GetADCbaseline(0));
1467 if (digits->HasData()) {
1471 // Create the sdigits for this chamber
1472 for (row = 0; row < nRowMax; row++ ) {
1473 for (col = 0; col < nColMax; col++ ) {
1476 //Float_t padgain = calGainFactorDetValue
1477 // * calGainFactorROC->GetValue(col,row);
1479 for (time = 0; time < nTimeTotal; time++) {
1481 Short_t adcVal = digits->GetData(row,col,time);
1482 Double_t signal = (Double_t) adcVal;
1483 // ADC -> signal in mV
1484 signal /= adcConvert;
1485 // Subtract baseline in mV
1487 // Signal in mV -> signal in #electrons
1490 //signal /= padgain; // Not needed for real data
1491 // Pad and time coupling
1494 sdigits->SetData(row,col,time,signal);
1495 tracks0->SetData(row,col,time,0);
1496 tracks1->SetData(row,col,time,0);
1497 tracks2->SetData(row,col,time,0);
1506 sdigits->Compress(0);
1507 tracks0->Compress();
1508 tracks1->Compress();
1509 tracks2->Compress();
1511 // No compress just remove
1512 manDig->RemoveDigits(det);
1513 manDig->RemoveDictionaries(det);
1521 //_____________________________________________________________________________
1522 Bool_t AliTRDdigitizer::SDigits2Digits()
1525 // Merges the input s-digits and converts them to normal digits
1528 if (!MergeSDigits()) {
1532 return ConvertSDigits();
1536 //_____________________________________________________________________________
1537 Bool_t AliTRDdigitizer::MergeSDigits()
1540 // Merges the input s-digits:
1541 // - The amplitude of the different inputs are summed up.
1542 // - Of the track IDs from the input dictionaries only one is
1543 // kept for each input. This works for maximal 3 different merged inputs.
1546 // Number of track dictionary arrays
1547 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1549 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1551 AliFatal("Could not get simulation parameters");
1555 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1557 AliFatal("Could not get calibration object");
1564 AliTRDarraySignal *digitsA;
1565 AliTRDarraySignal *digitsB;
1566 AliTRDarrayDictionary *dictionaryA[kNDict];
1567 AliTRDarrayDictionary *dictionaryB[kNDict];
1569 AliTRDdigitsManager *mergeSDigitsManager = 0x0;
1570 // Get the first s-digits
1571 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1572 if (!fSDigitsManager) {
1573 AliError("No SDigits manager");
1577 // Loop through the other sets of s-digits
1578 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(fSDigitsManager);
1580 if (mergeSDigitsManager) {
1581 AliDebug(1,Form("Merge %d input files.",fSDigitsManagerList->GetSize()));
1584 AliDebug(1,"Only one input file.");
1589 while (mergeSDigitsManager) {
1593 // Loop through the detectors
1594 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1596 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet);
1597 if (mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet) != nTimeTotal) {
1598 AliError(Form("Mismatch in the number of time bins [%d,%d] in detector %d"
1600 ,mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet)
1605 Int_t nRowMax = fGeo->GetPadPlane(iDet)->GetNrows();
1606 Int_t nColMax = fGeo->GetPadPlane(iDet)->GetNcols();
1608 // Loop through the pixels of one detector and add the signals
1609 digitsA = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(iDet);
1610 digitsB = (AliTRDarraySignal *) mergeSDigitsManager->GetSDigits(iDet);
1612 if (!digitsA->HasData()) continue;
1614 if (!digitsB->HasData()) continue;
1616 for (iDict = 0; iDict < kNDict; iDict++) {
1617 dictionaryA[iDict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(iDet,iDict);
1618 dictionaryB[iDict] = (AliTRDarrayDictionary *) mergeSDigitsManager->GetDictionary(iDet,iDict);
1619 dictionaryA[iDict]->Expand();
1620 dictionaryB[iDict]->Expand();
1623 // Merge only detectors that contain a signal
1624 Bool_t doMerge = kTRUE;
1625 if (fMergeSignalOnly) {
1626 if (digitsA->GetOverThreshold(0) == 0) {
1633 AliDebug(1,Form("Merge detector %d of input no.%d",iDet,iMerge+1));
1635 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1636 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1637 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1639 // Add the amplitudes of the summable digits
1640 Float_t ampA = digitsA->GetData(iRow,iCol,iTime);
1641 Float_t ampB = digitsB->GetData(iRow,iCol,iTime);
1643 digitsA->SetData(iRow,iCol,iTime,ampA);
1645 // Add the mask to the track id if defined.
1646 for (iDict = 0; iDict < kNDict; iDict++) {
1647 Int_t trackB = dictionaryB[iDict]->GetData(iRow,iCol,iTime);
1648 if ((fMasks) && (trackB > 0)) {
1649 for (jDict = 0; jDict < kNDict; jDict++) {
1650 Int_t trackA = dictionaryA[iDict]->GetData(iRow,iCol,iTime);
1652 trackA = trackB + fMasks[iMerge];
1653 dictionaryA[iDict]->SetData(iRow,iCol,iTime,trackA);
1654 } // if: track A == 0
1656 } // if: fMasks and trackB > 0
1665 mergeSDigitsManager->RemoveDigits(iDet);
1666 mergeSDigitsManager->RemoveDictionaries(iDet);
1669 digitsA->Compress(0);
1670 for (iDict = 0; iDict < kNDict; iDict++) {
1671 dictionaryA[iDict]->Compress();
1677 // The next set of s-digits
1678 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(mergeSDigitsManager);
1680 } // while: mergeDigitsManagers
1686 //_____________________________________________________________________________
1687 Bool_t AliTRDdigitizer::ConvertSDigits()
1690 // Converts s-digits to normal digits
1693 AliTRDarraySignal *digitsIn = 0x0;
1695 if (!fSDigitsManager->HasSDigits()) {
1696 AliError("No s-digits in digits manager");
1700 // Loop through the detectors
1701 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1703 // Get the merged s-digits (signals)
1704 digitsIn = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(det);
1705 if (!digitsIn->HasData()) {
1706 AliDebug(2,Form("No digits for det=%d",det));
1710 // Convert the merged sdigits to digits
1711 if (!Signal2ADC(det,digitsIn)) {
1715 // Copy the dictionary information to the output array
1716 if (!CopyDictionary(det)) {
1721 fSDigitsManager->RemoveDigits(det);
1722 fSDigitsManager->RemoveDictionaries(det);
1724 // Run digital processing
1725 RunDigitalProcessing(det);
1727 // Compress the arrays
1728 CompressOutputArrays(det);
1730 } // for: detector numbers
1732 if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
1733 if (trklLoader->Tree())
1734 trklLoader->WriteData("OVERWRITE");
1737 // Save the values for the raw data headers
1738 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
1739 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
1742 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
1744 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
1750 //_____________________________________________________________________________
1751 Bool_t AliTRDdigitizer::CopyDictionary(Int_t det)
1754 // Copies the dictionary information from the s-digits arrays
1755 // to the output arrays
1758 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1760 AliFatal("Could not get calibration object");
1764 AliDebug(1,Form("Start copying dictionaries for detector=%d",det));
1766 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1767 AliTRDarrayDictionary *dictionaryIn[kNDict];
1768 AliTRDarrayDictionary *dictionaryOut[kNDict];
1770 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1771 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1772 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1779 for (dict = 0; dict < kNDict; dict++) {
1781 dictionaryIn[dict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(det,dict);
1782 dictionaryIn[dict]->Expand();
1783 dictionaryOut[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1784 dictionaryOut[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
1786 for (row = 0; row < nRowMax; row++) {
1787 for (col = 0; col < nColMax; col++) {
1788 for (time = 0; time < nTimeTotal; time++) {
1789 Int_t track = dictionaryIn[dict]->GetData(row,col,time);
1790 dictionaryOut[dict]->SetData(row,col,time,track);
1795 } // for: dictionaries
1801 //_____________________________________________________________________________
1802 void AliTRDdigitizer::CompressOutputArrays(Int_t det)
1805 // Compress the output arrays
1808 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1809 AliTRDarrayDictionary *dictionary = 0x0;
1814 AliTRDarrayADC *digits = 0x0;
1815 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1820 AliTRDarraySignal *digits = 0x0;
1821 digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1822 digits->Compress(0);
1825 for (Int_t dict = 0; dict < kNDict; dict++) {
1826 dictionary = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1827 dictionary->Compress();
1834 //_____________________________________________________________________________
1835 Bool_t AliTRDdigitizer::WriteDigits() const
1838 // Writes out the TRD-digits and the dictionaries
1842 fRunLoader->CdGAFile();
1844 // Store the digits and the dictionary in the tree
1845 return fDigitsManager->WriteDigits();
1849 //_____________________________________________________________________________
1850 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1853 // Initializes the output branches
1859 AliError("Run Loader is NULL");
1863 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
1865 AliError("Can not get TRD loader from Run Loader");
1872 // If we produce SDigits
1873 tree = loader->TreeS();
1875 loader->MakeTree("S");
1876 tree = loader->TreeS();
1880 // If we produce Digits
1881 tree = loader->TreeD();
1883 loader->MakeTree("D");
1884 tree = loader->TreeD();
1887 fDigitsManager->SetEvent(iEvent);
1888 fDigitsManager->MakeBranch(tree);
1892 //_____________________________________________________________________________
1893 Int_t AliTRDdigitizer::Diffusion(Float_t vdrift, Double_t absdriftlength
1895 , Double_t &lRow, Double_t &lCol, Double_t &lTime)
1898 // Applies the diffusion smearing to the position of a single electron.
1899 // Depends on absolute drift length.
1902 Float_t diffL = 0.0;
1903 Float_t diffT = 0.0;
1905 if (AliTRDCommonParam::Instance()->GetDiffCoeff(diffL,diffT,vdrift)) {
1907 Float_t driftSqrt = TMath::Sqrt(absdriftlength);
1908 Float_t sigmaT = driftSqrt * diffT;
1909 Float_t sigmaL = driftSqrt * diffL;
1910 lRow = gRandom->Gaus(lRow ,sigmaT);
1911 if (AliTRDCommonParam::Instance()->ExBOn()) {
1912 lCol = gRandom->Gaus(lCol ,sigmaT * 1.0 / (1.0 + exbvalue*exbvalue));
1913 lTime = gRandom->Gaus(lTime,sigmaL * 1.0 / (1.0 + exbvalue*exbvalue));
1916 lCol = gRandom->Gaus(lCol ,sigmaT);
1917 lTime = gRandom->Gaus(lTime,sigmaL);
1931 //_____________________________________________________________________________
1932 void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
1935 // Run the digital processing in the TRAP
1938 AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();
1940 AliTRDarrayADC *digits = fDigitsManager->GetDigits(det);
1944 //Call the methods in the mcm class using the temporary array as input
1945 for(Int_t rob = 0; rob < digits->GetNrow() / 2; rob++)
1947 for(Int_t mcm = 0; mcm < 16; mcm++)
1949 fMcmSim->Init(det, rob, mcm);
1950 fMcmSim->SetDataByPad(digits, fDigitsManager);
1952 if (feeParam->GetTracklet()) {
1953 fMcmSim->Tracklet();
1954 fMcmSim->StoreTracklets();
1956 fMcmSim->ZSMapping();
1957 fMcmSim->WriteData(digits);