]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelection.cxx
- new class for TPC laser warmupt tracks rejection
[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       // Some "macros"
524       Bool_t mb1 = (fastOROffline > 0 || v0A || v0C) && (!v0BG);
525       Bool_t mb1prime = (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C) ) && (!v0BG);
526
527       // Background rejection
528       Bool_t bgID = kFALSE;
529       if (fBackgroundIdentification)
530         bgID = ! fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd));
531         
532       
533       /*Int_t ntrig = fastOROffline; // any 2 hits
534       if(v0A)              ntrig += 1;
535       if(v0C)              ntrig += 1; //v0C alone is enough
536       if(fmd)              ntrig += 1;
537       if(ssdClusters>1)    ntrig += 1;*/
538
539       
540       // <---
541
542       TString triggerLogicOnline  = fPSOADB->GetHardwareTrigger(triggerLogic);
543       TString triggerLogicOffline = fPSOADB->GetOfflineTrigger(triggerLogic);
544
545       AliDebug(AliLog::kDebug, Form("Triggers from OADB [0x%x][%d][%s][%s]",singleTriggerResult,AliOADBPhysicsSelection::GetActiveBit(singleTriggerResult),triggerLogicOffline.Data(),triggerLogicOnline.Data()));
546
547       // replay hardware trigger (should only remove events for MC)
548       Bool_t onlineTrigger  = EvaluateTriggerLogic(aEsd, triggerAnalysis, triggerLogicOnline, kFALSE);
549       // offline selection
550       Bool_t offlineTrigger = EvaluateTriggerLogic(aEsd, triggerAnalysis, triggerLogicOffline, kTRUE);
551       
552       // Printf("%s %s", triggerLogicOnline.Data(), triggerLogicOffline.Data());
553       // Printf("%d %d", onlineTrigger, offlineTrigger);
554            
555
556       // Fill trigger pattern histo
557       Int_t tpatt = 0;
558       if (fastORHW>0) tpatt+=1;
559       if (v0AHW)      tpatt+=2;
560       if (v0CHW)      tpatt+=4;
561       fHistTriggerPattern->Fill( tpatt );
562
563       // fill statistics and return decision
564       const Int_t nHistStat = 2;
565       for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
566         if (iHistStat == kStatIdxBin0 && !isBin0) continue; // skip the filling of bin0 stats if the event is not in the bin0
567       
568         fHistStatistics[iHistStat]->Fill(kStatTriggerClass, i);
569         if(iHistStat == kStatIdxAll) fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i); // Fill only for all (avoid double counting)
570
571         // We fill the rest only if hw trigger is ok
572         if (!onlineTrigger)
573           {
574             AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
575             continue;
576           } else {       
577           fHistStatistics[iHistStat]->Fill(kStatHWTrig, i);
578         }
579       
580
581         // v0 BG stats
582         if (v0ABG)
583           fHistStatistics[iHistStat]->Fill(kStatV0ABG, i);
584         if (v0CBG)
585           fHistStatistics[iHistStat]->Fill(kStatV0CBG, i);
586
587         // mb 1
588         if (mb1)
589           fHistStatistics[iHistStat]->Fill(kStatMB1, i);
590
591         if (mb1prime)
592           fHistStatistics[iHistStat]->Fill(kStatMB1Prime, i);
593
594         if (fmd)
595           fHistStatistics[iHistStat]->Fill(kStatFMD, i);
596
597         //if(ntrig >= 2 && !v0BG) 
598         //  fHistStatistics[iHistStat]->Fill(kStatAny2Hits, i);
599
600         if (fastOROffline > 0)
601           fHistStatistics[iHistStat]->Fill(kStatFO1, i);
602         if (fastOROffline > 1)
603           fHistStatistics[iHistStat]->Fill(kStatFO2, i);
604         if (fastOROfflineL1 > 1)
605           fHistStatistics[iHistStat]->Fill(kStatFO2L1, i);
606         
607         if (v0A)
608           fHistStatistics[iHistStat]->Fill(kStatV0A, i);
609         if (v0C)
610           fHistStatistics[iHistStat]->Fill(kStatV0C, i);
611           
612         if (zdcA)
613           fHistStatistics[iHistStat]->Fill(kStatZDCA, i);
614         if (zdcC)
615           fHistStatistics[iHistStat]->Fill(kStatZDCC, i);
616         if (zdcA && zdcC)
617           fHistStatistics[iHistStat]->Fill(kStatZDCAC, i);
618
619         if (zdcTime)
620           fHistStatistics[iHistStat]->Fill(kStatZDCTime, i);
621   
622         if (v0A && v0C && !v0BG && (!bgID && fIsPP))
623           fHistStatistics[iHistStat]->Fill(kStatV0, i);
624
625
626
627         // FIXME: check lines below
628         if ( offlineTrigger )
629           {
630             if (!v0BG || fSkipV0)
631               {
632                 if (!v0BG) fHistStatistics[iHistStat]->Fill(kStatOffline, i);
633       
634                 if (fBackgroundIdentification && bgID && fIsPP)
635                   {
636                     AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
637                     fHistStatistics[iHistStat]->Fill(kStatBG, i);
638                   }
639                 else
640                   {
641                     AliDebug(AliLog::kDebug, Form("Accepted event for histograms with trigger logic %d", triggerLogic));
642             
643                     fHistStatistics[iHistStat]->Fill(kStatAccepted, i);
644                     // if(iHistStat == kStatIdxAll) fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i); // Fill only for all (avoid double counting)
645                     if((i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection) && (iHistStat==kStatIdxAll))
646                       accept |= singleTriggerResult; // only set for "all" (should not really matter)
647                   }
648               }
649             else
650               AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
651           }
652         else
653           AliDebug(AliLog::kDebug, Form("Rejecting event because trigger logic %d is not fulfilled", triggerLogic));
654       }
655     }
656   }
657  
658   if (accept)
659     AliDebug(AliLog::kDebug, Form("Accepted event as collision candidate with bit mask %d", accept));
660   
661   return accept;
662 }
663
664
665 // const char * AliPhysicsSelection::GetFillingScheme(UInt_t runNumber)  {
666
667 //   if(fMC) return "MC";
668
669 //   if      (runNumber >= 104065 && runNumber <= 104160) {
670 //     return "4x4a";
671 //   } 
672 //   else if (runNumber >= 104315 && runNumber <= 104321) {
673 //     return "4x4a*";
674 //   }
675 //   else if (runNumber >= 104792 && runNumber <= 104803) {
676 //     return "4x4b";
677 //   }
678 //   else if (runNumber >= 104824 && runNumber <= 104892) {
679 //     return "4x4c";
680 //   }
681 //   else if (runNumber == 105143 || runNumber == 105160) {
682 //     return "16x16a";
683 //   }
684 //   else if (runNumber >= 105256 && runNumber <= 105268) {
685 //     return "4x4c";
686 //   } 
687 //   else if (runNumber >= 114786 && runNumber <= 116684) {
688 //     return "Single_2b_1_1_1";
689 //   }
690 //   else if (runNumber >= 117048 && runNumber <= 117120) {
691 //     return "Single_3b_2_2_2";
692 //   }
693 //   else if (runNumber >= 117220 && runNumber <= 119163) {
694 //     return "Single_2b_1_1_1";
695 //   }
696 //   else if (runNumber >= 119837 && runNumber <= 119862) {
697 //     return "Single_4b_2_2_2";
698 //   }
699 //   else if (runNumber >= 119902 && runNumber <= 120691) {
700 //     return "Single_6b_3_3_3";
701 //   }
702 //   else if (runNumber >= 120741 && runNumber <= 122375) {
703 //     return "Single_13b_8_8_8";
704 //   }
705 //   else if (runNumber >= 130148 && runNumber <= 130375) {
706 //     return "125n_48b_36_16_36";
707 //   } 
708 //   else if (runNumber >= 130601 && runNumber <= 130640) {
709 //     return "1000ns_50b_35_14_35";
710 //   }
711 //   else {
712 //     AliError(Form("Unknown filling scheme (run %d)", runNumber));
713 //   }
714
715 //   return "Unknown";
716 // }
717
718 // const char * AliPhysicsSelection::GetBXIDs(UInt_t runNumber, const char * trigger)  {
719
720 //   if (!fUseBXNumbers || fMC) return "";
721
722 //   if      (runNumber >= 104065 && runNumber <= 104160) {
723 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2128 #3019";
724 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #346 #3465";
725 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #1680";
726 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
727 //     //    else AliError(Form("Unknown trigger: %s", trigger));
728 //   } 
729 //   else if (runNumber >= 104315 && runNumber <= 104321) {
730 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2000 #2891";
731 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #218 #3337";
732 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1106 #1552";
733 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
734 //     //    else AliError(Form("Unknown trigger: %s", trigger));
735 //   }
736 //   else if (runNumber >= 104792 && runNumber <= 104803) {
737 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #2228 #3119";
738 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
739 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #769";
740 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
741 //     //    else AliError(Form("Unknown trigger: %s", trigger));
742 //   }
743 //   else if (runNumber >= 104824 && runNumber <= 104892) {
744 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3119 #769";
745 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2554 #446";
746 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1334 #2228";
747 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
748 //     //    else AliError(Form("Unknown trigger: %s", trigger));
749 //   }
750 //   else if (runNumber == 105143 || runNumber == 105160) {
751 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #1337 #1418 #2228 #2309 #3119 #3200 #446 #527";
752 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #1580  #1742  #1904  #2066  #2630  #2792  #2954  #3362";
753 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "  #845  #1007  #1169   #1577 #3359 #3521 #119  #281 ";
754 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
755 //     //    else AliError(Form("Unknown trigger: %s", trigger));
756 //   }
757 //   else if (runNumber >= 105256 && runNumber <= 105268) {
758 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #3019 #669";
759 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2454 #346";
760 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #1234 #2128";
761 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #790";
762 //     //    else AliError(Form("Unknown trigger: %s", trigger));
763 //   } else if (runNumber >= 114786 && runNumber <= 116684) { // 7 TeV 2010, assume always the same filling scheme
764 //     if     (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return " #346";
765 //     else if(!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return " #2131";
766 //     else if(!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return " #3019";
767 //     else if(!strcmp("CINT1-E-NOPF-ALL",trigger))     return " #1238";
768 //     //    else AliError(Form("Unknown trigger: %s", trigger));
769 //   }
770 //   else if (runNumber >= 117048 && runNumber <= 117120) {
771 //     //    return "Single_3b_2_2_2";
772 //    if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #1240 ";
773 //    else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131 ";
774 //    else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019 ";
775 //    else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238";
776 //    //   else AliError(Form("Unknown trigger: %s", trigger));
777
778 //   }
779 //   else if ((runNumber >= 117220 && runNumber <= 118555) || (runNumber >= 118784 && runNumber <= 119163)) 
780 //   {
781 //     //    return "Single_2b_1_1_1";
782 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346 ";
783 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131 ";
784 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019 ";
785 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238 ";
786 //     //    else AliError(Form("Unknown trigger: %s", trigger));                                                   
787 //   }
788 //   else if (runNumber >= 118556 && runNumber <= 118783) {
789 //     //    return "Single_2b_1_1_1"; 
790 //     // same as previous but was misaligned by 1 BX in fill 1069
791 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #345 ";
792 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2130 ";
793 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3018 ";
794 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1238 ";
795 //     //    else AliError(Form("Unknown trigger: %s", trigger));                                                   
796 //   }
797 //   else if (runNumber >= 119837 && runNumber <= 119862) {
798 //     //    return "Single_4b_2_2_2";
799 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #669  #3019 ";
800 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #346  #2454 ";
801 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #1234  #2128 ";
802 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1681 #3463";
803 //     //    else AliError(Form("Unknown trigger: %s", trigger));
804
805 //   }
806 //   else if (runNumber >= 119902 && runNumber <= 120691) {
807 //     //    return "Single_6b_3_3_3";
808 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #546  #746 ";
809 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #2131  #2331  #2531 ";
810 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019  #3219  #3419 ";
811 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1296 #1670";
812 //     //    else AliError(Form("Unknown trigger: %s", trigger));
813 //   }
814 //   else if (runNumber >= 120741 && runNumber <= 122375) {
815 //     //    return "Single_13b_8_8_8";
816 //     if      (!strcmp("CINT1B-ABCE-NOPF-ALL",trigger)) return "   #346  #446  #546  #646  #1240  #1340  #1440  #1540 ";
817 //     else if (!strcmp("CINT1A-ABCE-NOPF-ALL",trigger)) return "   #946  #2131  #2231  #2331  #2431 ";
818 //     else if (!strcmp("CINT1C-ABCE-NOPF-ALL",trigger)) return "   #3019  #3119  #3219  #3319  #3519 ";
819 //     else if (!strcmp("CINT1-E-NOPF-ALL",trigger)) return " #1835 #2726";
820 //     //    else AliError(Form("Unknown trigger: %s", trigger));
821     
822 //   } 
823 //   else if (runNumber >= 130148 && runNumber <= 130375) {
824 //     TString triggerString = trigger;
825 //     static TString returnString = " ";
826 //     returnString = "";
827 //     if (triggerString.Contains("B")) returnString += "   #346  #396  #446  #496  #546  #596  #646  #696  #1240  #1290  #1340  #1390  #1440  #1490  #1540  #1590 ";
828 //     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 ";
829 //     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 ";
830 //     // Printf("0x%x",returnString.Data());
831 //     // Printf("%s",returnString.Data());
832 //     return returnString.Data();
833 //   } 
834 //   else if (runNumber >= 130601 && runNumber <= 130640) {
835 //     TString triggerString = trigger;
836 //     static TString returnString = " ";
837 //     returnString = "";
838 //     if (triggerString.Contains("B")) returnString += "  #346  #386  #426  #466  #506  #546  #586  #1240  #1280  #1320  #1360  #1400  #1440  #1480 ";
839 //     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
840 //     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
841 //     return returnString.Data();
842 //   }
843
844 //   else {
845 //     AliWarning(Form("Unknown run %d, using all BXs!",runNumber));
846 //   }
847
848 //   return "";
849 // }
850
851 Bool_t AliPhysicsSelection::Initialize(const AliESDEvent* aEsd)
852 {
853   // initializes the object for the given ESD
854   
855   AliInfo(Form("Initializing for beam type: %s", aEsd->GetESDRun()->GetBeamType()));
856   fIsPP = kTRUE;
857   if (strcmp(aEsd->GetESDRun()->GetBeamType(), "Pb-Pb") == 0)
858     fIsPP = kFALSE;
859
860   return Initialize(aEsd->GetRunNumber());
861 }
862
863 Bool_t AliPhysicsSelection::Initialize(Int_t runNumber)
864 {
865   // initializes the object for the given run  
866
867
868   Bool_t oldStatus = TH1::AddDirectoryStatus();
869   TH1::AddDirectory(kFALSE);
870
871   /// Open OADB file and fetch OADB objects
872   TString oadbfilename = AliPhysicsSelection::GetOADBFileName();
873
874   TFile * foadb = new TFile (oadbfilename); 
875   if(!foadb->IsOpen()) AliFatal(Form("Cannot open OADB file %s", oadbfilename.Data()));
876
877
878   if(!fPSOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
879     AliInfo("Using Standard OADB");
880     AliOADBContainer * psContainer = (AliOADBContainer*) foadb->Get("physSel");
881     if (!psContainer) AliFatal("Cannot fetch OADB container for Physics selection");
882     fPSOADB = (AliOADBPhysicsSelection*) psContainer->GetObject(runNumber, fIsPP ? "oadbDefaultPP" : "oadbDefaultPbPb");
883     if (!fPSOADB) AliFatal(Form("Cannot find physics selection object for run %d", runNumber));
884   } else {
885     AliInfo("Using Custom OADB");
886   }
887   if(!fFillOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
888     AliOADBContainer * fillContainer = (AliOADBContainer*) foadb->Get("fillScheme");
889     if (!fillContainer) AliFatal("Cannot fetch OADB container for filling scheme");
890     fFillOADB = (AliOADBFillingScheme*) fillContainer->GetObject(runNumber, "Default");
891     if (!fFillOADB) AliFatal(Form("Cannot find  filling scheme object for run %d", runNumber));
892   }
893   if(!fTriggerOADB || !fUsingCustomClasses) { // if it's already set and custom class is required, we use the one provided by the user
894     AliOADBContainer * triggerContainer = (AliOADBContainer*) foadb->Get("trigAnalysis");
895     if (!triggerContainer) AliFatal("Cannot fetch OADB container for trigger analysis");
896     fTriggerOADB = (AliOADBTriggerAnalysis*) triggerContainer->GetObject(runNumber, "Default");
897     fTriggerOADB->Print();
898     if (!fTriggerOADB) AliFatal(Form("Cannot find  trigger analysis object for run %d", runNumber));
899   }
900
901   
902   if(!fBin0CallBack) 
903     AliError("Bin0 Callback not set: will not fill the statistics for the bin 0");
904
905   if (fMC) {
906     // override BX and bg options in case of MC
907     fComputeBG    = kFALSE;
908     fUseBXNumbers = kFALSE;
909   }
910
911   // FIXME: think how to implement this check with the OADB, trigger scheme is not a int number here 
912   // Int_t triggerScheme = GetTriggerScheme(runNumber);
913   // if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
914   //   AliFatal("Processing several runs with different trigger schemes is not supported");
915   
916   if(fComputeBG && fCurrentRun != -1 && fCurrentRun != runNumber) 
917     AliFatal("Cannot process several runs because BG computation is requested");
918
919   if(fComputeBG && !fUseBXNumbers) 
920     AliFatal("Cannot compute BG if BX numbers are not used");
921   
922   if(fUseBXNumbers && fFillingScheme != "" && fFillingScheme != fFillOADB->GetFillingSchemeName())
923     AliFatal("Cannot process runs with different filling scheme if usage of BX numbers is requested");
924
925   fFillingScheme      = fFillOADB->GetFillingSchemeName();
926
927   AliInfo(Form("Initializing for run %d", runNumber));
928   
929   // initialize first time?
930   if (fCurrentRun == -1)
931   {
932     const UInt_t ntriggerBits = fPSOADB->GetNTriggerBits(); // FIXME!
933     for(UInt_t ibit = 0; ibit < ntriggerBits; ibit++){
934       
935       TIterator * collIter = fPSOADB->GetCollTrigClass(ibit)->MakeIterator();
936       TIterator * bgIter   = fPSOADB->GetBGTrigClass(ibit)->MakeIterator();
937       TObjString * obj = 0;
938       while((obj = (TObjString*) collIter->Next())){
939         if (obj->String() != "") {
940           fCollTrigClasses.Add(new TObjString(GetTriggerString(obj)));
941         }
942       }
943       // BG classes only make sense for real data
944       if(!fMC) {
945         obj = 0 ;
946         while((obj = (TObjString*) bgIter->Next())){
947           if (obj->String() != "") {
948             fBGTrigClasses.Add(new TObjString(GetTriggerString(obj)));
949           }
950         }
951       }
952
953     }
954     // not sure how to handle this in the case of > x cuts
955     // // Book the token histo with the tokens actually used here!
956     // // 1. Cache the tokens
957     // for(UInt_t ibit = 0; ibit < ntriggerBits; ibit++){
958     //  EvaluateTriggerLogic(0,0, fPSOADB->GetHardwareTrigger(ibit), 0);
959     //  EvaluateTriggerLogic(0,0, fPSOADB->GetOfflineTrigger(ibit),  1);
960     // }
961     // // 2. Book the histogram
962     // if (!fHistStatisticsTokens) {
963     //  Int_t ntokens = fCashedTokens->GetEntries();
964     //  Int_t count = fCollTrigClasses->GetEntries() + fBGTrigClasses->GetEntries();
965     //  fHistStatisticsTokens = new TH2F("fHistStatisticsTokens", "fHistStatisticsTokens", ntokens + 2, 0.5, ntokens+2+0.5, count, -0.5, -0.5 + count);
966         
967     //  Int_t nrow = 0;
968     //  fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,  "Trigger class");
969     //  for(Int_t itoken = 0; itoken < ntoken; itoken++){
970     //    TParameter<Int_t> * param = fCashedTokens->At(itoken);
971     //    fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,  param->GetName());      
972     //  }
973         
974     //  fHistStatisticsTokens->GetXaxis()->SetBinLabel(nrow++,      "Accepted");
975
976     // }
977
978     
979     
980     // TODO: 
981     // Add a new statistics histo containing only the tokens actually used in the selection
982
983     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
984     
985     for (Int_t i=0; i<count; i++)
986     {
987       
988       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
989       triggerAnalysis->SetAnalyzeMC(fMC);
990       triggerAnalysis->EnableHistograms();
991       triggerAnalysis->SetSPDGFOThreshhold(1);
992       triggerAnalysis->SetDoFMD(kFALSE);
993       triggerAnalysis->SetCorrZDCCutParams(fTriggerOADB->GetZDCCutRefSumCorr(),
994                                            fTriggerOADB->GetZDCCutRefDeltaCorr(), 
995                                            fTriggerOADB->GetZDCCutSigmaSumCorr(),
996                                            fTriggerOADB->GetZDCCutSigmaDeltaCorr());
997       fTriggerAnalysis.Add(triggerAnalysis);
998     }
999
1000     
1001     // TODO: shall I really delete this?
1002     if (fHistStatistics[0])
1003       delete fHistStatistics[0];
1004     if (fHistStatistics[1])
1005       delete fHistStatistics[1];
1006   
1007     fHistStatistics[kStatIdxBin0] = BookHistStatistics("_Bin0");
1008     fHistStatistics[kStatIdxAll]  = BookHistStatistics("");
1009     
1010
1011     if (fHistBunchCrossing)
1012       delete fHistBunchCrossing;
1013   
1014     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
1015
1016     // TODO: remove fHistTriggerPattern
1017     if (fHistTriggerPattern)
1018       delete fHistTriggerPattern;
1019     
1020     const int ntrig=3;
1021     Int_t n = 1;
1022     const Int_t nbinTrig = TMath::Nint(TMath::Power(2,ntrig));
1023
1024     fHistTriggerPattern = new TH1F("fHistTriggerPattern", "Trigger pattern: FO + 2*v0A + 4*v0C", 
1025                                    nbinTrig, -0.5, nbinTrig-0.5);    
1026     fHistTriggerPattern->GetXaxis()->SetBinLabel(1,"NO TRIG");
1027     fHistTriggerPattern->GetXaxis()->SetBinLabel(2,"FO");
1028     fHistTriggerPattern->GetXaxis()->SetBinLabel(3,"v0A");
1029     fHistTriggerPattern->GetXaxis()->SetBinLabel(4,"FO & v0A");
1030     fHistTriggerPattern->GetXaxis()->SetBinLabel(5,"v0C");
1031     fHistTriggerPattern->GetXaxis()->SetBinLabel(6,"FO & v0C");
1032     fHistTriggerPattern->GetXaxis()->SetBinLabel(7,"v0A & v0C");
1033     fHistTriggerPattern->GetXaxis()->SetBinLabel(8,"FO & v0A & v0C");
1034
1035   
1036     n = 1;
1037     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1038     {
1039       
1040       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
1041       n++;
1042     }
1043     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1044     {
1045       
1046       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
1047       n++;
1048     }
1049
1050     
1051
1052   }
1053     
1054   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1055   for (Int_t i=0; i<count; i++)
1056   {
1057     
1058     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
1059   
1060     switch (runNumber)
1061     {
1062       case 104315:
1063       case 104316:
1064       case 104320:
1065       case 104321:
1066       case 104439:
1067         triggerAnalysis->SetV0TimeOffset(7.5);
1068         break;
1069       default:
1070         triggerAnalysis->SetV0TimeOffset(0);
1071     }
1072   }
1073     
1074   fCurrentRun = runNumber;
1075   
1076   TH1::AddDirectory(oldStatus);
1077   
1078   
1079   return kTRUE;
1080 }
1081
1082 TH2F * AliPhysicsSelection::BookHistStatistics(const char * tag) {
1083   // add 6 rows to count for the estimate of good, accidentals and
1084   // BG and the ratio of BG and accidentals to total +ratio goot to
1085   // first col + 2 for error on good.
1086   // TODO: Remember the the indexes of rows for the BG selection. Add new member fBGRows[] and use kStat as indexes
1087
1088   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1089 #ifdef VERBOSE_STAT
1090   Int_t extrarows = fComputeBG != 0 ? 11 : 0;
1091 #else
1092   Int_t extrarows = fComputeBG != 0 ? 6 : 0;
1093 #endif
1094   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);
1095
1096   h->GetXaxis()->SetBinLabel(kStatTriggerClass,  "Trigger class");
1097   h->GetXaxis()->SetBinLabel(kStatHWTrig,        "Hardware trigger");
1098   h->GetXaxis()->SetBinLabel(kStatFO1,           "FO >= 1");
1099   h->GetXaxis()->SetBinLabel(kStatFO2,           "FO >= 2");
1100   h->GetXaxis()->SetBinLabel(kStatFO2L1,         "FO (L1) >= 2");
1101   h->GetXaxis()->SetBinLabel(kStatV0A,           "V0A");
1102   h->GetXaxis()->SetBinLabel(kStatV0C,           "V0C");
1103   h->GetXaxis()->SetBinLabel(kStatFMD,           "FMD");
1104   h->GetXaxis()->SetBinLabel(kStatV0ABG,         "V0A BG");
1105   h->GetXaxis()->SetBinLabel(kStatV0CBG,         "V0C BG");
1106   h->GetXaxis()->SetBinLabel(kStatZDCA,          "ZDCA");
1107   h->GetXaxis()->SetBinLabel(kStatZDCC,          "ZDCC");
1108   h->GetXaxis()->SetBinLabel(kStatZDCAC,         "ZDCA & ZDCC");
1109   h->GetXaxis()->SetBinLabel(kStatZDCTime,       "ZDC Time Cut");
1110   h->GetXaxis()->SetBinLabel(kStatMB1,           "(FO >= 1 | V0A | V0C) & !V0 BG");
1111   h->GetXaxis()->SetBinLabel(kStatMB1Prime,      "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A &v0C) ) & !V0 BG");
1112   //h->GetXaxis()->SetBinLabel(kStatFO1AndV0,    "FO >= 1 & (V0A | V0C) & !V0 BG");
1113   h->GetXaxis()->SetBinLabel(kStatV0,            "V0A & V0C & !V0 BG & !BG ID");
1114   h->GetXaxis()->SetBinLabel(kStatOffline,       "Offline Trigger");
1115   //h->GetXaxis()->SetBinLabel(kStatAny2Hits,    "2 Hits & !V0 BG");
1116   h->GetXaxis()->SetBinLabel(kStatBG,            "Background identification");
1117   h->GetXaxis()->SetBinLabel(kStatAccepted,      "Accepted");
1118
1119   Int_t n = 1;
1120   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1121     {
1122       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
1123       n++;
1124     }
1125   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1126     {
1127       h->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
1128       n++;
1129     }
1130
1131   if(fComputeBG) {
1132     fBGStatOffset = n;
1133     h->GetYaxis()->SetBinLabel(n++, "All B");
1134     h->GetYaxis()->SetBinLabel(n++, "All A+C");
1135     h->GetYaxis()->SetBinLabel(n++, "All E");
1136     h->GetYaxis()->SetBinLabel(n++, Form("BG (A+C) (Mask [0x%x])", fComputeBG));
1137     h->GetYaxis()->SetBinLabel(n++, "ACC");
1138 #ifdef VERBOSE_STAT
1139     h->GetYaxis()->SetBinLabel(n++, "BG (A+C) %  (rel. to CINT1B)");
1140     h->GetYaxis()->SetBinLabel(n++, "ACC % (rel. to CINT1B)");
1141     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD %");
1142     h->GetYaxis()->SetBinLabel(n++, "GOOD % (rel. to 1st col)");
1143     h->GetYaxis()->SetBinLabel(n++, "ERR GOOD");
1144 #endif
1145     h->GetYaxis()->SetBinLabel(n++, "GOOD");
1146   }
1147
1148   return h;
1149 }
1150
1151 void AliPhysicsSelection::Print(const Option_t *option) const
1152 {
1153   // print the configuration
1154   TString msg;
1155   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
1156   msg += Form("Configuration initialized for run %d (MC: %d):\n", fCurrentRun, fMC);
1157   
1158   Printf("Collision trigger classes:");
1159   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
1160     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1161   
1162   Printf("Background trigger classes:");
1163   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
1164     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
1165
1166   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
1167   
1168   if (triggerAnalysis)
1169   {
1170     if (triggerAnalysis->GetV0TimeOffset() > 0)
1171       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
1172     
1173     Printf("\nTotal available events:");
1174     
1175     triggerAnalysis->PrintTriggerClasses();
1176   }
1177   
1178   if (fHistStatistics[kStatIdxAll])
1179   {
1180     for (Int_t i=0; i<fCollTrigClasses.GetEntries(); i++)
1181     {
1182       Printf("\nSelection statistics for collision trigger %s:", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1183       msg += Form("\nSelection statistics for collision trigger %s:\n", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
1184       
1185       Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(1, i+1));
1186       msg += Form("Total events with correct trigger class: %d\n", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(1, i+1));
1187       Printf("Selected collision candidates: %d", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(fHistStatistics[kStatIdxAll]->GetXaxis()->FindBin("Accepted"), i+1));
1188       msg += Form("Selected collision candidates: %d\n", (Int_t) fHistStatistics[kStatIdxAll]->GetBinContent(fHistStatistics[kStatIdxAll]->GetXaxis()->FindBin("Accepted"), i+1));
1189     }
1190   }
1191   
1192   if (fHistBunchCrossing)
1193   {
1194     Printf("\nBunch crossing statistics:");
1195     msg += "\nBunch crossing statistics:\n";
1196     
1197     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
1198     {
1199       TString str;
1200       str.Form("Trigger %s has events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
1201       
1202       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
1203         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
1204           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
1205              
1206       Printf("%s", str.Data());
1207       msg += str;
1208       msg += "\n";
1209     }
1210     
1211     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
1212     {
1213       Int_t countColl = 0;
1214       Int_t countBG = 0;
1215       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
1216       {
1217         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
1218         {
1219           if (fCollTrigClasses.FindObject(fHistBunchCrossing->GetYaxis()->GetBinLabel(i)))
1220             countColl++;
1221           if (fBGTrigClasses.FindObject(fHistBunchCrossing->GetYaxis()->GetBinLabel(i)))
1222             countBG++;
1223         }
1224       }
1225       if (countColl > 0 && countBG > 0)
1226         Printf("WARNING: Bunch crossing %d has collision and BG trigger classes active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
1227     }
1228   }
1229
1230   if (fUsingCustomClasses)        
1231     Printf("WARNING: Using custom trigger classes!");
1232   if (fSkipTriggerClassSelection) 
1233     Printf("WARNING: Skipping trigger class selection!");
1234   if (fSkipV0) 
1235     Printf("WARNING: Ignoring V0 information in selection");
1236   if(!fBin0CallBack) 
1237     Printf("WARNING: Callback not set: will not fill the statistics for the bin 0");
1238   TString opt(option);
1239   opt.ToUpper();
1240   if (opt == "STAT") {
1241      AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1242      if (mgr) mgr->AddStatisticsMsg(msg);
1243   }
1244 }
1245
1246 Long64_t AliPhysicsSelection::Merge(TCollection* list)
1247 {
1248   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1249   // PROOF).
1250   // Returns the number of merged objects (including this).
1251
1252   if (!list)
1253     return 0;
1254
1255   if (list->IsEmpty())
1256     return 1;
1257
1258   TIterator* iter = list->MakeIterator();
1259   TObject* obj;
1260   
1261   // collections of all histograms
1262   const Int_t nHists = 9;
1263   TList collections[nHists];
1264
1265   Int_t count = 0;
1266   while ((obj = iter->Next())) {
1267
1268     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
1269     if (entry == 0) 
1270       continue;
1271     // Update run number. If this one is not initialized (-1) take the one from 
1272     // the next physics selection to be merged with. In case of 2 different run
1273     // numbers issue a warning (should physics selections from different runs be 
1274     // merged together) A.G.
1275     Int_t currentRun = entry->GetCurrentRun();
1276     // Nothing to merge with since run number was not initialized.
1277     if (currentRun < 0) continue;
1278     if (fCurrentRun < 0) 
1279     {
1280       fCurrentRun = currentRun;
1281       fMC = entry->fMC;
1282       for (Int_t i=0; i<entry->fCollTrigClasses.GetEntries(); i++)
1283         fCollTrigClasses.Add(entry->fCollTrigClasses.At(i)->Clone());
1284       for (Int_t i=0; i<entry->fBGTrigClasses.GetEntries(); i++)
1285         fBGTrigClasses.Add(entry->fBGTrigClasses.At(i)->Clone());
1286     }
1287     if (fCurrentRun != currentRun)
1288        AliWarning(Form("Current run %d not matching the one to be merged with %d", fCurrentRun, currentRun));
1289
1290     // With the same strategy update fBGStatOffset
1291     Int_t bgstatoffset = entry->GetBGStatOffset();
1292     
1293     // Nothing to merge with since run number was not initialized.
1294     if (bgstatoffset < 0) continue;
1295     if (fBGStatOffset < 0) 
1296     {
1297       fBGStatOffset = bgstatoffset;
1298     }
1299     if (fBGStatOffset != bgstatoffset)
1300        AliWarning(Form("Current run %d not matching the one to be merged with %d", fBGStatOffset, bgstatoffset));
1301     
1302
1303     
1304     // Merge the OADBs (Take just the first instance you find
1305     if (!fPSOADB) {
1306       fPSOADB = (AliOADBPhysicsSelection*) entry->GetOADBPhysicsSelection()->Clone();
1307     }
1308     if (!fFillOADB){
1309       fFillOADB = (AliOADBFillingScheme*) entry->GetOADBFillingScheme()->Clone();
1310     }
1311
1312     if (entry->fTriggerAnalysis.GetEntries() > 0)
1313       collections[0].Add(&(entry->fTriggerAnalysis));
1314     if (entry->fHistStatistics[0])
1315       collections[1].Add(entry->fHistStatistics[0]);
1316     if (entry->fHistStatistics[1])
1317       collections[2].Add(entry->fHistStatistics[1]);
1318     // if (entry->fHistStatisticsTokens)
1319     //   collections[3].Add(entry->fHistStatisticsTokens);
1320     if (entry->fHistBunchCrossing)
1321       collections[3].Add(entry->fHistBunchCrossing);
1322     if (entry->fHistTriggerPattern)
1323       collections[4].Add(entry->fHistTriggerPattern);
1324     if (entry->fBackgroundIdentification)
1325       collections[5].Add(entry->fBackgroundIdentification);
1326
1327     count++;
1328   }
1329
1330   if (fTriggerAnalysis.GetEntries() == 0 && collections[0].GetEntries() > 0)
1331   {
1332     TList* firstList = (TList*) collections[0].First();
1333     for (Int_t i=0; i<firstList->GetEntries(); i++)
1334       fTriggerAnalysis.Add(firstList->At(i)->Clone());
1335     
1336     collections[0].RemoveAt(0);
1337   }
1338   fTriggerAnalysis.Merge(&collections[0]);
1339   
1340   // if this instance is empty (not initialized) nothing would be merged here --> copy first entry
1341   if (!fHistStatistics[0] && collections[1].GetEntries() > 0)
1342   {
1343     fHistStatistics[0] = (TH2F*) collections[1].First()->Clone();
1344     collections[1].RemoveAt(0);
1345   }
1346   if (fHistStatistics[0])
1347     fHistStatistics[0]->Merge(&collections[1]);
1348     
1349   if (!fHistStatistics[1] && collections[2].GetEntries() > 0)
1350   {
1351     fHistStatistics[1] = (TH2F*) collections[2].First()->Clone();
1352     collections[2].RemoveAt(0);
1353   }
1354   if (fHistStatistics[1])
1355     fHistStatistics[1]->Merge(&collections[2]);
1356     
1357   if (!fHistBunchCrossing && collections[3].GetEntries() > 0)
1358   {
1359     fHistBunchCrossing = (TH2F*) collections[3].First()->Clone();
1360     collections[3].RemoveAt(0);
1361   }
1362   if (fHistBunchCrossing)
1363     fHistBunchCrossing->Merge(&collections[3]);
1364   
1365   if (!fHistTriggerPattern && collections[4].GetEntries() > 0)
1366   {
1367     fHistTriggerPattern = (TH1F*) collections[4].First()->Clone();
1368     collections[4].RemoveAt(0);
1369   }
1370   if (fHistTriggerPattern)
1371     fHistTriggerPattern->Merge(&collections[4]);
1372   
1373   if (!fBackgroundIdentification && collections[5].GetEntries() > 0)
1374   {
1375     fBackgroundIdentification = (AliAnalysisCuts*) collections[5].First()->Clone();
1376     collections[5].RemoveAt(0);
1377   }
1378   if (fBackgroundIdentification)
1379     fBackgroundIdentification->Merge(&collections[5]);
1380   
1381   delete iter;
1382
1383   return count+1;
1384 }
1385
1386 void AliPhysicsSelection::SaveHistograms(const char* folder) 
1387 {
1388   // write histograms to current directory
1389   
1390   if (!fHistStatistics[0] || !fHistStatistics[1])
1391     return;
1392     
1393   if (folder)
1394   {
1395     gDirectory->mkdir(folder);
1396     gDirectory->cd(folder);
1397   }
1398   
1399
1400   // Fill the last rows of fHistStatistics before saving
1401   if (fComputeBG) {      
1402     AliInfo("BG estimate assumes that for a given run you only have A and C triggers separately or"
1403             " toghether as a AC class! Make sure this assumption holds in your case"); 
1404     
1405     // use an anum for the different trigger classes, to make loops easier to read
1406     enum {kClassB =0 , kClassA, kClassC, kClassAC, kClassE, kNClasses};
1407     const char * classFlags[] = {"B", "A", "C", "AC", "E"}; // labels
1408     
1409     UInt_t * rows[kNClasses] = {0}; // Array of matching rows
1410     Int_t nrows[kNClasses] = {0};
1411     // Get rows matching the requested trigger bits for all trigger classes
1412     for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1413       nrows[iTrigClass] = GetStatRow(classFlags[iTrigClass],fComputeBG,&rows[iTrigClass]);      
1414     }
1415     
1416     // 0. Determine the ratios of triggers E/B, A/B, C/B from the stat histogram
1417     // Those are used to rescale the different classes to the same number of bx ids
1418     // 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?         
1419     Int_t nBXIds[kNClasses] = {0};
1420       //      cout <<"Computing BG:" << endl;
1421     
1422     for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1423       for(Int_t irow = 0; irow < nrows[iTrigClass]; irow++) {
1424         if(irow==0) cout << "- Class " << classFlags[iTrigClass] << endl;   
1425         for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++) {
1426           if (fHistBunchCrossing->GetBinContent(j, rows[iTrigClass][irow]) > 0) {
1427             nBXIds[iTrigClass]++;        
1428           }
1429         }
1430         if(nBXIds[iTrigClass]>0) cout << "   Using row " << rows[iTrigClass][irow] <<  ": " 
1431                                       << fHistBunchCrossing->GetYaxis()->GetBinLabel(rows[iTrigClass][irow]) 
1432                                       << " (nBXID "<< nBXIds[iTrigClass] << ")"<< endl;
1433
1434       }
1435
1436     }
1437
1438     Float_t ratioToB[kNClasses];
1439     ratioToB[kClassE]  = nBXIds[kClassE]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassE]   : 0;
1440     ratioToB[kClassA]  = nBXIds[kClassA]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassA]   : 0;
1441     ratioToB[kClassC]  = nBXIds[kClassC]  >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassC]   : 0;
1442     ratioToB[kClassAC] = nBXIds[kClassAC] >0 ? Float_t(nBXIds[kClassB])/nBXIds[kClassAC]  : 0;
1443     Printf("Ratio between the BX ids in the different trigger classes:");
1444     Printf("  B/E  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassE], ratioToB[kClassE] );
1445     Printf("  B/A  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassA], ratioToB[kClassA] );
1446     Printf("  B/C  = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassC], ratioToB[kClassC] );
1447     Printf("  B/AC = %d/%d = %f", nBXIds[kClassB],nBXIds[kClassAC],ratioToB[kClassAC]);
1448     Int_t nHistStat = 2;
1449
1450     // 1. loop over all cols
1451     for(Int_t iHistStat = 0; iHistStat < nHistStat; iHistStat++){
1452       Int_t ncol = fHistStatistics[iHistStat]->GetNbinsX();
1453       Float_t good1 = 0;      
1454       for(Int_t icol = 1; icol <= ncol; icol++) {
1455         Int_t nEvents[kNClasses] = {0}; // number of events should be reset at every column
1456         // For all trigger classes, add up over row matching trigger mask (as selected before)
1457         for(Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1458           for(Int_t irow = 0; irow < nrows[iTrigClass]; irow++) {         
1459             nEvents[iTrigClass] += (Int_t) fHistStatistics[iHistStat]->GetBinContent(icol,rows[iTrigClass][irow]);      
1460           }
1461           //        cout << "Events " << classFlags[iTrigClass] << " ("<<icol<<") " << nEvents[iTrigClass] << endl;         
1462         }
1463         if (nEvents[kClassB]>0) {
1464           Float_t acc  = ratioToB[kClassE]*nEvents[kClassE]; 
1465           Double_t accErr = TMath::Sqrt(ratioToB[kClassE]*ratioToB[kClassE]*nEvents[kClassE]);
1466           //      Int_t bg   = cint1A + cint1C - 2*acc;
1467             
1468           // If intensity measurements are available, they already
1469           // contain the scaling for BX ratios, so we reset the
1470           // ratioToB entries
1471           if(icol == 1) {
1472             if(fBIFactorAC > 0 || fBIFactorA > 0 || fBIFactorC > 0) {
1473               if (fBIFactorAC <= 0 && (fBIFactorA <= 0 || fBIFactorC <= 0)) {
1474                 AliError("Not all intensities set!, assuming equal intensities");
1475                 fBIFactorA  = 1;
1476                 fBIFactorC  = 1;
1477                 fBIFactorAC = 1;
1478               } else {
1479                 AliInfo("Using ratio of number of bunch crossing embedded in the intensity measurements");
1480                 ratioToB[kClassA]  = ratioToB[kClassA]  >0 ? 1  : 0;
1481                 ratioToB[kClassC]  = ratioToB[kClassC]  >0 ? 1  : 0;
1482                 ratioToB[kClassAC] = ratioToB[kClassAC] >0 ? 1  : 0;
1483                 AliInfo(Form(" - BI Factor A:  %f",  fBIFactorA ));
1484                 AliInfo(Form(" - BI Factor C:  %f",  fBIFactorC ));
1485                 AliInfo(Form(" - BI Factor AC: %f",  fBIFactorAC ));
1486                 
1487               }
1488             } else {
1489               AliWarning("Intensities not set!, assuming equal intensities");
1490               fBIFactorA  = 1;
1491               fBIFactorC  = 1;
1492               fBIFactorAC = 1;
1493             }
1494           }
1495           // Assuming that for a given class the triggers are either recorded as A+C or AC
1496           Float_t bg  = nEvents[kClassAC] > 0 ?
1497             fBIFactorAC*(ratioToB[kClassAC]*nEvents[kClassAC] - 2*acc):
1498             fBIFactorA* (ratioToB[kClassA]*nEvents[kClassA]-acc) + 
1499             fBIFactorC* (ratioToB[kClassC]*nEvents[kClassC]-acc) ;
1500
1501           // cout << "-----------------------" << endl;
1502           // cout << "Factors: " << fBIFactorA << " " << fBIFactorC << " " << fBIFactorAC << endl;
1503           // cout << "Ratios: "  << ratioToB[kClassA] << " " << ratioToB[kClassC] << " " << ratioToB[kClassAC] << endl;
1504           // cout << "Evts:   "  << nEvents[kClassA] << " " << nEvents[kClassC] << " " << nEvents[kClassAC] << " " <<  nEvents[kClassB] << endl;
1505           // cout << "Acc: " << acc << endl;
1506           // cout << "BG: " << bg << endl;
1507           // cout  << "  " <<   fBIFactorA* (ratioToB[kClassA]*nEvents[kClassA]-acc) <<endl;
1508           // cout  << "  " <<   fBIFactorC* (ratioToB[kClassC]*nEvents[kClassC]-acc) <<endl;
1509           // cout  << "  " <<   fBIFactorAC*(ratioToB[kClassAC]*nEvents[kClassAC] - 2*acc) << endl;
1510           // cout << "-----------------------" << endl;
1511
1512           Float_t good = Float_t(nEvents[kClassB]) - bg - acc;
1513           if (icol ==1) good1 = good;
1514           //      Float_t errGood     = TMath::Sqrt(2*(nEvents[kClassA]+nEvents[kClassC]+nEvents[kClassE]));// Error on the number of goods assuming only bg fluctuates
1515           //      DeltaG^2 = B + FA^2 A + FC^2 C + Ratio^2 (FA+FC-1)^2 E.
1516           Float_t errGood     = nEvents[kClassAC] > 0 ? 
1517             TMath::Sqrt( nEvents[kClassB] +
1518                          fBIFactorAC*fBIFactorAC*ratioToB[kClassAC]*ratioToB[kClassAC]*nEvents[kClassAC] +
1519                          ratioToB[kClassE] * ratioToB[kClassE] * 
1520                          (fBIFactorAC - 1)*(fBIFactorAC - 1)*nEvents[kClassE]) :  
1521             TMath::Sqrt( nEvents[kClassB] + 
1522                          fBIFactorA*fBIFactorA*ratioToB[kClassA]*ratioToB[kClassA]*nEvents[kClassA] +
1523                          fBIFactorC*fBIFactorC*ratioToB[kClassC]*ratioToB[kClassC]*nEvents[kClassC] +
1524                          ratioToB[kClassE] * ratioToB[kClassE] * 
1525                          (fBIFactorA + fBIFactorC - 1)*(fBIFactorA + fBIFactorC - 1)*nEvents[kClassE]);
1526
1527           Float_t errBG = nEvents[kClassAC] > 0 ? 
1528             TMath::Sqrt(fBIFactorAC*fBIFactorAC*ratioToB[kClassAC]*ratioToB[kClassAC]*nEvents[kClassAC]+
1529                         4*ratioToB[kClassE]*ratioToB[kClassE]*(fBIFactorAC*fBIFactorAC)*nEvents[kClassE]) :
1530             TMath::Sqrt(fBIFactorA*fBIFactorA*ratioToB[kClassA]*ratioToB[kClassA]*nEvents[kClassA]+
1531                         fBIFactorC*fBIFactorC*ratioToB[kClassC]*ratioToB[kClassC]*nEvents[kClassC]+
1532                         ratioToB[kClassE]*ratioToB[kClassE]*(fBIFactorA+fBIFactorC)*(fBIFactorA+fBIFactorC)*nEvents[kClassE]);
1533         
1534         
1535           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllB, nEvents[kClassB]); 
1536           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllAC,nEvents[kClassA]+nEvents[kClassC]+nEvents[kClassAC]);      
1537           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAllE, nEvents[kClassE]); 
1538
1539           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowBG,bg);  
1540           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowBG,errBG);       
1541           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAcc,acc);        
1542           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowAcc,accErr);     
1543           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowGood,good);    
1544           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowGood,errGood);    
1545
1546 #ifdef VERBOSE_STAT
1547           //kStatRowBG=0,kStatRowAcc,kStatRowBGFrac,kStatRowAccFrac,kStatRowErrGoodFrac,kStatRowGoodFrac,kStatRowGood,kStatRowErrGood
1548           Float_t accFrac   = Float_t(acc) / nEvents[kClassB]  *100;
1549           Float_t errAccFrac= Float_t(accErr) / nEvents[kClassB]  *100;
1550           Float_t bgFrac    = Float_t(bg)  / nEvents[kClassB]  *100;
1551           Float_t goodFrac  = Float_t(good)  / good1 *100;
1552           Float_t errGoodFrac = errGood/good1 * 100;
1553           Float_t errFracBG = bg > 0 ? TMath::Sqrt((errBG/bg)*(errBG/bg) + 1/nEvents[kClassB])*bgFrac : 0;
1554           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowBGFrac,bgFrac);  
1555           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowBGFrac,errFracBG);       
1556           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowAccFrac,accFrac);    
1557           fHistStatistics[iHistStat]->SetBinError  (icol,fBGStatOffset+kStatRowAccFrac,errAccFrac);    
1558           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowGoodFrac,goodFrac);    
1559           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowErrGoodFrac,errGoodFrac);    
1560           fHistStatistics[iHistStat]->SetBinContent(icol,fBGStatOffset+kStatRowErrGood,errGood);    
1561 #endif
1562         }
1563       }
1564     }
1565     for (Int_t iTrigClass = 0; iTrigClass < kNClasses; iTrigClass++){
1566       delete [] rows[iTrigClass];
1567     }  
1568   }  
1569
1570   fHistStatistics[0]->Write();
1571   fHistStatistics[1]->Write();
1572   if(fHistBunchCrossing ) fHistBunchCrossing ->Write();
1573   if(fHistTriggerPattern) fHistTriggerPattern->Write();
1574   
1575   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
1576   for (Int_t i=0; i < count; i++)
1577     {
1578       TString triggerClass = "trigger_histograms_";
1579       if (i < fCollTrigClasses.GetEntries())
1580         triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
1581       else
1582         triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
1583       
1584       gDirectory->mkdir(triggerClass);
1585       gDirectory->cd(triggerClass);
1586       
1587       static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
1588       
1589       gDirectory->cd("..");
1590     }
1591   
1592   if (fBackgroundIdentification)
1593     {
1594       gDirectory->mkdir("background_identification");
1595       gDirectory->cd("background_identification");
1596       
1597       fBackgroundIdentification->GetOutput()->Write();
1598       
1599       gDirectory->cd("..");
1600     }
1601   
1602   if (folder)
1603     gDirectory->cd("..");
1604   
1605 }
1606
1607 Int_t AliPhysicsSelection::GetStatRow(const char * triggerBXClass, UInt_t offlineTriggerType, UInt_t ** rowIDs) const {
1608   // Puts inside the array rowIDs the row number for a given offline
1609   // trigger in a given bx class. Returns the total number of lines
1610   // matching the selection
1611   // triggerBXClass can be either "A", "AC", "B" or "E"
1612   // offlineTriggerType is one of the types defined in AliVEvent
1613   // User should delete rowIDs if no longer needed
1614
1615   if(!fHistStatistics[0]) {
1616     AliWarning("Not initialized, returning 0");
1617     return 0;
1618   }
1619   const Int_t nrows = fHistStatistics[0]->GetNbinsY();
1620
1621   // allocate memory for at maximum nrows
1622   Int_t nMatches = 0;
1623   (*rowIDs) = new UInt_t[nrows];
1624
1625
1626   // Loop over rows and find matching ones, using the PS OADB
1627   // FIXME: check BG estimates
1628   for(Int_t irow = 1; irow <= nrows; irow++){
1629     TString triggerClassCurrent = fHistStatistics[0]->GetYaxis()->GetBinLabel(irow);
1630     TPRegexp trigType (".*\\&(\\d+).*");
1631     TObjArray * arr = trigType.MatchS(triggerClassCurrent.Data());
1632     if(arr->GetEntries() > 1){
1633       UInt_t tType = ((TObjString*)arr->At(1))->GetString().Atoi();
1634       
1635       if (tType == offlineTriggerType && fPSOADB->GetBeamSide(triggerClassCurrent.Data()) == triggerBXClass) {
1636         (*rowIDs)[nMatches] = irow;
1637         nMatches++;
1638       }
1639
1640     }
1641     delete arr;
1642
1643   }
1644   
1645   return nMatches;
1646 }
1647
1648 void AliPhysicsSelection::SetBIFactors(const AliESDEvent * aESD) {
1649   // Set factors for realtive bunch intesities
1650   if(!aESD) { 
1651     AliFatal("ESD not given");
1652   }
1653   Int_t run = aESD->GetRunNumber();
1654   if (run > 105268) {
1655     // intensities stored in the ESDs
1656     const AliESDRun* esdRun = aESD->GetESDRun();
1657     Double_t intAB = esdRun->GetMeanIntensityIntecting(0);
1658     Double_t intCB = esdRun->GetMeanIntensityIntecting(1);
1659     Double_t intAA = esdRun->GetMeanIntensityNonIntecting(0);
1660     Double_t intCC = esdRun->GetMeanIntensityNonIntecting(1);
1661
1662     // cout << "INT " <<intAB <<endl;
1663     // cout << "INT " <<intCB <<endl;
1664     // cout << "INT " <<intAA <<endl;
1665     // cout << "INT " <<intCC <<endl;
1666
1667     if (intAB > 0 && intAA > 0) {
1668       fBIFactorA = intAB/intAA;
1669     } else {
1670       AliWarning("Cannot set fBIFactorA");
1671     }
1672     
1673     if (intCB > 0 && intCC > 0) {
1674       fBIFactorC = intCB/intCC;
1675     } else {
1676       AliWarning("Cannot set fBIFactorC");
1677     }
1678       
1679     if (intAB > 0 && intAA > 0 &&
1680         intCB > 0 && intCC > 0) {
1681       fBIFactorAC = (intAB+intCB)/(intAA+intCC);
1682     } else {
1683       AliWarning("Cannot set fBIFactorAC");
1684     }
1685         
1686   }  
1687   else {
1688     // First runs. Intensities hardcoded
1689     switch(run) {
1690     case 104155:
1691       fBIFactorA = 0.961912722908;
1692       fBIFactorC = 1.04992336081;
1693       break;
1694     case 104157:
1695       fBIFactorA = 0.947312854998;
1696       fBIFactorC = 1.01599706417;
1697       break;
1698     case 104159:
1699       fBIFactorA = 0.93659320151;
1700       fBIFactorC = 0.98580804207;
1701       break;
1702     case 104160:
1703       fBIFactorA = 0.929664189926;
1704       fBIFactorC = 0.963467679851;
1705       break;
1706     case 104315:
1707       fBIFactorA = 1.08939104979;
1708       fBIFactorC = 0.931113921925;
1709       break;
1710     case 104316:
1711       fBIFactorA = 1.08351880974;
1712       fBIFactorC = 0.916068345845;
1713       break;
1714     case 104320:
1715       fBIFactorA = 1.07669281245;
1716       fBIFactorC = 0.876818744763;
1717       break;
1718     case 104321:
1719       fBIFactorA = 1.00971079602;
1720       fBIFactorC = 0.773781299076;
1721       break;
1722     case 104792:
1723       fBIFactorA = 0.787215863962;
1724       fBIFactorC = 0.778253173071;
1725       break;
1726     case 104793:
1727       fBIFactorA = 0.692211363661;
1728       fBIFactorC = 0.733152456667;
1729       break;
1730     case 104799:
1731       fBIFactorA = 1.04027825161;
1732       fBIFactorC = 1.00530825942;
1733       break;
1734     case 104800:
1735       fBIFactorA = 1.05309910671;
1736       fBIFactorC = 1.00376801855;
1737       break;
1738     case 104801:
1739       fBIFactorA = 1.0531231922;
1740       fBIFactorC = 0.992439666758;
1741       break;
1742     case 104802:
1743       fBIFactorA = 1.04191478134;
1744       fBIFactorC = 0.979368585208;
1745       break;
1746     case 104803:
1747       fBIFactorA = 1.03121314094;
1748       fBIFactorC = 0.973379962609;
1749       break;
1750     case 104824:
1751       fBIFactorA = 0.969945926722;
1752       fBIFactorC = 0.39549745806;
1753       break;
1754     case 104825:
1755       fBIFactorA = 0.968627213937;
1756       fBIFactorC = 0.310100412205;
1757       break;
1758     case 104841:
1759       fBIFactorA = 0.991601393212;
1760       fBIFactorC = 0.83762204722;
1761       break;
1762     case 104845:
1763       fBIFactorA = 0.98040863886;
1764       fBIFactorC = 0.694824205793;
1765       break;
1766     case 104867:
1767       fBIFactorA = 1.10646173412;
1768       fBIFactorC = 0.841407246916;
1769       break;
1770     case 104876:
1771       fBIFactorA = 1.12063452421;
1772       fBIFactorC = 0.78726542895;
1773       break;
1774     case 104890:
1775       fBIFactorA = 1.02346137453;
1776       fBIFactorC = 1.03355663595;
1777       break;
1778     case 104892:
1779       fBIFactorA = 1.05406025913;
1780       fBIFactorC = 1.00029166135;
1781       break;
1782     case 105143:
1783       fBIFactorA = 0.947343384349;
1784       fBIFactorC = 0.972637444408;
1785       break;
1786     case 105160:
1787       fBIFactorA = 0.908854622177;
1788       fBIFactorC = 0.958851103977;
1789       break; 
1790     case 105256:
1791       fBIFactorA = 0.810076150206;
1792       fBIFactorC = 0.884663561883;
1793       break;
1794     case 105257:
1795       fBIFactorA = 0.80974912303;
1796       fBIFactorC = 0.878859123479;
1797       break;
1798     case 105268:
1799       fBIFactorA = 0.809052110679;
1800       fBIFactorC = 0.87233890989;
1801       break;
1802     default:
1803       fBIFactorA = 1;
1804       fBIFactorC = 1;
1805     }
1806   }
1807
1808 }
1809
1810 const char * AliPhysicsSelection::GetTriggerString(TObjString * obj) { 
1811   // Returns a formed objstring
1812   static TString retString;
1813   
1814   retString.Form("%s%s", 
1815                  obj->String().Data(), 
1816                  fUseBXNumbers ? fFillOADB->GetBXIDs(fPSOADB->GetBeamSide(obj->String().Data())) : ""  
1817                  );
1818
1819   if (fMC) {
1820     TPRegexp stripClasses("\\+\\S* ");
1821     stripClasses.Substitute(retString,"","g");
1822     stripClasses=TPRegexp("\\-\\S* ");
1823     stripClasses.Substitute(retString,"","g");
1824   }
1825
1826   return retString.Data();
1827 }
1828
1829 void AliPhysicsSelection::AddCollisionTriggerClass(const char* className){
1830   AliError("This method is deprecated! Will be removed soon! Please use SetCustomOADBObjects() instead!");
1831   if(!fPSOADB) {
1832     fPSOADB = new AliOADBPhysicsSelection("CustomPS");   
1833     fPSOADB->SetHardwareTrigger         ( 0,"SPDGFO >= 1 || V0A || V0C");
1834     fPSOADB->SetOfflineTrigger          ( 0,"(SPDGFO >= 1 || V0A || V0C) && !V0ABG && !V0CBG");
1835   }
1836
1837   // Strip all which is not needed, if BX ids are provided, they are still used here
1838   // offline trigger and logics are appended automagically, so we strip-em out if they are provided
1839   TString classNameStripped = className;
1840   if(classNameStripped.Index("*")>0)
1841     classNameStripped.Remove(classNameStripped.Index("*")); // keep only the class name (no bx, offline trigger...)   
1842   if(classNameStripped.Index("&")>0)
1843     classNameStripped.Remove(classNameStripped.Index("&")); // keep only the class name (no bx, offline trigger...)   
1844
1845   fPSOADB->AddCollisionTriggerClass   ( AliVEvent::kUserDefined,classNameStripped.Data(),"B",0);
1846   
1847   fUsingCustomClasses = kTRUE;
1848
1849
1850 }
1851
1852 void AliPhysicsSelection::AddBGTriggerClass(const char* className){
1853   // Add custom BG trigger class
1854  
1855   AliError("This method is deprecated! Will be removed soon! Please use SetCustomOADBObjects() instead!");
1856   if(!fPSOADB) {
1857     fPSOADB = new AliOADBPhysicsSelection("CustomPS");   
1858     fPSOADB->SetHardwareTrigger         ( 0,"SPDGFO >= 1 || V0A || V0C");
1859     fPSOADB->SetOfflineTrigger          ( 0,"(SPDGFO >= 1 || V0A || V0C) && !V0ABG && !V0CBG");
1860   }
1861
1862   // Strip all which is not needed, if BX ids are provided, they are still used here
1863   // offline trigger and logics are appended automagically, so we strip-em out if they are provided
1864   TString classNameStripped = className;
1865   if(classNameStripped.Index("*")>0)
1866     classNameStripped.Remove(classNameStripped.Index("*")); // keep only the class name (no bx, offline trigger...)   
1867   if(classNameStripped.Index("&")>0)
1868     classNameStripped.Remove(classNameStripped.Index("&")); // keep only the class name (no bx, offline trigger...)   
1869
1870   fPSOADB->AddBGTriggerClass   ( AliVEvent::kUserDefined,classNameStripped.Data(),"AC",0);
1871
1872   fUsingCustomClasses = kTRUE;
1873
1874 }