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 /* History of cvs commits:
20 * $Log: AliPHOSDigitizer.cxx,v $
21 * Revision 1.104 2007/12/18 09:08:18 hristov
22 * Splitting of the QA maker into simulation and reconstruction dependent parts (Yves)
24 * Revision 1.103 2007/11/07 11:25:06 schutz
25 * Comment out the QA checking before starting digitization
27 * Revision 1.102 2007/10/19 18:04:29 schutz
28 * The standalone QA data maker is called from AliSimulation and AliReconstruction outside the event loop; i.e. re-reading the data. The QA data making in the event loop has been commented out.
30 * Revision 1.101 2007/10/14 21:08:10 schutz
31 * Introduced the checking of QA results from previous step before entering the event loop
33 * Revision 1.100 2007/10/10 09:05:10 schutz
34 * Changing name QualAss to QA
36 * Revision 1.99 2007/09/30 17:08:20 schutz
37 * Introducing the notion of QA data acquisition cycle (needed by online)
39 * Revision 1.98 2007/09/26 14:22:17 cvetan
40 * Important changes to the reconstructor classes. Complete elimination of the run-loaders, which are now steered only from AliReconstruction. Removal of the corresponding Reconstruct() and FillESD() methods.
42 * Revision 1.97 2007/08/07 14:12:03 kharlov
43 * Quality assurance added (Yves Schutz)
45 * Revision 1.96 2007/04/28 10:43:36 policheh
46 * Dead channels simulation: digit energy sets to 0.
48 * Revision 1.95 2007/04/10 07:20:52 kharlov
49 * Decalibration should use the same CDB as calibration in AliPHOSClusterizerv1
51 * Revision 1.94 2007/02/01 10:34:47 hristov
52 * Removing warnings on Solaris x86
54 * Revision 1.93 2006/10/17 13:17:01 kharlov
55 * Replace AliInfo by AliDebug
57 * Revision 1.92 2006/08/28 10:01:56 kharlov
58 * Effective C++ warnings fixed (Timur Pocheptsov)
60 * Revision 1.91 2006/04/29 20:25:30 hristov
61 * Decalibration is implemented (Yu.Kharlov)
63 * Revision 1.90 2006/04/22 10:30:17 hristov
64 * Add fEnergy to AliPHOSDigit and operate with EMC amplitude in energy units (Yu.Kharlov)
66 * Revision 1.89 2006/04/11 15:22:59 hristov
67 * run number in query set to -1: forces AliCDBManager to use its run number (A.Colla)
69 * Revision 1.88 2006/03/13 14:05:43 kharlov
70 * Calibration objects for EMC and CPV
72 * Revision 1.87 2005/08/24 15:33:49 kharlov
73 * Calibration data for raw digits
75 * Revision 1.86 2005/07/12 20:07:35 hristov
76 * Changes needed to run simulation and reconstrruction in the same AliRoot session
78 * Revision 1.85 2005/05/28 14:19:04 schutz
79 * Compilation warnings fixed by T.P.
83 //_________________________________________________________________________
84 //*-- Author : Dmitri Peressounko (SUBATECH & Kurchatov Institute)
85 //////////////////////////////////////////////////////////////////////////////
86 // This TTask performs digitization of Summable digits (in the PHOS case it is just
87 // the sum of contributions from all primary particles into a given cell).
88 // In addition it performs mixing of summable digits from different events.
89 // The name of the TTask is also the title of the branch that will contain
90 // the created SDigits
91 // The title of the TTAsk is the name of the file that contains the hits from
92 // which the SDigits are created
94 // For each event two branches are created in TreeD:
95 // "PHOS" - list of digits
96 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
98 // Note, that one can set a title for new digits branch, and repeat digitization with
99 // another set of parameters.
102 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
103 // root[1] d->ExecuteTask()
104 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
105 // //Digitizes SDigitis in all events found in file galice.root
107 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
108 // // Will read sdigits from galice1.root
109 // root[3] d1->MixWith("galice2.root")
110 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
111 // // Reads another set of sdigits from galice2.root
112 // root[3] d1->MixWith("galice3.root")
113 // // Reads another set of sdigits from galice3.root
114 // root[4] d->ExecuteTask("deb timing")
115 // // Reads SDigits from files galice1.root, galice2.root ....
116 // // mixes them and stores produced Digits in file galice1.root
117 // // deb - prints number of produced digits
118 // // deb all - prints list of produced digits
119 // // timing - prints time used for digitization
122 // --- ROOT system ---
125 #include "TBenchmark.h"
129 // --- Standard library ---
131 // --- AliRoot header files ---
132 #include <TGeoManager.h>
134 #include "AliRunDigitizer.h"
135 #include "AliPHOSDigit.h"
136 #include "AliPHOSDigitizer.h"
137 #include "AliPHOSGeometry.h"
138 #include "AliPHOSTick.h"
139 #include "AliPHOSSimParam.h"
140 #include "AliPHOSCalibData.h"
141 #include "AliRunLoader.h"
142 #include "AliPHOSLoader.h"
144 ClassImp(AliPHOSDigitizer)
147 //____________________________________________________________________________
148 AliPHOSDigitizer::AliPHOSDigitizer() :
154 fInputFileNames(0x0),
157 fEventFolderName(""),
165 fManager = 0 ; // We work in the standalong mode
168 //____________________________________________________________________________
169 AliPHOSDigitizer::AliPHOSDigitizer(TString alirunFileName,
170 TString eventFolderName):
171 AliDigitizer("PHOS"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
172 fDefaultInit(kFALSE),
176 fInputFileNames(0x0),
179 fEventFolderName(eventFolderName),
188 fDefaultInit = kFALSE ;
189 fManager = 0 ; // We work in the standalong mode
190 fcdb = new AliPHOSCalibData(-1);
193 //____________________________________________________________________________
194 AliPHOSDigitizer::AliPHOSDigitizer(const AliPHOSDigitizer & d) :
196 fDefaultInit(d.fDefaultInit),
197 fDigitsInRun(d.fDigitsInRun),
200 fInputFileNames(0x0),//?
202 fEmcCrystals(d.fEmcCrystals),
203 fEventFolderName(d.fEventFolderName),
204 fFirstEvent(d.fFirstEvent),
205 fLastEvent(d.fLastEvent),
210 SetName(d.GetName()) ;
211 SetTitle(d.GetTitle()) ;
212 fcdb = new AliPHOSCalibData(-1);
215 //____________________________________________________________________________
216 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * rd) :
217 AliDigitizer(rd,"PHOS"+AliConfig::Instance()->GetDigitizerTaskName()),
218 fDefaultInit(kFALSE),
222 fInputFileNames(0x0),
225 fEventFolderName(fManager->GetInputFolderName(0)),
232 // ctor Init() is called by RunDigitizer
234 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
236 fDefaultInit = kFALSE ;
237 fcdb = new AliPHOSCalibData(-1);
240 //____________________________________________________________________________
241 AliPHOSDigitizer::~AliPHOSDigitizer()
244 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
246 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
247 phosLoader->CleanDigitizer() ;
250 delete [] fInputFileNames ;
251 delete [] fEventNames ;
253 if(fcdb){ delete fcdb ; fcdb=0;}
257 //____________________________________________________________________________
258 void AliPHOSDigitizer::Digitize(Int_t event)
261 // Makes the digitization of the collected summable digits.
262 // It first creates the array of all PHOS modules
263 // filled with noise (different for EMC, and CPV) and
264 // then adds contributions from SDigits.
265 // This design avoids scanning over the list of digits to add
266 // contribution to new SDigits only.
269 Bool_t toMakeNoise = kTRUE ; //Do not create noisy digits if merge with real data
272 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
273 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
275 Int_t readEvent = event ;
277 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
278 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
279 readEvent, GetTitle(), fEventFolderName.Data())) ;
280 rl->GetEvent(readEvent) ;
281 phosLoader->CleanSDigits() ;
282 phosLoader->LoadSDigits("READ") ;
285 TClonesArray * digits = phosLoader->Digits() ;
287 phosLoader->MakeDigitsArray() ;
288 digits = phosLoader->Digits() ;
294 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
295 //Making digits with noise, first EMC
296 //Check which PHOS modules are present
297 Bool_t isPresent[5] ;
300 for(Int_t i=0; i<5; i++){
301 volpath = "/ALIC_1/PHOS_";
303 if (gGeoManager->CheckPath(volpath.Data())) {
312 Int_t nEMC = nmod*geom->GetNPhi()*geom->GetNZ();
317 //check if CPV exists
318 Bool_t isCPVpresent=0 ;
319 for(Int_t i=1; i<=5 && !isCPVpresent; i++){
320 volpath = "/ALIC_1/PHOS_";
322 volpath += "/PCPV_1";
323 if (gGeoManager->CheckPath(volpath.Data()))
328 nCPV = nEMC + geom->GetNumberOfCPVPadsZ() * geom->GetNumberOfCPVPadsPhi() * nmod ;
334 digits->Expand(nCPV) ;
336 //take all the inputs to add together and load the SDigits
337 TObjArray * sdigArray = new TObjArray(fInput) ;
338 sdigArray->AddAt(phosLoader->SDigits(), 0) ;
340 for(Int_t i = 1 ; i < fInput ; i++){
341 TString tempo(fEventNames[i]) ;
343 AliRunLoader* rl2 = AliRunLoader::GetRunLoader(tempo) ;
345 rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ;
347 Fatal("AliPHOSDigitizer", "Could not find the Run Loader for %s - %s",fInputFileNames[i].Data(), tempo.Data()) ;
352 AliPHOSLoader * phosLoader2 = dynamic_cast<AliPHOSLoader*>(rl2->GetLoader("PHOSLoader"));
355 readEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
357 TClonesArray * digs ;
358 if(AliPHOSSimParam::GetInstance()->IsStreamDigits(i)){ //This is Digits Stream
359 AliInfo(Form("Adding event %d from input stream %d %s %s",
360 readEvent, i, fInputFileNames[i].Data(), tempo.Data())) ;
361 rl2->GetEvent(readEvent) ;
362 phosLoader2->CleanDigits() ;
363 phosLoader2->LoadDigits("READ") ;
364 digs = phosLoader2->Digits() ;
365 toMakeNoise=0 ; //Do not add noise, it is already in stream
368 AliInfo(Form("Adding event %d (SDigits) from input stream %d %s %s",
369 readEvent, i, fInputFileNames[i].Data(), tempo.Data())) ;
370 rl2->GetEvent(readEvent) ;
371 phosLoader2->CleanSDigits() ;
372 phosLoader2->LoadSDigits("READ") ;
373 digs = phosLoader2->SDigits() ;
375 sdigArray->AddAt(digs, i) ;
378 //Find the first crystal with signal
379 Int_t nextSig = 200000 ;
380 TClonesArray * sdigits ;
381 for(Int_t i = 0 ; i < fInput ; i++){
382 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
383 if ( !sdigits->GetEntriesFast() )
385 Int_t curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(0))->GetId() ;
386 if(curNext < nextSig)
390 TArrayI index(fInput) ;
391 index.Reset() ; //Set all indexes to zero
393 AliPHOSDigit * digit ;
394 AliPHOSDigit * curSDigit ;
396 // TClonesArray * ticks = new TClonesArray("AliPHOSTick",1000) ;
398 //Put Noise contribution
399 Double_t apdNoise = 0. ;
401 apdNoise = AliPHOSSimParam::GetInstance()->GetAPDNoise() ;
403 Int_t emcpermod=geom->GetNPhi()*geom->GetNZ();
405 for(Int_t imod=0; imod<5; imod++){
408 Int_t firstAbsId=imod*emcpermod+1 ;
409 Int_t lastAbsId =(imod+1)*emcpermod ;
410 for(absID = firstAbsId ; absID <= lastAbsId ; absID++){
411 Float_t noise = gRandom->Gaus(0.,apdNoise) ;
412 new((*digits)[idigit]) AliPHOSDigit( -1, absID, noise, TimeOfNoise() ) ;
413 //look if we have to add signal?
414 digit = dynamic_cast<AliPHOSDigit *>(digits->At(idigit)) ;
418 //Add SDigits from all inputs
420 // Int_t contrib = 0 ;
422 //New Timing model is necessary
423 // Float_t a = digit->GetEnergy() ;
424 // Float_t b = TMath::Abs( a / fTimeSignalLength) ;
425 // //Mark the beginning of the signal
426 // new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
427 // //Mark the end of the signal
428 // new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
430 // Calculate time as time of the largest digit
431 Float_t time = digit->GetTime() ;
432 Float_t eTime= digit->GetEnergy() ;
435 for(Int_t i = 0 ; i < fInput ; i++){
436 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
437 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
440 //May be several digits will contribute from the same input
441 while(curSDigit && curSDigit->GetId() == absID){
442 //Shift primary to separate primaries belonging different inputs
443 Int_t primaryoffset ;
445 primaryoffset = fManager->GetMask(i) ;
447 primaryoffset = 10000000*i ;
448 curSDigit->ShiftPrimary(primaryoffset) ;
450 //New Timing model is necessary
451 // a = curSDigit->GetEnergy() ;
452 // b = a /fTimeSignalLength ;
453 // new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
454 // new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
455 if(curSDigit->GetEnergy()>eTime){
456 eTime=curSDigit->GetEnergy() ;
457 time=curSDigit->GetTime() ;
460 *digit += *curSDigit ; //add energies
463 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
464 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
470 //calculate and set time
471 //New Timing model is necessary
472 // Float_t time = FrontEdgeTime(ticks) ;
473 digit->SetTime(time) ;
475 //Find next signal module
477 for(Int_t i = 0 ; i < fInput ; i++){
478 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
479 Int_t curNext = nextSig ;
480 if(sdigits->GetEntriesFast() > index[i] ){
481 curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(index[i]))->GetId() ;
483 if(curNext < nextSig) nextSig = curNext ;
489 //distretize energy if necessary
490 if(AliPHOSSimParam::GetInstance()->IsEDigitizationOn()){
491 Float_t adcW=AliPHOSSimParam::GetInstance()->GetADCchannelW() ;
492 for(Int_t i = 0 ; i < nEMC ; i++){
493 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
494 digit->SetEnergy(adcW*ceil(digit->GetEnergy()/adcW)) ;
498 //Apply Fast decalibration if necessary
499 if(AliPHOSSimParam::GetInstance()->GetFastDecalibration()>0.){
500 Float_t res=AliPHOSSimParam::GetInstance()->GetFastDecalibration() ;
501 for(Int_t i = 0 ; i < nEMC ; i++){
502 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
503 digit->SetEnergy(gRandom->Gaus(1.,res)*digit->GetEnergy()) ;
511 //Now CPV digits (different noise and no timing)
512 Int_t cpvpermod = geom->GetNumberOfCPVPadsZ() * geom->GetNumberOfCPVPadsPhi() ;
513 Int_t nEMCtotal=emcpermod*5 ;
514 Float_t cpvNoise = AliPHOSSimParam::GetInstance()->GetCPVNoise() ;
515 if(isCPVpresent){ //CPV is present in current geometry
516 for(Int_t imod=0; imod<5; imod++){ //module is present in current geometry
519 Int_t firstAbsId=nEMCtotal+imod*cpvpermod+1 ;
520 Int_t lastAbsId =nEMCtotal+(imod+1)*cpvpermod ;
521 for(absID = firstAbsId; absID <= lastAbsId; absID++){
522 Float_t noise = gRandom->Gaus(0., cpvNoise) ;
523 new((*digits)[idigit]) AliPHOSDigit( -1,absID,noise, TimeOfNoise() ) ;
525 //look if we have to add signal?
527 digit = dynamic_cast<AliPHOSDigit *>(digits->At(idigit-1)) ;
528 //Add SDigits from all inputs
529 for(Int_t i = 0 ; i < fInput ; i++){
530 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
531 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
535 //May be several digits will contribute from the same input
536 while(curSDigit && curSDigit->GetId() == absID){
537 //Shift primary to separate primaries belonging different inputs
538 Int_t primaryoffset ;
540 primaryoffset = fManager->GetMask(i) ;
542 primaryoffset = 10000000*i ;
543 curSDigit->ShiftPrimary(primaryoffset) ;
546 *digit += *curSDigit ;
548 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
549 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i]) ) ;
555 //Find next signal module
557 for(Int_t i = 0 ; i < fInput ; i++){
558 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
559 Int_t curNext = nextSig ;
560 if(sdigits->GetEntriesFast() > index[i] )
561 curNext = dynamic_cast<AliPHOSDigit *>( sdigits->At(index[i]) )->GetId() ;
562 if(curNext < nextSig) nextSig = curNext ;
570 delete sdigArray ; //We should not delete its contents
574 //set amplitudes in bad channels to zero
576 for(Int_t i = 0 ; i <digits->GetEntriesFast(); i++){
577 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
578 geom->AbsToRelNumbering(digit->GetId(),relId);
579 if(relId[1] == 0) // Emc
580 if(fcdb->IsBadChannelEmc(relId[0],relId[3],relId[2])) digit->SetEnergy(0.);
583 //remove digits below thresholds
584 Float_t emcThreshold = AliPHOSSimParam::GetInstance()->GetEmcDigitsThreshold() ;
585 for(Int_t i = 0 ; i < nEMC ; i++){
586 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
588 if(digit->GetEnergy() < emcThreshold){
589 digits->RemoveAt(i) ;
592 Float_t tres = TimeResolution(digit->GetEnergy()) ;
593 digit->SetTime(gRandom->Gaus(digit->GetTime(), tres) ) ;
596 Float_t cpvDigitThreshold = AliPHOSSimParam::GetInstance()->GetCpvDigitsThreshold() ;
597 for(Int_t i = nEMC; i < nCPV ; i++){
598 if( dynamic_cast<AliPHOSDigit*>(digits->At(i))->GetEnergy() < cpvDigitThreshold )
599 digits->RemoveAt(i) ;
604 Int_t ndigits = digits->GetEntriesFast() ;
605 digits->Expand(ndigits) ;
608 //Set indexes in list of digits and make true digitization of the energy
609 for (Int_t i = 0 ; i < ndigits ; i++) {
610 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
611 digit->SetIndexInList(i) ;
612 if(digit->GetId() > fEmcCrystals){ //digitize CPV only
613 digit->SetAmp(DigitizeCPV(digit->GetEnergy(),digit->GetId()) ) ;
619 //____________________________________________________________________________
620 void AliPHOSDigitizer::DecalibrateEMC(AliPHOSDigit *digit)
622 // Decalibrate EMC digit, i.e. change its energy by a factor read from CDB
624 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
626 //Determine rel.position of the cell absolute ID
628 geom->AbsToRelNumbering(digit->GetId(),relId);
629 Int_t module=relId[0];
631 Int_t column=relId[3];
632 Float_t decalibration = fcdb->GetADCchannelEmc(module,column,row);
633 Float_t energy = digit->GetEnergy() / decalibration;
634 digit->SetEnergy(energy);
636 //____________________________________________________________________________
637 Int_t AliPHOSDigitizer::DigitizeCPV(Float_t charge, Int_t absId)
639 // Returns digitized value of the CPV charge in a pad absId
641 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
643 //Determine rel.position of the cell absId
645 geom->AbsToRelNumbering(absId,relId);
646 Int_t module=relId[0];
648 Int_t column=relId[3];
652 if(absId > fEmcCrystals){ //digitize CPV only
654 //reading calibration data for cell absId.
655 Float_t adcPedestalCpv = fcdb->GetADCpedestalCpv(module,column,row);
656 Float_t adcChanelCpv = fcdb->GetADCchannelCpv( module,column,row);
658 channel = (Int_t) TMath::Ceil((charge - adcPedestalCpv)/adcChanelCpv) ;
659 Int_t nMax = AliPHOSSimParam::GetInstance()->GetNADCcpv() ;
660 if(channel > nMax ) channel = nMax ;
665 //____________________________________________________________________________
666 void AliPHOSDigitizer::Exec(Option_t *option)
668 // Steering method to process digitization for events
669 // in the range from fFirstEvent to fLastEvent.
670 // This range is optionally set by SetEventRange().
671 // if fLastEvent=-1, then process events until the end.
672 // by default fLastEvent = fFirstEvent (process only one event)
674 if (!fInit) { // to prevent overwrite existing file
675 AliError(Form("Give a version name different from %s",
676 fEventFolderName.Data() )) ;
680 if (strstr(option,"print")) {
685 if(strstr(option,"tim"))
686 gBenchmark->Start("PHOSDigitizer");
688 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
689 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
691 // Post Digitizer to the white board
692 phosLoader->PostDigitizer(this) ;
694 if (fLastEvent == -1)
695 fLastEvent = rl->GetNumberOfEvents() - 1 ;
697 fLastEvent = fFirstEvent ;
699 Int_t nEvents = fLastEvent - fFirstEvent + 1;
703 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
706 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
710 if(strstr(option,"deb"))
713 //increment the total number of Digits per run
714 fDigitsInRun += phosLoader->Digits()->GetEntriesFast() ;
717 phosLoader->CleanDigitizer();
719 if(strstr(option,"tim")){
720 gBenchmark->Stop("PHOSDigitizer");
722 message = " took %f seconds for Digitizing %f seconds per event\n" ;
723 AliInfo(Form( message.Data(),
724 gBenchmark->GetCpuTime("PHOSDigitizer"),
725 gBenchmark->GetCpuTime("PHOSDigitizer")/nEvents ));
728 //____________________________________________________________________________
729 Float_t AliPHOSDigitizer::TimeResolution(Float_t e){
730 //calculate TOF resolution using beam-test resutls
731 Float_t a=AliPHOSSimParam::GetInstance()->GetTOFa() ;
732 Float_t b=AliPHOSSimParam::GetInstance()->GetTOFb() ;
733 return TMath::Sqrt(a*a+b*b/e) ;
736 ////____________________________________________________________________________
737 //Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) const
739 // // Returns the shortest time among all time ticks
741 // ticks->Sort() ; //Sort in accordance with times of ticks
743 // AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
744 // Float_t time = ctick->CrossingTime(fTimeThreshold) ;
747 // while((t=(AliPHOSTick*) it.Next())){
748 // if(t->GetTime() < time) //This tick starts before crossing
753 // time = ctick->CrossingTime(fTimeThreshold) ;
758 //____________________________________________________________________________
759 Bool_t AliPHOSDigitizer::Init()
761 // Makes all memory allocations
764 AliPHOSGeometry *geom;
765 if (!(geom = AliPHOSGeometry::GetInstance()))
766 geom = AliPHOSGeometry::GetInstance("IHEP","");
767 // const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
769 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
772 fLastEvent = fFirstEvent ;
774 fInput = fManager->GetNinputs() ;
778 fInputFileNames = new TString[fInput] ;
779 fEventNames = new TString[fInput] ;
780 fInputFileNames[0] = GetTitle() ;
781 fEventNames[0] = fEventFolderName.Data() ;
783 for (index = 1 ; index < fInput ; index++) {
784 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
785 TString tempo = fManager->GetInputFolderName(index) ;
786 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added by fManager
789 //to prevent cleaning of this object while GetEvent is called
790 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
792 rl = AliRunLoader::Open(GetTitle(), fEventFolderName) ;
794 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
795 phosLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
800 //____________________________________________________________________________
801 void AliPHOSDigitizer::InitParameters()
803 // Set initial parameters Digitizer
806 SetEventRange(0,-1) ;
810 //__________________________________________________________________
811 void AliPHOSDigitizer::Print(const Option_t *)const
813 // Print Digitizer's parameters
814 AliInfo(Form("\n------------------- %s -------------", GetName() )) ;
815 if( strcmp(fEventFolderName.Data(), "") != 0 ){
816 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
820 nStreams = GetNInputStreams() ;
825 for (index = 0 ; index < nStreams ; index++) {
826 TString tempo(fEventNames[index]) ;
828 printf ("Adding SDigits from %s \n", fInputFileNames[index].Data()) ;
831 // printf("\nWith following parameters:\n") ;
832 // printf(" Electronics noise in EMC (fPinNoise) = %f GeV\n", fPinNoise ) ;
833 // printf(" Threshold in EMC (fEMCDigitThreshold) = %f GeV\n", fEMCDigitThreshold ) ;
834 // printf(" Noise in CPV (fCPVNoise) = %f aux units\n", fCPVNoise ) ;
835 // printf(" Threshold in CPV (fCPVDigitThreshold) = %f aux units\n",fCPVDigitThreshold ) ;
836 printf(" ---------------------------------------------------\n") ;
839 AliInfo(Form("AliPHOSDigitizer not initialized" )) ;
843 //__________________________________________________________________
844 void AliPHOSDigitizer::PrintDigits(Option_t * option)
846 // Print a table of digits
849 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
850 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
851 TClonesArray * digits = phosLoader->Digits() ;
852 const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
854 AliInfo(Form("%d", digits->GetEntriesFast())) ;
855 printf("\nevent %d", gAlice->GetEvNumber()) ;
856 printf("\n Number of entries in Digits list %d", digits->GetEntriesFast() ) ;
859 if(strstr(option,"all")||strstr(option,"EMC")){
861 AliPHOSDigit * digit;
862 printf("\nEMC digits (with primaries):\n") ;
863 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
864 Int_t maxEmc = geom->GetNModules()*geom->GetNCristalsInModule() ;
866 for (index = 0 ; (index < digits->GetEntriesFast()) &&
867 (dynamic_cast<AliPHOSDigit *>(digits->At(index))->GetId() <= maxEmc) ; index++) {
868 digit = (AliPHOSDigit * ) digits->At(index) ;
869 if(digit->GetNprimary() == 0)
871 // printf("%6d %8d %6.5e %4d %2d :",
872 // digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ; // YVK
873 printf("%6d %.4f %6.5e %4d %2d :",
874 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
876 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
877 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
883 if(strstr(option,"all")||strstr(option,"CPV")){
885 //loop over CPV digits
886 AliPHOSDigit * digit;
887 printf("\nCPV digits (with primaries):\n") ;
888 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
889 Int_t maxEmc = geom->GetNModules()*geom->GetNCristalsInModule() ;
891 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
892 digit = (AliPHOSDigit * ) digits->At(index) ;
893 if(digit->GetId() > maxEmc){
894 printf("%6d %8d %4d %2d :",
895 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
897 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
898 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
907 //__________________________________________________________________
908 Float_t AliPHOSDigitizer::TimeOfNoise(void) const
909 { // Calculates the time signal generated by noise
910 //PH Info("TimeOfNoise", "Change me") ;
911 return gRandom->Rndm() * 1.28E-5;
914 //__________________________________________________________________
915 void AliPHOSDigitizer::Unload()
919 for(i = 1 ; i < fInput ; i++){
920 TString tempo(fEventNames[i]) ;
922 AliRunLoader* rl = AliRunLoader::GetRunLoader(tempo) ;
923 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
924 phosLoader->UnloadSDigits() ;
927 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
928 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
929 phosLoader->UnloadDigits() ;
932 //____________________________________________________________________________
933 void AliPHOSDigitizer::WriteDigits()
936 // Makes TreeD in the output file.
937 // Check if branch already exists:
938 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
939 // already existing branches.
940 // else creates branch with Digits, named "PHOS", title "...",
941 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
942 // and names of files, from which digits are made.
944 AliRunLoader* rl = AliRunLoader::GetRunLoader(fEventFolderName) ;
945 AliPHOSLoader * phosLoader = dynamic_cast<AliPHOSLoader*>(rl->GetLoader("PHOSLoader"));
947 const TClonesArray * digits = phosLoader->Digits() ;
948 TTree * treeD = phosLoader->TreeD();
950 phosLoader->MakeTree("D");
951 treeD = phosLoader->TreeD();
954 // -- create Digits branch
955 Int_t bufferSize = 32000 ;
956 TBranch * digitsBranch = treeD->Branch("PHOS","TClonesArray",&digits,bufferSize);
957 digitsBranch->SetTitle(fEventFolderName);
958 digitsBranch->Fill() ;
960 phosLoader->WriteDigits("OVERWRITE");
961 phosLoader->WriteDigitizer("OVERWRITE");