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