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"
128 // --- Standard library ---
130 // --- AliRoot header files ---
132 #include "AliRunDigitizer.h"
133 #include "AliPHOSDigit.h"
134 #include "AliPHOSGetter.h"
135 #include "AliPHOSDigitizer.h"
136 #include "AliPHOSSDigitizer.h"
137 #include "AliPHOSGeometry.h"
138 #include "AliPHOSTick.h"
140 ClassImp(AliPHOSDigitizer)
143 //____________________________________________________________________________
144 AliPHOSDigitizer::AliPHOSDigitizer() :
150 fInputFileNames(0x0),
154 fEMCDigitThreshold(0.f),
156 fCPVDigitThreshold(0.f),
157 fTimeResolution(0.f),
159 fTimeSignalLength(0.f),
161 fADCpedestalEmc(0.f),
164 fADCpedestalCpv(0.f),
166 fEventFolderName(""),
174 fManager = 0 ; // We work in the standalong mode
177 //____________________________________________________________________________
178 AliPHOSDigitizer::AliPHOSDigitizer(TString alirunFileName,
179 TString eventFolderName):
180 AliDigitizer("PHOS"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
181 fDefaultInit(kFALSE),
185 fInputFileNames(0x0),
189 fEMCDigitThreshold(0.f),
191 fCPVDigitThreshold(0.f),
192 fTimeResolution(0.f),
194 fTimeSignalLength(0.f),
196 fADCpedestalEmc(0.f),
199 fADCpedestalCpv(0.f),
201 fEventFolderName(eventFolderName),
210 fDefaultInit = kFALSE ;
211 fManager = 0 ; // We work in the standalong mode
212 // //Initialize the quality assurance data maker only once
213 // fQADM = new AliPHOSQADataMaker() ;
214 // //FIXME: get the run number
217 // GetQADataMaker()->Init(AliQA::kDIGITS, run, fgkCycles) ;
218 // GetQADataMaker()->StartOfCycle(AliQA::kDIGITS) ;
221 //____________________________________________________________________________
222 AliPHOSDigitizer::AliPHOSDigitizer(const AliPHOSDigitizer & d) :
224 fDefaultInit(d.fDefaultInit),
225 fDigitsInRun(d.fDigitsInRun),
228 fInputFileNames(0x0),//?
230 fEmcCrystals(d.fEmcCrystals),
231 fPinNoise(d.fPinNoise),
232 fEMCDigitThreshold(d.fEMCDigitThreshold),
233 fCPVNoise(d.fCPVNoise),
234 fCPVDigitThreshold(d.fCPVDigitThreshold),
235 fTimeResolution(d.fTimeResolution),
236 fTimeThreshold(d.fTimeThreshold),
237 fTimeSignalLength(d.fTimeSignalLength),
238 fADCchanelEmc(d.fADCchanelEmc),
239 fADCpedestalEmc(d.fADCpedestalEmc),
240 fNADCemc(d.fNADCemc),
241 fADCchanelCpv(d.fADCchanelCpv),
242 fADCpedestalCpv(d.fADCpedestalCpv),
243 fNADCcpv(d.fNADCcpv),
244 fEventFolderName(d.fEventFolderName),
245 fFirstEvent(d.fFirstEvent),
246 fLastEvent(d.fLastEvent),
252 SetName(d.GetName()) ;
253 SetTitle(d.GetTitle()) ;
254 //Initialize the quality assurance data maker only once
255 //FIXME: get the run number
258 // GetQADataMaker()->Init(AliQA::kDIGITS, run, fgkCycles) ;
259 // GetQADataMaker()->StartOfCycle(AliQA::kDIGITS) ;
262 //____________________________________________________________________________
263 AliPHOSDigitizer::AliPHOSDigitizer(AliRunDigitizer * rd) :
264 AliDigitizer(rd,"PHOS"+AliConfig::Instance()->GetDigitizerTaskName()),
265 fDefaultInit(kFALSE),
269 fInputFileNames(0x0),
273 fEMCDigitThreshold(0.f),
275 fCPVDigitThreshold(0.f),
276 fTimeResolution(0.f),
278 fTimeSignalLength(0.f),
280 fADCpedestalEmc(0.f),
283 fADCpedestalCpv(0.f),
285 fEventFolderName(fManager->GetInputFolderName(0)),
292 // ctor Init() is called by RunDigitizer
294 SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
296 fDefaultInit = kFALSE ;
297 //Initialize the quality assurance data maker only once
298 // fQADM = new AliPHOSQADataMaker() ;
299 // //FIXME: get the run number
302 // GetQADataMaker()->Init(AliQA::kDIGITS, run) ;
303 // GetQADataMaker()->StartOfCycle(AliQA::kDIGITS) ;
306 //____________________________________________________________________________
307 AliPHOSDigitizer::~AliPHOSDigitizer()
309 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
311 // Clean Digitizer from the white board
312 gime->PhosLoader()->CleanDigitizer() ;
314 delete [] fInputFileNames ;
315 delete [] fEventNames ;
321 //____________________________________________________________________________
322 void AliPHOSDigitizer::Digitize(Int_t event)
325 // Makes the digitization of the collected summable digits.
326 // It first creates the array of all PHOS modules
327 // filled with noise (different for EMC, and CPV) and
328 // then adds contributions from SDigits.
329 // This design avoids scanning over the list of digits to add
330 // contribution to new SDigits only.
332 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
333 Int_t ReadEvent = event ;
335 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ;
336 AliDebug(1,Form("Adding event %d from input stream 0 %s %s",
337 ReadEvent, GetTitle(), fEventFolderName.Data())) ;
338 gime->Event(ReadEvent, "S") ;
339 TClonesArray * digits = gime->Digits() ;
342 const AliPHOSGeometry *geom = gime->PHOSGeometry() ;
343 //Making digits with noise, first EMC
344 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
349 nCPV = nEMC + geom->GetNumberOfCPVPadsZ() * geom->GetNumberOfCPVPadsPhi() * geom->GetNModules() ;
351 digits->Expand(nCPV) ;
353 // get first the sdigitizer from the tasks list
354 if ( !gime->SDigitizer() )
355 gime->LoadSDigitizer();
356 AliPHOSSDigitizer * sDigitizer = gime->SDigitizer();
359 AliFatal(Form("SDigitizer with name %s %s not found",
360 GetTitle(), fEventFolderName.Data() )) ;
362 //take all the inputs to add together and load the SDigits
363 TObjArray * sdigArray = new TObjArray(fInput) ;
364 sdigArray->AddAt(gime->SDigits(), 0) ;
366 for(i = 1 ; i < fInput ; i++){
367 TString tempo(fEventNames[i]) ;
369 AliPHOSGetter * gime1 = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
371 ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ;
372 AliInfo(Form("Adding event %d from input stream %d %s %s",
373 ReadEvent, i, fInputFileNames[i].Data(), tempo.Data())) ;
374 gime1->Event(ReadEvent,"S");
375 sdigArray->AddAt(gime1->SDigits(), i) ;
378 //Find the first crystal with signal
379 Int_t nextSig = 200000 ;
380 TClonesArray * sdigits ;
381 for(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 for(absID = 1 ; absID <= nEMC ; absID++){
400 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
401 // YVK: do not digitize amplitudes for EMC
402 // new((*digits)[absID-1]) AliPHOSDigit( -1, absID, sDigitizer->Digitize(noise), TimeOfNoise() ) ;
403 new((*digits)[absID-1]) AliPHOSDigit( -1, absID, noise, TimeOfNoise() ) ;
404 //look if we have to add signal?
405 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
408 //Add SDigits from all inputs
411 Float_t a = digit->GetEnergy() ;
412 Float_t b = TMath::Abs( a / fTimeSignalLength) ;
413 //Mark the beginning of the signal
414 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime(),0, b);
415 //Mark the end of the signal
416 new((*ticks)[contrib++]) AliPHOSTick(digit->GetTime()+fTimeSignalLength, -a, -b);
419 for(i = 0 ; i < fInput ; i++){
420 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
421 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
424 //May be several digits will contribute from the same input
425 while(curSDigit && curSDigit->GetId() == absID){
426 //Shift primary to separate primaries belonging different inputs
427 Int_t primaryoffset ;
429 primaryoffset = fManager->GetMask(i) ;
431 primaryoffset = 10000000*i ;
432 curSDigit->ShiftPrimary(primaryoffset) ;
434 a = curSDigit->GetEnergy() ;
435 b = a /fTimeSignalLength ;
436 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime(),0, b);
437 new((*ticks)[contrib++]) AliPHOSTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b);
439 *digit += *curSDigit ; //add energies
442 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
443 curSDigit = dynamic_cast<AliPHOSDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
449 //calculate and set time
450 Float_t time = FrontEdgeTime(ticks) ;
451 digit->SetTime(time) ;
453 //Find next signal module
455 for(i = 0 ; i < fInput ; i++){
456 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
457 Int_t curNext = nextSig ;
458 if(sdigits->GetEntriesFast() > index[i] ){
459 curNext = dynamic_cast<AliPHOSDigit *>(sdigits->At(index[i]))->GetId() ;
461 if(curNext < nextSig) nextSig = curNext ;
469 //Now CPV digits (different noise and no timing)
470 for(absID = nEMC+1; absID <= nCPV; absID++){
471 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
472 new((*digits)[absID-1]) AliPHOSDigit( -1,absID,sDigitizer->Digitize(noise), TimeOfNoise() ) ;
473 //look if we have to add signal?
475 digit = dynamic_cast<AliPHOSDigit *>(digits->At(absID-1)) ;
476 //Add SDigits from all inputs
477 for(i = 0 ; i < fInput ; i++){
478 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
479 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
483 //May be several digits will contribute from the same input
484 while(curSDigit && curSDigit->GetId() == absID){
485 //Shift primary to separate primaries belonging different inputs
486 Int_t primaryoffset ;
488 primaryoffset = fManager->GetMask(i) ;
490 primaryoffset = 10000000*i ;
491 curSDigit->ShiftPrimary(primaryoffset) ;
494 *digit += *curSDigit ;
496 if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
497 curSDigit = dynamic_cast<AliPHOSDigit*>( dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i]) ) ;
503 //Find next signal module
505 for(i = 0 ; i < fInput ; i++){
506 sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
507 Int_t curNext = nextSig ;
508 if(sdigits->GetEntriesFast() > index[i] )
509 curNext = dynamic_cast<AliPHOSDigit *>( sdigits->At(index[i]) )->GetId() ;
510 if(curNext < nextSig) nextSig = curNext ;
516 delete sdigArray ; //We should not delete its contents
518 //remove digits below thresholds
519 for(i = 0 ; i < nEMC ; i++){
520 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
521 //By default no decalibration should be applyed
522 // DecalibrateEMC(digit);
524 if(digit->GetEnergy() < fEMCDigitThreshold)
525 digits->RemoveAt(i) ;
527 digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
531 for(i = nEMC; i < nCPV ; i++)
532 // if( sDigitizer->Calibrate( dynamic_cast<AliPHOSDigit*>(digits->At(i))->GetAmp() ) < fCPVDigitThreshold )
533 if( dynamic_cast<AliPHOSDigit*>(digits->At(i))->GetEnergy() < fCPVDigitThreshold )
534 digits->RemoveAt(i) ;
538 Int_t ndigits = digits->GetEntriesFast() ;
539 digits->Expand(ndigits) ;
541 //Set indexes in list of digits and make true digitization of the energy
542 for (i = 0 ; i < ndigits ; i++) {
543 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
544 digit->SetIndexInList(i) ;
545 if(digit->GetId() > fEmcCrystals){ //digitize CPV only
546 digit->SetAmp(DigitizeCPV(digit->GetEnergy(),digit->GetId()) ) ;
552 //set amplitudes in bad channels to zero
553 for(i = 0 ; i <digits->GetEntries(); i++){
554 digit = dynamic_cast<AliPHOSDigit*>( digits->At(i) ) ;
555 gime->PHOSGeometry()->AbsToRelNumbering(digit->GetId(),relId);
556 if(relId[1] == 0) // Emc
557 if(gime->CalibData()->IsBadChannelEmc(relId[0],relId[3],relId[2])) digit->SetEnergy(0.);
562 //____________________________________________________________________________
563 void AliPHOSDigitizer::DecalibrateEMC(AliPHOSDigit *digit)
565 // Decalibrate EMC digit, i.e. change its energy by a factor read from CDB
567 AliPHOSGetter* gime = AliPHOSGetter::Instance();
569 if(!gime->CalibData()) {
570 AliPHOSCalibData* cdb = new AliPHOSCalibData(-1);
571 gime->SetCalibData(cdb);
574 //Determine rel.position of the cell absolute ID
576 gime->PHOSGeometry()->AbsToRelNumbering(digit->GetId(),relId);
577 Int_t module=relId[0];
579 Int_t column=relId[3];
580 Float_t decalibration = gime->CalibData()->GetADCchannelEmc(module,column,row);
581 Float_t energy = digit->GetEnergy() / decalibration;
582 digit->SetEnergy(energy);
584 //____________________________________________________________________________
585 Int_t AliPHOSDigitizer::DigitizeCPV(Float_t charge, Int_t absId)
587 // Returns digitized value of the CPV charge in a pad absId
589 AliPHOSGetter* gime = AliPHOSGetter::Instance();
591 if(!gime->CalibData()) {
592 AliPHOSCalibData* cdb = new AliPHOSCalibData(-1); // use AliCDBManager's run number
593 gime->SetCalibData(cdb);
596 //Determine rel.position of the cell absId
598 gime->PHOSGeometry()->AbsToRelNumbering(absId,relId);
599 Int_t module=relId[0];
601 Int_t column=relId[3];
605 if(absId > fEmcCrystals){ //digitize CPV only
607 //reading calibration data for cell absId.
608 //If no calibration DB found, accept default values.
610 if(gime->CalibData()) {
611 fADCpedestalCpv = gime->CalibData()->GetADCpedestalCpv(module,column,row);
612 fADCchanelCpv = gime->CalibData()->GetADCchannelCpv( module,column,row);
615 channel = (Int_t) TMath::Ceil((charge - fADCpedestalCpv)/fADCchanelCpv) ;
616 if(channel > fNADCcpv ) channel = fNADCcpv ;
621 //____________________________________________________________________________
622 void AliPHOSDigitizer::Exec(Option_t *option)
624 // Steering method to process digitization for events
625 // in the range from fFirstEvent to fLastEvent.
626 // This range is optionally set by SetEventRange().
627 // if fLastEvent=-1, then process events until the end.
628 // by default fLastEvent = fFirstEvent (process only one event)
630 if (!fInit) { // to prevent overwrite existing file
631 AliError(Form("Give a version name different from %s",
632 fEventFolderName.Data() )) ;
636 if (strstr(option,"print")) {
641 // // check the QA result for Hits and SDigits
642 // AliQA * qa = AliQA::Instance(AliQA::kPHOS) ;
643 // if ( qa->IsSet(AliQA::kPHOS, AliQA::kSIM, AliQA::kFATAL)) {
644 // AliFatal("QA status in Hits and/or SDIGITS was Fatal") ;
645 // } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kSIM, AliQA::kERROR)) {
646 // AliError("QA status in Hits and/or SDIGITS was Error") ;
647 // } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kSIM, AliQA::kWARNING) ) {
648 // AliWarning("QA status in Hits and/or SDIGITS was Warning") ;
649 // } else if ( qa->IsSet(AliQA::kPHOS, AliQA::kSIM, AliQA::kINFO) ) {
650 // AliInfo("QA status in Hits and/or SDIGITS was Info") ;
653 if(strstr(option,"tim"))
654 gBenchmark->Start("PHOSDigitizer");
656 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
658 // Post Digitizer to the white board
659 gime->PostDigitizer(this) ;
661 if (fLastEvent == -1)
662 fLastEvent = gime->MaxEvent() - 1 ;
664 fLastEvent = fFirstEvent ;
666 Int_t nEvents = fLastEvent - fFirstEvent + 1;
670 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
672 gime->Event(ievent,"S") ;
674 Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
676 // //makes the quality assurance data
677 // if (GetQADataMaker()->IsCycleDone() ) {
678 // GetQADataMaker()->EndOfCycle(AliQA::kDIGITS) ;
679 // GetQADataMaker()->StartOfCycle(AliQA::kDIGITS) ;
681 // GetQADataMaker()->Exec(AliQA::kDIGITS, gime->Digits()) ;
682 // GetQADataMaker()->Increment() ;
686 if(strstr(option,"deb"))
689 //increment the total number of Digits per run
690 fDigitsInRun += gime->Digits()->GetEntriesFast() ;
693 // //Write the quality assurance data only after the last event
694 // if ( fEventCounter == gime->MaxEvent() ) {
695 // GetQADataMaker()->EndOfCycle(AliQA::kDIGITS) ;
696 // GetQADataMaker()->Finish() ;
699 gime->PhosLoader()->CleanDigitizer();
701 if(strstr(option,"tim")){
702 gBenchmark->Stop("PHOSDigitizer");
704 message = " took %f seconds for Digitizing %f seconds per event\n" ;
705 AliInfo(Form( message.Data(),
706 gBenchmark->GetCpuTime("PHOSDigitizer"),
707 gBenchmark->GetCpuTime("PHOSDigitizer")/nEvents ));
711 //____________________________________________________________________________
712 Float_t AliPHOSDigitizer::FrontEdgeTime(TClonesArray * ticks) const
714 // Returns the shortest time among all time ticks
716 ticks->Sort() ; //Sort in accordance with times of ticks
718 AliPHOSTick * ctick = (AliPHOSTick *) it.Next() ;
719 Float_t time = ctick->CrossingTime(fTimeThreshold) ;
722 while((t=(AliPHOSTick*) it.Next())){
723 if(t->GetTime() < time) //This tick starts before crossing
728 time = ctick->CrossingTime(fTimeThreshold) ;
733 //____________________________________________________________________________
734 Bool_t AliPHOSDigitizer::Init()
736 // Makes all memory allocations
738 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
740 AliFatal(Form("Could not obtain the Getter object for file %s and event %s !",
741 GetTitle(), fEventFolderName.Data()));
745 const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
747 fEmcCrystals = geom->GetNModules() * geom->GetNCristalsInModule() ;
749 TString opt("Digits") ;
750 if(gime->VersionExists(opt) ) {
751 AliError(Form("Give a version name different from %s",
752 fEventFolderName.Data() )) ;
757 fLastEvent = fFirstEvent ;
759 fInput = fManager->GetNinputs() ;
763 fInputFileNames = new TString[fInput] ;
764 fEventNames = new TString[fInput] ;
765 fInputFileNames[0] = GetTitle() ;
766 fEventNames[0] = fEventFolderName.Data() ;
768 for (index = 1 ; index < fInput ; index++) {
769 fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0);
770 TString tempo = fManager->GetInputFolderName(index) ;
771 fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added by fManager
774 //to prevent cleaning of this object while GetEvent is called
775 gime->PhosLoader()->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
780 //____________________________________________________________________________
781 void AliPHOSDigitizer::InitParameters()
783 // Set initial parameters Digitizer
785 fPinNoise = 0.004 ; // [GeV]
786 fEMCDigitThreshold = 0.012 ; // [GeV]
787 fCPVNoise = 0.01; // [aux units]
788 fCPVDigitThreshold = 0.09 ; // [aux units]
789 fTimeResolution = 0.5e-9 ; // [sec]
790 fTimeSignalLength = 1.0e-9 ; // [sec]
792 fADCchanelEmc = 1.0; // Coefficient between real and measured energies in EMC
793 fADCpedestalEmc = 0. ; //
794 fNADCemc = (Int_t) TMath::Power(2,16) ; // number of channels in EMC ADC
796 fADCchanelCpv = 0.0012 ; // width of one ADC channel in CPV 'popugais'
797 fADCpedestalCpv = 0.012 ; //
798 fNADCcpv = (Int_t) TMath::Power(2,12); // number of channels in CPV ADC
800 // fTimeThreshold = 0.001*10000000 ; //Means 1 MeV in terms of SDigits amplitude
801 fTimeThreshold = 0.001 ; // [GeV]
802 SetEventRange(0,-1) ;
806 //__________________________________________________________________
807 void AliPHOSDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
809 // Allows to produce digits by superimposing background and signal event.
810 // It is assumed, that headers file with SIGNAL events is opened in
812 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
813 // Thus we avoid writing (changing) huge and expensive
814 // backgound files: all output will be writen into SIGNAL, i.e.
815 // opened in constructor file.
817 // One can open as many files to mix with as one needs.
818 // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
821 if( strcmp(fEventFolderName, "") == 0 )
825 Warning("MixWith", "Cannot use this method with AliRunDigitizer\n" ) ;
828 // looking for file which contains AliRun
829 if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
830 AliError(Form("File %s does not exist!", alirunFileName.Data())) ;
833 // looking for the file which contains SDigits
834 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
835 TString fileName( gime->GetSDigitsFileName() ) ;
836 if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
837 fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
838 if ( (gSystem->AccessPathName(fileName)) ) {
839 AliError(Form("The file %s does not exist!", fileName.Data())) ;
842 // need to increase the arrays
843 TString tempo = fInputFileNames[fInput-1] ;
844 delete [] fInputFileNames ;
845 fInputFileNames = new TString[fInput+1] ;
846 fInputFileNames[fInput-1] = tempo ;
848 tempo = fEventNames[fInput-1] ;
849 delete [] fEventNames ;
850 fEventNames = new TString[fInput+1] ;
851 fEventNames[fInput-1] = tempo ;
853 fInputFileNames[fInput] = alirunFileName ;
854 fEventNames[fInput] = eventFolderName ;
858 //__________________________________________________________________
859 void AliPHOSDigitizer::Print(const Option_t *)const
861 // Print Digitizer's parameters
862 AliInfo(Form("\n------------------- %s -------------", GetName() )) ;
863 if( strcmp(fEventFolderName.Data(), "") != 0 ){
864 printf(" Writing Digits to branch with title %s\n", fEventFolderName.Data()) ;
868 nStreams = GetNInputStreams() ;
873 for (index = 0 ; index < nStreams ; index++) {
874 TString tempo(fEventNames[index]) ;
876 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[index], tempo) ;
877 TString fileName( gime->GetSDigitsFileName() ) ;
878 if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name
879 fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index] + ".root" ;
880 printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ;
882 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
883 printf("\nWriting digits to %s", gime->GetDigitsFileName().Data()) ;
885 printf("\nWith following parameters:\n") ;
886 printf(" Electronics noise in EMC (fPinNoise) = %f GeV\n", fPinNoise ) ;
887 printf(" Threshold in EMC (fEMCDigitThreshold) = %f GeV\n", fEMCDigitThreshold ) ;
888 printf(" Noise in CPV (fCPVNoise) = %f aux units\n", fCPVNoise ) ;
889 printf(" Threshold in CPV (fCPVDigitThreshold) = %f aux units\n",fCPVDigitThreshold ) ;
890 printf(" ---------------------------------------------------\n") ;
893 AliInfo(Form("AliPHOSDigitizer not initialized" )) ;
897 //__________________________________________________________________
898 void AliPHOSDigitizer::PrintDigits(Option_t * option)
900 // Print a table of digits
902 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
903 TClonesArray * digits = gime->Digits() ;
905 AliInfo(Form("%d", digits->GetEntriesFast())) ;
906 printf("\nevent %d", gAlice->GetEvNumber()) ;
907 printf("\n Number of entries in Digits list %d", digits->GetEntriesFast() ) ;
910 if(strstr(option,"all")||strstr(option,"EMC")){
912 AliPHOSDigit * digit;
913 printf("\nEMC digits (with primaries):\n") ;
914 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
915 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
917 for (index = 0 ; (index < digits->GetEntriesFast()) &&
918 (dynamic_cast<AliPHOSDigit *>(digits->At(index))->GetId() <= maxEmc) ; index++) {
919 digit = (AliPHOSDigit * ) digits->At(index) ;
920 if(digit->GetNprimary() == 0)
922 // printf("%6d %8d %6.5e %4d %2d :",
923 // digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ; // YVK
924 printf("%6d %.4f %6.5e %4d %2d :",
925 digit->GetId(), digit->GetEnergy(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;
927 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
928 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
934 if(strstr(option,"all")||strstr(option,"CPV")){
936 //loop over CPV digits
937 AliPHOSDigit * digit;
938 printf("\nCPV digits (with primaries):\n") ;
939 printf("\n Id Amplitude Time Index Nprim: Primaries list \n") ;
940 Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
942 for (index = 0 ; index < digits->GetEntriesFast(); index++) {
943 digit = (AliPHOSDigit * ) digits->At(index) ;
944 if(digit->GetId() > maxEmc){
945 printf("%6d %8d %4d %2d :",
946 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;
948 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
949 printf("%d ",digit->GetPrimary(iprimary+1) ) ;
958 //__________________________________________________________________
959 Float_t AliPHOSDigitizer::TimeOfNoise(void) const
960 { // Calculates the time signal generated by noise
961 //PH Info("TimeOfNoise", "Change me") ;
962 return gRandom->Rndm() * 1.28E-5;
965 //__________________________________________________________________
966 void AliPHOSDigitizer::Unload()
970 for(i = 1 ; i < fInput ; i++){
971 TString tempo(fEventNames[i]) ;
973 AliPHOSGetter * gime = AliPHOSGetter::Instance(fInputFileNames[i], tempo) ;
974 gime->PhosLoader()->UnloadSDigits() ;
977 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName) ;
978 gime->PhosLoader()->UnloadDigits() ;
981 //____________________________________________________________________________
982 void AliPHOSDigitizer::WriteDigits()
985 // Makes TreeD in the output file.
986 // Check if branch already exists:
987 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
988 // already existing branches.
989 // else creates branch with Digits, named "PHOS", title "...",
990 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
991 // and names of files, from which digits are made.
993 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
994 const TClonesArray * digits = gime->Digits() ;
995 TTree * treeD = gime->TreeD();
997 // -- create Digits branch
998 Int_t bufferSize = 32000 ;
999 TBranch * digitsBranch = treeD->Branch("PHOS","TClonesArray",&digits,bufferSize);
1000 digitsBranch->SetTitle(fEventFolderName);
1001 digitsBranch->Fill() ;
1003 gime->WriteDigits("OVERWRITE");
1004 gime->WriteDigitizer("OVERWRITE");