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