]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelection.cxx
AliCaloPID: Correct matching rejection in case of recalculation in the analysis,...
[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 // The BX ids corresponding to real bunches crossings p2 are
42 // automatically selected. You cannot process runs with different
43 // filling schemes if this option is set. If you want to disable this,
44 // use: 
45 //   fPhysicsSelection->SetUseBXNumbers(0);
46 //
47 //
48 // If you are analizing muons and you want to keep the muon triggers
49 // besides the CINT1B you can set:
50 //   fPhysicsSelection->SetUseMuonTriggers();
51 //
52 //
53 // To compute the Background automatically using the control triggers
54 // use: 
55 //   fPhysicsSelection->SetComputeBG();
56 // this will show the value of the Beam Gas, accidentals and good
57 // events as additional rows in the statistic tables, but it will NOT
58 // subtract the background automatically.
59 // This option will only work for runs taken with the CINT1
60 // suite. This options enables automatically also the usage of BX
61 // numbers. You can only process one run at a time if you require this
62 // options, because it uses the bunch intensity estimated run by run.
63 // 
64 // The BG will usually be more important in the so-called "bin 0": the
65 // class can also compute the statistics table for events in this
66 // bin. Since the definition of bin 0 may in general change from
67 // analysis to analysis, the user needs to provide a callback
68 // implementing the definition of bin zero. The callback should be
69 // implemented as a method in the analysis task and should override
70 // the IsEventInBinZero method of AliAnalysisTaskSE, and should thus
71 // have the the following prototype:
72 //   Bool_t IsEventInBinZero(); 
73 // It should return true if the event is in the bin 0 and it is set by
74 // passing to the physics selection the NAME of the task where the
75 // callback is implemented: 
76 //   fPhysicsSelection->SetBin0Callback("MyTask").
77 //
78 //
79 // Usually the class selects the trigger scheme by itself depending on the run number.
80 // Nevertheless, you can do that manually by calling AddCollisionTriggerClass() and AddBGTriggerClass()
81 // Example:
82 // To define the class CINT1B-ABCE-NOPF-ALL as collision trigger (those will be accepted as  
83 // collision candidates when they pass the selection):
84 //   AddCollisionTriggerClass("+CINT1B-ABCE-NOPF-ALL #769 #3119");
85 // To select on bunch crossing IDs in addition, use:
86 //   AddCollisionTriggerClass("+CINT1B-ABCE-NOPF-ALL #769 #3119");
87 // To define the class CINT1A-ABCE-NOPF-ALL as a background trigger (those will only be counted
88 // for the control histograms):
89 //   AddBGTriggerClass("+CINT1A-ABCE-NOPF-ALL");
90 // You can also specify more than one trigger class in a string or you can require that some are *not*
91 // present. The following line would require CSMBA-ABCE-NOPF-ALL, but CSMBB-ABCE-NOPF-ALL is not allowed
92 // to be present:
93 //   AddBGTriggerClass("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL");
94 //
95 // The class also supports the triggers used in heavy ion runs
96 //
97 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN 
98 //           Michele Floris, CERN
99 //-------------------------------------------------------------------------
100
101 #include <Riostream.h>
102 #include <TH1F.h>
103 #include <TH2F.h>
104 #include <TList.h>
105 #include <TIterator.h>
106 #include <TDirectory.h>
107 #include <TObjArray.h>
108 #include <TPRegexp.h>
109 #include <TFormula.h>
110 #include <TParameter.h>
111 #include <TInterpreter.h>
112
113 #include <AliPhysicsSelection.h>
114
115 #include <AliTriggerAnalysis.h>
116 #include <AliLog.h>
117
118 #include <AliESDEvent.h>
119 #include <AliAnalysisTaskSE.h>
120 #include "AliAnalysisManager.h"
121 #include "TPRegexp.h"
122 #include "TFile.h"
123 #include <AliOADBContainer.h>
124 #include "AliOADBPhysicsSelection.h"
125 #include "AliOADBFillingScheme.h"
126 #include "AliOADBTriggerAnalysis.h"
127
128 ClassImp(AliPhysicsSelection)
129
130 AliPhysicsSelection::AliPhysicsSelection() :
131   AliAnalysisCuts("AliPhysicsSelection", "AliPhysicsSelection"),
132   fCurrentRun(-1),
133   fMC(kFALSE),
134   fCollTrigClasses(),
135   fBGTrigClasses(),
136   fTriggerAnalysis(),
137   fBackgroundIdentification(0),
138 //  fHistStatisticsTokens(0),
139   fHistBunchCrossing(0),
140   fHistTriggerPattern(0),
141   fSkipTriggerClassSelection(0),
142   fUsingCustomClasses(0),
143   fSkipV0(0),
144   fSkipZDCTime(0),
145   fBIFactorA(-1),
146   fBIFactorC(-1),
147   fBIFactorAC(-1), 
148   fComputeBG(0),
149   fBGStatOffset(-1),
150   fUseBXNumbers(1),
151   fUseMuonTriggers(0),
152   fFillingScheme(""),
153   fBin0CallBack(""),
154   fBin0CallBackPointer(0),
155   fIsPP(kFALSE),
156   fPSOADB(0),
157   fFillOADB(0),
158   fTriggerOADB(0),
159   fRegexp(0),
160   fCashedTokens(0)
161   
162 {
163   // constructor
164   
165   fCollTrigClasses.SetOwner(1);
166   fBGTrigClasses.SetOwner(1);
167   fTriggerAnalysis.SetOwner(1);
168   fHistStatistics[0] = 0;
169   fHistStatistics[1] = 0;
170   
171   fRegexp = new TPRegexp("([[:alpha:]]\\w*)");
172   fCashedTokens = new TList;
173   fCashedTokens->SetOwner();
174   
175   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
176 }
177     
178 AliPhysicsSelection::~AliPhysicsSelection()
179 {
180   // destructor
181
182   fCollTrigClasses.Delete();
183   fBGTrigClasses.Delete();
184   fTriggerAnalysis.Delete();
185
186   if (fHistStatistics[0])
187   {
188     delete fHistStatistics[0];
189     fHistStatistics[0] = 0;
190   }
191   if (fHistStatistics[1])
192   {
193     delete fHistStatistics[1];
194     fHistStatistics[1] = 0;
195   }
196   // if (fHistStatisticsTokens) {
197   //   delete fHistStatisticsTokens;
198   //   fHistStatisticsTokens = 0;
199   // }
200   if (fHistBunchCrossing)
201   {
202     delete fHistBunchCrossing;
203     fHistBunchCrossing = 0;
204   }
205   if (fHistTriggerPattern)
206   {
207     delete fHistTriggerPattern;
208     fHistTriggerPattern = 0;
209   }
210   if (fBackgroundIdentification)
211   { 
212     delete fBackgroundIdentification;
213     fBackgroundIdentification = 0;
214   }  
215
216   if (fPSOADB)
217   { 
218     delete fPSOADB;
219     fPSOADB = 0;    
220   }  
221   if (fFillOADB)
222   { 
223     delete fFillOADB;
224     fFillOADB = 0;
225   }  
226   if (fTriggerOADB)
227   { 
228     delete fTriggerOADB;
229     fTriggerOADB = 0;
230   }  
231   
232   if (fRegexp)
233   {
234     delete fRegexp;
235     fRegexp = 0;
236   }
237     
238   if (fCashedTokens)
239   {
240     delete fCashedTokens;
241     fCashedTokens = 0;
242   }
243
244
245 }
246
247 UInt_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger, Int_t& triggerLogic) const
248 {
249   // checks if the given trigger class(es) are found for the current event
250   // format of trigger: +TRIGGER1,TRIGGER1b,TRIGGER1c -TRIGGER2 [#XXX] [&YY] [*ZZ]
251   //   requires one out of TRIGGER1,TRIGGER1b,TRIGGER1c and rejects TRIGGER2
252   //   in bunch crossing XXX
253   //   if successful, YY is returned (for association between entry in fCollTrigClasses and AliVEvent::EOfflineTriggerTypes)
254   //   triggerLogic is filled with ZZ, defaults to kCINT1
255   
256   Bool_t foundBCRequirement = kFALSE;
257   Bool_t foundCorrectBC = kFALSE;
258   
259   UInt_t returnCode = AliVEvent::kUserDefined;
260   triggerLogic = kCINT1;
261   
262   AliDebug(AliLog::kDebug+1, Form("Processing event with triggers %s", aEsd->GetFiredTriggerClasses().Data()));
263   
264   TString str(trigger);
265   TObjArray* tokens = str.Tokenize(" ");
266   
267   for (Int_t i=0; i < tokens->GetEntries(); i++)
268   {
269     TString str2(((TObjString*) tokens->At(i))->String());
270     
271     if (str2[0] == '+' || str2[0] == '-')
272     {
273       Bool_t flag = (str2[0] == '+');
274       
275       str2.Remove(0, 1);
276       
277       TObjArray* tokens2 = str2.Tokenize(",");
278       
279       Bool_t foundTriggerClass = kFALSE;
280       for (Int_t j=0; j < tokens2->GetEntries(); j++)
281       {
282         TString str3(((TObjString*) tokens2->At(j))->String());
283         
284         if (flag && aEsd->IsTriggerClassFired(str3))
285           foundTriggerClass = kTRUE;
286         if (!flag && aEsd->IsTriggerClassFired(str3))
287         {
288           AliDebug(AliLog::kDebug+1, Form("Rejecting event because trigger class %s is present", str3.Data()));
289           delete tokens2;
290           delete tokens;
291           return kFALSE;
292         }
293       }
294       
295       delete tokens2;
296       
297       if (flag && !foundTriggerClass)
298       {
299         AliDebug(AliLog::kDebug+1, Form("Rejecting event because (none of the) trigger class(es) %s is present", str2.Data()));
300         delete tokens;
301         return kFALSE;
302       }
303     }
304     else if (str2[0] == '#')
305     {
306       foundBCRequirement = kTRUE;
307     
308       str2.Remove(0, 1);
309       
310       Int_t bcNumber = str2.Atoi();
311       AliDebug(AliLog::kDebug+1, Form("Checking for bunch crossing number %d", bcNumber));
312       
313       if (aEsd->GetBunchCrossNumber() == bcNumber)
314       {
315         foundCorrectBC = kTRUE;
316         AliDebug(AliLog::kDebug+1, Form("Found correct bunch crossing %d", bcNumber));
317       }
318     }
319     else if (str2[0] == '&')
320     {
321       str2.Remove(0, 1);
322       
323       returnCode = str2.Atoll();
324     }
325     else if (str2[0] == '*')
326     {
327       str2.Remove(0, 1);
328       
329       triggerLogic = str2.Atoi();
330     }
331     else
332       AliFatal(Form("Invalid trigger syntax: %s", trigger));
333   }
334   
335   delete tokens;
336   
337   if (foundBCRequirement && !foundCorrectBC)
338     return kFALSE;
339   
340   return returnCode;
341 }
342
343 //______________________________________________________________________________
344 TObject *AliPhysicsSelection::GetStatistics(const Option_t *option) const
345 {
346 // Get the statistics histograms ("ALL" and "BIN0" and "TOK")
347    TString opt(option);
348    opt.ToUpper();
349    Int_t ihist = 0;
350    if (opt == "ALL") ihist = kStatIdxAll;
351    if (opt == "BIN0") ihist = kStatIdxBin0;
352    //   if (opt == "TOK") return fHistStatisticsTokens;
353    return fHistStatistics[ihist];
354 }   
355
356 //______________________________________________________________________________
357 Bool_t AliPhysicsSelection::EvaluateTriggerLogic(const AliESDEvent* aEsd, AliTriggerAnalysis* triggerAnalysis, const char* triggerLogic, Bool_t offline)
358 {
359   // evaluates trigger logic. If called with no ESD pointer/triggerAnalysis pointer, it just caches the tokens
360   // Fills the statistics histogram, if booked at row i
361   TString trigger(triggerLogic);
362   
363   // add space after each token (to use ReplaceAll later)
364   fRegexp->Substitute(trigger, "$1 ", "g");
365   
366   while (1)
367   {
368     AliDebug(AliLog::kDebug, trigger.Data());
369     
370     TArrayI pos;
371     Int_t nMatches = fRegexp->Match(trigger, "", 0, 2, &pos);
372
373     if (nMatches <= 0)
374       break;
375       
376     TString token(trigger(pos[0], pos[1]-pos[0]+1));
377     
378     TParameter<Int_t>* param = (TParameter<Int_t>*) fCashedTokens->FindObject(token);
379     if (!param)
380     {
381       TInterpreter::EErrorCode error;
382       Int_t bit = gInterpreter->ProcessLine(Form("AliTriggerAnalysis::k%s;", token.Data()), &error);
383       
384       if (error > 0)
385         AliFatal(Form("Trigger token %s unknown", token.Data()));
386       
387       param = new TParameter<Int_t>(token, bit);
388       fCashedTokens->Add(param);
389       AliDebug(AliLog::kDebug, "Added token");
390     }
391     
392     Long64_t bit = param->GetVal();
393     
394     AliDebug(AliLog::kDebug, Form("Tok %d %d %s %lld", pos[0], pos[1], token.Data(), bit));
395     
396     if (offline) 
397       bit |= AliTriggerAnalysis::kOfflineFlag;
398     
399     if(aEsd && triggerAnalysis) {
400       trigger.ReplaceAll(token, Form("%d", triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) bit)));
401       //      if(fHistStatisticsTokens)               
402     }
403   }
404   
405   TFormula formula("formula", trigger);
406   Bool_t result = formula.Eval(0);
407   
408   AliDebug(AliLog::kDebug, Form("%s --> %d", trigger.Data(), result));
409   
410   return result;
411 }
412
413 //______________________________________________________________________________
414 UInt_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
415 {
416   // checks if the given event is a collision candidate
417   //
418   // returns a bit word describing the fired offline triggers (see AliVEvent::EOfflineTriggerTypes)
419   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
420   if (!mgr) {
421     AliError("Cannot get the analysis manager");
422     return 0;
423   } 
424   mgr->LoadBranch("AliESDHeader.");
425   mgr->LoadBranch("AliESDRun.");
426   mgr->LoadBranch("AliMultiplicity.");
427   mgr->LoadBranch("AliESDFMD.");
428   mgr->LoadBranch("AliESDVZERO.");
429   mgr->LoadBranch("AliESDZDC.");
430   mgr->LoadBranch("SPDVertex.");
431   mgr->LoadBranch("PrimaryVertex.");
432   if (fCurrentRun != aEsd->GetRunNumber()) {
433     if (!Initialize(aEsd))
434       AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
435     if(fComputeBG) SetBIFactors(aEsd); // is this safe here?
436   }
437   const AliESDHeader* esdHeader = aEsd->GetHeader();
438   if (!esdHeader)
439   {
440     AliError("ESD Header could not be retrieved");
441     return kFALSE;
442   }
443   
444   // check event type; should be PHYSICS = 7 for data and 0 for MC
445   if (!fMC)
446   {
447     if (esdHeader->GetEventType() != 7)
448       return kFALSE;
449   }
450   else
451   {
452     if (esdHeader->GetEventType() != 0)
453       AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
454   }
455   
456   UInt_t accept = 0;
457     
458   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
459   for (Int_t i=0; i < count; i++)
460   {
461     const char* triggerClass = 0;
462     if (i < fCollTrigClasses.GetEntries())
463       triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
464     else
465       triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
466   
467     AliDebug(AliLog::kDebug+1, Form("Processing trigger class %s", triggerClass));
468   
469     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
470   
471     triggerAnalysis->FillTriggerClasses(aEsd);
472     
473     Int_t triggerLogic = 0; 
474     UInt_t singleTriggerResult = CheckTriggerClass(aEsd, triggerClass, triggerLogic);
475     
476     if (singleTriggerResult)
477     {
478       triggerAnalysis->FillHistograms(aEsd);
479   
480       Bool_t isBin0 = kFALSE;
481       if (fBin0CallBack != "") {
482           isBin0 = ((AliAnalysisTaskSE*)mgr->GetTask(fBin0CallBack.Data()))->IsEventInBinZero();
483       } else if (fBin0CallBackPointer) {
484           isBin0 = (*fBin0CallBackPointer)(aEsd);
485       }
486       
487       // ---->
488       // FIXME: this stuff is still used to fill the stat
489       // tables. Decide wethere we are switching to a new stat table
490       // (with all AliTriggerAnalysis tokens? Only we the tokens
491       // actually used in the selection?) and clean up
492
493       // hardware trigger
494       Int_t fastORHW   = triggerAnalysis->EvaluateTrigger(aEsd, AliTriggerAnalysis::kSPDGFO); // SPD number of chips from trigger bits (!)
495       //      Int_t fastORHWL1 = triggerAnalysis->EvaluateTrigger(aEsd, AliTriggerAnalysis::kSPDGFOL1); // SPD number of chips from trigger bits in second layer (!)
496       Bool_t v0AHW     = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, AliTriggerAnalysis::kV0A);// should replay hw trigger
497       Bool_t v0CHW     = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, AliTriggerAnalysis::kV0C);// should replay hw trigger
498       
499       // offline trigger
500       Int_t fastOROffline   = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kSPDGFO)); // SPD number of chips from clusters (!)
501       Int_t fastOROfflineL1 = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kSPDGFOL1)); // SPD number of chips from clusters in second layer (!)
502       Bool_t v0A       = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kV0A));
503       Bool_t v0C       = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kV0C));
504       Bool_t v0ABG = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kV0ABG));
505       Bool_t v0CBG = fSkipV0 ? 0 : triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kV0CBG));
506       Bool_t v0BG = v0ABG || v0CBG;
507
508       // fmd
509       // Bool_t fmdA = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kFMDA));
510       // Bool_t fmdC = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kFMDC));
511       // Bool_t fmd  = fmdA || fmdC;
512     
513       // SSD
514       //Int_t ssdClusters = triggerAnalysis->SSDClusters(aEsd);
515       
516       // ZDC
517       // Bool_t zdcA = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kZDCA);
518       // Bool_t zdcC = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kZDCC);
519       Bool_t zdcA    = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kZDCTDCA));
520       Bool_t zdcC    = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kZDCTDCC));
521       Bool_t zdcTime = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kZDCTime));
522
523       Bool_t laserCut = triggerAnalysis->EvaluateTrigger(aEsd, (AliTriggerAnalysis::Trigger) (AliTriggerAnalysis::kOfflineFlag | AliTriggerAnalysis::kTPCLaserWarmUp));
524
525       // Some "macros"
526       Bool_t mb1 = (fastOROffline > 0 || v0A || v0C) && (!v0BG);
527       Bool_t mb1prime = (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C) ) && (!v0BG);
528
529       // Background rejection
530       Bool_t bgID = kFALSE;
531       if (fBackgroundIdentification)
532         bgID = ! fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd));
533         
534       
535       /*Int_t ntrig = fastOROffline; // any 2 hits
536       if(v0A)              ntrig += 1;
537       if(v0C)              ntrig += 1; //v0C alone is enough
538       if(fmd)              ntrig += 1;
539       if(ssdClusters>1)    ntrig += 1;*/
540
541       
542       // <---
543
544       TString triggerLogicOnline  = fPSOADB->GetHardwareTrigger(triggerLogic);
545       TString triggerLogicOffline = fPSOADB->GetOfflineTrigger(triggerLogic);
546
547       AliDebug(AliLog::kDebug, Form("Triggers from OADB [0x%x][%d][%s][%s]",singleTriggerResult,AliOADBPhysicsSelection::GetActiveBit(singleTriggerResult),triggerLogicOffline.Data(),triggerLogicOnline.Data()));
548
549       // replay hardware trigger (should only remove events for MC)
550       Bool_t onlineTrigger  = EvaluateTriggerLogic(aEsd, triggerAnalysis, triggerLogicOnline, kFALSE);
551       // offline selection
552       Bool_t offlineTrigger = EvaluateTriggerLogic(aEsd, triggerAnalysis, triggerLogicOffline, kTRUE);
553       
554       // Printf("%s %s", triggerLogicOnline.Data(), triggerLogicOffline.Data());
555       // Printf("%d %d", onlineTrigger, offlineTrigger);
556            
557
558       // Fill trigger pattern histo
559       Int_t tpatt = 0;
560       if (fastORHW>0) tpatt+=1;
561       if (v0AHW)      tpatt+=2;
562       if (v0CHW)      tpatt+=4;
563       fHistTriggerPattern->Fill( tpatt );
564
565       // fill statistics and return decision
566       const Int_t nHistStat = 2;
567       for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
568         if (iHistStat == kStatIdxBin0 && !isBin0) continue; // skip the filling of bin0 stats if the event is not in the bin0
569       
570         fHistStatistics[iHistStat]->Fill(kStatTriggerClass, i);
571         if(iHistStat == kStatIdxAll) fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i); // Fill only for all (avoid double counting)
572
573         // We fill the rest only if hw trigger is ok
574         if (!onlineTrigger)
575           {
576             AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
577             continue;
578           } else {       
579           fHistStatistics[iHistStat]->Fill(kStatHWTrig, i);
580         }
581       
582
583         // v0 BG stats
584         if (v0ABG)
585           fHistStatistics[iHistStat]->Fill(kStatV0ABG, i);
586         if (v0CBG)
587           fHistStatistics[iHistStat]->Fill(kStatV0CBG, i);
588
589         // mb 1
590         if (mb1)
591           fHistStatistics[iHistStat]->Fill(kStatMB1, i);
592
593         if (mb1prime)
594           fHistStatistics[iHistStat]->Fill(kStatMB1Prime, i);
595
596         if (laserCut)
597           fHistStatistics[iHistStat]->Fill(kStatLaserCut, i);
598
599         //if(ntrig >= 2 && !v0BG) 
600         //  fHistStatistics[iHistStat]->Fill(kStatAny2Hits, i);
601
602         if (fastOROffline > 0)
603           fHistStatistics[iHistStat]->Fill(kStatFO1, i);
604         if (fastOROffline > 1)
605           fHistStatistics[iHistStat]->Fill(kStatFO2, i);
606         if (fastOROfflineL1 > 1)
607           fHistStatistics[iHistStat]->Fill(kStatFO2L1, i);
608         
609         if (v0A)
610           fHistStatistics[iHistStat]->Fill(kStatV0A, i);
611         if (v0C)
612           fHistStatistics[iHistStat]->Fill(kStatV0C, i);
613           
614         if (zdcA)
615           fHistStatistics[iHistStat]->Fill(kStatZDCA, i);
616         if (zdcC)
617           fHistStatistics[iHistStat]->Fill(kStatZDCC, i);
618         if (zdcA && zdcC)
619           fHistStatistics[iHistStat]->Fill(kStatZDCAC, i);
620
621         if (zdcTime)
622           fHistStatistics[iHistStat]->Fill(kStatZDCTime, i);
623   
624         if (v0A && v0C && !v0BG && (!bgID && fIsPP))
625           fHistStatistics[iHistStat]->Fill(kStatV0, i);
626
627
628
629         // FIXME: check lines below
630         if ( offlineTrigger )
631           {
632             if (!v0BG || fSkipV0)
633               {
634                 if (!v0BG) fHistStatistics[iHistStat]->Fill(kStatOffline, i);
635       
636                 if (fBackgroundIdentification && bgID && fIsPP)
637                   {
638                     AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
639                     fHistStatistics[iHistStat]->Fill(kStatBG, i);
640                   }
641                 else
642                   {
643                     AliDebug(AliLog::kDebug, Form("Accepted event for histograms with trigger logic %d", triggerLogic));
644             
645                     fHistStatistics[iHistStat]->Fill(kStatAccepted, i);
646                     // if(iHistStat == kStatIdxAll) fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i); // Fill only for all (avoid double counting)
647                     if((i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection) && (iHistStat==kStatIdxAll))
648                       accept |= singleTriggerResult; // only set for "all" (should not really matter)
649                   }
650               }
651             else
652               AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
653           }
654         else
655           AliDebug(AliLog::kDebug, Form("Rejecting event because trigger logic %d is not fulfilled", triggerLogic));
656       }
657     }
658   }
659  
660   if (accept)
661     AliDebug(AliLog::kDebug, Form("Accepted event as collision candidate with bit mask %d", accept));
662   
663   return accept;
664 }
665
666
667 // const char * AliPhysicsSelection::GetFillingScheme(UInt_t runNumber)  {
668
669 //   if(fMC) return "MC";
670
671 //   if      (runNumber >= 104065 && runNumber <= 104160) {
672 //     return "4x4a";
673 //   } 
674 //   else if (runNumber >= 104315 && runNumber <= 104321) {
675 //     return "4x4a*";
676 //   }
677 //   else if (runNumber >= 104792 && runNumber <= 104803) {
678 //     return "4x4b";
679 //   }
680 //   else if (runNumber >= 104824 && runNumber <= 104892) {
681 //     return "4x4c";
682 //   }
683 //   else if (runNumber == 105143 || runNumber == 105160) {
684 //     return "16x16a";
685 //   }
686 //   else if (runNumber >= 105256 && runNumber <= 105268) {
687 //     return "4x4c";
688 //   } 
689 //   else if (runNumber >= 114786 && runNumber <= 116684) {
690 //     return "Single_2b_1_1_1";
691 //   }
692 //   else if (runNumber >= 117048 && runNumber <= 117120) {
693 //     return "Single_3b_2_2_2";
694 //   }
695 //   else if (runNumber >= 117220 && runNumber <= 119163) {
696 //     return "Single_2b_1_1_1";
697 //   }
698 //   else if (runNumber >= 119837 && runNumber <= 119862) {
699 //     return "Single_4b_2_2_2";
700 //   }
701 //   else if (runNumber >= 119902 && runNumber <= 120691) {
702 //     return "Single_6b_3_3_3";
703 //   }
704 //   else if (runNumber >= 120741 && runNumber <= 122375) {
705 //     return "Single_13b_8_8_8";
706 //   }
707 //   else if (runNumber >= 130148 && runNumber <= 130375) {
708 //     return "125n_48b_36_16_36";
709 //   } 
710 //   else if (runNumber >= 130601 && runNumber <= 130640) {
711 //     return "1000ns_50b_35_14_35";
712 //   }
713 //   else {
714 //     AliError(Form("Unknown filling scheme (run %d)", runNumber));
715 //   }
716
717 //   return "Unknown";
718 // }
719
720 // const char * AliPhysicsSelection::GetBXIDs(UInt_t runNumber, const char * trigger)  {
721
722 //   if (!fUseBXNumbers || fMC) return "";
723
724 //   if      (runNumber >= 104065 && runNumber <= 104160) {
725 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2128 #3019";
726 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #346 #3465";
727 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #1680";
728 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
729 //     //    else AliError(Form("Unknown trigger: %s", trigger));
730 //   } 
731 //   else if (runNumber >= 104315 && runNumber <= 104321) {
732 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2000 #2891";
733 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #218 #3337";
734 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1106 #1552";
735 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
736 //     //    else AliError(Form("Unknown trigger: %s", trigger));
737 //   }
738 //   else if (runNumber >= 104792 && runNumber <= 104803) {
739 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2228 #3119";
740 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
741 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #769";
742 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
743 //     //    else AliError(Form("Unknown trigger: %s", trigger));
744 //   }
745 //   else if (runNumber >= 104824 && runNumber <= 104892) {
746 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3119 #769";
747 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
748 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #2228";
749 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
750 //     //    else AliError(Form("Unknown trigger: %s", trigger));
751 //   }
752 //   else if (runNumber == 105143 || runNumber == 105160) {
753 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #1337 #1418 #2228 #2309 #3119 #3200 #446 #527";
754 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #1580  #1742  #1904  #2066  #2630  #2792  #2954  #3362";
755 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "  #845  #1007  #1169   #1577 #3359 #3521 #119  #281 ";
756 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
757 //     //    else AliError(Form("Unknown trigger: %s", trigger));
758 //   }
759 //   else if (runNumber >= 105256 && runNumber <= 105268) {
760 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3019 #669";
761 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2454 #346";
762 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #2128";
763 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
764 //     //    else AliError(Form("Unknown trigger: %s", trigger));
765 //   } else if (runNumber >= 114786 && runNumber <= 116684) { // 7 TeV 2010, assume always the same filling scheme
766 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #346";
767 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2131";
768 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #3019";
769 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #1238";
770 //     //    else AliError(Form("Unknown trigger: %s", trigger));
771 //   }
772 //   else if (runNumber >= 117048 && runNumber <= 117120) {
773 //     //    return "Single_3b_2_2_2";
774 //    if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #1240 ";
775 //    else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131 ";
776 //    else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019 ";
777 //    else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238";
778 //    //   else AliError(Form("Unknown trigger: %s", trigger));
779
780 //   }
781 //   else if ((runNumber >= 117220 && runNumber <= 118555) || (runNumber >= 118784 && runNumber <= 119163)) 
782 //   {
783 //     //    return "Single_2b_1_1_1";
784 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346 ";
785 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131 ";
786 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019 ";
787 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238 ";
788 //     //    else AliError(Form("Unknown trigger: %s", trigger));                                                   
789 //   }
790 //   else if (runNumber >= 118556 && runNumber <= 118783) {
791 //     //    return "Single_2b_1_1_1"; 
792 //     // same as previous but was misaligned by 1 BX in fill 1069
793 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #345 ";
794 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2130 ";
795 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3018 ";
796 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238 ";
797 //     //    else AliError(Form("Unknown trigger: %s", trigger));                                                   
798 //   }
799 //   else if (runNumber >= 119837 && runNumber <= 119862) {
800 //     //    return "Single_4b_2_2_2";
801 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #669  #3019 ";
802 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #346  #2454 ";
803 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #1234  #2128 ";
804 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1681 #3463";
805 //     //    else AliError(Form("Unknown trigger: %s", trigger));
806
807 //   }
808 //   else if (runNumber >= 119902 && runNumber <= 120691) {
809 //     //    return "Single_6b_3_3_3";
810 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #546  #746 ";
811 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131  #2331  #2531 ";
812 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019  #3219  #3419 ";
813 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1296 #1670";
814 //     //    else AliError(Form("Unknown trigger: %s", trigger));
815 //   }
816 //   else if (runNumber >= 120741 && runNumber <= 122375) {
817 //     //    return "Single_13b_8_8_8";
818 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #446  #546  #646  #1240  #1340  #1440  #1540 ";
819 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #946  #2131  #2231  #2331  #2431 ";
820 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019  #3119  #3219  #3319  #3519 ";
821 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1835 #2726";
822 //     //    else AliError(Form("Unknown trigger: %s", trigger));
823     
824 //   } 
825 //   else if (runNumber >= 130148 && runNumber <= 130375) {
826 //     TString triggerString = trigger;
827 //     static TString returnString = " ";
828 //     returnString = "";
829 //     if (triggerString.Contains("B")) returnString += "   #346  #396  #446  #496  #546  #596  #646  #696  #1240  #1290  #1340  #1390  #1440  #1490  #1540  #1590 ";
830 //     if (triggerString.Contains("A")) returnString += "   #755  #805  #855  #905  #955  #1005  #1799  #1849  #1899  #2131  #2181  #2231  #2281  #2331  #2381  #2431  #2481  #2531  #2581  #2631  #2846  #3016  #3066  #3116  #3166  #3216  #3266  #3316  #3366  #3425  #3475  #3525 ";
831 //     if (triggerString.Contains("C")) returnString += "   #3019  #3069  #3119  #3169  #3219  #3269  #3319  #3369  #14  #64  #114  #746  #796  #846  #908  #958  #1008  #1640  #1690  #1740  #2055  #2125  #2175  #2225  #2275  #2325  #2375  #2425  #2475  #2534  #2584  #2634 ";
832 //     // Printf("0x%x",returnString.Data());
833 //     // Printf("%s",returnString.Data());
834 //     return returnString.Data();
835 //   } 
836 //   else if (runNumber >= 130601 && runNumber <= 130640) {
837 //     TString triggerString = trigger;
838 //     static TString returnString = " ";
839 //     returnString = "";
840 //     if (triggerString.Contains("B")) returnString += "  #346  #386  #426  #466  #506  #546  #586  #1240  #1280  #1320  #1360  #1400  #1440  #1480 ";
841 //     if (triggerString.Contains("A")) returnString += "  #626  #666  #706  #746  #786  #826  #866  #1520  #1560  #1600  #1640  #1680  #1720  #1760  #2076  #2131  #2171  #2211  #2251  #2291  #2331  #2371  #2414  #2454  #2494  #2534  #2574  #2614  #2654  #2694  #2734  #2774  #2814 "; //#2854  #2894  #2934 not present in this run
842 //     if (triggerString.Contains("C")) returnString += "  #3019  #3059  #3099  #3139  #3179  #3219  #3259  #3299  #3339  #3379  #3419  #3459  #3499  #3539  #115  #629  #669  #709  #749  #789  #829  #869  #909  #949  #989  #1029  #1069  #1109  #1149  #1523  #1563  #1603  #1643 "; //#1683  #1723  #1763 not present in this run
843 //     return returnString.Data();
844 //   }
845
846 //   else {
847 //     AliWarning(Form("Unknown run %d, using all BXs!",runNumber));
848 //   }
849
850 //   return "";
851 // }
852
853 Bool_t AliPhysicsSelection::Initialize(const AliESDEvent* aEsd)
854 {
855   // initializes the object for the given ESD
856   
857   AliInfo(Form("Initializing for beam type: %s", aEsd->GetESDRun()->GetBeamType()));
858   fIsPP = kTRUE;
859   if (strcmp(aEsd->GetESDRun()->GetBeamType(), "Pb-Pb") == 0)
860     fIsPP = kFALSE;
861
862   return Initialize(aEsd->GetRunNumber());
863 }
864
865 Bool_t AliPhysicsSelection::Initialize(Int_t runNumber)
866 {
867   // initializes the object for the given run  
868
869
870   Bool_t oldStatus = TH1::AddDirectoryStatus();
871   TH1::AddDirectory(kFALSE);
872
873   /// Open OADB file and fetch OADB objects
874   TString oadbfilename = AliPhysicsSelection::GetOADBFileName();
875
876   TFile * foadb = TFile::Open(oadbfilename);
877   if(!foadb->IsOpen()) AliFatal(Form("Cannot open OADB file %s", oadbfilename.Data()));
878
879
880   if(!fPSOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
881     AliInfo("Using Standard OADB");
882     AliOADBContainer * psContainer = (AliOADBContainer*) foadb->Get("physSel");
883     if (!psContainer) AliFatal("Cannot fetch OADB container for Physics selection");
884     fPSOADB = (AliOADBPhysicsSelection*) psContainer->GetObject(runNumber, fIsPP ? "oadbDefaultPP" : "oadbDefaultPbPb");
885     if (!fPSOADB) AliFatal(Form("Cannot find physics selection object for run %d", runNumber));
886   } else {
887     AliInfo("Using Custom OADB");
888   }
889   if(!fFillOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
890     AliOADBContainer * fillContainer = (AliOADBContainer*) foadb->Get("fillScheme");
891     if (!fillContainer) AliFatal("Cannot fetch OADB container for filling scheme");
892     fFillOADB = (AliOADBFillingScheme*) fillContainer->GetObject(runNumber, "Default");
893     if (!fFillOADB) AliFatal(Form("Cannot find  filling scheme object for run %d", runNumber));
894   }
895   if(!fTriggerOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
896     AliOADBContainer * triggerContainer = (AliOADBContainer*) foadb->Get("trigAnalysis");
897     if (!triggerContainer) AliFatal("Cannot fetch OADB container for trigger analysis");
898     fTriggerOADB = (AliOADBTriggerAnalysis*) triggerContainer->GetObject(runNumber, "Default");
899     fTriggerOADB->Print();
900     if (!fTriggerOADB) AliFatal(Form("Cannot find  trigger analysis object for run %d", runNumber));
901   }
902
903   
904   if(!fBin0CallBack) 
905     AliError("Bin0 Callback not set: will not fill the statistics for the bin 0");
906
907   if (fMC) {
908     // override BX and bg options in case of MC
909     fComputeBG    = kFALSE;
910     fUseBXNumbers = kFALSE;
911   }
912
913   // FIXME: think how to implement this check with the OADB, trigger scheme is not a int number here 
914   // Int_t triggerScheme = GetTriggerScheme(runNumber);
915   // if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
916   //   AliFatal("Processing several runs with different trigger schemes is not supported");
917   
918   if(fComputeBG && fCurrentRun != -1 && fCurrentRun != runNumber) 
919     AliFatal("Cannot process several runs because BG computation is requested");
920
921   if(fComputeBG && !fUseBXNumbers) 
922     AliFatal("Cannot compute BG if BX numbers are not used");
923   
924   if(fUseBXNumbers && fFillingScheme != "" && fFillingScheme != fFillOADB->GetFillingSchemeName())
925     AliFatal("Cannot process runs with different filling scheme if usage of BX numbers is requested");
926
927   fFillingScheme      = fFillOADB->GetFillingSchemeName();
928
929   AliInfo(Form("Initializing for run %d", runNumber));
930   
931   // initialize first time?
932   if (fCurrentRun == -1)
933   {
934     const UInt_t ntriggerBits = fPSOADB->GetNTriggerBits(); // FIXME!
935     for(UInt_t ibit = 0; ibit < ntriggerBits; ibit++){
936       
937       TIterator * collIter = fPSOADB->GetCollTrigClass(ibit)->MakeIterator();
938       TIterator * bgIter   = fPSOADB->GetBGTrigClass(ibit)->MakeIterator();
939       TObjString * obj = 0;
940       while((obj = (TObjString*) collIter->Next())){
941         if (obj->String() != "") {
942           fCollTrigClasses.Add(new TObjString(GetTriggerString(obj)));
943         }
944       }
945       // BG classes only make sense for real data
946       if(!fMC) {
947         obj = 0 ;
948         while((obj = (TObjString*) bgIter->Next())){
949           if (obj->String() != "") {
950             fBGTrigClasses.Add(new TObjString(GetTriggerString(obj)));
951           }
952         }
953       }
954
955     }
956     // not sure how to handle this in the case of > x cuts
957     // // Book the token histo with the tokens actually used here!
958     // // 1. Cache the tokens
959     // for(UInt_t ibit = 0; ibit < ntriggerBits; ibit++){
960     //  EvaluateTriggerLogic(0,0, fPSOADB->GetHardwareTrigger(ibit), 0);
961     //  EvaluateTriggerLogic(0,0, fPSOADB->GetOfflineTrigger(ibit),  1);
962     // }
963     // // 2. Book the histogram
964     // if (!fHistStatisticsTokens) {
965     //  Int_t ntokens = fCashedTokens->GetEntries();
966     //  Int_t count = fCollTrigClasses->GetEntries() + fBGTrigClasses->GetEntries();
967     //  fHistStatisticsTokens = new TH2F("fHistStatisticsTokens", "fHistStatisticsTokens", ntokens + 2, 0.5, ntokens+2+0.5, count, -0.5, -0.5 + count);
968         
969     //  Int_t nrow = 0;
970     //  fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,  "Trigger class");
971     //  for(Int_t itoken = 0; itoken < ntoken; itoken++){
972     //    TParameter<Int_t> * param = fCashedTokens->At(itoken);
973     //    fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,  param->GetName());      
974     //  }
975         
976     //  fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,      "Accepted");
977
978     // }
979
980     
981     
982     // TODO: 
983     // Add a new statistics histo containing only the tokens actually used in the selection
984
985     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
986     
987     for (Int_t i=0; i<count; i++)
988     {
989       
990       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
991       triggerAnalysis->SetAnalyzeMC(fMC);
992       triggerAnalysis->EnableHistograms();
993       triggerAnalysis->SetSPDGFOThreshhold(1);
994       triggerAnalysis->SetDoFMD(kFALSE);
995       triggerAnalysis->SetCorrZDCCutParams(fTriggerOADB->GetZDCCutRefSumCorr(),
996                                            fTriggerOADB->GetZDCCutRefDeltaCorr(), 
997                                            fTriggerOADB->GetZDCCutSigmaSumCorr(),
998                                            fTriggerOADB->GetZDCCutSigmaDeltaCorr());
999       fTriggerAnalysis.Add(triggerAnalysis);
1000     }
1001
1002     
1003     // TODO: shall I really delete this?
1004     if (fHistStatistics[0])
1005       delete fHistStatistics[0];
1006     if (fHistStatistics[1])
1007       delete fHistStatistics[1];
1008   
1009     fHistStatistics[kStatIdxBin0] = BookHistStatistics("_Bin0");
1010     fHistStatistics[kStatIdxAll]  = BookHistStatistics("");
1011     
1012
1013     if (fHistBunchCrossing)
1014       delete fHistBunchCrossing;
1015   
1016     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
1017
1018     // TODO: remove fHistTriggerPattern
1019     if (fHistTriggerPattern)
1020       delete fHistTriggerPattern;
1021     
1022     const int ntrig=3;
1023     Int_t n = 1;
1024     const Int_t nbinTrig = TMath::Nint(TMath::Power(2,ntrig));
1025
1026     fHistTriggerPattern = new TH1F("fHistTriggerPattern", "Trigger pattern: FO + 2*v0A + 4*v0C", 
1027                                    nbinTrig, -0.5, nbinTrig-0.5);    
1028     fHistTriggerPattern->GetXaxis()->SetBinLabel(1,"NO TRIG");
1029     fHistTriggerPattern->GetXaxis()->SetBinLabel(2,"FO");
1030     fHistTriggerPattern->GetXaxis()->SetBinLabel(3,"v0A");
1031     fHistTriggerPattern->GetXaxis()->SetBinLabel(4,"FO & v0A");
1032     fHistTriggerPattern->GetXaxis()->SetBinLabel(5,"v0C");
1033     fHistTriggerPattern->GetXaxis()->SetBinLabel(6,"FO & v0C");
1034     fHistTriggerPattern->GetXaxis()->SetBinLabel(7,"v0A & v0C");
1035     fHistTriggerPattern->GetXaxis()->SetBinLabel(8,"FO & v0A & v0C");
1036
1037   
1038     n = 1;
1039     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1040     {
1041       
1042       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
1043       n++;
1044     }
1045     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1046     {
1047       
1048       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
1049       n++;
1050     }
1051
1052     
1053
1054   }
1055     
1056   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1057   for (Int_t i=0; i<count; i++)
1058   {
1059     
1060     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
1061   
1062     switch (runNumber)
1063     {
1064       case 104315:
1065       case 104316:
1066       case 104320:
1067       case 104321:
1068       case 104439:
1069         triggerAnalysis->SetV0TimeOffset(7.5);
1070         break;
1071       default:
1072         triggerAnalysis->SetV0TimeOffset(0);
1073     }
1074   }
1075     
1076   fCurrentRun = runNumber;
1077   
1078   TH1::AddDirectory(oldStatus);
1079   
1080   
1081   return kTRUE;
1082 }
1083
1084 TH2F * AliPhysicsSelection::BookHistStatistics(const char * tag) {
1085   // add 6 rows to count for the estimate of good, accidentals and
1086   // BG and the ratio of BG and accidentals to total +ratio goot to
1087   // first col + 2 for error on good.
1088   // TODO: Remember the the indexes of rows for the BG selection. Add new member fBGRows[] and use kStat as indexes
1089
1090   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1091 #ifdef VERBOSE_STAT
1092   Int_t extrarows = fComputeBG != 0 ? 11 : 0;
1093 #else
1094   Int_t extrarows = fComputeBG != 0 ? 6 : 0;
1095 #endif
1096   TH2F * h = new TH2F(Form("fHistStatistics%s",tag), Form("fHistStatistics - %s ;;",tag), kStatAccepted, 0.5, kStatAccepted+0.5, count+extrarows, -0.5, -0.5 + count+extrarows);
1097
1098   h->GetXaxis()->SetBinLabel(kStatTriggerClass,  "Trigger class");
1099   h->GetXaxis()->SetBinLabel(kStatHWTrig,        "Hardware trigger");
1100   h->GetXaxis()->SetBinLabel(kStatFO1,           "FO >= 1");
1101   h->GetXaxis()->SetBinLabel(kStatFO2,           "FO >= 2");
1102   h->GetXaxis()->SetBinLabel(kStatFO2L1,         "FO (L1) >= 2");
1103   h->GetXaxis()->SetBinLabel(kStatV0A,           "V0A");
1104   h->GetXaxis()->SetBinLabel(kStatV0C,           "V0C");
1105   h->GetXaxis()->SetBinLabel(kStatLaserCut,      "TPC Laser Wup Cut");
1106   h->GetXaxis()->SetBinLabel(kStatV0ABG,         "V0A BG");
1107   h->GetXaxis()->SetBinLabel(kStatV0CBG,         "V0C BG");
1108   h->GetXaxis()->SetBinLabel(kStatZDCA,          "ZDCA");
1109   h->GetXaxis()->SetBinLabel(kStatZDCC,          "ZDCC");
1110   h->GetXaxis()->SetBinLabel(kStatZDCAC,         "ZDCA & ZDCC");
1111   h->GetXaxis()->SetBinLabel(kStatZDCTime,       "ZDC Time Cut");
1112   h->GetXaxis()->SetBinLabel(kStatMB1,           "(FO >= 1 | V0A | V0C) & !V0 BG");
1113   h->GetXaxis()->SetBinLabel(kStatMB1Prime,      "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A &v0C) ) & !V0 BG");
1114   //h->GetXaxis()->SetBinLabel(kStatFO1AndV0,    "FO >= 1 & (V0A | V0C) & !V0 BG");
1115   h->GetXaxis()->SetBinLabel(kStatV0,            "V0A & V0C & !V0 BG & !BG ID");
1116   h->GetXaxis()->SetBinLabel(kStatOffline,       "Offline Trigger");
1117   //h->GetXaxis()->SetBinLabel(kStatAny2Hits,    "2 Hits & !V0 BG");
1118   h->GetXaxis()->SetBinLabel(kStatBG,            "Background identification");
1119   h->GetXaxis()->SetBinLabel(kStatAccepted,      "Accepted");
1120
1121   Int_t n = 1;
1122   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1123     {
1124       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
1125       n++;
1126     }
1127   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1128     {
1129       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
1130       n++;
1131     }
1132
1133   if(fComputeBG) {
1134     fBGStatOffset = n;
1135     h->GetYaxis()->SetBinLabel(n++, "All B");
1136     h->GetYaxis()->SetBinLabel(n++, "All A+C");
1137     h->GetYaxis()->SetBinLabel(n++, "All E");
1138     h->GetYaxis()->SetBinLabel(n++, Form("BG (A+C) (Mask [0x%x])", fComputeBG));
1139     h->GetYaxis()->SetBinLabel(n++, "ACC");
1140 #ifdef VERBOSE_STAT
1141     h->GetYaxis()->SetBinLabel(n++, "BG (A+C) %  (rel. to CINT1B)");
1142     h->GetYaxis()->SetBinLabel(n++, "ACC % (rel. to CINT1B)");
1143     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD %");
1144     h->GetYaxis()->SetBinLabel(n++, "GOOD % (rel. to 1st col)");
1145     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD");
1146 #endif
1147     h->GetYaxis()->SetBinLabel(n++, "GOOD");
1148   }
1149
1150   return h;
1151 }
1152
1153 void AliPhysicsSelection::Print(const Option_t *option) const
1154 {
1155   // print the configuration
1156   TString msg;
1157   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
1158   msg += Form("Configuration initialized for run %d (MC: %d):\n", fCurrentRun, fMC);
1159   
1160   Printf("Collision trigger classes:");
1161   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1162     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1163   
1164   Printf("Background trigger classes:");
1165   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1166     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
1167
1168   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
1169   
1170   if (triggerAnalysis)
1171   {
1172     if (triggerAnalysis->GetV0TimeOffset() > 0)
1173       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
1174     
1175     Printf("\nTotal available events:");
1176     
1177     triggerAnalysis->PrintTriggerClasses();
1178   }
1179   
1180   if (fHistStatistics[kStatIdxAll])
1181   {
1182     for (Int_t i=0; i<fCollTrigClasses.GetEntries(); i++)
1183     {
1184       Printf("\nSelection statistics for collision trigger %s:", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1185       msg += Form("\nSelection statistics for collision trigger %s:\n", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1186       
1187       Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(1, i+1));
1188       msg += Form("Total events with correct trigger class: %d\n", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(1, i+1));
1189       Printf("Selected collision candidates: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(fHistStatistics[kStatIdxAll]->GetXaxis()->FindBin("Accepted"), i+1));
1190       msg += Form("Selected collision candidates: %d\n", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(fHistStatistics[kStatIdxAll]->GetXaxis()->FindBin("Accepted"), i+1));
1191     }
1192   }
1193   
1194   if (fHistBunchCrossing)
1195   {
1196     Printf("\nBunch crossing statistics:");
1197     msg += "\nBunch crossing statistics:\n";
1198     
1199     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
1200     {
1201       TString str;
1202       str.Form("Trigger %s has events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
1203       
1204       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
1205         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
1206           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
1207              
1208       Printf("%s", str.Data());
1209       msg += str;
1210       msg += "\n";
1211     }
1212     
1213     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
1214     {
1215       Int_t countColl = 0;
1216       Int_t countBG = 0;
1217       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
1218       {
1219         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
1220         {
1221           if (fCollTrigClasses.FindObject(fHistBunchCrossing->GetYaxis()->GetBinLabel(i)))
1222             countColl++;
1223           if (fBGTrigClasses.FindObject(fHistBunchCrossing->GetYaxis()->GetBinLabel(i)))
1224             countBG++;
1225         }
1226       }
1227       if (countColl > 0 && countBG > 0)
1228         Printf("WARNING: Bunch crossing %d has collision and BG trigger classes active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
1229     }
1230   }
1231
1232   if (fUsingCustomClasses)        
1233     Printf("WARNING: Using custom trigger classes!");
1234   if (fSkipTriggerClassSelection) 
1235     Printf("WARNING: Skipping trigger class selection!");
1236   if (fSkipV0) 
1237     Printf("WARNING: Ignoring V0 information in selection");
1238   if(!fBin0CallBack) 
1239     Printf("WARNING: Callback not set: will not fill the statistics for the bin 0");
1240   TString opt(option);
1241   opt.ToUpper();
1242   if (opt == "STAT") {
1243      AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1244      if (mgr) mgr->AddStatisticsMsg(msg);
1245   }
1246 }
1247
1248 Long64_t AliPhysicsSelection::Merge(TCollection* list)
1249 {
1250   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1251   // PROOF).
1252   // Returns the number of merged objects (including this).
1253
1254   if (!list)
1255     return 0;
1256
1257   if (list->IsEmpty())
1258     return 1;
1259
1260   TIterator* iter = list->MakeIterator();
1261   TObject* obj;
1262   
1263   // collections of all histograms
1264   const Int_t nHists = 9;
1265   TList collections[nHists];
1266
1267   Int_t count = 0;
1268   while ((obj = iter->Next())) {
1269
1270     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
1271     if (entry == 0) 
1272       continue;
1273     // Update run number. If this one is not initialized (-1) take the one from 
1274     // the next physics selection to be merged with. In case of 2 different run
1275     // numbers issue a warning (should physics selections from different runs be 
1276     // merged together) A.G.
1277     Int_t currentRun = entry->GetCurrentRun();
1278     // Nothing to merge with since run number was not initialized.
1279     if (currentRun < 0) continue;
1280     if (fCurrentRun < 0) 
1281     {
1282       fCurrentRun = currentRun;
1283       fMC = entry->fMC;
1284       for (Int_t i=0; i<entry->fCollTrigClasses.GetEntries(); i++)
1285         fCollTrigClasses.Add(entry->fCollTrigClasses.At(i)->Clone());
1286       for (Int_t i=0; i<entry->fBGTrigClasses.GetEntries(); i++)
1287         fBGTrigClasses.Add(entry->fBGTrigClasses.At(i)->Clone());
1288     }
1289     if (fCurrentRun != currentRun)
1290        AliWarning(Form("Current run %d not matching the one to be merged with %d", fCurrentRun, currentRun));
1291
1292     // With the same strategy update fBGStatOffset
1293     Int_t bgstatoffset = entry->GetBGStatOffset();
1294     
1295     // Nothing to merge with since BG was not initialized.
1296     if (!(bgstatoffset < 0)){
1297       if (fBGStatOffset < 0) 
1298         {
1299           fBGStatOffset = bgstatoffset;
1300         }
1301     }
1302     if (fBGStatOffset != bgstatoffset)
1303       AliWarning(Form("Current run %d not matching the one to be merged with %d", fBGStatOffset, bgstatoffset));
1304     
1305
1306     
1307     // Merge the OADBs (Take just the first instance you find
1308     if (!fPSOADB && entry->GetOADBPhysicsSelection()) {
1309       fPSOADB = (AliOADBPhysicsSelection*) entry->GetOADBPhysicsSelection()->Clone();
1310     }
1311     if (!fFillOADB && entry->GetOADBFillingScheme()){
1312       fFillOADB = (AliOADBFillingScheme*) entry->GetOADBFillingScheme()->Clone();
1313     }
1314
1315     if (entry->fTriggerAnalysis.GetEntries() > 0)
1316       collections[0].Add(&(entry->fTriggerAnalysis));
1317     if (entry->fHistStatistics[0])
1318       collections[1].Add(entry->fHistStatistics[0]);
1319     if (entry->fHistStatistics[1])
1320       collections[2].Add(entry->fHistStatistics[1]);
1321     // if (entry->fHistStatisticsTokens)
1322     //   collections[3].Add(entry->fHistStatisticsTokens);
1323     if (entry->fHistBunchCrossing)
1324       collections[3].Add(entry->fHistBunchCrossing);
1325     if (entry->fHistTriggerPattern)
1326       collections[4].Add(entry->fHistTriggerPattern);
1327     if (entry->fBackgroundIdentification)
1328       collections[5].Add(entry->fBackgroundIdentification);
1329
1330     count++;
1331   }
1332
1333   if (fTriggerAnalysis.GetEntries() == 0 && collections[0].GetEntries() > 0)
1334   {
1335     TList* firstList = (TList*) collections[0].First();
1336     for (Int_t i=0; i<firstList->GetEntries(); i++)
1337       fTriggerAnalysis.Add(firstList->At(i)->Clone());
1338     
1339     collections[0].RemoveAt(0);
1340   }
1341   fTriggerAnalysis.Merge(&collections[0]);
1342   
1343   // if this instance is empty (not initialized) nothing would be merged here --> copy first entry
1344   if (!fHistStatistics[0] && collections[1].GetEntries() > 0)
1345   {
1346     fHistStatistics[0] = (TH2F*) collections[1].First()->Clone();
1347     collections[1].RemoveAt(0);
1348   }
1349   if (fHistStatistics[0])
1350     fHistStatistics[0]->Merge(&collections[1]);
1351     
1352   if (!fHistStatistics[1] && collections[2].GetEntries() > 0)
1353   {
1354     fHistStatistics[1] = (TH2F*) collections[2].First()->Clone();
1355     collections[2].RemoveAt(0);
1356   }
1357   if (fHistStatistics[1])
1358     fHistStatistics[1]->Merge(&collections[2]);
1359     
1360   if (!fHistBunchCrossing && collections[3].GetEntries() > 0)
1361   {
1362     fHistBunchCrossing = (TH2F*) collections[3].First()->Clone();
1363     collections[3].RemoveAt(0);
1364   }
1365   if (fHistBunchCrossing)
1366     fHistBunchCrossing->Merge(&collections[3]);
1367   
1368   if (!fHistTriggerPattern && collections[4].GetEntries() > 0)
1369   {
1370     fHistTriggerPattern = (TH1F*) collections[4].First()->Clone();
1371     collections[4].RemoveAt(0);
1372   }
1373   if (fHistTriggerPattern)
1374     fHistTriggerPattern->Merge(&collections[4]);
1375   
1376   if (!fBackgroundIdentification && collections[5].GetEntries() > 0)
1377   {
1378     fBackgroundIdentification = (AliAnalysisCuts*) collections[5].First()->Clone();
1379     collections[5].RemoveAt(0);
1380   }
1381   if (fBackgroundIdentification)
1382     fBackgroundIdentification->Merge(&collections[5]);
1383   
1384   delete iter;
1385
1386   return count+1;
1387 }
1388
1389 void AliPhysicsSelection::SaveHistograms(const char* folder) 
1390 {
1391   // write histograms to current directory
1392   
1393   if (!fHistStatistics[0] || !fHistStatistics[1])
1394     return;
1395     
1396   if (folder)
1397   {
1398     gDirectory->mkdir(folder);
1399     gDirectory->cd(folder);
1400   }
1401   
1402
1403   // Fill the last rows of fHistStatistics before saving
1404   if (fComputeBG) {      
1405     AliInfo("BG estimate assumes that for a given run you only have A and C triggers separately or"
1406             " toghether as a AC class! Make sure this assumption holds in your case"); 
1407     
1408     // use an anum for the different trigger classes, to make loops easier to read
1409     enum {kClassB =0 , kClassA, kClassC, kClassAC, kClassE, kNClasses};
1410     const char * classFlags[] = {"B", "A", "C", "AC", "E"}; // labels
1411     
1412     UInt_t * rows[kNClasses] = {0}; // Array of matching rows
1413     Int_t nrows[kNClasses] = {0};
1414     // Get rows matching the requested trigger bits for all trigger classes
1415     for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1416       nrows[iTrigClass] = GetStatRow(classFlags[iTrigClass],fComputeBG,&rows[iTrigClass]);      
1417     }
1418     
1419     // 0. Determine the ratios of triggers E/B, A/B, C/B from the stat histogram
1420     // Those are used to rescale the different classes to the same number of bx ids
1421     // TODO: pass names of the rows for B, CA and E and look names of the rows. How do I handle the case in which both AC are in the same row?         
1422     Int_t nBXIds[kNClasses] = {0};
1423       //      cout <<"Computing BG:" << endl;
1424     
1425     for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1426       for(Int_t irow = 0; irow < nrows[iTrigClass]; irow++) {
1427         if(irow==0) cout << "- Class " << classFlags[iTrigClass] << endl;   
1428         for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++) {
1429           if (fHistBunchCrossing->GetBinContent(j, rows[iTrigClass][irow]) > 0) {
1430             nBXIds[iTrigClass]++;        
1431           }
1432         }
1433         if(nBXIds[iTrigClass]>0) cout << "   Using row " << rows[iTrigClass][irow] <<  ": " 
1434                                       << fHistBunchCrossing->GetYaxis()->GetBinLabel(rows[iTrigClass][irow]) 
1435                                       << " (nBXID "<< nBXIds[iTrigClass] << ")"<< endl;
1436
1437       }
1438
1439     }
1440
1441     Float_t ratioToB[kNClasses];
1442     ratioToB[kClassE]  = nBXIds[kClassE]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassE]   : 0;
1443     ratioToB[kClassA]  = nBXIds[kClassA]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassA]   : 0;
1444     ratioToB[kClassC]  = nBXIds[kClassC]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassC]   : 0;
1445     ratioToB[kClassAC] = nBXIds[kClassAC] >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassAC]  : 0;
1446     Printf("Ratio between the BX ids in the different trigger classes:");
1447     Printf("  B/E  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassE], ratioToB[kClassE] );
1448     Printf("  B/A  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassA], ratioToB[kClassA] );
1449     Printf("  B/C  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassC], ratioToB[kClassC] );
1450     Printf("  B/AC = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassAC],ratioToB[kClassAC]);
1451     Int_t nHistStat = 2;
1452
1453     // 1. loop over all cols
1454     for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
1455       Int_t ncol = fHistStatistics[iHistStat]->GetNbinsX();
1456       Float_t good1 = 0;      
1457       for(Int_t icol = 1; icol <= ncol; icol++) {
1458         Int_t nEvents[kNClasses] = {0}; // number of events should be reset at every column
1459         // For all trigger classes, add up over row matching trigger mask (as selected before)
1460         for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1461           for(Int_t irow = 0; irow < nrows[iTrigClass]; irow++) {         
1462             nEvents[iTrigClass] += (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,rows[iTrigClass][irow]);      
1463           }
1464           //        cout << "Events " << classFlags[iTrigClass] << " ("<<icol<<") " << nEvents[iTrigClass] << endl;         
1465         }
1466         if (nEvents[kClassB]>0) {
1467           Float_t acc  = ratioToB[kClassE]*nEvents[kClassE]; 
1468           Double_t accErr = TMath::Sqrt(ratioToB[kClassE]*ratioToB[kClassE]*nEvents[kClassE]);
1469           //      Int_t bg   = cint1A + cint1C - 2*acc;
1470             
1471           // If intensity measurements are available, they already
1472           // contain the scaling for BX ratios, so we reset the
1473           // ratioToB entries
1474           if(icol == 1) {
1475             if(fBIFactorAC > 0 || fBIFactorA > 0 || fBIFactorC > 0) {
1476               if (fBIFactorAC <= 0 && (fBIFactorA <= 0 || fBIFactorC <= 0)) {
1477                 AliError("Not all intensities set!, assuming equal intensities");
1478                 fBIFactorA  = 1;
1479                 fBIFactorC  = 1;
1480                 fBIFactorAC = 1;
1481               } else {
1482                 AliInfo("Using ratio of number of bunch crossing embedded in the intensity measurements");
1483                 ratioToB[kClassA]  = ratioToB[kClassA]  >0 ? 1  : 0;
1484                 ratioToB[kClassC]  = ratioToB[kClassC]  >0 ? 1  : 0;
1485                 ratioToB[kClassAC] = ratioToB[kClassAC] >0 ? 1  : 0;
1486                 AliInfo(Form(" - BI Factor A:  %f",  fBIFactorA ));
1487                 AliInfo(Form(" - BI Factor C:  %f",  fBIFactorC ));
1488                 AliInfo(Form(" - BI Factor AC: %f",  fBIFactorAC ));
1489                 
1490               }
1491             } else {
1492               AliWarning("Intensities not set!, assuming equal intensities");
1493               fBIFactorA  = 1;
1494               fBIFactorC  = 1;
1495               fBIFactorAC = 1;
1496             }
1497           }
1498           // Assuming that for a given class the triggers are either recorded as A+C or AC
1499           Float_t bg  = nEvents[kClassAC] > 0 ?
1500             fBIFactorAC*(ratioToB[kClassAC]*nEvents[kClassAC] - 2*acc):
1501             fBIFactorA* (ratioToB[kClassA]*nEvents[kClassA]-acc) + 
1502             fBIFactorC* (ratioToB[kClassC]*nEvents[kClassC]-acc) ;
1503
1504           // cout << "-----------------------" << endl;
1505           // cout << "Factors: " << fBIFactorA << " " << fBIFactorC << " " << fBIFactorAC << endl;
1506           // cout << "Ratios: "  << ratioToB[kClassA] << " " << ratioToB[kClassC] << " " << ratioToB[kClassAC] << endl;
1507           // cout << "Evts:   "  << nEvents[kClassA] << " " << nEvents[kClassC] << " " << nEvents[kClassAC] << " " <<  nEvents[kClassB] << endl;
1508           // cout << "Acc: " << acc << endl;
1509           // cout << "BG: " << bg << endl;
1510           // cout  << "  " <<   fBIFactorA* (ratioToB[kClassA]*nEvents[kClassA]-acc) <<endl;
1511           // cout  << "  " <<   fBIFactorC* (ratioToB[kClassC]*nEvents[kClassC]-acc) <<endl;
1512           // cout  << "  " <<   fBIFactorAC*(ratioToB[kClassAC]*nEvents[kClassAC] - 2*acc) << endl;
1513           // cout << "-----------------------" << endl;
1514
1515           Float_t good = Float_t(nEvents[kClassB]) - bg - acc;
1516           if (icol ==1) good1 = good;
1517           //      Float_t errGood     = TMath::Sqrt(2*(nEvents[kClassA]+nEvents[kClassC]+nEvents[kClassE]));// Error on the number of goods assuming only bg fluctuates
1518           //      DeltaG^2 = B + FA^2 A + FC^2 C + Ratio^2 (FA+FC-1)^2 E.
1519           Float_t errGood     = nEvents[kClassAC] > 0 ? 
1520             TMath::Sqrt( nEvents[kClassB] +
1521                          fBIFactorAC*fBIFactorAC*ratioToB[kClassAC]*ratioToB[kClassAC]*nEvents[kClassAC] +
1522                          ratioToB[kClassE] * ratioToB[kClassE] * 
1523                          (fBIFactorAC - 1)*(fBIFactorAC - 1)*nEvents[kClassE]) :  
1524             TMath::Sqrt( nEvents[kClassB] + 
1525                          fBIFactorA*fBIFactorA*ratioToB[kClassA]*ratioToB[kClassA]*nEvents[kClassA] +
1526                          fBIFactorC*fBIFactorC*ratioToB[kClassC]*ratioToB[kClassC]*nEvents[kClassC] +
1527                          ratioToB[kClassE] * ratioToB[kClassE] * 
1528                          (fBIFactorA + fBIFactorC - 1)*(fBIFactorA + fBIFactorC - 1)*nEvents[kClassE]);
1529
1530           Float_t errBG = nEvents[kClassAC] > 0 ? 
1531             TMath::Sqrt(fBIFactorAC*fBIFactorAC*ratioToB[kClassAC]*ratioToB[kClassAC]*nEvents[kClassAC]+
1532                         4*ratioToB[kClassE]*ratioToB[kClassE]*(fBIFactorAC*fBIFactorAC)*nEvents[kClassE]) :
1533             TMath::Sqrt(fBIFactorA*fBIFactorA*ratioToB[kClassA]*ratioToB[kClassA]*nEvents[kClassA]+
1534                         fBIFactorC*fBIFactorC*ratioToB[kClassC]*ratioToB[kClassC]*nEvents[kClassC]+
1535                         ratioToB[kClassE]*ratioToB[kClassE]*(fBIFactorA+fBIFactorC)*(fBIFactorA+fBIFactorC)*nEvents[kClassE]);
1536         
1537         
1538           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllB, nEvents[kClassB]); 
1539           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllAC,nEvents[kClassA]+nEvents[kClassC]+nEvents[kClassAC]);      
1540           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllE, nEvents[kClassE]); 
1541
1542           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowBG,bg);  
1543           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowBG,errBG);       
1544           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAcc,acc);        
1545           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowAcc,accErr);     
1546           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowGood,good);    
1547           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowGood,errGood);    
1548
1549 #ifdef VERBOSE_STAT
1550           //kStatRowBG=0,kStatRowAcc,kStatRowBGFrac,kStatRowAccFrac,kStatRowErrGoodFrac,kStatRowGoodFrac,kStatRowGood,kStatRowErrGood
1551           Float_t accFrac   = Float_t(acc) / nEvents[kClassB]  *100;
1552           Float_t errAccFrac= Float_t(accErr) / nEvents[kClassB]  *100;
1553           Float_t bgFrac    = Float_t(bg)  / nEvents[kClassB]  *100;
1554           Float_t goodFrac  = Float_t(good)  / good1 *100;
1555           Float_t errGoodFrac = errGood/good1 * 100;
1556           Float_t errFracBG = bg > 0 ? TMath::Sqrt((errBG/bg)*(errBG/bg) + 1/nEvents[kClassB])*bgFrac : 0;
1557           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowBGFrac,bgFrac);  
1558           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowBGFrac,errFracBG);       
1559           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAccFrac,accFrac);    
1560           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowAccFrac,errAccFrac);    
1561           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowGoodFrac,goodFrac);    
1562           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowErrGoodFrac,errGoodFrac);    
1563           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowErrGood,errGood);    
1564 #endif
1565         }
1566       }
1567     }
1568     for (Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1569       delete [] rows[iTrigClass];
1570     }  
1571   }  
1572
1573   fHistStatistics[0]->Write();
1574   fHistStatistics[1]->Write();
1575   if(fHistBunchCrossing ) fHistBunchCrossing ->Write();
1576   if(fHistTriggerPattern) fHistTriggerPattern->Write();
1577   
1578   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1579   for (Int_t i=0; i < count; i++)
1580     {
1581       TString triggerClass = "trigger_histograms_";
1582       if (i < fCollTrigClasses.GetEntries())
1583         triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
1584       else
1585         triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
1586       
1587       gDirectory->mkdir(triggerClass);
1588       gDirectory->cd(triggerClass);
1589       
1590       static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
1591       
1592       gDirectory->cd("..");
1593     }
1594   
1595   if (fBackgroundIdentification)
1596     {
1597       gDirectory->mkdir("background_identification");
1598       gDirectory->cd("background_identification");
1599       
1600       fBackgroundIdentification->GetOutput()->Write();
1601       
1602       gDirectory->cd("..");
1603     }
1604   
1605   if (folder)
1606     gDirectory->cd("..");
1607   
1608 }
1609
1610 Int_t AliPhysicsSelection::GetStatRow(const char * triggerBXClass, UInt_t offlineTriggerType, UInt_t ** rowIDs) const {
1611   // Puts inside the array rowIDs the row number for a given offline
1612   // trigger in a given bx class. Returns the total number of lines
1613   // matching the selection
1614   // triggerBXClass can be either "A", "AC", "B" or "E"
1615   // offlineTriggerType is one of the types defined in AliVEvent
1616   // User should delete rowIDs if no longer needed
1617
1618   if(!fHistStatistics[0]) {
1619     AliWarning("Not initialized, returning 0");
1620     return 0;
1621   }
1622   const Int_t nrows = fHistStatistics[0]->GetNbinsY();
1623
1624   // allocate memory for at maximum nrows
1625   Int_t nMatches = 0;
1626   (*rowIDs) = new UInt_t[nrows];
1627
1628
1629   // Loop over rows and find matching ones, using the PS OADB
1630   // FIXME: check BG estimates
1631   for(Int_t irow = 1; irow <= nrows; irow++){
1632     TString triggerClassCurrent = fHistStatistics[0]->GetYaxis()->GetBinLabel(irow);
1633     TPRegexp trigType (".*\\&(\\d+).*");
1634     TObjArray * arr = trigType.MatchS(triggerClassCurrent.Data());
1635     if(arr->GetEntries() > 1){
1636       UInt_t tType = ((TObjString*)arr->At(1))->GetString().Atoi();
1637       
1638       if (tType == offlineTriggerType && fPSOADB->GetBeamSide(triggerClassCurrent.Data()) == triggerBXClass) {
1639         (*rowIDs)[nMatches] = irow;
1640         nMatches++;
1641       }
1642
1643     }
1644     delete arr;
1645
1646   }
1647   
1648   return nMatches;
1649 }
1650
1651 void AliPhysicsSelection::SetBIFactors(const AliESDEvent * aESD) {
1652   // Set factors for realtive bunch intesities
1653   if(!aESD) { 
1654     AliFatal("ESD not given");
1655   }
1656   Int_t run = aESD->GetRunNumber();
1657   if (run > 105268) {
1658     // intensities stored in the ESDs
1659     const AliESDRun* esdRun = aESD->GetESDRun();
1660     Double_t intAB = esdRun->GetMeanIntensityIntecting(0);
1661     Double_t intCB = esdRun->GetMeanIntensityIntecting(1);
1662     Double_t intAA = esdRun->GetMeanIntensityNonIntecting(0);
1663     Double_t intCC = esdRun->GetMeanIntensityNonIntecting(1);
1664
1665     // cout << "INT " <<intAB <<endl;
1666     // cout << "INT " <<intCB <<endl;
1667     // cout << "INT " <<intAA <<endl;
1668     // cout << "INT " <<intCC <<endl;
1669
1670     if (intAB > 0 && intAA > 0) {
1671       fBIFactorA = intAB/intAA;
1672     } else {
1673       AliWarning("Cannot set fBIFactorA");
1674     }
1675     
1676     if (intCB > 0 && intCC > 0) {
1677       fBIFactorC = intCB/intCC;
1678     } else {
1679       AliWarning("Cannot set fBIFactorC");
1680     }
1681       
1682     if (intAB > 0 && intAA > 0 &&
1683         intCB > 0 && intCC > 0) {
1684       fBIFactorAC = (intAB+intCB)/(intAA+intCC);
1685     } else {
1686       AliWarning("Cannot set fBIFactorAC");
1687     }
1688         
1689   }  
1690   else {
1691     // First runs. Intensities hardcoded
1692     switch(run) {
1693     case 104155:
1694       fBIFactorA = 0.961912722908;
1695       fBIFactorC = 1.04992336081;
1696       break;
1697     case 104157:
1698       fBIFactorA = 0.947312854998;
1699       fBIFactorC = 1.01599706417;
1700       break;
1701     case 104159:
1702       fBIFactorA = 0.93659320151;
1703       fBIFactorC = 0.98580804207;
1704       break;
1705     case 104160:
1706       fBIFactorA = 0.929664189926;
1707       fBIFactorC = 0.963467679851;
1708       break;
1709     case 104315:
1710       fBIFactorA = 1.08939104979;
1711       fBIFactorC = 0.931113921925;
1712       break;
1713     case 104316:
1714       fBIFactorA = 1.08351880974;
1715       fBIFactorC = 0.916068345845;
1716       break;
1717     case 104320:
1718       fBIFactorA = 1.07669281245;
1719       fBIFactorC = 0.876818744763;
1720       break;
1721     case 104321:
1722       fBIFactorA = 1.00971079602;
1723       fBIFactorC = 0.773781299076;
1724       break;
1725     case 104792:
1726       fBIFactorA = 0.787215863962;
1727       fBIFactorC = 0.778253173071;
1728       break;
1729     case 104793:
1730       fBIFactorA = 0.692211363661;
1731       fBIFactorC = 0.733152456667;
1732       break;
1733     case 104799:
1734       fBIFactorA = 1.04027825161;
1735       fBIFactorC = 1.00530825942;
1736       break;
1737     case 104800:
1738       fBIFactorA = 1.05309910671;
1739       fBIFactorC = 1.00376801855;
1740       break;
1741     case 104801:
1742       fBIFactorA = 1.0531231922;
1743       fBIFactorC = 0.992439666758;
1744       break;
1745     case 104802:
1746       fBIFactorA = 1.04191478134;
1747       fBIFactorC = 0.979368585208;
1748       break;
1749     case 104803:
1750       fBIFactorA = 1.03121314094;
1751       fBIFactorC = 0.973379962609;
1752       break;
1753     case 104824:
1754       fBIFactorA = 0.969945926722;
1755       fBIFactorC = 0.39549745806;
1756       break;
1757     case 104825:
1758       fBIFactorA = 0.968627213937;
1759       fBIFactorC = 0.310100412205;
1760       break;
1761     case 104841:
1762       fBIFactorA = 0.991601393212;
1763       fBIFactorC = 0.83762204722;
1764       break;
1765     case 104845:
1766       fBIFactorA = 0.98040863886;
1767       fBIFactorC = 0.694824205793;
1768       break;
1769     case 104867:
1770       fBIFactorA = 1.10646173412;
1771       fBIFactorC = 0.841407246916;
1772       break;
1773     case 104876:
1774       fBIFactorA = 1.12063452421;
1775       fBIFactorC = 0.78726542895;
1776       break;
1777     case 104890:
1778       fBIFactorA = 1.02346137453;
1779       fBIFactorC = 1.03355663595;
1780       break;
1781     case 104892:
1782       fBIFactorA = 1.05406025913;
1783       fBIFactorC = 1.00029166135;
1784       break;
1785     case 105143:
1786       fBIFactorA = 0.947343384349;
1787       fBIFactorC = 0.972637444408;
1788       break;
1789     case 105160:
1790       fBIFactorA = 0.908854622177;
1791       fBIFactorC = 0.958851103977;
1792       break; 
1793     case 105256:
1794       fBIFactorA = 0.810076150206;
1795       fBIFactorC = 0.884663561883;
1796       break;
1797     case 105257:
1798       fBIFactorA = 0.80974912303;
1799       fBIFactorC = 0.878859123479;
1800       break;
1801     case 105268:
1802       fBIFactorA = 0.809052110679;
1803       fBIFactorC = 0.87233890989;
1804       break;
1805     default:
1806       fBIFactorA = 1;
1807       fBIFactorC = 1;
1808     }
1809   }
1810
1811 }
1812
1813 const char * AliPhysicsSelection::GetTriggerString(TObjString * obj) { 
1814   // Returns a formed objstring
1815   static TString retString;
1816   
1817   retString.Form("%s%s", 
1818                  obj->String().Data(), 
1819                  fUseBXNumbers ? fFillOADB->GetBXIDs(fPSOADB->GetBeamSide(obj->String().Data())) : ""  
1820                  );
1821
1822   if (fMC) {
1823     TPRegexp stripClasses("\\+\\S* ");
1824     stripClasses.Substitute(retString,"","g");
1825     stripClasses=TPRegexp("\\-\\S* ");
1826     stripClasses.Substitute(retString,"","g");
1827   }
1828
1829   return retString.Data();
1830 }
1831
1832 void AliPhysicsSelection::AddCollisionTriggerClass(const char* className){
1833   AliError("This method is deprecated! Will be removed soon! Please use SetCustomOADBObjects() instead!");
1834   if(!fPSOADB) {
1835     fPSOADB = new AliOADBPhysicsSelection("CustomPS");   
1836     fPSOADB->SetHardwareTrigger         ( 0,"SPDGFO >= 1 || V0A || V0C");
1837     fPSOADB->SetOfflineTrigger          ( 0,"(SPDGFO >= 1 || V0A || V0C) && !V0ABG && !V0CBG");
1838   }
1839
1840   // Strip all which is not needed, if BX ids are provided, they are still used here
1841   // offline trigger and logics are appended automagically, so we strip-em out if they are provided
1842   TString classNameStripped = className;
1843   if(classNameStripped.Index("*")>0)
1844     classNameStripped.Remove(classNameStripped.Index("*")); // keep only the class name (no bx, offline trigger...)   
1845   if(classNameStripped.Index("&")>0)
1846     classNameStripped.Remove(classNameStripped.Index("&")); // keep only the class name (no bx, offline trigger...)   
1847
1848   fPSOADB->AddCollisionTriggerClass   ( AliVEvent::kUserDefined,classNameStripped.Data(),"B",0);
1849   
1850   fUsingCustomClasses = kTRUE;
1851
1852
1853 }
1854
1855 void AliPhysicsSelection::AddBGTriggerClass(const char* className){
1856   // Add custom BG trigger class
1857  
1858   AliError("This method is deprecated! Will be removed soon! Please use SetCustomOADBObjects() instead!");
1859   if(!fPSOADB) {
1860     fPSOADB = new AliOADBPhysicsSelection("CustomPS");   
1861     fPSOADB->SetHardwareTrigger         ( 0,"SPDGFO >= 1 || V0A || V0C");
1862     fPSOADB->SetOfflineTrigger          ( 0,"(SPDGFO >= 1 || V0A || V0C) && !V0ABG && !V0CBG");
1863   }
1864
1865   // Strip all which is not needed, if BX ids are provided, they are still used here
1866   // offline trigger and logics are appended automagically, so we strip-em out if they are provided
1867   TString classNameStripped = className;
1868   if(classNameStripped.Index("*")>0)
1869     classNameStripped.Remove(classNameStripped.Index("*")); // keep only the class name (no bx, offline trigger...)   
1870   if(classNameStripped.Index("&")>0)
1871     classNameStripped.Remove(classNameStripped.Index("&")); // keep only the class name (no bx, offline trigger...)   
1872
1873   fPSOADB->AddBGTriggerClass   ( AliVEvent::kUserDefined,classNameStripped.Data(),"AC",0);
1874
1875   fUsingCustomClasses = kTRUE;
1876
1877 }