]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFDigitizer.cxx
TTask and TFolder structures implemented
[u/mrichter/AliRoot.git] / TOF / AliTOFDigitizer.cxx
CommitLineData
68861244 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16
17//_________________________________________________________________________
18// This is a TTask that makes TOF-Digits out of TOF-Hits.
19// A TOF Digit is essentially a TOF hit with the smearing for strip
20// time resolution and simulation of the ADC correlation signal.
21// Digits are written to TreeD in branch "TOF".
22// AliTOFDigitizer with all current parameters used for digitization
23// (time resolution and ADC parameter) is written
24// to TreeD branch "AliTOFDigitizer".
25// Both branches have the same title. If necessary one can produce
26// another set of Digits with different parameters. Two versions
27// can be distunguished using titles of the branches.
28// User case:
29// root [0] AliTOFDigitizer * s = new AliTOFDigitizer("galice.root")
30// Warning in <AliITSgeomSPD425Long::Default Creator>: Detector size may not be write.
31// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
32// root [1] s->ExecuteTask()
33// // Makes Digits for all events stored in galice.root
34// root [2] s->SetTimeRes(100.)
35// // One can change parameters of digitization
36// root [3] s->SetDigitsBranch("Time Resolution 100. ps")
37// // and write them into the new branch
38// root [4] s->ExecuteTask("deb all tim")
39// // available parameters:
40// deb - print number of produced Digits
41// deb all - print number and list of produced Digits
42// tim - print benchmarking information
43//
44// -- Author : F. Pierella (Bologna University) pierella@bo.infn.it
45//////////////////////////////////////////////////////////////////////////////
46
47#include "TFile.h"
48#include "TTask.h"
49#include "TTree.h"
50#include "TSystem.h"
51#include "TROOT.h"
52#include "TFolder.h"
53#include "TBenchmark.h"
54#include "TRandom.h"
55
56#include <iomanip.h>
57
58#include "AliRun.h"
59#include "AliTOFdigit.h"
60#include "AliTOFhit.h"
61#include "AliTOFDigitizer.h"
62
63
64ClassImp(AliTOFDigitizer)
65
66
67//____________________________________________________________________________
68 AliTOFDigitizer::AliTOFDigitizer():TTask("AliTOFDigitizer","")
69{
70 // ctor
71 fTimeRes = 100; // ps
72 fChrgRes = 100.; // pC
73 fNevents = 0 ;
74 fDigits = 0 ;
75 fHits = 0 ;
76 fIsInitialized = kFALSE ;
77
78}
79
80//____________________________________________________________________________
81AliTOFDigitizer::AliTOFDigitizer(const char* headerFile, const char *digitsTitle):TTask("AliTOFDigitizer","")
82{
83 // ctor
84 fTimeRes = 100; // ps
85 fChrgRes = 100.; // pC
86 fNevents = 0 ;
87 fDigitsTitle = digitsTitle ;
88 fHeadersFile = headerFile ;
89 fIsInitialized = kFALSE ;
90 Init();
91}
92
93//____________________________________________________________________________
94AliTOFDigitizer::~AliTOFDigitizer()
95{
96 // dtor
97 if(fDigits)
98 delete fDigits ;
99 if(fHits)
100 delete fHits ;
101}
102//____________________________________________________________________________
103void AliTOFDigitizer::Init()
104{
105 // Initialization: open root-file, allocate arrays for hits and digits,
106 // attach task Digitizer to the list of TOF tasks
107 //
108 // Initialization can not be done in the default constructor
109
110 if(!fIsInitialized){
111
112 if(fHeadersFile.IsNull())
113 fHeadersFile="galice.root" ;
114
115 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
116
117 //if file was not opened yet, read gAlice
118 if(file == 0){
119 if(fHeadersFile.Contains("rfio"))
120 file = TFile::Open(fHeadersFile,"update") ;
121 else
122 file = new TFile(fHeadersFile.Data(),"update") ;
123 gAlice = (AliRun *) file->Get("gAlice") ;
124 }
125 fHits = new TClonesArray("AliTOFhit" , 405);
126 fDigits = new TClonesArray("AliTOFdigit", 405);
127
128 //add Task to //FPAlice/tasks/(S)Digitizer/TOF
129 TFolder * alice = (TFolder*)gROOT->GetListOfBrowsables()->FindObject("FPAlice") ;
130 TTask * aliceSD = (TTask*)alice->FindObject("tasks/(S)Digitizer") ;
131 TTask * tofD = (TTask*)aliceSD->GetListOfTasks()->FindObject("TOF") ;
132 tofD->Add(this) ;
133
134 fIsInitialized = kTRUE ;
135 }
136}
137//____________________________________________________________________________
138void AliTOFDigitizer::Exec(Option_t *option)
139{
140 Int_t tracks[3]; // track info
141 Int_t vol[5]; // dummy location for digit
142 Float_t digit[2]; // TOF digit variables
143
144 TRandom* rnd = new TRandom();
145
146 // Collects all hits in the same active volume into digit
147
148 if(!fIsInitialized)
149 Init() ;
150
151 if(strstr(option,"tim"))
152 gBenchmark->Start("TOFDigitizer");
153
154 fNevents = (Int_t) gAlice->TreeE()->GetEntries() ;
155
156 Int_t ievent ;
157 // start loop on events
158 for(ievent = 0; ievent < fNevents; ievent++){
159 gAlice->GetEvent(ievent) ;
160 gAlice->SetEvent(ievent) ;
161
162 if(gAlice->TreeH()==0){
163 cout << "AliTOFDigitizer: There is no Hit Tree" << endl;
164 return ;
165 }
166
167 //set address of the hits
168 TBranch * branch = gAlice->TreeH()->GetBranch("TOF");
169 if (branch)
170 branch->SetAddress(&fHits);
171 else{
172 cout << "ERROR in AliTOFDigitizer: "<< endl ;
173 cout << " no branch TOF in TreeH"<< endl ;
174 cout << " do nothing " << endl ;
175 return ;
176 }
177
178 fDigits->Clear();
179 Int_t ndigits = 0 ;
180
181 //Now made Digits from hits, for TOF it is the same a part for the tof smearing
182 // and some missing MC variables
183 Int_t itrack ;
184 for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
185
186 //=========== Get the TOF branch from Hits Tree for the Primary track itrack
187 branch->GetEntry(itrack,0);
188
189 Int_t i;
190 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
191
192 AliTOFhit * hit = (AliTOFhit*)fHits->At(i) ;
193
194 vol[0] = hit->GetSector();
195 vol[1] = hit->GetPlate();
196 vol[2] = hit->GetPadx();
197 vol[3] = hit->GetPadz();
198 vol[4] = hit->GetStrip();
199 // 95% of efficiency to be inserted here
200 // edge effect to be inserted here
201 // cross talk to be inserted here
202 // simulation of the detector response
203 Float_t idealtime = hit->GetTof(); // unit s
204 idealtime *= 1.E+12; // conversion from s to ps
205 // fTimeRes is given usually in ps
206 Float_t tdctime = rnd->Gaus(idealtime, fTimeRes);
207 digit[0] = tdctime;
208 // typical Landau Distribution to be inserted here
209 Float_t idealcharge = hit->GetEdep();
210 Float_t adccharge = rnd->Gaus(idealcharge, fChrgRes);
211 digit[1] = adccharge;
212 // to be added a check for overlaps
213 new((*fDigits)[ndigits]) AliTOFdigit(tracks, vol, digit);
214 ndigits++ ;
215 }
216
217 } // end loop over tracks
218
219 delete rnd;
220 rnd = 0;
221
222 ndigits = fDigits->GetEntriesFast() ;
223 printf("AliTOFDigitizer: Total number of digits %d\n",ndigits);
224
225 if(gAlice->TreeD() == 0)
226 gAlice->MakeTree("D") ;
227
228 //check, if this branch already exits
229 TBranch * digitsBranch = 0;
230 TBranch * digitizerBranch = 0;
231
232 TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
233 Int_t ibranch;
234 Bool_t tofNotFound = kTRUE ;
235 Bool_t digitizerNotFound = kTRUE ;
236
237 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
238
239 if(tofNotFound){
240 digitsBranch=(TBranch *) branches->At(ibranch) ;
241 if( (strcmp("TOF",digitsBranch->GetName())==0 ) &&
242 (fDigitsTitle.CompareTo(digitsBranch->GetTitle()) == 0) )
243 tofNotFound = kFALSE ;
244 }
245 if(digitizerNotFound){
246 digitizerBranch = (TBranch *) branches->At(ibranch) ;
247 if( (strcmp(digitizerBranch->GetName(),"AliTOFDigitizer") == 0)&&
248 (fDigitsTitle.CompareTo(digitizerBranch->GetTitle()) == 0) )
249 digitizerNotFound = kFALSE ;
250 }
251 }
252
253 if(!(digitizerNotFound && tofNotFound)){
254 cout << "AliTOFdigitizer error:" << endl ;
255 cout << "Can not overwrite existing branches: do not write" << endl ;
256 return ;
257 }
258
259 //Make (if necessary) branches
260 char * file =0;
261 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
262 file = new char[strlen(gAlice->GetBaseFile())+20] ;
263 sprintf(file,"%s/TOF.Digits.root",gAlice->GetBaseFile()) ;
264 }
265
266 TDirectory *cwd = gDirectory;
267
268 //First list of digits
269 Int_t bufferSize = 32000 ;
270 digitsBranch = gAlice->TreeD()->Branch("TOF",&fDigits,bufferSize);
271 digitsBranch->SetTitle(fDigitsTitle.Data());
272 if (file) {
273 digitsBranch->SetFile(file);
274 TIter next( digitsBranch->GetListOfBranches());
275 TBranch * subbr;
276 while ((subbr=(TBranch*)next())) {
277 subbr->SetFile(file);
278 }
279 cwd->cd();
280 }
281
282 //second - Digitizer
283 Int_t splitlevel = 0 ;
284 AliTOFDigitizer * digtz = this ;
285 digitizerBranch = gAlice->TreeD()->Branch("AliTOFDigitizer","AliTOFDigitizer",
286 &digtz,bufferSize,splitlevel);
287 digitizerBranch->SetTitle(fDigitsTitle.Data());
288 if (file) {
289 digitizerBranch->SetFile(file);
290 TIter next( digitizerBranch->GetListOfBranches());
291 TBranch * subbr ;
292 while ((subbr=(TBranch*)next())) {
293 subbr->SetFile(file);
294 }
295 cwd->cd();
296 delete file;
297 }
298
299 digitsBranch->Fill();
300 digitizerBranch->Fill();
301 gAlice->TreeD()->Write(0,TObject::kOverwrite) ;
302
303 if(strstr(option,"deb"))
304 PrintDigits(option) ;
305
306 }
307
308 if(strstr(option,"tim")){
309 gBenchmark->Stop("TOFDigitizer");
310 cout << "AliTOFDigitizer:" << endl ;
311 cout << " took " << gBenchmark->GetCpuTime("TOFDigitizer") << " seconds for Digitizing "
312 << gBenchmark->GetCpuTime("TOFDigitizer")/fNevents << " seconds per event " << endl ;
313 cout << endl ;
314 }
315
316
317}
318//__________________________________________________________________
319void AliTOFDigitizer::SetDigitsBranch(const char* title )
320{
321 // Setting title to branch Digits
322 if(!fDigitsTitle.IsNull())
323 cout << "AliTOFdigitizer: changing Digits file from " <<fDigitsTitle.Data() << " to " << title << endl ;
324 fDigitsTitle=title ;
325}
326//__________________________________________________________________
327void AliTOFDigitizer::Print(Option_t* option)const
328{
329 // Prints parameters of Digitizer
330 cout << "------------------- "<< GetName() << " -------------" << endl ;
331 cout << " Writing Digits to branch with title " << fDigitsTitle.Data() << endl ;
332 cout << " with digitization parameters Time resolution = " << fTimeRes << endl ;
333 cout << " Adc smearing parameter = " << fChrgRes << endl ;
334 cout << "---------------------------------------------------"<<endl ;
335
336}
337//__________________________________________________________________
338Bool_t AliTOFDigitizer::operator==( AliTOFDigitizer const &digtz )const
339{
340 // Equal operator.
341 // Dititizers are equal if their time resolution and Adc
342 // smearing parameter are equal
343
344 if( (fTimeRes==digtz.fTimeRes)&&(fChrgRes==digtz.fChrgRes))
345 return kTRUE ;
346 else
347 return kFALSE ;
348}
349//__________________________________________________________________
350void AliTOFDigitizer::PrintDigits(Option_t * option)
351{
352 // Prints list of digits produced in the current pass of AliTOFDigitizer
353
354 cout << "AliTOFDigitizer: " << endl ;
355 cout << " Number of entries in Digits list " << fDigits->GetEntriesFast() << endl ;
356 cout << endl ;
357
358 if(strstr(option,"all")){// print all digits
359
360 //loop over digits
361 AliTOFdigit * digit;
362 cout << "Digit # " << " Time of flight(ps) " <<
363 " ADC(pC)" << " Sector #" << " Plate #" <<
364 " Strip #" << " Pad # (x) " << " Pad # (z) " << endl;
365 Int_t index ;
366 for (index = 0 ; index < fDigits->GetEntries() ; index++) {
367 digit = (AliTOFdigit * ) fDigits->At(index) ;
368
369 // set the width of the output with the setw(int width) manipulator
370
371 cout << setw(7) << index << " "
372 << setw(13) << digit->GetTdc() << " "
373 << setw(13) << digit->GetAdc() << " "
374 << setw(6) << digit->GetSector() << " "
375 << setw(6) << digit->GetPlate() << " "
376 << setw(6) << digit->GetStrip() << " "
377 << setw(7) << digit->GetPadx() << " "
378 << setw(7) << digit->GetPadz() << " " << endl;
379
380 } // end loop on digits
381
382 } // close if "all" option selected
383}