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 //_________________________________________________________________________
19 //*-- Author : Dmitri Peressounko (SUBATECH & Kurchatov Institute)
20 //////////////////////////////////////////////////////////////////////////////
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).
23 // In addition it performs mixing of summable digits from different events.
25 // For each event two branches are created in TreeD:
26 // "PHOS" - list of digits
27 // "AliPHOSDigitizer" - AliPHOSDigitizer with all parameters used in digitization
29 // Note, that one can set a title for new digits branch, and repeat digitization with
30 // another set of parameters.
33 // root[0] AliPHOSDigitizer * d = new AliPHOSDigitizer() ;
34 // root[1] d->ExecuteTask()
35 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 // //Digitizes SDigitis in all events found in file galice.root
38 // root[2] AliPHOSDigitizer * d1 = new AliPHOSDigitizer("galice1.root") ;
39 // // Will read sdigits from galice1.root
40 // root[3] d1->MixWith("galice2.root")
41 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 // // Reads another set of sdigits from galice2.root
43 // root[3] d1->MixWith("galice3.root")
44 // // Reads another set of sdigits from galice3.root
45 // root[4] d->ExecuteTask("deb timing")
46 // // Reads SDigits from files galice1.root, galice2.root ....
47 // // mixes them and stores produced Digits in file galice1.root
48 // // deb - prints number of produced digits
49 // // deb all - prints list of produced digits
50 // // timing - prints time used for digitization
53 // --- ROOT system ---
59 #include "TObjString.h"
60 #include "TBenchmark.h"
61 // --- Standard library ---
64 // --- AliRoot header files ---
67 #include "AliPHOSDigit.h"
69 #include "AliPHOSDigitizer.h"
70 #include "AliPHOSSDigitizer.h"
71 #include "AliPHOSGeometry.h"
73 ClassImp(AliPHOSDigitizer)
76 //____________________________________________________________________________
77 AliPHOSDigitizer::AliPHOSDigitizer():TTask("AliPHOSDigitizer","")
85 fEMCDigitThreshold = 0.01 ;
87 fCPVDigitThreshold = 0.09 ;
88 fPPSDNoise = 0.0000001;
89 fPPSDDigitThreshold = 0.0000002 ;
90 fInitialized = kFALSE ;
98 //____________________________________________________________________________
99 void AliPHOSDigitizer::Init()
101 // Makes all memory allocations
102 // Adds Digitizer task to the folder of PHOS tasks
106 if (fHeaderFiles == 0) {
107 fHeaderFiles = new TClonesArray("TObjString",1) ;
108 new((*fHeaderFiles)[0]) TObjString("galice.root") ;
111 //Test, if this file already open
113 TFile *file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(0))->GetString() ) ;
116 if(((TObjString *) fHeaderFiles->At(0))->GetString().Contains("rfio"))
117 file = TFile::Open(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
119 file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString(),"update") ;
120 gAlice = (AliRun *) file->Get("gAlice") ; //If not read yet
123 file = new TFile(((TObjString *) fHeaderFiles->At(0))->GetString()) ;
127 if (fSDigitsTitles == 0) {
128 fSDigitsTitles = new TClonesArray("TObjString",1);
129 new((*fSDigitsTitles)[0]) TObjString("") ;
132 fSDigits = new TClonesArray("TClonesArray",1) ;
133 new((*fSDigits)[0]) TClonesArray("AliPHOSDigit",1000) ;
139 fDigits = new TClonesArray("AliPHOSDigit",200000) ;
141 fIevent = new TArrayI(1) ;
142 fIevent->AddAt(-1,0 ) ;
143 fIeventMax = new TArrayI(1) ;
145 fIeventMax->AddAt((Int_t) gAlice->TreeE()->GetEntries(), 0 );
147 //add Task to //YSAlice/tasks/(S)Diditizer/PHOS
148 TFolder * alice = (TFolder*)gROOT->GetListOfBrowsables()->FindObject("YSAlice") ;
149 TTask * aliceSD = (TTask*)alice->FindObject("tasks/(S)Digitizer") ;
150 TTask * phosSD = (TTask*)aliceSD->GetListOfTasks()->FindObject("PHOS") ;
153 fInitialized = kTRUE ;
157 //____________________________________________________________________________
158 AliPHOSDigitizer::AliPHOSDigitizer(const char *headerFile,const char *sDigitsTitle):
159 TTask("AliPHOSDigitizer","")
162 fHeaderFiles = new TClonesArray("TObjString",1) ;
163 new((*fHeaderFiles)[0]) TObjString(headerFile) ;
165 fSDigitsTitles = new TClonesArray("TObjString",1); // Title name of the SDigits branch
166 new((*fSDigitsTitles)[0]) TObjString(sDigitsTitle) ;
170 fEMCDigitThreshold = 0.01 ;
172 fCPVDigitThreshold = 0.09 ;
173 fPPSDNoise = 0.0000001;
174 fPPSDDigitThreshold = 0.0000002 ;
175 fInitialized = kFALSE ;
181 //____________________________________________________________________________
182 AliPHOSDigitizer::~AliPHOSDigitizer()
186 if(fHeaderFiles) delete fHeaderFiles ;
187 if(fSDigitsTitles) delete fSDigitsTitles ;
188 if(fSDigits) delete fSDigits ;
189 if(fDigits) delete fDigits ;
191 //____________________________________________________________________________
192 void AliPHOSDigitizer::Reset()
194 // sets current event number to the first simulated event
200 for(inputs = 0; inputs < fNinputs ;inputs++)
201 fIevent->AddAt(-1, inputs ) ;
204 //____________________________________________________________________________
205 Bool_t AliPHOSDigitizer::Combinator()
207 // Makes all desirable combinations of Signal+Background,
208 // returns kFALSE when all combinations are made.
209 // May be useful to introduce options like "One-to-One", "All-to-One" and "All-to-All" ?
210 //realizing "One-to-One" option...
216 Bool_t endNotReached = kTRUE ;
218 for(inputs = 0; (inputs < fNinputs) && endNotReached ;inputs++){
219 if(fIevent->At(inputs)+1 < fIeventMax->At(inputs))
220 fIevent->AddAt(fIevent->At(inputs)+1, inputs ) ;
223 endNotReached = kFALSE ;
224 else //for inputs other than base one start from the beginning
225 fIevent->AddAt(0, inputs ) ;
228 return endNotReached ;
232 //____________________________________________________________________________
233 void AliPHOSDigitizer::Digitize(Option_t *option)
236 // Makes the digitization of the collected summable digits.
237 // It first creates the array of all PHOS modules
238 // filled with noise (different for EMC, CPV and PPSD) and
239 // then adds contributions from SDigits.
240 // This design avoids scanning over the list of digits to add
241 // contribution to new SDigits only.
248 AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;
249 AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance( phos->GetGeometry()->GetName(), phos->GetGeometry()->GetTitle() );
251 //Making digits with noise, first EMC
252 Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
257 TString name = geom->GetName() ;
259 if ( name == "IHEP" || name == "MIXT" )
260 nCPV =nEMC + geom->GetNumberOfCPVPadsZ()*geom->GetNumberOfCPVPadsPhi()*
261 geom->GetNCPVModules()*geom->GetNumberOfCPVLayers() ;
265 if ( name == "GPS2" || name == "MIXT" )
266 nPPSD =nCPV+2*geom->GetNPPSDModules()*geom->GetNumberOfModulesPhi()*geom->GetNumberOfModulesZ()*
267 geom->GetNumberOfPadsPhi()*geom->GetNumberOfPadsZ() ;
272 fDigits->Expand(nPPSD) ;
275 for(absID = 1; absID <= nEMC; absID++){
276 Float_t noise = gRandom->Gaus(0., fPinNoise) ;
277 new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
280 for(absID = nEMC+1; absID <= nCPV; absID++){
281 Float_t noise = gRandom->Gaus(0., fCPVNoise) ;
282 new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
285 for(absID = nCPV+1; absID <= nPPSD; absID++){
286 Float_t noise = gRandom->Gaus(0., fPPSDNoise) ;
287 new((*fDigits)[absID-1]) AliPHOSDigit( -1,absID,fSDigitizer->Digitize(noise) ) ;
291 // Now look throught (unsorted) list of SDigits and add corresponding digits
292 AliPHOSDigit *curSDigit ;
293 AliPHOSDigit *digit ;
296 for(inputs = 0; inputs< fNinputs ; inputs++){ //loop over (possible) merge sources
298 TClonesArray * sdigits= (TClonesArray *)fSDigits->At(inputs) ;
301 Int_t nSDigits = sdigits->GetEntries() ;
302 for(isdigit=0;isdigit< nSDigits; isdigit++){
303 curSDigit = (AliPHOSDigit *)sdigits->At(isdigit) ;
304 if(inputs) //Shift primaries for non-background sdigits
305 curSDigit->ShiftPrimary(inputs) ;
306 digit = (AliPHOSDigit *)fDigits->At(curSDigit->GetId() - 1);
307 *digit = *digit + *curSDigit ;
312 //remove digits below thresholds
313 for(absID = 0; absID < nEMC ; absID++)
314 if(fSDigitizer->Calibrate(((AliPHOSDigit*)fDigits->At(absID))->GetAmp()) < fEMCDigitThreshold)
315 fDigits->RemoveAt(absID) ;
316 for(absID = nEMC; absID < nCPV ; absID++)
317 if(fSDigitizer->Calibrate(((AliPHOSDigit*)fDigits->At(absID))->GetAmp()) < fCPVDigitThreshold)
318 fDigits->RemoveAt(absID) ;
319 for(absID = nCPV; absID < nPPSD ; absID++)
320 if(fSDigitizer->Calibrate(((AliPHOSDigit *)fDigits->At(absID))->GetAmp()) < fPPSDDigitThreshold)
321 fDigits->RemoveAt(absID) ;
323 fDigits->Compress() ;
325 Int_t ndigits = fDigits->GetEntriesFast() ;
327 fDigits->Expand(ndigits) ;
330 //Set indexes in list of digits
332 for (i = 0 ; i < ndigits ; i++) {
333 AliPHOSDigit * digit = (AliPHOSDigit *) fDigits->At(i) ;
334 digit->SetIndexInList(i) ;
337 //____________________________________________________________________________
338 void AliPHOSDigitizer::WriteDigits()
341 // Makes TreeD in the output file.
342 // Check if branch already exists:
343 // if yes, exit without writing: ROOT TTree does not support overwriting/updating of
344 // already existing branches.
345 // else creates branch with Digits, named "PHOS", title "...",
346 // and branch "AliPHOSDigitizer", with the same title to keep all the parameters
347 // and names of files, from which digits are made.
349 gAlice->GetEvent(fIevent->At(0)) ; // Suitable only for One-To-One mixing
350 gAlice->SetEvent(fIevent->At(0)) ; // for all-to-all will produce a lot of branches in TreeD
352 if(gAlice->TreeD()==0)
353 gAlice->MakeTree("D") ;
355 //Check, if this branch already exits?
356 TBranch * digitsBranch = 0;
357 TBranch * digitizerBranch = 0;
359 TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
361 Bool_t phosNotFound = kTRUE ;
362 Bool_t digitizerNotFound = kTRUE ;
364 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
367 digitsBranch=(TBranch *) branches->At(ibranch) ;
368 if( (strcmp("PHOS",digitsBranch->GetName())==0 ) &&
369 (fDigitsTitle.CompareTo(digitsBranch->GetTitle()) == 0) )
370 phosNotFound = kFALSE ;
372 if(digitizerNotFound){
373 digitizerBranch = (TBranch *) branches->At(ibranch) ;
374 if( (strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) &&
375 (fDigitsTitle.CompareTo(digitizerBranch->GetTitle()) == 0))
376 digitizerNotFound = kFALSE ;
381 if(!(digitizerNotFound && phosNotFound)){
382 cout << "AliPHOSDigitizer error: " << endl ;
383 cout << " can not update/overwrite existing branches "<< endl ;
384 cout << " do not write " << endl ;
388 // create new branches
390 //First generate file name
392 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
393 file = new char[strlen(gAlice->GetBaseFile())+20] ;
394 sprintf(file,"%s/PHOS.Digits.root",gAlice->GetBaseFile()) ;
397 TDirectory *cwd = gDirectory;
399 //First create list of sdigits
400 Int_t bufferSize = 32000 ;
401 digitsBranch = gAlice->TreeD()->Branch("PHOS",&fDigits,bufferSize);
402 digitsBranch->SetTitle(fDigitsTitle.Data());
404 digitsBranch->SetFile(file);
405 TIter next( digitsBranch->GetListOfBranches());
407 while ((sbr=(TBranch*)next())) {
413 //second - create Digitizer
414 Int_t splitlevel = 0 ;
415 AliPHOSDigitizer * d = this ;
416 digitizerBranch = gAlice->TreeD()->Branch("AliPHOSDigitizer","AliPHOSDigitizer",
417 &d,bufferSize,splitlevel);
418 digitizerBranch->SetTitle(fDigitsTitle.Data());
420 digitizerBranch->SetFile(file);
421 TIter next( digitizerBranch->GetListOfBranches());
423 while ((sbr=(TBranch*)next())) {
429 digitsBranch->Fill() ;
430 digitizerBranch->Fill() ;
432 gAlice->TreeD()->Write(0,kOverwrite) ;
434 //remove fSDigitizer before new event.
443 //____________________________________________________________________________
444 void AliPHOSDigitizer::Exec(Option_t *option)
448 if(!fInitialized) Init() ;
450 if(strstr(option,"tim"))
451 gBenchmark->Start("PHOSDigitizer");
453 //reset events numbers to start from the beginnig
458 if(!ReadSDigits()) //read sdigits event(s) evaluated by Combinator() from file(s)
461 Digitize(option) ; //Add prepared SDigits to digits and add the noise
464 if(strstr(option,"deb"))
469 if(strstr(option,"tim")){
470 gBenchmark->Stop("PHOSDigitizer");
471 cout << "AliPHOSDigitizer:" << endl ;
472 cout << " took " << gBenchmark->GetCpuTime("PHOSDigitizer") << " seconds for SDigitizing "
473 << gBenchmark->GetCpuTime("PHOSDigitizer")/(fIeventMax->At(0)) << " seconds per event " << endl ;
479 //__________________________________________________________________
480 Bool_t AliPHOSDigitizer::ReadSDigits()
482 // Reads summable digits from the opened files for the particular set of events given by fIevent
484 if(!fInitialized) Init() ;
488 for(inputs = fNinputs-1; inputs >= 0; inputs --){
490 Int_t event = fIevent->At(inputs) ;
492 TFile * file = (TFile*) gROOT->GetFile(((TObjString *) fHeaderFiles->At(inputs))->GetString() ) ;
495 // Get SDigits Tree header from file
497 sprintf(treeName,"TreeS%d",event);
498 TTree * treeS = (TTree*)file->Get(treeName);
501 cout << "Error at AliPHOSDigitizer: no "<<treeName << " in file " << file->GetName() << endl ;
502 cout << "Do nothing " << endl ;
506 TBranch * sdigitsBranch = 0;
507 TBranch * sdigitizerBranch = 0;
509 TObjArray * branches = treeS->GetListOfBranches() ;
511 Bool_t phosNotFound = kTRUE ;
512 Bool_t sdigitizerNotFound = kTRUE ;
514 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
517 sdigitsBranch=(TBranch *) branches->At(ibranch) ;
518 if(( strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
519 ((TObjString*) fSDigitsTitles->At(inputs))->GetString().CompareTo(sdigitsBranch->GetTitle())== 0 )
520 phosNotFound = kFALSE ;
524 if(sdigitizerNotFound){
525 sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
526 if(( strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0) &&
527 ((TObjString*) fSDigitsTitles->At(inputs))->GetString().CompareTo(sdigitizerBranch->GetTitle())== 0 )
528 sdigitizerNotFound = kFALSE ;
533 if(sdigitizerNotFound || phosNotFound){
534 cout << "Can't find Branch with sdigits or SDigitizer in the file " ;
535 if( ((TObjString*)fSDigitsTitles->At(inputs))->GetString().IsNull() )
536 cout << file->GetName() << endl ;
538 cout << ((TObjString*)fSDigitsTitles->At(inputs))->GetString().Data() << endl ;
539 cout << "Do nothing" <<endl ;
543 TClonesArray * sdigits = (TClonesArray*) fSDigits->At(inputs) ;
544 sdigitsBranch->SetAddress(&sdigits) ;
546 AliPHOSSDigitizer *sDigitizer = new AliPHOSSDigitizer();
547 sdigitizerBranch->SetAddress(&sDigitizer) ;
549 sdigitsBranch->GetEntry(0) ;
550 sdigitizerBranch->GetEntry(0) ;
553 fSDigitizer = sDigitizer ;
555 if(!((*fSDigitizer)==(*sDigitizer)) ){
556 cout << "AliPHOSDigitizer ERROR:" << endl ;
557 cout << " you are using sdigits made with different SDigitizers" << endl ;
558 cout << "fSD " << fSDigitizer << " SD" << sDigitizer << endl ;
559 fSDigitizer->Print("") ;
560 sDigitizer->Print("") ;
561 cout << "Do Nothing " << endl ;
566 fPedestal = fSDigitizer->GetPedestalParameter() ;
567 fSlope = fSDigitizer->GetCalibrationParameter() ;
572 //__________________________________________________________________
573 void AliPHOSDigitizer::MixWith(char* HeaderFile, char* sDigitsTitle)
575 // Alows to produce digits by superimposing background and signal event.
576 // It is assumed, that headers file with SIGNAL events is opened in
578 // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed
579 // Thus we avoid writing (changing) huge and expensive
580 // backgound files: all output will be writen into SIGNAL, i.e.
581 // opened in constructor file.
583 // One can open as many files to mix with as one needs.
591 cout << "Specify at least header file to merge"<< endl ;
596 for(inputs = 0; inputs < fNinputs ; inputs++){
597 if(strcmp(((TObjString *)fHeaderFiles->At(inputs))->GetString(),HeaderFile) == 0 ){
598 if(sDigitsTitle == 0){
599 if(((TObjString*)fSDigitsTitles->At(inputs))->GetString().CompareTo("") == 0){
600 cout << "Entry already exists, do not add" << endl ;
605 if(((TObjString*)fSDigitsTitles->At(inputs))->GetString().CompareTo(sDigitsTitle)){
606 cout << "Entry already exists, do not add" << endl ;
612 fHeaderFiles->Expand(fNinputs+1) ;
613 new((*fHeaderFiles)[fNinputs]) TObjString(HeaderFile) ;
616 TFile * file = new TFile(((TObjString *) fHeaderFiles->At(fNinputs))->GetString()) ;
620 fSDigitsTitles->Expand(fNinputs+1) ;
621 new((*fSDigitsTitles)[fNinputs]) TObjString(sDigitsTitle) ;
623 fSDigits->Expand(fNinputs+1) ;
624 new((*fSDigits)[fNinputs]) TClonesArray("AliPHOSDigit",1000) ;
626 fIevent->Set(fNinputs+1) ;
627 fIevent->AddAt(-1, fNinputs) ;
629 fIeventMax->Set(fNinputs+1) ;
631 TTree * te = (TTree *) file->Get("TE") ;
632 fIeventMax->AddAt((Int_t) te->GetEntries(), fNinputs );
637 //__________________________________________________________________
638 void AliPHOSDigitizer::Print(Option_t* option)const {
639 // Print Digitizer's parameters
642 cout << "------------------- "<< GetName() << " -------------" << endl ;
643 cout << "Digitizing sDigits from file(s): " <<endl ;
645 for(input = 0; input < fNinputs ; input++) {
646 cout << " " << ((TObjString *) fHeaderFiles->At(input))->GetString() <<
647 " Branch title:" << ((TObjString *) fSDigitsTitles->At(input))->GetString() << endl ;
650 cout << "Writing digits to " << ((TObjString *) fHeaderFiles->At(0))->GetString() << endl ;
653 cout << "With following parameters: " << endl ;
654 cout << " Electronics noise in EMC (fPinNoise) = " << fPinNoise << endl ;
655 cout << " Threshold in EMC (fEMCDigitThreshold) = " << fEMCDigitThreshold << endl ; ;
656 cout << " Noise in CPV (fCPVNoise) = " << fCPVNoise << endl ;
657 cout << " Threshold in CPV (fCPVDigitThreshold) = " << fCPVDigitThreshold << endl ;
658 cout << " Noise in PPSD (fPPSDNoise) = " << fPPSDNoise << endl ;
659 cout << " Threshold in PPSD (fPPSDDigitThreshold) = " << fPPSDDigitThreshold << endl ;
660 cout << "---------------------------------------------------" << endl ;
663 cout << "AliPHOSDigitizer not initialized " << endl ;
666 //__________________________________________________________________
667 void AliPHOSDigitizer::PrintDigits(Option_t * option){
668 // Print a table of digits
670 cout << "AliPHOSDigitiser:"<< endl ;
671 cout << " Number of entries in Digits list " << fDigits->GetEntriesFast() << endl ;
673 if(strstr(option,"all")){
676 AliPHOSDigit * digit;
677 cout << "Digit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
679 for (index = 0 ; index < fDigits->GetEntries() ; index++) {
680 digit = (AliPHOSDigit * ) fDigits->At(index) ;
681 cout << setw(8) << digit->GetId() << " " << setw(3) << digit->GetAmp() << " "
682 << setw(6) << digit->GetIndexInList() << " "
683 << setw(5) << digit->GetNprimary() <<" ";
686 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
687 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";
693 //__________________________________________________________________
694 void AliPHOSDigitizer::SetSDigitsBranch(const char* title)
696 // we set title (comment) of the SDigits branch in the first! header file
697 if(!fInitialized) Init() ;
699 ((TObjString*) fSDigitsTitles->At(0) )->SetString((char*)title) ;
702 //__________________________________________________________________
703 void AliPHOSDigitizer::SetDigitsBranch(const char* title)
705 //Sets the title (comment) of the branch to which Digits branch
706 if(!fInitialized) Init() ;
708 fDigitsTitle = title ;
711 //__________________________________________________________________