]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelection.cxx
AliAnalysisManager::StartAnalysis return now a code:
[u/mrichter/AliRoot.git] / ANALYSIS / AliPhysicsSelection.cxx
1 /* $Id: AliPhysicsSelection.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 //-------------------------------------------------------------------------
19 //                      Implementation of   Class AliPhysicsSelection
20 // This class selects collision candidates from data runs, applying selection cuts on triggers 
21 // and background rejection based on the content of the ESD
22 //
23 // Usage:
24 //
25 // Create the object:
26 //   fPhysicsSelection = new AliPhysicsSelection;
27 //
28 // For MC data, call
29 //   fPhysicsSelection->SetAnalyzeMC()
30 //
31 // To check if an event is a collision candidate, use:
32 //   fPhysicsSelection->IsCollisionCandidate(fESD)
33 //
34 // After processing save the resulting histograms to a file with (a folder physics_selection 
35 //   will be created that contains the histograms):
36 //   fPhysicsSelection->SaveHistograms("physics_selection")
37 //
38 // To print statistics after processing use:
39 //   fPhysicsSelection->Print();
40 //
41 // Usually the class selects the trigger scheme by itself depending on the run number.
42 // Nevertheless, you can do that manually by calling AddCollisionTriggerClass() and AddBGTriggerClass()
43 // Example:
44 // To define the class CINT1B-ABCE-NOPF-ALL as collision trigger (those will be accepted as  
45 // collision candidates when they pass the selection):
46 //   AddCollisionTriggerClass("+CINT1B-ABCE-NOPF-ALL #769 #3119");
47 // To select on bunch crossing IDs in addition, use:
48 //   AddCollisionTriggerClass("+CINT1B-ABCE-NOPF-ALL #769 #3119");
49 // To define the class CINT1A-ABCE-NOPF-ALL as a background trigger (those will only be counted
50 // for the control histograms):
51 //   AddBGTriggerClass("+CINT1A-ABCE-NOPF-ALL");
52 // You can also specify more than one trigger class in a string or you can require that some are *not*
53 // present. The following line would require CSMBA-ABCE-NOPF-ALL, but CSMBB-ABCE-NOPF-ALL is not allowed
54 // to be present:
55 //   AddBGTriggerClass("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL");
56 //
57 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN
58 //-------------------------------------------------------------------------
59
60 #include <Riostream.h>
61 #include <TH1F.h>
62 #include <TH2F.h>
63 #include <TList.h>
64 #include <TIterator.h>
65 #include <TDirectory.h>
66 #include <TObjArray.h>
67
68 #include <AliPhysicsSelection.h>
69
70 #include <AliTriggerAnalysis.h>
71 #include <AliLog.h>
72
73 #include <AliESDEvent.h>
74
75 ClassImp(AliPhysicsSelection)
76
77 AliPhysicsSelection::AliPhysicsSelection() :
78   AliAnalysisCuts("AliPhysicsSelection", "AliPhysicsSelection"),
79   fCurrentRun(-1),
80   fMC(kFALSE),
81   fCollTrigClasses(),
82   fBGTrigClasses(),
83   fTriggerAnalysis(),
84   fBackgroundIdentification(0),
85   fHistStatistics(0),
86   fHistBunchCrossing(0),
87   fSkipTriggerClassSelection(0),
88   fUsingCustomClasses(0)
89 {
90   // constructor
91   
92   fCollTrigClasses.SetOwner(1);
93   fBGTrigClasses.SetOwner(1);
94   fTriggerAnalysis.SetOwner(1);
95   
96   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
97 }
98     
99 AliPhysicsSelection::~AliPhysicsSelection()
100 {
101   // destructor
102   
103   fCollTrigClasses.Delete();
104   fBGTrigClasses.Delete();
105   fTriggerAnalysis.Delete();
106
107   if (fHistStatistics)
108   {
109     delete fHistStatistics;
110     fHistStatistics = 0;
111   }
112   
113   if (fHistBunchCrossing)
114   {
115     delete fHistBunchCrossing;
116     fHistBunchCrossing = 0;
117   }
118 }
119
120 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
121 {
122   // checks if the given trigger class(es) are found for the current event
123   // format of trigger: +TRIGGER1 -TRIGGER2
124   //   requires TRIGGER1 and rejects TRIGGER2
125   
126   Bool_t foundBCRequirement = kFALSE;
127   Bool_t foundCorrectBC = kFALSE;
128   
129   TString str(trigger);
130   TObjArray* tokens = str.Tokenize(" ");
131   
132   for (Int_t i=0; i < tokens->GetEntries(); i++)
133   {
134     TString str2(((TObjString*) tokens->At(i))->String());
135     
136     if (str2[0] == '+' || str2[0] == '-')
137     {
138       Bool_t flag = (str2[0] == '+');
139       
140       str2.Remove(0, 1);
141       
142       if (flag && !aEsd->IsTriggerClassFired(str2))
143       {
144         AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
145         delete tokens;
146         return kFALSE;
147       }
148       if (!flag && aEsd->IsTriggerClassFired(str2))
149       {
150         AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
151         delete tokens;
152         return kFALSE;
153       }
154     }
155     else if (str2[0] == '#')
156     {
157       foundBCRequirement = kTRUE;
158     
159       str2.Remove(0, 1);
160       
161       Int_t bcNumber = str2.Atoi();
162       AliDebug(AliLog::kDebug, Form("Checking for bunch crossing number %d", bcNumber));
163       
164       if (aEsd->GetBunchCrossNumber() == bcNumber)
165       {
166         foundCorrectBC = kTRUE;
167         AliDebug(AliLog::kDebug, Form("Found correct bunch crossing %d", bcNumber));
168       }
169     }
170     else
171       AliFatal(Form("Invalid trigger syntax: %s", trigger));
172   }
173   
174   delete tokens;
175   
176   if (foundBCRequirement && !foundCorrectBC)
177     return kFALSE;
178   
179   return kTRUE;
180 }
181     
182 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
183 {
184   // checks if the given event is a collision candidate
185   
186   if (fCurrentRun != aEsd->GetRunNumber())
187     if (!Initialize(aEsd->GetRunNumber()))
188       AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
189     
190   const AliESDHeader* esdHeader = aEsd->GetHeader();
191   if (!esdHeader)
192   {
193     AliError("ESD Header could not be retrieved");
194     return kFALSE;
195   }
196   
197   // check event type; should be PHYSICS = 7 for data and 0 for MC
198   if (!fMC)
199   {
200     if (esdHeader->GetEventType() != 7)
201       return kFALSE;
202   }
203   else
204   {
205     if (esdHeader->GetEventType() != 0)
206       AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
207   }
208   
209   Bool_t accept = kFALSE;
210     
211   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
212   for (Int_t i=0; i < count; i++)
213   {
214     const char* triggerClass = 0;
215     if (i < fCollTrigClasses.GetEntries())
216       triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
217     else
218       triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
219   
220     AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
221   
222     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
223   
224     triggerAnalysis->FillTriggerClasses(aEsd);
225     
226     if (CheckTriggerClass(aEsd, triggerClass))
227     {
228       triggerAnalysis->FillHistograms(aEsd);
229   
230       fHistStatistics->Fill(1, i);
231       
232       // hardware trigger (should only remove events for MC)
233       // replay CINT1B hardware trigger
234       // TODO this has to depend on the actual hardware trigger (and that depends on the run...)
235       Int_t fastORHW = triggerAnalysis->SPDFiredChips(aEsd, 1); // SPD number of chips from trigger bits (!)
236       Bool_t v0AOnline = (triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kASide, kTRUE) == AliTriggerAnalysis::kV0BB);
237       Bool_t v0COnline = (triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kCSide, kTRUE) == AliTriggerAnalysis::kV0BB);
238       
239       if (fastORHW == 0 && !v0AOnline && !v0COnline)
240       {
241         AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
242         continue;
243       }
244       
245       fHistStatistics->Fill(2, i);
246     
247       // offline trigger
248       AliTriggerAnalysis::V0Decision v0ADecision = triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kASide, kFALSE);
249       AliTriggerAnalysis::V0Decision v0CDecision = triggerAnalysis->V0Trigger(aEsd, AliTriggerAnalysis::kCSide, kFALSE);
250       
251       Bool_t v0AOffline = (v0ADecision == AliTriggerAnalysis::kV0BB);
252       Bool_t v0COffline = (v0CDecision == AliTriggerAnalysis::kV0BB);
253         
254       Int_t fastOROffline = triggerAnalysis->SPDFiredChips(aEsd, 0); // SPD number of chips from clusters (!)
255       Bool_t v0ABG = (v0ADecision == AliTriggerAnalysis::kV0BG);
256       Bool_t v0CBG = (v0CDecision == AliTriggerAnalysis::kV0BG);
257       Bool_t v0BG = v0ABG || v0CBG;
258     
259       if (fastOROffline > 0)
260         fHistStatistics->Fill(3, i);
261       if (fastOROffline > 1)
262         fHistStatistics->Fill(4, i);
263         
264       if (v0AOffline)
265         fHistStatistics->Fill(5, i);
266       if (v0COffline)
267         fHistStatistics->Fill(6, i);
268       if (v0ABG)
269         fHistStatistics->Fill(7, i);
270       if (v0CBG)
271         fHistStatistics->Fill(8, i);
272         
273       if (v0ADecision == AliTriggerAnalysis::kV0Fake)
274         fHistStatistics->Fill(9, i);
275       if (v0CDecision == AliTriggerAnalysis::kV0Fake)
276         fHistStatistics->Fill(10, i);
277       
278       if (fastOROffline > 1 && !v0BG)
279         fHistStatistics->Fill(11, i);
280         
281       if ((fastOROffline > 0 || v0AOffline || v0COffline) && !v0BG)
282         fHistStatistics->Fill(12, i);
283     
284       if (fastOROffline > 0 && (v0AOffline || v0COffline) && !v0BG)
285         fHistStatistics->Fill(13, i);
286   
287       if (v0AOffline && v0COffline && !v0BG)
288         fHistStatistics->Fill(14, i);
289       
290       if (fastOROffline > 1 || (fastOROffline > 0 && (v0AOffline || v0COffline)) || (v0AOffline && v0COffline))
291       {
292         if (!v0BG)
293         {
294           fHistStatistics->Fill(15, i);
295       
296           if (fBackgroundIdentification && !fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd)))
297           {
298             AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
299             fHistStatistics->Fill(16, i);
300           }
301           else
302           {
303             AliDebug(AliLog::kDebug, "Accepted event for histograms");
304             
305             fHistStatistics->Fill(17, i);
306             fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i);
307             if (i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection)
308               accept = kTRUE;
309           }
310         }
311         else
312           AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
313       }
314       else
315         AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
316     }
317   }
318  
319   if (accept)
320     AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
321   
322   return accept;
323 }
324
325 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t runNumber)
326 {
327   // returns the current trigger scheme (classes that are accepted/rejected)
328   
329   if (fMC)
330     return 0;
331     
332   // TODO dependent on run number
333   
334   switch (runNumber)
335   {
336     // CSMBB triggers
337     case 104044:
338     case 105054:
339     case 105057:
340       return 2;
341   }
342   
343   // default: CINT1 suite
344   return 1;
345 }  
346     
347 Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
348 {
349   // initializes the object for the given run
350   // TODO having the run number here and parameters hardcoded is clearly temporary, a way needs to be found to have a CDB-like configuration also for analysis
351   
352   Bool_t oldStatus = TH1::AddDirectoryStatus();
353   TH1::AddDirectory(kFALSE);
354   
355   Int_t triggerScheme = GetTriggerScheme(runNumber);
356   
357   if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
358     AliFatal("Processing several runs with different trigger schemes is not supported");
359   
360   AliInfo(Form("Initializing for run %d", runNumber));
361   
362   // initialize first time?
363   if (fCurrentRun == -1)
364   {
365     if (fUsingCustomClasses) {
366       AliInfo("Using user-provided trigger classes");
367     } else {
368       switch (triggerScheme)
369       {
370       case 0:
371         fCollTrigClasses.Add(new TObjString(""));
372         break;
373         
374       case 1:
375         fCollTrigClasses.Add(new TObjString("+CINT1B-ABCE-NOPF-ALL"));
376         fBGTrigClasses.Add(new TObjString("+CINT1A-ABCE-NOPF-ALL"));
377         fBGTrigClasses.Add(new TObjString("+CINT1C-ABCE-NOPF-ALL"));
378         fBGTrigClasses.Add(new TObjString("+CINT1-E-NOPF-ALL"));
379         break;
380         
381       case 2:
382         fCollTrigClasses.Add(new TObjString("+CSMBB-ABCE-NOPF-ALL"));
383         fBGTrigClasses.Add(new TObjString("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
384         fBGTrigClasses.Add(new TObjString("+CSMBC-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
385         break;
386         
387       default:
388         AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
389       }
390     }
391     
392     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
393     
394     for (Int_t i=0; i<count; i++)
395     {
396       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
397       triggerAnalysis->SetAnalyzeMC(fMC);
398       triggerAnalysis->EnableHistograms();
399       triggerAnalysis->SetSPDGFOThreshhold(1);
400       fTriggerAnalysis.Add(triggerAnalysis);
401     }
402       
403     if (fHistStatistics)
404       delete fHistStatistics;
405   
406     fHistStatistics = new TH2F("fHistStatistics", "fHistStatistics;;", 17, 0.5, 17.5, count, -0.5, -0.5 + count);
407       
408     Int_t n = 1;
409     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Trigger class");
410     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Hardware trigger");
411     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1");
412     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2");
413     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A");
414     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C");
415     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A BG");
416     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C BG");
417     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A Fake");
418     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C Fake");
419     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2 &!V0 BG");
420     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 1 | V0A | V0C) & !V0 BG");
421     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1 & (V0A | V0C) & !V0 BG");
422     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A & V0C & !V0 BG");
423     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A & V0C)) & !V0 BG");
424     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification");
425     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted");
426     
427     if (fHistBunchCrossing)
428       delete fHistBunchCrossing;
429   
430     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
431       
432     n = 1;
433     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
434     {
435       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
436       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
437       n++;
438     }
439     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
440     {
441       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
442       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
443       n++;
444     }
445   }
446     
447   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
448   for (Int_t i=0; i<count; i++)
449   {
450     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
451   
452     switch (runNumber)
453     {
454       case 104315:
455       case 104316:
456       case 104320:
457       case 104321:
458       case 104439:
459         triggerAnalysis->SetV0TimeOffset(7.5);
460         break;
461       default:
462         triggerAnalysis->SetV0TimeOffset(0);
463     }
464   }
465     
466   fCurrentRun = runNumber;
467   
468   TH1::AddDirectory(oldStatus);
469   
470   return kTRUE;
471 }
472
473 void AliPhysicsSelection::Print(Option_t* /* option */) const
474 {
475   // print the configuration
476   
477   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
478   
479   Printf("Collision trigger classes:");
480   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
481     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
482   
483   Printf("Background trigger classes:");
484   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
485     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
486
487   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
488   
489   if (triggerAnalysis)
490   {
491     if (triggerAnalysis->GetV0TimeOffset() > 0)
492       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
493     
494     Printf("\nTotal available events:");
495     
496     triggerAnalysis->PrintTriggerClasses();
497   }
498   
499   if (fHistStatistics && fCollTrigClasses.GetEntries() > 0)
500   {
501     Printf("\nSelection statistics for first collision trigger (%s):", ((TObjString*) fCollTrigClasses.First())->String().Data());
502     
503     Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics->GetBinContent(1, 1));
504     Printf("Selected collision candidates: %d", (Int_t) fHistStatistics->GetBinContent(fHistStatistics->GetXaxis()->FindBin("Accepted"), 1));
505   }
506   
507   if (fHistBunchCrossing)
508   {
509     Printf("\nBunch crossing statistics:");
510     
511     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
512     {
513       TString str;
514       str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
515       
516       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
517         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
518           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
519        
520       Printf("%s", str.Data());
521     }
522     
523     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
524     {
525       Int_t count = 0;
526       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
527       {
528         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
529           count++;
530       }
531       if (count > 1)
532         Printf("WARNING: Bunch crossing %d has more than one trigger class active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
533     }
534   }
535
536   if (fUsingCustomClasses)        
537     Printf("WARNING: Using custom trigger classes!");
538   if (fSkipTriggerClassSelection) 
539     Printf("WARNING: Skipping trigger class selection!");
540 }
541
542 Long64_t AliPhysicsSelection::Merge(TCollection* list)
543 {
544   // Merge a list of AliMultiplicityCorrection objects with this (needed for
545   // PROOF).
546   // Returns the number of merged objects (including this).
547
548   if (!list)
549     return 0;
550
551   if (list->IsEmpty())
552     return 1;
553
554   TIterator* iter = list->MakeIterator();
555   TObject* obj;
556   
557   // collections of all histograms
558   const Int_t nHists = 9;
559   TList collections[nHists];
560
561   Int_t count = 0;
562   while ((obj = iter->Next())) {
563
564     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
565     if (entry == 0) 
566       continue;
567       
568     collections[0].Add(&(entry->fTriggerAnalysis));
569     if (entry->fHistStatistics)
570       collections[1].Add(entry->fHistStatistics);
571     if (entry->fHistBunchCrossing)
572       collections[2].Add(entry->fHistBunchCrossing);
573     if (entry->fBackgroundIdentification)
574       collections[3].Add(entry->fBackgroundIdentification);
575
576     count++;
577   }
578
579   fTriggerAnalysis.Merge(&collections[0]);
580   if (fHistStatistics)
581     fHistStatistics->Merge(&collections[1]);
582   if (fHistBunchCrossing)
583     fHistBunchCrossing->Merge(&collections[2]);
584   if (fBackgroundIdentification)
585     fBackgroundIdentification->Merge(&collections[3]);
586   
587   delete iter;
588
589   return count+1;
590 }
591
592 void AliPhysicsSelection::SaveHistograms(const char* folder) const
593 {
594   // write histograms to current directory
595   
596   if (!fHistStatistics)
597     return;
598     
599   if (folder)
600   {
601     gDirectory->mkdir(folder);
602     gDirectory->cd(folder);
603   }
604   
605   fHistStatistics->Write();
606   fHistBunchCrossing->Write();
607   
608   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
609   for (Int_t i=0; i < count; i++)
610   {
611     TString triggerClass = "trigger_histograms_";
612     if (i < fCollTrigClasses.GetEntries())
613       triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
614     else
615       triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
616   
617     gDirectory->mkdir(triggerClass);
618     gDirectory->cd(triggerClass);
619   
620     static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
621     
622     gDirectory->cd("..");
623   }
624  
625   if (fBackgroundIdentification)
626   {
627     gDirectory->mkdir("background_identification");
628     gDirectory->cd("background_identification");
629       
630     fBackgroundIdentification->GetOutput()->Write();
631       
632     gDirectory->cd("..");
633   }
634   
635   if (folder)
636     gDirectory->cd("..");
637 }