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 **************************************************************************/
16 //_________________________________________________________________________
17 // This is a TTask that constructs SDigits out of Hits
18 // A Summable Digits is the sum of all hits in a pad
21 //-- Author: F. Pierella
22 //////////////////////////////////////////////////////////////////////////////
30 #include "AliTOFHitMap.h"
31 #include "AliTOFSDigit.h"
32 #include "AliTOFhit.h"
38 #include "AliTOFSDigitizer.h"
40 #include "AliDetector.h"
53 ClassImp(AliTOFSDigitizer)
55 //____________________________________________________________________________
56 AliTOFSDigitizer::AliTOFSDigitizer():TTask("AliTOFSDigitizer","")
65 //____________________________________________________________________________
66 AliTOFSDigitizer::AliTOFSDigitizer(char* HeaderFile,char *SdigitsFile ):TTask("AliTOFSDigitizer","")
68 fNevents = 0 ; // Number of events to digitize, 0 means all evens in current file
69 // add Task to //root/Tasks folder
70 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
71 roottasks->Add(this) ;
74 //____________________________________________________________________________
75 AliTOFSDigitizer::~AliTOFSDigitizer()
80 //____________________________________________________________________________
81 void AliTOFSDigitizer::Exec(Option_t *option) {
84 AliTOF *TOF = (AliTOF *) gAlice->GetDetector ("TOF");
87 Error("AliTOFSDigitizer","TOF not found");
92 fNevents = (Int_t) gAlice->TreeE()->GetEntries();
94 for (Int_t ievent = 0; ievent < fNevents; ievent++) {
95 gAlice->GetEvent(ievent);
96 TTree *TH = gAlice->TreeH ();
99 if (gAlice->TreeS () == 0)
100 gAlice->MakeTree ("S");
105 sprintf (branchname, "%s", TOF->GetName ());
106 //Make branch for digits
107 TOF->MakeBranch ("S");
109 //Now made SDigits from hits
111 Int_t vol[5]; // location for a digit
112 Float_t digit[2]; // TOF digit variables
115 TClonesArray *TOFhits = TOF->Hits();
118 AliTOFHitMap *hitMap = new AliTOFHitMap(TOF->SDigits());
120 Int_t ntracks = static_cast<Int_t>(TH->GetEntries());
121 for (Int_t track = 0; track < ntracks; track++)
125 particle = gAlice->Particle(track);
126 Int_t nhits = TOFhits->GetEntriesFast();
128 for (Int_t hit = 0; hit < nhits; hit++)
130 tofHit = (AliTOFhit *) TOFhits->UncheckedAt(hit);
131 vol[0] = tofHit->GetSector();
132 vol[1] = tofHit->GetPlate();
133 vol[2] = tofHit->GetStrip();
134 vol[3] = tofHit->GetPadx();
135 vol[4] = tofHit->GetPadz();
137 // 95% of efficiency to be inserted here
138 // edge effect to be inserted here
139 // cross talk to be inserted here
141 Float_t idealtime = tofHit->GetTof(); // unit s
142 idealtime *= 1.E+12; // conversion from s to ps
143 // fTimeRes is given usually in ps
144 Float_t tdctime = gRandom->Gaus(idealtime, TOF->GetTimeRes());
147 // typical Landau Distribution to be inserted here
148 // instead of Gaussian Distribution
149 Float_t idealcharge = tofHit->GetEdep();
150 Float_t adccharge = gRandom->Gaus(idealcharge, TOF->GetChrgRes());
151 digit[1] = adccharge;
152 Int_t tracknum = tofHit->GetTrack();
154 // check if two digit are on the same pad; in that case we sum
155 // the two or more digits
156 if (hitMap->TestHit(vol) != kEmpty) {
157 AliTOFSDigit *sdig = static_cast<AliTOFSDigit*>(hitMap->GetHit(vol));
158 sdig->Update(tdctime,adccharge,tracknum);
160 TOF->AddSDigit(tracknum, vol, digit);
163 } // end loop on hits for the current track
164 } // end loop on ntracks
168 gAlice->TreeS()->Reset();
169 gAlice->TreeS()->Fill();
170 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
176 //__________________________________________________________________
177 void AliTOFSDigitizer::SetSDigitsFile(char * file ){
178 if(!fSDigitsFile.IsNull())
179 cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ;
182 //__________________________________________________________________
183 void AliTOFSDigitizer::Print(Option_t* option)const
185 cout << "------------------- "<< GetName() << " -------------" << endl ;
186 if(fSDigitsFile.IsNull())
187 cout << " Writing SDigitis to file galice.root "<< endl ;
189 cout << " Writing SDigitis to file " << (char*) fSDigitsFile.Data() << endl ;