]>
Commit | Line | Data |
---|---|---|
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 | /* $Id$ */ | |
18 | ||
19 | //_________________________________________________________________________ | |
20 | // This is a TTask that makes SDigits out of Hits | |
21 | // The name of the TTask is also the title of the branch that will contain | |
22 | // the created SDigits | |
23 | // The title of the TTAsk is the name of the file that contains the hits from | |
24 | // which the SDigits are created | |
25 | // A Summable Digits is the sum of all hits originating | |
26 | // from one primary in one active cell | |
27 | // A threshold for assignment of the primary to SDigit is applied | |
28 | // SDigits are written to TreeS, branch "PHOS" | |
29 | // AliPHOSSDigitizer with all current parameters is written | |
30 | // to TreeS branch "AliPHOSSDigitizer". | |
31 | // Both branches have the same title. If necessary one can produce | |
32 | // another set of SDigits with different parameters. Two versions | |
33 | // can be distunguished using titles of the branches. | |
34 | // User case: | |
35 | // root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root") | |
36 | // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated | |
37 | // root [1] s->ExecuteTask() | |
38 | // // Makes SDigitis for all events stored in galice.root | |
39 | // root [2] s->SetPedestalParameter(0.001) | |
40 | // // One can change parameters of digitization | |
41 | // root [3] s->SetSDigitsBranch("Pedestal 0.001") | |
42 | // // and write them into the new branch | |
43 | // root [4] s->ExecuteTask("deb all tim") | |
44 | // // available parameters: | |
45 | // deb - print # of produced SDigitis | |
46 | // deb all - print # and list of produced SDigits | |
47 | // tim - print benchmarking information | |
48 | // | |
49 | //*-- Author : Dmitri Peressounko (SUBATECH & KI) | |
50 | ////////////////////////////////////////////////////////////////////////////// | |
51 | ||
52 | ||
53 | // --- ROOT system --- | |
54 | #include "TFile.h" | |
55 | #include "TTask.h" | |
56 | #include "TTree.h" | |
57 | #include "TSystem.h" | |
58 | #include "TROOT.h" | |
59 | #include "TFolder.h" | |
60 | #include "TBenchmark.h" | |
61 | #include "TGeometry.h" | |
62 | ||
63 | // --- Standard library --- | |
64 | #include <iomanip.h> | |
65 | ||
66 | // --- AliRoot header files --- | |
67 | #include "AliRun.h" | |
68 | #include "AliHeader.h" | |
69 | #include "AliPHOSDigit.h" | |
70 | #include "AliPHOSGeometry.h" | |
71 | #include "AliPHOSGetter.h" | |
72 | #include "AliPHOSHit.h" | |
73 | #include "AliPHOSSDigitizer.h" | |
74 | ||
75 | ||
76 | ClassImp(AliPHOSSDigitizer) | |
77 | ||
78 | ||
79 | //____________________________________________________________________________ | |
80 | AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") | |
81 | { | |
82 | // ctor | |
83 | fA = 0; | |
84 | fB = 10000000.; | |
85 | fPrimThreshold = 0.01 ; | |
86 | fSDigitsInRun = 0 ; | |
87 | fSplitFile = 0 ; | |
88 | } | |
89 | ||
90 | //____________________________________________________________________________ | |
91 | AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile) | |
92 | { | |
93 | // ctor | |
94 | fA = 0; | |
95 | fB = 10000000.; | |
96 | fPrimThreshold = 0.01 ; | |
97 | fSDigitsInRun = 0 ; | |
98 | fSplitFile = 0 ; | |
99 | Init(); | |
100 | } | |
101 | ||
102 | //____________________________________________________________________________ | |
103 | AliPHOSSDigitizer::~AliPHOSSDigitizer() | |
104 | { | |
105 | if (fSplitFile) | |
106 | if ( fSplitFile->IsOpen() ) | |
107 | fSplitFile->Close() ; | |
108 | AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; | |
109 | gime->CloseFile() ; | |
110 | } | |
111 | ||
112 | //____________________________________________________________________________ | |
113 | void AliPHOSSDigitizer::Init() | |
114 | { | |
115 | // Initialization: open root-file, allocate arrays for hits and sdigits, | |
116 | // attach task SDigitizer to the list of PHOS tasks | |
117 | // | |
118 | // Initialization can not be done in the default constructor | |
119 | //============================================================= YS | |
120 | // The initialisation is now done by AliPHOSGetter | |
121 | ||
122 | if( strcmp(GetTitle(), "") == 0 ) | |
123 | SetTitle("galice.root") ; | |
124 | ||
125 | AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ; | |
126 | if ( gime == 0 ) { | |
127 | cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ; | |
128 | return ; | |
129 | } | |
130 | ||
131 | gime->PostSDigits( GetName(), GetTitle() ) ; | |
132 | ||
133 | TString sdname(GetName() ); | |
134 | sdname.Append(":") ; | |
135 | sdname.Append(GetTitle() ) ; | |
136 | SetName(sdname) ; | |
137 | gime->PostSDigitizer(this) ; | |
138 | } | |
139 | ||
140 | //____________________________________________________________________________ | |
141 | void AliPHOSSDigitizer::Exec(Option_t *option) | |
142 | { | |
143 | // Collects all hits in the same active volume into digit | |
144 | if( strcmp(GetName(), "") == 0 ) | |
145 | Init() ; | |
146 | ||
147 | if (strstr(option, "print") ) { | |
148 | Print("") ; | |
149 | return ; | |
150 | } | |
151 | ||
152 | if(strstr(option,"tim")) | |
153 | gBenchmark->Start("PHOSSDigitizer"); | |
154 | ||
155 | //Check, if this branch already exits | |
156 | gAlice->GetEvent(0) ; | |
157 | ||
158 | TString sdname(GetName()) ; | |
159 | sdname.Remove(sdname.Index(GetTitle())-1) ; | |
160 | ||
161 | if(gAlice->TreeS() ) { | |
162 | TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ; | |
163 | TIter next(lob) ; | |
164 | TBranch * branch = 0 ; | |
165 | Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; | |
166 | ||
167 | while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) { | |
168 | TString thisName( GetName() ) ; | |
169 | TString branchName( branch->GetTitle() ) ; | |
170 | branchName.Append(":") ; | |
171 | if ( (strcmp(branch->GetName(), "PHOS")==0) && thisName.BeginsWith(branchName) ) | |
172 | phosfound = kTRUE ; | |
173 | else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && thisName.BeginsWith(branchName) ) | |
174 | sdigitizerfound = kTRUE ; | |
175 | } | |
176 | ||
177 | if ( phosfound || sdigitizerfound ) { | |
178 | cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() | |
179 | << " already exits" << endl ; | |
180 | return ; | |
181 | } | |
182 | } | |
183 | ||
184 | AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; | |
185 | ||
186 | Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; | |
187 | ||
188 | Int_t ievent ; | |
189 | for(ievent = 0; ievent < nevents; ievent++){ | |
190 | gime->Event(ievent,"H") ; | |
191 | const TClonesArray * hits = gime->Hits() ; | |
192 | TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; | |
193 | sdigits->Clear(); | |
194 | Int_t nSdigits = 0 ; | |
195 | ||
196 | ||
197 | //Now make SDigits from hits, for PHOS it is the same, so just copy | |
198 | ||
199 | Int_t Nprim = (Int_t) (gAlice->TreeH())->GetEntries(); | |
200 | // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK; | |
201 | Int_t iprim; | |
202 | for (iprim=0; iprim<Nprim; iprim++) { | |
203 | //=========== Get the PHOS branch from Hits Tree for the Primary iprim | |
204 | gime->Track(iprim) ; | |
205 | Int_t i; | |
206 | for ( i = 0 ; i < hits->GetEntries() ; i++ ) { | |
207 | AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ; | |
208 | // Assign primary number only if contribution is significant | |
209 | ||
210 | if( hit->GetEnergy() > fPrimThreshold) | |
211 | new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(), | |
212 | Digitize(hit->GetEnergy()), hit->GetTime()) ; | |
213 | else | |
214 | new((*sdigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(), | |
215 | Digitize(hit->GetEnergy()), hit->GetTime()) ; | |
216 | nSdigits++ ; | |
217 | ||
218 | } | |
219 | ||
220 | } // loop over iprim | |
221 | ||
222 | sdigits->Sort() ; | |
223 | ||
224 | nSdigits = sdigits->GetEntriesFast() ; | |
225 | fSDigitsInRun += nSdigits ; | |
226 | sdigits->Expand(nSdigits) ; | |
227 | ||
228 | Int_t i ; | |
229 | for (i = 0 ; i < nSdigits ; i++) { | |
230 | AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; | |
231 | digit->SetIndexInList(i) ; | |
232 | } | |
233 | ||
234 | if(gAlice->TreeS() == 0) | |
235 | gAlice->MakeTree("S", fSplitFile); | |
236 | ||
237 | //First list of sdigits | |
238 | Int_t bufferSize = 32000 ; | |
239 | TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize); | |
240 | sdigitsBranch->SetTitle(sdname); | |
241 | ||
242 | //Next - SDigitizer | |
243 | Int_t splitlevel = 0 ; | |
244 | AliPHOSSDigitizer * sd = this ; | |
245 | TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer", | |
246 | &sd,bufferSize,splitlevel); | |
247 | sdigitizerBranch->SetTitle(sdname); | |
248 | ||
249 | sdigitsBranch->Fill() ; | |
250 | sdigitizerBranch->Fill() ; | |
251 | gAlice->TreeS()->AutoSave() ; | |
252 | ||
253 | if(strstr(option,"deb")) | |
254 | PrintSDigits(option) ; | |
255 | ||
256 | } | |
257 | ||
258 | if (fSplitFile) | |
259 | if ( fSplitFile->IsOpen() ) | |
260 | fSplitFile->Close() ; | |
261 | ||
262 | ||
263 | if(strstr(option,"tim")){ | |
264 | gBenchmark->Stop("PHOSSDigitizer"); | |
265 | cout << "AliPHOSSDigitizer:" << endl ; | |
266 | cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " | |
267 | << gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ; | |
268 | cout << endl ; | |
269 | } | |
270 | ||
271 | ||
272 | } | |
273 | ||
274 | //__________________________________________________________________ | |
275 | void AliPHOSSDigitizer::SetSDigitsBranch(const char * title ) | |
276 | { | |
277 | // Setting title to branch SDigits | |
278 | ||
279 | TString stitle(title) ; | |
280 | ||
281 | // check if branch with title already exists | |
282 | TBranch * sdigitsBranch = | |
283 | static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; | |
284 | TBranch * sdigitizerBranch = | |
285 | static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ; | |
286 | const char * sdigitsTitle = sdigitsBranch ->GetTitle() ; | |
287 | const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ; | |
288 | if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){ | |
289 | cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ; | |
290 | return ; | |
291 | } | |
292 | ||
293 | cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ; | |
294 | ||
295 | SetName(title) ; | |
296 | ||
297 | // Post to the WhiteBoard | |
298 | AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; | |
299 | gime->PostSDigits( title, GetTitle()) ; | |
300 | } | |
301 | ||
302 | //__________________________________________________________________ | |
303 | void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName) | |
304 | { | |
305 | // Diverts the SDigits in a file separate from the hits file | |
306 | ||
307 | TDirectory * cwd = gDirectory ; | |
308 | ||
309 | if ( !(gAlice->GetTreeSFileName() == splitFileName) ) { | |
310 | if (gAlice->GetTreeSFile() ) | |
311 | gAlice->GetTreeSFile()->Close() ; | |
312 | } | |
313 | ||
314 | fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data()); | |
315 | fSplitFile->cd() ; | |
316 | if ( !fSplitFile->Get("gAlice") ) | |
317 | gAlice->Write(); | |
318 | ||
319 | TTree *treeE = gAlice->TreeE(); | |
320 | if (!treeE) { | |
321 | cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl; | |
322 | abort() ; | |
323 | } | |
324 | ||
325 | // copy TreeE | |
326 | if ( !fSplitFile->Get("TreeE") ) { | |
327 | AliHeader *header = new AliHeader(); | |
328 | treeE->SetBranchAddress("Header", &header); | |
329 | treeE->SetBranchStatus("*",1); | |
330 | TTree *treeENew = treeE->CloneTree(); | |
331 | treeENew->Write(); | |
332 | } | |
333 | ||
334 | // copy AliceGeom | |
335 | if ( !fSplitFile->Get("AliceGeom") ) { | |
336 | TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom")); | |
337 | if (!AliceGeom) { | |
338 | cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl; | |
339 | abort() ; | |
340 | } | |
341 | AliceGeom->Write(); | |
342 | } | |
343 | ||
344 | gAlice->MakeTree("S",fSplitFile); | |
345 | cwd->cd() ; | |
346 | } | |
347 | ||
348 | //__________________________________________________________________ | |
349 | void AliPHOSSDigitizer::Print(Option_t* option)const | |
350 | { | |
351 | // Prints parameters of SDigitizer | |
352 | cout << "------------------- "<< GetName() << " -------------" << endl ; | |
353 | cout << " Writing SDigits to branch with title " << GetName() << endl ; | |
354 | cout << " with digitization parameters A = " << fA << endl ; | |
355 | cout << " B = " << fB << endl ; | |
356 | cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ; | |
357 | cout << "---------------------------------------------------"<<endl ; | |
358 | ||
359 | } | |
360 | //__________________________________________________________________ | |
361 | Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const | |
362 | { | |
363 | // Equal operator. | |
364 | // SDititizers are equal if their pedestal, slope and threshold are equal | |
365 | ||
366 | if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold)) | |
367 | return kTRUE ; | |
368 | else | |
369 | return kFALSE ; | |
370 | } | |
371 | //__________________________________________________________________ | |
372 | //__________________________________________________________________ | |
373 | void AliPHOSSDigitizer::PrintSDigits(Option_t * option) | |
374 | { | |
375 | // Prints list of digits produced in the current pass of AliPHOSDigitizer | |
376 | ||
377 | ||
378 | AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; | |
379 | TString sdname(GetName()) ; | |
380 | sdname.Remove(sdname.Index(GetTitle())-1) ; | |
381 | const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; | |
382 | ||
383 | cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ; | |
384 | cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ; | |
385 | cout << endl ; | |
386 | if(strstr(option,"all")||strstr(option,"EMC")){ | |
387 | ||
388 | //loop over digits | |
389 | AliPHOSDigit * digit; | |
390 | cout << "EMC sdigits " << endl ; | |
391 | cout << "Digit Id Amplitude Index Nprim Primaries list " << endl; | |
392 | Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ; | |
393 | Int_t index ; | |
394 | for (index = 0 ; (index < sdigits->GetEntriesFast()) && | |
395 | (((AliPHOSDigit * ) sdigits->At(index))->GetId() <= maxEmc) ; index++) { | |
396 | digit = (AliPHOSDigit * ) sdigits->At(index) ; | |
397 | if(digit->GetNprimary() == 0) continue; | |
398 | cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " " | |
399 | << setw(6) << digit->GetIndexInList() << " " | |
400 | << setw(5) << digit->GetNprimary() <<" "; | |
401 | ||
402 | Int_t iprimary; | |
403 | for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) | |
404 | cout << setw(5) << digit->GetPrimary(iprimary+1) << " "; | |
405 | cout << endl; | |
406 | } | |
407 | cout << endl; | |
408 | } | |
409 | ||
410 | if(strstr(option,"all")||strstr(option,"CPV")){ | |
411 | ||
412 | //loop over CPV digits | |
413 | AliPHOSDigit * digit; | |
414 | cout << "CPV sdigits " << endl ; | |
415 | cout << "Digit Id Amplitude Index Nprim Primaries list " << endl; | |
416 | Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ; | |
417 | Int_t index ; | |
418 | for (index = 0 ; index < sdigits->GetEntriesFast(); index++) { | |
419 | digit = (AliPHOSDigit * ) sdigits->At(index) ; | |
420 | if(digit->GetId() > maxEmc){ | |
421 | cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " " | |
422 | << setw(6) << digit->GetIndexInList() << " " | |
423 | << setw(5) << digit->GetNprimary() <<" "; | |
424 | ||
425 | Int_t iprimary; | |
426 | for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) | |
427 | cout << setw(5) << digit->GetPrimary(iprimary+1) << " "; | |
428 | cout << endl; | |
429 | } | |
430 | } | |
431 | } | |
432 | ||
433 | } | |
434 | ||
435 | //____________________________________________________________________________ | |
436 | void AliPHOSSDigitizer::UseHitsFrom(const char * filename) | |
437 | { | |
438 | SetTitle(filename) ; | |
439 | Init() ; | |
440 | } |