]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSDigitizer.cxx
ReadField() used in Constructor.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSDigitizer.cxx
CommitLineData
990119d6 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/* $Id$ */
17
18//_________________________________________________________________________
990119d6 19//*-- Author : Dmitri Peressounko (SUBATECH & Kurchatov Institute)
20//////////////////////////////////////////////////////////////////////////////
a4e98857 21// This TTask performs digitization of Summable digits (in the PHOS case it is just
22// the sum of contributions from all primary particles into a given cell).
990119d6 23// In addition it performs mixing of summable digits from different events.
7b7c1533 24// The name of the TTask is also the title of the branch that will contain
25// the created SDigits
26// The title of the TTAsk is the name of the file that contains the hits from
27// which the SDigits are created
bca3b32a 28//
29// For each event two branches are created in TreeD:
30// "PHOS" - list of digits
31// "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
32//
a4e98857 33// Note, that one can set a title for new digits branch, and repeat digitization with
bca3b32a 34// another set of parameters.
35//
a4e98857 36// Use case:
990119d6 37// root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
38// root[1] d->ExecuteTask()
39// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
40// //Digitizes SDigitis in all events found in file galice.root
bca3b32a 41//
8cb3533f 42// root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
43// // Will read sdigits from galice1.root
44// root[3] d1->MixWith("galice2.root")
990119d6 45// Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
a4e98857 46// // Reads another set of sdigits from galice2.root
8cb3533f 47// root[3] d1->MixWith("galice3.root")
a4e98857 48// // Reads another set of sdigits from galice3.root
8cb3533f 49// root[4] d->ExecuteTask("deb timing")
50// // Reads SDigits from files galice1.root, galice2.root ....
51// // mixes them and stores produced Digits in file galice1.root
52// // deb - prints number of produced digits
53// // deb all - prints list of produced digits
54// // timing - prints time used for digitization
990119d6 55//
990119d6 56
57// --- ROOT system ---
8cb3533f 58#include "TFile.h"
990119d6 59#include "TTree.h"
60#include "TSystem.h"
8cb3533f 61#include "TROOT.h"
62#include "TFolder.h"
63#include "TObjString.h"
64#include "TBenchmark.h"
ba54256b 65
990119d6 66// --- Standard library ---
8cb3533f 67#include <iomanip.h>
990119d6 68
69// --- AliRoot header files ---
70
8cb3533f 71#include "AliRun.h"
3f81a70b 72#include "AliRunDigitizer.h"
990119d6 73#include "AliPHOSDigit.h"
dd5c4038 74#include "AliPHOS.h"
7b7c1533 75#include "AliPHOSGetter.h"
990119d6 76#include "AliPHOSDigitizer.h"
77#include "AliPHOSSDigitizer.h"
8cb3533f 78#include "AliPHOSGeometry.h"
7437a0f7 79#include "AliPHOSTick.h"
990119d6 80
81ClassImp(AliPHOSDigitizer)
82
83
84//____________________________________________________________________________
3f81a70b 85 AliPHOSDigitizer::AliPHOSDigitizer()
990119d6 86{
87 // ctor
88
2bd5457f 89 fPinNoise = 0.01 ;
90 fEMCDigitThreshold = 0.01 ;
91 fCPVNoise = 0.01;
92 fCPVDigitThreshold = 0.09 ;
aaf8a71c 93 fTimeResolution = 0.5e-9 ;
7437a0f7 94 fTimeSignalLength = 1.0e-9 ;
2b60655b 95 fDigitsInRun = 0 ;
3758d9fc 96 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
97 fADCpedestalEmc = 0.005 ; //
98 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
99
100 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
101 fADCpedestalCpv = 0.012 ; //
102 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
103
104 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
9891b76e 105 fManager = 0 ; // We work in the standalong mode
3f81a70b 106}
107
108//____________________________________________________________________________
109AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char * name)
110{
111 // ctor
112 SetName(name) ;
113 SetTitle(headerFile) ;
9891b76e 114 fManager = 0 ; // We work in the standalong mode
3f81a70b 115 Init() ;
116
990119d6 117}
118
990119d6 119//____________________________________________________________________________
9891b76e 120AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * ard):AliDigitizer(ard)
990119d6 121{
122 // ctor
93aeb834 123 SetName(""); //Will call init in the digitizing
124 SetTitle("aliroot") ;
990119d6 125}
126
127//____________________________________________________________________________
128 AliPHOSDigitizer::~AliPHOSDigitizer()
129{
130 // dtor
8cb3533f 131
8cb3533f 132
990119d6 133}
134
135//____________________________________________________________________________
7b7c1533 136void AliPHOSDigitizer::Digitize(const Int_t event)
a4e98857 137{
138
139 // Makes the digitization of the collected summable digits.
140 // It first creates the array of all PHOS modules
141 // filled with noise (different for EMC, CPV and PPSD) and
142 // then adds contributions from SDigits.
143 // This design avoids scanning over the list of digits to add
144 // contribution to new SDigits only.
990119d6 145
7b7c1533 146 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
3f81a70b 147 TClonesArray * digits = gime->Digits(GetName()) ;
7b7c1533 148
149 digits->Clear() ;
990119d6 150
7b7c1533 151 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
990119d6 152
153 //Making digits with noise, first EMC
154 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
155
156 Int_t nCPV ;
990119d6 157 Int_t absID ;
158 TString name = geom->GetName() ;
159
9688c1dd 160 nCPV = nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
161 geom->GetNModules() ;
8cb3533f 162
9688c1dd 163 digits->Expand(nCPV) ;
8cb3533f 164
7b7c1533 165 // get first the sdigitizer from the tasks list (must have same name as the digitizer)
166 const AliPHOSSDigitizer * sDigitizer = gime->SDigitizer(GetName());
167 if ( !sDigitizer) {
168 cerr << "ERROR: AliPHOSDigitizer::Digitize -> SDigitizer with name " << GetName() << " not found " << endl ;
169 abort() ;
170 }
9688c1dd 171
7b7c1533 172 // loop through the sdigits posted to the White Board and add them to the noise
7a9d98f9 173 TCollection * folderslist = gime->SDigitsFolder()->GetListOfFolders() ;
7b7c1533 174 TIter next(folderslist) ;
175 TFolder * folder = 0 ;
3f81a70b 176 TClonesArray * sdigits = 0 ;
177 Int_t input = 0 ;
9688c1dd 178 TObjArray * sdigArray = new TObjArray(2) ;
179 while ( (folder = (TFolder*)next()) )
3f81a70b 180 if ( (sdigits = (TClonesArray*)folder->FindObject(GetName()) ) ) {
54b82aa4 181 TString fileName(folder->GetName()) ;
182 fileName.ReplaceAll("_","/") ;
3f81a70b 183 cout << "INFO: AliPHOSDigitizer::Digitize -> Adding SDigits "
54b82aa4 184 << GetName() << " from " << fileName << endl ;
9688c1dd 185 sdigArray->AddAt(sdigits, input) ;
186 input++ ;
187 }
188
7a9d98f9 189 //Find the first crystall with signal
9688c1dd 190 Int_t nextSig = 200000 ;
191 Int_t i;
192 for(i=0; i<input; i++){
193 sdigits = (TClonesArray *)sdigArray->At(i) ;
a6eedfad 194 if ( !sdigits->GetEntriesFast() )
7a9d98f9 195 continue ;
9688c1dd 196 Int_t curNext = ((AliPHOSDigit *)sdigits->At(0))->GetId() ;
7a9d98f9 197 if(curNext < nextSig)
198 nextSig = curNext ;
9688c1dd 199 }
7437a0f7 200
9688c1dd 201 TArrayI index(input) ;
202 index.Reset() ; //Set all indexes to zero
203
204 AliPHOSDigit * digit ;
205 AliPHOSDigit * curSDigit ;
206
7437a0f7 207 TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
208
9688c1dd 209 //Put Noise contribution
210 for(absID = 1; absID <= nEMC; absID++){
211 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
212 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
213 //look if we have to add signal?
214 if(absID==nextSig){
215 //Add SDigits from all inputs
216 digit = (AliPHOSDigit *) digits->At(absID-1) ;
7437a0f7 217
218 ticks->Clear() ;
9688c1dd 219 Int_t contrib = 0 ;
7437a0f7 220 Float_t a = digit->GetAmp() ;
221 Float_t b = TMath::Abs( a /fTimeSignalLength) ;
222 //Mark the beginnign of the signal
223 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
224 //Mark the end of the ignal
225 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
226
9688c1dd 227 //loop over inputs
228 for(i=0; i<input; i++){
5c9dfca0 229 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
230 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
231 else
232 curSDigit = 0 ;
9688c1dd 233 //May be several digits will contribute from the same input
234 while(curSDigit && curSDigit->GetId() == absID){
235 //Shift primary to separate primaries belonging different inputs
236 Int_t primaryoffset ;
9891b76e 237 if(fManager)
238 primaryoffset = fManager->GetMask(i) ;
9688c1dd 239 else
21c293b7 240 primaryoffset = 10000000*i ;
241 curSDigit->ShiftPrimary(primaryoffset) ;
7437a0f7 242
243 a = curSDigit->GetAmp() ;
244 b = a /fTimeSignalLength ;
245 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
246 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
9688c1dd 247
248 *digit = *digit + *curSDigit ; //add energies
7437a0f7 249
9688c1dd 250 index[i]++ ;
5c9dfca0 251 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
252 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
253 else
254 curSDigit = 0 ;
9688c1dd 255 }
256 }
257
258 //calculate and set time
7437a0f7 259 Float_t time = FrontEdgeTime(ticks) ;
9688c1dd 260 digit->SetTime(time) ;
7437a0f7 261
9688c1dd 262 //Find next signal module
7437a0f7 263 nextSig = 200000 ;
9688c1dd 264 for(i=0; i<input; i++){
265 sdigits = ((TClonesArray *)sdigArray->At(i)) ;
5c9dfca0 266 Int_t curNext = nextSig ;
267 if(sdigits->GetEntriesFast() > index[i] ){
268 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
269
270 }
9688c1dd 271 if(curNext < nextSig) nextSig = curNext ;
7b7c1533 272 }
273 }
990119d6 274 }
3f81a70b 275
7437a0f7 276 ticks->Delete() ;
277 delete ticks ;
9688c1dd 278
279
280 //Now CPV digits (different noise and no timing)
281 for(absID = nEMC+1; absID <= nCPV; absID++){
282 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
283 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
284 //look if we have to add signal?
285 if(absID==nextSig){
286 digit = (AliPHOSDigit *) digits->At(absID-1) ;
287 //Add SDigits from all inputs
288 for(i=0; i<input; i++){
5c9dfca0 289 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
290 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
291 else
292 curSDigit = 0 ;
a6eedfad 293
9688c1dd 294 //May be several digits will contribute from the same input
295 while(curSDigit && curSDigit->GetId() == absID){
296 //Shift primary to separate primaries belonging different inputs
297 Int_t primaryoffset ;
9891b76e 298 if(fManager)
299 primaryoffset = fManager->GetMask(i) ;
9688c1dd 300 else
301 primaryoffset = 10000000*i ;
302 curSDigit->ShiftPrimary(primaryoffset) ;
303
304 //add energies
305 *digit = *digit + *curSDigit ;
306 index[i]++ ;
5c9dfca0 307 if(((TClonesArray *)sdigArray->At(i))->GetEntriesFast() > index[i] )
308 curSDigit = (AliPHOSDigit*)((TClonesArray *)sdigArray->At(i))->At(index[i]) ;
309 else
310 curSDigit = 0 ;
9688c1dd 311 }
312 }
a6eedfad 313
9688c1dd 314 //Find next signal module
a6eedfad 315 nextSig = 200000 ;
9688c1dd 316 for(i=0; i<input; i++){
317 sdigits = (TClonesArray *)sdigArray->At(i) ;
5c9dfca0 318 Int_t curNext = nextSig ;
319 if(sdigits->GetEntriesFast() > index[i] )
320 curNext = ((AliPHOSDigit *) sdigits->At(index[i]))->GetId() ;
9688c1dd 321 if(curNext < nextSig) nextSig = curNext ;
322 }
323
324 }
325 }
9688c1dd 326 delete sdigArray ; //We should not delete its contents
327
328
329
990119d6 330 //remove digits below thresholds
a6eedfad 331 for(i = 0; i < nEMC ; i++){
332 digit = (AliPHOSDigit*) digits->At(i) ;
aaf8a71c 333 if(sDigitizer->Calibrate( digit->GetAmp() ) < fEMCDigitThreshold)
a6eedfad 334 digits->RemoveAt(i) ;
aaf8a71c 335 else
336 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
337 }
338
a6eedfad 339
340 for(i = nEMC; i < nCPV ; i++)
341 if(sDigitizer->Calibrate(((AliPHOSDigit*)digits->At(i))->GetAmp()) < fCPVDigitThreshold)
342 digits->RemoveAt(i) ;
9688c1dd 343
7b7c1533 344 digits->Compress() ;
990119d6 345
7b7c1533 346 Int_t ndigits = digits->GetEntriesFast() ;
7b7c1533 347 digits->Expand(ndigits) ;
990119d6 348
3758d9fc 349 //Set indexes in list of digits and make true digitization of the energy
990119d6 350 for (i = 0 ; i < ndigits ; i++) {
7b7c1533 351 AliPHOSDigit * digit = (AliPHOSDigit *) digits->At(i) ;
990119d6 352 digit->SetIndexInList(i) ;
f672adc8 353 Float_t energy = sDigitizer->Calibrate(digit->GetAmp()) ;
3758d9fc 354 digit->SetAmp(DigitizeEnergy(energy,digit->GetId()) ) ;
990119d6 355 }
990119d6 356
7b7c1533 357}
3758d9fc 358//____________________________________________________________________________
359Int_t AliPHOSDigitizer::DigitizeEnergy(Float_t energy, Int_t absId)
360{
361 Int_t chanel ;
362 if(absId <= fEmcCrystals){ //digitize as EMC
363 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalEmc)/fADCchanelEmc) ;
364 if(chanel > fNADCemc ) chanel = fNADCemc ;
365 }
366 else{ //Digitize as CPV
367 chanel = (Int_t) TMath::Ceil((energy - fADCpedestalCpv)/fADCchanelCpv) ;
368 if(chanel > fNADCcpv ) chanel = fNADCcpv ;
369 }
370 return chanel ;
371}
7b7c1533 372//____________________________________________________________________________
373void AliPHOSDigitizer::Exec(Option_t *option)
374{
375 // Managing method
990119d6 376
f672adc8 377 if(strcmp(GetName(), "") == 0 )
7b7c1533 378 Init() ;
8cb3533f 379
7b7c1533 380 if (strstr(option,"print")) {
381 Print("");
382 return ;
8cb3533f 383 }
990119d6 384
7b7c1533 385 if(strstr(option,"tim"))
386 gBenchmark->Start("PHOSDigitizer");
8cb3533f 387
3f81a70b 388 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
389
390 Int_t nevents ;
391
392 TTree * treeD ;
393
9891b76e 394 if(fManager){
395 treeD = fManager->GetTreeD() ;
3f81a70b 396 nevents = 1 ; // Will process only one event
397 }
398 else {
399 gAlice->GetEvent(0) ;
400 nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
401 treeD=gAlice->TreeD() ;
402 }
403 if(treeD == 0 ){
404 cerr << " AliPHOSDigitizer :: Can not find TreeD " << endl ;
405 return ;
406 }
407
7b7c1533 408 //Check, if this branch already exits
3f81a70b 409 TObjArray * lob = (TObjArray*)treeD->GetListOfBranches() ;
7b7c1533 410 TIter next(lob) ;
411 TBranch * branch = 0 ;
412 Bool_t phosfound = kFALSE, digitizerfound = kFALSE ;
8cb3533f 413
7b7c1533 414 while ( (branch = (TBranch*)next()) && (!phosfound || !digitizerfound) ) {
3f81a70b 415 if ( (strcmp(branch->GetName(), "PHOS")==0) &&
416 (strcmp(branch->GetTitle(), GetName())==0) )
7b7c1533 417 phosfound = kTRUE ;
8cb3533f 418
3f81a70b 419 else if ( (strcmp(branch->GetName(), "AliPHOSDigitizer")==0) &&
420 (strcmp(branch->GetTitle(), GetName())==0) )
7b7c1533 421 digitizerfound = kTRUE ;
8cb3533f 422 }
990119d6 423
3f81a70b 424 if ( phosfound ) {
425 cerr << "WARNING: AliPHOSDigitizer -> Digits branch with name " << GetName()
426 << " already exits" << endl ;
427 return ;
428 }
429 if ( digitizerfound ) {
430 cerr << "WARNING: AliPHOSDigitizer -> Digitizer branch with name " << GetName()
7b7c1533 431 << " already exits" << endl ;
432 return ;
433 }
eec3ac52 434
7b7c1533 435 Int_t ievent ;
bca3b32a 436
7b7c1533 437 for(ievent = 0; ievent < nevents; ievent++){
3f81a70b 438
9891b76e 439 if(fManager){
3f81a70b 440 Int_t input ;
9891b76e 441 for(input = 0 ; input < fManager->GetNinputs(); input ++){
442 TTree * treeS = fManager->GetInputTreeS(input) ;
3f81a70b 443 if(!treeS){
444 cerr << "AliPHOSDigitizer -> No Input " << endl ;
445 return ;
446 }
447 gime->ReadTreeS(treeS,input) ;
448 }
449 }
450 else
451 gime->Event(ievent,"S") ;
990119d6 452
7b7c1533 453 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
3f81a70b 454
7b7c1533 455 WriteDigits(ievent) ;
3f81a70b 456
01a599c9 457 if(strstr(option,"deb"))
458 PrintDigits(option);
94de8339 459
460 //increment the total number of Digits per run
461 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
01a599c9 462 }
3f81a70b 463
8cb3533f 464 if(strstr(option,"tim")){
465 gBenchmark->Stop("PHOSDigitizer");
466 cout << "AliPHOSDigitizer:" << endl ;
ba54256b 467 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for Digitizing "
468 << gBenchmark->GetCpuTime("PHOSDigitizer")/nevents << " seconds per event " << endl ;
8cb3533f 469 cout << endl ;
470 }
990119d6 471
472}
473
9688c1dd 474//____________________________________________________________________________
7437a0f7 475Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks)
9688c1dd 476{ //
7437a0f7 477 ticks->Sort() ; //Sort in accordance with times of ticks
478 TIter it(ticks) ;
479 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
480 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
481
482 AliPHOSTick * t ;
483 while((t=(AliPHOSTick*) it.Next())){
484 if(t->GetTime() < time) //This tick starts before crossing
485 *ctick+=*t ;
486 else
487 return time ;
488
489 time = ctick->CrossingTime(fTimeThreshold) ;
9688c1dd 490 }
491 return time ;
9688c1dd 492}
7b7c1533 493//____________________________________________________________________________
3f81a70b 494Bool_t AliPHOSDigitizer::Init()
a4e98857 495{
3758d9fc 496 fPinNoise = 0.01 ;
497 fEMCDigitThreshold = 0.01 ;
498 fCPVNoise = 0.01;
499 fCPVDigitThreshold = 0.09 ;
aaf8a71c 500 fTimeResolution = 0.5e-9 ;
3758d9fc 501 fTimeSignalLength = 1.0e-9 ;
502 fDigitsInRun = 0 ;
503 fADCchanelEmc = 0.0015; // width of one ADC channel in GeV
504 fADCpedestalEmc = 0.005 ; //
505 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
506
507 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
508 fADCpedestalCpv = 0.012 ; //
509 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
510
511 fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
8cb3533f 512
9891b76e 513 if(fManager)
f672adc8 514 SetTitle("aliroot") ;
4a72b8de 515 else if (strcmp(GetTitle(),"")==0)
53c09043 516 SetTitle("galice.root") ;
f672adc8 517
518 if( strcmp(GetName(), "") == 0 )
519 SetName("Default") ;
7b7c1533 520
7b7c1533 521 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;
522 if ( gime == 0 ) {
523 cerr << "ERROR: AliPHOSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
3f81a70b 524 return kFALSE;
7b7c1533 525 }
3758d9fc 526
527 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
528 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
93aeb834 529
aad24ee7 530 // Post Digits to the white board
3f81a70b 531 gime->PostDigits(GetName() ) ;
7b7c1533 532
aad24ee7 533 // Post Digitizer to the white board
93aeb834 534 gime->PostDigitizer(this) ;
3f81a70b 535
536 //Mark that we will use current header file
9891b76e 537 if(!fManager){
3f81a70b 538 gime->PostSDigits(GetName(),GetTitle()) ;
539 gime->PostSDigitizer(GetName(),GetTitle()) ;
540 }
541 return kTRUE ;
990119d6 542}
7b7c1533 543
990119d6 544//__________________________________________________________________
7b7c1533 545void AliPHOSDigitizer::MixWith(const char* headerFile)
a4e98857 546{
ba54256b 547 // Allows to produce digits by superimposing background and signal event.
bca3b32a 548 // It is assumed, that headers file with SIGNAL events is opened in
a4e98857 549 // the constructor.
550 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
551 // Thus we avoid writing (changing) huge and expensive
bca3b32a 552 // backgound files: all output will be writen into SIGNAL, i.e.
553 // opened in constructor file.
554 //
a4e98857 555 // One can open as many files to mix with as one needs.
7b7c1533 556 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
557 // can be mixed.
bca3b32a 558
7b7c1533 559 if( strcmp(GetName(), "") == 0 )
8cb3533f 560 Init() ;
561
9891b76e 562 if(fManager){
3f81a70b 563 cout << "Can not use this method under AliRunDigitizer " << endl ;
564 return ;
565 }
7b7c1533 566
567 // check if the specified SDigits do not already exist on the White Board:
aad24ee7 568 // //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
990119d6 569
aad24ee7 570 TString path = "Folders/RunMC/Event/Data/PHOS/SDigits" ;
7b7c1533 571 path += headerFile ;
572 path += "/" ;
3f81a70b 573 path += GetName() ;
7b7c1533 574 if ( gROOT->FindObjectAny(path.Data()) ) {
575 cerr << "WARNING: AliPHOSDigitizer::MixWith -> Entry already exists, do not add" << endl ;
576 return;
990119d6 577 }
3f81a70b 578
579 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
580 gime->PostSDigits(GetName(),headerFile) ;
581
7b7c1533 582 // check if the requested file is already open or exist and if SDigits Branch exist
583 TFile * file = (TFile*)gROOT->FindObject(headerFile);
584 if ( !file ) {
585 file = new TFile(headerFile, "READ") ;
586 if (!file) {
587 cerr << "ERROR: AliPHOSDigitizer::MixWith -> File " << headerFile << " does not exist!" << endl ;
588 return ;
589 }
590 }
3f81a70b 591
990119d6 592}
7b7c1533 593
990119d6 594//__________________________________________________________________
dd5c4038 595void AliPHOSDigitizer::Print(Option_t* option)const {
596 // Print Digitizer's parameters
7b7c1533 597 if( strcmp(GetName(), "") != 0 ){
8cb3533f 598
599 cout << "------------------- "<< GetName() << " -------------" << endl ;
600 cout << "Digitizing sDigits from file(s): " <<endl ;
7b7c1533 601
aad24ee7 602 TCollection * folderslist = ((TFolder*)gROOT->FindObjectAny("Folders/RunMC/Event/Data/PHOS/SDigits"))->GetListOfFolders() ;
7b7c1533 603 TIter next(folderslist) ;
604 TFolder * folder = 0 ;
605
606 while ( (folder = (TFolder*)next()) ) {
607 if ( folder->FindObject(GetName()) )
608 cout << "Adding SDigits " << GetName() << " from " << folder->GetName() << endl ;
8cb3533f 609 }
610 cout << endl ;
7b7c1533 611 cout << "Writing digits to " << GetTitle() << endl ;
8cb3533f 612
613 cout << endl ;
614 cout << "With following parameters: " << endl ;
615 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
616 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
617 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
618 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
8cb3533f 619 cout << "---------------------------------------------------" << endl ;
990119d6 620 }
8cb3533f 621 else
622 cout << "AliPHOSDigitizer not initialized " << endl ;
623
624}
9688c1dd 625
8cb3533f 626//__________________________________________________________________
dd5c4038 627void AliPHOSDigitizer::PrintDigits(Option_t * option){
628 // Print a table of digits
a4e98857 629
7b7c1533 630 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
631 TClonesArray * digits = gime->Digits() ;
632
01a599c9 633 cout << "AliPHOSDigitiser: event " << gAlice->GetEvNumber() << endl ;
7b7c1533 634 cout << " Number of entries in Digits list " << digits->GetEntriesFast() << endl ;
990119d6 635 cout << endl ;
9688c1dd 636 if(strstr(option,"all")||strstr(option,"EMC")){
8cb3533f 637
638 //loop over digits
639 AliPHOSDigit * digit;
a6eedfad 640 cout << "EMC digits (with primaries): " << endl ;
641 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
9688c1dd 642 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
8cb3533f 643 Int_t index ;
a6eedfad 644 for (index = 0 ; (index < digits->GetEntriesFast()) &&
9688c1dd 645 (((AliPHOSDigit * ) digits->At(index))->GetId() <= maxEmc) ; index++) {
7b7c1533 646 digit = (AliPHOSDigit * ) digits->At(index) ;
9688c1dd 647 if(digit->GetNprimary() == 0) continue;
648 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
649 << setw(6) << digit->GetIndexInList() << " "
650 << setw(5) << digit->GetNprimary() <<" ";
8cb3533f 651
652 Int_t iprimary;
653 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
9688c1dd 654 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
8cb3533f 655 cout << endl;
9688c1dd 656 }
657 cout << endl;
658 }
659
660 if(strstr(option,"all")||strstr(option,"CPV")){
8cb3533f 661
9688c1dd 662 //loop over CPV digits
663 AliPHOSDigit * digit;
a6eedfad 664 cout << "CPV digits: " << endl ;
665 cout << "Digit Id Amplitude Index Nprim Primaries list " << endl;
9688c1dd 666 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
667 Int_t index ;
a6eedfad 668 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
9688c1dd 669 digit = (AliPHOSDigit * ) digits->At(index) ;
670 if(digit->GetId() > maxEmc){
671 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " "
672 << setw(6) << digit->GetIndexInList() << " "
673 << setw(5) << digit->GetNprimary() <<" ";
674
675 Int_t iprimary;
676 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
677 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
678 cout << endl;
679 }
680 }
8cb3533f 681 }
9688c1dd 682
8cb3533f 683}
7b7c1533 684
8cb3533f 685//__________________________________________________________________
a4e98857 686void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
687{
bca3b32a 688 // we set title (comment) of the SDigits branch in the first! header file
7b7c1533 689 if( strcmp(GetName(), "") == 0 )
690 Init() ;
990119d6 691
7b7c1533 692 AliPHOSGetter::GetInstance()->SDigits()->SetName(title) ;
693
9688c1dd 694}
695//__________________________________________________________________
696Float_t AliPHOSDigitizer::TimeOfNoise(void)
697{ // Calculates the time signal generated by noise
698 //to be rewritten, now returns just big number
699 return 1. ;
700
8cb3533f 701}
7b7c1533 702//____________________________________________________________________________
703void AliPHOSDigitizer::Reset()
704{
705 // sets current event number to the first simulated event
706
707 if( strcmp(GetName(), "") == 0 )
708 Init() ;
709
710 // Int_t inputs ;
711// for(inputs = 0; inputs < fNinputs ;inputs++)
712// fIevent->AddAt(-1, inputs ) ;
713
714}
715
716//____________________________________________________________________________
717void AliPHOSDigitizer::WriteDigits(Int_t event)
718{
719
720 // Makes TreeD in the output file.
721 // Check if branch already exists:
722 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
723 // already existing branches.
724 // else creates branch with Digits, named "PHOS", title "...",
725 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
726 // and names of files, from which digits are made.
727
728 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
aad24ee7 729 const TClonesArray * digits = gime->Digits(GetName()) ;
730 TTree * treeD ;
7b7c1533 731
9891b76e 732 if(fManager)
733 treeD = fManager->GetTreeD() ;
3f81a70b 734 else
735 treeD = gAlice->TreeD();
736
7b7c1533 737 // create new branches
738 // -- generate file name if necessary
739 char * file =0;
740 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
741 file = new char[strlen(gAlice->GetBaseFile())+20] ;
742 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
743 }
744
745 TDirectory *cwd = gDirectory;
7b7c1533 746 // -- create Digits branch
747 Int_t bufferSize = 32000 ;
ba54256b 748 TBranch * digitsBranch = treeD->Branch("PHOS",&digits,bufferSize);
7b7c1533 749 digitsBranch->SetTitle(GetName());
750 if (file) {
751 digitsBranch->SetFile(file);
752 TIter next( digitsBranch->GetListOfBranches());
753 TBranch * sbr ;
754 while ((sbr=(TBranch*)next())) {
755 sbr->SetFile(file);
756 }
757 cwd->cd();
758 }
759
760 // -- Create Digitizer branch
761 Int_t splitlevel = 0 ;
762 AliPHOSDigitizer * d = gime->Digitizer(GetName()) ;
ba54256b 763 TBranch * digitizerBranch = treeD->Branch("AliPHOSDigitizer", "AliPHOSDigitizer", &d,bufferSize,splitlevel);
7b7c1533 764 digitizerBranch->SetTitle(GetName());
765 if (file) {
766 digitizerBranch->SetFile(file);
767 TIter next( digitizerBranch->GetListOfBranches());
768 TBranch * sbr;
769 while ((sbr=(TBranch*)next())) {
770 sbr->SetFile(file);
771 }
772 cwd->cd();
773 }
3f81a70b 774
6b3966b7 775 treeD->Fill() ;
ba54256b 776 treeD->Write(0,kOverwrite) ;
3f81a70b 777
7b7c1533 778}