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