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 "AliRunDigitizer.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"
74 ClassImp(AliTRDdigitizer)
76 //_____________________________________________________________________________
77 AliTRDdigitizer::AliTRDdigitizer()
82 ,fSDigitsManagerList(0)
85 ,fMcmSim(new AliTRDmcmSim)
90 ,fMergeSignalOnly(kFALSE)
93 // AliTRDdigitizer default constructor
98 //_____________________________________________________________________________
99 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
100 :AliDigitizer(name,title)
104 ,fSDigitsManagerList(0)
107 ,fMcmSim(new AliTRDmcmSim)
112 ,fMergeSignalOnly(kFALSE)
115 // AliTRDdigitizer constructor
120 //_____________________________________________________________________________
121 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
122 , const Text_t *name, const Text_t *title)
123 :AliDigitizer(manager,name,title)
127 ,fSDigitsManagerList(0)
130 ,fMcmSim(new AliTRDmcmSim)
135 ,fMergeSignalOnly(kFALSE)
138 // AliTRDdigitizer constructor
143 //_____________________________________________________________________________
144 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
145 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
149 ,fSDigitsManagerList(0)
152 ,fMcmSim(new AliTRDmcmSim)
157 ,fMergeSignalOnly(kFALSE)
160 // AliTRDdigitizer constructor
165 //_____________________________________________________________________________
166 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
171 ,fSDigitsManagerList(0)
174 ,fMcmSim(new AliTRDmcmSim)
177 ,fCompress(d.fCompress)
178 ,fSDigits(d.fSDigits)
179 ,fMergeSignalOnly(d.fMergeSignalOnly)
182 // AliTRDdigitizer copy constructor
187 //_____________________________________________________________________________
188 AliTRDdigitizer::~AliTRDdigitizer()
191 // AliTRDdigitizer destructor
194 if (fDigitsManager) {
195 delete fDigitsManager;
199 if (fSDigitsManager) {
200 // s-digitsmanager will be deleted via list
204 if (fSDigitsManagerList) {
205 fSDigitsManagerList->Delete();
206 delete fSDigitsManagerList;
207 fSDigitsManagerList = 0;
227 //_____________________________________________________________________________
228 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
231 // Assignment operator
235 ((AliTRDdigitizer &) d).Copy(*this);
242 //_____________________________________________________________________________
243 void AliTRDdigitizer::Copy(TObject &d) const
249 ((AliTRDdigitizer &) d).fRunLoader = 0;
250 ((AliTRDdigitizer &) d).fDigitsManager = 0;
251 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
252 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
253 ((AliTRDdigitizer &) d).fTRD = 0;
254 ((AliTRDdigitizer &) d).fGeo = 0;
255 ((AliTRDdigitizer &) d).fEvent = 0;
256 ((AliTRDdigitizer &) d).fMasks = 0;
257 ((AliTRDdigitizer &) d).fCompress = fCompress;
258 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
259 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
263 //_____________________________________________________________________________
264 void AliTRDdigitizer::Exec(const Option_t * const option)
267 // Executes the merging
272 AliTRDdigitsManager *sdigitsManager;
274 TString optionString = option;
275 if (optionString.Contains("deb")) {
276 AliLog::SetClassDebugLevel("AliTRDdigitizer",1);
277 AliInfo("Called with debug option");
280 // The AliRoot file is already connected by the manager
281 AliRunLoader *inrl = 0x0;
284 AliDebug(1,"AliRun object found on file.");
287 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
289 gAlice = inrl->GetAliRun();
291 AliError("Could not find AliRun object.");
296 Int_t nInput = fManager->GetNinputs();
297 fMasks = new Int_t[nInput];
298 for (iInput = 0; iInput < nInput; iInput++) {
299 fMasks[iInput] = fManager->GetMask(iInput);
306 AliRunLoader *orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
308 if (InitDetector()) {
310 AliLoader *ogime = orl->GetLoader("TRDLoader");
314 // If we produce SDigits
315 tree = ogime->TreeS();
317 ogime->MakeTree("S");
318 tree = ogime->TreeS();
322 // If we produce Digits
323 tree = ogime->TreeD();
325 ogime->MakeTree("D");
326 tree = ogime->TreeD();
334 for (iInput = 0; iInput < nInput; iInput++) {
336 AliDebug(1,Form("Add input stream %d",iInput));
338 // Check if the input tree exists
339 inrl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
340 AliLoader *gime = inrl->GetLoader("TRDLoader");
342 TTree *treees = gime->TreeS();
344 if (gime->LoadSDigits()) {
345 AliError(Form("Error Occured while loading S. Digits for input %d.",iInput));
348 treees = gime->TreeS();
352 AliError(Form("Input stream %d does not exist",iInput));
356 // Read the s-digits via digits manager
357 sdigitsManager = new AliTRDdigitsManager();
358 sdigitsManager->SetSDigits(kTRUE);
360 AliRunLoader *rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(iInput));
361 AliLoader *gimme = rl->GetLoader("TRDLoader");
364 gimme->LoadSDigits();
367 sdigitsManager->ReadDigits(gimme->TreeS());
369 // Add the s-digits to the input list
370 AddSDigitsManager(sdigitsManager);
374 // Convert the s-digits to normal digits
375 AliDebug(1,"Do the conversion");
379 AliDebug(1,"Write the digits");
380 fDigitsManager->WriteDigits();
386 DeleteSDigitsManager();
392 //_____________________________________________________________________________
393 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
396 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
398 // Connect the AliRoot file containing Geometry, Kine, and Hits
401 TString evfoldname = AliConfig::GetDefaultEventFolderName();
403 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
405 fRunLoader = AliRunLoader::Open(file,evfoldname,"UPDATE");
408 AliError(Form("Can not open session for file %s.",file));
412 if (!fRunLoader->GetAliRun()) {
413 fRunLoader->LoadgAlice();
415 gAlice = fRunLoader->GetAliRun();
418 AliDebug(1,"AliRun object found on file.");
421 AliError("Could not find AliRun object.");
427 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
429 AliError("Can not get TRD loader from Run Loader");
433 if (InitDetector()) {
436 // If we produce SDigits
437 tree = loader->TreeS();
439 loader->MakeTree("S");
440 tree = loader->TreeS();
444 // If we produce Digits
445 tree = loader->TreeD();
447 loader->MakeTree("D");
448 tree = loader->TreeD();
451 return MakeBranch(tree);
459 //_____________________________________________________________________________
460 Bool_t AliTRDdigitizer::Open(AliRunLoader * const runLoader, Int_t nEvent)
463 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
465 // Connect the AliRoot file containing Geometry, Kine, and Hits
468 fRunLoader = runLoader;
470 AliError("RunLoader does not exist");
474 if (!fRunLoader->GetAliRun()) {
475 fRunLoader->LoadgAlice();
477 gAlice = fRunLoader->GetAliRun();
480 AliDebug(1,"AliRun object found on file.");
483 AliError("Could not find AliRun object.");
489 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
491 AliError("Can not get TRD loader from Run Loader");
495 if (InitDetector()) {
498 // If we produce SDigits
499 tree = loader->TreeS();
501 loader->MakeTree("S");
502 tree = loader->TreeS();
506 // If we produce Digits
507 tree = loader->TreeD();
509 loader->MakeTree("D");
510 tree = loader->TreeD();
513 return MakeBranch(tree);
521 //_____________________________________________________________________________
522 Bool_t AliTRDdigitizer::InitDetector()
525 // Sets the pointer to the TRD detector and the geometry
528 // Get the pointer to the detector class and check for version 1
529 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
531 AliFatal("No TRD module found");
534 if (fTRD->IsVersion() != 1) {
535 AliFatal("TRD must be version 1 (slow simulator)");
540 fGeo = new AliTRDgeometry();
542 // Create a digits manager
543 if (fDigitsManager) {
544 delete fDigitsManager;
546 fDigitsManager = new AliTRDdigitsManager();
547 fDigitsManager->SetSDigits(fSDigits);
548 fDigitsManager->CreateArrays();
549 fDigitsManager->SetEvent(fEvent);
551 // The list for the input s-digits manager to be merged
552 if (fSDigitsManagerList) {
553 fSDigitsManagerList->Delete();
556 fSDigitsManagerList = new TList();
563 //_____________________________________________________________________________
564 Bool_t AliTRDdigitizer::MakeBranch(TTree *tree) const
567 // Create the branches for the digits array
570 return fDigitsManager->MakeBranch(tree);
574 //_____________________________________________________________________________
575 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
578 // Add a digits manager for s-digits to the input list.
581 fSDigitsManagerList->Add(man);
585 //_____________________________________________________________________________
586 void AliTRDdigitizer::DeleteSDigitsManager()
589 // Removes digits manager from the input list.
592 fSDigitsManagerList->Delete();
596 //_____________________________________________________________________________
597 Bool_t AliTRDdigitizer::MakeDigits()
603 AliDebug(1,"Start creating digits");
606 AliError("No geometry defined");
610 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
612 AliFatal("Could not get calibration object");
616 const Int_t kNdet = AliTRDgeometry::Ndet();
618 Float_t **hits = new Float_t*[kNdet];
619 Int_t *nhit = new Int_t[kNdet];
621 AliTRDarraySignal *signals = 0x0;
623 // Check the number of time bins from simParam against OCDB,
624 // if OCDB value is not supposed to be used.
625 // As default, the value from OCDB is taken
626 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
627 if (calibration->GetNumberOfTimeBinsDCS() != AliTRDSimParam::Instance()->GetNTimeBins()) {
628 AliWarning(Form("Number of time bins is different to OCDB value [SIM=%d, OCDB=%d]"
629 ,AliTRDSimParam::Instance()->GetNTimeBins()
630 ,calibration->GetNumberOfTimeBinsDCS()));
632 // Save the values for the raw data headers
633 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
636 // Save the values for the raw data headers
637 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(calibration->GetNumberOfTimeBinsDCS());
640 // Save the values for the raw data headers
641 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
643 // Sort all hits according to detector number
644 if (!SortHits(hits,nhit)) {
645 AliError("Sorting hits failed");
651 // Loop through all detectors
652 for (Int_t det = 0; det < kNdet; det++) {
654 // Detectors that are switched off, not installed, etc.
655 if (( calibration->IsChamberInstalled(det)) &&
656 (!calibration->IsChamberMasked(det)) &&
657 ( fGeo->ChamberInGeometry(det)) &&
660 signals = new AliTRDarraySignal();
662 // Convert the hits of the current detector to detector signals
663 if (!ConvertHits(det,hits[det],nhit[det],signals)) {
664 AliError(Form("Conversion of hits failed for detector=%d",det));
671 // Convert the detector signals to digits or s-digits
672 if (!ConvertSignals(det,signals)) {
673 AliError(Form("Conversion of signals failed for detector=%d",det));
681 // Delete the signals array
685 } // if: detector status
692 if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
693 if (trklLoader->Tree())
694 trklLoader->WriteData("OVERWRITE");
705 //_____________________________________________________________________________
706 Bool_t AliTRDdigitizer::SortHits(Float_t **hits, Int_t *nhit)
709 // Read all the hits and sorts them according to detector number
710 // in the output array <hits>.
713 AliDebug(1,"Start sorting hits");
715 const Int_t kNdet = AliTRDgeometry::Ndet();
716 // Size of the hit vector
717 const Int_t kNhit = 6;
722 Int_t *lhit = new Int_t[kNdet];
724 for (Int_t det = 0; det < kNdet; det++) {
730 AliLoader *gimme = fRunLoader->GetLoader("TRDLoader");
731 if (!gimme->TreeH()) {
734 TTree *hitTree = gimme->TreeH();
735 if (hitTree == 0x0) {
736 AliError("Can not get TreeH");
740 fTRD->SetTreeAddress();
742 // Get the number of entries in the hit tree
743 // (Number of primary particles creating a hit somewhere)
744 Int_t nTrk = (Int_t) hitTree->GetEntries();
745 AliDebug(1,Form("Found %d tracks",nTrk));
747 // Loop through all the tracks in the tree
748 for (Int_t iTrk = 0; iTrk < nTrk; iTrk++) {
750 gAlice->GetMCApp()->ResetHits();
751 hitTree->GetEvent(iTrk);
754 AliError(Form("No hits array for track = %d",iTrk));
758 // Number of hits for this track
759 nhitTrk = fTRD->Hits()->GetEntriesFast();
762 // Loop through the TRD hits
763 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
768 // Don't analyze test hits
769 if (((Int_t) hit->GetCharge()) != 0) {
771 Int_t trk = hit->Track();
772 Int_t det = hit->GetDetector();
773 Int_t q = hit->GetCharge();
774 Float_t x = hit->X();
775 Float_t y = hit->Y();
776 Float_t z = hit->Z();
777 Float_t time = hit->GetTime();
779 if (nhit[det] == lhit[det]) {
780 // Inititialization of new detector
781 xyz = new Float_t[kNhit*(nhitTrk+lhit[det])];
783 memcpy(xyz,hits[det],sizeof(Float_t)*kNhit*lhit[det]);
786 lhit[det] += nhitTrk;
792 xyz[nhit[det]*kNhit+0] = x;
793 xyz[nhit[det]*kNhit+1] = y;
794 xyz[nhit[det]*kNhit+2] = z;
795 xyz[nhit[det]*kNhit+3] = q;
796 xyz[nhit[det]*kNhit+4] = trk;
797 xyz[nhit[det]*kNhit+5] = time;
802 hit = (AliTRDhit *) fTRD->NextHit();
804 } // for: hits of one track
814 //_____________________________________________________________________________
815 Bool_t AliTRDdigitizer::ConvertHits(Int_t det
816 , const Float_t * const hits
818 , AliTRDarraySignal *signals)
821 // Converts the detectorwise sorted hits to detector signals
824 AliDebug(1,Form("Start converting hits for detector=%d (nhits=%d)",det,nhit));
826 // Number of pads included in the pad response
827 const Int_t kNpad = 3;
828 // Number of track dictionary arrays
829 const Int_t kNdict = AliTRDdigitsManager::kNDict;
830 // Size of the hit vector
831 const Int_t kNhit = 6;
833 // Width of the amplification region
834 const Float_t kAmWidth = AliTRDgeometry::AmThick();
835 // Width of the drift region
836 const Float_t kDrWidth = AliTRDgeometry::DrThick();
837 // Drift + amplification region
838 const Float_t kDrMin = - 0.5 * kAmWidth;
839 const Float_t kDrMax = kDrWidth + 0.5 * kAmWidth;
843 Int_t timeBinTRFend = 1;
847 Double_t padSignal[kNpad];
848 Double_t signalOld[kNpad];
850 AliTRDarrayDictionary *dictionary[kNdict];
852 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
853 AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
854 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
857 AliFatal("Could not get common parameterss");
861 AliFatal("Could not get simulation parameters");
865 AliFatal("Could not get calibration object");
869 // Get the detector wise calibration objects
870 AliTRDCalROC *calVdriftROC = 0;
871 Float_t calVdriftDetValue = 0.0;
872 const AliTRDCalDet *calVdriftDet = calibration->GetVdriftDet();
873 AliTRDCalROC *calT0ROC = 0;
874 Float_t calT0DetValue = 0.0;
875 const AliTRDCalDet *calT0Det = calibration->GetT0Det();
877 if (simParam->TRFOn()) {
878 timeBinTRFend = ((Int_t) (simParam->GetTRFhi()
879 * commonParam->GetSamplingFrequency())) - 1;
882 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
883 Float_t samplingRate = commonParam->GetSamplingFrequency();
884 Float_t elAttachProp = simParam->GetElAttachProp() / 100.0;
886 AliTRDpadPlane *padPlane = fGeo->GetPadPlane(det);
887 Int_t layer = fGeo->GetLayer(det); //update
888 Float_t row0 = padPlane->GetRow0ROC();
889 Int_t nRowMax = padPlane->GetNrows();
890 Int_t nColMax = padPlane->GetNcols();
892 // Create a new array for the signals
893 signals->Allocate(nRowMax,nColMax,nTimeTotal);
895 // Create a new array for the dictionary
896 for (dict = 0; dict < kNdict; dict++) {
897 dictionary[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
898 dictionary[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
901 // Loop through the hits in this detector
902 for (Int_t hit = 0; hit < nhit; hit++) {
904 pos[0] = hits[hit*kNhit+0];
905 pos[1] = hits[hit*kNhit+1];
906 pos[2] = hits[hit*kNhit+2];
907 Float_t q = hits[hit*kNhit+3];
908 Float_t hittime = hits[hit*kNhit+5];
909 Int_t track = ((Int_t) hits[hit*kNhit+4]);
913 // Find the current volume with the geo manager
914 gGeoManager->SetCurrentPoint(pos);
915 gGeoManager->FindNode();
916 if (strstr(gGeoManager->GetPath(),"/UK")) {
920 // Get the calibration objects
921 calVdriftROC = calibration->GetVdriftROC(det);
922 calVdriftDetValue = calVdriftDet->GetValue(det);
923 calT0ROC = calibration->GetT0ROC(det);
924 calT0DetValue = calT0Det->GetValue(det);
926 // Go to the local coordinate system:
927 // loc[0] - col direction in amplification or driftvolume
928 // loc[1] - row direction in amplification or driftvolume
929 // loc[2] - time direction in amplification or driftvolume
930 gGeoManager->MasterToLocal(pos,loc);
932 // Relative to middle of amplification region
933 loc[2] = loc[2] - kDrWidth/2.0 - kAmWidth/2.0;
936 // The driftlength [cm] (w/o diffusion yet !).
937 // It is negative if the hit is between pad plane and anode wires.
938 Double_t driftlength = -1.0 * loc[2];
940 // Stupid patch to take care of TR photons that are absorbed
941 // outside the chamber volume. A real fix would actually need
942 // a more clever implementation of the TR hit generation
944 if ((loc[1] < padPlane->GetRowEndROC()) ||
945 (loc[1] > padPlane->GetRow0ROC())) {
948 if ((driftlength < kDrMin) ||
949 (driftlength > kDrMax)) {
954 // Get row and col of unsmeared electron to retrieve drift velocity
955 // The pad row (z-direction)
956 Int_t rowE = padPlane->GetPadRowNumberROC(loc[1]);
960 Double_t rowOffset = padPlane->GetPadRowOffsetROC(rowE,loc[1]);
961 // The pad column (rphi-direction)
962 Double_t offsetTilt = padPlane->GetTiltOffset(rowOffset);
963 Int_t colE = padPlane->GetPadColNumber(loc[0]+offsetTilt);
967 Double_t colOffset = 0.0;
969 // Normalized drift length
970 Float_t driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
971 Double_t absdriftlength = TMath::Abs(driftlength);
972 if (commonParam->ExBOn()) {
973 absdriftlength /= TMath::Sqrt(GetLorentzFactor(driftvelocity));
976 // Loop over all electrons of this hit
977 // TR photons produce hits with negative charge
978 Int_t nEl = ((Int_t) TMath::Abs(q));
979 for (Int_t iEl = 0; iEl < nEl; iEl++) {
981 // Now the real local coordinate system of the ROC
982 // column direction: locC
983 // row direction: locR
984 // time direction: locT
985 // locR and locC are identical to the coordinates of the corresponding
986 // volumina of the drift or amplification region.
987 // locT is defined relative to the wire plane (i.e. middle of amplification
988 // region), meaning locT = 0, and is negative for hits coming from the
990 Double_t locC = loc[0];
991 Double_t locR = loc[1];
992 Double_t locT = loc[2];
994 // Electron attachment
995 if (simParam->ElAttachOn()) {
996 if (gRandom->Rndm() < (absdriftlength * elAttachProp)) {
1001 // Apply the diffusion smearing
1002 if (simParam->DiffusionOn()) {
1003 if (!(Diffusion(driftvelocity,absdriftlength,locR,locC,locT))) {
1008 // Apply E x B effects (depends on drift direction)
1009 if (commonParam->ExBOn()) {
1010 if (!(ExB(driftvelocity,driftlength,locC))) {
1015 // The electron position after diffusion and ExB in pad coordinates.
1016 // The pad row (z-direction)
1017 rowE = padPlane->GetPadRowNumberROC(locR);
1018 if (rowE < 0) continue;
1019 rowOffset = padPlane->GetPadRowOffsetROC(rowE,locR);
1021 // The pad column (rphi-direction)
1022 offsetTilt = padPlane->GetTiltOffset(rowOffset);
1023 colE = padPlane->GetPadColNumber(locC+offsetTilt);
1024 if (colE < 0) continue;
1025 colOffset = padPlane->GetPadColOffset(colE,locC+offsetTilt);
1027 // Also re-retrieve drift velocity because col and row may have changed
1028 driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
1029 Float_t t0 = calT0DetValue + calT0ROC->GetValue(colE,rowE);
1031 // Convert the position to drift time [mus], using either constant drift velocity or
1032 // time structure of drift cells (non-isochronity, GARFIELD calculation).
1033 // Also add absolute time of hits to take pile-up events into account properly
1035 if (simParam->TimeStructOn()) {
1036 // Get z-position with respect to anode wire
1037 Double_t zz = row0 - locR + padPlane->GetAnodeWireOffset();
1038 zz -= ((Int_t)(2 * zz)) / 2.0;
1042 // Use drift time map (GARFIELD)
1043 drifttime = commonParam->TimeStruct(driftvelocity,0.5*kAmWidth-1.0*locT,zz)
1047 // Use constant drift velocity
1048 drifttime = TMath::Abs(locT) / driftvelocity
1052 // Apply the gas gain including fluctuations
1053 Double_t ggRndm = 0.0;
1055 ggRndm = gRandom->Rndm();
1056 } while (ggRndm <= 0);
1057 Double_t signal = -(simParam->GetGasGain()) * TMath::Log(ggRndm);
1059 // Apply the pad response
1060 if (simParam->PRFOn()) {
1061 // The distance of the electron to the center of the pad
1062 // in units of pad width
1063 Double_t dist = (colOffset - 0.5*padPlane->GetColSize(colE))
1064 / padPlane->GetColSize(colE);
1065 // This is a fixed parametrization, i.e. not dependent on
1066 // calibration values !
1067 if (!(calibration->PadResponse(signal,dist,layer,padSignal))) continue;
1071 padSignal[1] = signal;
1075 // The time bin (always positive), with t0 distortion
1076 Double_t timeBinIdeal = drifttime * samplingRate + t0;
1078 if (TMath::Abs(timeBinIdeal) > 2*nTimeTotal) {
1079 timeBinIdeal = 2 * nTimeTotal;
1081 Int_t timeBinTruncated = ((Int_t) timeBinIdeal);
1082 // The distance of the position to the middle of the timebin
1083 Double_t timeOffset = ((Float_t) timeBinTruncated
1084 + 0.5 - timeBinIdeal) / samplingRate;
1086 // Sample the time response inside the drift region
1087 // + additional time bins before and after.
1088 // The sampling is done always in the middle of the time bin
1089 for (Int_t iTimeBin = TMath::Max(timeBinTruncated,0)
1090 ;iTimeBin < TMath::Min(timeBinTruncated+timeBinTRFend,nTimeTotal)
1093 // Apply the time response
1094 Double_t timeResponse = 1.0;
1095 Double_t crossTalk = 0.0;
1096 Double_t time = (iTimeBin - timeBinTruncated) / samplingRate + timeOffset;
1098 if (simParam->TRFOn()) {
1099 timeResponse = simParam->TimeResponse(time);
1101 if (simParam->CTOn()) {
1102 crossTalk = simParam->CrossTalk(time);
1109 for (iPad = 0; iPad < kNpad; iPad++) {
1111 Int_t colPos = colE + iPad - 1;
1112 if (colPos < 0) continue;
1113 if (colPos >= nColMax) break;
1116 signalOld[iPad] = signals->GetData(rowE,colPos,iTimeBin);
1118 if (colPos != colE) {
1119 // Cross talk added to non-central pads
1120 signalOld[iPad] += padSignal[iPad]
1121 * (timeResponse + crossTalk);
1124 // W/o cross talk at central pad
1125 signalOld[iPad] += padSignal[iPad]
1129 signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);
1131 // Store the track index in the dictionary
1132 // Note: We store index+1 in order to allow the array to be compressed
1133 // Note2: Taking out the +1 in track
1134 if (signalOld[iPad] > 0.0) {
1135 for (dict = 0; dict < kNdict; dict++) {
1136 Int_t oldTrack = dictionary[dict]->GetData(rowE,colPos,iTimeBin);
1137 if (oldTrack == track) break;
1138 if (oldTrack == -1 ) {
1139 dictionary[dict]->SetData(rowE,colPos,iTimeBin,track);
1147 } // Loop: time bins
1149 } // Loop: electrons of a single hit
1153 AliDebug(2,Form("Finished analyzing %d hits",nhit));
1159 //_____________________________________________________________________________
1160 Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
1163 // Convert signals to digits
1166 AliDebug(1,Form("Start converting the signals for detector %d",det));
1169 // Convert the signal array to s-digits
1170 if (!Signal2SDigits(det,signals)) {
1175 // Convert the signal array to digits
1176 if (!Signal2ADC(det,signals)) {
1179 // Run digital processing for digits
1180 RunDigitalProcessing(det);
1183 // Compress the arrays
1184 CompressOutputArrays(det);
1190 //_____________________________________________________________________________
1191 Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
1194 // Converts the sampled electron signals to ADC values for a given chamber
1197 AliDebug(1,Form("Start converting signals to ADC values for detector=%d",det));
1199 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1201 AliFatal("Could not get calibration object");
1205 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1207 AliFatal("Could not get simulation parameters");
1211 // Converts number of electrons to fC
1212 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1215 Double_t coupling = simParam->GetPadCoupling()
1216 * simParam->GetTimeCoupling();
1217 // Electronics conversion factor
1218 Double_t convert = kEl2fC
1219 * simParam->GetChipGain();
1220 // ADC conversion factor
1221 Double_t adcConvert = simParam->GetADCoutRange()
1222 / simParam->GetADCinRange();
1223 // The electronics baseline in mV
1224 Double_t baseline = simParam->GetADCbaseline()
1226 // The electronics baseline in electrons
1227 Double_t baselineEl = baseline
1234 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1235 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1236 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1237 if (fSDigitsManager->GetDigitsParam()->GetNTimeBins(det)) {
1238 nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1241 AliFatal("Could not get number of time bins");
1245 // The gainfactor calibration objects
1246 const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1247 AliTRDCalROC *calGainFactorROC = 0;
1248 Float_t calGainFactorDetValue = 0.0;
1250 AliTRDarrayADC *digits = 0x0;
1253 AliError(Form("Signals array for detector %d does not exist\n",det));
1256 if (signals->HasData()) {
1257 // Expand the container if neccessary
1261 // Create missing containers
1262 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1265 // Get the container for the digits of this detector
1266 if (fDigitsManager->HasSDigits()) {
1267 AliError("Digits manager has s-digits");
1271 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1272 // Allocate memory space for the digits buffer
1273 if (!digits->HasData()) {
1274 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1277 // Get the calibration objects
1278 calGainFactorROC = calibration->GetGainFactorROC(det);
1279 calGainFactorDetValue = calGainFactorDet->GetValue(det);
1281 // Create the digits for this chamber
1282 for (row = 0; row < nRowMax; row++ ) {
1283 for (col = 0; col < nColMax; col++ ) {
1285 // Check whether pad is masked
1286 // Bridged pads are not considered yet!!!
1287 if (calibration->IsPadMasked(det,col,row) || calibration->IsPadNotConnected(det,col,row)) {
1292 Float_t padgain = calGainFactorDetValue
1293 * calGainFactorROC->GetValue(col,row);
1295 AliError(Form("Not a valid gain %f, %d %d %d",padgain,det,col,row));
1298 for (time = 0; time < nTimeTotal; time++) {
1300 // Get the signal amplitude
1301 Float_t signalAmp = signals->GetData(row,col,time);
1302 // Pad and time coupling
1303 signalAmp *= coupling;
1305 signalAmp *= padgain;
1307 // Add the noise, starting from minus ADC baseline in electrons
1308 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,simParam->GetNoise())
1312 signalAmp *= convert;
1313 // Add ADC baseline in mV
1314 signalAmp += baseline;
1316 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1317 // signal is larger than fADCinRange
1319 if (signalAmp >= simParam->GetADCinRange()) {
1320 adc = ((Short_t) simParam->GetADCoutRange());
1323 adc = TMath::Nint(signalAmp * adcConvert);
1326 // Saving all digits
1327 digits->SetData(row,col,time,adc);
1338 //_____________________________________________________________________________
1339 Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
1342 // Converts the sampled electron signals to s-digits
1345 AliDebug(1,Form("Start converting signals to s-digits for detector=%d",det));
1347 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1349 AliFatal("Could not get calibration object");
1357 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1358 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1359 Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1361 // Get the container for the digits of this detector
1362 if (!fDigitsManager->HasSDigits()) {
1363 AliError("Digits manager has no s-digits");
1367 AliTRDarraySignal *digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1368 // Allocate memory space for the digits buffer
1369 if (!digits->HasData()) {
1370 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1373 // Create the sdigits for this chamber
1374 for (row = 0; row < nRowMax; row++ ) {
1375 for (col = 0; col < nColMax; col++ ) {
1376 for (time = 0; time < nTimeTotal; time++) {
1377 digits->SetData(row,col,time,signals->GetData(row,col,time));
1386 //_____________________________________________________________________________
1387 Bool_t AliTRDdigitizer::Digits2SDigits(AliTRDdigitsManager * const manDig
1388 , AliTRDdigitsManager * const manSDig)
1391 // Converts digits into s-digits. Needed for embedding into real data.
1394 AliDebug(1,"Start converting digits to s-digits");
1397 fGeo = new AliTRDgeometry();
1400 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1402 AliFatal("Could not get calibration object");
1406 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1408 AliFatal("Could not get simulation parameters");
1412 // Converts number of electrons to fC
1413 const Double_t kEl2fC = 1.602e-19 * 1.0e15;
1416 Double_t coupling = simParam->GetPadCoupling()
1417 * simParam->GetTimeCoupling();
1418 // Electronics conversion factor
1419 Double_t convert = kEl2fC
1420 * simParam->GetChipGain();
1421 // ADC conversion factor
1422 Double_t adcConvert = simParam->GetADCoutRange()
1423 / simParam->GetADCinRange();
1424 // The electronics baseline in mV
1425 Double_t baseline = simParam->GetADCbaseline()
1427 // The electronics baseline in electrons
1428 //Double_t baselineEl = baseline
1431 // The gainfactor calibration objects
1432 //const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
1433 //AliTRDCalROC *calGainFactorROC = 0;
1434 //Float_t calGainFactorDetValue = 0.0;
1440 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1442 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1443 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1444 Int_t nTimeTotal = manDig->GetDigitsParam()->GetNTimeBins(det);
1446 // Get the calibration objects
1447 //calGainFactorROC = calibration->GetGainFactorROC(det);
1448 //calGainFactorDetValue = calGainFactorDet->GetValue(det);
1451 AliTRDarrayADC *digits = (AliTRDarrayADC *) manDig->GetDigits(det);
1453 if (!manSDig->HasSDigits()) {
1454 AliError("SDigits manager has no s-digits");
1458 AliTRDarraySignal *sdigits = (AliTRDarraySignal *) manSDig->GetSDigits(det);
1459 AliTRDarrayDictionary *tracks0 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,0);
1460 AliTRDarrayDictionary *tracks1 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,1);
1461 AliTRDarrayDictionary *tracks2 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,2);
1462 // Allocate memory space for the digits buffer
1463 sdigits->Allocate(nRowMax,nColMax,nTimeTotal);
1464 tracks0->Allocate(nRowMax,nColMax,nTimeTotal);
1465 tracks1->Allocate(nRowMax,nColMax,nTimeTotal);
1466 tracks2->Allocate(nRowMax,nColMax,nTimeTotal);
1468 // Keep the digits param
1469 manSDig->GetDigitsParam()->SetNTimeBinsAll(manDig->GetDigitsParam()->GetNTimeBins(0));
1470 manSDig->GetDigitsParam()->SetADCbaselineAll(manDig->GetDigitsParam()->GetADCbaseline(0));
1472 if (digits->HasData()) {
1476 // Create the sdigits for this chamber
1477 for (row = 0; row < nRowMax; row++ ) {
1478 for (col = 0; col < nColMax; col++ ) {
1481 //Float_t padgain = calGainFactorDetValue
1482 // * calGainFactorROC->GetValue(col,row);
1484 for (time = 0; time < nTimeTotal; time++) {
1486 Short_t adcVal = digits->GetData(row,col,time);
1487 Double_t signal = (Double_t) adcVal;
1488 // ADC -> signal in mV
1489 signal /= adcConvert;
1490 // Subtract baseline in mV
1492 // Signal in mV -> signal in #electrons
1495 //signal /= padgain; // Not needed for real data
1496 // Pad and time coupling
1499 sdigits->SetData(row,col,time,signal);
1500 tracks0->SetData(row,col,time,0);
1501 tracks1->SetData(row,col,time,0);
1502 tracks2->SetData(row,col,time,0);
1511 sdigits->Compress(0);
1512 tracks0->Compress();
1513 tracks1->Compress();
1514 tracks2->Compress();
1516 // No compress just remove
1517 manDig->RemoveDigits(det);
1518 manDig->RemoveDictionaries(det);
1526 //_____________________________________________________________________________
1527 Bool_t AliTRDdigitizer::SDigits2Digits()
1530 // Merges the input s-digits and converts them to normal digits
1533 if (!MergeSDigits()) {
1537 return ConvertSDigits();
1541 //_____________________________________________________________________________
1542 Bool_t AliTRDdigitizer::MergeSDigits()
1545 // Merges the input s-digits:
1546 // - The amplitude of the different inputs are summed up.
1547 // - Of the track IDs from the input dictionaries only one is
1548 // kept for each input. This works for maximal 3 different merged inputs.
1551 // Number of track dictionary arrays
1552 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1554 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
1556 AliFatal("Could not get simulation parameters");
1560 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1562 AliFatal("Could not get calibration object");
1569 AliTRDarraySignal *digitsA;
1570 AliTRDarraySignal *digitsB;
1571 AliTRDarrayDictionary *dictionaryA[kNDict];
1572 AliTRDarrayDictionary *dictionaryB[kNDict];
1574 AliTRDdigitsManager *mergeSDigitsManager = 0x0;
1575 // Get the first s-digits
1576 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1577 if (!fSDigitsManager) {
1578 AliError("No SDigits manager");
1582 // Loop through the other sets of s-digits
1583 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(fSDigitsManager);
1585 if (mergeSDigitsManager) {
1586 AliDebug(1,Form("Merge %d input files.",fSDigitsManagerList->GetSize()));
1589 AliDebug(1,"Only one input file.");
1594 while (mergeSDigitsManager) {
1598 // Loop through the detectors
1599 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1601 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet);
1602 if (mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet) != nTimeTotal) {
1603 AliError(Form("Mismatch in the number of time bins [%d,%d] in detector %d"
1605 ,mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet)
1610 Int_t nRowMax = fGeo->GetPadPlane(iDet)->GetNrows();
1611 Int_t nColMax = fGeo->GetPadPlane(iDet)->GetNcols();
1613 // Loop through the pixels of one detector and add the signals
1614 digitsA = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(iDet);
1615 digitsB = (AliTRDarraySignal *) mergeSDigitsManager->GetSDigits(iDet);
1617 if (!digitsA->HasData()) continue;
1619 if (!digitsB->HasData()) continue;
1621 for (iDict = 0; iDict < kNDict; iDict++) {
1622 dictionaryA[iDict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(iDet,iDict);
1623 dictionaryB[iDict] = (AliTRDarrayDictionary *) mergeSDigitsManager->GetDictionary(iDet,iDict);
1624 dictionaryA[iDict]->Expand();
1625 dictionaryB[iDict]->Expand();
1628 // Merge only detectors that contain a signal
1629 Bool_t doMerge = kTRUE;
1630 if (fMergeSignalOnly) {
1631 if (digitsA->GetOverThreshold(0) == 0) {
1638 AliDebug(1,Form("Merge detector %d of input no.%d",iDet,iMerge+1));
1640 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1641 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1642 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1644 // Add the amplitudes of the summable digits
1645 Float_t ampA = digitsA->GetData(iRow,iCol,iTime);
1646 Float_t ampB = digitsB->GetData(iRow,iCol,iTime);
1648 digitsA->SetData(iRow,iCol,iTime,ampA);
1650 // Add the mask to the track id if defined.
1651 for (iDict = 0; iDict < kNDict; iDict++) {
1652 Int_t trackB = dictionaryB[iDict]->GetData(iRow,iCol,iTime);
1653 if ((fMasks) && (trackB > 0)) {
1654 for (jDict = 0; jDict < kNDict; jDict++) {
1655 Int_t trackA = dictionaryA[iDict]->GetData(iRow,iCol,iTime);
1657 trackA = trackB + fMasks[iMerge];
1658 dictionaryA[iDict]->SetData(iRow,iCol,iTime,trackA);
1659 } // if: track A == 0
1661 } // if: fMasks and trackB > 0
1670 mergeSDigitsManager->RemoveDigits(iDet);
1671 mergeSDigitsManager->RemoveDictionaries(iDet);
1674 digitsA->Compress(0);
1675 for (iDict = 0; iDict < kNDict; iDict++) {
1676 dictionaryA[iDict]->Compress();
1682 // The next set of s-digits
1683 mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(mergeSDigitsManager);
1685 } // while: mergeDigitsManagers
1691 //_____________________________________________________________________________
1692 Bool_t AliTRDdigitizer::ConvertSDigits()
1695 // Converts s-digits to normal digits
1698 AliTRDarraySignal *digitsIn = 0x0;
1700 if (!fSDigitsManager->HasSDigits()) {
1701 AliError("No s-digits in digits manager");
1705 // Loop through the detectors
1706 for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1708 // Get the merged s-digits (signals)
1709 digitsIn = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(det);
1710 if (!digitsIn->HasData()) {
1711 AliDebug(2,Form("No digits for det=%d",det));
1715 // Convert the merged sdigits to digits
1716 if (!Signal2ADC(det,digitsIn)) {
1720 // Copy the dictionary information to the output array
1721 if (!CopyDictionary(det)) {
1726 fSDigitsManager->RemoveDigits(det);
1727 fSDigitsManager->RemoveDictionaries(det);
1729 // Run digital processing
1730 RunDigitalProcessing(det);
1732 // Compress the arrays
1733 CompressOutputArrays(det);
1735 } // for: detector numbers
1737 if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
1738 if (trklLoader->Tree())
1739 trklLoader->WriteData("OVERWRITE");
1742 // Save the values for the raw data headers
1743 if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
1744 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
1747 fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
1749 fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
1755 //_____________________________________________________________________________
1756 Bool_t AliTRDdigitizer::CopyDictionary(Int_t det)
1759 // Copies the dictionary information from the s-digits arrays
1760 // to the output arrays
1763 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1765 AliFatal("Could not get calibration object");
1769 AliDebug(1,Form("Start copying dictionaries for detector=%d",det));
1771 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1772 AliTRDarrayDictionary *dictionaryIn[kNDict];
1773 AliTRDarrayDictionary *dictionaryOut[kNDict];
1775 Int_t nRowMax = fGeo->GetPadPlane(det)->GetNrows();
1776 Int_t nColMax = fGeo->GetPadPlane(det)->GetNcols();
1777 Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1784 for (dict = 0; dict < kNDict; dict++) {
1786 dictionaryIn[dict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(det,dict);
1787 dictionaryIn[dict]->Expand();
1788 dictionaryOut[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1789 dictionaryOut[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
1791 for (row = 0; row < nRowMax; row++) {
1792 for (col = 0; col < nColMax; col++) {
1793 for (time = 0; time < nTimeTotal; time++) {
1794 Int_t track = dictionaryIn[dict]->GetData(row,col,time);
1795 dictionaryOut[dict]->SetData(row,col,time,track);
1800 } // for: dictionaries
1806 //_____________________________________________________________________________
1807 void AliTRDdigitizer::CompressOutputArrays(Int_t det)
1810 // Compress the output arrays
1813 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1814 AliTRDarrayDictionary *dictionary = 0x0;
1819 AliTRDarrayADC *digits = 0x0;
1820 digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1825 AliTRDarraySignal *digits = 0x0;
1826 digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1827 digits->Compress(0);
1830 for (Int_t dict = 0; dict < kNDict; dict++) {
1831 dictionary = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1832 dictionary->Compress();
1839 //_____________________________________________________________________________
1840 Bool_t AliTRDdigitizer::WriteDigits() const
1843 // Writes out the TRD-digits and the dictionaries
1847 fRunLoader->CdGAFile();
1849 // Store the digits and the dictionary in the tree
1850 return fDigitsManager->WriteDigits();
1854 //_____________________________________________________________________________
1855 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1858 // Initializes the output branches
1864 AliError("Run Loader is NULL");
1868 AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
1870 AliError("Can not get TRD loader from Run Loader");
1877 // If we produce SDigits
1878 tree = loader->TreeS();
1880 loader->MakeTree("S");
1881 tree = loader->TreeS();
1885 // If we produce Digits
1886 tree = loader->TreeD();
1888 loader->MakeTree("D");
1889 tree = loader->TreeD();
1892 fDigitsManager->SetEvent(iEvent);
1893 fDigitsManager->MakeBranch(tree);
1897 //_____________________________________________________________________________
1898 Int_t AliTRDdigitizer::Diffusion(Float_t vdrift, Double_t absdriftlength
1899 , Double_t &lRow, Double_t &lCol, Double_t &lTime)
1902 // Applies the diffusion smearing to the position of a single electron.
1903 // Depends on absolute drift length.
1906 Float_t diffL = 0.0;
1907 Float_t diffT = 0.0;
1909 if (AliTRDCommonParam::Instance()->GetDiffCoeff(diffL,diffT,vdrift)) {
1911 Float_t driftSqrt = TMath::Sqrt(absdriftlength);
1912 Float_t sigmaT = driftSqrt * diffT;
1913 Float_t sigmaL = driftSqrt * diffL;
1914 lRow = gRandom->Gaus(lRow ,sigmaT);
1915 lCol = gRandom->Gaus(lCol ,sigmaT * GetLorentzFactor(vdrift));
1916 lTime = gRandom->Gaus(lTime,sigmaL * GetLorentzFactor(vdrift));
1929 //_____________________________________________________________________________
1930 Float_t AliTRDdigitizer::GetLorentzFactor(Float_t vd)
1933 // Returns the Lorentz factor
1936 Double_t omegaTau = AliTRDCommonParam::Instance()->GetOmegaTau(vd);
1937 Double_t lorentzFactor = 1.0;
1938 if (AliTRDCommonParam::Instance()->ExBOn()) {
1939 lorentzFactor = 1.0 / (1.0 + omegaTau*omegaTau);
1942 return lorentzFactor;
1946 //_____________________________________________________________________________
1947 Int_t AliTRDdigitizer::ExB(Float_t vdrift, Double_t driftlength, Double_t &lCol)
1950 // Applies E x B effects to the position of a single electron.
1951 // Depends on signed drift length.
1955 + AliTRDCommonParam::Instance()->GetOmegaTau(vdrift)
1962 //_____________________________________________________________________________
1963 void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
1966 // Run the digital processing in the TRAP
1969 AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();
1971 AliTRDarrayADC *digits = fDigitsManager->GetDigits(det);
1975 //Call the methods in the mcm class using the temporary array as input
1976 for(Int_t rob = 0; rob < digits->GetNrow() / 2; rob++)
1978 for(Int_t mcm = 0; mcm < 16; mcm++)
1980 fMcmSim->Init(det, rob, mcm);
1981 fMcmSim->SetDataByPad(digits, fDigitsManager);
1983 if (feeParam->GetTracklet()) {
1984 fMcmSim->Tracklet();
1985 fMcmSim->StoreTracklets();
1987 fMcmSim->ZSMapping();
1988 fMcmSim->WriteData(digits);