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 if (fSDigitsManagerList) {
201 fSDigitsManagerList->Delete();
202 delete fSDigitsManagerList;
204 fSDigitsManagerList = 0;
217 //_____________________________________________________________________________
218 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
221 // Assignment operator
225 ((AliTRDdigitizer &) d).Copy(*this);
232 //_____________________________________________________________________________
233 void AliTRDdigitizer::Copy(TObject &d) const
239 ((AliTRDdigitizer &) d).fRunLoader = 0;
240 ((AliTRDdigitizer &) d).fDigitsManager = 0;
241 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
242 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
243 ((AliTRDdigitizer &) d).fTRD = 0;
244 ((AliTRDdigitizer &) d).fGeo = 0;
245 ((AliTRDdigitizer &) d).fEvent = 0;
246 ((AliTRDdigitizer &) d).fMasks = 0;
247 ((AliTRDdigitizer &) d).fCompress = fCompress;
248 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
249 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
253 //_____________________________________________________________________________
254 void AliTRDdigitizer::Digitize(const Option_t* option)
257 // Executes the merging
262 AliTRDdigitsManager *sdigitsManager;
264 TString optionString = option;
265 if (optionString.Contains("deb")) {
266 AliLog::SetClassDebugLevel("AliTRDdigitizer",1);
267 AliInfo("Called with debug option");
270 // The AliRoot file is already connected by the manager
271 AliRunLoader *inrl = 0x0;
274 AliDebug(1,"AliRun object found on file.");
277 inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0));
279 gAlice = inrl->GetAliRun();
281 AliError("Could not find AliRun object.");
286 Int_t nInput = fDigInput->GetNinputs();
287 fMasks = new Int_t[nInput];
288 for (iInput = 0; iInput < nInput; iInput++) {
289 fMasks[iInput] = fDigInput->GetMask(iInput);
296 AliRunLoader *orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
298 if (InitDetector()) {
300 AliLoader *ogime = orl->GetLoader("TRDLoader");
304 // If we produce SDigits
305 tree = ogime->TreeS();
307 ogime->MakeTree("S");
308 tree = ogime->TreeS();
312 // If we produce Digits
313 tree = ogime->TreeD();
315 ogime->MakeTree("D");
316 tree = ogime->TreeD();
324 for (iInput = 0; iInput < nInput; iInput++) {
326 AliDebug(1,Form("Add input stream %d",iInput));
328 // Check if the input tree exists
329 inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
330 AliLoader *gime = inrl->GetLoader("TRDLoader");
332 TTree *treees = gime->TreeS();
334 if (gime->LoadSDigits()) {
335 AliError(Form("Error Occured while loading S. Digits for input %d.",iInput));
338 treees = gime->TreeS();
342 AliError(Form("Input stream %d does not exist",iInput));
346 // Read the s-digits via digits manager
347 sdigitsManager = new AliTRDdigitsManager();
348 sdigitsManager->SetSDigits(kTRUE);
350 AliRunLoader *rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
351 AliLoader *gimme = rl->GetLoader("TRDLoader");
354 gimme->LoadSDigits();
357 sdigitsManager->ReadDigits(gimme->TreeS());
359 // Add the s-digits to the input list
360 AddSDigitsManager(sdigitsManager);
364 // Convert the s-digits to normal digits
365 AliDebug(1,"Do the conversion");
369 AliDebug(1,"Write the digits");
370 fDigitsManager->WriteDigits();
376 DeleteSDigitsManager();
382 //_____________________________________________________________________________
383 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
386 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
388 // Connect the AliRoot file containing Geometry, Kine, and Hits
391 TString evfoldname = AliConfig::GetDefaultEventFolderName();
393 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
395 fRunLoader = AliRunLoader::Open(file,evfoldname,"UPDATE");
398 AliError(Form("Can not open session for file %s.",file));
402 if (!fRunLoader->GetAliRun()) {
403 fRunLoader->LoadgAlice();
405 gAlice = fRunLoader->GetAliRun();
408 AliDebug(1,"AliRun object found on file.");
411 AliError("Could not find AliRun object.");
417 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
419 AliError("Can not get TRD loader from Run Loader");
423 if (InitDetector()) {
426 // If we produce SDigits
427 tree = loader->TreeS();
429 loader->MakeTree("S");
430 tree = loader->TreeS();
434 // If we produce Digits
435 tree = loader->TreeD();
437 loader->MakeTree("D");
438 tree = loader->TreeD();
441 return MakeBranch(tree);
449 //_____________________________________________________________________________
450 Bool_t AliTRDdigitizer::Open(AliRunLoader * const runLoader, Int_t nEvent)
453 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
455 // Connect the AliRoot file containing Geometry, Kine, and Hits
458 fRunLoader = runLoader;
460 AliError("RunLoader does not exist");
464 if (!fRunLoader->GetAliRun()) {
465 fRunLoader->LoadgAlice();
467 gAlice = fRunLoader->GetAliRun();
470 AliDebug(1,"AliRun object found on file.");
473 AliError("Could not find AliRun object.");
479 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
481 AliError("Can not get TRD loader from Run Loader");
485 if (InitDetector()) {
488 // If we produce SDigits
489 tree = loader->TreeS();
491 loader->MakeTree("S");
492 tree = loader->TreeS();
496 // If we produce Digits
497 tree = loader->TreeD();
499 loader->MakeTree("D");
500 tree = loader->TreeD();
503 return MakeBranch(tree);
511 //_____________________________________________________________________________
512 Bool_t AliTRDdigitizer::InitDetector()
515 // Sets the pointer to the TRD detector and the geometry
518 // Get the pointer to the detector class and check for version 1
519 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
521 AliFatal("No TRD module found");
524 if (fTRD->IsVersion() != 1) {
525 AliFatal("TRD must be version 1 (slow simulator)");
530 fGeo = new AliTRDgeometry();
532 // Create a digits manager
533 if (fDigitsManager) {
534 delete fDigitsManager;
536 fDigitsManager = new AliTRDdigitsManager();
537 fDigitsManager->SetSDigits(fSDigits);
538 fDigitsManager->CreateArrays();
539 fDigitsManager->SetEvent(fEvent);
541 // The list for the input s-digits manager to be merged
542 if (fSDigitsManagerList) {
543 fSDigitsManagerList->Delete();
546 fSDigitsManagerList = new TList();
553 //_____________________________________________________________________________
554 Bool_t AliTRDdigitizer::MakeBranch(TTree *tree) const
557 // Create the branches for the digits array
560 return fDigitsManager->MakeBranch(tree);
564 //_____________________________________________________________________________
565 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
568 // Add a digits manager for s-digits to the input list.
571 fSDigitsManagerList->Add(man);
575 //_____________________________________________________________________________
576 void AliTRDdigitizer::DeleteSDigitsManager()
579 // Removes digits manager from the input list.
582 fSDigitsManagerList->Delete();
586 //_____________________________________________________________________________
587 Bool_t AliTRDdigitizer::MakeDigits()
593 AliDebug(1,"Start creating digits");
596 AliError("No geometry defined");
600 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
602 AliFatal("Could not get calibration object");
606 const Int_t kNdet = AliTRDgeometry::Ndet();
608 Float_t **hits = new Float_t*[kNdet];
609 Int_t *nhit = new Int_t[kNdet];
610 memset(nhit,0,kNdet*sizeof(Int_t));
612 AliTRDarraySignal *signals = 0x0;
614 // Check the number of time bins from simParam against OCDB,
615 // if OCDB value is not supposed to be used.
616 // As default, the value from OCDB is taken
617 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
618 if (calibration->GetNumberOfTimeBinsDCS() != AliTRDSimParam::Instance()->GetNTimeBins()) {
619 AliWarning(Form("Number of time bins is different to OCDB value [SIM=%d, OCDB=%d]"
620 ,AliTRDSimParam::Instance()->GetNTimeBins()
621 ,calibration->GetNumberOfTimeBinsDCS()));
623 // Save the values for the raw data headers
624 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
627 // Save the values for the raw data headers
628 if (calibration->GetNumberOfTimeBinsDCS() == -1) {
629 AliFatal("No useful DCS information available for this run!");
631 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(calibration->GetNumberOfTimeBinsDCS());
634 // Save the values for the raw data headers
635 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
637 // Sort all hits according to detector number
638 if (!SortHits(hits,nhit)) {
639 AliError("Sorting hits failed");
645 // Loop through all detectors
646 for (Int_t det = 0; det < kNdet; det++) {
648 // Detectors that are switched off, not installed, etc.
649 if ((!calibration->IsChamberNoData(det)) &&
650 ( fGeo->ChamberInGeometry(det)) &&
653 signals = new AliTRDarraySignal();
655 // Convert the hits of the current detector to detector signals
656 if (!ConvertHits(det,hits[det],nhit[det],signals)) {
657 AliError(Form("Conversion of hits failed for detector=%d",det));
665 // Convert the detector signals to digits or s-digits
666 if (!ConvertSignals(det,signals)) {
667 AliError(Form("Conversion of signals failed for detector=%d",det));
675 // Delete the signals array
679 } // if: detector status
686 if (AliDataLoader *trklLoader
687 = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
688 if (trklLoader->Tree())
689 trklLoader->WriteData("OVERWRITE");
700 //_____________________________________________________________________________
701 Bool_t AliTRDdigitizer::SortHits(Float_t **hits, Int_t *nhit)
704 // Read all the hits and sorts them according to detector number
705 // in the output array <hits>.
708 AliDebug(1,"Start sorting hits");
710 const Int_t kNdet = AliTRDgeometry::Ndet();
711 // Size of the hit vector
712 const Int_t kNhit = 6;
717 Int_t *lhit = new Int_t[kNdet];
718 memset(lhit,0,kNdet*sizeof(Int_t));
720 for (Int_t det = 0; det < kNdet; det++) {
724 AliLoader *gimme = fRunLoader->GetLoader("TRDLoader");
725 if (!gimme->TreeH()) {
728 TTree *hitTree = gimme->TreeH();
729 if (hitTree == 0x0) {
730 AliError("Can not get TreeH");
734 fTRD->SetTreeAddress();
736 // Get the number of entries in the hit tree
737 // (Number of primary particles creating a hit somewhere)
738 Int_t nTrk = (Int_t) hitTree->GetEntries();
739 AliDebug(1,Form("Found %d tracks",nTrk));
741 // Loop through all the tracks in the tree
742 for (Int_t iTrk = 0; iTrk < nTrk; iTrk++) {
744 gAlice->GetMCApp()->ResetHits();
745 hitTree->GetEvent(iTrk);
748 AliError(Form("No hits array for track = %d",iTrk));
752 // Number of hits for this track
753 nhitTrk = fTRD->Hits()->GetEntriesFast();
756 // Loop through the TRD hits
757 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
762 // Don't analyze test hits
763 if (((Int_t) hit->GetCharge()) != 0) {
765 Int_t trk = hit->Track();
766 Int_t det = hit->GetDetector();
767 Int_t q = hit->GetCharge();
768 Float_t x = hit->X();
769 Float_t y = hit->Y();
770 Float_t z = hit->Z();
771 Float_t time = hit->GetTime();
773 if (nhit[det] == lhit[det]) {
774 // Inititialization of new detector
775 xyz = new Float_t[kNhit*(nhitTrk+lhit[det])];
777 memcpy(xyz,hits[det],sizeof(Float_t)*kNhit*lhit[det]);
780 lhit[det] += nhitTrk;
786 xyz[nhit[det]*kNhit+0] = x;
787 xyz[nhit[det]*kNhit+1] = y;
788 xyz[nhit[det]*kNhit+2] = z;
789 xyz[nhit[det]*kNhit+3] = q;
790 xyz[nhit[det]*kNhit+4] = trk;
791 xyz[nhit[det]*kNhit+5] = time;
796 hit = (AliTRDhit *) fTRD->NextHit();
798 } // for: hits of one track
808 //_____________________________________________________________________________
809 Bool_t AliTRDdigitizer::ConvertHits(Int_t det
810 , const Float_t * const hits
812 , AliTRDarraySignal *signals)
815 // Converts the detectorwise sorted hits to detector signals
818 AliDebug(1,Form("Start converting hits for detector=%d (nhits=%d)",det,nhit));
820 // Number of pads included in the pad response
821 const Int_t kNpad = 3;
822 // Number of track dictionary arrays
823 const Int_t kNdict = AliTRDdigitsManager::kNDict;
824 // Size of the hit vector
825 const Int_t kNhit = 6;
827 // Width of the amplification region
828 const Float_t kAmWidth = AliTRDgeometry::AmThick();
829 // Width of the drift region
830 const Float_t kDrWidth = AliTRDgeometry::DrThick();
831 // Drift + amplification region
832 const Float_t kDrMin = - 0.5 * kAmWidth;
833 const Float_t kDrMax = kDrWidth + 0.5 * kAmWidth;
837 Int_t timeBinTRFend = 1;
841 Double_t padSignal[kNpad];
842 Double_t signalOld[kNpad];
844 AliTRDarrayDictionary *dictionary[kNdict];
846 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
847 AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
848 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
851 AliFatal("Could not get common parameterss");
855 AliFatal("Could not get simulation parameters");
859 AliFatal("Could not get calibration object");
863 // Get the detector wise calibration objects
864 AliTRDCalROC *calVdriftROC = 0;
865 Float_t calVdriftDetValue = 0.0;
866 const AliTRDCalDet *calVdriftDet = calibration->GetVdriftDet();
867 AliTRDCalROC *calT0ROC = 0;
868 Float_t calT0DetValue = 0.0;
869 const AliTRDCalDet *calT0Det = calibration->GetT0Det();
870 Double_t calExBDetValue = 0.0;
871 const AliTRDCalDet *calExBDet = calibration->GetExBDet();
873 if (simParam->TRFOn()) {
874 timeBinTRFend = ((Int_t) (simParam->GetTRFhi()
875 * commonParam->GetSamplingFrequency())) - 1;
878 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
879 Float_t samplingRate = commonParam->GetSamplingFrequency();
880 Float_t elAttachProp = simParam->GetElAttachProp() / 100.0;
882 AliTRDpadPlane *padPlane = fGeo->GetPadPlane(det);
883 Int_t layer = fGeo->GetLayer(det); //update
884 Float_t row0 = padPlane->GetRow0ROC();
885 Int_t nRowMax = padPlane->GetNrows();
886 Int_t nColMax = padPlane->GetNcols();
888 // Create a new array for the signals
889 signals->Allocate(nRowMax,nColMax,nTimeTotal);
891 // Create a new array for the dictionary
892 for (dict = 0; dict < kNdict; dict++) {
893 dictionary[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
894 dictionary[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
897 // Loop through the hits in this detector
898 for (Int_t hit = 0; hit < nhit; hit++) {
900 pos[0] = hits[hit*kNhit+0];
901 pos[1] = hits[hit*kNhit+1];
902 pos[2] = hits[hit*kNhit+2];
903 Float_t q = hits[hit*kNhit+3];
904 Float_t hittime = hits[hit*kNhit+5];
905 Int_t track = ((Int_t) hits[hit*kNhit+4]);
909 // Find the current volume with the geo manager
910 gGeoManager->SetCurrentPoint(pos);
911 gGeoManager->FindNode();
912 if (strstr(gGeoManager->GetPath(),"/UK")) {
916 // Get the calibration objects
917 calVdriftROC = calibration->GetVdriftROC(det);
918 calVdriftDetValue = calVdriftDet->GetValue(det);
919 calT0ROC = calibration->GetT0ROC(det);
920 calT0DetValue = calT0Det->GetValue(det);
921 calExBDetValue = calExBDet->GetValue(det);
923 // Go to the local coordinate system:
924 // loc[0] - col direction in amplification or driftvolume
925 // loc[1] - row direction in amplification or driftvolume
926 // loc[2] - time direction in amplification or driftvolume
927 gGeoManager->MasterToLocal(pos,loc);
929 // Relative to middle of amplification region
930 loc[2] = loc[2] - kDrWidth/2.0 - kAmWidth/2.0;
933 // The driftlength [cm] (w/o diffusion yet !).
934 // It is negative if the hit is between pad plane and anode wires.
935 Double_t driftlength = -1.0 * loc[2];
937 // Stupid patch to take care of TR photons that are absorbed
938 // outside the chamber volume. A real fix would actually need
939 // a more clever implementation of the TR hit generation
941 if ((loc[1] < padPlane->GetRowEndROC()) ||
942 (loc[1] > padPlane->GetRow0ROC())) {
945 if ((driftlength < kDrMin) ||
946 (driftlength > kDrMax)) {
951 // Get row and col of unsmeared electron to retrieve drift velocity
952 // The pad row (z-direction)
953 Int_t rowE = padPlane->GetPadRowNumberROC(loc[1]);
957 Double_t rowOffset = padPlane->GetPadRowOffsetROC(rowE,loc[1]);
958 // The pad column (rphi-direction)
959 Double_t offsetTilt = padPlane->GetTiltOffset(rowOffset);
960 Int_t colE = padPlane->GetPadColNumber(loc[0]+offsetTilt);
964 Double_t colOffset = 0.0;
966 // Normalized drift length
967 Float_t driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
968 Double_t absdriftlength = TMath::Abs(driftlength);
969 if (commonParam->ExBOn()) {
970 absdriftlength /= TMath::Sqrt(1.0 / (1.0 + calExBDetValue*calExBDetValue));
973 // Loop over all electrons of this hit
974 // TR photons produce hits with negative charge
975 Int_t nEl = ((Int_t) TMath::Abs(q));
976 for (Int_t iEl = 0; iEl < nEl; iEl++) {
978 // Now the real local coordinate system of the ROC
979 // column direction: locC
980 // row direction: locR
981 // time direction: locT
982 // locR and locC are identical to the coordinates of the corresponding
983 // volumina of the drift or amplification region.
984 // locT is defined relative to the wire plane (i.e. middle of amplification
985 // region), meaning locT = 0, and is negative for hits coming from the
987 Double_t locC = loc[0];
988 Double_t locR = loc[1];
989 Double_t locT = loc[2];
991 // Electron attachment
992 if (simParam->ElAttachOn()) {
993 if (gRandom->Rndm() < (absdriftlength * elAttachProp)) {
998 // Apply the diffusion smearing
999 if (simParam->DiffusionOn()) {
1000 if (!(Diffusion(driftvelocity,absdriftlength,calExBDetValue,locR,locC,locT))) {
1005 // Apply E x B effects (depends on drift direction)
1006 if (commonParam->ExBOn()) {
1007 locC = locC + calExBDetValue * driftlength;
1010 // The electron position after diffusion and ExB in pad coordinates.
1011 // The pad row (z-direction)
1012 rowE = padPlane->GetPadRowNumberROC(locR);
1013 if (rowE < 0) continue;
1014 rowOffset = padPlane->GetPadRowOffsetROC(rowE,locR);
1016 // The pad column (rphi-direction)
1017 offsetTilt = padPlane->GetTiltOffset(rowOffset);
1018 colE = padPlane->GetPadColNumber(locC+offsetTilt);
1019 if (colE < 0) continue;
1020 colOffset = padPlane->GetPadColOffset(colE,locC+offsetTilt);
1022 // Also re-retrieve drift velocity because col and row may have changed
1023 driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
1024 Float_t t0 = calT0DetValue + calT0ROC->GetValue(colE,rowE);
1026 // Convert the position to drift time [mus], using either constant drift velocity or
1027 // time structure of drift cells (non-isochronity, GARFIELD calculation).
1028 // Also add absolute time of hits to take pile-up events into account properly
1030 if (simParam->TimeStructOn()) {
1031 // Get z-position with respect to anode wire
1032 Double_t zz = row0 - locR + padPlane->GetAnodeWireOffset();
1033 zz -= ((Int_t)(2 * zz)) / 2.0;
1037 // Use drift time map (GARFIELD)
1038 drifttime = commonParam->TimeStruct(driftvelocity,0.5*kAmWidth-1.0*locT,zz)
1042 // Use constant drift velocity
1043 drifttime = TMath::Abs(locT) / driftvelocity
1047 // Apply the gas gain including fluctuations
1048 Double_t ggRndm = 0.0;
1050 ggRndm = gRandom->Rndm();
1051 } while (ggRndm <= 0);
1052 Double_t signal = -(simParam->GetGasGain()) * TMath::Log(ggRndm);
1054 // Apply the pad response
1055 if (simParam->PRFOn()) {
1056 // The distance of the electron to the center of the pad
1057 // in units of pad width
1058 Double_t dist = (colOffset - 0.5*padPlane->GetColSize(colE))
1059 / padPlane->GetColSize(colE);
1060 // This is a fixed parametrization, i.e. not dependent on
1061 // calibration values !
1062 if (!(calibration->PadResponse(signal,dist,layer,padSignal))) continue;
1066 padSignal[1] = signal;
1070 // The time bin (always positive), with t0 distortion
1071 Double_t timeBinIdeal = drifttime * samplingRate + t0;
1073 if (TMath::Abs(timeBinIdeal) > 2*nTimeTotal) {
1074 timeBinIdeal = 2 * nTimeTotal;
1076 Int_t timeBinTruncated = ((Int_t) timeBinIdeal);
1077 // The distance of the position to the middle of the timebin
1078 Double_t timeOffset = ((Float_t) timeBinTruncated
1079 + 0.5 - timeBinIdeal) / samplingRate;
1081 // Sample the time response inside the drift region
1082 // + additional time bins before and after.
1083 // The sampling is done always in the middle of the time bin
1084 for (Int_t iTimeBin = TMath::Max(timeBinTruncated,0)
1085 ;iTimeBin < TMath::Min(timeBinTruncated+timeBinTRFend,nTimeTotal)
1088 // Apply the time response
1089 Double_t timeResponse = 1.0;
1090 Double_t crossTalk = 0.0;
1091 Double_t time = (iTimeBin - timeBinTruncated) / samplingRate + timeOffset;
1093 if (simParam->TRFOn()) {
1094 timeResponse = simParam->TimeResponse(time);
1096 if (simParam->CTOn()) {
1097 crossTalk = simParam->CrossTalk(time);
1104 for (iPad = 0; iPad < kNpad; iPad++) {
1106 Int_t colPos = colE + iPad - 1;
1107 if (colPos < 0) continue;
1108 if (colPos >= nColMax) break;
1111 signalOld[iPad] = signals->GetData(rowE,colPos,iTimeBin);
1113 if (colPos != colE) {
1114 // Cross talk added to non-central pads
1115 signalOld[iPad] += padSignal[iPad]
1116 * (timeResponse + crossTalk);
1119 // W/o cross talk at central pad
1120 signalOld[iPad] += padSignal[iPad]
1124 signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);
1126 // Store the track index in the dictionary
1127 // Note: We store index+1 in order to allow the array to be compressed
1128 // Note2: Taking out the +1 in track
1129 if (signalOld[iPad] > 0.0) {
1130 for (dict = 0; dict < kNdict; dict++) {
1131 Int_t oldTrack = dictionary[dict]->GetData(rowE,colPos,iTimeBin);
1132 if (oldTrack == track) break;
1133 if (oldTrack == -1 ) {
1134 dictionary[dict]->SetData(rowE,colPos,iTimeBin,track);
1142 } // Loop: time bins
1144 } // Loop: electrons of a single hit
1148 AliDebug(2,Form("Finished analyzing %d hits",nhit));
1154 //_____________________________________________________________________________
1155 Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
1158 // Convert signals to digits
1161 AliDebug(1,Form("Start converting the signals for detector %d",det));
1164 // Convert the signal array to s-digits
1165 if (!Signal2SDigits(det,signals)) {
1170 // Convert the signal array to digits
1171 if (!Signal2ADC(det,signals)) {
1174 // Run digital processing for digits
1175 RunDigitalProcessing(det);
1178 // Compress the arrays
1179 CompressOutputArrays(det);
1185 //_____________________________________________________________________________
1186 Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
1189 // Converts the sampled electron signals to ADC values for a given chamber
1192 AliDebug(1,Form("Start converting signals to ADC values for detector=%d",det));
1194 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1196 AliFatal("Could not get calibration object");
1200 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1202 AliFatal("Could not get simulation parameters");
1206 // Converts number of electrons to fC
1207 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1210 Double_t coupling = simParam->GetPadCoupling()
1211 * simParam->GetTimeCoupling();
1212 // Electronics conversion factor
1213 Double_t convert = kEl2fC
1214 * simParam->GetChipGain();
1215 // ADC conversion factor
1216 Double_t adcConvert = simParam->GetADCoutRange()
1217 / simParam->GetADCinRange();
1218 // The electronics baseline in mV
1219 Double_t baseline = simParam->GetADCbaseline()
1221 // The electronics baseline in electrons
1222 Double_t baselineEl = baseline
1229 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1230 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1231 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1232 if (fSDigitsManager->GetDigitsParam()->GetNTimeBins(det)) {
1233 nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1236 AliFatal("Could not get number of time bins");
1240 // The gain factor calibration objects
1241 const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1242 AliTRDCalROC *calGainFactorROC = 0x0;
1243 Float_t calGainFactorDetValue = 0.0;
1245 AliTRDarrayADC *digits = 0x0;
1248 AliError(Form("Signals array for detector %d does not exist\n",det));
1251 if (signals->HasData()) {
1252 // Expand the container if neccessary
1256 // Create missing containers
1257 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1260 // Get the container for the digits of this detector
1261 if (fDigitsManager->HasSDigits()) {
1262 AliError("Digits manager has s-digits");
1266 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1267 // Allocate memory space for the digits buffer
1268 if (!digits->HasData()) {
1269 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1272 // Get the calibration objects
1273 calGainFactorROC = calibration->GetGainFactorROC(det);
1274 calGainFactorDetValue = calGainFactorDet->GetValue(det);
1276 // Create the digits for this chamber
1277 for (row = 0; row < nRowMax; row++ ) {
1278 for (col = 0; col < nColMax; col++ ) {
1280 // Check whether pad is masked
1281 // Bridged pads are not considered yet!!!
1282 if (calibration->IsPadMasked(det,col,row) ||
1283 calibration->IsPadNotConnected(det,col,row)) {
1288 Float_t padgain = calGainFactorDetValue
1289 * calGainFactorROC->GetValue(col,row);
1291 AliError(Form("Not a valid gain %f, %d %d %d",padgain,det,col,row));
1294 for (time = 0; time < nTimeTotal; time++) {
1296 // Get the signal amplitude
1297 Float_t signalAmp = signals->GetData(row,col,time);
1298 // Pad and time coupling
1299 signalAmp *= coupling;
1301 signalAmp *= padgain;
1303 // Add the noise, starting from minus ADC baseline in electrons
1304 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,simParam->GetNoise())
1308 signalAmp *= convert;
1309 // Add ADC baseline in mV
1310 signalAmp += baseline;
1312 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1313 // signal is larger than fADCinRange
1315 if (signalAmp >= simParam->GetADCinRange()) {
1316 adc = ((Short_t) simParam->GetADCoutRange());
1319 adc = TMath::Nint(signalAmp * adcConvert);
1322 // Saving all digits
1323 digits->SetData(row,col,time,adc);
1334 //_____________________________________________________________________________
1335 Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
1338 // Converts the sampled electron signals to s-digits
1341 AliDebug(1,Form("Start converting signals to s-digits for detector=%d",det));
1343 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1345 AliFatal("Could not get calibration object");
1353 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1354 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1355 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1357 // Get the container for the digits of this detector
1358 if (!fDigitsManager->HasSDigits()) {
1359 AliError("Digits manager has no s-digits");
1363 AliTRDarraySignal *digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1364 // Allocate memory space for the digits buffer
1365 if (!digits->HasData()) {
1366 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1369 // Create the sdigits for this chamber
1370 for (row = 0; row < nRowMax; row++ ) {
1371 for (col = 0; col < nColMax; col++ ) {
1372 for (time = 0; time < nTimeTotal; time++) {
1373 digits->SetData(row,col,time,signals->GetData(row,col,time));
1382 //_____________________________________________________________________________
1383 Bool_t AliTRDdigitizer::Digits2SDigits(AliTRDdigitsManager * const manDig
1384 , AliTRDdigitsManager * const manSDig)
1387 // Converts digits into s-digits. Needed for embedding into real data.
1390 AliDebug(1,"Start converting digits to s-digits");
1393 fGeo = new AliTRDgeometry();
1396 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1398 AliFatal("Could not get calibration object");
1402 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1404 AliFatal("Could not get simulation parameters");
1408 // Converts number of electrons to fC
1409 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1412 Double_t coupling = simParam->GetPadCoupling()
1413 * simParam->GetTimeCoupling();
1414 // Electronics conversion factor
1415 Double_t convert = kEl2fC
1416 * simParam->GetChipGain();
1417 // ADC conversion factor
1418 Double_t adcConvert = simParam->GetADCoutRange()
1419 / simParam->GetADCinRange();
1420 // The electronics baseline in mV
1421 Double_t baseline = simParam->GetADCbaseline()
1423 // The electronics baseline in electrons
1424 //Double_t baselineEl = baseline
1427 // The gainfactor calibration objects
1428 // Not used since these digits are supposed to be from real raw data
1429 //const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1430 //AliTRDCalROC *calGainFactorROC = 0;
1431 //Float_t calGainFactorDetValue = 0.0;
1437 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1439 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1440 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1441 Int_t nTimeTotal = manDig->GetDigitsParam()->GetNTimeBins(det);
1443 // Get the calibration objects
1444 //calGainFactorROC = calibration->GetGainFactorROC(det);
1445 //calGainFactorDetValue = calGainFactorDet->GetValue(det);
1448 AliTRDarrayADC *digits = (AliTRDarrayADC *) manDig->GetDigits(det);
1450 if (!manSDig->HasSDigits()) {
1451 AliError("SDigits manager has no s-digits");
1455 AliTRDarraySignal *sdigits = (AliTRDarraySignal *) manSDig->GetSDigits(det);
1456 AliTRDarrayDictionary *tracks0 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,0);
1457 AliTRDarrayDictionary *tracks1 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,1);
1458 AliTRDarrayDictionary *tracks2 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,2);
1459 // Allocate memory space for the digits buffer
1460 sdigits->Allocate(nRowMax,nColMax,nTimeTotal);
1461 tracks0->Allocate(nRowMax,nColMax,nTimeTotal);
1462 tracks1->Allocate(nRowMax,nColMax,nTimeTotal);
1463 tracks2->Allocate(nRowMax,nColMax,nTimeTotal);
1465 // Keep the digits param
1466 manSDig->GetDigitsParam()->SetNTimeBinsAll(manDig->GetDigitsParam()->GetNTimeBins(0));
1467 manSDig->GetDigitsParam()->SetADCbaselineAll(manDig->GetDigitsParam()->GetADCbaseline(0));
1469 if (digits->HasData()) {
1473 // Create the sdigits for this chamber
1474 for (row = 0; row < nRowMax; row++ ) {
1475 for (col = 0; col < nColMax; col++ ) {
1478 //Float_t padgain = calGainFactorDetValue
1479 // * calGainFactorROC->GetValue(col,row);
1481 for (time = 0; time < nTimeTotal; time++) {
1483 Short_t adcVal = digits->GetData(row,col,time);
1484 Double_t signal = (Double_t) adcVal;
1485 // ADC -> signal in mV
1486 signal /= adcConvert;
1487 // Subtract baseline in mV
1489 // Signal in mV -> signal in #electrons
1492 //signal /= padgain; // Not needed for real data
1493 // Pad and time coupling
1496 sdigits->SetData(row,col,time,signal);
1497 tracks0->SetData(row,col,time,0);
1498 tracks1->SetData(row,col,time,0);
1499 tracks2->SetData(row,col,time,0);
1508 sdigits->Compress(0);
1509 tracks0->Compress();
1510 tracks1->Compress();
1511 tracks2->Compress();
1513 // No compress just remove
1514 manDig->RemoveDigits(det);
1515 manDig->RemoveDictionaries(det);
1523 //_____________________________________________________________________________
1524 Bool_t AliTRDdigitizer::SDigits2Digits()
1527 // Merges the input s-digits and converts them to normal digits
1530 if (!MergeSDigits()) {
1534 return ConvertSDigits();
1538 //_____________________________________________________________________________
1539 Bool_t AliTRDdigitizer::MergeSDigits()
1542 // Merges the input s-digits:
1543 // - The amplitude of the different inputs are summed up.
1544 // - Of the track IDs from the input dictionaries only one is
1545 // kept for each input. This works for maximal 3 different merged inputs.
1548 // Number of track dictionary arrays
1549 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1551 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1553 AliFatal("Could not get simulation parameters");
1557 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1559 AliFatal("Could not get calibration object");
1566 AliTRDarraySignal *digitsA;
1567 AliTRDarraySignal *digitsB;
1568 AliTRDarrayDictionary *dictionaryA[kNDict];
1569 AliTRDarrayDictionary *dictionaryB[kNDict];
1571 AliTRDdigitsManager *mergeSDigitsManager = 0x0;
1572 // Get the first s-digits
1573 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1574 if (!fSDigitsManager) {
1575 AliError("No SDigits manager");
1579 // Loop through the other sets of s-digits
1580 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(fSDigitsManager);
1582 if (mergeSDigitsManager) {
1583 AliDebug(1,Form("Merge %d input files.",fSDigitsManagerList->GetSize()));
1586 AliDebug(1,"Only one input file.");
1591 while (mergeSDigitsManager) {
1595 // Loop through the detectors
1596 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1598 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet);
1599 if (mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet) != nTimeTotal) {
1600 AliError(Form("Mismatch in the number of time bins [%d,%d] in detector %d"
1602 ,mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet)
1607 Int_t nRowMax = fGeo->GetPadPlane(iDet)->GetNrows();
1608 Int_t nColMax = fGeo->GetPadPlane(iDet)->GetNcols();
1610 // Loop through the pixels of one detector and add the signals
1611 digitsA = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(iDet);
1612 digitsB = (AliTRDarraySignal *) mergeSDigitsManager->GetSDigits(iDet);
1614 if (!digitsA->HasData()) continue;
1616 if (!digitsB->HasData()) continue;
1618 for (iDict = 0; iDict < kNDict; iDict++) {
1619 dictionaryA[iDict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(iDet,iDict);
1620 dictionaryB[iDict] = (AliTRDarrayDictionary *) mergeSDigitsManager->GetDictionary(iDet,iDict);
1621 dictionaryA[iDict]->Expand();
1622 dictionaryB[iDict]->Expand();
1625 // Merge only detectors that contain a signal
1626 Bool_t doMerge = kTRUE;
1627 if (fMergeSignalOnly) {
1628 if (digitsA->GetOverThreshold(0) == 0) {
1635 AliDebug(1,Form("Merge detector %d of input no.%d",iDet,iMerge+1));
1637 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1638 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1639 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1641 // Add the amplitudes of the summable digits
1642 Float_t ampA = digitsA->GetData(iRow,iCol,iTime);
1643 Float_t ampB = digitsB->GetData(iRow,iCol,iTime);
1645 digitsA->SetData(iRow,iCol,iTime,ampA);
1647 // Add the mask to the track id if defined.
1648 for (iDict = 0; iDict < kNDict; iDict++) {
1649 Int_t trackB = dictionaryB[iDict]->GetData(iRow,iCol,iTime);
1650 if ((fMasks) && (trackB > 0)) {
1651 for (jDict = 0; jDict < kNDict; jDict++) {
1652 Int_t trackA = dictionaryA[iDict]->GetData(iRow,iCol,iTime);
1654 trackA = trackB + fMasks[iMerge];
1655 dictionaryA[iDict]->SetData(iRow,iCol,iTime,trackA);
1656 } // if: track A == 0
1658 } // if: fMasks and trackB > 0
1667 mergeSDigitsManager->RemoveDigits(iDet);
1668 mergeSDigitsManager->RemoveDictionaries(iDet);
1671 digitsA->Compress(0);
1672 for (iDict = 0; iDict < kNDict; iDict++) {
1673 dictionaryA[iDict]->Compress();
1679 // The next set of s-digits
1680 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(mergeSDigitsManager);
1682 } // while: mergeDigitsManagers
1688 //_____________________________________________________________________________
1689 Bool_t AliTRDdigitizer::ConvertSDigits()
1692 // Converts s-digits to normal digits
1695 AliTRDarraySignal *digitsIn = 0x0;
1697 if (!fSDigitsManager->HasSDigits()) {
1698 AliError("No s-digits in digits manager");
1702 // Loop through the detectors
1703 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1705 // Get the merged s-digits (signals)
1706 digitsIn = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(det);
1707 if (!digitsIn->HasData()) {
1708 AliDebug(2,Form("No digits for det=%d",det));
1712 // Convert the merged sdigits to digits
1713 if (!Signal2ADC(det,digitsIn)) {
1717 // Copy the dictionary information to the output array
1718 if (!CopyDictionary(det)) {
1723 fSDigitsManager->RemoveDigits(det);
1724 fSDigitsManager->RemoveDictionaries(det);
1726 // Run digital processing
1727 RunDigitalProcessing(det);
1729 // Compress the arrays
1730 CompressOutputArrays(det);
1732 } // for: detector numbers
1734 if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
1735 if (trklLoader->Tree())
1736 trklLoader->WriteData("OVERWRITE");
1739 // Save the values for the raw data headers
1740 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
1741 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
1744 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
1746 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
1752 //_____________________________________________________________________________
1753 Bool_t AliTRDdigitizer::CopyDictionary(Int_t det)
1756 // Copies the dictionary information from the s-digits arrays
1757 // to the output arrays
1760 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1762 AliFatal("Could not get calibration object");
1766 AliDebug(1,Form("Start copying dictionaries for detector=%d",det));
1768 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1769 AliTRDarrayDictionary *dictionaryIn[kNDict];
1770 AliTRDarrayDictionary *dictionaryOut[kNDict];
1772 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1773 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1774 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1781 for (dict = 0; dict < kNDict; dict++) {
1783 dictionaryIn[dict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(det,dict);
1784 dictionaryIn[dict]->Expand();
1785 dictionaryOut[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1786 dictionaryOut[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
1788 for (row = 0; row < nRowMax; row++) {
1789 for (col = 0; col < nColMax; col++) {
1790 for (time = 0; time < nTimeTotal; time++) {
1791 Int_t track = dictionaryIn[dict]->GetData(row,col,time);
1792 dictionaryOut[dict]->SetData(row,col,time,track);
1797 } // for: dictionaries
1803 //_____________________________________________________________________________
1804 void AliTRDdigitizer::CompressOutputArrays(Int_t det)
1807 // Compress the output arrays
1810 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1811 AliTRDarrayDictionary *dictionary = 0x0;
1816 AliTRDarrayADC *digits = 0x0;
1817 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1822 AliTRDarraySignal *digits = 0x0;
1823 digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1824 digits->Compress(0);
1827 for (Int_t dict = 0; dict < kNDict; dict++) {
1828 dictionary = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1829 dictionary->Compress();
1836 //_____________________________________________________________________________
1837 Bool_t AliTRDdigitizer::WriteDigits() const
1840 // Writes out the TRD-digits and the dictionaries
1844 fRunLoader->CdGAFile();
1846 // Store the digits and the dictionary in the tree
1847 return fDigitsManager->WriteDigits();
1851 //_____________________________________________________________________________
1852 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1855 // Initializes the output branches
1861 AliError("Run Loader is NULL");
1865 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
1867 AliError("Can not get TRD loader from Run Loader");
1874 // If we produce SDigits
1875 tree = loader->TreeS();
1877 loader->MakeTree("S");
1878 tree = loader->TreeS();
1882 // If we produce Digits
1883 tree = loader->TreeD();
1885 loader->MakeTree("D");
1886 tree = loader->TreeD();
1889 fDigitsManager->SetEvent(iEvent);
1890 fDigitsManager->MakeBranch(tree);
1894 //_____________________________________________________________________________
1895 Int_t AliTRDdigitizer::Diffusion(Float_t vdrift, Double_t absdriftlength
1897 , Double_t &lRow, Double_t &lCol, Double_t &lTime)
1900 // Applies the diffusion smearing to the position of a single electron.
1901 // Depends on absolute drift length.
1904 Float_t diffL = 0.0;
1905 Float_t diffT = 0.0;
1907 if (AliTRDCommonParam::Instance()->GetDiffCoeff(diffL,diffT,vdrift)) {
1909 Float_t driftSqrt = TMath::Sqrt(absdriftlength);
1910 Float_t sigmaT = driftSqrt * diffT;
1911 Float_t sigmaL = driftSqrt * diffL;
1912 lRow = gRandom->Gaus(lRow ,sigmaT);
1913 if (AliTRDCommonParam::Instance()->ExBOn()) {
1914 lCol = gRandom->Gaus(lCol ,sigmaT * 1.0 / (1.0 + exbvalue*exbvalue));
1915 lTime = gRandom->Gaus(lTime,sigmaL * 1.0 / (1.0 + exbvalue*exbvalue));
1918 lCol = gRandom->Gaus(lCol ,sigmaT);
1919 lTime = gRandom->Gaus(lTime,sigmaL);
1933 //_____________________________________________________________________________
1934 void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
1937 // Run the digital processing in the TRAP
1940 AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();
1942 AliTRDarrayADC *digits = fDigitsManager->GetDigits(det);
1946 //Call the methods in the mcm class using the temporary array as input
1947 // process the data in the same order as in hardware
1948 for (Int_t side = 0; side <= 1; side++) {
1949 for(Int_t rob = side; rob < digits->GetNrow() / 2; rob += 2) {
1950 for(Int_t mcm = 0; mcm < 16; mcm++) {
1951 fMcmSim->Init(det, rob, mcm);
1952 fMcmSim->SetDataByPad(digits, fDigitsManager);
1954 if (feeParam->GetTracklet()) {
1955 fMcmSim->Tracklet();
1956 fMcmSim->StoreTracklets();
1958 fMcmSim->ZSMapping();
1959 fMcmSim->WriteData(digits);