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